Actual source code: snesrichardson.c
petsc-3.6.4 2016-04-12
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(PetscOptions *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;
124: if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
125: SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
126: }
128: snes->reason = SNES_CONVERGED_ITERATING;
130: maxits = snes->max_its; /* maximum number of iterations */
131: X = snes->vec_sol; /* X^n */
132: Y = snes->vec_sol_update; /* \tilde X */
133: F = snes->vec_func; /* residual vector */
135: PetscObjectSAWsTakeAccess((PetscObject)snes);
136: snes->iter = 0;
137: snes->norm = 0.;
138: PetscObjectSAWsGrantAccess((PetscObject)snes);
140: if (snes->pc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
141: SNESApplyNPC(snes,X,NULL,F);
142: SNESGetConvergedReason(snes->pc,&reason);
143: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
144: snes->reason = SNES_DIVERGED_INNER;
145: return(0);
146: }
147: VecNorm(F,NORM_2,&fnorm);
148: } else {
149: if (!snes->vec_func_init_set) {
150: SNESComputeFunction(snes,X,F);
151: } else snes->vec_func_init_set = PETSC_FALSE;
153: VecNorm(F,NORM_2,&fnorm);
154: SNESCheckFunctionNorm(snes,fnorm);
155: }
156: if (snes->pc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
157: SNESApplyNPC(snes,X,F,Y);
158: SNESGetConvergedReason(snes->pc,&reason);
159: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
160: snes->reason = SNES_DIVERGED_INNER;
161: return(0);
162: }
163: } else {
164: VecCopy(F,Y);
165: }
167: PetscObjectSAWsTakeAccess((PetscObject)snes);
168: snes->norm = fnorm;
169: PetscObjectSAWsGrantAccess((PetscObject)snes);
170: SNESLogConvergenceHistory(snes,fnorm,0);
171: SNESMonitor(snes,0,fnorm);
173: /* test convergence */
174: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
175: if (snes->reason) return(0);
177: /* Call general purpose update function */
178: if (snes->ops->update) {
179: (*snes->ops->update)(snes, snes->iter);
180: }
182: /* set parameter for default relative tolerance convergence test */
183: snes->ttol = fnorm*snes->rtol;
184: /* test convergence */
185: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
186: if (snes->reason) return(0);
188: for (i = 1; i < maxits+1; i++) {
189: SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
190: SNESLineSearchGetReason(snes->linesearch, &lsresult);
191: SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
192: if (lsresult) {
193: if (++snes->numFailures >= snes->maxFailures) {
194: snes->reason = SNES_DIVERGED_LINE_SEARCH;
195: break;
196: }
197: }
198: if (snes->nfuncs >= snes->max_funcs) {
199: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
200: break;
201: }
203: /* Monitor convergence */
204: PetscObjectSAWsTakeAccess((PetscObject)snes);
205: snes->iter = i;
206: snes->norm = fnorm;
207: PetscObjectSAWsGrantAccess((PetscObject)snes);
208: SNESLogConvergenceHistory(snes,snes->norm,0);
209: SNESMonitor(snes,snes->iter,snes->norm);
210: /* Test for convergence */
211: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
212: if (snes->reason) break;
214: /* Call general purpose update function */
215: if (snes->ops->update) {
216: (*snes->ops->update)(snes, snes->iter);
217: }
219: if (snes->pc) {
220: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
221: SNESApplyNPC(snes,X,NULL,Y);
222: VecNorm(F,NORM_2,&fnorm);
223: VecCopy(Y,F);
224: } else {
225: SNESApplyNPC(snes,X,F,Y);
226: }
227: SNESGetConvergedReason(snes->pc,&reason);
228: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
229: snes->reason = SNES_DIVERGED_INNER;
230: return(0);
231: }
232: } else {
233: VecCopy(F,Y);
234: }
235: }
236: if (i == maxits+1) {
237: PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
238: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
239: }
240: return(0);
241: }
243: /*MC
244: SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
246: Level: beginner
248: Options Database:
249: + -snes_linesearch_type <l2,cp,basic> Line search type.
250: - -snes_linesearch_damping<1.0> Damping for the line search.
252: Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
253: (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If
254: an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
255: solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
256: where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
258: The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic
259: linesearch, one may have to scale the update with -snes_linesearch_damping
261: This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
263: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
264: M*/
267: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
268: {
269: PetscErrorCode ierr;
270: SNES_NRichardson *neP;
273: snes->ops->destroy = SNESDestroy_NRichardson;
274: snes->ops->setup = SNESSetUp_NRichardson;
275: snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
276: snes->ops->view = SNESView_NRichardson;
277: snes->ops->solve = SNESSolve_NRichardson;
278: snes->ops->reset = SNESReset_NRichardson;
280: snes->usesksp = PETSC_FALSE;
281: snes->usespc = PETSC_TRUE;
283: snes->pcside = PC_LEFT;
285: PetscNewLog(snes,&neP);
286: snes->data = (void*) neP;
288: if (!snes->tolerancesset) {
289: snes->max_funcs = 30000;
290: snes->max_its = 10000;
291: snes->stol = 1e-20;
292: }
293: return(0);
294: }