Actual source code: snesrichardson.c
petsc-3.10.5 2019-03-28
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: {
57: SNESLineSearch linesearch;
60: PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");
61: PetscOptionsTail();
62: if (!snes->linesearch) {
63: SNESGetLineSearch(snes, &linesearch);
64: SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);
65: }
66: return(0);
67: }
69: /*
70: SNESView_NRichardson - Prints info from the SNESRichardson data structure.
72: Input Parameters:
73: + SNES - the SNES context
74: - viewer - visualization context
76: Application Interface Routine: SNESView()
77: */
78: static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
79: {
80: PetscBool iascii;
84: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
85: if (iascii) {
86: }
87: return(0);
88: }
90: /*
91: SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
93: Input Parameters:
94: . snes - the SNES context
96: Output Parameter:
97: . outits - number of iterations until termination
99: Application Interface Routine: SNESSolve()
100: */
101: PetscErrorCode SNESSolve_NRichardson(SNES snes)
102: {
103: Vec X, Y, F;
104: PetscReal xnorm, fnorm, ynorm;
105: PetscInt maxits, i;
106: PetscErrorCode ierr;
107: SNESLineSearchReason lsresult;
108: SNESConvergedReason reason;
111: 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);
113: snes->reason = SNES_CONVERGED_ITERATING;
115: maxits = snes->max_its; /* maximum number of iterations */
116: X = snes->vec_sol; /* X^n */
117: Y = snes->vec_sol_update; /* \tilde X */
118: F = snes->vec_func; /* residual vector */
120: PetscObjectSAWsTakeAccess((PetscObject)snes);
121: snes->iter = 0;
122: snes->norm = 0.;
123: PetscObjectSAWsGrantAccess((PetscObject)snes);
125: if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
126: SNESApplyNPC(snes,X,NULL,F);
127: SNESGetConvergedReason(snes->npc,&reason);
128: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
129: snes->reason = SNES_DIVERGED_INNER;
130: return(0);
131: }
132: VecNorm(F,NORM_2,&fnorm);
133: } else {
134: if (!snes->vec_func_init_set) {
135: SNESComputeFunction(snes,X,F);
136: } else snes->vec_func_init_set = PETSC_FALSE;
138: VecNorm(F,NORM_2,&fnorm);
139: SNESCheckFunctionNorm(snes,fnorm);
140: }
141: if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
142: SNESApplyNPC(snes,X,F,Y);
143: SNESGetConvergedReason(snes->npc,&reason);
144: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
145: snes->reason = SNES_DIVERGED_INNER;
146: return(0);
147: }
148: } else {
149: VecCopy(F,Y);
150: }
152: PetscObjectSAWsTakeAccess((PetscObject)snes);
153: snes->norm = fnorm;
154: PetscObjectSAWsGrantAccess((PetscObject)snes);
155: SNESLogConvergenceHistory(snes,fnorm,0);
156: SNESMonitor(snes,0,fnorm);
158: /* test convergence */
159: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
160: if (snes->reason) return(0);
162: /* Call general purpose update function */
163: if (snes->ops->update) {
164: (*snes->ops->update)(snes, snes->iter);
165: }
167: /* set parameter for default relative tolerance convergence test */
168: snes->ttol = fnorm*snes->rtol;
169: /* test convergence */
170: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
171: if (snes->reason) return(0);
173: for (i = 1; i < maxits+1; i++) {
174: SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);
175: SNESLineSearchGetReason(snes->linesearch, &lsresult);
176: SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);
177: if (lsresult) {
178: if (++snes->numFailures >= snes->maxFailures) {
179: snes->reason = SNES_DIVERGED_LINE_SEARCH;
180: break;
181: }
182: }
183: if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
184: snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
185: break;
186: }
188: /* Monitor convergence */
189: PetscObjectSAWsTakeAccess((PetscObject)snes);
190: snes->iter = i;
191: snes->norm = fnorm;
192: PetscObjectSAWsGrantAccess((PetscObject)snes);
193: SNESLogConvergenceHistory(snes,snes->norm,0);
194: SNESMonitor(snes,snes->iter,snes->norm);
195: /* Test for convergence */
196: (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);
197: if (snes->reason) break;
199: /* Call general purpose update function */
200: if (snes->ops->update) {
201: (*snes->ops->update)(snes, snes->iter);
202: }
204: if (snes->npc) {
205: if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
206: SNESApplyNPC(snes,X,NULL,Y);
207: VecNorm(F,NORM_2,&fnorm);
208: VecCopy(Y,F);
209: } else {
210: SNESApplyNPC(snes,X,F,Y);
211: }
212: SNESGetConvergedReason(snes->npc,&reason);
213: if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
214: snes->reason = SNES_DIVERGED_INNER;
215: return(0);
216: }
217: } else {
218: VecCopy(F,Y);
219: }
220: }
221: if (i == maxits+1) {
222: PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);
223: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
224: }
225: return(0);
226: }
228: /*MC
229: SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
231: Level: beginner
233: Options Database:
234: + -snes_linesearch_type <l2,cp,basic> Line search type.
235: - -snes_linesearch_damping<1.0> Damping for the line search.
237: Notes:
238: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
239: (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search. If
240: an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
241: solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
242: where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
244: The update, especially without inner nonlinear preconditioner, may be ill-scaled. If using the basic
245: linesearch, one may have to scale the update with -snes_linesearch_damping
247: This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
249: Only supports left non-linear preconditioning.
251: .seealso: SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
252: M*/
253: PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
254: {
255: PetscErrorCode ierr;
256: SNES_NRichardson *neP;
259: snes->ops->destroy = SNESDestroy_NRichardson;
260: snes->ops->setup = SNESSetUp_NRichardson;
261: snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
262: snes->ops->view = SNESView_NRichardson;
263: snes->ops->solve = SNESSolve_NRichardson;
264: snes->ops->reset = SNESReset_NRichardson;
266: snes->usesksp = PETSC_FALSE;
267: snes->usesnpc = PETSC_TRUE;
269: snes->npcside= PC_LEFT;
271: snes->alwayscomputesfinalresidual = PETSC_TRUE;
273: PetscNewLog(snes,&neP);
274: snes->data = (void*) neP;
276: if (!snes->tolerancesset) {
277: snes->max_funcs = 30000;
278: snes->max_its = 10000;
279: snes->stol = 1e-20;
280: }
281: return(0);
282: }