xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 98921bda46e76d7aaed9e0138c5ff9d0ce93f355)
142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith 
3d5c3842bSBarry Smith PetscErrorCode SNESReset_NRichardson(SNES snes)
4028b6e76SBarry Smith {
5028b6e76SBarry Smith   PetscFunctionBegin;
6028b6e76SBarry Smith   PetscFunctionReturn(0);
7028b6e76SBarry Smith }
8028b6e76SBarry Smith 
9028b6e76SBarry Smith /*
10d5c3842bSBarry Smith   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
11028b6e76SBarry Smith 
12028b6e76SBarry Smith   Input Parameter:
13028b6e76SBarry Smith . snes - the SNES context
14028b6e76SBarry Smith 
15028b6e76SBarry Smith   Application Interface Routine: SNESDestroy()
16028b6e76SBarry Smith */
17d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes)
18028b6e76SBarry Smith {
19028b6e76SBarry Smith   PetscErrorCode ierr;
20028b6e76SBarry Smith 
21028b6e76SBarry Smith   PetscFunctionBegin;
22d5c3842bSBarry Smith   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
23028b6e76SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
24028b6e76SBarry Smith   PetscFunctionReturn(0);
25028b6e76SBarry Smith }
26028b6e76SBarry Smith 
27028b6e76SBarry Smith /*
28d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
29d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
30028b6e76SBarry Smith 
31028b6e76SBarry Smith    Input Parameters:
32028b6e76SBarry Smith +  snes - the SNES context
33028b6e76SBarry Smith -  x - the solution vector
34028b6e76SBarry Smith 
35028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
36028b6e76SBarry Smith  */
37d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes)
38028b6e76SBarry Smith {
39028b6e76SBarry Smith   PetscFunctionBegin;
402da392ccSBarry Smith   if (snes->npcside== PC_RIGHT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");
416c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
42028b6e76SBarry Smith   PetscFunctionReturn(0);
43028b6e76SBarry Smith }
44028b6e76SBarry Smith 
45028b6e76SBarry Smith /*
4604d7464bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.
47028b6e76SBarry Smith 
48028b6e76SBarry Smith   Input Parameter:
49028b6e76SBarry Smith . snes - the SNES context
50028b6e76SBarry Smith 
51028b6e76SBarry Smith   Application Interface Routine: SNESSetFromOptions()
52028b6e76SBarry Smith */
534416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(PetscOptionItems *PetscOptionsObject,SNES snes)
54028b6e76SBarry Smith {
55028b6e76SBarry Smith   PetscErrorCode ierr;
565fd66863SKarl Rupp 
57028b6e76SBarry Smith   PetscFunctionBegin;
58e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");CHKERRQ(ierr);
59028b6e76SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
60028b6e76SBarry Smith   PetscFunctionReturn(0);
61028b6e76SBarry Smith }
62028b6e76SBarry Smith 
63028b6e76SBarry Smith /*
64d5c3842bSBarry Smith   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
65028b6e76SBarry Smith 
66028b6e76SBarry Smith   Input Parameters:
67028b6e76SBarry Smith + SNES - the SNES context
68028b6e76SBarry Smith - viewer - visualization context
69028b6e76SBarry Smith 
70028b6e76SBarry Smith   Application Interface Routine: SNESView()
71028b6e76SBarry Smith */
72d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
73028b6e76SBarry Smith {
74028b6e76SBarry Smith   PetscBool      iascii;
75028b6e76SBarry Smith   PetscErrorCode ierr;
76028b6e76SBarry Smith 
77028b6e76SBarry Smith   PetscFunctionBegin;
78251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
79028b6e76SBarry Smith   if (iascii) {
80028b6e76SBarry Smith   }
81028b6e76SBarry Smith   PetscFunctionReturn(0);
82028b6e76SBarry Smith }
83028b6e76SBarry Smith 
84028b6e76SBarry Smith /*
85d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
86028b6e76SBarry Smith 
87028b6e76SBarry Smith   Input Parameters:
88028b6e76SBarry Smith . snes - the SNES context
89028b6e76SBarry Smith 
90028b6e76SBarry Smith   Output Parameter:
91028b6e76SBarry Smith . outits - number of iterations until termination
92028b6e76SBarry Smith 
93028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
94028b6e76SBarry Smith */
95d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes)
96028b6e76SBarry Smith {
97bf7f4e0aSPeter Brune   Vec                  X, Y, F;
98e7058c64SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
99028b6e76SBarry Smith   PetscInt             maxits, i;
100028b6e76SBarry Smith   PetscErrorCode       ierr;
101422a814eSBarry Smith   SNESLineSearchReason lsresult;
102028b6e76SBarry Smith   SNESConvergedReason  reason;
103028b6e76SBarry Smith 
104028b6e76SBarry Smith   PetscFunctionBegin;
105*98921bdaSJacob Faibussowitsch   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
106c579b300SPatrick Farrell 
107028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
108028b6e76SBarry Smith 
109028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
110028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
111028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
112028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
113028b6e76SBarry Smith 
114e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
115028b6e76SBarry Smith   snes->iter = 0;
116028b6e76SBarry Smith   snes->norm = 0.;
117e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
118b7281c8aSPeter Brune 
119efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
120be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
121efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
122b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
123b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
124b7281c8aSPeter Brune       PetscFunctionReturn(0);
125b7281c8aSPeter Brune     }
126b7281c8aSPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
127b7281c8aSPeter Brune   } else {
128e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
129028b6e76SBarry Smith       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1301aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
131c1c75074SPeter Brune 
132b7281c8aSPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
133422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
134b7281c8aSPeter Brune   }
135efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
136be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
137efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
138b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
139b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
140b7281c8aSPeter Brune         PetscFunctionReturn(0);
141b7281c8aSPeter Brune       }
142b7281c8aSPeter Brune   } else {
143b7281c8aSPeter Brune     ierr = VecCopy(F,Y);CHKERRQ(ierr);
144b7281c8aSPeter Brune   }
145e4ed7901SPeter Brune 
146e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
147028b6e76SBarry Smith   snes->norm = fnorm;
148e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
149a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
150028b6e76SBarry Smith   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
151028b6e76SBarry Smith 
152028b6e76SBarry Smith   /* test convergence */
153028b6e76SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
154028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
155028b6e76SBarry Smith 
156028b6e76SBarry Smith   /* Call general purpose update function */
157028b6e76SBarry Smith   if (snes->ops->update) {
158028b6e76SBarry Smith     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
159028b6e76SBarry Smith   }
160b7281c8aSPeter Brune 
161b7281c8aSPeter Brune   /* set parameter for default relative tolerance convergence test */
162b7281c8aSPeter Brune   snes->ttol = fnorm*snes->rtol;
163b7281c8aSPeter Brune   /* test convergence */
164b7281c8aSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
165b7281c8aSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
166b7281c8aSPeter Brune 
167b7281c8aSPeter Brune   for (i = 1; i < maxits+1; i++) {
168f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
169422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
170f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
171422a814eSBarry Smith     if (lsresult) {
172028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
173028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
174028b6e76SBarry Smith         break;
175028b6e76SBarry Smith       }
176028b6e76SBarry Smith     }
177e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
178028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
179028b6e76SBarry Smith       break;
180028b6e76SBarry Smith     }
181b7281c8aSPeter Brune 
182028b6e76SBarry Smith     /* Monitor convergence */
183e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
184b7281c8aSPeter Brune     snes->iter = i;
185028b6e76SBarry Smith     snes->norm = fnorm;
186c1e67a49SFande Kong     snes->xnorm = xnorm;
187c1e67a49SFande Kong     snes->ynorm = ynorm;
188e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
189a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
190028b6e76SBarry Smith     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
191028b6e76SBarry Smith     /* Test for convergence */
192e7058c64SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
193028b6e76SBarry Smith     if (snes->reason) break;
194b7281c8aSPeter Brune 
195b7281c8aSPeter Brune     /* Call general purpose update function */
196b7281c8aSPeter Brune     if (snes->ops->update) {
197b7281c8aSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
198028b6e76SBarry Smith     }
199b7281c8aSPeter Brune 
200efd4aadfSBarry Smith     if (snes->npc) {
201b7281c8aSPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
202be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,NULL,Y);CHKERRQ(ierr);
203b7281c8aSPeter Brune         ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
204b7281c8aSPeter Brune         ierr = VecCopy(Y,F);CHKERRQ(ierr);
205b7281c8aSPeter Brune       } else {
206be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
207b7281c8aSPeter Brune       }
208efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
209b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
210b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
211b7281c8aSPeter Brune         PetscFunctionReturn(0);
212b7281c8aSPeter Brune       }
213b7281c8aSPeter Brune     } else {
214b7281c8aSPeter Brune       ierr = VecCopy(F,Y);CHKERRQ(ierr);
215b7281c8aSPeter Brune     }
216b7281c8aSPeter Brune   }
217b7281c8aSPeter Brune   if (i == maxits+1) {
218028b6e76SBarry Smith     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
219028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
220028b6e76SBarry Smith   }
221028b6e76SBarry Smith   PetscFunctionReturn(0);
222028b6e76SBarry Smith }
223028b6e76SBarry Smith 
224028b6e76SBarry Smith /*MC
225d5c3842bSBarry Smith   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
226028b6e76SBarry Smith 
227028b6e76SBarry Smith   Level: beginner
228028b6e76SBarry Smith 
229028b6e76SBarry Smith   Options Database:
23072835e02SPeter Brune +   -snes_linesearch_type <l2,cp,basic> Line search type.
23172835e02SPeter Brune -   -snes_linesearch_damping<1.0> Damping for the line search.
232028b6e76SBarry Smith 
23395452b02SPatrick Sanan   Notes:
23495452b02SPatrick Sanan     If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
235af60355fSPeter Brune             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
236be95d8f1SBarry Smith             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
237af60355fSPeter Brune             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
238af60355fSPeter Brune             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
239028b6e76SBarry Smith 
24072835e02SPeter Brune             The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
24172835e02SPeter Brune             linesearch, one may have to scale the update with -snes_linesearch_damping
24272835e02SPeter Brune 
243028b6e76SBarry Smith      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
244028b6e76SBarry Smith 
2452d547940SBarry Smith      Only supports left non-linear preconditioning.
2462d547940SBarry Smith 
24704d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
248028b6e76SBarry Smith M*/
2498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
250028b6e76SBarry Smith {
25192c02d66SPeter Brune   PetscErrorCode   ierr;
252bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
253d8d34be6SBarry Smith   SNESLineSearch   linesearch;
2545fd66863SKarl Rupp 
255028b6e76SBarry Smith   PetscFunctionBegin;
256d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
257d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
258d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
259d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
260d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
261d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
262028b6e76SBarry Smith 
263028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
264efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
265028b6e76SBarry Smith 
266efd4aadfSBarry Smith   snes->npcside= PC_LEFT;
267c6b63b32SPeter Brune 
268d8d34be6SBarry Smith   ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
269ec786807SJed Brown   if (!((PetscObject)linesearch)->type_name) {
270d8d34be6SBarry Smith     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
271ec786807SJed Brown   }
272d8d34be6SBarry Smith 
2734fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2744fc747eaSLawrence Mitchell 
275b00a9115SJed Brown   ierr       = PetscNewLog(snes,&neP);CHKERRQ(ierr);
276bf7f4e0aSPeter Brune   snes->data = (void*) neP;
277bf7f4e0aSPeter Brune 
27888976e71SPeter Brune   if (!snes->tolerancesset) {
2790e444f03SPeter Brune     snes->max_funcs = 30000;
2800e444f03SPeter Brune     snes->max_its   = 10000;
281c522fa08SPeter Brune     snes->stol      = 1e-20;
28288976e71SPeter Brune   }
283028b6e76SBarry Smith   PetscFunctionReturn(0);
284028b6e76SBarry Smith }
285