Actual source code: qn.c
petsc-3.7.1 2016-05-15
1: #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2: #include <petscdm.h>
4: #define H(i,j) qn->dXdFmat[i*qn->m + j]
6: const char *const SNESQNScaleTypes[] = {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
7: const char *const SNESQNRestartTypes[] = {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
8: const char *const SNESQNTypes[] = {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0};
10: typedef struct {
11: Vec *U; /* Stored past states (vary from method to method) */
12: Vec *V; /* Stored past states (vary from method to method) */
13: PetscInt m; /* The number of kept previous steps */
14: PetscReal *lambda; /* The line search history of the method */
15: PetscReal *norm; /* norms of the steps */
16: PetscScalar *alpha, *beta;
17: PetscScalar *dXtdF, *dFtdX, *YtdX;
18: PetscBool singlereduction; /* Aggregated reduction implementation */
19: PetscScalar *dXdFmat; /* A matrix of values for dX_i dot dF_j */
20: PetscViewer monitor;
21: PetscReal powell_gamma; /* Powell angle restart condition */
22: PetscReal powell_downhill; /* Powell descent restart condition */
23: PetscReal scaling; /* scaling of H0 */
24: SNESQNType type; /* the type of quasi-newton method used */
25: SNESQNScaleType scale_type; /* the type of scaling used */
26: SNESQNRestartType restart_type; /* determine the frequency and type of restart conditions */
27: } SNES_QN;
31: PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D)
32: {
33: PetscErrorCode ierr;
34: SNES_QN *qn = (SNES_QN*)snes->data;
35: Vec W = snes->work[3];
36: Vec *U = qn->U;
37: PetscInt m = qn->m;
38: PetscInt k,i,j,lits,l = m;
39: PetscReal unorm,a,b;
40: PetscReal *lambda=qn->lambda;
41: PetscScalar gdot;
42: PetscReal udot;
45: if (it < m) l = it;
46: if (it > 0) {
47: k = (it-1)%l;
48: SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);
49: VecCopy(Xold,U[k]);
50: VecAXPY(U[k],-1.0,X);
51: if (qn->monitor) {
52: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
53: PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);
54: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
55: }
56: }
57: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
58: KSPSolve(snes->ksp,D,W);
59: SNESCheckKSPSolve(snes);
60: KSPGetIterationNumber(snes->ksp,&lits);
61: snes->linear_its += lits;
62: VecCopy(W,Y);
63: } else {
64: VecCopy(D,Y);
65: VecScale(Y,qn->scaling);
66: }
68: /* inward recursion starting at the first update and working forward */
69: for (i = 0; i < l-1; i++) {
70: j = (it+i-l)%l;
71: k = (it+i-l+1)%l;
72: VecNorm(U[j],NORM_2,&unorm);
73: VecDot(U[j],Y,&gdot);
74: unorm *= unorm;
75: udot = PetscRealPart(gdot);
76: a = (lambda[j]/lambda[k]);
77: b = -(1.-lambda[j]);
78: a *= udot/unorm;
79: b *= udot/unorm;
80: VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);
82: if (qn->monitor) {
83: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
84: PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,PetscRealPart(gdot));
85: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
86: }
87: }
88: if (l > 0) {
89: k = (it-1)%l;
90: VecDot(U[k],Y,&gdot);
91: VecNorm(U[k],NORM_2,&unorm);
92: unorm *= unorm;
93: udot = PetscRealPart(gdot);
94: a = unorm/(unorm-lambda[k]*udot);
95: b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
96: if (qn->monitor) {
97: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
98: PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);
99: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
100: }
101: VecAXPBY(Y,b,a,U[k]);
102: }
103: l = m;
104: if (it+1<m)l=it+1;
105: k = it%l;
106: if (qn->monitor) {
107: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
108: PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);
109: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
110: }
111: return(0);
112: }
116: PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
117: {
119: SNES_QN *qn = (SNES_QN*)snes->data;
120: Vec W = snes->work[3];
121: Vec *U = qn->U;
122: Vec *T = qn->V;
124: /* ksp thing for Jacobian scaling */
125: PetscInt h,k,j,i,lits;
126: PetscInt m = qn->m;
127: PetscScalar gdot,udot;
128: PetscInt l = m;
131: if (it < m) l = it;
132: VecCopy(D,Y);
133: if (l > 0) {
134: k = (it-1)%l;
135: SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);
136: VecCopy(Dold,U[k]);
137: VecAXPY(U[k],-1.0,D);
138: VecCopy(Xold,T[k]);
139: VecAXPY(T[k],-1.0,X);
140: }
142: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
143: KSPSolve(snes->ksp,Y,W);
144: SNESCheckKSPSolve(snes);
145: KSPGetIterationNumber(snes->ksp,&lits);
146: snes->linear_its += lits;
147: VecCopy(W,Y);
148: } else {
149: VecScale(Y,qn->scaling);
150: }
152: /* inward recursion starting at the first update and working forward */
153: if (l > 0) {
154: for (i = 0; i < l-1; i++) {
155: j = (it+i-l)%l;
156: k = (it+i-l+1)%l;
157: h = (it-1)%l;
158: VecDotBegin(U[j],U[h],&gdot);
159: VecDotBegin(U[j],U[j],&udot);
160: VecDotEnd(U[j],U[h],&gdot);
161: VecDotEnd(U[j],U[j],&udot);
162: VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);
163: VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);
164: if (qn->monitor) {
165: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
166: PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));
167: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
168: }
169: }
170: }
171: return(0);
172: }
176: PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
177: {
179: SNES_QN *qn = (SNES_QN*)snes->data;
180: Vec W = snes->work[3];
181: Vec *dX = qn->U;
182: Vec *dF = qn->V;
183: PetscScalar *alpha = qn->alpha;
184: PetscScalar *beta = qn->beta;
185: PetscScalar *dXtdF = qn->dXtdF;
186: PetscScalar *dFtdX = qn->dFtdX;
187: PetscScalar *YtdX = qn->YtdX;
189: /* ksp thing for Jacobian scaling */
190: PetscInt k,i,j,g,lits;
191: PetscInt m = qn->m;
192: PetscScalar t;
193: PetscInt l = m;
196: if (it < m) l = it;
197: VecCopy(D,Y);
198: if (it > 0) {
199: k = (it - 1) % l;
200: VecCopy(D,dF[k]);
201: VecAXPY(dF[k], -1.0, Dold);
202: VecCopy(X, dX[k]);
203: VecAXPY(dX[k], -1.0, Xold);
204: if (qn->singlereduction) {
205: VecMDotBegin(dF[k],l,dX,dXtdF);
206: VecMDotBegin(dX[k],l,dF,dFtdX);
207: VecMDotBegin(Y,l,dX,YtdX);
208: VecMDotEnd(dF[k],l,dX,dXtdF);
209: VecMDotEnd(dX[k],l,dF,dFtdX);
210: VecMDotEnd(Y,l,dX,YtdX);
211: for (j = 0; j < l; j++) {
212: H(k, j) = dFtdX[j];
213: H(j, k) = dXtdF[j];
214: }
215: /* copy back over to make the computation of alpha and beta easier */
216: for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
217: } else {
218: VecDot(dX[k], dF[k], &dXtdF[k]);
219: }
220: if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
221: SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);
222: }
223: }
225: /* outward recursion starting at iteration k's update and working back */
226: for (i=0; i<l; i++) {
227: k = (it-i-1)%l;
228: if (qn->singlereduction) {
229: /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
230: t = YtdX[k];
231: for (j=0; j<i; j++) {
232: g = (it-j-1)%l;
233: t -= alpha[g]*H(k, g);
234: }
235: alpha[k] = t/H(k,k);
236: } else {
237: VecDot(dX[k],Y,&t);
238: alpha[k] = t/dXtdF[k];
239: }
240: if (qn->monitor) {
241: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
242: PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha: %14.12e\n", it, k, PetscRealPart(alpha[k]));
243: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
244: }
245: VecAXPY(Y,-alpha[k],dF[k]);
246: }
248: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
249: KSPSolve(snes->ksp,Y,W);
250: SNESCheckKSPSolve(snes);
251: KSPGetIterationNumber(snes->ksp,&lits);
252: snes->linear_its += lits;
253: VecCopy(W, Y);
254: } else {
255: VecScale(Y, qn->scaling);
256: }
257: if (qn->singlereduction) {
258: VecMDot(Y,l,dF,YtdX);
259: }
260: /* inward recursion starting at the first update and working forward */
261: for (i = 0; i < l; i++) {
262: k = (it + i - l) % l;
263: if (qn->singlereduction) {
264: t = YtdX[k];
265: for (j = 0; j < i; j++) {
266: g = (it + j - l) % l;
267: t += (alpha[g] - beta[g])*H(g, k);
268: }
269: beta[k] = t / H(k, k);
270: } else {
271: VecDot(dF[k], Y, &t);
272: beta[k] = t / dXtdF[k];
273: }
274: VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
275: if (qn->monitor) {
276: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
277: PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));
278: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
279: }
280: }
281: return(0);
282: }
286: static PetscErrorCode SNESSolve_QN(SNES snes)
287: {
288: PetscErrorCode ierr;
289: SNES_QN *qn = (SNES_QN*) snes->data;
290: Vec X,Xold;
291: Vec F,W;
292: Vec Y,D,Dold;
293: PetscInt i, i_r;
294: PetscReal fnorm,xnorm,ynorm,gnorm;
295: SNESLineSearchReason lssucceed;
296: PetscBool powell,periodic;
297: PetscScalar DolddotD,DolddotDold;
298: SNESConvergedReason reason;
300: /* basically just a regular newton's method except for the application of the Jacobian */
303: if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
305: PetscCitationsRegister(SNESCitation,&SNEScite);
306: F = snes->vec_func; /* residual vector */
307: Y = snes->vec_sol_update; /* search direction generated by J^-1D*/
308: W = snes->work[3];
309: X = snes->vec_sol; /* solution vector */
310: Xold = snes->work[0];
312: /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
313: D = snes->work[1];
314: Dold = snes->work[2];
316: snes->reason = SNES_CONVERGED_ITERATING;
318: PetscObjectSAWsTakeAccess((PetscObject)snes);
319: snes->iter = 0;
320: snes->norm = 0.;
321: PetscObjectSAWsGrantAccess((PetscObject)snes);
323: if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
324: SNESApplyNPC(snes,X,NULL,F);
325: SNESGetConvergedReason(snes->pc,&reason);
326: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
327: snes->reason = SNES_DIVERGED_INNER;
328: return(0);
329: }
330: VecNorm(F,NORM_2,&fnorm);
331: } else {
332: if (!snes->vec_func_init_set) {
333: SNESComputeFunction(snes,X,F);
334: } else snes->vec_func_init_set = PETSC_FALSE;
336: VecNorm(F,NORM_2,&fnorm);
337: SNESCheckFunctionNorm(snes,fnorm);
338: }
339: if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
340: SNESApplyNPC(snes,X,F,D);
341: SNESGetConvergedReason(snes->pc,&reason);
342: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
343: snes->reason = SNES_DIVERGED_INNER;
344: return(0);
345: }
346: } else {
347: VecCopy(F,D);
348: }
350: PetscObjectSAWsTakeAccess((PetscObject)snes);
351: snes->norm = fnorm;
352: PetscObjectSAWsGrantAccess((PetscObject)snes);
353: SNESLogConvergenceHistory(snes,fnorm,0);
354: SNESMonitor(snes,0,fnorm);
356: /* test convergence */
357: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
358: if (snes->reason) return(0);
360: if (snes->pc && snes->pcside == PC_RIGHT) {
361: PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);
362: SNESSolve(snes->pc,snes->vec_rhs,X);
363: PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);
364: SNESGetConvergedReason(snes->pc,&reason);
365: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
366: snes->reason = SNES_DIVERGED_INNER;
367: return(0);
368: }
369: SNESGetNPCFunction(snes,F,&fnorm);
370: VecCopy(F,D);
371: }
373: /* scale the initial update */
374: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
375: SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
376: KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);
377: }
379: for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
380: if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
381: PetscScalar ff,xf;
382: VecCopy(Dold,Y);
383: VecCopy(Xold,W);
384: VecAXPY(Y,-1.0,D);
385: VecAXPY(W,-1.0,X);
386: VecDotBegin(Y,Y,&ff);
387: VecDotBegin(W,Y,&xf);
388: VecDotEnd(Y,Y,&ff);
389: VecDotEnd(W,Y,&xf);
390: qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
391: }
392: switch (qn->type) {
393: case SNES_QN_BADBROYDEN:
394: SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);
395: break;
396: case SNES_QN_BROYDEN:
397: SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);
398: break;
399: case SNES_QN_LBFGS:
400: SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);
401: break;
402: }
403: /* line search for lambda */
404: ynorm = 1; gnorm = fnorm;
405: VecCopy(D, Dold);
406: VecCopy(X, Xold);
407: SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
408: if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
409: SNESLineSearchGetReason(snes->linesearch, &lssucceed);
410: SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
411: if (lssucceed) {
412: if (++snes->numFailures >= snes->maxFailures) {
413: snes->reason = SNES_DIVERGED_LINE_SEARCH;
414: break;
415: }
416: }
417: if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
418: SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);
419: }
421: /* convergence monitoring */
422: PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);
424: if (snes->pc && snes->pcside == PC_RIGHT) {
425: PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);
426: SNESSolve(snes->pc,snes->vec_rhs,X);
427: PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);
428: SNESGetConvergedReason(snes->pc,&reason);
429: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
430: snes->reason = SNES_DIVERGED_INNER;
431: return(0);
432: }
433: SNESGetNPCFunction(snes,F,&fnorm);
434: }
436: SNESSetIterationNumber(snes, i+1);
437: snes->norm = fnorm;
439: SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
440: SNESMonitor(snes,snes->iter,snes->norm);
441: /* set parameter for default relative tolerance convergence test */
442: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
443: if (snes->reason) return(0);
444: if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
445: SNESApplyNPC(snes,X,F,D);
446: SNESGetConvergedReason(snes->pc,&reason);
447: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
448: snes->reason = SNES_DIVERGED_INNER;
449: return(0);
450: }
451: } else {
452: VecCopy(F, D);
453: }
454: powell = PETSC_FALSE;
455: if (qn->restart_type == SNES_QN_RESTART_POWELL) {
456: /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
457: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
458: MatMult(snes->jacobian_pre,Dold,W);
459: } else {
460: VecCopy(Dold,W);
461: }
462: VecDotBegin(W, Dold, &DolddotDold);
463: VecDotBegin(W, D, &DolddotD);
464: VecDotEnd(W, Dold, &DolddotDold);
465: VecDotEnd(W, D, &DolddotD);
466: if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
467: }
468: periodic = PETSC_FALSE;
469: if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
470: if (i_r>qn->m-1) periodic = PETSC_TRUE;
471: }
472: /* restart if either powell or periodic restart is satisfied. */
473: if (powell || periodic) {
474: if (qn->monitor) {
475: PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
476: PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);
477: PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
478: }
479: i_r = -1;
480: /* general purpose update */
481: if (snes->ops->update) {
482: (*snes->ops->update)(snes, snes->iter);
483: }
484: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
485: SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
486: }
487: }
488: /* general purpose update */
489: if (snes->ops->update) {
490: (*snes->ops->update)(snes, snes->iter);
491: }
492: }
493: if (i == snes->max_its) {
494: PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);
495: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
496: }
497: return(0);
498: }
502: static PetscErrorCode SNESSetUp_QN(SNES snes)
503: {
504: SNES_QN *qn = (SNES_QN*)snes->data;
506: DM dm;
510: if (!snes->vec_sol) {
511: SNESGetDM(snes,&dm);
512: DMCreateGlobalVector(dm,&snes->vec_sol);
513: }
515: VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);
516: if (qn->type != SNES_QN_BROYDEN) VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);
517: PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);
519: if (qn->singlereduction) {
520: PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);
521: }
522: SNESSetWorkVecs(snes,4);
523: /* set method defaults */
524: if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
525: if (qn->type == SNES_QN_BADBROYDEN) {
526: qn->scale_type = SNES_QN_SCALE_NONE;
527: } else {
528: qn->scale_type = SNES_QN_SCALE_SHANNO;
529: }
530: }
531: if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
532: if (qn->type == SNES_QN_LBFGS) {
533: qn->restart_type = SNES_QN_RESTART_POWELL;
534: } else {
535: qn->restart_type = SNES_QN_RESTART_PERIODIC;
536: }
537: }
539: if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
540: SNESSetUpMatrices(snes);
541: }
542: if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
543: return(0);
544: }
548: static PetscErrorCode SNESReset_QN(SNES snes)
549: {
551: SNES_QN *qn;
554: if (snes->data) {
555: qn = (SNES_QN*)snes->data;
556: if (qn->U) {
557: VecDestroyVecs(qn->m, &qn->U);
558: }
559: if (qn->V) {
560: VecDestroyVecs(qn->m, &qn->V);
561: }
562: if (qn->singlereduction) {
563: PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);
564: }
565: PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);
566: }
567: return(0);
568: }
572: static PetscErrorCode SNESDestroy_QN(SNES snes)
573: {
577: SNESReset_QN(snes);
578: PetscFree(snes->data);
579: PetscObjectComposeFunction((PetscObject)snes,"",NULL);
580: return(0);
581: }
585: static PetscErrorCode SNESSetFromOptions_QN(PetscOptionItems *PetscOptionsObject,SNES snes)
586: {
588: PetscErrorCode ierr;
589: SNES_QN *qn = (SNES_QN*)snes->data;
590: PetscBool monflg = PETSC_FALSE,flg;
591: SNESLineSearch linesearch;
592: SNESQNRestartType rtype = qn->restart_type;
593: SNESQNScaleType stype = qn->scale_type;
594: SNESQNType qtype = qn->type;
597: PetscOptionsHead(PetscOptionsObject,"SNES QN options");
598: PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);
599: PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance", "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);
600: PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance", "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);
601: PetscOptionsBool("-snes_qn_monitor", "Monitor for the QN methods", "SNESQN", monflg, &monflg, NULL);
602: PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions", "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);
603: PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);
604: if (flg) SNESQNSetScaleType(snes,stype);
606: PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);
607: if (flg) SNESQNSetRestartType(snes,rtype);
609: PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);
610: if (flg) {SNESQNSetType(snes,qtype);}
611: PetscOptionsTail();
612: if (!snes->linesearch) {
613: SNESGetLineSearch(snes, &linesearch);
614: if (qn->type == SNES_QN_LBFGS) {
615: SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
616: } else if (qn->type == SNES_QN_BROYDEN) {
617: SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);
618: } else {
619: SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
620: }
621: }
622: if (monflg) {
623: qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
624: }
625: return(0);
626: }
630: static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
631: {
632: SNES_QN *qn = (SNES_QN*)snes->data;
633: PetscBool iascii;
637: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
638: if (iascii) {
639: PetscViewerASCIIPrintf(viewer," QN type is %s, restart type is %s, scale type is %s\n",SNESQNTypes[qn->type],SNESQNRestartTypes[qn->restart_type],SNESQNScaleTypes[qn->scale_type]);
640: PetscViewerASCIIPrintf(viewer," Stored subspace size: %d\n", qn->m);
641: if (qn->singlereduction) {
642: PetscViewerASCIIPrintf(viewer," Using the single reduction variant.\n");
643: }
644: }
645: return(0);
646: }
650: /*@
651: SNESQNSetRestartType - Sets the restart type for SNESQN.
653: Logically Collective on SNES
655: Input Parameters:
656: + snes - the iterative context
657: - rtype - restart type
659: Options Database:
660: + -snes_qn_restart_type <powell,periodic,none> - set the restart type
661: - -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic
663: Level: intermediate
665: SNESQNRestartTypes:
666: + SNES_QN_RESTART_NONE - never restart
667: . SNES_QN_RESTART_POWELL - restart based upon descent criteria
668: - SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations
670: .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType
671: @*/
672: PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
673: {
678: PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));
679: return(0);
680: }
684: /*@
685: SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.
687: Logically Collective on SNES
689: Input Parameters:
690: + snes - the iterative context
691: - stype - scale type
693: Options Database:
694: . -snes_qn_scale_type <shanno,none,linesearch,jacobian>
696: Level: intermediate
698: SNESQNScaleTypes:
699: + SNES_QN_SCALE_NONE - don't scale the problem
700: . SNES_QN_SCALE_SHANNO - use shanno scaling
701: . SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
702: - SNES_QN_SCALE_JACOBIAN - scale by solving a linear system coming from the Jacobian you provided with SNESSetJacobian() computed at the first iteration
703: of QN and at ever restart.
705: .keywords: scaling type
707: .seealso: SNES, SNESQN, SNESLineSearch, SNESQNScaleType, SNESetJacobian()
708: @*/
710: PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
711: {
716: PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));
717: return(0);
718: }
722: PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
723: {
724: SNES_QN *qn = (SNES_QN*)snes->data;
727: qn->scale_type = stype;
728: return(0);
729: }
733: PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
734: {
735: SNES_QN *qn = (SNES_QN*)snes->data;
738: qn->restart_type = rtype;
739: return(0);
740: }
744: /*@
745: SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.
747: Logically Collective on SNES
749: Input Parameters:
750: + snes - the iterative context
751: - qtype - variant type
753: Options Database:
754: . -snes_qn_type <lbfgs,broyden,badbroyden>
756: Level: beginner
758: SNESQNTypes:
759: + SNES_QN_LBFGS - LBFGS variant
760: . SNES_QN_BROYDEN - Broyden variant
761: - SNES_QN_BADBROYDEN - Bad Broyden variant
763: .keywords: SNES, SNESQN, type, set, SNESQNType
764: @*/
766: PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
767: {
772: PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));
773: return(0);
774: }
778: PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
779: {
780: SNES_QN *qn = (SNES_QN*)snes->data;
783: qn->type = qtype;
784: return(0);
785: }
787: /* -------------------------------------------------------------------------- */
788: /*MC
789: SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.
791: Options Database:
793: + -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
794: + -snes_qn_restart_type <powell,periodic,none> - set the restart type
795: . -snes_qn_powell_angle - Angle condition for restart.
796: . -snes_qn_powell_descent - Descent condition for restart.
797: . -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
798: . -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian
799: . -snes_linesearch_type <cp, l2, basic> - Type of line search.
800: - -snes_qn_monitor - Monitors the quasi-newton Jacobian.
802: Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
803: previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
804: updates.
806: When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied. The first of
807: these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
808: iteration as the current iteration's values when constructing the approximate Jacobian. The second, composed,
809: perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.
811: References:
812: + 1. - Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.
813: . 2. - R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi Newton Matrices and their use in Limited Memory Methods,
814: Technical Report, Northwestern University, June 1992.
815: . 3. - Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
816: Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.
817: - 4. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
818: SIAM Review, 57(4), 2015
820: Level: beginner
822: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR
824: M*/
827: PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
828: {
830: SNES_QN *qn;
833: snes->ops->setup = SNESSetUp_QN;
834: snes->ops->solve = SNESSolve_QN;
835: snes->ops->destroy = SNESDestroy_QN;
836: snes->ops->setfromoptions = SNESSetFromOptions_QN;
837: snes->ops->view = SNESView_QN;
838: snes->ops->reset = SNESReset_QN;
840: snes->pcside = PC_LEFT;
842: snes->usespc = PETSC_TRUE;
843: snes->usesksp = PETSC_FALSE;
845: if (!snes->tolerancesset) {
846: snes->max_funcs = 30000;
847: snes->max_its = 10000;
848: }
850: PetscNewLog(snes,&qn);
851: snes->data = (void*) qn;
852: qn->m = 10;
853: qn->scaling = 1.0;
854: qn->U = NULL;
855: qn->V = NULL;
856: qn->dXtdF = NULL;
857: qn->dFtdX = NULL;
858: qn->dXdFmat = NULL;
859: qn->monitor = NULL;
860: qn->singlereduction = PETSC_TRUE;
861: qn->powell_gamma = 0.9999;
862: qn->scale_type = SNES_QN_SCALE_DEFAULT;
863: qn->restart_type = SNES_QN_RESTART_DEFAULT;
864: qn->type = SNES_QN_LBFGS;
866: PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);
867: PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);
868: PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);
869: return(0);
870: }