xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
142f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
2028b6e76SBarry Smith 
3d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NRichardson(SNES snes)
4d71ae5a4SJacob Faibussowitsch {
5028b6e76SBarry Smith   PetscFunctionBegin;
6*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 */
17d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NRichardson(SNES snes)
18d71ae5a4SJacob Faibussowitsch {
19028b6e76SBarry Smith   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(SNESReset_NRichardson(snes));
219566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
22*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23028b6e76SBarry Smith }
24028b6e76SBarry Smith 
25028b6e76SBarry Smith /*
26d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
27d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
28028b6e76SBarry Smith 
29028b6e76SBarry Smith    Input Parameters:
30028b6e76SBarry Smith +  snes - the SNES context
31028b6e76SBarry Smith -  x - the solution vector
32028b6e76SBarry Smith 
33028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
34028b6e76SBarry Smith  */
35d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NRichardson(SNES snes)
36d71ae5a4SJacob Faibussowitsch {
37028b6e76SBarry Smith   PetscFunctionBegin;
3808401ef6SPierre Jolivet   PetscCheck(snes->npcside != PC_RIGHT, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "NRichardson only supports left preconditioning");
396c67d002SPeter Brune   if (snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_UNPRECONDITIONED;
40*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41028b6e76SBarry Smith }
42028b6e76SBarry Smith 
43028b6e76SBarry Smith /*
4404d7464bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESNEWTONLS method.
45028b6e76SBarry Smith 
46028b6e76SBarry Smith   Input Parameter:
47028b6e76SBarry Smith . snes - the SNES context
48028b6e76SBarry Smith 
49028b6e76SBarry Smith   Application Interface Routine: SNESSetFromOptions()
50028b6e76SBarry Smith */
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes, PetscOptionItems *PetscOptionsObject)
52d71ae5a4SJacob Faibussowitsch {
53028b6e76SBarry Smith   PetscFunctionBegin;
54d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SNES Richardson options");
55d0609cedSBarry Smith   PetscOptionsHeadEnd();
56*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57028b6e76SBarry Smith }
58028b6e76SBarry Smith 
59028b6e76SBarry Smith /*
60d5c3842bSBarry Smith   SNESView_NRichardson - Prints info from the SNESRichardson data structure.
61028b6e76SBarry Smith 
62028b6e76SBarry Smith   Input Parameters:
63028b6e76SBarry Smith + SNES - the SNES context
64028b6e76SBarry Smith - viewer - visualization context
65028b6e76SBarry Smith 
66028b6e76SBarry Smith   Application Interface Routine: SNESView()
67028b6e76SBarry Smith */
68d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NRichardson(SNES snes, PetscViewer viewer)
69d71ae5a4SJacob Faibussowitsch {
70028b6e76SBarry Smith   PetscBool iascii;
71028b6e76SBarry Smith 
72028b6e76SBarry Smith   PetscFunctionBegin;
739566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
749371c9d4SSatish Balay   if (iascii) { }
75*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76028b6e76SBarry Smith }
77028b6e76SBarry Smith 
78028b6e76SBarry Smith /*
79d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
80028b6e76SBarry Smith 
81028b6e76SBarry Smith   Input Parameters:
82028b6e76SBarry Smith . snes - the SNES context
83028b6e76SBarry Smith 
84028b6e76SBarry Smith   Output Parameter:
85028b6e76SBarry Smith . outits - number of iterations until termination
86028b6e76SBarry Smith 
87028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
88028b6e76SBarry Smith */
89d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_NRichardson(SNES snes)
90d71ae5a4SJacob Faibussowitsch {
91bf7f4e0aSPeter Brune   Vec                  X, Y, F;
92e7058c64SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
93028b6e76SBarry Smith   PetscInt             maxits, i;
94422a814eSBarry Smith   SNESLineSearchReason lsresult;
95028b6e76SBarry Smith   SNESConvergedReason  reason;
96028b6e76SBarry Smith 
97028b6e76SBarry Smith   PetscFunctionBegin;
980b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
99c579b300SPatrick Farrell 
100028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
101028b6e76SBarry Smith 
102028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
103028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
104028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
105028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
106028b6e76SBarry Smith 
1079566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
108028b6e76SBarry Smith   snes->iter = 0;
109028b6e76SBarry Smith   snes->norm = 0.;
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
111b7281c8aSPeter Brune 
112efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1139566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1149566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
115b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
116b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
117*3ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
118b7281c8aSPeter Brune     }
1199566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
120b7281c8aSPeter Brune   } else {
121e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1229566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1231aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
124c1c75074SPeter Brune 
1259566063dSJacob Faibussowitsch     PetscCall(VecNorm(F, NORM_2, &fnorm));
126422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
127b7281c8aSPeter Brune   }
128efd4aadfSBarry Smith   if (snes->npc && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
1299566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, F, Y));
1309566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
131b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
132b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
133*3ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
134b7281c8aSPeter Brune     }
135b7281c8aSPeter Brune   } else {
1369566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, Y));
137b7281c8aSPeter Brune   }
138e4ed7901SPeter Brune 
1399566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
140028b6e76SBarry Smith   snes->norm = fnorm;
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1429566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1439566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
144028b6e76SBarry Smith 
145028b6e76SBarry Smith   /* test convergence */
146dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
147*3ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
148028b6e76SBarry Smith 
149028b6e76SBarry Smith   /* Call general purpose update function */
150dbbe0bcdSBarry Smith   PetscTryTypeMethod(snes, update, snes->iter);
151b7281c8aSPeter Brune 
152b7281c8aSPeter Brune   /* set parameter for default relative tolerance convergence test */
153b7281c8aSPeter Brune   snes->ttol = fnorm * snes->rtol;
154b7281c8aSPeter Brune   /* test convergence */
155dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
156*3ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
157b7281c8aSPeter Brune 
158b7281c8aSPeter Brune   for (i = 1; i < maxits + 1; i++) {
1599566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Y));
1609566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult));
1619566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm));
162422a814eSBarry Smith     if (lsresult) {
163028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
164028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
165028b6e76SBarry Smith         break;
166028b6e76SBarry Smith       }
167028b6e76SBarry Smith     }
168e71169deSBarry Smith     if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
169028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
170028b6e76SBarry Smith       break;
171028b6e76SBarry Smith     }
172b7281c8aSPeter Brune 
173028b6e76SBarry Smith     /* Monitor convergence */
1749566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
175b7281c8aSPeter Brune     snes->iter  = i;
176028b6e76SBarry Smith     snes->norm  = fnorm;
177c1e67a49SFande Kong     snes->xnorm = xnorm;
178c1e67a49SFande Kong     snes->ynorm = ynorm;
1799566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1809566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
1819566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
182028b6e76SBarry Smith     /* Test for convergence */
183dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
184028b6e76SBarry Smith     if (snes->reason) break;
185b7281c8aSPeter Brune 
186b7281c8aSPeter Brune     /* Call general purpose update function */
187dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
188b7281c8aSPeter Brune 
189efd4aadfSBarry Smith     if (snes->npc) {
190b7281c8aSPeter Brune       if (snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1919566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, NULL, Y));
1929566063dSJacob Faibussowitsch         PetscCall(VecNorm(F, NORM_2, &fnorm));
1939566063dSJacob Faibussowitsch         PetscCall(VecCopy(Y, F));
194b7281c8aSPeter Brune       } else {
1959566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, Y));
196b7281c8aSPeter Brune       }
1979566063dSJacob Faibussowitsch       PetscCall(SNESGetConvergedReason(snes->npc, &reason));
198b7281c8aSPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
199b7281c8aSPeter Brune         snes->reason = SNES_DIVERGED_INNER;
200*3ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
201b7281c8aSPeter Brune       }
202b7281c8aSPeter Brune     } else {
2039566063dSJacob Faibussowitsch       PetscCall(VecCopy(F, Y));
204b7281c8aSPeter Brune     }
205b7281c8aSPeter Brune   }
206b7281c8aSPeter Brune   if (i == maxits + 1) {
20763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
208028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
209028b6e76SBarry Smith   }
210*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
211028b6e76SBarry Smith }
212028b6e76SBarry Smith 
213028b6e76SBarry Smith /*MC
214d5c3842bSBarry Smith    SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
215028b6e76SBarry Smith 
216f6dfbefdSBarry Smith    Options Database Keys:
21767b8a455SSatish Balay +  -snes_linesearch_type <l2,cp,basic> - Line search type.
21867b8a455SSatish Balay -  -snes_linesearch_damping<1.0> - Damping for the line search.
219028b6e76SBarry Smith 
220f6dfbefdSBarry Smith    Level: beginner
221f6dfbefdSBarry Smith 
22295452b02SPatrick Sanan    Notes:
22395452b02SPatrick Sanan    If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
224f6dfbefdSBarry Smith    (F(x^n) - b) where lambda is obtained either `SNESLineSearchSetDamping()`, -snes_damping or a line search.  If
225f6dfbefdSBarry Smith    an inner nonlinear preconditioner is provided (either with -npc_snes_type or `SNESSetNPC())` then the inner
226af60355fSPeter Brune    solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
227af60355fSPeter Brune    where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
228028b6e76SBarry Smith 
22972835e02SPeter Brune    The update, especially without inner nonlinear preconditioner, may be ill-scaled.  If using the basic
23072835e02SPeter Brune    linesearch, one may have to scale the update with -snes_linesearch_damping
23172835e02SPeter Brune 
232028b6e76SBarry Smith    This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
233028b6e76SBarry Smith 
2342d547940SBarry Smith    Only supports left non-linear preconditioning.
2352d547940SBarry Smith 
236db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONLS`, `SNESNEWTONTR`, `SNESNGMRES`, `SNESQN`, `SNESNCG`
237028b6e76SBarry Smith M*/
238d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NRichardson(SNES snes)
239d71ae5a4SJacob Faibussowitsch {
240bf7f4e0aSPeter Brune   SNES_NRichardson *neP;
241d8d34be6SBarry Smith   SNESLineSearch    linesearch;
2425fd66863SKarl Rupp 
243028b6e76SBarry Smith   PetscFunctionBegin;
244d5c3842bSBarry Smith   snes->ops->destroy        = SNESDestroy_NRichardson;
245d5c3842bSBarry Smith   snes->ops->setup          = SNESSetUp_NRichardson;
246d5c3842bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NRichardson;
247d5c3842bSBarry Smith   snes->ops->view           = SNESView_NRichardson;
248d5c3842bSBarry Smith   snes->ops->solve          = SNESSolve_NRichardson;
249d5c3842bSBarry Smith   snes->ops->reset          = SNESReset_NRichardson;
250028b6e76SBarry Smith 
251028b6e76SBarry Smith   snes->usesksp = PETSC_FALSE;
252efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
253028b6e76SBarry Smith 
254efd4aadfSBarry Smith   snes->npcside = PC_LEFT;
255c6b63b32SPeter Brune 
2569566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
25748a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2));
258d8d34be6SBarry Smith 
2594fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
2604fc747eaSLawrence Mitchell 
2614dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
262bf7f4e0aSPeter Brune   snes->data = (void *)neP;
263bf7f4e0aSPeter Brune 
26488976e71SPeter Brune   if (!snes->tolerancesset) {
2650e444f03SPeter Brune     snes->max_funcs = 30000;
2660e444f03SPeter Brune     snes->max_its   = 10000;
267c522fa08SPeter Brune     snes->stol      = 1e-20;
26888976e71SPeter Brune   }
269*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
270028b6e76SBarry Smith }
271