Actual source code: qn.c

petsc-3.6.0 2015-06-09
Report Typos and Errors
  1: #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/

  3: #define H(i,j)  qn->dXdFmat[i*qn->m + j]

  5: const char *const SNESQNScaleTypes[] =        {"DEFAULT","NONE","SHANNO","LINESEARCH","JACOBIAN","SNESQNScaleType","SNES_QN_SCALING_",0};
  6: const char *const SNESQNRestartTypes[] =      {"DEFAULT","NONE","POWELL","PERIODIC","SNESQNRestartType","SNES_QN_RESTART_",0};
  7: const char *const SNESQNTypes[] =             {"LBFGS","BROYDEN","BADBROYDEN","SNESQNType","SNES_QN_",0};

  9: typedef struct {
 10:   Vec               *U;                   /* Stored past states (vary from method to method) */
 11:   Vec               *V;                   /* Stored past states (vary from method to method) */
 12:   PetscInt          m;                    /* The number of kept previous steps */
 13:   PetscReal         *lambda;              /* The line search history of the method */
 14:   PetscReal         *norm;                /* norms of the steps */
 15:   PetscScalar       *alpha, *beta;
 16:   PetscScalar       *dXtdF, *dFtdX, *YtdX;
 17:   PetscBool         singlereduction;      /* Aggregated reduction implementation */
 18:   PetscScalar       *dXdFmat;             /* A matrix of values for dX_i dot dF_j */
 19:   PetscViewer       monitor;
 20:   PetscReal         powell_gamma;         /* Powell angle restart condition */
 21:   PetscReal         powell_downhill;      /* Powell descent restart condition */
 22:   PetscReal         scaling;              /* scaling of H0 */
 23:   SNESQNType        type;                 /* the type of quasi-newton method used */
 24:   SNESQNScaleType   scale_type;           /* the type of scaling used */
 25:   SNESQNRestartType restart_type;         /* determine the frequency and type of restart conditions */
 26: } SNES_QN;

 30: PetscErrorCode SNESQNApply_Broyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D)
 31: {
 32:   PetscErrorCode     ierr;
 33:   SNES_QN            *qn = (SNES_QN*)snes->data;
 34:   Vec                W   = snes->work[3];
 35:   Vec                *U  = qn->U;
 36:   PetscInt           m = qn->m;
 37:   PetscInt           k,i,j,lits,l = m;
 38:   PetscReal          unorm,a,b;
 39:   PetscReal          *lambda=qn->lambda;
 40:   PetscScalar        gdot;
 41:   PetscReal          udot;

 44:   if (it < m) l = it;
 45:   if (it > 0) {
 46:     k = (it-1)%l;
 47:     SNESLineSearchGetLambda(snes->linesearch,&lambda[k]);
 48:     VecCopy(Xold,U[k]);
 49:     VecAXPY(U[k],-1.0,X);
 50:     if (qn->monitor) {
 51:       PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
 52:       PetscViewerASCIIPrintf(qn->monitor, "scaling vector %d of %d by lambda: %14.12e \n",k,l,lambda[k]);
 53:       PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
 54:     }
 55:   }
 56:   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
 57:     KSPSolve(snes->ksp,D,W);
 58:     SNESCheckKSPSolve(snes);
 59:     KSPGetIterationNumber(snes->ksp,&lits);
 60:     snes->linear_its += lits;
 61:     VecCopy(W,Y);
 62:   } else {
 63:     VecCopy(D,Y);
 64:     VecScale(Y,qn->scaling);
 65:   }

 67:   /* inward recursion starting at the first update and working forward */
 68:   for (i = 0; i < l-1; i++) {
 69:     j = (it+i-l)%l;
 70:     k = (it+i-l+1)%l;
 71:     VecNorm(U[j],NORM_2,&unorm);
 72:     VecDot(U[j],Y,&gdot);
 73:     unorm *= unorm;
 74:     udot = PetscRealPart(gdot);
 75:     a = (lambda[j]/lambda[k]);
 76:     b = -(1.-lambda[j]);
 77:     a *= udot/unorm;
 78:     b *= udot/unorm;
 79:     VecAXPBYPCZ(Y,a,b,1.,U[k],U[j]);

 81:     if (qn->monitor) {
 82:       PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
 83:       PetscViewerASCIIPrintf(qn->monitor, "using vector %d and %d, gdot: %14.12e\n",k,j,PetscRealPart(gdot));
 84:       PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
 85:     }
 86:   }
 87:   if (l > 0) {
 88:     k = (it-1)%l;
 89:     VecDot(U[k],Y,&gdot);
 90:     VecNorm(U[k],NORM_2,&unorm);
 91:     unorm *= unorm;
 92:     udot = PetscRealPart(gdot);
 93:     a = unorm/(unorm-lambda[k]*udot);
 94:     b = -(1.-lambda[k])*udot/(unorm-lambda[k]*udot);
 95:     if (qn->monitor) {
 96:       PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
 97:       PetscViewerASCIIPrintf(qn->monitor, "using vector %d: a: %14.12e b: %14.12e \n",k,a,b);
 98:       PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
 99:     }
100:     VecAXPBY(Y,b,a,U[k]);
101:   }
102:   l = m;
103:   if (it+1<m)l=it+1;
104:   k = it%l;
105:   if (qn->monitor) {
106:     PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
107:     PetscViewerASCIIPrintf(qn->monitor, "setting vector %d of %d\n",k,l);
108:     PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
109:   }
110:   return(0);
111: }

