xref: /petsc/src/snes/impls/richardson/snesrichardson.c (revision ea630c6e93734680a55c9e053bce70cd1b499b87)
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;
52028b6e76SBarry Smith   ierr = SNESDefaultGetWork(snes,1);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) {
94*ea630c6eSPeter 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__ "SNESLineSearchNo_NRichardson"
101d5c3842bSBarry Smith PetscErrorCode SNESLineSearchNo_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
102028b6e76SBarry Smith {
103028b6e76SBarry Smith   PetscErrorCode   ierr;
104028b6e76SBarry Smith 
105028b6e76SBarry Smith   PetscFunctionBegin;
106*ea630c6eSPeter Brune   ierr = VecAXPY(X, snes->damping, Y);CHKERRQ(ierr);
107028b6e76SBarry Smith   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
108028b6e76SBarry Smith   ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
109028b6e76SBarry Smith   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated norm");
110028b6e76SBarry Smith   PetscFunctionReturn(0);
111028b6e76SBarry Smith }
112028b6e76SBarry Smith 
113028b6e76SBarry Smith #undef __FUNCT__
114d5c3842bSBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms_NRichardson"
115d5c3842bSBarry Smith PetscErrorCode SNESLineSearchNoNorms_NRichardson(SNES snes,void *lsctx,Vec X,Vec F,Vec Y,PetscReal fnorm,PetscReal dummyXnorm,Vec dummyG,Vec W,PetscReal *dummyYnorm,PetscReal *gnorm,PetscBool *flag)
116028b6e76SBarry Smith {
117028b6e76SBarry Smith   PetscErrorCode   ierr;
118028b6e76SBarry Smith 
119028b6e76SBarry Smith   PetscFunctionBegin;
120*ea630c6eSPeter Brune   ierr = VecAXPY(X, snes->damping, Y);CHKERRQ(ierr);
121028b6e76SBarry Smith   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
122028b6e76SBarry Smith   PetscFunctionReturn(0);
123028b6e76SBarry Smith }
124028b6e76SBarry Smith 
125028b6e76SBarry Smith #undef __FUNCT__
126d5c3842bSBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic_NRichardson"
127d5c3842bSBarry 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)
128028b6e76SBarry Smith {
129028b6e76SBarry Smith   PetscInt         i;
130028b6e76SBarry Smith   PetscReal        alphas[3] = {0.0, 0.5, 1.0};
131028b6e76SBarry Smith   PetscReal        norms[3];
132028b6e76SBarry Smith   PetscReal        alpha,a,b;
133028b6e76SBarry Smith   PetscErrorCode   ierr;
134028b6e76SBarry Smith 
135028b6e76SBarry Smith   PetscFunctionBegin;
136028b6e76SBarry Smith   norms[0]  = fnorm;
137028b6e76SBarry Smith   for(i=1; i < 3; ++i) {
138028b6e76SBarry Smith     ierr = VecWAXPY(W, alphas[i], Y, X);CHKERRQ(ierr);     /* W =  X^n - \alpha Y */
139028b6e76SBarry Smith     ierr = SNESComputeFunction(snes, W, F);CHKERRQ(ierr);
140028b6e76SBarry Smith     ierr = VecNorm(F, NORM_2, &norms[i]);CHKERRQ(ierr);
141028b6e76SBarry Smith   }
142028b6e76SBarry Smith   for(i = 0; i < 3; ++i) {
143028b6e76SBarry Smith     norms[i] = PetscSqr(norms[i]);
144028b6e76SBarry Smith   }
145028b6e76SBarry Smith   /* Fit a quadratic:
146028b6e76SBarry Smith        If we have x_{0,1,2} = 0, x_1, x_2 which generate norms y_{0,1,2}
147028b6e76SBarry 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)
148028b6e76SBarry 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)
149028b6e76SBarry Smith        c = y_0
150028b6e76SBarry Smith        x_min = -b/2a
151028b6e76SBarry Smith 
152028b6e76SBarry Smith        If we let x_{0,1,2} = 0, 0.5, 1.0
153028b6e76SBarry Smith        a = 2 y_2 - 4 y_1 + 2 y_0
154028b6e76SBarry Smith        b =  -y_2 + 4 y_1 - 3 y_0
155028b6e76SBarry Smith        c =   y_0
156028b6e76SBarry Smith   */
157028b6e76SBarry 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]));
158028b6e76SBarry 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]);
159028b6e76SBarry Smith   /* Check for positive a (concave up) */
160028b6e76SBarry Smith   if (a >= 0.0) {
161028b6e76SBarry Smith     alpha = -b/(2.0*a);
162028b6e76SBarry Smith     alpha = PetscMin(alpha, alphas[2]);
163028b6e76SBarry Smith     alpha = PetscMax(alpha, alphas[0]);
164028b6e76SBarry Smith   } else {
165028b6e76SBarry Smith     alpha = 1.0;
166028b6e76SBarry Smith   }
167*ea630c6eSPeter Brune   if (snes->ls_monitor) {
168*ea630c6eSPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
169*ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: norms[0] = %g, norms[1] = %g, norms[2] = %g alpha %g\n", sqrt(norms[0]),sqrt(norms[1]),sqrt(norms[2]),alpha);CHKERRQ(ierr);
170*ea630c6eSPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
171028b6e76SBarry Smith   }
172028b6e76SBarry Smith   ierr = VecAXPY(X, alpha, Y);CHKERRQ(ierr);
173028b6e76SBarry Smith   if (alpha != 1.0) {
174028b6e76SBarry Smith     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
175028b6e76SBarry Smith     ierr = VecNorm(F, NORM_2, gnorm);CHKERRQ(ierr);
176028b6e76SBarry Smith   } else {
177028b6e76SBarry Smith     *gnorm = PetscSqrtReal(norms[2]);
178028b6e76SBarry Smith   }
179028b6e76SBarry Smith   if (alpha == 0.0) *flag = PETSC_FALSE;
180028b6e76SBarry Smith   else              *flag = PETSC_TRUE;
181028b6e76SBarry Smith   PetscFunctionReturn(0);
182028b6e76SBarry Smith }
183028b6e76SBarry Smith 
184028b6e76SBarry Smith /*
185d5c3842bSBarry Smith   SNESSolve_NRichardson - Solves a nonlinear system with the Richardson method.
186028b6e76SBarry Smith 
187028b6e76SBarry Smith   Input Parameters:
188028b6e76SBarry Smith . snes - the SNES context
189028b6e76SBarry Smith 
190028b6e76SBarry Smith   Output Parameter:
191028b6e76SBarry Smith . outits - number of iterations until termination
192028b6e76SBarry Smith 
193028b6e76SBarry Smith   Application Interface Routine: SNESSolve()
194028b6e76SBarry Smith */
195028b6e76SBarry Smith #undef __FUNCT__
196d5c3842bSBarry Smith #define __FUNCT__ "SNESSolve_NRichardson"
197d5c3842bSBarry Smith PetscErrorCode SNESSolve_NRichardson(SNES snes)
198028b6e76SBarry Smith {
199028b6e76SBarry Smith   Vec                 X, Y, F, W;
200028b6e76SBarry Smith   PetscReal           fnorm;
201028b6e76SBarry Smith   PetscInt            maxits, i;
202028b6e76SBarry Smith   PetscErrorCode      ierr;
203028b6e76SBarry Smith   SNESConvergedReason reason;
204028b6e76SBarry Smith 
205028b6e76SBarry Smith   PetscFunctionBegin;
206028b6e76SBarry Smith   snes->reason = SNES_CONVERGED_ITERATING;
207028b6e76SBarry Smith 
208028b6e76SBarry Smith   maxits = snes->max_its;	     /* maximum number of iterations */
209028b6e76SBarry Smith   X      = snes->vec_sol;	     /* X^n */
210028b6e76SBarry Smith   Y      = snes->vec_sol_update; /* \tilde X */
211028b6e76SBarry Smith   F      = snes->vec_func;       /* residual vector */
212028b6e76SBarry Smith   W      = snes->work[0];        /* work vector */
213028b6e76SBarry Smith 
214028b6e76SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
215028b6e76SBarry Smith   snes->iter = 0;
216028b6e76SBarry Smith   snes->norm = 0.;
217028b6e76SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
218028b6e76SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
219028b6e76SBarry Smith   if (snes->domainerror) {
220028b6e76SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
221028b6e76SBarry Smith     PetscFunctionReturn(0);
222028b6e76SBarry Smith   }
223028b6e76SBarry Smith   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
224028b6e76SBarry Smith   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
225028b6e76SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
226028b6e76SBarry Smith   snes->norm = fnorm;
227028b6e76SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
228028b6e76SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
229028b6e76SBarry Smith   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
230028b6e76SBarry Smith 
231028b6e76SBarry Smith   /* set parameter for default relative tolerance convergence test */
232028b6e76SBarry Smith   snes->ttol = fnorm*snes->rtol;
233028b6e76SBarry Smith   /* test convergence */
234028b6e76SBarry Smith   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
235028b6e76SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
236028b6e76SBarry Smith 
237028b6e76SBarry Smith   for(i = 0; i < maxits; i++) {
238028b6e76SBarry Smith     PetscBool  lsSuccess = PETSC_TRUE;
239028b6e76SBarry Smith     PetscReal  dummyNorm;
240028b6e76SBarry Smith 
241028b6e76SBarry Smith     /* Call general purpose update function */
242028b6e76SBarry Smith     if (snes->ops->update) {
243028b6e76SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
244028b6e76SBarry Smith     }
245028b6e76SBarry Smith     if (!snes->pc) {
246028b6e76SBarry Smith       ierr = VecCopy(F,Y);CHKERRQ(ierr);
247028b6e76SBarry Smith       ierr = VecScale(Y,-1.0);CHKERRQ(ierr);
248028b6e76SBarry Smith     } else {
249028b6e76SBarry Smith       ierr = VecCopy(X,Y);CHKERRQ(ierr);
250c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, snes->vec_rhs, Y);CHKERRQ(ierr);
251028b6e76SBarry Smith       ierr = SNESGetConvergedReason(snes->pc,&reason);CHKERRQ(ierr);
25234b225ddSBarry Smith       if (reason < 0  && reason != SNES_DIVERGED_MAX_IT) {
253028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_INNER;
254028b6e76SBarry Smith         PetscFunctionReturn(0);
255028b6e76SBarry Smith       }
256028b6e76SBarry Smith       ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
257028b6e76SBarry Smith     }
258028b6e76SBarry Smith 
259*ea630c6eSPeter Brune       ierr = (*snes->ops->linesearch)(snes, snes->lsP, X, F, Y, fnorm, 0.0, 0, W, &dummyNorm, &fnorm, &lsSuccess);CHKERRQ(ierr);
260028b6e76SBarry Smith     if (!lsSuccess) {
261028b6e76SBarry Smith       if (++snes->numFailures >= snes->maxFailures) {
262028b6e76SBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
263028b6e76SBarry Smith         break;
264028b6e76SBarry Smith       }
265028b6e76SBarry Smith     }
266028b6e76SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
267028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
268028b6e76SBarry Smith       break;
269028b6e76SBarry Smith     }
270028b6e76SBarry Smith     if (snes->domainerror) {
271028b6e76SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
272028b6e76SBarry Smith       PetscFunctionReturn(0);
273028b6e76SBarry Smith     }
274028b6e76SBarry Smith     /* Monitor convergence */
275028b6e76SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
276028b6e76SBarry Smith     snes->iter = i+1;
277028b6e76SBarry Smith     snes->norm = fnorm;
278028b6e76SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
279028b6e76SBarry Smith     SNESLogConvHistory(snes,snes->norm,0);
280028b6e76SBarry Smith     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
281028b6e76SBarry Smith     /* Test for convergence */
282028b6e76SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
283028b6e76SBarry Smith     if (snes->reason) break;
284028b6e76SBarry Smith   }
285028b6e76SBarry Smith   if (i == maxits) {
286028b6e76SBarry Smith     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
287028b6e76SBarry Smith     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
288028b6e76SBarry Smith   }
289028b6e76SBarry Smith   PetscFunctionReturn(0);
290028b6e76SBarry Smith }
291028b6e76SBarry Smith 
292028b6e76SBarry Smith 
293028b6e76SBarry Smith /*MC
294d5c3842bSBarry Smith   SNESNRICHARDSON - Richardson nonlinear solver that uses successive substitutions, also sometimes known as Picard iteration.
295028b6e76SBarry Smith 
296028b6e76SBarry Smith   Level: beginner
297028b6e76SBarry Smith 
298028b6e76SBarry Smith   Options Database:
299028b6e76SBarry Smith +   -snes_ls_damping - damping factor to apply to F(x) (used only if -snes_ls is basic or basicnonorms)
300028b6e76SBarry Smith -   -snes_ls <basic,basicnormnorms,quadratic>
301028b6e76SBarry Smith 
302d5c3842bSBarry Smith   Notes: If no inner nonlinear preconditioner is provided then solves F(x) - b = 0 using x^{n+1} = x^{n} - lambda (F(x^n) - b)
303d5c3842bSBarry Smith             where lambda is obtained either SNESLineSearchSetDamping(), -snes_damping or a line search.
304d5c3842bSBarry Smith          If an inner nonlinear preconditioner is provided (either with -npc_snes_type or SNESSetPC()) then the inner solver is called an initial solution x^n and the
305d5c3842bSBarry Smith             nonlinear Richardson uses x^{n+1} = x^{n} + lambda d^{n} where d^{n} = \hat{x}^{n} - x^{n} where \hat{x}^{n} is the solution returned from the inner solver.
306028b6e76SBarry Smith 
307028b6e76SBarry Smith      This uses no derivative information thus will be much slower then Newton's method obtained with -snes_type ls
308028b6e76SBarry Smith 
309d5c3842bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESLS, SNESTR, SNESNGMRES, SNESNQN
310028b6e76SBarry Smith M*/
311028b6e76SBarry Smith EXTERN_C_BEGIN
312028b6e76SBarry Smith #undef __FUNCT__
313d5c3842bSBarry Smith #define __FUNCT__ "SNESCreate_NRichardson"
314d5c3842bSBarry Smith PetscErrorCode  SNESCreate_NRichardson(SNES snes)
315028b6e76SBarry Smith {
316028b6e76SBarry Smith   PetscFunctionBegin;
317d5c3842bSBarry Smith   snes->ops->destroy	     = SNESDestroy_NRichardson;
318d5c3842bSBarry Smith   snes->ops->setup	     = SNESSetUp_NRichardson;
319d5c3842bSBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_NRichardson;
320d5c3842bSBarry Smith   snes->ops->view            = SNESView_NRichardson;
321d5c3842bSBarry Smith   snes->ops->solve	     = SNESSolve_NRichardson;
322d5c3842bSBarry Smith   snes->ops->reset           = SNESReset_NRichardson;
323028b6e76SBarry Smith 
324028b6e76SBarry Smith   snes->usesksp              = PETSC_FALSE;
32542f4f86dSBarry Smith   snes->usespc               = PETSC_TRUE;
326028b6e76SBarry Smith 
327*ea630c6eSPeter Brune   snes->ops->linesearch         = SNESLineSearchQuadratic_NRichardson;
328*ea630c6eSPeter Brune   snes->ops->linesearchno       = SNESLineSearchNo_NRichardson;
329*ea630c6eSPeter Brune   snes->ops->linesearchnonorms  = SNESLineSearchNoNorms_NRichardson;
330*ea630c6eSPeter Brune   snes->ops->linesearchquadratic= SNESLineSearchQuadratic_NRichardson;
331*ea630c6eSPeter Brune   snes->ops->linesearchcubic    = PETSC_NULL;
332*ea630c6eSPeter Brune   snes->ops->postcheckstep      = PETSC_NULL;
333*ea630c6eSPeter Brune   snes->postcheck               = PETSC_NULL;
334*ea630c6eSPeter Brune   snes->ops->precheckstep       = PETSC_NULL;
335*ea630c6eSPeter Brune   snes->precheck                = PETSC_NULL;
336028b6e76SBarry Smith 
337028b6e76SBarry Smith   PetscFunctionReturn(0);
338028b6e76SBarry Smith }
339028b6e76SBarry Smith EXTERN_C_END
340