Actual source code: snesncg.c

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

  6: PetscErrorCode SNESReset_NCG(SNES snes)
  7: {
  9:   return(0);
 10: }

 12: #define SNESLINESEARCHNCGLINEAR "ncglinear"

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

 17:   Input Parameter:
 18: . snes - the SNES context

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

 29:   PetscFree(snes->data);
 30:   return(0);
 31: }

 33: /*
 34:    SNESSetUp_NCG - Sets up the internal data structures for the later use
 35:    of the SNESNCG nonlinear solver.

 37:    Input Parameters:
 38: +  snes - the SNES context
 39: -  x - the solution vector

 41:    Application Interface Routine: SNESSetUp()
 42:  */

 44: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch);

 48: PetscErrorCode SNESSetUp_NCG(SNES snes)
 49: {

 53:   SNESSetWorkVecs(snes,2);
 54:   if (snes->pcside == PC_RIGHT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNESNCG only supports left preconditioning");
 55:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 56:   return(0);
 57: }
 58: /*
 59:   SNESSetFromOptions_NCG - Sets various parameters for the SNESNCG method.

 61:   Input Parameter:
 62: . snes - the SNES context

 64:   Application Interface Routine: SNESSetFromOptions()
 65: */
 68: static PetscErrorCode SNESSetFromOptions_NCG(SNES snes)
 69: {
 70:   SNES_NCG       *ncg = (SNES_NCG*)snes->data;
 72:   PetscBool      debug;
 73:   SNESLineSearch linesearch;
 74:   SNESNCGType    ncgtype=ncg->type;

 77:   PetscOptionsHead("SNES NCG options");
 78:   PetscOptionsBool("-snes_ncg_monitor","Monitor NCG iterations","SNES",ncg->monitor ? PETSC_TRUE : PETSC_FALSE, &debug, NULL);
 79:   PetscOptionsEnum("-snes_ncg_type","NCG Beta type used","SNESNCGSetType",SNESNCGTypes,(PetscEnum)ncg->type,(PetscEnum*)&ncgtype,NULL);
 80:   SNESNCGSetType(snes, ncgtype);
 81:   if (debug) {
 82:     ncg->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)snes));
 83:   }
 84:   PetscOptionsTail();
 85:   if (!snes->linesearch) {
 86:     SNESGetLineSearch(snes, &linesearch);
 87:     if (!snes->pc) {
 88:       SNESLineSearchSetType(linesearch, SNESLINESEARCHCP);
 89:     } else {
 90:       SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
 91:     }
 92:   }
 93:   SNESLineSearchRegister(SNESLINESEARCHNCGLINEAR, SNESLineSearchCreate_NCGLinear);
 94:   return(0);
 95: }

 97: /*
 98:   SNESView_NCG - Prints info from the SNESNCG data structure.

100:   Input Parameters:
101: + SNES - the SNES context
102: - viewer - visualization context

104:   Application Interface Routine: SNESView()
105: */
108: static PetscErrorCode SNESView_NCG(SNES snes, PetscViewer viewer)
109: {
110:   PetscBool      iascii;

114:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
115:   if (iascii) {
116:   }
117:   return(0);
118: }

122: PetscErrorCode SNESLineSearchApply_NCGLinear(SNESLineSearch linesearch)
123: {
124:   PetscScalar    alpha, ptAp;
125:   Vec            X, Y, F, W;
126:   SNES           snes;
128:   PetscReal      *fnorm, *xnorm, *ynorm;

131:   SNESLineSearchGetSNES(linesearch, &snes);
132:   X     = linesearch->vec_sol;
133:   W     = linesearch->vec_sol_new;
134:   F     = linesearch->vec_func;
135:   Y     = linesearch->vec_update;
136:   fnorm = &linesearch->fnorm;
137:   xnorm = &linesearch->xnorm;
138:   ynorm = &linesearch->ynorm;

140:   if (!snes->jacobian) {
141:     SNESSetUpMatrices(snes);
142:   }

144:   /*

146:    The exact step size for unpreconditioned linear CG is just:
147:    alpha = (r, r) / (p, Ap) = (f, f) / (y, Jy)
148:    */
149:   SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre);
150:   VecDot(F, F, &alpha);
151:   MatMult(snes->jacobian, Y, W);
152:   VecDot(Y, W, &ptAp);
153:   alpha = alpha / ptAp;
154:   VecAXPY(X,-alpha,Y);
155:   SNESComputeFunction(snes,X,F);

157:   VecNorm(F,NORM_2,fnorm);
158:   VecNorm(X,NORM_2,xnorm);
159:   VecNorm(Y,NORM_2,ynorm);
160:   return(0);
161: }

165: PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NCGLinear(SNESLineSearch linesearch)
166: {
168:   linesearch->ops->apply          = SNESLineSearchApply_NCGLinear;
169:   linesearch->ops->destroy        = NULL;
170:   linesearch->ops->setfromoptions = NULL;
171:   linesearch->ops->reset          = NULL;
172:   linesearch->ops->view           = NULL;
173:   linesearch->ops->setup          = NULL;
174:   return(0);
175: }

