Actual source code: snesncg.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1:  #include <../src/snes/impls/ncg/snesncgimpl.h>
  2: const char *const SNESNCGTypes[] = {"FR","PRP","HS","DY","CD","SNESNCGType","SNES_NCG_",0};

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

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

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

 16:   Application Interface Routine: SNESDestroy()
 17: */
 18: static PetscErrorCode SNESDestroy_NCG(SNES snes)
 19: {

 23:   PetscFree(snes->data);
 24:   return(0);
 25: }

 27: /*
 28:    SNESSetUp_NCG - Sets up the internal data structures for the later use
 29:    of the SNESNCG nonlinear solver.

 31:    Input Parameters:
 32: +  snes - the SNES context
 33: -  x - the solution vector

 35:    Application Interface Routine: SNESSetUp()
 36:  */

 38: static PetscErrorCode SNESSetUp_NCG(SNES snes)
 39: {

 43:   SNESSetWorkVecs(snes,2);
 44:   if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
 45:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 46:   return(0);
 47: }

 49: static PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
 50: {
 51:   PetscScalar    alpha, ptAp;
 52:   Vec            X, Y, F, W;
 53:   SNES           snes;
 55:   PetscReal      *fnorm, *xnorm, *ynorm;

 58:   SNESLineSearchGetSNES(linesearch, &snes);
 59:   X     = linesearch->vec_sol;
 60:   W     = linesearch->vec_sol_new;
 61:   F     = linesearch->vec_func;
 62:   Y     = linesearch->vec_update;
 63:   fnorm = &linesearch->fnorm;
 64:   xnorm = &linesearch->xnorm;
 65:   ynorm = &linesearch->ynorm;

 67:   if (!snes->jacobian) {
 68:     SNESSetUpMatrices(snes);
 69:   }

 71:   /*

 73:    The exact step size for unpreconditioned linear CG is just:
 74:    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
 75:    */
 76:   SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);
 77:   VecDot(F, F, &alpha);
 78:   MatMult(snes->jacobian, Y, W);
 79:   VecDot(Y, W, &ptAp);
 80:   alpha = alpha / ptAp;
 81:   VecAXPY(X,-alpha,Y);
 82:   SNESComputeFunction(snes,X,F);

 84:   VecNorm(F,NORM_2,fnorm);
 85:   VecNorm(X,NORM_2,xnorm);
 86:   VecNorm(Y,NORM_2,ynorm);
 87:   return(0);
 88: }

 90: /*MC
 91:    SNESLINESEARCHNCGLINEAR - Special line search only for SNESNCG

 93:    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.
 94:    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.

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

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

100:    Level: advanced

102: .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
103: M*/

105: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
106: {
108:   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
109:   linesearch->ops->destroy        = NULL;
110:   linesearch->ops->setfromoptions = NULL;
111:   linesearch->ops->reset          = NULL;
112:   linesearch->ops->view           = NULL;
113:   linesearch->ops->setup          = NULL;
114:   return(0);
115: }

117: /*
118:   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.

120:   Input Parameter:
121: . snes - the SNES context

123:   Application Interface Routine: SNESSetFromOptions()
124: */
125: static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
126: {
127:   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
129:   PetscBool      debug = PETSC_FALSE;
130:   SNESNCGType    ncgtype=ncg->type;
131:   SNESLineSearch linesearch;

134:   PetscOptionsHead(PetscOptionsObject,"SNES NCG options");
135:   PetscOptionsBool("-snes_ncg_monitor","Monitor the beta values used in the NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
136:   if (debug) {
137:     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
138:   }
139:   PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
140:   SNESNCGSetType(snes, ncgtype);
141:   PetscOptionsTail();
142:   if (!snes->linesearch) {
143:     SNESGetLineSearch(snes, &linesearch);
144:     if (!snes->npc) {
145:       SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
146:     } else {
147:       SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
148:     }
149:   }
150:   return(0);
151: }

153: /*
154:   SNESView_NCG - Prints info from the SNESNCG data structure.

156:   Input Parameters:
157: + SNES - the SNES context
158: - viewer - visualization context

160:   Application Interface Routine: SNESView()
161: */
162: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
163: {
164:   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
165:   PetscBool      iascii;

169:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
170:   if (iascii) {
171:     PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]);
172:   }
173:   return(0);
174: }

