Actual source code: snesncg.c

petsc-3.9.4 2018-09-11
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: PetscErrorCode SNESReset_NCG(SNES snes)
  5: {
  7:   return(0);
  8: }

 10: #define SNESLINESEARCHNCGLINEAR "ncglinear"

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

 15:   Input Parameter:
 16: . snes - the SNES context

 18:   Application Interface Routine: SNESDestroy()
 19: */
 20: PetscErrorCode SNESDestroy_NCG(SNES snes)
 21: {

 25:   PetscFree(snes->data);
 26:   return(0);
 27: }

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

 33:    Input Parameters:
 34: +  snes - the SNES context
 35: -  x - the solution vector

 37:    Application Interface Routine: SNESSetUp()
 38:  */

 40: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch);

 42: PetscErrorCode SNESSetUp_NCG(SNES snes)
 43: {

 47:   SNESSetWorkVecs(snes,2);
 48:   if (snes->npcside== PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
 49:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 50:   return(0);
 51: }
 52: /*
 53:   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.

 55:   Input Parameter:
 56: . snes - the SNES context

 58:   Application Interface Routine: SNESSetFromOptions()
 59: */
 60: static PetscErrorCode SNESSetFromOptions_NCG(PetscOptionItems *PetscOptionsObject,SNES snes)
 61: {
 62:   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
 64:   PetscBool      debug = PETSC_FALSE;
 65:   SNESLineSearch linesearch;
 66:   SNESNCGType    ncgtype=ncg->type;

 69:   PetscOptionsHead(PetscOptionsObject,"SNES NCG options");
 70:   PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
 71:   if (debug) {
 72:     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 73:   }
 74:   PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
 75:   SNESNCGSetType(snes, ncgtype);
 76:   PetscOptionsTail();
 77:   if (!snes->linesearch) {
 78:     SNESGetLineSearch(snes, &linesearch);
 79:     if (!snes->npc) {
 80:       SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
 81:     } else {
 82:       SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
 83:     }
 84:   }
 85:   SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);
 86:   return(0);
 87: }

 89: /*
 90:   SNESView_NCG - Prints info from the SNESNCG data structure.

 92:   Input Parameters:
 93: + SNES - the SNES context
 94: - viewer - visualization context

 96:   Application Interface Routine: SNESView()
 97: */
 98: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
 99: {
100:   SNES_NCG      *ncg = (SNES_NCG *) snes->data;
101:   PetscBool      iascii;

105:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
106:   if (iascii) {
107:     PetscViewerASCIIPrintf(viewer, "  type: %s\n", SNESNCGTypes[ncg->type]);
108:   }
109:   return(0);
110: }

112: PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
113: {
114:   PetscScalar    alpha, ptAp;
115:   Vec            X, Y, F, W;
116:   SNES           snes;
118:   PetscReal      *fnorm, *xnorm, *ynorm;

121:   SNESLineSearchGetSNES(linesearch, &snes);
122:   X     = linesearch->vec_sol;
123:   W     = linesearch->vec_sol_new;
124:   F     = linesearch->vec_func;
125:   Y     = linesearch->vec_update;
126:   fnorm = &linesearch->fnorm;
127:   xnorm = &linesearch->xnorm;
128:   ynorm = &linesearch->ynorm;

130:   if (!snes->jacobian) {
131:     SNESSetUpMatrices(snes);
132:   }

134:   /*

136:    The exact step size for unpreconditioned linear CG is just:
137:    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
138:    */
139:   SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);
140:   VecDot(F, F, &alpha);
141:   MatMult(snes->jacobian, Y, W);
142:   VecDot(Y, W, &ptAp);
143:   alpha = alpha / ptAp;
144:   VecAXPY(X,-alpha,Y);
145:   SNESComputeFunction(snes,X,F);

147:   VecNorm(F,NORM_2,fnorm);
148:   VecNorm(X,NORM_2,xnorm);
149:   VecNorm(Y,NORM_2,ynorm);
150:   return(0);
151: }

153: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
154: {
156:   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
157:   linesearch->ops->destroy        = NULL;
158:   linesearch->ops->setfromoptions = NULL;
159:   linesearch->ops->reset          = NULL;
160:   linesearch->ops->view           = NULL;
161:   linesearch->ops->setup          = NULL;
162:   return(0);
163: }

165: /*

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

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

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

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

190:     Logically Collective on SNES

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

196:     Options Database:
197: .   -snes_ncg_type<prp,fr,hs,dy,cd>

199:     Level: intermediate

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

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

211: .keywords: SNES, SNESNCG, selection, type, set
212: @*/
213: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
214: {

219:   PetscTryMethod(snes,"SNESNCGSetType_C",(SNES,SNESNCGType),(snes,btype));
220:   return(0);
221: }

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

228:   ncg->type = btype;
229:   return(0);
230: }

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

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

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

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

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

258:   PetscCitationsRegister(SNESCitation,&SNEScite);
259:   snes->reason = SNES_CONVERGED_ITERATING;

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

268:   SNESGetLineSearch(snes, &linesearch);

270:   PetscObjectSAWsTakeAccess((PetscObject)snes);
271:   snes->iter = 0;
272:   snes->norm = 0.;
273:   PetscObjectSAWsGrantAccess((PetscObject)snes);

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

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

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


310:   PetscObjectSAWsTakeAccess((PetscObject)snes);
311:   snes->norm = fnorm;
312:   PetscObjectSAWsGrantAccess((PetscObject)snes);
313:   SNESLogConvergenceHistory(snes,fnorm,0);
314:   SNESMonitor(snes,0,fnorm);

316:   /* test convergence */
317:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
318:   if (snes->reason) return(0);

320:   /* Call general purpose update function */
321:   if (snes->ops->update) {
322:     (*snes->ops->update)(snes, snes->iter);
323:   }

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

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

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

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

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



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

441:   Level: beginner

443:   Options Database:
444: +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter, default is prp.
445: .   -snes_linesearch_type <cp,l2,basic> - Line search type.
446: -   -snes_ncg_monitor - Print relevant information about the ncg iteration.

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

452:           Only supports left non-linear preconditioning.

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


459: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
460: M*/
461: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
462: {
464:   SNES_NCG       * neP;

467:   snes->ops->destroy        = SNESDestroy_NCG;
468:   snes->ops->setup          = SNESSetUp_NCG;
469:   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
470:   snes->ops->view           = SNESView_NCG;
471:   snes->ops->solve          = SNESSolve_NCG;
472:   snes->ops->reset          = SNESReset_NCG;

474:   snes->usesksp = PETSC_FALSE;
475:   snes->usesnpc = PETSC_TRUE;
476:   snes->npcside = PC_LEFT;

478:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

480:   if (!snes->tolerancesset) {
481:     snes->max_funcs = 30000;
482:     snes->max_its   = 10000;
483:     snes->stol      = 1e-20;
484:   }

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