115: PetscErrorCode SNESQNApply_BadBroyden(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
116: {
118:   SNES_QN        *qn = (SNES_QN*)snes->data;
119:   Vec            W   = snes->work[3];
120:   Vec            *U  = qn->U;
121:   Vec            *T  = qn->V;

123:   /* ksp thing for Jacobian scaling */
124:   PetscInt           h,k,j,i,lits;
125:   PetscInt           m = qn->m;
126:   PetscScalar        gdot,udot;
127:   PetscInt           l = m;

130:   if (it < m) l = it;
131:   VecCopy(D,Y);
132:   if (l > 0) {
133:     k    = (it-1)%l;
134:     SNESLineSearchGetLambda(snes->linesearch,&qn->lambda[k]);
135:     VecCopy(Dold,U[k]);
136:     VecAXPY(U[k],-1.0,D);
137:     VecCopy(Xold,T[k]);
138:     VecAXPY(T[k],-1.0,X);
139:   }

141:   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
142:     KSPSolve(snes->ksp,Y,W);
143:     SNESCheckKSPSolve(snes);
144:     KSPGetIterationNumber(snes->ksp,&lits);
145:     snes->linear_its += lits;
146:     VecCopy(W,Y);
147:   } else {
148:     VecScale(Y,qn->scaling);
149:   }

151:   /* inward recursion starting at the first update and working forward */
152:   if (l > 0) {
153:     for (i = 0; i < l-1; i++) {
154:       j    = (it+i-l)%l;
155:       k    = (it+i-l+1)%l;
156:       h    = (it-1)%l;
157:       VecDotBegin(U[j],U[h],&gdot);
158:       VecDotBegin(U[j],U[j],&udot);
159:       VecDotEnd(U[j],U[h],&gdot);
160:       VecDotEnd(U[j],U[j],&udot);
161:       VecAXPY(Y,PetscRealPart(gdot)/PetscRealPart(udot),T[k]);
162:       VecAXPY(Y,-(1.-qn->lambda[k])*PetscRealPart(gdot)/PetscRealPart(udot),T[j]);
163:       if (qn->monitor) {
164:         PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
165:         PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d gdot: %14.12e\n", it, k, PetscRealPart(gdot));
166:         PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
167:       }
168:     }
169:   }
170:   return(0);
171: }