176: /*

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

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

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

184:  */
185: PetscErrorCode SNESNCGComputeYtJtF_Private(SNES snes, Vec X, Vec F, Vec Y, Vec W, Vec G, PetscScalar * ytJtf)
186: {
188:   PetscScalar    ftf, ftg, fty, h;

191:   VecDot(F, F, &ftf);
192:   VecDot(F, Y, &fty);
193:   h      = 1e-5*fty / fty;
194:   VecCopy(X, W);
195:   VecAXPY(W, -h, Y);          /* this is arbitrary */
196:   SNESComputeFunction(snes, W, G);
197:   VecDot(G, F, &ftg);
198:   *ytJtf = (ftg - ftf) / h;
199:   return(0);
200: }

202: /*@
203:     SNESNCGSetType - Sets the conjugate update type for SNESNCG.

205:     Logically Collective on SNES

207:     Input Parameters:
208: +   snes - the iterative context
209: -   btype - update type

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

214:     Level: intermediate

216:     SNESNCGSelectTypes:
217: +   SNES_NCG_FR - Fletcher-Reeves update
218: .   SNES_NCG_PRP - Polak-Ribiere-Polyak update
219: .   SNES_NCG_HS - Hestenes-Steifel update
220: .   SNES_NCG_DY - Dai-Yuan update
221: -   SNES_NCG_CD - Conjugate Descent update

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

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

229: @*/
230: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
231: {

236:   PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
237:   return(0);
238: }

240: static PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
241: {
242:   SNES_NCG *ncg = (SNES_NCG*)snes->data;

245:   ncg->type = btype;
246:   return(0);
247: }

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

252:   Input Parameters:
253: . snes - the SNES context

255:   Output Parameter:
256: . outits - number of iterations until termination

258:   Application Interface Routine: SNESSolve()
259: */
260: static PetscErrorCode SNESSolve_NCG(SNES snes)
261: {
262:   SNES_NCG             *ncg = (SNES_NCG*)snes->data;
263:   Vec                  X,dX,lX,F,dXold;
264:   PetscReal            fnorm, ynorm, xnorm, beta = 0.0;
265:   PetscScalar          dXdotdX, dXolddotdXold, dXdotdXold, lXdotdX, lXdotdXold;
266:   PetscInt             maxits, i;
267:   PetscErrorCode       ierr;
268:   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
269:   SNESLineSearch       linesearch;
270:   SNESConvergedReason  reason;

273:   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);

275:   PetscCitationsRegister(SNESCitation,&SNEScite);
276:   snes->reason = SNES_CONVERGED_ITERATING;

278:   maxits = snes->max_its;            /* maximum number of iterations */
279:   X      = snes->vec_sol;            /* X^n */
280:   dXold  = snes->work[0];            /* The previous iterate of X */
281:   dX     = snes->work[1];            /* the preconditioned direction */
282:   lX     = snes->vec_sol_update;     /* the conjugate direction */
283:   F      = snes->vec_func;           /* residual vector */

285:   SNESGetLineSearch(snes, &linesearch);

287:   PetscObjectSAWsTakeAccess((PetscObject)snes);
288:   snes->iter = 0;
289:   snes->norm = 0.;
290:   PetscObjectSAWsGrantAccess((PetscObject)snes);

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

