xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision e4ed7901070ce1ca36eb7d68dd07557542e66e31)
1 #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2 
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "SNESReset_NRichardson"
6 PetscErrorCode SNESReset_NRichardson(SNES snes)
7 {
8 
9   PetscFunctionBegin;
10   PetscFunctionReturn(0);
11 }
12 
13 /*
14   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
15 
16   Input Parameter:
17 . snes - the SNES context
18 
19   Application Interface Routine: SNESDestroy()
20 */
21 #undef __FUNCT__
22 #define __FUNCT__ "SNESDestroy_NRichardson"
23 PetscErrorCode SNESDestroy_NRichardson(SNES snes)
24 {
25   PetscErrorCode   ierr;
26 
27   PetscFunctionBegin;
28   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
29   ierr = PetscFree(snes->data);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 /*
34    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
35    of the SNESNRICHARDSON nonlinear solver.
36 
37    Input Parameters:
38 +  snes - the SNES context
39 -  x - the solution vector
40 
41    Application Interface Routine: SNESSetUp()
42  */
43 #undef __FUNCT__
44 #define __FUNCT__ "SNESSetUp_NRichardson"
45 PetscErrorCode SNESSetUp_NRichardson(SNES snes)
46 {
47   PetscFunctionBegin;
48   PetscFunctionReturn(0);
49 }
50 
51 /*
52   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method.
53 
54   Input Parameter:
55 . snes - the SNES context
56 
57   Application Interface Routine: SNESSetFromOptions()
58 */
59 #undef __FUNCT__
60 #define __FUNCT__ "SNESSetFromOptions_NRichardson"
61 static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
62 {
63   PetscErrorCode ierr;
64   SNESLineSearch linesearch;
65   PetscFunctionBegin;
66   ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr);
67   ierr = PetscOptionsTail();CHKERRQ(ierr);
68   if (!snes->linesearch) {
69     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
70     ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
71   }
72   PetscFunctionReturn(0);
73 }
74 
75 /*
76   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
77 
78   Input Parameters:
79 + SNES - the SNES context
80 - viewer - visualization context
81 
82   Application Interface Routine: SNESView()
83 */
84 #undef __FUNCT__
85 #define __FUNCT__ "SNESView_NRichardson"
86 static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
87 {
88   PetscBool        iascii;
89   PetscErrorCode   ierr;
90 
91   PetscFunctionBegin;
92   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
93   if (iascii) {
94   }
95   PetscFunctionReturn(0);
96 }
97 
98 /*
99   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
100 
101   Input Parameters:
102 . snes - the SNES context
103 
104   Output Parameter:
105 . outits - number of iterations until termination
106 
107   Application Interface Routine: SNESSolve()
108 */
109 #undef __FUNCT__
110 #define __FUNCT__ "SNESSolve_NRichardson"
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;
119 
120   PetscFunctionBegin;
121   snes->reason = SNES_CONVERGED_ITERATING;
122 
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 */
127 
128   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
129   snes->iter = 0;
130   snes->norm = 0.;
131   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
132   if (!snes->vec_func_init_set) {
133     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
134     if (snes->domainerror) {
135       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
136       PetscFunctionReturn(0);
137     }
138   } else {
139     snes->vec_func_init_set = PETSC_FALSE;
140   }
141   if (!snes->norm_init_set) {
142     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* 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   }
148 
149   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
150   snes->norm = fnorm;
151   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
152   SNESLogConvHistory(snes,fnorm,0);
153   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
154 
155   /* set parameter for default relative tolerance convergence test */
156   snes->ttol = fnorm*snes->rtol;
157   /* test convergence */
158   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
159   if (snes->reason) PetscFunctionReturn(0);
160 
161   for(i = 0; i < maxits; i++) {
162     lsSuccess = PETSC_TRUE;
163     /* Call general purpose update function */
164     if (snes->ops->update) {
165       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
166     }
167     else if (snes->pc) {
168       ierr = VecCopy(X,Y);CHKERRQ(ierr);
169       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
170       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
171       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
172         snes->reason = SNES_DIVERGED_INNER;
173         PetscFunctionReturn(0);
174       }
175       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
176     } else {
177       ierr = VecCopy(F,Y);CHKERRQ(ierr);
178     }
179     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
180     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
181     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr);
182     if (!lsSuccess) {
183       if (++snes->numFailures >= snes->maxFailures) {
184         snes->reason = SNES_DIVERGED_LINE_SEARCH;
185         break;
186       }
187     }
188     if (snes->nfuncs >= snes->max_funcs) {
189       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
190       break;
191     }
192     if (snes->domainerror) {
193       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
194       PetscFunctionReturn(0);
195     }
196     /* Monitor convergence */
197     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
198     snes->iter = i+1;
199     snes->norm = fnorm;
200     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
201     SNESLogConvHistory(snes,snes->norm,0);
202     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
203     /* Test for convergence */
204     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
205     if (snes->reason) break;
206   }
207   if (i == maxits) {
208     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
209     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
210   }
211   PetscFunctionReturn(0);
212 }
213 
214 /*MC
215   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
216 
217   Level: beginner
218 
219   Options Database:
220 .   -snes_ls <basic,basicnormnorms,quadratic,critical> Line search type.
221 
222   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
223             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
224             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner
225             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
226             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
227 
228      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
229 
230 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESQN, SNESNCG
231 M*/
232 EXTERN_C_BEGIN
233 #undef __FUNCT__
234 #define __FUNCT__ "SNESCreate_NRichardson"
235 PetscErrorCode  SNESCreate_NRichardson(SNES snes)
236 {
237   PetscErrorCode          ierr;
238   SNES_NRichardson        *neP;
239   PetscFunctionBegin;
240   snes->ops->destroy         = SNESDestroy_NRichardson;
241   snes->ops->setup           = SNESSetUp_NRichardson;
242   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
243   snes->ops->view            = SNESView_NRichardson;
244   snes->ops->solve           = SNESSolve_NRichardson;
245   snes->ops->reset           = SNESReset_NRichardson;
246 
247   snes->usesksp              = PETSC_FALSE;
248   snes->usespc               = PETSC_TRUE;
249 
250   ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr);
251   snes->data = (void*) neP;
252 
253   if (!snes->tolerancesset) {
254     snes->max_funcs = 30000;
255     snes->max_its   = 10000;
256   }
257 
258   PetscFunctionReturn(0);
259 }
260 EXTERN_C_END
261