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 SNESDestroy_NCG(SNES snes)
  5: {
  6:   PetscFunctionBegin;
  7:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", NULL));
  8:   PetscCall(PetscFree(snes->data));
  9:   PetscFunctionReturn(PETSC_SUCCESS);
 10: }

 12: static PetscErrorCode SNESSetUp_NCG(SNES snes)
 13: {
 14:   PetscFunctionBegin;
 15:   PetscCall(SNESSetWorkVecs(snes, 2));
 16:   PetscCheck(snes->npcside != PC_RIGHT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
 17:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 18:   PetscFunctionReturn(PETSC_SUCCESS);
 19: }

 21: static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
 22: {
 23:   PetscScalar alpha, ptAp;
 24:   Vec         X, Y, F, W;
 25:   SNES        snes;
 26:   PetscReal  *fnorm, *xnorm, *ynorm;

 28:   PetscFunctionBegin;
 29:   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
 30:   X     = linesearch->vec_sol;
 31:   W     = linesearch->vec_sol_new;
 32:   F     = linesearch->vec_func;
 33:   Y     = linesearch->vec_update;
 34:   fnorm = &linesearch->fnorm;
 35:   xnorm = &linesearch->xnorm;
 36:   ynorm = &linesearch->ynorm;

 38:   if (!snes->jacobian) PetscCall(SNESSetUpMatrices(snes));

 40:   /*
 41:    The exact step size for unpreconditioned linear CG is just:
 42:    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
 43:    */
 44:   PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
 45:   PetscCall(VecDot(F, F, &alpha));
 46:   PetscCall(MatMult(snes->jacobian, Y, W));
 47:   PetscCall(VecDot(Y, W, &ptAp));
 48:   alpha = alpha / ptAp;
 49:   PetscCall(VecAXPY(X, -alpha, Y));
 50:   PetscCall(SNESComputeFunction(snes, X, F));

 52:   PetscCall(VecNorm(F, NORM_2, fnorm));
 53:   PetscCall(VecNorm(X, NORM_2, xnorm));
 54:   PetscCall(VecNorm(Y, NORM_2, ynorm));
 55:   PetscFunctionReturn(PETSC_SUCCESS);
 56: }

 58: /*MC
 59:    SNESLINESEARCHNCGLINEAR - Special line search only for the nonlinear CG solver `SNESNCG`

 61:    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.
 62:    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.

 64:    Level: advanced

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

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

 71: .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
 72: M*/

 74: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
 75: {
 76:   PetscFunctionBegin;
 77:   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
 78:   linesearch->ops->destroy        = NULL;
 79:   linesearch->ops->setfromoptions = NULL;
 80:   linesearch->ops->reset          = NULL;
 81:   linesearch->ops->view           = NULL;
 82:   linesearch->ops->setup          = NULL;
 83:   PetscFunctionReturn(PETSC_SUCCESS);
 84: }

 86: static PetscErrorCode SNESSetFromOptions_NCG(SNES snes, PetscOptionItems PetscOptionsObject)
 87: {
 88:   SNES_NCG      *ncg     = (SNES_NCG *)snes->data;
 89:   PetscBool      debug   = PETSC_FALSE;
 90:   SNESNCGType    ncgtype = ncg->type;
 91:   SNESLineSearch linesearch;

 93:   PetscFunctionBegin;
 94:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES NCG options");
 95:   PetscCall(PetscOptionsBool("-snes_ncg_monitor", "Monitor the beta values used in the NCG iterations", "SNES", ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL));
 96:   if (debug) ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 97:   PetscCall(PetscOptionsEnum("-snes_ncg_type", "NCG Beta type used", "SNESNCGSetType", SNESNCGTypes, (PetscEnum)ncg->type, (PetscEnum *)&ncgtype, NULL));
 98:   PetscCall(SNESNCGSetType(snes, ncgtype));
 99:   PetscOptionsHeadEnd();
100:   if (!snes->linesearch) {
101:     PetscCall(SNESGetLineSearch(snes, &linesearch));
102:     if (!((PetscObject)linesearch)->type_name) {
103:       if (!snes->npc) {
104:         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHCP));
105:       } else {
106:         PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHSECANT));
107:       }
108:     }
109:   }
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
114: {
115:   SNES_NCG *ncg = (SNES_NCG *)snes->data;
116:   PetscBool isascii;

118:   PetscFunctionBegin;
119:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
120:   if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]));
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }

