Actual source code: snesncg.c

  1: #include <../src/snes/impls/ncg/snesncgimpl.h>
  2: const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",NULL};

  4: static PetscErrorCode SNESReset_NCG(SNES snes)
  5: {
  6:   return 0;
  7: }

  9: /*
 10:   SNESDestroy_NCG - Destroys the private SNES_NCG context that was created with SNESCreate_NCG().

 12:   Input Parameter:
 13: . snes - the SNES context

 15:   Application Interface Routine: SNESDestroy()
 16: */
 17: static PetscErrorCode SNESDestroy_NCG(SNES snes)
 18: {
 19:   PetscFree(snes->data);
 20:   return 0;
 21: }

 23: /*
 24:    SNESSetUp_NCG - Sets up the internal data structures for the later use
 25:    of the SNESNCG nonlinear solver.

 27:    Input Parameters:
 28: +  snes - the SNES context
 29: -  x - the solution vector

 31:    Application Interface Routine: SNESSetUp()
 32:  */

 34: static PetscErrorCode SNESSetUp_NCG(SNES snes)
 35: {
 36:   SNESSetWorkVecs(snes,2);
 38:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 39:   return 0;
 40: }

 42: static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
 43: {
 44:   PetscScalar    alpha, ptAp;
 45:   Vec            X, Y, F, W;
 46:   SNES           snes;
 47:   PetscReal      *fnorm, *xnorm, *ynorm;

 49:   SNESLineSearchGetSNES(linesearch, &snes);
 50:   X     = linesearch->vec_sol;
 51:   W     = linesearch->vec_sol_new;
 52:   F     = linesearch->vec_func;
 53:   Y     = linesearch->vec_update;
 54:   fnorm = &linesearch->fnorm;
 55:   xnorm = &linesearch->xnorm;
 56:   ynorm = &linesearch->ynorm;

 58:   if (!snes->jacobian) {
 59:     SNESSetUpMatrices(snes);
 60:   }

 62:   /*

 64:    The exact step size for unpreconditioned linear CG is just:
 65:    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
 66:    */
 67:   SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);
 68:   VecDot(F, F, &alpha);
 69:   MatMult(snes->jacobian, Y, W);
 70:   VecDot(Y, W, &ptAp);
 71:   alpha = alpha / ptAp;
 72:   VecAXPY(X,-alpha,Y);
 73:   SNESComputeFunction(snes,X,F);

 75:   VecNorm(F,NORM_2,fnorm);
 76:   VecNorm(X,NORM_2,xnorm);
 77:   VecNorm(Y,NORM_2,ynorm);
 78:   return 0;
 79: }

 81: /*MC
 82:    SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG

 84:    This line search uses the length "as if" the problem is linear (that is what is computed by the linear CG method) using the Jacobian of the function.
 85:    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy) where r (f) is the current residual (function value), p (y) is the current search direction.

 87:    Notes: This requires a Jacobian-vector product but does not require the solution of a linear system with the Jacobian

 89:    This is a "odd-ball" line search, we don't know if it is in the literature or used in practice by anyone.

 91:    Level: advanced

 93: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
 94: M*/

 96: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
 97: {
 98:   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
 99:   linesearch->ops->destroy        = NULL;
100:   linesearch->ops->setfromoptions = NULL;
101:   linesearch->ops->reset          = NULL;
102:   linesearch->ops->view           = NULL;
103:   linesearch->ops->setup          = NULL;
104:   return 0;
105: }

107: /*
108:   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.

110:   Input Parameter:
111: . snes - the SNES context

113:   Application Interface Routine: SNESSetFromOptions()
114: */
115: static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
116: {
117:   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
118:   PetscBool      debug = PETSC_FALSE;
119:   SNESNCGType    ncgtype=ncg->type;
120:   SNESLineSearch linesearch;

122:   PetscOptionsHead(PetscOptionsObject,"SNES NCG options");
123:   PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
124:   if (debug) {
125:     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
126:   }
127:   PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
128:   SNESNCGSetType(snes, ncgtype);
129:   PetscOptionsTail();
130:   if (!snes->linesearch) {
131:     SNESGetLineSearch(snes, &linesearch);
132:     if (!((PetscObject)linesearch)->type_name) {
133:       if (!snes->npc) {
134:         SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
135:       } else {
136:         SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
137:       }
138:     }
139:   }
140:   return 0;
141: }

