Actual source code: snesrichardson.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1:  #include <../src/snes/impls/richardson/snesrichardsonimpl.h>


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

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

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

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

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

 28: /*
 29:    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
 30:    of the SNESNRICHARDSON nonlinear solver.

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

 36:    Application Interface Routine: SNESSetUp()
 37:  */
 38: PetscErrorCode SNESSetUp_NRichardson(SNES snes)
 39: {
 41:   if (snes->npcside== PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");}
 42:   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
 43:   return(0);
 44: }

 46: /*
 47:   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.

 49:   Input Parameter:
 50: . snes - the SNES context

 52:   Application Interface Routine: SNESSetFromOptions()
 53: */
 54: static PetscErrorCode SNESSetFromOptions_NRichardson(PetscOptionItems *PetscOptionsObject,SNES snes)
 55: {

 59:   PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");
 60:   PetscOptionsTail();
 61:   return(0);
 62: }

 64: /*
 65:   SNESView_NRichardson - Prints info from the SNESRichardson data structure.

 67:   Input Parameters:
 68: + SNES - the SNES context
 69: - viewer - visualization context

 71:   Application Interface Routine: SNESView()
 72: */
 73: static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
 74: {
 75:   PetscBool      iascii;

 79:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 80:   if (iascii) {
 81:   }
 82:   return(0);
 83: }

 85: /*
 86:   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.

 88:   Input Parameters:
 89: . snes - the SNES context

 91:   Output Parameter:
 92: . outits - number of iterations until termination

 94:   Application Interface Routine: SNESSolve()
 95: */
 96: PetscErrorCode SNESSolve_NRichardson(SNES snes)
 97: {
 98:   Vec                  X, Y, F;
 99:   PetscReal            xnorm, fnorm, ynorm;
100:   PetscInt             maxits, i;
101:   PetscErrorCode       ierr;
102:   SNESLineSearchReason lsresult;
103:   SNESConvergedReason  reason;

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

108:   snes->reason = SNES_CONVERGED_ITERATING;

110:   maxits = snes->max_its;        /* maximum number of iterations */
111:   X      = snes->vec_sol;        /* X^n */
112:   Y      = snes->vec_sol_update; /* \tilde X */
113:   F      = snes->vec_func;       /* residual vector */

115:   PetscObjectSAWsTakeAccess((PetscObject)snes);
116:   snes->iter = 0;
117:   snes->norm = 0.;
118:   PetscObjectSAWsGrantAccess((PetscObject)snes);

120:   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
121:     SNESApplyNPC(snes,X,NULL,F);
122:     SNESGetConvergedReason(snes->npc,&reason);
123:     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
124:       snes->reason = SNES_DIVERGED_INNER;
125:       return(0);
126:     }
127:     VecNorm(F,NORM_2,&fnorm);
128:   } else {
129:     if (!snes->vec_func_init_set) {
130:       SNESComputeFunction(snes,X,F);
131:     } else snes->vec_func_init_set = PETSC_FALSE;

133:     VecNorm(F,NORM_2,&fnorm);
134:     SNESCheckFunctionNorm(snes,fnorm);
135:   }
136:   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
137:       SNESApplyNPC(snes,X,F,Y);
138:       SNESGetConvergedReason(snes->npc,&reason);
139:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
140:         snes->reason = SNES_DIVERGED_INNER;
141:         return(0);
142:       }
143:   } else {
144:     VecCopy(F,Y);
145:   }

147:   PetscObjectSAWsTakeAccess((PetscObject)snes);
148:   snes->norm = fnorm;
149:   PetscObjectSAWsGrantAccess((PetscObject)snes);
150:   SNESLogConvergenceHistory(snes,fnorm,0);
151:   SNESMonitor(snes,0,fnorm);

153:   /* test convergence */
154:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
155:   if (snes->reason) return(0);

157:   /* Call general purpose update function */
158:   if (snes->ops->update) {
159:     (*snes->ops->update)(snes, snes->iter);
160:   }

162:   /* set parameter for default relative tolerance convergence test */
163:   snes->ttol = fnorm*snes->rtol;
164:   /* test convergence */
165:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
166:   if (snes->reason) return(0);

168:   for (i = 1; i < maxits+1; i++) {
169:     SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
170:     SNESLineSearchGetReason(snes->linesearch, &lsresult);
171:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
172:     if (lsresult) {
173:       if (++snes->numFailures >= snes->maxFailures) {
174:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
175:         break;
176:       }
177:     }
178:     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
179:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
180:       break;
181:     }

183:     /* Monitor convergence */
184:     PetscObjectSAWsTakeAccess((PetscObject)snes);
185:     snes->iter = i;
186:     snes->norm = fnorm;
187:     snes->xnorm = xnorm;
188:     snes->ynorm = ynorm;
189:     PetscObjectSAWsGrantAccess((PetscObject)snes);
190:     SNESLogConvergenceHistory(snes,snes->norm,0);
191:     SNESMonitor(snes,snes->iter,snes->norm);
192:     /* Test for convergence */
193:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
194:     if (snes->reason) break;

196:     /* Call general purpose update function */
197:     if (snes->ops->update) {
198:       (*snes->ops->update)(snes, snes->iter);
199:     }

201:     if (snes->npc) {
202:       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
203:         SNESApplyNPC(snes,X,NULL,Y);
204:         VecNorm(F,NORM_2,&fnorm);
205:         VecCopy(Y,F);
206:       } else {
207:         SNESApplyNPC(snes,X,F,Y);
208:       }
209:       SNESGetConvergedReason(snes->npc,&reason);
210:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
211:         snes->reason = SNES_DIVERGED_INNER;
212:         return(0);
213:       }
214:     } else {
215:       VecCopy(F,Y);
216:     }
217:   }
218:   if (i == maxits+1) {
219:     PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
220:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
221:   }
222:   return(0);
223: }

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

228:   Level: beginner

230:   Options Database:
231: +   -snes_linesearch_type <l2,cp,basic> Line search type.
232: -   -snes_linesearch_damping<1.0> Damping for the line search.

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

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

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

246:      Only supports left non-linear preconditioning.

248: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
249: M*/
250: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
251: {
252:   PetscErrorCode   ierr;
253:   SNES_NRichardson *neP;
254:   SNESLineSearch   linesearch;

257:   snes->ops->destroy        = SNESDestroy_NRichardson;
258:   snes->ops->setup          = SNESSetUp_NRichardson;
259:   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
260:   snes->ops->view           = SNESView_NRichardson;
261:   snes->ops->solve          = SNESSolve_NRichardson;
262:   snes->ops->reset          = SNESReset_NRichardson;

264:   snes->usesksp = PETSC_FALSE;
265:   snes->usesnpc = PETSC_TRUE;

267:   snes->npcside= PC_LEFT;

269:   SNESGetLineSearch(snes, &linesearch);
270:   if (!((PetscObject)linesearch)->type_name) {
271:     SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
272:   }

274:   snes->alwayscomputesfinalresidual = PETSC_TRUE;

276:   PetscNewLog(snes,&neP);
277:   snes->data = (void*) neP;

279:   if (!snes->tolerancesset) {
280:     snes->max_funcs = 30000;
281:     snes->max_its   = 10000;
282:     snes->stol      = 1e-20;
283:   }
284:   return(0);
285: }