Actual source code: snesrichardson.c

petsc-3.3-p7 2013-05-11
  1: #include <../src/snes/impls/richardson/snesrichardsonimpl.h>


  6: PetscErrorCode SNESReset_NRichardson(SNES snes)
  7: {

 10:   return(0);
 11: }

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

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

 19:   Application Interface Routine: SNESDestroy()
 20: */
 23: PetscErrorCode SNESDestroy_NRichardson(SNES snes)
 24: {
 25:   PetscErrorCode   ierr;

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

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

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

 41:    Application Interface Routine: SNESSetUp()
 42:  */
 45: PetscErrorCode SNESSetUp_NRichardson(SNES snes)
 46: {
 48:   return(0);
 49: }

 51: /*
 52:   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method.

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

 57:   Application Interface Routine: SNESSetFromOptions()
 58: */
 61: static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
 62: {
 64:   SNESLineSearch linesearch;
 66:   PetscOptionsHead("SNES Richardson options");
 67:   PetscOptionsTail();
 68:   if (!snes->linesearch) {
 69:     SNESGetSNESLineSearch(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;
 89:   PetscErrorCode   ierr;

 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:   PetscObjectTakeAccess(snes);
129:   snes->iter = 0;
130:   snes->norm = 0.;
131:   PetscObjectGrantAccess(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 {
139:     snes->vec_func_init_set = PETSC_FALSE;
140:   }
141:   if (!snes->norm_init_set) {
142:     VecNorm(F,NORM_2,&fnorm); /* fnorm <- ||F||  */
143:     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
144:   } else {
145:     fnorm = snes->norm_init;
146:     snes->norm_init_set = PETSC_FALSE;
147:   }

149:   PetscObjectTakeAccess(snes);
150:   snes->norm = fnorm;
151:   PetscObjectGrantAccess(snes);
152:   SNESLogConvHistory(snes,fnorm,0);
153:   SNESMonitor(snes,0,fnorm);

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

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

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

219:   Level: beginner

221:   Options Database:
222: +   -snes_linesearch_type <l2,cp,basic> Line search type.
223: -   -snes_linesearch_damping<1.0> Damping for the line search.

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

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

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

236: .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESQN, SNESNCG
237: M*/
238: EXTERN_C_BEGIN
241: PetscErrorCode  SNESCreate_NRichardson(SNES snes)
242: {
243:   PetscErrorCode          ierr;
244:   SNES_NRichardson        *neP;
246:   snes->ops->destroy         = SNESDestroy_NRichardson;
247:   snes->ops->setup           = SNESSetUp_NRichardson;
248:   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
249:   snes->ops->view            = SNESView_NRichardson;
250:   snes->ops->solve           = SNESSolve_NRichardson;
251:   snes->ops->reset           = SNESReset_NRichardson;

253:   snes->usesksp              = PETSC_FALSE;
254:   snes->usespc               = PETSC_TRUE;

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

259:   if (!snes->tolerancesset) {
260:     snes->max_funcs = 30000;
261:     snes->max_its   = 10000;
262:     snes->stol      = 1e-20;
263:   }

265:   return(0);
266: }
267: EXTERN_C_END