Actual source code: snesncg.c

petsc-master 2019-07-16
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:   SNESNCGType    ncgtype=ncg->type;
 66:   SNESLineSearch linesearch;

 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: @*/
212: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
213: {

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

222: PetscErrorCode SNESNCGSetType_NCG(SNES snes, SNESNCGType btype)
223: {
224:   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: 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:   PetscErrorCode       ierr;
250:   SNESLineSearchReason lsresult = SNES_LINESEARCH_SUCCEEDED;
251:   SNESLineSearch       linesearch;
252:   SNESConvergedReason  reason;

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

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

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

267:   SNESGetLineSearch(snes, &linesearch);

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

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

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

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


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

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

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

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

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

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

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

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



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

442:   Level: beginner

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

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

454:           Only supports left non-linear preconditioning.

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


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

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

476:   snes->usesksp = PETSC_FALSE;
477:   snes->usesnpc = PETSC_TRUE;
478:   snes->npcside = PC_LEFT;

480:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

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

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