Actual source code: snesrichardson.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: #include <../src/snes/impls/richardson/snesrichardsonimpl.h>


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

 12: /*
 13:   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().

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

 18:   Application Interface Routine: SNESDestroy()
 19: */
 22: PetscErrorCode SNESDestroy_NRichardson(SNES snes)
 23: {

 27:   SNESReset_NRichardson(snes);
 28:   PetscFree(snes->data);
 29:   return(0);
 30: }

 32: /*
 33:    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
 34:    of the SNESNRICHARDSON nonlinear solver.

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

 40:    Application Interface Routine: SNESSetUp()
 41:  */
 44: PetscErrorCode SNESSetUp_NRichardson(SNES snes)
 45: {
 47:   if (snes->pcside == PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");}
 48:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 49:   return(0);
 50: }

 52: /*
 53:   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.

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

 58:   Application Interface Routine: SNESSetFromOptions()
 59: */
 62: static PetscErrorCode SNESSetFromOptions_NRichardson(PetscOptionItems *PetscOptionsObject,SNES snes)
 63: {
 65:   SNESLineSearch linesearch;

 68:   PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");
 69:   PetscOptionsTail();
 70:   if (!snes->linesearch) {
 71:     SNESGetLineSearch(snes, &linesearch);
 72:     SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
 73:   }
 74:   return(0);
 75: }

 77: /*
 78:   SNESView_NRichardson - Prints info from the SNESRichardson data structure.

 80:   Input Parameters:
 81: + SNES - the SNES context
 82: - viewer - visualization context

 84:   Application Interface Routine: SNESView()
 85: */
 88: static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
 89: {
 90:   PetscBool      iascii;

 94:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 95:   if (iascii) {
 96:   }
 97:   return(0);
 98: }

100: /*
101:   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.

103:   Input Parameters:
104: . snes - the SNES context

106:   Output Parameter:
107: . outits - number of iterations until termination

109:   Application Interface Routine: SNESSolve()
110: */
113: PetscErrorCode SNESSolve_NRichardson(SNES snes)
114: {
115:   Vec                  X, Y, F;
116:   PetscReal            xnorm, fnorm, ynorm;
117:   PetscInt             maxits, i;
118:   PetscErrorCode       ierr;
119:   SNESLineSearchReason lsresult;
120:   SNESConvergedReason  reason;

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

125:   snes->reason = SNES_CONVERGED_ITERATING;

127:   maxits = snes->max_its;        /* maximum number of iterations */
128:   X      = snes->vec_sol;        /* X^n */
129:   Y      = snes->vec_sol_update; /* \tilde X */
130:   F      = snes->vec_func;       /* residual vector */

132:   PetscObjectSAWsTakeAccess((PetscObject)snes);
133:   snes->iter = 0;
134:   snes->norm = 0.;
135:   PetscObjectSAWsGrantAccess((PetscObject)snes);

137:   if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
138:     SNESApplyNPC(snes,X,NULL,F);
139:     SNESGetConvergedReason(snes->pc,&reason);
140:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
141:       snes->reason = SNES_DIVERGED_INNER;
142:       return(0);
143:     }
144:     VecNorm(F,NORM_2,&fnorm);
145:   } else {
146:     if (!snes->vec_func_init_set) {
147:       SNESComputeFunction(snes,X,F);
148:     } else snes->vec_func_init_set = PETSC_FALSE;

150:     VecNorm(F,NORM_2,&fnorm);
151:     SNESCheckFunctionNorm(snes,fnorm);
152:   }
153:   if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
154:       SNESApplyNPC(snes,X,F,Y);
155:       SNESGetConvergedReason(snes->pc,&reason);
156:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
157:         snes->reason = SNES_DIVERGED_INNER;
158:         return(0);
159:       }
160:   } else {
161:     VecCopy(F,Y);
162:   }

164:   PetscObjectSAWsTakeAccess((PetscObject)snes);
165:   snes->norm = fnorm;
166:   PetscObjectSAWsGrantAccess((PetscObject)snes);
167:   SNESLogConvergenceHistory(snes,fnorm,0);
168:   SNESMonitor(snes,0,fnorm);

170:   /* test convergence */
171:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
172:   if (snes->reason) return(0);

174:   /* Call general purpose update function */
175:   if (snes->ops->update) {
176:     (*snes->ops->update)(snes, snes->iter);
177:   }

179:   /* set parameter for default relative tolerance convergence test */
180:   snes->ttol = fnorm*snes->rtol;
181:   /* test convergence */
182:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
183:   if (snes->reason) return(0);

185:   for (i = 1; i < maxits+1; i++) {
186:     SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
187:     SNESLineSearchGetReason(snes->linesearch, &lsresult);
188:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
189:     if (lsresult) {
190:       if (++snes->numFailures >= snes->maxFailures) {
191:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
192:         break;
193:       }
194:     }
195:     if (snes->nfuncs >= snes->max_funcs) {
196:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
197:       break;
198:     }

200:     /* Monitor convergence */
201:     PetscObjectSAWsTakeAccess((PetscObject)snes);
202:     snes->iter = i;
203:     snes->norm = fnorm;
204:     PetscObjectSAWsGrantAccess((PetscObject)snes);
205:     SNESLogConvergenceHistory(snes,snes->norm,0);
206:     SNESMonitor(snes,snes->iter,snes->norm);
207:     /* Test for convergence */
208:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
209:     if (snes->reason) break;

211:     /* Call general purpose update function */
212:     if (snes->ops->update) {
213:       (*snes->ops->update)(snes, snes->iter);
214:     }

216:     if (snes->pc) {
217:       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
218:         SNESApplyNPC(snes,X,NULL,Y);
219:         VecNorm(F,NORM_2,&fnorm);
220:         VecCopy(Y,F);
221:       } else {
222:         SNESApplyNPC(snes,X,F,Y);
223:       }
224:       SNESGetConvergedReason(snes->pc,&reason);
225:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
226:         snes->reason = SNES_DIVERGED_INNER;
227:         return(0);
228:       }
229:     } else {
230:       VecCopy(F,Y);
231:     }
232:   }
233:   if (i == maxits+1) {
234:     PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
235:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
236:   }
237:   return(0);
238: }

240: /*MC
241:   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.

243:   Level: beginner

245:   Options Database:
246: +   -snes_linesearch_type <l2,cp,basic> Line search type.
247: -   -snes_linesearch_damping<1.0> Damping for the line search.

249:   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
250:             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
251:             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
252:             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
253:             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.

255:             The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
256:             linesearch, one may have to scale the update with -snes_linesearch_damping

258:      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls

260: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
261: M*/
264: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
265: {
266:   PetscErrorCode   ierr;
267:   SNES_NRichardson *neP;

270:   snes->ops->destroy        = SNESDestroy_NRichardson;
271:   snes->ops->setup          = SNESSetUp_NRichardson;
272:   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
273:   snes->ops->view           = SNESView_NRichardson;
274:   snes->ops->solve          = SNESSolve_NRichardson;
275:   snes->ops->reset          = SNESReset_NRichardson;

277:   snes->usesksp = PETSC_FALSE;
278:   snes->usespc  = PETSC_TRUE;

280:   snes->pcside = PC_LEFT;

282:   PetscNewLog(snes,&neP);
283:   snes->data = (void*) neP;

285:   if (!snes->tolerancesset) {
286:     snes->max_funcs = 30000;
287:     snes->max_its   = 10000;
288:     snes->stol      = 1e-20;
289:   }
290:   return(0);
291: }