175: PetscErrorCode SNESQNApply_LBFGS(SNES snes,PetscInt it,Vec Y,Vec X,Vec Xold,Vec D,Vec Dold)
176: {
178:   SNES_QN        *qn    = (SNES_QN*)snes->data;
179:   Vec            W      = snes->work[3];
180:   Vec            *dX    = qn->U;
181:   Vec            *dF    = qn->V;
182:   PetscScalar    *alpha = qn->alpha;
183:   PetscScalar    *beta  = qn->beta;
184:   PetscScalar    *dXtdF = qn->dXtdF;
185:   PetscScalar    *dFtdX = qn->dFtdX;
186:   PetscScalar    *YtdX  = qn->YtdX;

188:   /* ksp thing for Jacobian scaling */
189:   PetscInt           k,i,j,g,lits;
190:   PetscInt           m = qn->m;
191:   PetscScalar        t;
192:   PetscInt           l = m;

195:   if (it < m) l = it;
196:   VecCopy(D,Y);
197:   if (it > 0) {
198:     k    = (it - 1) % l;
199:     VecCopy(D,dF[k]);
200:     VecAXPY(dF[k], -1.0, Dold);
201:     VecCopy(X, dX[k]);
202:     VecAXPY(dX[k], -1.0, Xold);
203:     if (qn->singlereduction) {
204:       VecMDotBegin(dF[k],l,dX,dXtdF);
205:       VecMDotBegin(dX[k],l,dF,dFtdX);
206:       VecMDotBegin(Y,l,dX,YtdX);
207:       VecMDotEnd(dF[k],l,dX,dXtdF);
208:       VecMDotEnd(dX[k],l,dF,dFtdX);
209:       VecMDotEnd(Y,l,dX,YtdX);
210:       for (j = 0; j < l; j++) {
211:         H(k, j) = dFtdX[j];
212:         H(j, k) = dXtdF[j];
213:       }
214:       /* copy back over to make the computation of alpha and beta easier */
215:       for (j = 0; j < l; j++) dXtdF[j] = H(j, j);
216:     } else {
217:       VecDot(dX[k], dF[k], &dXtdF[k]);
218:     }
219:     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
220:       SNESLineSearchGetLambda(snes->linesearch,&qn->scaling);
221:     }
222:   }

224:   /* outward recursion starting at iteration k's update and working back */
225:   for (i=0; i<l; i++) {
226:     k = (it-i-1)%l;
227:     if (qn->singlereduction) {
228:       /* construct t = dX[k] dot Y as Y_0 dot dX[k] + sum(-alpha[j]dX[k]dF[j]) */
229:       t = YtdX[k];
230:       for (j=0; j<i; j++) {
231:         g  = (it-j-1)%l;
232:         t -= alpha[g]*H(k, g);
233:       }
234:       alpha[k] = t/H(k,k);
235:     } else {
236:       VecDot(dX[k],Y,&t);
237:       alpha[k] = t/dXtdF[k];
238:     }
239:     if (qn->monitor) {
240:       PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
241:       PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha:        %14.12e\n", it, k, PetscRealPart(alpha[k]));
242:       PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
243:     }
244:     VecAXPY(Y,-alpha[k],dF[k]);
245:   }

247:   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
248:     KSPSolve(snes->ksp,Y,W);
249:     SNESCheckKSPSolve(snes);
250:     KSPGetIterationNumber(snes->ksp,&lits);
251:     snes->linear_its += lits;
252:     VecCopy(W, Y);
253:   } else {
254:     VecScale(Y, qn->scaling);
255:   }
256:   if (qn->singlereduction) {
257:     VecMDot(Y,l,dF,YtdX);
258:   }
259:   /* inward recursion starting at the first update and working forward */
260:   for (i = 0; i < l; i++) {
261:     k = (it + i - l) % l;
262:     if (qn->singlereduction) {
263:       t = YtdX[k];
264:       for (j = 0; j < i; j++) {
265:         g  = (it + j - l) % l;
266:         t += (alpha[g] - beta[g])*H(g, k);
267:       }
268:       beta[k] = t / H(k, k);
269:     } else {
270:       VecDot(dF[k], Y, &t);
271:       beta[k] = t / dXtdF[k];
272:     }
273:     VecAXPY(Y, (alpha[k] - beta[k]), dX[k]);
274:     if (qn->monitor) {
275:       PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
276:       PetscViewerASCIIPrintf(qn->monitor, "it: %d k: %d alpha - beta: %14.12e\n", it, k, PetscRealPart(alpha[k] - beta[k]));
277:       PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
278:     }
279:   }
280:   return(0);
281: }

285: static PetscErrorCode SNESSolve_QN(SNES snes)
286: {
287:   PetscErrorCode       ierr;
288:   SNES_QN              *qn = (SNES_QN*) snes->data;
289:   Vec                  X,Xold;
290:   Vec                  F,W;
291:   Vec                  Y,D,Dold;
292:   PetscInt             i, i_r;
293:   PetscReal            fnorm,xnorm,ynorm,gnorm;
294:   SNESLineSearchReason lssucceed;
295:   PetscBool            powell,periodic;
296:   PetscScalar          DolddotD,DolddotDold;
297:   SNESConvergedReason  reason;

299:   /* basically just a regular newton's method except for the application of the Jacobian */


303:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
304:     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
305:   }