124: /*@
125:   SNESNCGSetType - Sets the conjugate update type for nonlinear CG `SNESNCG`.

127:   Logically Collective

129:   Input Parameters:
130: + snes  - the iterative context
131: - btype - update type, see `SNESNCGType`

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

136:   Level: intermediate

138:   Notes:
139:   `SNES_NCG_PRP` is the default, and the only one that tolerates generalized search directions.

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

144: .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESNCGType`, `SNES_NCG_FR`, `SNES_NCG_PRP`, `SNES_NCG_HS`, `SNES_NCG_DY`, `SNES_NCG_CD`
145: @*/
146: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
147: {
148:   PetscFunctionBegin;
150:   PetscTryMethod(snes, "SNESNCGSetType_C", (SNES, SNESNCGType), (snes, btype));
151:   PetscFunctionReturn(PETSC_SUCCESS);
152: }

154: static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
155: {
156:   SNES_NCG *ncg = (SNES_NCG *)snes->data;

158:   PetscFunctionBegin;
159:   ncg->type = btype;
160:   PetscFunctionReturn(PETSC_SUCCESS);
161: }

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

166:   Input Parameter:
167: . snes - the `SNES` context

169:   Output Parameter:
170: . outits - number of iterations until termination

172:   Application Interface Routine: SNESSolve()
173: */
174: static PetscErrorCode SNESSolve_NCG(SNES snes)
175: {
176:   SNES_NCG            *ncg = (SNES_NCG *)snes->data;
177:   Vec                  X, dX, lX, F, dXold;
178:   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
179:   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
180:   PetscInt             maxits, i;
181:   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
182:   SNESLineSearch       linesearch;
183:   SNESConvergedReason  reason;

185:   PetscFunctionBegin;
186:   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

188:   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
189:   snes->reason = SNES_CONVERGED_ITERATING;

191:   maxits = snes->max_its;        /* maximum number of iterations */
192:   X      = snes->vec_sol;        /* X^n */
193:   dXold  = snes->work[0];        /* The previous iterate of X */
194:   dX     = snes->work[1];        /* the preconditioned direction */
195:   lX     = snes->vec_sol_update; /* the conjugate direction */
196:   F      = snes->vec_func;       /* residual vector */

198:   PetscCall(SNESGetLineSearch(snes, &linesearch));

200:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
201:   snes->iter = 0;
202:   snes->norm = 0.;
203:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

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

207:   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
208:     PetscCall(SNESApplyNPC(snes, X, NULL, dX));
209:     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
210:     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
211:       snes->reason = SNES_DIVERGED_INNER;
212:       PetscFunctionReturn(PETSC_SUCCESS);
213:     }
214:     PetscCall(VecCopy(dX, F));
215:     PetscCall(VecNorm(F, NORM_2, &fnorm));
216:   } else {
217:     if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F));
218:     else snes->vec_func_init_set = PETSC_FALSE;

220:     /* convergence test */
221:     PetscCall(VecNorm(F, NORM_2, &fnorm));
222:     SNESCheckFunctionNorm(snes, fnorm);
223:     PetscCall(VecCopy(F, dX));
224:   }
225:   if (snes->npc) {
226:     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
227:       PetscCall(SNESApplyNPC(snes, X, F, dX));
228:       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
229:       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
230:         snes->reason = SNES_DIVERGED_INNER;
231:         PetscFunctionReturn(PETSC_SUCCESS);
232:       }
233:     }
234:   }
235:   PetscCall(VecCopy(dX, lX));
236:   PetscCall(VecDot(dX, dX, &dXdotdX));

238:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
239:   snes->norm = fnorm;
240:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
241:   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));

243:   /* test convergence */
244:   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
245:   PetscCall(SNESMonitor(snes, 0, fnorm));
246:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);

248:   /* Call general purpose update function */
249:   PetscTryTypeMethod(snes, update, snes->iter);

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

253:   for (i = 1; i < maxits + 1; i++) {
254:     /* some update types require the old update direction or conjugate direction */
255:     if (ncg->type != SNES_NCG_FR) PetscCall(VecCopy(dX, dXold));
256:     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, lX));
257:     PetscCall(SNESLineSearchGetReason(linesearch, &lsresult));
258:     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
259:     if (lsresult) {
260:       if (++snes->numFailures >= snes->maxFailures) {
261:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
262:         PetscFunctionReturn(PETSC_SUCCESS);
263:       }
264:     }
265:     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
266:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
267:       PetscFunctionReturn(PETSC_SUCCESS);
268:     }
269:     /* Monitor convergence */
270:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
271:     snes->iter  = i;
272:     snes->norm  = fnorm;
273:     snes->xnorm = xnorm;
274:     snes->ynorm = ynorm;
275:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
276:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));

278:     /* Test for convergence */
279:     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
280:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
281:     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);