294:   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
295:     SNESApplyNPC(snes,X,NULL,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:     VecCopy(dX,F);
302:     VecNorm(F,NORM_2,&fnorm);
303:   } else {
304:     if (!snes->vec_func_init_set) {
305:       SNESComputeFunction(snes,X,F);
306:     } else snes->vec_func_init_set = PETSC_FALSE;

308:     /* convergence test */
309:     VecNorm(F,NORM_2,&fnorm);
310:     SNESCheckFunctionNorm(snes,fnorm);
311:     VecCopy(F,dX);
312:   }
313:   if (snes->npc) {
314:     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
315:       SNESApplyNPC(snes,X,F,dX);
316:       SNESGetConvergedReason(snes->npc,&reason);
317:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
318:         snes->reason = SNES_DIVERGED_INNER;
319:         return(0);
320:       }
321:     }
322:   }
323:   VecCopy(dX,lX);
324:   VecDot(dX, dX, &dXdotdX);


327:   PetscObjectSAWsTakeAccess((PetscObject)snes);
328:   snes->norm = fnorm;
329:   PetscObjectSAWsGrantAccess((PetscObject)snes);
330:   SNESLogConvergenceHistory(snes,fnorm,0);
331:   SNESMonitor(snes,0,fnorm);

333:   /* test convergence */
334:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
335:   if (snes->reason) return(0);

337:   /* Call general purpose update function */
338:   if (snes->ops->update) {
339:     (*snes->ops->update)(snes, snes->iter);
340:   }

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

344:   for (i = 1; i < maxits + 1; i++) {
345:     /* some update types require the old update direction or conjugate direction */
346:     if (ncg->type != SNES_NCG_FR) {
347:       VecCopy(dX, dXold);
348:     }
349:     SNESLineSearchApply(linesearch,X,F,&fnorm,lX);
350:     SNESLineSearchGetReason(linesearch, &lsresult);
351:     SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
352:     if (lsresult) {
353:       if (++snes->numFailures >= snes->maxFailures) {
354:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
355:         return(0);
356:       }
357:     }
358:     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
359:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
360:       return(0);
361:     }
362:     /* Monitor convergence */
363:     PetscObjectSAWsTakeAccess((PetscObject)snes);
364:     snes->iter = i;
365:     snes->norm = fnorm;
366:     snes->xnorm = xnorm;
367:     snes->ynorm = ynorm;
368:     PetscObjectSAWsGrantAccess((PetscObject)snes);
369:     SNESLogConvergenceHistory(snes,snes->norm,0);
370:     SNESMonitor(snes,snes->iter,snes->norm);

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

376:     /* Call general purpose update function */
377:     if (snes->ops->update) {
378:       (*snes->ops->update)(snes, snes->iter);
379:     }
380:     if (snes->npc) {
381:       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
382:         SNESApplyNPC(snes,X,NULL,dX);
383:         SNESGetConvergedReason(snes->npc,&reason);
384:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
385:           snes->reason = SNES_DIVERGED_INNER;
386:           return(0);
387:         }
388:         VecCopy(dX,F);
389:       } else {
390:         SNESApplyNPC(snes,X,F,dX);
391:         SNESGetConvergedReason(snes->npc,&reason);
392:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
393:           snes->reason = SNES_DIVERGED_INNER;
394:           return(0);
395:         }
396:       }
397:     } else {
398:       VecCopy(F,dX);
399:     }