307:   PetscCitationsRegister(SNESCitation,&SNEScite);
308:   F    = snes->vec_func;                /* residual vector */
309:   Y    = snes->vec_sol_update;          /* search direction generated by J^-1D*/
310:   W    = snes->work[3];
311:   X    = snes->vec_sol;                 /* solution vector */
312:   Xold = snes->work[0];

314:   /* directions generated by the preconditioned problem with F_pre = F or x - M(x, b) */
315:   D    = snes->work[1];
316:   Dold = snes->work[2];

318:   snes->reason = SNES_CONVERGED_ITERATING;

320:   PetscObjectSAWsTakeAccess((PetscObject)snes);
321:   snes->iter = 0;
322:   snes->norm = 0.;
323:   PetscObjectSAWsGrantAccess((PetscObject)snes);

325:   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
326:     SNESApplyNPC(snes,X,NULL,F);
327:     SNESGetConvergedReason(snes->pc,&reason);
328:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
329:       snes->reason = SNES_DIVERGED_INNER;
330:       return(0);
331:     }
332:     VecNorm(F,NORM_2,&fnorm);
333:   } else {
334:     if (!snes->vec_func_init_set) {
335:       SNESComputeFunction(snes,X,F);
336:     } else snes->vec_func_init_set = PETSC_FALSE;

338:     VecNorm(F,NORM_2,&fnorm);
339:     SNESCheckFunctionNorm(snes,fnorm);
340:   }
341:   if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
342:       SNESApplyNPC(snes,X,F,D);
343:       SNESGetConvergedReason(snes->pc,&reason);
344:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
345:         snes->reason = SNES_DIVERGED_INNER;
346:         return(0);
347:       }
348:   } else {
349:     VecCopy(F,D);
350:   }

352:   PetscObjectSAWsTakeAccess((PetscObject)snes);
353:   snes->norm = fnorm;
354:   PetscObjectSAWsGrantAccess((PetscObject)snes);
355:   SNESLogConvergenceHistory(snes,fnorm,0);
356:   SNESMonitor(snes,0,fnorm);

358:   /* test convergence */
359:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
360:   if (snes->reason) return(0);

362:   if (snes->pc && snes->pcside == PC_RIGHT) {
363:     PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);
364:     SNESSolve(snes->pc,snes->vec_rhs,X);
365:     PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);
366:     SNESGetConvergedReason(snes->pc,&reason);
367:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
368:       snes->reason = SNES_DIVERGED_INNER;
369:       return(0);
370:     }
371:     SNESGetNPCFunction(snes,F,&fnorm);
372:     VecCopy(F,D);
373:   }

375:   /* scale the initial update */
376:   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
377:     SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
378:   }

380:   for (i = 0, i_r = 0; i < snes->max_its; i++, i_r++) {
381:     if (qn->scale_type == SNES_QN_SCALE_SHANNO && i_r > 0) {
382:       PetscScalar ff,xf;
383:       VecCopy(Dold,Y);
384:       VecCopy(Xold,W);
385:       VecAXPY(Y,-1.0,D);
386:       VecAXPY(W,-1.0,X);
387:       VecDotBegin(Y,Y,&ff);
388:       VecDotBegin(W,Y,&xf);
389:       VecDotEnd(Y,Y,&ff);
390:       VecDotEnd(W,Y,&xf);
391:       qn->scaling = PetscRealPart(xf)/PetscRealPart(ff);
392:     }
393:     switch (qn->type) {
394:     case SNES_QN_BADBROYDEN:
395:       SNESQNApply_BadBroyden(snes,i_r,Y,X,Xold,D,Dold);
396:       break;
397:     case SNES_QN_BROYDEN:
398:       SNESQNApply_Broyden(snes,i_r,Y,X,Xold,D);
399:       break;
400:     case SNES_QN_LBFGS:
401:       SNESQNApply_LBFGS(snes,i_r,Y,X,Xold,D,Dold);
402:       break;
403:     }
404:     /* line search for lambda */
405:     ynorm = 1; gnorm = fnorm;
406:     VecCopy(D, Dold);
407:     VecCopy(X, Xold);
408:     SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
409:     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
410:     SNESLineSearchGetReason(snes->linesearch, &lssucceed);
411:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
412:     if (lssucceed) {
413:       if (++snes->numFailures >= snes->maxFailures) {
414:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
415:         break;
416:       }
417:     }
418:     if (qn->scale_type == SNES_QN_SCALE_LINESEARCH) {
419:       SNESLineSearchGetLambda(snes->linesearch, &qn->scaling);
420:     }

422:     /* convergence monitoring */
423:     PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);