283:     /* Call general purpose update function */
284:     PetscTryTypeMethod(snes, update, snes->iter);
285:     if (snes->npc) {
286:       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
287:         PetscCall(SNESApplyNPC(snes, X, NULL, dX));
288:         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
289:         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
290:           snes->reason = SNES_DIVERGED_INNER;
291:           PetscFunctionReturn(PETSC_SUCCESS);
292:         }
293:         PetscCall(VecCopy(dX, F));
294:       } else {
295:         PetscCall(SNESApplyNPC(snes, X, F, dX));
296:         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
297:         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
298:           snes->reason = SNES_DIVERGED_INNER;
299:           PetscFunctionReturn(PETSC_SUCCESS);
300:         }
301:       }
302:     } else {
303:       PetscCall(VecCopy(F, dX));
304:     }

306:     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
307:     switch (ncg->type) {
308:     case SNES_NCG_FR: /* Fletcher-Reeves */
309:       dXolddotdXold = dXdotdX;
310:       PetscCall(VecDot(dX, dX, &dXdotdX));
311:       beta = PetscRealPart(dXdotdX / dXolddotdXold);
312:       break;
313:     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
314:       dXolddotdXold = dXdotdX;
315:       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
316:       PetscCall(VecDotBegin(dXold, dX, &dXdotdXold));
317:       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
318:       PetscCall(VecDotEnd(dXold, dX, &dXdotdXold));
319:       beta = PetscRealPart((dXdotdX - dXdotdXold) / dXolddotdXold);
320:       if (beta < 0.0) beta = 0.0; /* restart */
321:       break;
322:     case SNES_NCG_HS: /* Hestenes-Stiefel */
323:       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
324:       PetscCall(VecDotBegin(dX, dXold, &dXdotdXold));
325:       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
326:       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
327:       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
328:       PetscCall(VecDotEnd(dX, dXold, &dXdotdXold));
329:       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
330:       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
331:       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
332:       break;
333:     case SNES_NCG_DY: /* Dai-Yuan */
334:       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
335:       PetscCall(VecDotBegin(lX, dX, &lXdotdX));
336:       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
337:       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
338:       PetscCall(VecDotEnd(lX, dX, &lXdotdX));
339:       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
340:       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
341:       break;
342:     case SNES_NCG_CD: /* Conjugate Descent */
343:       PetscCall(VecDotBegin(dX, dX, &dXdotdX));
344:       PetscCall(VecDotBegin(lX, dXold, &lXdotdXold));
345:       PetscCall(VecDotEnd(dX, dX, &dXdotdX));
346:       PetscCall(VecDotEnd(lX, dXold, &lXdotdXold));
347:       beta = PetscRealPart(dXdotdX / lXdotdXold);
348:       break;
349:     }
350:     if (ncg->monitor) PetscCall(PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta));
351:     PetscCall(VecAYPX(lX, beta, dX));
352:   }
353:   PetscFunctionReturn(PETSC_SUCCESS);
354: }

356: /*MC
357:   SNESNCG - Nonlinear Conjugate-Gradient method for the solution of nonlinear systems {cite}`bruneknepleysmithtu15`.

359:   Level: beginner

361:   Options Database Keys:
362: +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is `prp`.
363: .   -snes_linesearch_type <cp,l2,basic>  - Line search type.
364: -   -snes_ncg_monitor                    - Print the beta values nonlinear Conjugate-Gradient used in the  iteration, .

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

371:    Only supports left non-linear preconditioning.

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

376: .seealso: [](ch_snes), `SNES`, `SNESNCG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESLINESEARCHNCGLINEAR`, `SNESNCGSetType()`, `SNESLineSearchSetType()`
377: M*/
378: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
379: {
380:   SNES_NCG *neP;

382:   PetscFunctionBegin;
383:   snes->ops->destroy        = SNESDestroy_NCG;
384:   snes->ops->setup          = SNESSetUp_NCG;
385:   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
386:   snes->ops->view           = SNESView_NCG;
387:   snes->ops->solve          = SNESSolve_NCG;

389:   snes->usesksp = PETSC_FALSE;
390:   snes->usesnpc = PETSC_TRUE;
391:   snes->npcside = PC_LEFT;

393:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

395:   PetscCall(SNESParametersInitialize(snes));
396:   PetscObjectParameterSetDefault(snes, max_funcs, 30000);
397:   PetscObjectParameterSetDefault(snes, max_its, 10000);
398:   PetscObjectParameterSetDefault(snes, stol, 1e-20);

400:   PetscCall(PetscNew(&neP));
401:   snes->data   = (void *)neP;
402:   neP->monitor = NULL;
403:   neP->type    = SNES_NCG_PRP;
404:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNCGSetType_C", SNESNCGSetType_NCG));
405:   PetscFunctionReturn(PETSC_SUCCESS);
406: }