xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision c6b63b32180a8dd433035330cddea33048fd2f9d)
142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith 
3bf7f4e0aSPeter Brune 
4028b6e76SBarry Smith #undef __FUNCT__
5d5c3842bSBarry Smith #define __FUNCT__ "SNESReset_NRichardson"
6d5c3842bSBarry Smith PetscErrorCode SNESReset_NRichardson(SNES snes)
7028b6e76SBarry Smith {
8028b6e76SBarry Smith   PetscFunctionBegin;
9028b6e76SBarry Smith   PetscFunctionReturn(0);
10028b6e76SBarry Smith }
11028b6e76SBarry Smith 
12028b6e76SBarry Smith /*
13d5c3842bSBarry Smith   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
14028b6e76SBarry Smith 
15028b6e76SBarry Smith   Input Parameter:
16028b6e76SBarry Smith . snes - the SNES context
17028b6e76SBarry Smith 
18028b6e76SBarry Smith   Application Interface Routine: SNESDestroy()
19028b6e76SBarry Smith */
20028b6e76SBarry Smith #undef __FUNCT__
21d5c3842bSBarry Smith #define __FUNCT__ "SNESDestroy_NRichardson"
22d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes)
23028b6e76SBarry Smith {
24028b6e76SBarry Smith   PetscErrorCode ierr;
25028b6e76SBarry Smith 
26028b6e76SBarry Smith   PetscFunctionBegin;
27d5c3842bSBarry Smith   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
28028b6e76SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
29028b6e76SBarry Smith   PetscFunctionReturn(0);
30028b6e76SBarry Smith }
31028b6e76SBarry Smith 
32028b6e76SBarry Smith /*
33d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
34d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
35028b6e76SBarry Smith 
36028b6e76SBarry Smith    Input Parameters:
37028b6e76SBarry Smith +  snes - the SNES context
38028b6e76SBarry Smith -  x - the solution vector
39028b6e76SBarry Smith 
40028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
41028b6e76SBarry Smith  */
42028b6e76SBarry Smith #undef __FUNCT__
43d5c3842bSBarry Smith #define __FUNCT__ "SNESSetUp_NRichardson"
44d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes)
45028b6e76SBarry Smith {
46028b6e76SBarry Smith   PetscFunctionBegin;
47*c6b63b32SPeter Brune   if (snes->pcside == PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");}
48028b6e76SBarry Smith   PetscFunctionReturn(0);
49028b6e76SBarry Smith }
50028b6e76SBarry Smith 
51028b6e76SBarry Smith /*
5204d7464bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.
53028b6e76SBarry Smith 
54028b6e76SBarry Smith   Input Parameter:
55028b6e76SBarry Smith . snes - the SNES context
56028b6e76SBarry Smith 
57028b6e76SBarry Smith   Application Interface Routine: SNESSetFromOptions()
58028b6e76SBarry Smith */
59028b6e76SBarry Smith #undef __FUNCT__
60d5c3842bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NRichardson"
61d5c3842bSBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
62028b6e76SBarry Smith {
63028b6e76SBarry Smith   PetscErrorCode ierr;
64f1c6b773SPeter Brune   SNESLineSearch linesearch;
655fd66863SKarl Rupp 
66028b6e76SBarry Smith   PetscFunctionBegin;
6742f4f86dSBarry Smith   ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr);
68028b6e76SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
699e764e56SPeter Brune   if (!snes->linesearch) {
707601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
711a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
729e764e56SPeter Brune   }
73028b6e76SBarry Smith   PetscFunctionReturn(0);
74028b6e76SBarry Smith }
75028b6e76SBarry Smith 
76028b6e76SBarry Smith /*
77d5c3842bSBarry Smith   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
78028b6e76SBarry Smith 
79028b6e76SBarry Smith   Input Parameters:
80028b6e76SBarry Smith + SNES - the SNES context
81028b6e76SBarry Smith - viewer - visualization context
82028b6e76SBarry Smith 
83028b6e76SBarry Smith   Application Interface Routine: SNESView()
84028b6e76SBarry Smith */
85028b6e76SBarry Smith #undef __FUNCT__
86d5c3842bSBarry Smith #define __FUNCT__ "SNESView_NRichardson"
87d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
88028b6e76SBarry Smith {
89028b6e76SBarry Smith   PetscBool      iascii;
90028b6e76SBarry Smith   PetscErrorCode ierr;
91028b6e76SBarry Smith 
92028b6e76SBarry Smith   PetscFunctionBegin;
93251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
94028b6e76SBarry Smith   if (iascii) {
95028b6e76SBarry Smith   }
96028b6e76SBarry Smith   PetscFunctionReturn(0);
97028b6e76SBarry Smith }
98028b6e76SBarry Smith 
99028b6e76SBarry Smith /*
100d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
101028b6e76SBarry Smith 
102028b6e76SBarry Smith   Input Parameters:
103028b6e76SBarry Smith . snes - the SNES context
104028b6e76SBarry Smith 
105028b6e76SBarry Smith   Output Parameter:
106028b6e76SBarry Smith . outits - number of iterations until termination
107028b6e76SBarry Smith 
108028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
109028b6e76SBarry Smith */
110028b6e76SBarry Smith #undef __FUNCT__
111d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson"
112d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes)
113028b6e76SBarry Smith {
114bf7f4e0aSPeter Brune   Vec                 X, Y, F;
115e7058c64SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
116028b6e76SBarry Smith   PetscInt            maxits, i;
117028b6e76SBarry Smith   PetscErrorCode      ierr;
118bf7f4e0aSPeter Brune   PetscBool           lsSuccess;
119028b6e76SBarry Smith   SNESConvergedReason reason;
120028b6e76SBarry Smith 
121028b6e76SBarry Smith   PetscFunctionBegin;
122028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
123028b6e76SBarry Smith 
124028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
125028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
126028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
127028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
128028b6e76SBarry Smith 
129ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
130028b6e76SBarry Smith   snes->iter = 0;
131028b6e76SBarry Smith   snes->norm = 0.;
132ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
133e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
134028b6e76SBarry Smith     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
135028b6e76SBarry Smith     if (snes->domainerror) {
136028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
137028b6e76SBarry Smith       PetscFunctionReturn(0);
138028b6e76SBarry Smith     }
1391aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
1401aa26658SKarl Rupp 
141e4ed7901SPeter Brune   if (!snes->norm_init_set) {
142028b6e76SBarry Smith     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
143189a9710SBarry Smith     if (PetscIsInfOrNanReal(fnorm)) {
144189a9710SBarry Smith       snes->reason = SNES_DIVERGED_FNORM_NAN;
145189a9710SBarry Smith       PetscFunctionReturn(0);
146189a9710SBarry Smith     }
147e4ed7901SPeter Brune   } else {
148e4ed7901SPeter Brune     fnorm               = snes->norm_init;
149e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
150e4ed7901SPeter Brune   }
151e4ed7901SPeter Brune 
152ce8c27fbSBarry Smith   ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
153028b6e76SBarry Smith   snes->norm = fnorm;
154ce8c27fbSBarry Smith   ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
155a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
156028b6e76SBarry Smith   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
157028b6e76SBarry Smith 
158028b6e76SBarry Smith   /* set parameter for default relative tolerance convergence test */
159028b6e76SBarry Smith   snes->ttol = fnorm*snes->rtol;
160028b6e76SBarry Smith   /* test convergence */
161028b6e76SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
162028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
163028b6e76SBarry Smith 
164028b6e76SBarry Smith   for (i = 0; i < maxits; i++) {
165bf7f4e0aSPeter Brune     lsSuccess = PETSC_TRUE;
166028b6e76SBarry Smith     /* Call general purpose update function */
167028b6e76SBarry Smith     if (snes->ops->update) {
168028b6e76SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
169028b6e76SBarry Smith     }
170*c6b63b32SPeter Brune     if (snes->pc) {
171028b6e76SBarry Smith       ierr = VecCopy(X,Y);CHKERRQ(ierr);
172*c6b63b32SPeter Brune       if (snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
173217b9c2eSPeter Brune         ierr = SNESSetInitialFunction(snes->pc, F);CHKERRQ(ierr);
174217b9c2eSPeter Brune         ierr = SNESSetInitialFunctionNorm(snes->pc, fnorm);CHKERRQ(ierr);
175*c6b63b32SPeter Brune       }
17690cf3547SPeter Brune       ierr = PetscLogEventBegin(SNES_NPCSolve,snes->pc,Y,snes->vec_rhs,0);CHKERRQ(ierr);
177c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
17890cf3547SPeter Brune       ierr = PetscLogEventEnd(SNES_NPCSolve,snes->pc,Y,snes->vec_rhs,0);CHKERRQ(ierr);
179028b6e76SBarry Smith       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
18034b225ddSBarry Smith       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
181028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_INNER;
182028b6e76SBarry Smith         PetscFunctionReturn(0);
183028b6e76SBarry Smith       }
184e5583369SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
185a3685310SPeter Brune     } else {
186a3685310SPeter Brune       ierr = VecCopy(F,Y);CHKERRQ(ierr);
187028b6e76SBarry Smith     }
188f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
189f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
190f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr);
191028b6e76SBarry Smith     if (!lsSuccess) {
192028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
193028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
194028b6e76SBarry Smith         break;
195028b6e76SBarry Smith       }
196028b6e76SBarry Smith     }
197028b6e76SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
198028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
199028b6e76SBarry Smith       break;
200028b6e76SBarry Smith     }
201028b6e76SBarry Smith     if (snes->domainerror) {
202028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
203028b6e76SBarry Smith       PetscFunctionReturn(0);
204028b6e76SBarry Smith     }
205028b6e76SBarry Smith     /* Monitor convergence */
206ce8c27fbSBarry Smith     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
207028b6e76SBarry Smith     snes->iter = i+1;
208028b6e76SBarry Smith     snes->norm = fnorm;
209ce8c27fbSBarry Smith     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
210a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
211028b6e76SBarry Smith     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
212028b6e76SBarry Smith     /* Test for convergence */
213e7058c64SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
214028b6e76SBarry Smith     if (snes->reason) break;
215028b6e76SBarry Smith   }
216028b6e76SBarry Smith   if (i == maxits) {
217028b6e76SBarry Smith     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
218028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
219028b6e76SBarry Smith   }
220028b6e76SBarry Smith   PetscFunctionReturn(0);
221028b6e76SBarry Smith }
222028b6e76SBarry Smith 
223028b6e76SBarry Smith /*MC
224d5c3842bSBarry Smith   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
225028b6e76SBarry Smith 
226028b6e76SBarry Smith   Level: beginner
227028b6e76SBarry Smith 
228028b6e76SBarry Smith   Options Database:
22972835e02SPeter Brune +   -snes_linesearch_type <l2,cp,basic> Line search type.
23072835e02SPeter Brune -   -snes_linesearch_damping<1.0> Damping for the line search.
231028b6e76SBarry Smith 
232af60355fSPeter Brune   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
233af60355fSPeter Brune             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
234af60355fSPeter Brune             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner
235af60355fSPeter Brune             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
236af60355fSPeter Brune             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
237028b6e76SBarry Smith 
23872835e02SPeter Brune             The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
23972835e02SPeter Brune             linesearch, one may have to scale the update with -snes_linesearch_damping
24072835e02SPeter Brune 
241028b6e76SBarry Smith      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
242028b6e76SBarry Smith 
24304d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
244028b6e76SBarry Smith M*/
245028b6e76SBarry Smith #undef __FUNCT__
246d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson"
2478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
248028b6e76SBarry Smith {
24992c02d66SPeter Brune   PetscErrorCode   ierr;
250bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
2515fd66863SKarl Rupp 
252028b6e76SBarry Smith   PetscFunctionBegin;
253d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
254d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
255d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
256d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
257d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
258d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
259028b6e76SBarry Smith 
260028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
26142f4f86dSBarry Smith   snes->usespc  = PETSC_TRUE;
262028b6e76SBarry Smith 
263*c6b63b32SPeter Brune   snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
264*c6b63b32SPeter Brune   snes->pcside = PC_LEFT;
265*c6b63b32SPeter Brune 
266bf7f4e0aSPeter Brune   ierr       = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr);
267bf7f4e0aSPeter Brune   snes->data = (void*) neP;
268bf7f4e0aSPeter Brune 
26988976e71SPeter Brune   if (!snes->tolerancesset) {
2700e444f03SPeter Brune     snes->max_funcs = 30000;
2710e444f03SPeter Brune     snes->max_its   = 10000;
272c522fa08SPeter Brune     snes->stol      = 1e-20;
27388976e71SPeter Brune   }
274028b6e76SBarry Smith   PetscFunctionReturn(0);
275028b6e76SBarry Smith }
276