Actual source code: lmvmmat.c
petsc-3.5.0 2014-06-30
1: #include <../src/tao/matrix/lmvmmat.h> /*I "lmvmmat.h" */
2: #include <petsctao.h> /*I "petsctao.h" */
3: #include <petscksp.h>
5: #define TaoMid(a,b,c) (((a) < (b)) ? \
6: (((b) < (c)) ? (b) : \
7: (((a) < (c)) ? (c) : (a))) : \
8: (((a) < (c)) ? (a) : \
9: (((b) < (c)) ? (c) : (b))))
11: /* These lists are used for setting options */
12: static const char *Scale_Table[64] = {"none","scalar","broyden"};
14: static const char *Rescale_Table[64] = {"none","scalar","gl"};
16: static const char *Limit_Table[64] = {"none","average","relative","absolute"};
20: /*@C
21: MatCreateLMVM - Creates a limited memory matrix for lmvm algorithms.
23: Collective on A
25: Input Parameters:
26: + comm - MPI Communicator
27: . n - local size of vectors
28: - N - global size of vectors
30: Output Parameters:
31: . A - New LMVM matrix
33: Level: developer
35: @*/
36: extern PetscErrorCode MatCreateLMVM(MPI_Comm comm, PetscInt n, PetscInt N, Mat *A)
37: {
38: MatLMVMCtx *ctx;
40: PetscInt nhistory;
43: /* create data structure and populate with default values */
44: PetscNew(&ctx);
45: ctx->lm=5;
46: ctx->eps=0.0;
47: ctx->limitType=MatLMVM_Limit_None;
48: ctx->scaleType=MatLMVM_Scale_Broyden;
49: ctx->rScaleType = MatLMVM_Rescale_Scalar;
50: ctx->s_alpha = 1.0;
51: ctx->r_alpha = 1.0;
52: ctx->r_beta = 0.5;
53: ctx->mu = 1.0;
54: ctx->nu = 100.0;
55:
56: ctx->phi = 0.125;
57:
58: ctx->scalar_history = 1;
59: ctx->rescale_history = 1;
60:
61: ctx->delta_min = 1e-7;
62: ctx->delta_max = 100.0;
64: /* Begin configuration */
65: PetscOptionsInt("-tao_lmm_vectors", "vectors to use for approximation", "", ctx->lm, &ctx->lm, 0);
66: PetscOptionsReal("-tao_lmm_limit_mu", "mu limiting factor", "", ctx->mu, &ctx->mu, 0);
67: PetscOptionsReal("-tao_lmm_limit_nu", "nu limiting factor", "", ctx->nu, &ctx->nu, 0);
68: PetscOptionsReal("-tao_lmm_broyden_phi", "phi factor for Broyden scaling", "", ctx->phi, &ctx->phi, 0);
69: PetscOptionsReal("-tao_lmm_scalar_alpha", "alpha factor for scalar scaling", "",ctx->s_alpha, &ctx->s_alpha, 0);
70: PetscOptionsReal("-tao_lmm_rescale_alpha", "alpha factor for rescaling diagonal", "", ctx->r_alpha, &ctx->r_alpha, 0);
71: PetscOptionsReal("-tao_lmm_rescale_beta", "beta factor for rescaling diagonal", "", ctx->r_beta, &ctx->r_beta, 0);
72: PetscOptionsInt("-tao_lmm_scalar_history", "amount of history for scalar scaling", "", ctx->scalar_history, &ctx->scalar_history, 0);
73: PetscOptionsInt("-tao_lmm_rescale_history", "amount of history for rescaling diagonal", "", ctx->rescale_history, &ctx->rescale_history, 0);
74: PetscOptionsReal("-tao_lmm_eps", "rejection tolerance", "", ctx->eps, &ctx->eps, 0);
75: PetscOptionsEList("-tao_lmm_scale_type", "scale type", "", Scale_Table, MatLMVM_Scale_Types, Scale_Table[ctx->scaleType], &ctx->scaleType, 0);
76: PetscOptionsEList("-tao_lmm_rescale_type", "rescale type", "", Rescale_Table, MatLMVM_Rescale_Types, Rescale_Table[ctx->rScaleType], &ctx->rScaleType, 0);
77: PetscOptionsEList("-tao_lmm_limit_type", "limit type", "", Limit_Table, MatLMVM_Limit_Types, Limit_Table[ctx->limitType], &ctx->limitType, 0);
78: PetscOptionsReal("-tao_lmm_delta_min", "minimum delta value", "", ctx->delta_min, &ctx->delta_min, 0);
79: PetscOptionsReal("-tao_lmm_delta_max", "maximum delta value", "", ctx->delta_max, &ctx->delta_max, 0);
80:
81: /* Complete configuration */
82: ctx->rescale_history = PetscMin(ctx->rescale_history, ctx->lm);
83: PetscMalloc1(ctx->lm+1,&ctx->rho);
84: PetscMalloc1(ctx->lm+1,&ctx->beta);
86: nhistory = PetscMax(ctx->scalar_history,1);
87: PetscMalloc1(nhistory,&ctx->yy_history);
88: PetscMalloc1(nhistory,&ctx->ys_history);
89: PetscMalloc1(nhistory,&ctx->ss_history);
91: nhistory = PetscMax(ctx->rescale_history,1);
92: PetscMalloc1(nhistory,&ctx->yy_rhistory);
93: PetscMalloc1(nhistory,&ctx->ys_rhistory);
94: PetscMalloc1(nhistory,&ctx->ss_rhistory);
96: /* Finish initializations */
97: ctx->lmnow = 0;
98: ctx->iter = 0;
99: ctx->nupdates = 0;
100: ctx->nrejects = 0;
101: ctx->delta = 1.0;
103: ctx->Gprev = 0;
104: ctx->Xprev = 0;
106: ctx->scale = 0;
107: ctx->useScale = PETSC_FALSE;
109: ctx->H0 = 0;
110: ctx->useDefaultH0=PETSC_TRUE;
112: MatCreateShell(comm, n, n, N, N, ctx, A);
113: MatShellSetOperation(*A,MATOP_DESTROY,(void(*)(void))MatDestroy_LMVM);
114: MatShellSetOperation(*A,MATOP_VIEW,(void(*)(void))MatView_LMVM);
115: return(0);
116: }
120: extern PetscErrorCode MatLMVMSolve(Mat A, Vec b, Vec x)
121: {
122: PetscReal sq, yq, dd;
123: PetscInt ll;
124: PetscBool scaled;
125: MatLMVMCtx *shell;
132: MatShellGetContext(A,(void**)&shell);
133: if (shell->lmnow < 1) {
134: shell->rho[0] = 1.0;
135: }
137: VecCopy(b,x);
138: for (ll = 0; ll < shell->lmnow; ++ll) {
139: VecDot(x,shell->S[ll],&sq);
140: shell->beta[ll] = sq * shell->rho[ll];
141: VecAXPY(x,-shell->beta[ll],shell->Y[ll]);
142: }
144: scaled = PETSC_FALSE;
145: if (!scaled && !shell->useDefaultH0 && shell->H0) {
146: MatSolve(shell->H0,x,shell->U);
147: VecDot(x,shell->U,&dd);
148: if ((dd > 0.0) && !PetscIsInfOrNanReal(dd)) {
149: /* Accept Hessian solve */
150: VecCopy(shell->U,x);
151: scaled = PETSC_TRUE;
152: }
153: }
155: if (!scaled && shell->useScale) {
156: VecPointwiseMult(shell->U,x,shell->scale);
157: VecDot(x,shell->U,&dd);
158: if ((dd > 0.0) && !PetscIsInfOrNanReal(dd)) {
159: /* Accept scaling */
160: VecCopy(shell->U,x);
161: scaled = PETSC_TRUE;
162: }
163: }
165: if (!scaled) {
166: switch(shell->scaleType) {
167: case MatLMVM_Scale_None:
168: break;
170: case MatLMVM_Scale_Scalar:
171: VecScale(x,shell->sigma);
172: break;
174: case MatLMVM_Scale_Broyden:
175: VecPointwiseMult(x,x,shell->D);
176: break;
177: }
178: }
179: for (ll = shell->lmnow-1; ll >= 0; --ll) {
180: VecDot(x,shell->Y[ll],&yq);
181: VecAXPY(x,shell->beta[ll]-yq*shell->rho[ll],shell->S[ll]);
182: }
183: return(0);
184: }
188: extern PetscErrorCode MatView_LMVM(Mat A, PetscViewer pv)
189: {
190: PetscBool isascii;
192: MatLMVMCtx *lmP;
195: MatShellGetContext(A,(void**)&lmP);
196: PetscObjectTypeCompare((PetscObject)pv,PETSCVIEWERASCII,&isascii);
197: if (isascii) {
198: PetscViewerASCIIPrintf(pv,"LMVM Matrix\n");
199: PetscViewerASCIIPrintf(pv," Number of vectors: %D\n",lmP->lm);
200: PetscViewerASCIIPrintf(pv," scale type: %s\n",Scale_Table[lmP->scaleType]);
201: PetscViewerASCIIPrintf(pv," rescale type: %s\n",Rescale_Table[lmP->rScaleType]);
202: PetscViewerASCIIPrintf(pv," limit type: %s\n",Limit_Table[lmP->limitType]);
203: PetscViewerASCIIPrintf(pv," updates: %D\n",lmP->nupdates);
204: PetscViewerASCIIPrintf(pv," rejects: %D\n",lmP->nrejects);
205: }
206: return(0);
207: }
211: extern PetscErrorCode MatDestroy_LMVM(Mat M)
212: {
213: MatLMVMCtx *ctx;
217: MatShellGetContext(M,(void**)&ctx);
218: if (ctx->allocated) {
219: if (ctx->Xprev) {
220: PetscObjectDereference((PetscObject)ctx->Xprev);
221: }
222: if (ctx->Gprev) {
223: PetscObjectDereference((PetscObject)ctx->Gprev);
224: }
226: VecDestroyVecs(ctx->lm+1,&ctx->S);
227: VecDestroyVecs(ctx->lm+1,&ctx->Y);
228: VecDestroy(&ctx->D);
229: VecDestroy(&ctx->U);
230: VecDestroy(&ctx->V);
231: VecDestroy(&ctx->W);
232: VecDestroy(&ctx->P);
233: VecDestroy(&ctx->Q);
234: if (ctx->scale) {
235: VecDestroy(&ctx->scale);
236: }
237: }
238: PetscFree(ctx->rho);
239: PetscFree(ctx->beta);
240: PetscFree(ctx->yy_history);
241: PetscFree(ctx->ys_history);
242: PetscFree(ctx->ss_history);
243: PetscFree(ctx->yy_rhistory);
244: PetscFree(ctx->ys_rhistory);
245: PetscFree(ctx->ss_rhistory);
246: PetscFree(ctx);
247: return(0);
248: }
252: extern PetscErrorCode MatLMVMReset(Mat M)
253: {
255: MatLMVMCtx *ctx;
256: PetscInt i;
259: MatShellGetContext(M,(void**)&ctx);
260: if (ctx->Gprev) {
261: PetscObjectDereference((PetscObject)ctx->Gprev);
262: }
263: if (ctx->Xprev) {
264: PetscObjectDereference((PetscObject)ctx->Xprev);
265: }
266: ctx->Gprev = ctx->Y[ctx->lm];
267: ctx->Xprev = ctx->S[ctx->lm];
268: PetscObjectReference((PetscObject)ctx->Gprev);
269: PetscObjectReference((PetscObject)ctx->Xprev);
270: for (i=0; i<ctx->lm; ++i) {
271: ctx->rho[i] = 0.0;
272: }
273: ctx->rho[0] = 1.0;
275: /* Set the scaling and diagonal scaling matrix */
276: switch(ctx->scaleType) {
277: case MatLMVM_Scale_None:
278: ctx->sigma = 1.0;
279: break;
280: case MatLMVM_Scale_Scalar:
281: ctx->sigma = ctx->delta;
282: break;
283: case MatLMVM_Scale_Broyden:
284: VecSet(ctx->D,ctx->delta);
285: break;
286: }
288: ctx->iter=0;
289: ctx->nupdates=0;
290: ctx->lmnow=0;
291: return(0);
292: }
296: extern PetscErrorCode MatLMVMUpdate(Mat M, Vec x, Vec g)
297: {
298: MatLMVMCtx *ctx;
299: PetscReal rhotemp, rhotol;
300: PetscReal y0temp, s0temp;
301: PetscReal yDy, yDs, sDs;
302: PetscReal sigmanew, denom;
304: PetscInt i;
305: PetscBool same;
306: PetscReal yy_sum=0.0, ys_sum=0.0, ss_sum=0.0;
311: PetscObjectTypeCompare((PetscObject)M,MATSHELL,&same);
312: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix M is not type MatLMVM");
313: MatShellGetContext(M,(void**)&ctx);
314: if (!ctx->allocated) {
315: MatLMVMAllocateVectors(M, x);
316: }
318: if (0 == ctx->iter) {
319: MatLMVMReset(M);
320: } else {
321: VecAYPX(ctx->Gprev,-1.0,g);
322: VecAYPX(ctx->Xprev,-1.0,x);
324: VecDot(ctx->Gprev,ctx->Xprev,&rhotemp);
325: VecDot(ctx->Gprev,ctx->Gprev,&y0temp);
327: rhotol = ctx->eps * y0temp;
328: if (rhotemp > rhotol) {
329: ++ctx->nupdates;
331: ctx->lmnow = PetscMin(ctx->lmnow+1, ctx->lm);
332: ierr=PetscObjectDereference((PetscObject)ctx->S[ctx->lm]);
333: ierr=PetscObjectDereference((PetscObject)ctx->Y[ctx->lm]);
334: for (i = ctx->lm-1; i >= 0; --i) {
335: ctx->S[i+1] = ctx->S[i];
336: ctx->Y[i+1] = ctx->Y[i];
337: ctx->rho[i+1] = ctx->rho[i];
338: }
339: ctx->S[0] = ctx->Xprev;
340: ctx->Y[0] = ctx->Gprev;
341: PetscObjectReference((PetscObject)ctx->S[0]);
342: PetscObjectReference((PetscObject)ctx->Y[0]);
343: ctx->rho[0] = 1.0 / rhotemp;
345: /* Compute the scaling */
346: switch(ctx->scaleType) {
347: case MatLMVM_Scale_None:
348: break;
350: case MatLMVM_Scale_Scalar:
351: /* Compute s^T s */
352: VecDot(ctx->Xprev,ctx->Xprev,&s0temp);
354: /* Scalar is positive; safeguards are not required. */
356: /* Save information for scalar scaling */
357: ctx->yy_history[(ctx->nupdates - 1) % ctx->scalar_history] = y0temp;
358: ctx->ys_history[(ctx->nupdates - 1) % ctx->scalar_history] = rhotemp;
359: ctx->ss_history[(ctx->nupdates - 1) % ctx->scalar_history] = s0temp;
361: /* Compute summations for scalar scaling */
362: yy_sum = 0; /* No safeguard required; y^T y > 0 */
363: ys_sum = 0; /* No safeguard required; y^T s > 0 */
364: ss_sum = 0; /* No safeguard required; s^T s > 0 */
365: for (i = 0; i < PetscMin(ctx->nupdates, ctx->scalar_history); ++i) {
366: yy_sum += ctx->yy_history[i];
367: ys_sum += ctx->ys_history[i];
368: ss_sum += ctx->ss_history[i];
369: }
371: if (0.0 == ctx->s_alpha) {
372: /* Safeguard ys_sum */
373: if (0.0 == ys_sum) {
374: ys_sum = TAO_ZERO_SAFEGUARD;
375: }
377: sigmanew = ss_sum / ys_sum;
378: } else if (1.0 == ctx->s_alpha) {
379: /* Safeguard yy_sum */
380: if (0.0 == yy_sum) {
381: yy_sum = TAO_ZERO_SAFEGUARD;
382: }
384: sigmanew = ys_sum / yy_sum;
385: } else {
386: denom = 2*ctx->s_alpha*yy_sum;
388: /* Safeguard denom */
389: if (0.0 == denom) {
390: denom = TAO_ZERO_SAFEGUARD;
391: }
393: sigmanew = ((2*ctx->s_alpha-1)*ys_sum + PetscSqrtScalar((2*ctx->s_alpha-1)*(2*ctx->s_alpha-1)*ys_sum*ys_sum - 4*(ctx->s_alpha)*(ctx->s_alpha-1)*yy_sum*ss_sum)) / denom;
394: }
396: switch(ctx->limitType) {
397: case MatLMVM_Limit_Average:
398: if (1.0 == ctx->mu) {
399: ctx->sigma = sigmanew;
400: } else if (ctx->mu) {
401: ctx->sigma = ctx->mu * sigmanew + (1.0 - ctx->mu) * ctx->sigma;
402: }
403: break;
405: case MatLMVM_Limit_Relative:
406: if (ctx->mu) {
407: ctx->sigma = TaoMid((1.0 - ctx->mu) * ctx->sigma, sigmanew, (1.0 + ctx->mu) * ctx->sigma);
408: }
409: break;
411: case MatLMVM_Limit_Absolute:
412: if (ctx->nu) {
413: ctx->sigma = TaoMid(ctx->sigma - ctx->nu, sigmanew, ctx->sigma + ctx->nu);
414: }
415: break;
417: default:
418: ctx->sigma = sigmanew;
419: break;
420: }
421: break;
423: case MatLMVM_Scale_Broyden:
424: /* Original version */
425: /* Combine DFP and BFGS */
427: /* This code appears to be numerically unstable. We use the */
428: /* original version because this was used to generate all of */
429: /* the data and because it may be the least unstable of the */
430: /* bunch. */
432: /* P = Q = inv(D); */
433: VecCopy(ctx->D,ctx->P);
434: VecReciprocal(ctx->P);
435: VecCopy(ctx->P,ctx->Q);
437: /* V = y*y */
438: VecPointwiseMult(ctx->V,ctx->Gprev,ctx->Gprev);
440: /* W = inv(D)*s */
441: VecPointwiseMult(ctx->W,ctx->Xprev,ctx->P);
442: VecDot(ctx->W,ctx->Xprev,&sDs);
444: /* Safeguard rhotemp and sDs */
445: if (0.0 == rhotemp) {
446: rhotemp = TAO_ZERO_SAFEGUARD;
447: }
449: if (0.0 == sDs) {
450: sDs = TAO_ZERO_SAFEGUARD;
451: }
453: if (1.0 != ctx->phi) {
454: /* BFGS portion of the update */
455: /* U = (inv(D)*s)*(inv(D)*s) */
456: VecPointwiseMult(ctx->U,ctx->W,ctx->W);
458: /* Assemble */
459: VecAXPY(ctx->P,1.0/rhotemp,ctx->V);
460: VecAXPY(ctx->P,-1.0/sDs,ctx->U);
461: }
463: if (0.0 != ctx->phi) {
464: /* DFP portion of the update */
465: /* U = inv(D)*s*y */
466: VecPointwiseMult(ctx->U, ctx->W, ctx->Gprev);
468: /* Assemble */
469: VecAXPY(ctx->Q,1.0/rhotemp + sDs/(rhotemp*rhotemp), ctx->V);
470: VecAXPY(ctx->Q,-2.0/rhotemp,ctx->U);
471: }
473: if (0.0 == ctx->phi) {
474: VecCopy(ctx->P,ctx->U);
475: } else if (1.0 == ctx->phi) {
476: VecCopy(ctx->Q,ctx->U);
477: } else {
478: /* Broyden update U=(1-phi)*P + phi*Q */
479: VecCopy(ctx->Q,ctx->U);
480: VecAXPBY(ctx->U,1.0-ctx->phi, ctx->phi, ctx->P);
481: }
483: /* Obtain inverse and ensure positive definite */
484: VecReciprocal(ctx->U);
485: VecAbs(ctx->U);
487: switch(ctx->rScaleType) {
488: case MatLMVM_Rescale_None:
489: break;
491: case MatLMVM_Rescale_Scalar:
492: case MatLMVM_Rescale_GL:
493: if (ctx->rScaleType == MatLMVM_Rescale_GL) {
494: /* Gilbert and Lemarachal use the old diagonal */
495: VecCopy(ctx->D,ctx->P);
496: } else {
497: /* The default version uses the current diagonal */
498: VecCopy(ctx->U,ctx->P);
499: }
501: /* Compute s^T s */
502: VecDot(ctx->Xprev,ctx->Xprev,&s0temp);
504: /* Save information for special cases of scalar rescaling */
505: ctx->yy_rhistory[(ctx->nupdates - 1) % ctx->rescale_history] = y0temp;
506: ctx->ys_rhistory[(ctx->nupdates - 1) % ctx->rescale_history] = rhotemp;
507: ctx->ss_rhistory[(ctx->nupdates - 1) % ctx->rescale_history] = s0temp;
509: if (0.5 == ctx->r_beta) {
510: if (1 == PetscMin(ctx->nupdates, ctx->rescale_history)) {
511: VecPointwiseMult(ctx->V,ctx->Y[0],ctx->P);
512: VecDot(ctx->V,ctx->Y[0],&yy_sum);
514: VecPointwiseDivide(ctx->W,ctx->S[0],ctx->P);
515: VecDot(ctx->W,ctx->S[0],&ss_sum);
517: ys_sum = ctx->ys_rhistory[0];
518: } else {
519: VecCopy(ctx->P,ctx->Q);
520: VecReciprocal(ctx->Q);
522: /* Compute summations for scalar scaling */
523: yy_sum = 0; /* No safeguard required */
524: ys_sum = 0; /* No safeguard required */
525: ss_sum = 0; /* No safeguard required */
526: for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) {
527: VecPointwiseMult(ctx->V,ctx->Y[i],ctx->P);
528: VecDot(ctx->V,ctx->Y[i],&yDy);
529: yy_sum += yDy;
531: VecPointwiseMult(ctx->W,ctx->S[i],ctx->Q);
532: VecDot(ctx->W,ctx->S[i],&sDs);
533: ss_sum += sDs;
534: ys_sum += ctx->ys_rhistory[i];
535: }
536: }
537: } else if (0.0 == ctx->r_beta) {
538: if (1 == PetscMin(ctx->nupdates, ctx->rescale_history)) {
539: /* Compute summations for scalar scaling */
540: VecPointwiseDivide(ctx->W,ctx->S[0],ctx->P);
542: VecDot(ctx->W, ctx->Y[0], &ys_sum);
543: VecDot(ctx->W, ctx->W, &ss_sum);
544: yy_sum += ctx->yy_rhistory[0];
545: } else {
546: VecCopy(ctx->Q, ctx->P);
547: VecReciprocal(ctx->Q);
549: /* Compute summations for scalar scaling */
550: yy_sum = 0; /* No safeguard required */
551: ys_sum = 0; /* No safeguard required */
552: ss_sum = 0; /* No safeguard required */
553: for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) {
554: VecPointwiseMult(ctx->W, ctx->S[i], ctx->Q);
555: VecDot(ctx->W, ctx->Y[i], &yDs);
556: ys_sum += yDs;
558: VecDot(ctx->W, ctx->W, &sDs);
559: ss_sum += sDs;
561: yy_sum += ctx->yy_rhistory[i];
562: }
563: }
564: } else if (1.0 == ctx->r_beta) {
565: /* Compute summations for scalar scaling */
566: yy_sum = 0; /* No safeguard required */
567: ys_sum = 0; /* No safeguard required */
568: ss_sum = 0; /* No safeguard required */
569: for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) {
570: VecPointwiseMult(ctx->V, ctx->Y[i], ctx->P);
571: VecDot(ctx->V, ctx->S[i], &yDs);
572: ys_sum += yDs;
574: VecDot(ctx->V, ctx->V, &yDy);
575: yy_sum += yDy;
577: ss_sum += ctx->ss_rhistory[i];
578: }
579: } else {
580: VecCopy(ctx->Q, ctx->P);
582: VecPow(ctx->P, ctx->r_beta);
583: VecPointwiseDivide(ctx->Q, ctx->P, ctx->Q);
585: /* Compute summations for scalar scaling */
586: yy_sum = 0; /* No safeguard required */
587: ys_sum = 0; /* No safeguard required */
588: ss_sum = 0; /* No safeguard required */
589: for (i = 0; i < PetscMin(ctx->nupdates, ctx->rescale_history); ++i) {
590: VecPointwiseMult(ctx->V, ctx->P, ctx->Y[i]);
591: VecPointwiseMult(ctx->W, ctx->Q, ctx->S[i]);
593: VecDot(ctx->V, ctx->V, &yDy);
594: VecDot(ctx->V, ctx->W, &yDs);
595: VecDot(ctx->W, ctx->W, &sDs);
597: yy_sum += yDy;
598: ys_sum += yDs;
599: ss_sum += sDs;
600: }
601: }
603: if (0.0 == ctx->r_alpha) {
604: /* Safeguard ys_sum */
605: if (0.0 == ys_sum) {
606: ys_sum = TAO_ZERO_SAFEGUARD;
607: }
609: sigmanew = ss_sum / ys_sum;
610: } else if (1.0 == ctx->r_alpha) {
611: /* Safeguard yy_sum */
612: if (0.0 == yy_sum) {
613: ys_sum = TAO_ZERO_SAFEGUARD;
614: }
616: sigmanew = ys_sum / yy_sum;
617: } else {
618: denom = 2*ctx->r_alpha*yy_sum;
620: /* Safeguard denom */
621: if (0.0 == denom) {
622: denom = TAO_ZERO_SAFEGUARD;
623: }
625: sigmanew = ((2*ctx->r_alpha-1)*ys_sum + PetscSqrtScalar((2*ctx->r_alpha-1)*(2*ctx->r_alpha-1)*ys_sum*ys_sum - 4*ctx->r_alpha*(ctx->r_alpha-1)*yy_sum*ss_sum)) / denom;
626: }
628: /* If Q has small values, then Q^(r_beta - 1) */
629: /* can have very large values. Hence, ys_sum */
630: /* and ss_sum can be infinity. In this case, */
631: /* sigmanew can either be not-a-number or infinity. */
633: if (PetscIsInfOrNanReal(sigmanew)) {
634: /* sigmanew is not-a-number; skip rescaling */
635: } else if (!sigmanew) {
636: /* sigmanew is zero; this is a bad case; skip rescaling */
637: } else {
638: /* sigmanew is positive */
639: VecScale(ctx->U, sigmanew);
640: }
641: break;
642: }
644: /* Modify for previous information */
645: switch(ctx->limitType) {
646: case MatLMVM_Limit_Average:
647: if (1.0 == ctx->mu) {
648: VecCopy(ctx->D, ctx->U);
649: } else if (ctx->mu) {
650: VecAXPBY(ctx->D,ctx->mu, 1.0-ctx->mu,ctx->U);
651: }
652: break;
654: case MatLMVM_Limit_Relative:
655: if (ctx->mu) {
656: /* P = (1-mu) * D */
657: VecAXPBY(ctx->P, 1.0-ctx->mu, 0.0, ctx->D);
658: /* Q = (1+mu) * D */
659: VecAXPBY(ctx->Q, 1.0+ctx->mu, 0.0, ctx->D);
660: VecMedian(ctx->P, ctx->U, ctx->Q, ctx->D);
661: }
662: break;
664: case MatLMVM_Limit_Absolute:
665: if (ctx->nu) {
666: VecCopy(ctx->P, ctx->D);
667: VecShift(ctx->P, -ctx->nu);
668: VecCopy(ctx->D, ctx->Q);
669: VecShift(ctx->Q, ctx->nu);
670: VecMedian(ctx->P, ctx->U, ctx->Q, ctx->P);
671: }
672: break;
674: default:
675: VecCopy(ctx->U, ctx->D);
676: break;
677: }
678: break;
679: }
680: PetscObjectDereference((PetscObject)ctx->Xprev);
681: PetscObjectDereference((PetscObject)ctx->Gprev);
682: ctx->Xprev = ctx->S[ctx->lm];
683: ctx->Gprev = ctx->Y[ctx->lm];
684: PetscObjectReference((PetscObject)ctx->S[ctx->lm]);
685: PetscObjectReference((PetscObject)ctx->Y[ctx->lm]);
687: } else {
688: ++ctx->nrejects;
689: }
690: }
692: ++ctx->iter;
693: VecCopy(x, ctx->Xprev);
694: VecCopy(g, ctx->Gprev);
695: return(0);
696: }
700: extern PetscErrorCode MatLMVMSetDelta(Mat m, PetscReal d)
701: {
702: MatLMVMCtx *ctx;
704: PetscBool same;
708: PetscObjectTypeCompare((PetscObject)m,MATSHELL,&same);
709: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix m is not type MatLMVM");
710: MatShellGetContext(m,(void**)&ctx);
711: ctx->delta = PetscAbsReal(d);
712: ctx->delta = PetscMax(ctx->delta_min, ctx->delta);
713: ctx->delta = PetscMin(ctx->delta_max, ctx->delta);
714: return(0);
715: }
719: extern PetscErrorCode MatLMVMSetScale(Mat m, Vec s)
720: {
721: MatLMVMCtx *ctx;
723: PetscBool same;
727: PetscObjectTypeCompare((PetscObject)m,MATSHELL,&same);
728: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix m is not type MatLMVM");
729: MatShellGetContext(m,(void**)&ctx);
731: VecDestroy(&ctx->scale);
732: if (s) {
733: VecDuplicate(s,&ctx->scale);
734: }
735: return(0);
736: }
740: extern PetscErrorCode MatLMVMGetRejects(Mat m, PetscInt *nrejects)
741: {
742: MatLMVMCtx *ctx;
744: PetscBool same;
748: PetscObjectTypeCompare((PetscObject)m,MATSHELL,&same);
749: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix m is not type MatLMVM");
750: MatShellGetContext(m,(void**)&ctx);
751: *nrejects = ctx->nrejects;
752: return(0);
753: }
757: extern PetscErrorCode MatLMVMSetH0(Mat m, Mat A)
758: {
760: return(0);
761: }
765: extern PetscErrorCode MatLMVMGetX0(Mat m, Vec x)
766: {
768: return(0);
769: }
773: extern PetscErrorCode MatLMVMSetPrev(Mat M, Vec x, Vec g)
774: {
775: MatLMVMCtx *ctx;
777: PetscBool same;
782: PetscObjectTypeCompare((PetscObject)M,MATSHELL,&same);
783: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix M is not type MatLMVM");
784: MatShellGetContext(M,(void**)&ctx);
785: if (ctx->nupdates == 0) {
786: MatLMVMUpdate(M,x,g);
787: } else {
788: VecCopy(x,ctx->Xprev);
789: VecCopy(g,ctx->Gprev);
790: /* TODO scaling specific terms */
791: }
792: return(0);
793: }
797: extern PetscErrorCode MatLMVMRefine(Mat coarse, Mat op, Vec fineX, Vec fineG)
798: {
800: PetscBool same;
807: PetscObjectTypeCompare((PetscObject)coarse,MATSHELL,&same);
808: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix m is not type MatLMVM");
809: PetscObjectTypeCompare((PetscObject)op,MATSHELL,&same);
810: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix m is not type MatLMVM");
811: return(0);
812: }
816: extern PetscErrorCode MatLMVMAllocateVectors(Mat m, Vec v)
817: {
819: MatLMVMCtx *ctx;
820: PetscBool same;
825: PetscObjectTypeCompare((PetscObject)m,MATSHELL,&same);
826: if (!same) SETERRQ(PETSC_COMM_SELF,1,"Matrix m is not type MatLMVM");
827: MatShellGetContext(m,(void**)&ctx);
829: /* Perform allocations */
830: VecDuplicateVecs(v,ctx->lm+1,&ctx->S);
831: VecDuplicateVecs(v,ctx->lm+1,&ctx->Y);
832: VecDuplicate(v,&ctx->D);
833: VecDuplicate(v,&ctx->U);
834: VecDuplicate(v,&ctx->V);
835: VecDuplicate(v,&ctx->W);
836: VecDuplicate(v,&ctx->P);
837: VecDuplicate(v,&ctx->Q);
838: ctx->allocated = PETSC_TRUE;
839: return(0);
840: }