xref: /petsc/src/snes/impls/tr/tr.c (revision 58c0e5077dcf40d6a880c19c87f1075ae1d22c8e)
1 
2 #include <../src/snes/impls/tr/trimpl.h>                /*I   "petscsnes.h"   I*/
3 
4 typedef struct {
5   SNES           snes;
6   /*  Information on the regular SNES convergence test; which may have been user provided */
7   PetscErrorCode (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
8   PetscErrorCode (*convdestroy)(void*);
9   void           *convctx;
10 } SNES_TR_KSPConverged_Ctx;
11 
12 static PetscErrorCode SNESTR_KSPConverged_Private(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *cctx)
13 {
14   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
15   SNES                     snes = ctx->snes;
16   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
17   Vec                      x;
18   PetscReal                nrm;
19   PetscErrorCode           ierr;
20 
21   PetscFunctionBegin;
22   ierr = (*ctx->convtest)(ksp,n,rnorm,reason,ctx->convctx);CHKERRQ(ierr);
23   if (*reason) {
24     ierr = PetscInfo2(snes,"Default or user provided convergence test KSP iterations=%D, rnorm=%g\n",n,(double)rnorm);CHKERRQ(ierr);
25   }
26   /* Determine norm of solution */
27   ierr = KSPBuildSolution(ksp,0,&x);CHKERRQ(ierr);
28   ierr = VecNorm(x,NORM_2,&nrm);CHKERRQ(ierr);
29   if (nrm >= neP->delta) {
30     ierr    = PetscInfo2(snes,"Ending linear iteration early, delta=%g, length=%g\n",(double)neP->delta,(double)nrm);CHKERRQ(ierr);
31     *reason = KSP_CONVERGED_STEP_LENGTH;
32   }
33   PetscFunctionReturn(0);
34 }
35 
36 static PetscErrorCode SNESTR_KSPConverged_Destroy(void *cctx)
37 {
38   SNES_TR_KSPConverged_Ctx *ctx = (SNES_TR_KSPConverged_Ctx*)cctx;
39   PetscErrorCode           ierr;
40 
41   PetscFunctionBegin;
42   ierr = (*ctx->convdestroy)(ctx->convctx);CHKERRQ(ierr);
43   ierr = PetscFree(ctx);CHKERRQ(ierr);
44   PetscFunctionReturn(0);
45 }
46 
47 /* ---------------------------------------------------------------- */
48 /*
49    SNESTR_Converged_Private -test convergence JUST for
50    the trust region tolerance.
51 
52 */
53 static PetscErrorCode SNESTR_Converged_Private(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
54 {
55   SNES_NEWTONTR  *neP = (SNES_NEWTONTR*)snes->data;
56   PetscErrorCode ierr;
57 
58   PetscFunctionBegin;
59   *reason = SNES_CONVERGED_ITERATING;
60   if (neP->delta < xnorm * snes->deltatol) {
61     ierr    = PetscInfo3(snes,"Converged due to trust region param %g<%g*%g\n",(double)neP->delta,(double)xnorm,(double)snes->deltatol);CHKERRQ(ierr);
62     *reason = SNES_DIVERGED_TR_DELTA;
63   } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
64     ierr    = PetscInfo1(snes,"Exceeded maximum number of function evaluations: %D\n",snes->max_funcs);CHKERRQ(ierr);
65     *reason = SNES_DIVERGED_FUNCTION_COUNT;
66   }
67   PetscFunctionReturn(0);
68 }
69 
70 /*@C
71    SNESNewtonTRSetPostCheck - Sets a user function that is called after the search step has been determined but before the next
72        function evaluation. Allows the user a chance to change or override the decision of the line search routine
73 
74    Logically Collective on snes
75 
76    Input Parameters:
77 +  snes - the nonlinear solver object
78 .  func - [optional] function evaluation routine, see SNESNewtonTRPostCheck()  for the calling sequence
79 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
80 
81    Level: intermediate
82 
83    Note: This function is called BEFORE the function evalaulation within the SNESNEWTONTR solver while the function set in
84    SNESLineSearchSetPostCheck() is called AFTER the function evaluation.
85 
86 .seealso: SNESNewtonTRPostCheck(), SNESNewtonTRGetPostCheck()
87 @*/
88 PetscErrorCode  SNESNewtonTRSetPostCheck(SNES snes, PetscErrorCode (*func)(SNES,Vec,Vec,Vec,PetscBool*,void*),void *ctx)
89 {
90   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
91 
92   PetscFunctionBegin;
93   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
94   if (func) tr->postcheck    = func;
95   if (ctx)  tr->postcheckctx = ctx;
96   PetscFunctionReturn(0);
97 }
98 
99 /*@C
100    SNESNewtonTRGetPostCheck - Gets the post-check function
101 
102    Not collective
103 
104    Input Parameter:
105 .  snes - the nonlinear solver context
106 
107    Output Parameters:
108 +  func - [optional] function evaluation routine, see for the calling sequence SNESNewtonTRPostCheck()
109 -  ctx  - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
110 
111    Level: intermediate
112 
113 .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRPostCheck()
114 @*/
115 PetscErrorCode  SNESNewtonTRGetPostCheck(SNES snes, PetscErrorCode (**func)(SNES,Vec,Vec,Vec,PetscBool*,void*),void **ctx)
116 {
117   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
118 
119   PetscFunctionBegin;
120   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
121   if (func) *func = tr->postcheck;
122   if (ctx)  *ctx  = tr->postcheckctx;
123   PetscFunctionReturn(0);
124 }
125 
126 /*@C
127    SNESNewtonTRPostCheck - Called after the step has been determined in SNESNEWTONTR but before the function evaluation
128 
129    Logically Collective on snes
130 
131    Input Parameters:
132 +  snes - the solver.  X - The last solution
133 .  Y - The full step direction
134 -  W - The updated solution, W = X + tao*Y for some tao
135 
136    Output Parameters:
137 .  changed_W - Indicator if the new candidate solution W has been changed.
138 
139    Notes:
140      If Y is changed then W is recomputed as X + tao*Y where tao is the current update value in SNESNEWTONTR
141 
142    Level: developer
143 
144 .seealso: SNESNewtonTRSetPostCheck(), SNESNewtonTRGetPostCheck()
145 @*/
146 static PetscErrorCode SNESNewtonTRPostCheck(SNES snes,Vec X,Vec Y,Vec W,PetscBool *changed_W)
147 {
148   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
149   PetscErrorCode ierr;
150 
151   PetscFunctionBegin;
152   *changed_W = PETSC_FALSE;
153   if (tr->postcheck) {
154     ierr = (*tr->postcheck)(snes,X,Y,W,changed_W,tr->postcheckctx);CHKERRQ(ierr);
155     PetscValidLogicalCollectiveBool(snes,*changed_W,6);
156   }
157   PetscFunctionReturn(0);
158 }
159 
160 /*
161    SNESSolve_NEWTONTR - Implements Newton's Method with a very simple trust
162    region approach for solving systems of nonlinear equations.
163 
164 
165 */
166 static PetscErrorCode SNESSolve_NEWTONTR(SNES snes)
167 {
168   SNES_NEWTONTR            *neP = (SNES_NEWTONTR*)snes->data;
169   Vec                      X,F,Y,G,Ytmp;
170   PetscErrorCode           ierr;
171   PetscInt                 maxits,i,lits;
172   PetscReal                rho,fnorm,gnorm,gpnorm,xnorm=0,delta,nrm,ynorm,norm1;
173   PetscScalar              cnorm;
174   KSP                      ksp;
175   SNESConvergedReason      reason = SNES_CONVERGED_ITERATING;
176   PetscBool                breakout = PETSC_FALSE;
177   SNES_TR_KSPConverged_Ctx *ctx;
178   PetscErrorCode           (*convtest)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*);
179 
180   PetscFunctionBegin;
181   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);
182 
183   maxits = snes->max_its;               /* maximum number of iterations */
184   X      = snes->vec_sol;               /* solution vector */
185   F      = snes->vec_func;              /* residual vector */
186   Y      = snes->work[0];               /* work vectors */
187   G      = snes->work[1];
188   Ytmp   = snes->work[2];
189 
190   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
191   snes->iter = 0;
192   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
193 
194   /* Set the linear stopping criteria to use the More' trick. */
195   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
196   ierr = KSPGetConvergenceTest(ksp,&convtest,NULL,NULL);CHKERRQ(ierr);
197   if (convtest != SNESTR_KSPConverged_Private) {
198     ierr                  = PetscNew(&ctx);CHKERRQ(ierr);
199     ctx->snes             = snes;
200     ierr                  = KSPGetAndClearConvergenceTest(ksp,&ctx->convtest,&ctx->convctx,&ctx->convdestroy);CHKERRQ(ierr);
201     ierr                  = KSPSetConvergenceTest(ksp,SNESTR_KSPConverged_Private,ctx,SNESTR_KSPConverged_Destroy);CHKERRQ(ierr);
202     ierr                  = PetscInfo(snes,"Using Krylov convergence test SNESTR_KSPConverged_Private\n");CHKERRQ(ierr);
203   }
204 
205   if (!snes->vec_func_init_set) {
206     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);          /* F(X) */
207   } else snes->vec_func_init_set = PETSC_FALSE;
208 
209   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
210   SNESCheckFunctionNorm(snes,fnorm);
211   ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);             /* fnorm <- || F || */
212   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
213   snes->norm = fnorm;
214   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
215   delta      = xnorm ? neP->delta0*xnorm : neP->delta0;
216   neP->delta = delta;
217   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
218   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
219 
220   /* test convergence */
221   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
222   if (snes->reason) PetscFunctionReturn(0);
223 
224 
225   for (i=0; i<maxits; i++) {
226 
227     /* Call general purpose update function */
228     if (snes->ops->update) {
229       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
230     }
231 
232     /* Solve J Y = F, where J is Jacobian matrix */
233     ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
234     SNESCheckJacobianDomainerror(snes);
235     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr);
236     ierr = KSPSolve(snes->ksp,F,Ytmp);CHKERRQ(ierr);
237     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
238 
239     snes->linear_its += lits;
240 
241     ierr  = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
242     ierr  = VecNorm(Ytmp,NORM_2,&nrm);CHKERRQ(ierr);
243     norm1 = nrm;
244     while (1) {
245       PetscBool changed_w;
246       ierr = VecCopy(Ytmp,Y);CHKERRQ(ierr);
247       nrm  = norm1;
248 
249       /* Scale Y if need be and predict new value of F norm */
250       if (nrm >= delta) {
251         nrm    = delta/nrm;
252         gpnorm = (1.0 - nrm)*fnorm;
253         cnorm  = nrm;
254         ierr   = PetscInfo1(snes,"Scaling direction by %g\n",(double)nrm);CHKERRQ(ierr);
255         ierr   = VecScale(Y,cnorm);CHKERRQ(ierr);
256         nrm    = gpnorm;
257         ynorm  = delta;
258       } else {
259         gpnorm = 0.0;
260         ierr   = PetscInfo(snes,"Direction is in Trust Region\n");CHKERRQ(ierr);
261         ynorm  = nrm;
262       }
263       ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
264       ierr = VecAYPX(Y,-1.0,X);CHKERRQ(ierr);            /* Y <- X - Y */
265       ierr = SNESNewtonTRPostCheck(snes,X,Ytmp,Y,&changed_w);CHKERRQ(ierr);
266       ierr = SNESComputeFunction(snes,Y,G);CHKERRQ(ierr); /*  F(X) */
267       ierr = VecNorm(G,NORM_2,&gnorm);CHKERRQ(ierr);      /* gnorm <- || g || */
268       SNESCheckFunctionNorm(snes,gnorm);
269       if (fnorm == gpnorm) rho = 0.0;
270       else rho = (fnorm*fnorm - gnorm*gnorm)/(fnorm*fnorm - gpnorm*gpnorm);
271 
272       /* Update size of trust region */
273       if      (rho < neP->mu)  delta *= neP->delta1;
274       else if (rho < neP->eta) delta *= neP->delta2;
275       else                     delta *= neP->delta3;
276       ierr = PetscInfo3(snes,"fnorm=%g, gnorm=%g, ynorm=%g\n",(double)fnorm,(double)gnorm,(double)ynorm);CHKERRQ(ierr);
277       ierr = PetscInfo3(snes,"gpred=%g, rho=%g, delta=%g\n",(double)gpnorm,(double)rho,(double)delta);CHKERRQ(ierr);
278 
279       neP->delta = delta;
280       if (rho > neP->sigma) break;
281       ierr = PetscInfo(snes,"Trying again in smaller region\n");CHKERRQ(ierr);
282       /* check to see if progress is hopeless */
283       neP->itflag = PETSC_FALSE;
284       ierr        = SNESTR_Converged_Private(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
285       if (!reason) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);}
286       if (reason == SNES_CONVERGED_SNORM_RELATIVE) reason = SNES_DIVERGED_INNER;
287       if (reason) {
288         /* We're not progressing, so return with the current iterate */
289         ierr     = SNESMonitor(snes,i+1,fnorm);CHKERRQ(ierr);
290         breakout = PETSC_TRUE;
291         break;
292       }
293       snes->numFailures++;
294     }
295     if (!breakout) {
296       /* Update function and solution vectors */
297       fnorm = gnorm;
298       ierr  = VecCopy(G,F);CHKERRQ(ierr);
299       ierr  = VecCopy(Y,X);CHKERRQ(ierr);
300       /* Monitor convergence */
301       ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
302       snes->iter = i+1;
303       snes->norm = fnorm;
304       snes->xnorm = xnorm;
305       snes->ynorm = ynorm;
306       ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
307       ierr       = SNESLogConvergenceHistory(snes,snes->norm,lits);CHKERRQ(ierr);
308       ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
309       /* Test for convergence, xnorm = || X || */
310       neP->itflag = PETSC_TRUE;
311       if (snes->ops->converged != SNESConvergedSkip) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
312       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
313       if (reason) break;
314     } else break;
315   }
316   if (i == maxits) {
317     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
318     if (!reason) reason = SNES_DIVERGED_MAX_IT;
319   }
320   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
321   snes->reason = reason;
322   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 /*------------------------------------------------------------*/
326 static PetscErrorCode SNESSetUp_NEWTONTR(SNES snes)
327 {
328   PetscErrorCode ierr;
329 
330   PetscFunctionBegin;
331   ierr = SNESSetWorkVecs(snes,3);CHKERRQ(ierr);
332   ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
333   PetscFunctionReturn(0);
334 }
335 
336 PetscErrorCode SNESReset_NEWTONTR(SNES snes)
337 {
338 
339   PetscFunctionBegin;
340   PetscFunctionReturn(0);
341 }
342 
343 static PetscErrorCode SNESDestroy_NEWTONTR(SNES snes)
344 {
345   PetscErrorCode ierr;
346 
347   PetscFunctionBegin;
348   ierr = SNESReset_NEWTONTR(snes);CHKERRQ(ierr);
349   ierr = PetscFree(snes->data);CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352 /*------------------------------------------------------------*/
353 
354 static PetscErrorCode SNESSetFromOptions_NEWTONTR(PetscOptionItems *PetscOptionsObject,SNES snes)
355 {
356   SNES_NEWTONTR  *ctx = (SNES_NEWTONTR*)snes->data;
357   PetscErrorCode ierr;
358 
359   PetscFunctionBegin;
360   ierr = PetscOptionsHead(PetscOptionsObject,"SNES trust region options for nonlinear equations");CHKERRQ(ierr);
361   ierr = PetscOptionsReal("-snes_trtol","Trust region tolerance","SNESSetTrustRegionTolerance",snes->deltatol,&snes->deltatol,NULL);CHKERRQ(ierr);
362   ierr = PetscOptionsReal("-snes_tr_mu","mu","None",ctx->mu,&ctx->mu,NULL);CHKERRQ(ierr);
363   ierr = PetscOptionsReal("-snes_tr_eta","eta","None",ctx->eta,&ctx->eta,NULL);CHKERRQ(ierr);
364   ierr = PetscOptionsReal("-snes_tr_sigma","sigma","None",ctx->sigma,&ctx->sigma,NULL);CHKERRQ(ierr);
365   ierr = PetscOptionsReal("-snes_tr_delta0","delta0","None",ctx->delta0,&ctx->delta0,NULL);CHKERRQ(ierr);
366   ierr = PetscOptionsReal("-snes_tr_delta1","delta1","None",ctx->delta1,&ctx->delta1,NULL);CHKERRQ(ierr);
367   ierr = PetscOptionsReal("-snes_tr_delta2","delta2","None",ctx->delta2,&ctx->delta2,NULL);CHKERRQ(ierr);
368   ierr = PetscOptionsReal("-snes_tr_delta3","delta3","None",ctx->delta3,&ctx->delta3,NULL);CHKERRQ(ierr);
369   ierr = PetscOptionsTail();CHKERRQ(ierr);
370   PetscFunctionReturn(0);
371 }
372 
373 static PetscErrorCode SNESView_NEWTONTR(SNES snes,PetscViewer viewer)
374 {
375   SNES_NEWTONTR  *tr = (SNES_NEWTONTR*)snes->data;
376   PetscErrorCode ierr;
377   PetscBool      iascii;
378 
379   PetscFunctionBegin;
380   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
381   if (iascii) {
382     ierr = PetscViewerASCIIPrintf(viewer,"  Trust region tolerance (-snes_trtol)\n",(double)snes->deltatol);CHKERRQ(ierr);
383     ierr = PetscViewerASCIIPrintf(viewer,"  mu=%g, eta=%g, sigma=%g\n",(double)tr->mu,(double)tr->eta,(double)tr->sigma);CHKERRQ(ierr);
384     ierr = PetscViewerASCIIPrintf(viewer,"  delta0=%g, delta1=%g, delta2=%g, delta3=%g\n",(double)tr->delta0,(double)tr->delta1,(double)tr->delta2,(double)tr->delta3);CHKERRQ(ierr);
385   }
386   PetscFunctionReturn(0);
387 }
388 /* ------------------------------------------------------------ */
389 /*MC
390       SNESNEWTONTR - Newton based nonlinear solver that uses a trust region
391 
392    Options Database:
393 +    -snes_trtol <tol> - trust region tolerance
394 .    -snes_tr_mu <mu> - trust region parameter
395 .    -snes_tr_eta <eta> - trust region parameter
396 .    -snes_tr_sigma <sigma> - trust region parameter
397 .    -snes_tr_delta0 <delta0> -  initial size of the trust region is delta0*norm2(x)
398 .    -snes_tr_delta1 <delta1> - trust region parameter
399 .    -snes_tr_delta2 <delta2> - trust region parameter
400 -    -snes_tr_delta3 <delta3> - trust region parameter
401 
402    The basic algorithm is taken from "The Minpack Project", by More',
403    Sorensen, Garbow, Hillstrom, pages 88-111 of "Sources and Development
404    of Mathematical Software", Wayne Cowell, editor.
405 
406    Level: intermediate
407 
408 .seealso:  SNESCreate(), SNES, SNESSetType(), SNESNEWTONLS, SNESSetTrustRegionTolerance()
409 
410 M*/
411 PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONTR(SNES snes)
412 {
413   SNES_NEWTONTR  *neP;
414   PetscErrorCode ierr;
415 
416   PetscFunctionBegin;
417   snes->ops->setup          = SNESSetUp_NEWTONTR;
418   snes->ops->solve          = SNESSolve_NEWTONTR;
419   snes->ops->destroy        = SNESDestroy_NEWTONTR;
420   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONTR;
421   snes->ops->view           = SNESView_NEWTONTR;
422   snes->ops->reset          = SNESReset_NEWTONTR;
423 
424   snes->usesksp = PETSC_TRUE;
425   snes->usesnpc = PETSC_FALSE;
426 
427   snes->alwayscomputesfinalresidual = PETSC_TRUE;
428 
429   ierr        = PetscNewLog(snes,&neP);CHKERRQ(ierr);
430   snes->data  = (void*)neP;
431   neP->mu     = 0.25;
432   neP->eta    = 0.75;
433   neP->delta  = 0.0;
434   neP->delta0 = 0.2;
435   neP->delta1 = 0.3;
436   neP->delta2 = 0.75;
437   neP->delta3 = 2.0;
438   neP->sigma  = 0.0001;
439   neP->itflag = PETSC_FALSE;
440   neP->rnorm0 = 0.0;
441   neP->ttol   = 0.0;
442   PetscFunctionReturn(0);
443 }
444 
445