425:     if (snes->pc && snes->pcside == PC_RIGHT) {
426:       PetscLogEventBegin(SNES_NPCSolve,snes->pc,X,0,0);
427:       SNESSolve(snes->pc,snes->vec_rhs,X);
428:       PetscLogEventEnd(SNES_NPCSolve,snes->pc,X,0,0);
429:       SNESGetConvergedReason(snes->pc,&reason);
430:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
431:         snes->reason = SNES_DIVERGED_INNER;
432:         return(0);
433:       }
434:       SNESGetNPCFunction(snes,F,&fnorm);
435:     }

437:     SNESSetIterationNumber(snes, i+1);
438:     snes->norm = fnorm;

440:     SNESLogConvergenceHistory(snes,snes->norm,snes->iter);
441:     SNESMonitor(snes,snes->iter,snes->norm);
442:     /* set parameter for default relative tolerance convergence test */
443:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
444:     if (snes->reason) return(0);
445:     if (snes->pc && snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
446:       SNESApplyNPC(snes,X,F,D);
447:       SNESGetConvergedReason(snes->pc,&reason);
448:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
449:         snes->reason = SNES_DIVERGED_INNER;
450:         return(0);
451:       }
452:     } else {
453:       VecCopy(F, D);
454:     }
455:     powell = PETSC_FALSE;
456:     if (qn->restart_type == SNES_QN_RESTART_POWELL) {
457:       /* check restart by Powell's Criterion: |F^T H_0 Fold| > 0.2 * |Fold^T H_0 Fold| */
458:       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
459:         MatMult(snes->jacobian_pre,Dold,W);
460:       } else {
461:         VecCopy(Dold,W);
462:       }
463:       VecDotBegin(W, Dold, &DolddotDold);
464:       VecDotBegin(W, D, &DolddotD);
465:       VecDotEnd(W, Dold, &DolddotDold);
466:       VecDotEnd(W, D, &DolddotD);
467:       if (PetscAbs(PetscRealPart(DolddotD)) > qn->powell_gamma*PetscAbs(PetscRealPart(DolddotDold))) powell = PETSC_TRUE;
468:     }
469:     periodic = PETSC_FALSE;
470:     if (qn->restart_type == SNES_QN_RESTART_PERIODIC) {
471:       if (i_r>qn->m-1) periodic = PETSC_TRUE;
472:     }
473:     /* restart if either powell or periodic restart is satisfied. */
474:     if (powell || periodic) {
475:       if (qn->monitor) {
476:         PetscViewerASCIIAddTab(qn->monitor,((PetscObject)snes)->tablevel+2);
477:         PetscViewerASCIIPrintf(qn->monitor, "restart! |%14.12e| > %4.2f*|%14.12e| or i_r = %d\n", PetscRealPart(DolddotD), qn->powell_gamma, PetscRealPart(DolddotDold), i_r);
478:         PetscViewerASCIISubtractTab(qn->monitor,((PetscObject)snes)->tablevel+2);
479:       }
480:       i_r = -1;
481:       /* general purpose update */
482:       if (snes->ops->update) {
483:         (*snes->ops->update)(snes, snes->iter);
484:       }
485:       if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
486:         SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
487:       }
488:     }
489:     /* general purpose update */
490:     if (snes->ops->update) {
491:       (*snes->ops->update)(snes, snes->iter);
492:     }
493:   }
494:   if (i == snes->max_its) {
495:     PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", snes->max_its);
496:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
497:   }
498:   return(0);
499: }