143: /*
144:   SNESView_NCG - Prints info from the SNESNCG data structure.

146:   Input Parameters:
147: + SNES - the SNES context
148: - viewer - visualization context

150:   Application Interface Routine: SNESView()
151: */
152: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
153: {
154:   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
155:   PetscBool      iascii;

157:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
158:   if (iascii) {
159:     PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]);
160:   }
161:   return 0;
162: }

164: /*

166:  Assuming F = SNESComputeFunction(X) compute Y^tJ^tF using a simple secant approximation of the jacobian.

168:  This routine is not currently used. I don't know what its intended purpose is.

170:  Note it has a hardwired differencing parameter of 1e-5

172:  */
173: PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
174: {
175:   PetscScalar    ftf, ftg, fty, h;

177:   VecDot(F, F, &ftf);
178:   VecDot(F, Y, &fty);
179:   h      = 1e-5*fty / fty;
180:   VecCopy(X, W);
181:   VecAXPY(W, -h, Y);          /* this is arbitrary */
182:   SNESComputeFunction(snes, W, G);
183:   VecDot(G, F, &ftg);
184:   *ytJtf = (ftg - ftf) / h;
185:   return 0;
186: }

188: /*@
189:     SNESNCGSetType - Sets the conjugate update type for SNESNCG.

191:     Logically Collective on SNES

193:     Input Parameters:
194: +   snes - the iterative context
195: -   btype - update type

197:     Options Database:
198: .   -snes_ncg_type <prp,fr,hs,dy,cd> -strategy for selecting algorithm for computing beta

200:     Level: intermediate

202:     SNESNCGSelectTypes:
203: +   SNES_NCG_FR - Fletcher-Reeves update
204: .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
205: .   SNES_NCG_HS - Hestenes-Steifel update
206: .   SNES_NCG_DY - Dai-Yuan update
207: -   SNES_NCG_CD - Conjugate Descent update

209:    Notes:
210:    SNES_NCG_PRP is the default, and the only one that tolerates generalized search directions.

212:    It is not clear what "generalized search directions" means, does it mean use with a nonlinear preconditioner,
213:    that is using -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC()?

215: @*/
216: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
217: {
219:   PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
220:   return 0;
221: }

223: static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
224: {
225:   SNES_NCG *ncg = (SNES_NCG*)snes->data;

227:   ncg->type = btype;
228:   return 0;
229: }

231: /*
232:   SNESSolve_NCG - Solves a nonlinear system with the Nonlinear Conjugate Gradient method.

234:   Input Parameters:
235: . snes - the SNES context

237:   Output Parameter:
238: . outits - number of iterations until termination

240:   Application Interface Routine: SNESSolve()
241: */
242: static PetscErrorCode SNESSolve_NCG(SNES snes)
243: {
244:   SNES_NCG             *ncg = (SNES_NCG*)snes->data;
245:   Vec                  X,dX,lX,F,dXold;
246:   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
247:   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
248:   PetscInt             maxits, i;
249:   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
250:   SNESLineSearch       linesearch;
251:   SNESConvergedReason  reason;


255:   PetscCitationsRegister(SNESCitation,&SNEScite);
256:   snes->reason = SNES_CONVERGED_ITERATING;

258:   maxits = snes->max_its;            /* maximum number of iterations */
259:   X      = snes->vec_sol;            /* X^n */
260:   dXold  = snes->work[0];            /* The previous iterate of X */
261:   dX     = snes->work[1];            /* the preconditioned direction */
262:   lX     = snes->vec_sol_update;     /* the conjugate direction */
263:   F      = snes->vec_func;           /* residual vector */

265:   SNESGetLineSearch(snes, &linesearch);

267:   PetscObjectSAWsTakeAccess((PetscObject)snes);
268:   snes->iter = 0;
269:   snes->norm = 0.;
270:   PetscObjectSAWsGrantAccess((PetscObject)snes);

272:   /* compute the initial function and preconditioned update dX */

274:   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
275:     SNESApplyNPC(snes,X,NULL,dX);
276:     SNESGetConvergedReason(snes->npc,&reason);
277:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
278:       snes->reason = SNES_DIVERGED_INNER;
279:       return 0;
280:     }
281:     VecCopy(dX,F);
282:     VecNorm(F,NORM_2,&fnorm);
283:   } else {
284:     if (!snes->vec_func_init_set) {
285:       SNESComputeFunction(snes,X,F);
286:     } else snes->vec_func_init_set = PETSC_FALSE;

288:     /* convergence test */
289:     VecNorm(F,NORM_2,&fnorm);
290:     SNESCheckFunctionNorm(snes,fnorm);
291:     VecCopy(F,dX);
292:   }
293:   if (snes->npc) {
294:     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
295:       SNESApplyNPC(snes,X,F,dX);
296:       SNESGetConvergedReason(snes->npc,&reason);
297:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
298:         snes->reason = SNES_DIVERGED_INNER;
299:         return 0;
300:       }
301:     }
302:   }
303:   VecCopy(dX,lX);
304:   VecDot(dX, dX, &dXdotdX);

