Actual source code: snesrichardson.c

petsc-3.4.5 2014-06-29
  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:   return(0);
 48: }

 50: /*
 51:   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.

 53:   Input Parameter:
 54: . snes - the SNES context

 56:   Application Interface Routine: SNESSetFromOptions()
 57: */
 60: static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
 61: {
 63:   SNESLineSearch linesearch;

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

 75: /*
 76:   SNESView_NRichardson - Prints info from the SNESRichardson data structure.

 78:   Input Parameters:
 79: + SNES - the SNES context
 80: - viewer - visualization context

 82:   Application Interface Routine: SNESView()
 83: */
 86: static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
 87: {
 88:   PetscBool      iascii;

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

 98: /*
 99:   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.

101:   Input Parameters:
102: . snes - the SNES context

104:   Output Parameter:
105: . outits - number of iterations until termination

107:   Application Interface Routine: SNESSolve()
108: */
111: PetscErrorCode SNESSolve_NRichardson(SNES snes)
112: {
113:   Vec                 X, Y, F;
114:   PetscReal           xnorm, fnorm, ynorm;
115:   PetscInt            maxits, i;
116:   PetscErrorCode      ierr;
117:   PetscBool           lsSuccess;
118:   SNESConvergedReason reason;

121:   snes->reason = SNES_CONVERGED_ITERATING;

123:   maxits = snes->max_its;        /* maximum number of iterations */
124:   X      = snes->vec_sol;        /* X^n */
125:   Y      = snes->vec_sol_update; /* \tilde X */
126:   F      = snes->vec_func;       /* residual vector */

128:   PetscObjectAMSTakeAccess((PetscObject)snes);
129:   snes->iter = 0;
130:   snes->norm = 0.;
131:   PetscObjectAMSGrantAccess((PetscObject)snes);
132:   if (!snes->vec_func_init_set) {
133:     SNESComputeFunction(snes,X,F);
134:     if (snes->domainerror) {
135:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
136:       return(0);
137:     }
138:   } else snes->vec_func_init_set = PETSC_FALSE;

140:   if (!snes->norm_init_set) {
141:     VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F||  */
142:     if (PetscIsInfOrNanReal(fnorm)) {
143:       snes->reason = SNES_DIVERGED_FNORM_NAN;
144:       return(0);
145:     }
146:   } else {
147:     fnorm               = snes->norm_init;
148:     snes->norm_init_set = PETSC_FALSE;
149:   }

151:   PetscObjectAMSTakeAccess((PetscObject)snes);
152:   snes->norm = fnorm;
153:   PetscObjectAMSGrantAccess((PetscObject)snes);
154:   SNESLogConvergenceHistory(snes,fnorm,0);
155:   SNESMonitor(snes,0,fnorm);

157:   /* set parameter for default relative tolerance convergence test */
158:   snes->ttol = fnorm*snes->rtol;
159:   /* test convergence */
160:   (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
161:   if (snes->reason) return(0);

163:   for (i = 0; i < maxits; i++) {
164:     lsSuccess = PETSC_TRUE;
165:     /* Call general purpose update function */
166:     if (snes->ops->update) {
167:       (*snes->ops->update)(snes, snes->iter);
168:     }
169:     if (snes->pc && snes->pcside == PC_RIGHT) {
170:       VecCopy(X,Y);
171:       SNESSetInitialFunction(snes->pc, F);
172:       SNESSetInitialFunctionNorm(snes->pc, fnorm);
173:       PetscLogEventBegin(SNES_NPCSolve,snes->pc,Y,snes->vec_rhs,0);
174:       SNESSolve(snes->pc, snes->vec_rhs, Y);
175:       PetscLogEventEnd(SNES_NPCSolve,snes->pc,Y,snes->vec_rhs,0);
176:       SNESGetConvergedReason(snes->pc,&reason);
177:       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
178:         snes->reason = SNES_DIVERGED_INNER;
179:         return(0);
180:       }
181:       VecAYPX(Y,-1.0,X);
182:     } else {
183:       VecCopy(F,Y);
184:     }
185:     SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
186:     SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
187:     SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);
188:     if (!lsSuccess) {
189:       if (++snes->numFailures >= snes->maxFailures) {
190:         snes->reason = SNES_DIVERGED_LINE_SEARCH;
191:         break;
192:       }
193:     }
194:     if (snes->nfuncs >= snes->max_funcs) {
195:       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
196:       break;
197:     }
198:     if (snes->domainerror) {
199:       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
200:       return(0);
201:     }
202:     /* Monitor convergence */
203:     PetscObjectAMSTakeAccess((PetscObject)snes);
204:     snes->iter = i+1;
205:     snes->norm = fnorm;
206:     PetscObjectAMSGrantAccess((PetscObject)snes);
207:     SNESLogConvergenceHistory(snes,snes->norm,0);
208:     SNESMonitor(snes,snes->iter,snes->norm);
209:     /* Test for convergence */
210:     (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
211:     if (snes->reason) break;
212:   }
213:   if (i == maxits) {
214:     PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
215:     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
216:   }
217:   return(0);
218: }

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

223:   Level: beginner

225:   Options Database:
226: +   -snes_linesearch_type <l2,cp,basic> Line search type.
227: -   -snes_linesearch_damping<1.0> Damping for the line search.

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

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

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

240: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
241: M*/
244: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
245: {
246:   PetscErrorCode   ierr;
247:   SNES_NRichardson *neP;

250:   snes->ops->destroy        = SNESDestroy_NRichardson;
251:   snes->ops->setup          = SNESSetUp_NRichardson;
252:   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
253:   snes->ops->view           = SNESView_NRichardson;
254:   snes->ops->solve          = SNESSolve_NRichardson;
255:   snes->ops->reset          = SNESReset_NRichardson;

257:   snes->usesksp = PETSC_FALSE;
258:   snes->usespc  = PETSC_TRUE;

260:   PetscNewLog(snes, SNES_NRichardson, &neP);
261:   snes->data = (void*) neP;

263:   if (!snes->tolerancesset) {
264:     snes->max_funcs = 30000;
265:     snes->max_its   = 10000;
266:     snes->stol      = 1e-20;
267:   }
268:   return(0);
269: }