503: static PetscErrorCode SNESSetUp_QN(SNES snes)
504: {
505:   SNES_QN        *qn = (SNES_QN*)snes->data;

509:   VecDuplicateVecs(snes->vec_sol, qn->m, &qn->U);
510:   if (qn->type != SNES_QN_BROYDEN) VecDuplicateVecs(snes->vec_sol, qn->m, &qn->V);
511:   PetscMalloc4(qn->m,&qn->alpha,qn->m,&qn->beta,qn->m,&qn->dXtdF,qn->m,&qn->lambda);

513:   if (qn->singlereduction) {
514:     PetscMalloc3(qn->m*qn->m,&qn->dXdFmat,qn->m,&qn->dFtdX,qn->m,&qn->YtdX);
515:   }
516:   SNESSetWorkVecs(snes,4);
517:   /* set method defaults */
518:   if (qn->scale_type == SNES_QN_SCALE_DEFAULT) {
519:     if (qn->type == SNES_QN_BADBROYDEN) {
520:       qn->scale_type = SNES_QN_SCALE_NONE;
521:     } else {
522:       qn->scale_type = SNES_QN_SCALE_SHANNO;
523:     }
524:   }
525:   if (qn->restart_type == SNES_QN_RESTART_DEFAULT) {
526:     if (qn->type == SNES_QN_LBFGS) {
527:       qn->restart_type = SNES_QN_RESTART_POWELL;
528:     } else {
529:       qn->restart_type = SNES_QN_RESTART_PERIODIC;
530:     }
531:   }

533:   if (qn->scale_type == SNES_QN_SCALE_JACOBIAN) {
534:     SNESSetUpMatrices(snes);
535:   }
536:   if (snes->pcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) {snes->functype = SNES_FUNCTION_UNPRECONDITIONED;}
537:   return(0);
538: }

542: static PetscErrorCode SNESReset_QN(SNES snes)
543: {
545:   SNES_QN        *qn;

548:   if (snes->data) {
549:     qn = (SNES_QN*)snes->data;
550:     if (qn->U) {
551:       VecDestroyVecs(qn->m, &qn->U);
552:     }
553:     if (qn->V) {
554:       VecDestroyVecs(qn->m, &qn->V);
555:     }
556:     if (qn->singlereduction) {
557:       PetscFree3(qn->dXdFmat, qn->dFtdX, qn->YtdX);
558:     }
559:     PetscFree4(qn->alpha,qn->beta,qn->dXtdF,qn->lambda);
560:   }
561:   return(0);
562: }

566: static PetscErrorCode SNESDestroy_QN(SNES snes)
567: {

571:   SNESReset_QN(snes);
572:   PetscFree(snes->data);
573:   PetscObjectComposeFunction((PetscObject)snes,"",NULL);
574:   return(0);
575: }

579: static PetscErrorCode SNESSetFromOptions_QN(PetscOptions *PetscOptionsObject,SNES snes)
580: {

582:   PetscErrorCode    ierr;
583:   SNES_QN           *qn    = (SNES_QN*)snes->data;
584:   PetscBool         monflg = PETSC_FALSE,flg;
585:   SNESLineSearch    linesearch;
586:   SNESQNRestartType rtype = qn->restart_type;
587:   SNESQNScaleType   stype = qn->scale_type;
588:   SNESQNType        qtype = qn->type;

591:   PetscOptionsHead(PetscOptionsObject,"SNES QN options");
592:   PetscOptionsInt("-snes_qn_m","Number of past states saved for L-BFGS methods","SNESQN",qn->m,&qn->m,NULL);
593:   PetscOptionsReal("-snes_qn_powell_gamma","Powell angle tolerance",          "SNESQN", qn->powell_gamma, &qn->powell_gamma, NULL);
594:   PetscOptionsReal("-snes_qn_powell_downhill","Powell descent tolerance",        "SNESQN", qn->powell_downhill, &qn->powell_downhill, NULL);
595:   PetscOptionsBool("-snes_qn_monitor",         "Monitor for the QN methods",      "SNESQN", monflg, &monflg, NULL);
596:   PetscOptionsBool("-snes_qn_single_reduction", "Aggregate reductions",           "SNESQN", qn->singlereduction, &qn->singlereduction, NULL);
597:   PetscOptionsEnum("-snes_qn_scale_type","Scaling type","SNESQNSetScaleType",SNESQNScaleTypes,(PetscEnum)stype,(PetscEnum*)&stype,&flg);
598:   if (flg) SNESQNSetScaleType(snes,stype);

600:   PetscOptionsEnum("-snes_qn_restart_type","Restart type","SNESQNSetRestartType",SNESQNRestartTypes,(PetscEnum)rtype,(PetscEnum*)&rtype,&flg);
601:   if (flg) SNESQNSetRestartType(snes,rtype);

603:   PetscOptionsEnum("-snes_qn_type","Quasi-Newton update type","",SNESQNTypes,(PetscEnum)qtype,(PetscEnum*)&qtype,&flg);
604:   if (flg) {SNESQNSetType(snes,qtype);}
605:   PetscOptionsTail();
606:   if (!snes->linesearch) {
607:     SNESGetLineSearch(snes, &linesearch);
608:     if (qn->type == SNES_QN_LBFGS) {
609:       SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
610:     } else if (qn->type == SNES_QN_BROYDEN) {
611:       SNESLineSearchSetType(linesearch, SNESLINESEARCHBASIC);
612:     } else {
613:       SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
614:     }
615:   }
616:   if (monflg) {
617:     qn->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
618:   }
619:   return(0);
620: }