306:   PetscObjectSAWsTakeAccess((PetscObject)snes);
307:   snes->norm = fnorm;
308:   PetscObjectSAWsGrantAccess((PetscObject)snes);
309:   SNESLogConvergenceHistory(snes,fnorm,0);
310:   SNESMonitor(snes,0,fnorm);

312:   /* test convergence */
313:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
314:   if (snes->reason) return 0;

316:   /* Call general purpose update function */
317:   if (snes->ops->update) {
318:     (*snes->ops->update)(snes, snes->iter);
319:   }

321:   /* first update -- just use the (preconditioned) residual direction for the initial conjugate direction */

323:   for (i = 1; i < maxits + 1; i++) {
324:     /* some update types require the old update direction or conjugate direction */
325:     if (ncg->type != SNES_NCG_FR) {
326:       VecCopy(dX, dXold);
327:     }
328:     SNESLineSearchApply(linesearch,X,F,&fnorm,lX);
329:     SNESLineSearchGetReason(linesearch, &lsresult);
330:     SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
331:     if (lsresult) {
332:       if (++snes->numFailures >= snes->maxFailures) {
333:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
334:         return 0;
335:       }
336:     }
337:     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
338:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
339:       return 0;
340:     }
341:     /* Monitor convergence */
342:     PetscObjectSAWsTakeAccess((PetscObject)snes);
343:     snes->iter = i;
344:     snes->norm = fnorm;
345:     snes->xnorm = xnorm;
346:     snes->ynorm = ynorm;
347:     PetscObjectSAWsGrantAccess((PetscObject)snes);
348:     SNESLogConvergenceHistory(snes,snes->norm,0);
349:     SNESMonitor(snes,snes->iter,snes->norm);

351:     /* Test for convergence */
352:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
353:     if (snes->reason) return 0;

355:     /* Call general purpose update function */
356:     if (snes->ops->update) {
357:       (*snes->ops->update)(snes, snes->iter);
358:     }
359:     if (snes->npc) {
360:       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
361:         SNESApplyNPC(snes,X,NULL,dX);
362:         SNESGetConvergedReason(snes->npc,&reason);
363:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
364:           snes->reason = SNES_DIVERGED_INNER;
365:           return 0;
366:         }
367:         VecCopy(dX,F);
368:       } else {
369:         SNESApplyNPC(snes,X,F,dX);
370:         SNESGetConvergedReason(snes->npc,&reason);
371:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
372:           snes->reason = SNES_DIVERGED_INNER;
373:           return 0;
374:         }
375:       }
376:     } else {
377:       VecCopy(F,dX);
378:     }