179: /*

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

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

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

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

206:     Logically Collective on SNES

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

212:     Options Database:
213: .   -snes_ncg_type<prp,fr,hs,dy,cd>

215:     Level: intermediate

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

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

227: .keywords: SNES, SNESNCG, selection, type, set
228: @*/
229: PetscErrorCode SNESNCGSetType(SNES snes, SNESNCGType btype)
230: {

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

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

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

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

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

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

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

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

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

286:   SNESGetLineSearch(snes, &linesearch);

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

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

295:   if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
296:     SNESApplyNPC(snes,X,NULL,dX);
297:     SNESGetConvergedReason(snes->pc,&reason);
298:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
299:       snes->reason = SNES_DIVERGED_INNER;
300:       return(0);
301:     }
302:     VecCopy(dX,F);
303:     VecNorm(F,NORM_2,&fnorm);
304:   } else {
305:     if (!snes->vec_func_init_set) {
306:       SNESComputeFunction(snes,X,F);
307:       if (snes->domainerror) {
308:         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
309:         return(0);
310:       }
311:     } else snes->vec_func_init_set = PETSC_FALSE;

313:     /* convergence test */
314:     VecNorm(F,NORM_2,&fnorm);
315:     if (PetscIsInfOrNanReal(fnorm)) {
316:       snes->reason = SNES_DIVERGED_FNORM_NAN;
317:       return(0);
318:     }

320:     VecCopy(F,dX);
321:   }
322:   if (snes->pc) {
323:     if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
324:       SNESApplyNPC(snes,X,F,dX);
325:       SNESGetConvergedReason(snes->pc,&reason);
326:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
327:         snes->reason = SNES_DIVERGED_INNER;
328:         return(0);
329:       }
330:     }
331:   }
332:   VecCopy(dX,lX);
333:   VecDot(dX, dX, &dXdotdX);


336:   PetscObjectSAWsTakeAccess((PetscObject)snes);
337:   snes->norm = fnorm;
338:   PetscObjectSAWsGrantAccess((PetscObject)snes);
339:   SNESLogConvergenceHistory(snes,fnorm,0);
340:   SNESMonitor(snes,0,fnorm);

342:   /* test convergence */
343:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
344:   if (snes->reason) return(0);

346:   /* Call general purpose update function */
347:   if (snes->ops->update) {
348:     (*snes->ops->update)(snes, snes->iter);
349:   }

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

353:   for (i = 1; i < maxits + 1; i++) {
354:     lsSuccess = PETSC_TRUE;
355:     /* some update types require the old update direction or conjugate direction */
356:     if (ncg->type != SNES_NCG_FR) {
357:       VecCopy(dX, dXold);
358:     }
359:     SNESLineSearchApply(linesearch,X,F,&fnorm,lX);
360:     SNESLineSearchGetSuccess(linesearch, &lsSuccess);
361:     if (!lsSuccess) {
362:       if (++snes->numFailures >= snes->maxFailures) {
363:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
364:         return(0);
365:       }
366:     }
367:     if (snes->nfuncs >= snes->max_funcs) {
368:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
369:       return(0);
370:     }
371:     if (snes->domainerror) {
372:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
373:       return(0);
374:     }
375:     SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);
376:     /* Monitor convergence */
377:     PetscObjectSAWsTakeAccess((PetscObject)snes);
378:     snes->iter = i;
379:     snes->norm = fnorm;
380:     PetscObjectSAWsGrantAccess((PetscObject)snes);
381:     SNESLogConvergenceHistory(snes,snes->norm,0);
382:     SNESMonitor(snes,snes->iter,snes->norm);

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

388:     /* Call general purpose update function */
389:     if (snes->ops->update) {
390:       (*snes->ops->update)(snes, snes->iter);
391:     }
392:     if (snes->pc) {
393:       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
394:         SNESApplyNPC(snes,X,NULL,dX);
395:         SNESGetConvergedReason(snes->pc,&reason);
396:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
397:           snes->reason = SNES_DIVERGED_INNER;
398:           return(0);
399:         }
400:         VecCopy(dX,F);
401:       } else {
402:         SNESApplyNPC(snes,X,F,dX);
403:         SNESGetConvergedReason(snes->pc,&reason);
404:         if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
405:           snes->reason = SNES_DIVERGED_INNER;
406:           return(0);
407:         }
408:       }
409:     } else {
410:       VecCopy(F,dX);
411:     }

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



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

472:   Level: beginner

474:   Options Database:
475: +   -snes_ncg_type <fr, prp, dy, hs, cd> - Choice of conjugate-gradient update parameter.
476: .   -snes_linesearch_type <cp,l2,basic> - Line search type.
477: -   -snes_ncg_monitor - Print relevant information about the ncg iteration.

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


484: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESNQN
485: M*/
488: PETSC_EXTERN PetscErrorCode SNESCreate_NCG(SNES snes)
489: {
491:   SNES_NCG       * neP;

494:   snes->ops->destroy        = SNESDestroy_NCG;
495:   snes->ops->setup          = SNESSetUp_NCG;
496:   snes->ops->setfromoptions = SNESSetFromOptions_NCG;
497:   snes->ops->view           = SNESView_NCG;
498:   snes->ops->solve          = SNESSolve_NCG;
499:   snes->ops->reset          = SNESReset_NCG;

501:   snes->usesksp = PETSC_FALSE;
502:   snes->usespc  = PETSC_TRUE;
503:   snes->pcside  = PC_LEFT;

505:   if (!snes->tolerancesset) {
506:     snes->max_funcs = 30000;
507:     snes->max_its   = 10000;
508:     snes->stol      = 1e-20;
509:   }

511:   PetscNewLog(snes,&neP);
512:   snes->data   = (void*) neP;
513:   neP->monitor = NULL;
514:   neP->type    = SNES_NCG_PRP;
515:   PetscObjectComposeFunction((PetscObject)snes,"SNESNCGSetType_C", SNESNCGSetType_NCG);
516:   return(0);
517: }