401:     /* compute the conjugate direction lX = dX + beta*lX with beta = ((dX, dX) / (dX_old, dX_old) (Fletcher-Reeves update)*/
402:     switch (ncg->type) {
403:     case SNES_NCG_FR: /* Fletcher-Reeves */
404:       dXolddotdXold= dXdotdX;
405:       VecDot(dX, dX, &dXdotdX);
406:       beta         = PetscRealPart(dXdotdX / dXolddotdXold);
407:       break;
408:     case SNES_NCG_PRP: /* Polak-Ribiere-Poylak */
409:       dXolddotdXold= dXdotdX;
410:       VecDotBegin(dX, dX, &dXdotdXold);
411:       VecDotBegin(dXold, dX, &dXdotdXold);
412:       VecDotEnd(dX, dX, &dXdotdX);
413:       VecDotEnd(dXold, dX, &dXdotdXold);
414:       beta         = PetscRealPart(((dXdotdX - dXdotdXold) / dXolddotdXold));
415:       if (beta < 0.0) beta = 0.0; /* restart */
416:       break;
417:     case SNES_NCG_HS: /* Hestenes-Stiefel */
418:       VecDotBegin(dX, dX, &dXdotdX);
419:       VecDotBegin(dX, dXold, &dXdotdXold);
420:       VecDotBegin(lX, dX, &lXdotdX);
421:       VecDotBegin(lX, dXold, &lXdotdXold);
422:       VecDotEnd(dX, dX, &dXdotdX);
423:       VecDotEnd(dX, dXold, &dXdotdXold);
424:       VecDotEnd(lX, dX, &lXdotdX);
425:       VecDotEnd(lX, dXold, &lXdotdXold);
426:       beta = PetscRealPart((dXdotdX - dXdotdXold) / (lXdotdX - lXdotdXold));
427:       break;
428:     case SNES_NCG_DY: /* Dai-Yuan */
429:       VecDotBegin(dX, dX, &dXdotdX);
430:       VecDotBegin(lX, dX, &lXdotdX);
431:       VecDotBegin(lX, dXold, &lXdotdXold);
432:       VecDotEnd(dX, dX, &dXdotdX);
433:       VecDotEnd(lX, dX, &lXdotdX);
434:       VecDotEnd(lX, dXold, &lXdotdXold);
435:       beta = PetscRealPart(dXdotdX / (lXdotdXold - lXdotdX));
436:       break;
437:     case SNES_NCG_CD: /* Conjugate Descent */
438:       VecDotBegin(dX, dX, &dXdotdX);
439:       VecDotBegin(lX, dXold, &lXdotdXold);
440:       VecDotEnd(dX, dX, &dXdotdX);
441:       VecDotEnd(lX, dXold, &lXdotdXold);
442:       beta = PetscRealPart(dXdotdX / lXdotdXold);
443:       break;
444:     }
445:     if (ncg->monitor) {
446:       PetscViewerASCIIPrintf(ncg->monitor, "beta = %e\n", (double)beta);
447:     }
448:     VecAYPX(lX, beta, dX);
449:   }
450:   PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
451:   if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
452:   return(0);
453: }

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

458:   Level: beginner

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

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

470:           Only supports left non-linear preconditioning.

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

475:    References:
476: .  1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu,"Composing Scalable Nonlinear Algebraic Solvers",
477:    SIAM Review, 57(4), 2015


480: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESLINESEARCHNCGLINEAR, SNESNCGSetType(), SNESNCGGetType(), SNESLineSearchSetType()
481: M*/
482: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
483: {
485:   SNES_NCG       * neP;

488:   snes->ops->destroy        = SNESDestroy_NCG;
489:   snes->ops->setup          = SNESSetUp_NCG;
490:   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
491:   snes->ops->view           = SNESView_NCG;
492:   snes->ops->solve          = SNESSolve_NCG;
493:   snes->ops->reset          = SNESReset_NCG;

495:   snes->usesksp = PETSC_FALSE;
496:   snes->usesnpc = PETSC_TRUE;
497:   snes->npcside = PC_LEFT;

499:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

501:   if (!snes->tolerancesset) {
502:     snes->max_funcs = 30000;
503:     snes->max_its   = 10000;
504:     snes->stol      = 1e-20;
505:   }

507:   PetscNewLog(snes,&neP);
508:   snes->data   = (void*) neP;
509:   neP->monitor = NULL;
510:   neP->type    = SNES_NCG_PRP;
511:   PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);
512:   return(0);
513: }