380:     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
381:     switch (ncg->type) {
382:     case SNES_NCG_FR: /* Fletcher-Reeves */
383:       dXolddotdXold= dXdotdX;
384:       VecDot(dX, dX, &dXdotdX);
385:       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
386:       break;
387:     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
388:       dXolddotdXold= dXdotdX;
389:       VecDotBegin(dX, dX, &dXdotdXold);
390:       VecDotBegin(dXold, dX, &dXdotdXold);
391:       VecDotEnd(dX, dX, &dXdotdX);
392:       VecDotEnd(dXold, dX, &dXdotdXold);
393:       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
394:       if (beta < 0.0) beta = 0.0; /* restart */
395:       break;
396:     case SNES_NCG_HS: /* Hestenes-Stiefel */
397:       VecDotBegin(dX, dX, &dXdotdX);
398:       VecDotBegin(dX, dXold, &dXdotdXold);
399:       VecDotBegin(lX, dX, &lXdotdX);
400:       VecDotBegin(lX, dXold, &lXdotdXold);
401:       VecDotEnd(dX, dX, &dXdotdX);
402:       VecDotEnd(dX, dXold, &dXdotdXold);
403:       VecDotEnd(lX, dX, &lXdotdX);
404:       VecDotEnd(lX, dXold, &lXdotdXold);
405:       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
406:       break;
407:     case SNES_NCG_DY: /* Dai-Yuan */
408:       VecDotBegin(dX, dX, &dXdotdX);
409:       VecDotBegin(lX, dX, &lXdotdX);
410:       VecDotBegin(lX, dXold, &lXdotdXold);
411:       VecDotEnd(dX, dX, &dXdotdX);
412:       VecDotEnd(lX, dX, &lXdotdX);
413:       VecDotEnd(lX, dXold, &lXdotdXold);
414:       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
415:       break;
416:     case SNES_NCG_CD: /* Conjugate Descent */
417:       VecDotBegin(dX, dX, &dXdotdX);
418:       VecDotBegin(lX, dXold, &lXdotdXold);
419:       VecDotEnd(dX, dX, &dXdotdX);
420:       VecDotEnd(lX, dXold, &lXdotdXold);
421:       beta = PetscRealPart(dXdotdX / lXdotdXold);
422:       break;
423:     }
424:     if (ncg->monitor) {
425:       PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);
426:     }
427:     VecAYPX(lX, beta, dX);
428:   }
429:   PetscInfo(snes, "Maximum number of iterations has been reached: %D\n", maxits);
430:   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
431:   return 0;
432: }

434: /*MC
435:   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of general nonlinear systems.

437:   Level: beginner

439:   Options Database:
440: +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
441: .   -snes_linesearch_type <cp,l2,basic> - Line search type.
442: -   -snes_ncg_monitor - Print the beta values used in the ncg iteration, .

444:    Notes:
445:     This solves the nonlinear system of equations F(x) = 0 using the nonlinear generalization of the conjugate
446:           gradient method.  This may be used with a nonlinear preconditioner used to pick the new search directions, but otherwise
447:           chooses the initial search direction as F(x) for the initial guess x.

449:           Only supports left non-linear preconditioning.

451:     Default line search is SNESLINESEARCHCP, unless a nonlinear preconditioner is used with -npc_snes_type <type>, SNESSetNPC(), or SNESGetNPC() then
452:     SNESLINESEARCHL2 is used. Also supports the special purpose line search SNESLINESEARCHNCGLINEAR

454:    References:
455: .  * -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
456:    SIAM Review, 57(4), 2015

458: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType()
459: M*/
460: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
461: {
462:   SNES_NCG       * neP;

464:   snes->ops->destroy        = SNESDestroy_NCG;
465:   snes->ops->setup          = SNESSetUp_NCG;
466:   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
467:   snes->ops->view           = SNESView_NCG;
468:   snes->ops->solve          = SNESSolve_NCG;
469:   snes->ops->reset          = SNESReset_NCG;

471:   snes->usesksp = PETSC_FALSE;
472:   snes->usesnpc = PETSC_TRUE;
473:   snes->npcside = PC_LEFT;

475:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

477:   if (!snes->tolerancesset) {
478:     snes->max_funcs = 30000;
479:     snes->max_its   = 10000;
480:     snes->stol      = 1e-20;
481:   }

483:   PetscNewLog(snes,&neP);
484:   snes->data   = (void*) neP;
485:   neP->monitor = NULL;
486:   neP->type    = SNES_NCG_PRP;
487:   PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);
488:   return 0;
489: }