624: static PetscErrorCode SNESView_QN(SNES snes, PetscViewer viewer)
625: {
626:   SNES_QN        *qn    = (SNES_QN*)snes->data;
627:   PetscBool      iascii;

631:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
632:   if (iascii) {
633:     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]);
634:     PetscViewerASCIIPrintf(viewer,"  Stored subspace size: %d\n", qn->m);
635:     if (qn->singlereduction) {
636:       PetscViewerASCIIPrintf(viewer,"  Using the single reduction variant.\n");
637:     }
638:   }
639:   return(0);
640: }

644: /*@
645:     SNESQNSetRestartType - Sets the restart type for SNESQN.

647:     Logically Collective on SNES

649:     Input Parameters:
650: +   snes - the iterative context
651: -   rtype - restart type

653:     Options Database:
654: +   -snes_qn_restart_type <powell,periodic,none> - set the restart type
655: -   -snes_qn_m <m> - sets the number of stored updates and the restart period for periodic

657:     Level: intermediate

659:     SNESQNRestartTypes:
660: +   SNES_QN_RESTART_NONE - never restart
661: .   SNES_QN_RESTART_POWELL - restart based upon descent criteria
662: -   SNES_QN_RESTART_PERIODIC - restart after a fixed number of iterations

664: .keywords: SNES, SNESQN, restart, type, set SNESLineSearch, SNESQNRestartType
665: @*/
666: PetscErrorCode SNESQNSetRestartType(SNES snes, SNESQNRestartType rtype)
667: {

672:   PetscTryMethod(snes,"SNESQNSetRestartType_C",(SNES,SNESQNRestartType),(snes,rtype));
673:   return(0);
674: }

678: /*@
679:     SNESQNSetScaleType - Sets the scaling type for the inner inverse Jacobian in SNESQN.

681:     Logically Collective on SNES

683:     Input Parameters:
684: +   snes - the iterative context
685: -   stype - scale type

687:     Options Database:
688: .   -snes_qn_scale_type <shanno,none,linesearch,jacobian>

690:     Level: intermediate

692:     SNESQNScaleTypes:
693: +   SNES_QN_SCALE_NONE - don't scale the problem
694: .   SNES_QN_SCALE_SHANNO - use shanno scaling
695: .   SNES_QN_SCALE_LINESEARCH - scale based upon line search lambda
696: -   SNES_QN_SCALE_JACOBIAN - scale by inverting a previously computed Jacobian.

698: .keywords: SNES, SNESQN, scaling, type, set SNESLineSearch, SNESQNScaleType
699: @*/

701: PetscErrorCode SNESQNSetScaleType(SNES snes, SNESQNScaleType stype)
702: {

707:   PetscTryMethod(snes,"SNESQNSetScaleType_C",(SNES,SNESQNScaleType),(snes,stype));
708:   return(0);
709: }

713: PetscErrorCode SNESQNSetScaleType_QN(SNES snes, SNESQNScaleType stype)
714: {
715:   SNES_QN *qn = (SNES_QN*)snes->data;

718:   qn->scale_type = stype;
719:   return(0);
720: }

724: PetscErrorCode SNESQNSetRestartType_QN(SNES snes, SNESQNRestartType rtype)
725: {
726:   SNES_QN *qn = (SNES_QN*)snes->data;

729:   qn->restart_type = rtype;
730:   return(0);
731: }

735: /*@
736:     SNESQNSetType - Sets the quasi-Newton variant to be used in SNESQN.

738:     Logically Collective on SNES

740:     Input Parameters:
741: +   snes - the iterative context
742: -   qtype - variant type

744:     Options Database:
745: .   -snes_qn_type <lbfgs,broyden,badbroyden>

747:     Level: beginner

749:     SNESQNTypes:
750: +   SNES_QN_LBFGS - LBFGS variant
751: .   SNES_QN_BROYDEN - Broyden variant
752: -   SNES_QN_BADBROYDEN - Bad Broyden variant

754: .keywords: SNES, SNESQN, type, set, SNESQNType
755: @*/

