xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision efd4aadf157bf1ba2d80c2be092fcf4247860003)
142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith 
3bf7f4e0aSPeter Brune 
4d5c3842bSBarry Smith PetscErrorCode SNESReset_NRichardson(SNES snes)
5028b6e76SBarry Smith {
6028b6e76SBarry Smith   PetscFunctionBegin;
7028b6e76SBarry Smith   PetscFunctionReturn(0);
8028b6e76SBarry Smith }
9028b6e76SBarry Smith 
10028b6e76SBarry Smith /*
11d5c3842bSBarry Smith   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
12028b6e76SBarry Smith 
13028b6e76SBarry Smith   Input Parameter:
14028b6e76SBarry Smith . snes - the SNES context
15028b6e76SBarry Smith 
16028b6e76SBarry Smith   Application Interface Routine: SNESDestroy()
17028b6e76SBarry Smith */
18d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes)
19028b6e76SBarry Smith {
20028b6e76SBarry Smith   PetscErrorCode ierr;
21028b6e76SBarry Smith 
22028b6e76SBarry Smith   PetscFunctionBegin;
23d5c3842bSBarry Smith   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
24028b6e76SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
25028b6e76SBarry Smith   PetscFunctionReturn(0);
26028b6e76SBarry Smith }
27028b6e76SBarry Smith 
28028b6e76SBarry Smith /*
29d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
30d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
31028b6e76SBarry Smith 
32028b6e76SBarry Smith    Input Parameters:
33028b6e76SBarry Smith +  snes - the SNES context
34028b6e76SBarry Smith -  x - the solution vector
35028b6e76SBarry Smith 
36028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
37028b6e76SBarry Smith  */
38d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes)
39028b6e76SBarry Smith {
40028b6e76SBarry Smith   PetscFunctionBegin;
41*efd4aadfSBarry Smith   if (snes->npcside== PC_RIGHT) {SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"NRichardson only supports left preconditioning");}
426c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
43028b6e76SBarry Smith   PetscFunctionReturn(0);
44028b6e76SBarry Smith }
45028b6e76SBarry Smith 
46028b6e76SBarry Smith /*
4704d7464bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.
48028b6e76SBarry Smith 
49028b6e76SBarry Smith   Input Parameter:
50028b6e76SBarry Smith . snes - the SNES context
51028b6e76SBarry Smith 
52028b6e76SBarry Smith   Application Interface Routine: SNESSetFromOptions()
53028b6e76SBarry Smith */
544416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(PetscOptionItems *PetscOptionsObject,SNES snes)
55028b6e76SBarry Smith {
56028b6e76SBarry Smith   PetscErrorCode ierr;
57f1c6b773SPeter Brune   SNESLineSearch linesearch;
585fd66863SKarl Rupp 
59028b6e76SBarry Smith   PetscFunctionBegin;
60e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNES Richardson options");CHKERRQ(ierr);
61028b6e76SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
629e764e56SPeter Brune   if (!snes->linesearch) {
637601faf0SJed Brown     ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
641a4f838cSPeter Brune     ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
659e764e56SPeter Brune   }
66028b6e76SBarry Smith   PetscFunctionReturn(0);
67028b6e76SBarry Smith }
68028b6e76SBarry Smith 
69028b6e76SBarry Smith /*
70d5c3842bSBarry Smith   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
71028b6e76SBarry Smith 
72028b6e76SBarry Smith   Input Parameters:
73028b6e76SBarry Smith + SNES - the SNES context
74028b6e76SBarry Smith - viewer - visualization context
75028b6e76SBarry Smith 
76028b6e76SBarry Smith   Application Interface Routine: SNESView()
77028b6e76SBarry Smith */
78d5c3842bSBarry Smith static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
79028b6e76SBarry Smith {
80028b6e76SBarry Smith   PetscBool      iascii;
81028b6e76SBarry Smith   PetscErrorCode ierr;
82028b6e76SBarry Smith 
83028b6e76SBarry Smith   PetscFunctionBegin;
84251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
85028b6e76SBarry Smith   if (iascii) {
86028b6e76SBarry Smith   }
87028b6e76SBarry Smith   PetscFunctionReturn(0);
88028b6e76SBarry Smith }
89028b6e76SBarry Smith 
90028b6e76SBarry Smith /*
91d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
92028b6e76SBarry Smith 
93028b6e76SBarry Smith   Input Parameters:
94028b6e76SBarry Smith . snes - the SNES context
95028b6e76SBarry Smith 
96028b6e76SBarry Smith   Output Parameter:
97028b6e76SBarry Smith . outits - number of iterations until termination
98028b6e76SBarry Smith 
99028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
100028b6e76SBarry Smith */
101d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes)
102028b6e76SBarry Smith {
103bf7f4e0aSPeter Brune   Vec                  X, Y, F;
104e7058c64SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
105028b6e76SBarry Smith   PetscInt             maxits, i;
106028b6e76SBarry Smith   PetscErrorCode       ierr;
107422a814eSBarry Smith   SNESLineSearchReason lsresult;
108028b6e76SBarry Smith   SNESConvergedReason  reason;
109028b6e76SBarry Smith 
110028b6e76SBarry Smith   PetscFunctionBegin;
1116c4ed002SBarry Smith   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);
112c579b300SPatrick Farrell 
113028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
114028b6e76SBarry Smith 
115028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
116028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
117028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
118028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
119028b6e76SBarry Smith 
120e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
121028b6e76SBarry Smith   snes->iter = 0;
122028b6e76SBarry Smith   snes->norm = 0.;
123e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
124b7281c8aSPeter Brune 
125*efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
126be95d8f1SBarry Smith     ierr = SNESApplyNPC(snes,X,NULL,F);CHKERRQ(ierr);
127*efd4aadfSBarry Smith     ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
128b7281c8aSPeter Brune     if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
129b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
130b7281c8aSPeter Brune       PetscFunctionReturn(0);
131b7281c8aSPeter Brune     }
132b7281c8aSPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
133b7281c8aSPeter Brune   } else {
134e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
135028b6e76SBarry Smith       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1361aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
137c1c75074SPeter Brune 
138b7281c8aSPeter Brune     ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
139422a814eSBarry Smith     SNESCheckFunctionNorm(snes,fnorm);
140b7281c8aSPeter Brune   }
141*efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
142be95d8f1SBarry Smith       ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
143*efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
144b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
145b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
146b7281c8aSPeter Brune         PetscFunctionReturn(0);
147b7281c8aSPeter Brune       }
148b7281c8aSPeter Brune   } else {
149b7281c8aSPeter Brune     ierr = VecCopy(F,Y);CHKERRQ(ierr);
150b7281c8aSPeter Brune   }
151e4ed7901SPeter Brune 
152e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
153028b6e76SBarry Smith   snes->norm = fnorm;
154e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((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   /* test convergence */
159028b6e76SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
160028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
161028b6e76SBarry Smith 
162028b6e76SBarry Smith   /* Call general purpose update function */
163028b6e76SBarry Smith   if (snes->ops->update) {
164028b6e76SBarry Smith     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
165028b6e76SBarry Smith   }
166b7281c8aSPeter Brune 
167b7281c8aSPeter Brune   /* set parameter for default relative tolerance convergence test */
168b7281c8aSPeter Brune   snes->ttol = fnorm*snes->rtol;
169b7281c8aSPeter Brune   /* test convergence */
170b7281c8aSPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
171b7281c8aSPeter Brune   if (snes->reason) PetscFunctionReturn(0);
172b7281c8aSPeter Brune 
173b7281c8aSPeter Brune   for (i = 1; i < maxits+1; i++) {
174f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y);CHKERRQ(ierr);
175422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
176f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
177422a814eSBarry Smith     if (lsresult) {
178028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
179028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
180028b6e76SBarry Smith         break;
181028b6e76SBarry Smith       }
182028b6e76SBarry Smith     }
183028b6e76SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
184028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
185028b6e76SBarry Smith       break;
186028b6e76SBarry Smith     }
187b7281c8aSPeter Brune 
188028b6e76SBarry Smith     /* Monitor convergence */
189e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
190b7281c8aSPeter Brune     snes->iter = i;
191028b6e76SBarry Smith     snes->norm = fnorm;
192e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
193a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
194028b6e76SBarry Smith     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
195028b6e76SBarry Smith     /* Test for convergence */
196e7058c64SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
197028b6e76SBarry Smith     if (snes->reason) break;
198b7281c8aSPeter Brune 
199b7281c8aSPeter Brune     /* Call general purpose update function */
200b7281c8aSPeter Brune     if (snes->ops->update) {
201b7281c8aSPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
202028b6e76SBarry Smith     }
203b7281c8aSPeter Brune 
204*efd4aadfSBarry Smith     if (snes->npc) {
205b7281c8aSPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
206be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,NULL,Y);CHKERRQ(ierr);
207b7281c8aSPeter Brune         ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);
208b7281c8aSPeter Brune         ierr = VecCopy(Y,F);CHKERRQ(ierr);
209b7281c8aSPeter Brune       } else {
210be95d8f1SBarry Smith         ierr = SNESApplyNPC(snes,X,F,Y);CHKERRQ(ierr);
211b7281c8aSPeter Brune       }
212*efd4aadfSBarry Smith       ierr = SNESGetConvergedReason(snes->npc,&reason);CHKERRQ(ierr);
213b7281c8aSPeter Brune       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
214b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
215b7281c8aSPeter Brune         PetscFunctionReturn(0);
216b7281c8aSPeter Brune       }
217b7281c8aSPeter Brune     } else {
218b7281c8aSPeter Brune       ierr = VecCopy(F,Y);CHKERRQ(ierr);
219b7281c8aSPeter Brune     }
220b7281c8aSPeter Brune   }
221b7281c8aSPeter Brune   if (i == maxits+1) {
222028b6e76SBarry Smith     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
223028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
224028b6e76SBarry Smith   }
225028b6e76SBarry Smith   PetscFunctionReturn(0);
226028b6e76SBarry Smith }
227028b6e76SBarry Smith 
228028b6e76SBarry Smith /*MC
229d5c3842bSBarry Smith   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
230028b6e76SBarry Smith 
231028b6e76SBarry Smith   Level: beginner
232028b6e76SBarry Smith 
233028b6e76SBarry Smith   Options Database:
23472835e02SPeter Brune +   -snes_linesearch_type <l2,cp,basic> Line search type.
23572835e02SPeter Brune -   -snes_linesearch_damping<1.0> Damping for the line search.
236028b6e76SBarry Smith 
237af60355fSPeter Brune   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
238af60355fSPeter Brune             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
239be95d8f1SBarry Smith             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetNPC()) then the inner
240af60355fSPeter Brune             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
241af60355fSPeter Brune             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
242028b6e76SBarry Smith 
24372835e02SPeter Brune             The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
24472835e02SPeter Brune             linesearch, one may have to scale the update with -snes_linesearch_damping
24572835e02SPeter Brune 
246028b6e76SBarry Smith      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
247028b6e76SBarry Smith 
2482d547940SBarry Smith      Only supports left non-linear preconditioning.
2492d547940SBarry Smith 
25004d7464bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESNEWTONTR, SNESNGMRES, SNESQN, SNESNCG
251028b6e76SBarry Smith M*/
2528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
253028b6e76SBarry Smith {
25492c02d66SPeter Brune   PetscErrorCode   ierr;
255bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
2565fd66863SKarl Rupp 
257028b6e76SBarry Smith   PetscFunctionBegin;
258d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
259d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
260d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
261d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
262d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
263d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
264028b6e76SBarry Smith 
265028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
266*efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
267028b6e76SBarry Smith 
268*efd4aadfSBarry Smith   snes->npcside= PC_LEFT;
269c6b63b32SPeter Brune 
2704fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2714fc747eaSLawrence Mitchell 
272b00a9115SJed Brown   ierr       = PetscNewLog(snes,&neP);CHKERRQ(ierr);
273bf7f4e0aSPeter Brune   snes->data = (void*) neP;
274bf7f4e0aSPeter Brune 
27588976e71SPeter Brune   if (!snes->tolerancesset) {
2760e444f03SPeter Brune     snes->max_funcs = 30000;
2770e444f03SPeter Brune     snes->max_its   = 10000;
278c522fa08SPeter Brune     snes->stol      = 1e-20;
27988976e71SPeter Brune   }
280028b6e76SBarry Smith   PetscFunctionReturn(0);
281028b6e76SBarry Smith }
282