xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 88976e7148bb891b38778b9a5d5ae4c2b2cf44bf)
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 
9028b6e76SBarry Smith   PetscFunctionBegin;
10028b6e76SBarry Smith   PetscFunctionReturn(0);
11028b6e76SBarry Smith }
12028b6e76SBarry Smith 
13028b6e76SBarry Smith /*
14d5c3842bSBarry Smith   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
15028b6e76SBarry Smith 
16028b6e76SBarry Smith   Input Parameter:
17028b6e76SBarry Smith . snes - the SNES context
18028b6e76SBarry Smith 
19028b6e76SBarry Smith   Application Interface Routine: SNESDestroy()
20028b6e76SBarry Smith */
21028b6e76SBarry Smith #undef __FUNCT__
22d5c3842bSBarry Smith #define __FUNCT__ "SNESDestroy_NRichardson"
23d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes)
24028b6e76SBarry Smith {
25028b6e76SBarry Smith   PetscErrorCode   ierr;
26028b6e76SBarry Smith 
27028b6e76SBarry Smith   PetscFunctionBegin;
28d5c3842bSBarry Smith   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
29028b6e76SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
30028b6e76SBarry Smith   PetscFunctionReturn(0);
31028b6e76SBarry Smith }
32028b6e76SBarry Smith 
33028b6e76SBarry Smith /*
34d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
35d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
36028b6e76SBarry Smith 
37028b6e76SBarry Smith    Input Parameters:
38028b6e76SBarry Smith +  snes - the SNES context
39028b6e76SBarry Smith -  x - the solution vector
40028b6e76SBarry Smith 
41028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
42028b6e76SBarry Smith  */
43028b6e76SBarry Smith #undef __FUNCT__
44d5c3842bSBarry Smith #define __FUNCT__ "SNESSetUp_NRichardson"
45d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes)
46028b6e76SBarry Smith {
47028b6e76SBarry Smith   PetscFunctionBegin;
48028b6e76SBarry Smith   PetscFunctionReturn(0);
49028b6e76SBarry Smith }
50028b6e76SBarry Smith 
51028b6e76SBarry Smith /*
52d5c3842bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS 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;
65028b6e76SBarry Smith   PetscFunctionBegin;
6642f4f86dSBarry Smith   ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr);
67028b6e76SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
689e764e56SPeter Brune   if (!snes->linesearch) {
69f1c6b773SPeter Brune     ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
70f1c6b773SPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
719e764e56SPeter Brune   }
72028b6e76SBarry Smith   PetscFunctionReturn(0);
73028b6e76SBarry Smith }
74028b6e76SBarry Smith 
75028b6e76SBarry Smith /*
76d5c3842bSBarry Smith   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
77028b6e76SBarry Smith 
78028b6e76SBarry Smith   Input Parameters:
79028b6e76SBarry Smith + SNES - the SNES context
80028b6e76SBarry Smith - viewer - visualization context
81028b6e76SBarry Smith 
82028b6e76SBarry Smith   Application Interface Routine: SNESView()
83028b6e76SBarry Smith */
84028b6e76SBarry Smith #undef __FUNCT__
85d5c3842bSBarry Smith #define __FUNCT__ "SNESView_NRichardson"
86d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
87028b6e76SBarry Smith {
88028b6e76SBarry Smith   PetscBool        iascii;
89028b6e76SBarry Smith   PetscErrorCode   ierr;
90028b6e76SBarry Smith 
91028b6e76SBarry Smith   PetscFunctionBegin;
92028b6e76SBarry Smith   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
93028b6e76SBarry Smith   if (iascii) {
94028b6e76SBarry Smith   }
95028b6e76SBarry Smith   PetscFunctionReturn(0);
96028b6e76SBarry Smith }
97028b6e76SBarry Smith 
98028b6e76SBarry Smith /*
99d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
100028b6e76SBarry Smith 
101028b6e76SBarry Smith   Input Parameters:
102028b6e76SBarry Smith . snes - the SNES context
103028b6e76SBarry Smith 
104028b6e76SBarry Smith   Output Parameter:
105028b6e76SBarry Smith . outits - number of iterations until termination
106028b6e76SBarry Smith 
107028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
108028b6e76SBarry Smith */
109028b6e76SBarry Smith #undef __FUNCT__
110d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson"
111d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes)
112028b6e76SBarry Smith {
113bf7f4e0aSPeter Brune   Vec                 X, Y, F;
114e7058c64SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
115028b6e76SBarry Smith   PetscInt            maxits, i;
116028b6e76SBarry Smith   PetscErrorCode      ierr;
117bf7f4e0aSPeter Brune   PetscBool           lsSuccess;
118028b6e76SBarry Smith   SNESConvergedReason reason;
119028b6e76SBarry Smith 
120028b6e76SBarry Smith   PetscFunctionBegin;
121028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
122028b6e76SBarry Smith 
123028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
124028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
125028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
126028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
127028b6e76SBarry Smith 
128028b6e76SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
129028b6e76SBarry Smith   snes->iter = 0;
130028b6e76SBarry Smith   snes->norm = 0.;
131028b6e76SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
132028b6e76SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
133028b6e76SBarry Smith   if (snes->domainerror) {
134028b6e76SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
135028b6e76SBarry Smith     PetscFunctionReturn(0);
136028b6e76SBarry Smith   }
137028b6e76SBarry Smith   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
138028b6e76SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
139028b6e76SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
140028b6e76SBarry Smith   snes->norm = fnorm;
141028b6e76SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
142028b6e76SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
143028b6e76SBarry Smith   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
144028b6e76SBarry Smith 
145028b6e76SBarry Smith   /* set parameter for default relative tolerance convergence test */
146028b6e76SBarry Smith   snes->ttol = fnorm*snes->rtol;
147028b6e76SBarry Smith   /* test convergence */
148028b6e76SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
149028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
150028b6e76SBarry Smith 
151028b6e76SBarry Smith   for(i = 0; i < maxits; i++) {
152bf7f4e0aSPeter Brune     lsSuccess = PETSC_TRUE;
153028b6e76SBarry Smith     /* Call general purpose update function */
154028b6e76SBarry Smith     if (snes->ops->update) {
155028b6e76SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
156028b6e76SBarry Smith     }
15783947a81SPeter Brune     else if (snes->pc) {
158028b6e76SBarry Smith       ierr = VecCopy(X,Y);CHKERRQ(ierr);
159c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
160028b6e76SBarry Smith       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
16134b225ddSBarry Smith       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
162028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_INNER;
163028b6e76SBarry Smith         PetscFunctionReturn(0);
164028b6e76SBarry Smith       }
165e5583369SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
166a3685310SPeter Brune     } else {
167a3685310SPeter Brune       ierr = VecCopy(F,Y);CHKERRQ(ierr);
168028b6e76SBarry Smith     }
169f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
170f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
171f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lsSuccess);CHKERRQ(ierr);
172028b6e76SBarry Smith     if (!lsSuccess) {
173028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
174028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
175028b6e76SBarry Smith         break;
176028b6e76SBarry Smith       }
177028b6e76SBarry Smith     }
178028b6e76SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
179028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
180028b6e76SBarry Smith       break;
181028b6e76SBarry Smith     }
182028b6e76SBarry Smith     if (snes->domainerror) {
183028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
184028b6e76SBarry Smith       PetscFunctionReturn(0);
185028b6e76SBarry Smith     }
186028b6e76SBarry Smith     /* Monitor convergence */
187028b6e76SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
188028b6e76SBarry Smith     snes->iter = i+1;
189028b6e76SBarry Smith     snes->norm = fnorm;
190028b6e76SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
191028b6e76SBarry Smith     SNESLogConvHistory(snes,snes->norm,0);
192028b6e76SBarry Smith     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
193028b6e76SBarry Smith     /* Test for convergence */
194e7058c64SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
195028b6e76SBarry Smith     if (snes->reason) break;
196028b6e76SBarry Smith   }
197028b6e76SBarry Smith   if (i == maxits) {
198028b6e76SBarry Smith     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
199028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
200028b6e76SBarry Smith   }
201028b6e76SBarry Smith   PetscFunctionReturn(0);
202028b6e76SBarry Smith }
203028b6e76SBarry Smith 
204028b6e76SBarry Smith /*MC
205d5c3842bSBarry Smith   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
206028b6e76SBarry Smith 
207028b6e76SBarry Smith   Level: beginner
208028b6e76SBarry Smith 
209028b6e76SBarry Smith   Options Database:
210f3aaf736SPeter Brune .   -snes_ls <basic,basicnormnorms,quadratic,critical> Line search type.
211028b6e76SBarry Smith 
212af60355fSPeter Brune   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
213af60355fSPeter Brune             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
214af60355fSPeter Brune             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner
215af60355fSPeter Brune             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
216af60355fSPeter Brune             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
217028b6e76SBarry Smith 
218028b6e76SBarry Smith      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
219028b6e76SBarry Smith 
2201867fe5bSPeter Brune .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESQN, SNESNCG
221028b6e76SBarry Smith M*/
222028b6e76SBarry Smith EXTERN_C_BEGIN
223028b6e76SBarry Smith #undef __FUNCT__
224d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson"
225d5c3842bSBarry Smith PetscErrorCode  SNESCreate_NRichardson(SNES snes)
226028b6e76SBarry Smith {
22792c02d66SPeter Brune   PetscErrorCode          ierr;
228bf7f4e0aSPeter Brune   SNES_NRichardson        *neP;
229028b6e76SBarry Smith   PetscFunctionBegin;
230d5c3842bSBarry Smith   snes->ops->destroy         = SNESDestroy_NRichardson;
231d5c3842bSBarry Smith   snes->ops->setup           = SNESSetUp_NRichardson;
232d5c3842bSBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
233d5c3842bSBarry Smith   snes->ops->view            = SNESView_NRichardson;
234d5c3842bSBarry Smith   snes->ops->solve           = SNESSolve_NRichardson;
235d5c3842bSBarry Smith   snes->ops->reset           = SNESReset_NRichardson;
236028b6e76SBarry Smith 
237028b6e76SBarry Smith   snes->usesksp              = PETSC_FALSE;
23842f4f86dSBarry Smith   snes->usespc               = PETSC_TRUE;
239028b6e76SBarry Smith 
240bf7f4e0aSPeter Brune   ierr = PetscNewLog(snes, SNES_NRichardson, &neP);CHKERRQ(ierr);
241bf7f4e0aSPeter Brune   snes->data = (void*) neP;
242bf7f4e0aSPeter Brune 
243*88976e71SPeter Brune   if (!snes->tolerancesset) {
2440e444f03SPeter Brune     snes->max_funcs = 30000;
2450e444f03SPeter Brune     snes->max_its   = 10000;
246*88976e71SPeter Brune   }
2470e444f03SPeter Brune 
248028b6e76SBarry Smith   PetscFunctionReturn(0);
249028b6e76SBarry Smith }
250028b6e76SBarry Smith EXTERN_C_END
251