xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision e5583369c2ebdb51c6868e1dcb6336efbdb688fd)
1028b6e76SBarry Smith 
242f4f86dSBarry Smith #include <../src/snes/impls/richardson/snesrichardsonimpl.h>
3028b6e76SBarry Smith 
4028b6e76SBarry Smith #undef __FUNCT__
5d5c3842bSBarry Smith #define __FUNCT__ "SNESReset_NRichardson"
6d5c3842bSBarry Smith PetscErrorCode SNESReset_NRichardson(SNES snes)
7028b6e76SBarry Smith {
8028b6e76SBarry Smith   PetscErrorCode ierr;
9028b6e76SBarry Smith 
10028b6e76SBarry Smith   PetscFunctionBegin;
11028b6e76SBarry Smith   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
12028b6e76SBarry Smith   PetscFunctionReturn(0);
13028b6e76SBarry Smith }
14028b6e76SBarry Smith 
15028b6e76SBarry Smith /*
16d5c3842bSBarry Smith   SNESDestroy_NRichardson - Destroys the private SNES_NRichardson context that was created with SNESCreate_NRichardson().
17028b6e76SBarry Smith 
18028b6e76SBarry Smith   Input Parameter:
19028b6e76SBarry Smith . snes - the SNES context
20028b6e76SBarry Smith 
21028b6e76SBarry Smith   Application Interface Routine: SNESDestroy()
22028b6e76SBarry Smith */
23028b6e76SBarry Smith #undef __FUNCT__
24d5c3842bSBarry Smith #define __FUNCT__ "SNESDestroy_NRichardson"
25d5c3842bSBarry Smith PetscErrorCode SNESDestroy_NRichardson(SNES snes)
26028b6e76SBarry Smith {
27028b6e76SBarry Smith   PetscErrorCode   ierr;
28028b6e76SBarry Smith 
29028b6e76SBarry Smith   PetscFunctionBegin;
30d5c3842bSBarry Smith   ierr = SNESReset_NRichardson(snes);CHKERRQ(ierr);
31028b6e76SBarry Smith   ierr = PetscFree(snes->data);CHKERRQ(ierr);
32028b6e76SBarry Smith   PetscFunctionReturn(0);
33028b6e76SBarry Smith }
34028b6e76SBarry Smith 
35028b6e76SBarry Smith /*
36d5c3842bSBarry Smith    SNESSetUp_NRichardson - Sets up the internal data structures for the later use
37d5c3842bSBarry Smith    of the SNESNRICHARDSON nonlinear solver.
38028b6e76SBarry Smith 
39028b6e76SBarry Smith    Input Parameters:
40028b6e76SBarry Smith +  snes - the SNES context
41028b6e76SBarry Smith -  x - the solution vector
42028b6e76SBarry Smith 
43028b6e76SBarry Smith    Application Interface Routine: SNESSetUp()
44028b6e76SBarry Smith  */
45028b6e76SBarry Smith #undef __FUNCT__
46d5c3842bSBarry Smith #define __FUNCT__ "SNESSetUp_NRichardson"
47d5c3842bSBarry Smith PetscErrorCode SNESSetUp_NRichardson(SNES snes)
48028b6e76SBarry Smith {
49028b6e76SBarry Smith   PetscErrorCode ierr;
50028b6e76SBarry Smith 
51028b6e76SBarry Smith   PetscFunctionBegin;
52f9a8e2e6SPeter Brune   ierr = SNESDefaultGetWork(snes,2);CHKERRQ(ierr);
53028b6e76SBarry Smith   PetscFunctionReturn(0);
54028b6e76SBarry Smith }
55028b6e76SBarry Smith 
56028b6e76SBarry Smith /*
57d5c3842bSBarry Smith   SNESSetFromOptions_NRichardson - Sets various parameters for the SNESLS method.
58028b6e76SBarry Smith 
59028b6e76SBarry Smith   Input Parameter:
60028b6e76SBarry Smith . snes - the SNES context
61028b6e76SBarry Smith 
62028b6e76SBarry Smith   Application Interface Routine: SNESSetFromOptions()
63028b6e76SBarry Smith */
64028b6e76SBarry Smith #undef __FUNCT__
65d5c3842bSBarry Smith #define __FUNCT__ "SNESSetFromOptions_NRichardson"
66d5c3842bSBarry Smith static PetscErrorCode SNESSetFromOptions_NRichardson(SNES snes)
67028b6e76SBarry Smith {
68028b6e76SBarry Smith   PetscErrorCode ierr;
69028b6e76SBarry Smith   PetscFunctionBegin;
7042f4f86dSBarry Smith     ierr = PetscOptionsHead("SNES Richardson options");CHKERRQ(ierr);
71028b6e76SBarry Smith     ierr = PetscOptionsTail();CHKERRQ(ierr);
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) {
94ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  richardson variant: %s\n", SNESLineSearchTypeName(snes->ls_type));CHKERRQ(ierr);
95028b6e76SBarry Smith   }
96028b6e76SBarry Smith   PetscFunctionReturn(0);
97028b6e76SBarry Smith }
98028b6e76SBarry Smith 
99028b6e76SBarry Smith #undef __FUNCT__
100d5c3842bSBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic_NRichardson"
101d5c3842bSBarry Smith PetscErrorCode SNESLineSearchQuadratic_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec G,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
102028b6e76SBarry Smith {
103028b6e76SBarry Smith   PetscInt         i;
104028b6e76SBarry Smith   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
105028b6e76SBarry Smith   PetscReal        norms[3];
106028b6e76SBarry Smith   PetscReal        alpha,a,b;
107028b6e76SBarry Smith   PetscErrorCode   ierr;
108028b6e76SBarry Smith 
109028b6e76SBarry Smith   PetscFunctionBegin;
110028b6e76SBarry Smith   norms[0]  = fnorm;
111028b6e76SBarry Smith   for(i=1; i < 3; ++i) {
112eb1825c3SPeter Brune     ierr = VecWAXPY(W, -alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
113f9a8e2e6SPeter Brune     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
114f9a8e2e6SPeter Brune     ierr = VecNorm(G, NORM_2, &norms[i]);CHKERRQ(ierr);
115028b6e76SBarry Smith   }
116028b6e76SBarry Smith   for(i = 0; i < 3; ++i) {
117028b6e76SBarry Smith     norms[i] = PetscSqr(norms[i]);
118028b6e76SBarry Smith   }
119028b6e76SBarry Smith   /* Fit a quadratic:
120028b6e76SBarry Smith        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
121028b6e76SBarry Smith        a = (x_1 y_2 - x_2 y_1 + (x_2 - x_1) y_0)/(x^2_2 x_1 - x_2 x^2_1)
122028b6e76SBarry Smith        b = (x^2_1 y_2 - x^2_2 y_1 + (x^2_2 - x^2_1) y_0)/(x_2 x^2_1 - x^2_2 x_1)
123028b6e76SBarry Smith        c = y_0
124028b6e76SBarry Smith        x_min = -b/2a
125028b6e76SBarry Smith 
126028b6e76SBarry Smith        If we let x_{0,1,2} = 0, 0.5, 1.0
127028b6e76SBarry Smith        a = 2 y_2 - 4 y_1 + 2 y_0
128028b6e76SBarry Smith        b =  -y_2 + 4 y_1 - 3 y_0
129028b6e76SBarry Smith        c =   y_0
130028b6e76SBarry Smith   */
131028b6e76SBarry Smith   a = (alphas[1]*norms[2] - alphas[2]*norms[1] + (alphas[2] - alphas[1])*norms[0])/(PetscSqr(alphas[2])*alphas[1] - alphas[2]*PetscSqr(alphas[1]));
132028b6e76SBarry Smith   b = (PetscSqr(alphas[1])*norms[2] - PetscSqr(alphas[2])*norms[1] + (PetscSqr(alphas[2]) - PetscSqr(alphas[1]))*norms[0])/(alphas[2]*PetscSqr(alphas[1]) - PetscSqr(alphas[2])*alphas[1]);
133028b6e76SBarry Smith   /* Check for positive a (concave up) */
134028b6e76SBarry Smith   if (a >= 0.0) {
135028b6e76SBarry Smith     alpha = -b/(2.0*a);
136028b6e76SBarry Smith     alpha = PetscMin(alpha, alphas[2]);
137028b6e76SBarry Smith     alpha = PetscMax(alpha, alphas[0]);
138028b6e76SBarry Smith   } else {
139028b6e76SBarry Smith     alpha = 1.0;
140028b6e76SBarry Smith   }
141ea630c6eSPeter Brune   if (snes->ls_monitor) {
142ea630c6eSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
143*cf22de31SPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n",
144*cf22de31SPeter Brune                                   PetscSqrtReal(norms[0]),PetscSqrtReal(norms[1]),PetscSqrtReal(norms[2]),alpha);CHKERRQ(ierr);
145ea630c6eSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
146028b6e76SBarry Smith   }
147f9a8e2e6SPeter Brune   ierr = VecCopy(X, W);CHKERRQ(ierr);
148eb1825c3SPeter Brune   ierr = VecAXPY(W, -alpha, Y);CHKERRQ(ierr);
149028b6e76SBarry Smith   if (alpha != 1.0) {
150f9a8e2e6SPeter Brune     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
151f9a8e2e6SPeter Brune     ierr = VecNorm(G, NORM_2, gnorm);CHKERRQ(ierr);
152028b6e76SBarry Smith   } else {
153028b6e76SBarry Smith     *gnorm = PetscSqrtReal(norms[2]);
154028b6e76SBarry Smith   }
155028b6e76SBarry Smith   if (alpha == 0.0) *flag = PETSC_FALSE;
156028b6e76SBarry Smith   else              *flag = PETSC_TRUE;
157028b6e76SBarry Smith   PetscFunctionReturn(0);
158028b6e76SBarry Smith }
159028b6e76SBarry Smith 
160028b6e76SBarry Smith /*
161d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
162028b6e76SBarry Smith 
163028b6e76SBarry Smith   Input Parameters:
164028b6e76SBarry Smith . snes - the SNES context
165028b6e76SBarry Smith 
166028b6e76SBarry Smith   Output Parameter:
167028b6e76SBarry Smith . outits - number of iterations until termination
168028b6e76SBarry Smith 
169028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
170028b6e76SBarry Smith */
171028b6e76SBarry Smith #undef __FUNCT__
172d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson"
173d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes)
174028b6e76SBarry Smith {
175f9a8e2e6SPeter Brune   Vec                 X, Y, F, W, G;
176028b6e76SBarry Smith   PetscReal           fnorm;
177028b6e76SBarry Smith   PetscInt            maxits, i;
178028b6e76SBarry Smith   PetscErrorCode      ierr;
179028b6e76SBarry Smith   SNESConvergedReason reason;
180028b6e76SBarry Smith 
181028b6e76SBarry Smith   PetscFunctionBegin;
182028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
183028b6e76SBarry Smith 
184028b6e76SBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
185028b6e76SBarry Smith   X      = snes->vec_sol;        /* X^n */
186028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
187028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
188028b6e76SBarry Smith   W      = snes->work[0];        /* work vector */
189f9a8e2e6SPeter Brune   G      = snes->work[1];        /* line search function vector */
190028b6e76SBarry Smith 
191028b6e76SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
192028b6e76SBarry Smith   snes->iter = 0;
193028b6e76SBarry Smith   snes->norm = 0.;
194028b6e76SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
195028b6e76SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
196028b6e76SBarry Smith   if (snes->domainerror) {
197028b6e76SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
198028b6e76SBarry Smith     PetscFunctionReturn(0);
199028b6e76SBarry Smith   }
200028b6e76SBarry Smith   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
201028b6e76SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
202028b6e76SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
203028b6e76SBarry Smith   snes->norm = fnorm;
204028b6e76SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
205028b6e76SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
206028b6e76SBarry Smith   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
207028b6e76SBarry Smith 
208028b6e76SBarry Smith   /* set parameter for default relative tolerance convergence test */
209028b6e76SBarry Smith   snes->ttol = fnorm*snes->rtol;
210028b6e76SBarry Smith   /* test convergence */
211028b6e76SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
212028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
213028b6e76SBarry Smith 
214028b6e76SBarry Smith   for(i = 0; i < maxits; i++) {
215028b6e76SBarry Smith     PetscBool  lsSuccess = PETSC_TRUE;
216028b6e76SBarry Smith     PetscReal  dummyNorm;
217028b6e76SBarry Smith 
218028b6e76SBarry Smith     /* Call general purpose update function */
219028b6e76SBarry Smith     if (snes->ops->update) {
220028b6e76SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
221028b6e76SBarry Smith     }
222a3685310SPeter Brune     if (snes->usegs && snes->ops->computegs) {
223a3685310SPeter Brune       ierr = VecCopy(X, Y);CHKERRQ(ierr);
224a3685310SPeter Brune       ierr = SNESComputeGS(snes, snes->vec_rhs, Y);CHKERRQ(ierr);
225e5583369SPeter Brune       ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr);
226a3685310SPeter Brune     } else if (snes->pc) {
227028b6e76SBarry Smith       ierr = VecCopy(X,Y);CHKERRQ(ierr);
228c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
229028b6e76SBarry Smith       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
23034b225ddSBarry Smith       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
231028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_INNER;
232028b6e76SBarry Smith         PetscFunctionReturn(0);
233028b6e76SBarry Smith       }
234e5583369SPeter Brune       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);
235a3685310SPeter Brune     } else {
236a3685310SPeter Brune       ierr = VecCopy(F,Y);CHKERRQ(ierr);
237028b6e76SBarry Smith     }
238f9a8e2e6SPeter Brune     ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, Y, fnorm, 0.0, G, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
239028b6e76SBarry Smith     if (!lsSuccess) {
240028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
241028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
242028b6e76SBarry Smith         break;
243028b6e76SBarry Smith       }
244028b6e76SBarry Smith     }
245028b6e76SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
246028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
247028b6e76SBarry Smith       break;
248028b6e76SBarry Smith     }
249028b6e76SBarry Smith     if (snes->domainerror) {
250028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
251028b6e76SBarry Smith       PetscFunctionReturn(0);
252028b6e76SBarry Smith     }
253f9a8e2e6SPeter Brune     ierr = VecCopy(G, F);CHKERRQ(ierr);
254f9a8e2e6SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
255028b6e76SBarry Smith     /* Monitor convergence */
256028b6e76SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
257028b6e76SBarry Smith     snes->iter = i+1;
258028b6e76SBarry Smith     snes->norm = fnorm;
259028b6e76SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
260028b6e76SBarry Smith     SNESLogConvHistory(snes,snes->norm,0);
261028b6e76SBarry Smith     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
262028b6e76SBarry Smith     /* Test for convergence */
263028b6e76SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
264028b6e76SBarry Smith     if (snes->reason) break;
265028b6e76SBarry Smith   }
266028b6e76SBarry Smith   if (i == maxits) {
267028b6e76SBarry Smith     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
268028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
269028b6e76SBarry Smith   }
270028b6e76SBarry Smith   PetscFunctionReturn(0);
271028b6e76SBarry Smith }
272028b6e76SBarry Smith 
273028b6e76SBarry Smith 
27492c02d66SPeter Brune EXTERN_C_BEGIN
27592c02d66SPeter Brune #undef __FUNCT__
27692c02d66SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_NRichardson"
27792c02d66SPeter Brune PetscErrorCode  SNESLineSearchSetType_NRichardson(SNES snes, SNESLineSearchType type)
27892c02d66SPeter Brune {
27992c02d66SPeter Brune   PetscErrorCode ierr;
28092c02d66SPeter Brune   PetscFunctionBegin;
28192c02d66SPeter Brune 
28292c02d66SPeter Brune   switch (type) {
28392c02d66SPeter Brune   case SNES_LS_BASIC:
284f9a8e2e6SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
28592c02d66SPeter Brune     break;
28692c02d66SPeter Brune   case SNES_LS_BASIC_NONORMS:
287f9a8e2e6SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
28892c02d66SPeter Brune     break;
28992c02d66SPeter Brune   case SNES_LS_QUADRATIC:
29092c02d66SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_NRichardson,PETSC_NULL);CHKERRQ(ierr);
29192c02d66SPeter Brune     break;
292f9a8e2e6SPeter Brune   case SNES_LS_SECANT:
293af60355fSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
294f9a8e2e6SPeter Brune     break;
29592c02d66SPeter Brune   default:
29692c02d66SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type");
29792c02d66SPeter Brune     break;
29892c02d66SPeter Brune   }
29992c02d66SPeter Brune   snes->ls_type = type;
30092c02d66SPeter Brune   PetscFunctionReturn(0);
30192c02d66SPeter Brune }
30292c02d66SPeter Brune EXTERN_C_END
30392c02d66SPeter Brune 
304028b6e76SBarry Smith /*MC
305d5c3842bSBarry Smith   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
306028b6e76SBarry Smith 
307028b6e76SBarry Smith   Level: beginner
308028b6e76SBarry Smith 
309028b6e76SBarry Smith   Options Database:
310028b6e76SBarry Smith +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
311028b6e76SBarry Smith -   -snes_ls <basic,basicnormnorms,quadratic>
312028b6e76SBarry Smith 
313af60355fSPeter Brune   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda
314af60355fSPeter Brune             (F(x^n) - b) where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.  If
315af60355fSPeter Brune             an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner
316af60355fSPeter Brune             solver is called an initial solution x^n and the nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n}
317af60355fSPeter Brune             where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
318028b6e76SBarry Smith 
319028b6e76SBarry Smith      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
320028b6e76SBarry Smith 
321d5c3842bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
322028b6e76SBarry Smith M*/
323028b6e76SBarry Smith EXTERN_C_BEGIN
324028b6e76SBarry Smith #undef __FUNCT__
325d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson"
326d5c3842bSBarry Smith PetscErrorCode  SNESCreate_NRichardson(SNES snes)
327028b6e76SBarry Smith {
32892c02d66SPeter Brune   PetscErrorCode ierr;
329028b6e76SBarry Smith   PetscFunctionBegin;
330d5c3842bSBarry Smith   snes->ops->destroy         = SNESDestroy_NRichardson;
331d5c3842bSBarry Smith   snes->ops->setup           = SNESSetUp_NRichardson;
332d5c3842bSBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
333d5c3842bSBarry Smith   snes->ops->view            = SNESView_NRichardson;
334d5c3842bSBarry Smith   snes->ops->solve           = SNESSolve_NRichardson;
335d5c3842bSBarry Smith   snes->ops->reset           = SNESReset_NRichardson;
336028b6e76SBarry Smith 
337028b6e76SBarry Smith   snes->usesksp              = PETSC_FALSE;
33842f4f86dSBarry Smith   snes->usespc               = PETSC_TRUE;
339028b6e76SBarry Smith 
34092c02d66SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_NRichardson",SNESLineSearchSetType_NRichardson);CHKERRQ(ierr);
341af60355fSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_SECANT);CHKERRQ(ierr);
342028b6e76SBarry Smith 
343028b6e76SBarry Smith   PetscFunctionReturn(0);
344028b6e76SBarry Smith }
345028b6e76SBarry Smith EXTERN_C_END
346