757: PetscErrorCode SNESQNSetType(SNES snes, SNESQNType qtype)
758: {

763:   PetscTryMethod(snes,"SNESQNSetType_C",(SNES,SNESQNType),(snes,qtype));
764:   return(0);
765: }

769: PetscErrorCode SNESQNSetType_QN(SNES snes, SNESQNType qtype)
770: {
771:   SNES_QN *qn = (SNES_QN*)snes->data;

774:   qn->type = qtype;
775:   return(0);
776: }

778: /* -------------------------------------------------------------------------- */
779: /*MC
780:       SNESQN - Limited-Memory Quasi-Newton methods for the solution of nonlinear systems.

782:       Options Database:

784: +     -snes_qn_m <m> - Number of past states saved for the L-Broyden methods.
785: +     -snes_qn_restart_type <powell,periodic,none> - set the restart type
786: .     -snes_qn_powell_angle - Angle condition for restart.
787: .     -snes_qn_powell_descent - Descent condition for restart.
788: .     -snes_qn_type <lbfgs,broyden,badbroyden> - QN type
789: .     -snes_qn_scale_type <shanno,none,linesearch,jacobian> - scaling performed on inner Jacobian
790: .     -snes_linesearch_type <cp, l2, basic> - Type of line search.
791: -     -snes_qn_monitor - Monitors the quasi-newton Jacobian.

793:       Notes: This implements the L-BFGS, Broyden, and "Bad" Broyden algorithms for the solution of F(x) = b using
794:       previous change in F(x) and x to form the approximate inverse Jacobian using a series of multiplicative rank-one
795:       updates.

797:       When using a nonlinear preconditioner, one has two options as to how the preconditioner is applied.  The first of
798:       these options, sequential, uses the preconditioner to generate a new solution and function and uses those at this
799:       iteration as the current iteration's values when constructing the approximate Jacobian.  The second, composed,
800:       perturbs the problem the Jacobian represents to be P(x, b) - x = 0, where P(x, b) is the preconditioner.

802:       References:

804:       Kelley, C.T., Iterative Methods for Linear and Nonlinear Equations, Chapter 8, SIAM, 1995.

806:       R. Byrd, J. Nocedal, R. Schnabel, Representations of Quasi-Newton Matrices and their use in Limited Memory Methods,
807:       Technical Report, Northwestern University, June 1992.

809:       Peter N. Brown, Alan C. Hindmarsh, Homer F. Walker, Experiments with Quasi-Newton Methods in Solving Stiff ODE
810:       Systems, SIAM J. Sci. Stat. Comput. Vol 6(2), April 1985.


813:       Level: beginner

815: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR

817: M*/
820: PETSC_EXTERN PetscErrorCode SNESCreate_QN(SNES snes)
821: {
823:   SNES_QN        *qn;

826:   snes->ops->setup          = SNESSetUp_QN;
827:   snes->ops->solve          = SNESSolve_QN;
828:   snes->ops->destroy        = SNESDestroy_QN;
829:   snes->ops->setfromoptions = SNESSetFromOptions_QN;
830:   snes->ops->view           = SNESView_QN;
831:   snes->ops->reset          = SNESReset_QN;

833:   snes->pcside = PC_LEFT;

835:   snes->usespc  = PETSC_TRUE;
836:   snes->usesksp = PETSC_FALSE;

838:   if (!snes->tolerancesset) {
839:     snes->max_funcs = 30000;
840:     snes->max_its   = 10000;
841:   }

843:   PetscNewLog(snes,&qn);
844:   snes->data          = (void*) qn;
845:   qn->m               = 10;
846:   qn->scaling         = 1.0;
847:   qn->U               = NULL;
848:   qn->V               = NULL;
849:   qn->dXtdF           = NULL;
850:   qn->dFtdX           = NULL;
851:   qn->dXdFmat         = NULL;
852:   qn->monitor         = NULL;
853:   qn->singlereduction = PETSC_TRUE;
854:   qn->powell_gamma    = 0.9999;
855:   qn->scale_type      = SNES_QN_SCALE_DEFAULT;
856:   qn->restart_type    = SNES_QN_RESTART_DEFAULT;
857:   qn->type            = SNES_QN_LBFGS;

859:   PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetScaleType_C",SNESQNSetScaleType_QN);
860:   PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetRestartType_C",SNESQNSetRestartType_QN);
861:   PetscObjectComposeFunction((PetscObject)snes,"SNESQNSetType_C",SNESQNSetType_QN);
862:   return(0);
863: }