xref: /petsc/src/snes/impls/ls/ls.c (revision a5892487d3a72298dc86de4f4aaa34ef49eae3db)
15e76c431SBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h>
35e42470aSBarry Smith 
48a5d9ee4SBarry Smith /*
520f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
620f52e01SBarry Smith     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
736669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
820f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
98a5d9ee4SBarry Smith */
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal      a1;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
16ace3abfcSBarry Smith   PetscBool      hastranspose;
178a5d9ee4SBarry Smith 
188a5d9ee4SBarry Smith   PetscFunctionBegin;
198a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2036669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2136669109SBarry Smith   if (hastranspose) {
228a5d9ee4SBarry Smith     /* Compute || J^T F|| */
238a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
248a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
258f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
268a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2736669109SBarry Smith   } else {
2836669109SBarry Smith     Vec         work;
29ea709b57SSatish Balay     PetscScalar result;
3036669109SBarry Smith     PetscReal   wnorm;
3136669109SBarry Smith 
32abb0e124SMatthew Knepley     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
3336669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
376bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
3836669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
398f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
4036669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4136669109SBarry Smith   }
428a5d9ee4SBarry Smith   PetscFunctionReturn(0);
438a5d9ee4SBarry Smith }
448a5d9ee4SBarry Smith 
4574637425SBarry Smith /*
465ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4774637425SBarry Smith */
484a2ae208SSatish Balay #undef __FUNCT__
494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
501e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
5174637425SBarry Smith {
5274637425SBarry Smith   PetscReal      a1,a2;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
54ace3abfcSBarry Smith   PetscBool      hastranspose;
5574637425SBarry Smith 
5674637425SBarry Smith   PetscFunctionBegin;
5774637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5874637425SBarry Smith   if (hastranspose) {
5974637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6079f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6174637425SBarry Smith 
6274637425SBarry Smith     /* Compute || J^T W|| */
6374637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6474637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6574637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
6675567043SBarry Smith     if (a1 != 0.0) {
678f1a2a5eSBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
6885471664SBarry Smith     }
6974637425SBarry Smith   }
7074637425SBarry Smith   PetscFunctionReturn(0);
7174637425SBarry Smith }
7274637425SBarry Smith 
73*a5892487SPeter Brune #undef __FUNCT__
74*a5892487SPeter Brune #define __FUNCT__ "SNESLineSearchCubic_LS"
75*a5892487SPeter Brune /*@C
76*a5892487SPeter Brune    SNESLineSearchCubic - Performs a cubic line search (default line search method).
77*a5892487SPeter Brune 
78*a5892487SPeter Brune    Collective on SNES
79*a5892487SPeter Brune 
80*a5892487SPeter Brune    Input Parameters:
81*a5892487SPeter Brune +  snes - nonlinear context
82*a5892487SPeter Brune .  lsctx - optional context for line search (not used here)
83*a5892487SPeter Brune .  x - current iterate
84*a5892487SPeter Brune .  f - residual evaluated at x
85*a5892487SPeter Brune .  y - search direction
86*a5892487SPeter Brune .  fnorm - 2-norm of f
87*a5892487SPeter Brune -  xnorm - norm of x if known, otherwise 0
88*a5892487SPeter Brune 
89*a5892487SPeter Brune    Output Parameters:
90*a5892487SPeter Brune +  g - residual evaluated at new iterate y
91*a5892487SPeter Brune .  w - new iterate
92*a5892487SPeter Brune .  gnorm - 2-norm of g
93*a5892487SPeter Brune .  ynorm - 2-norm of search length
94*a5892487SPeter Brune -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
95*a5892487SPeter Brune 
96*a5892487SPeter Brune    Options Database Key:
97*a5892487SPeter Brune +  -snes_ls cubic - Activates SNESLineSearchCubic()
98*a5892487SPeter Brune .   -snes_ls_alpha <alpha> - Sets alpha
99*a5892487SPeter Brune .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
100*a5892487SPeter Brune -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
101*a5892487SPeter Brune 
102*a5892487SPeter Brune 
103*a5892487SPeter Brune    Notes:
104*a5892487SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
105*a5892487SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
106*a5892487SPeter Brune 
107*a5892487SPeter Brune    Level: advanced
108*a5892487SPeter Brune 
109*a5892487SPeter Brune .keywords: SNES, nonlinear, line search, cubic
110*a5892487SPeter Brune 
111*a5892487SPeter Brune .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
112*a5892487SPeter Brune @*/
113*a5892487SPeter Brune PetscErrorCode  SNESLineSearchCubic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
114*a5892487SPeter Brune {
115*a5892487SPeter Brune   /*
116*a5892487SPeter Brune      Note that for line search purposes we work with with the related
117*a5892487SPeter Brune      minimization problem:
118*a5892487SPeter Brune         min  z(x):  R^n -> R,
119*a5892487SPeter Brune      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
120*a5892487SPeter Brune    */
121*a5892487SPeter Brune 
122*a5892487SPeter Brune   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
123*a5892487SPeter Brune   PetscReal      minlambda,lambda,lambdatemp;
124*a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
125*a5892487SPeter Brune   PetscScalar    cinitslope;
126*a5892487SPeter Brune #endif
127*a5892487SPeter Brune   PetscErrorCode ierr;
128*a5892487SPeter Brune   PetscInt       count;
129*a5892487SPeter Brune   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
130*a5892487SPeter Brune   MPI_Comm       comm;
131*a5892487SPeter Brune 
132*a5892487SPeter Brune   PetscFunctionBegin;
133*a5892487SPeter Brune   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
134*a5892487SPeter Brune   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
135*a5892487SPeter Brune   *flag   = PETSC_TRUE;
136*a5892487SPeter Brune 
137*a5892487SPeter Brune   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
138*a5892487SPeter Brune   if (*ynorm == 0.0) {
139*a5892487SPeter Brune     if (snes->ls_monitor) {
140*a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
141*a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
142*a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
143*a5892487SPeter Brune     }
144*a5892487SPeter Brune     *gnorm = fnorm;
145*a5892487SPeter Brune     ierr   = VecCopy(x,w);CHKERRQ(ierr);
146*a5892487SPeter Brune     ierr   = VecCopy(f,g);CHKERRQ(ierr);
147*a5892487SPeter Brune     *flag  = PETSC_FALSE;
148*a5892487SPeter Brune     goto theend1;
149*a5892487SPeter Brune   }
150*a5892487SPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
151*a5892487SPeter Brune     if (snes->ls_monitor) {
152*a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
153*a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(snes->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
154*a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
155*a5892487SPeter Brune     }
156*a5892487SPeter Brune     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
157*a5892487SPeter Brune     *ynorm = snes->maxstep;
158*a5892487SPeter Brune   }
159*a5892487SPeter Brune   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
160*a5892487SPeter Brune   minlambda = snes->steptol/rellength;
161*a5892487SPeter Brune   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
162*a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
163*a5892487SPeter Brune   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
164*a5892487SPeter Brune   initslope = PetscRealPart(cinitslope);
165*a5892487SPeter Brune #else
166*a5892487SPeter Brune   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
167*a5892487SPeter Brune #endif
168*a5892487SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
169*a5892487SPeter Brune   if (initslope == 0.0) initslope = -1.0;
170*a5892487SPeter Brune 
171*a5892487SPeter Brune   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
172*a5892487SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
173*a5892487SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
174*a5892487SPeter Brune     *flag = PETSC_FALSE;
175*a5892487SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
176*a5892487SPeter Brune     goto theend1;
177*a5892487SPeter Brune   }
178*a5892487SPeter Brune   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
179*a5892487SPeter Brune   if (snes->domainerror) {
180*a5892487SPeter Brune     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
181*a5892487SPeter Brune     PetscFunctionReturn(0);
182*a5892487SPeter Brune   }
183*a5892487SPeter Brune   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
184*a5892487SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
185*a5892487SPeter Brune   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
186*a5892487SPeter Brune   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
187*a5892487SPeter Brune     if (snes->ls_monitor) {
188*a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
189*a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
190*a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
191*a5892487SPeter Brune     }
192*a5892487SPeter Brune     goto theend1;
193*a5892487SPeter Brune   }
194*a5892487SPeter Brune 
195*a5892487SPeter Brune   /* Fit points with quadratic */
196*a5892487SPeter Brune   lambda     = 1.0;
197*a5892487SPeter Brune   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
198*a5892487SPeter Brune   lambdaprev = lambda;
199*a5892487SPeter Brune   gnormprev  = *gnorm;
200*a5892487SPeter Brune   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
201*a5892487SPeter Brune   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
202*a5892487SPeter Brune   else                         lambda = lambdatemp;
203*a5892487SPeter Brune 
204*a5892487SPeter Brune   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
205*a5892487SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
206*a5892487SPeter Brune     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
207*a5892487SPeter Brune     *flag = PETSC_FALSE;
208*a5892487SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
209*a5892487SPeter Brune     goto theend1;
210*a5892487SPeter Brune   }
211*a5892487SPeter Brune   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
212*a5892487SPeter Brune   if (snes->domainerror) {
213*a5892487SPeter Brune     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
214*a5892487SPeter Brune     PetscFunctionReturn(0);
215*a5892487SPeter Brune   }
216*a5892487SPeter Brune   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
217*a5892487SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
218*a5892487SPeter Brune   if (snes->ls_monitor) {
219*a5892487SPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
220*a5892487SPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr);
221*a5892487SPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
222*a5892487SPeter Brune   }
223*a5892487SPeter Brune   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
224*a5892487SPeter Brune     if (snes->ls_monitor) {
225*a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
226*a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
227*a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
228*a5892487SPeter Brune     }
229*a5892487SPeter Brune     goto theend1;
230*a5892487SPeter Brune   }
231*a5892487SPeter Brune 
232*a5892487SPeter Brune   /* Fit points with cubic */
233*a5892487SPeter Brune   count = 1;
234*a5892487SPeter Brune   while (PETSC_TRUE) {
235*a5892487SPeter Brune     if (lambda <= minlambda) {
236*a5892487SPeter Brune       if (snes->ls_monitor) {
237*a5892487SPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
238*a5892487SPeter Brune 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
239*a5892487SPeter Brune 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr);
240*a5892487SPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
241*a5892487SPeter Brune       }
242*a5892487SPeter Brune       *flag = PETSC_FALSE;
243*a5892487SPeter Brune       break;
244*a5892487SPeter Brune     }
245*a5892487SPeter Brune     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
246*a5892487SPeter Brune     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
247*a5892487SPeter Brune     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
248*a5892487SPeter Brune     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
249*a5892487SPeter Brune     d  = b*b - 3*a*initslope;
250*a5892487SPeter Brune     if (d < 0.0) d = 0.0;
251*a5892487SPeter Brune     if (a == 0.0) {
252*a5892487SPeter Brune       lambdatemp = -initslope/(2.0*b);
253*a5892487SPeter Brune     } else {
254*a5892487SPeter Brune       lambdatemp = (-b + sqrt(d))/(3.0*a);
255*a5892487SPeter Brune     }
256*a5892487SPeter Brune     lambdaprev = lambda;
257*a5892487SPeter Brune     gnormprev  = *gnorm;
258*a5892487SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
259*a5892487SPeter Brune     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
260*a5892487SPeter Brune     else                         lambda     = lambdatemp;
261*a5892487SPeter Brune     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
262*a5892487SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
263*a5892487SPeter Brune       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
264*a5892487SPeter Brune       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
265*a5892487SPeter Brune       *flag = PETSC_FALSE;
266*a5892487SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
267*a5892487SPeter Brune       break;
268*a5892487SPeter Brune     }
269*a5892487SPeter Brune     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
270*a5892487SPeter Brune     if (snes->domainerror) {
271*a5892487SPeter Brune       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
272*a5892487SPeter Brune       PetscFunctionReturn(0);
273*a5892487SPeter Brune     }
274*a5892487SPeter Brune     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
275*a5892487SPeter Brune     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
276*a5892487SPeter Brune     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* is reduction enough? */
277*a5892487SPeter Brune       if (snes->ls_monitor) {
278*a5892487SPeter Brune 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
279*a5892487SPeter Brune       }
280*a5892487SPeter Brune       break;
281*a5892487SPeter Brune     } else {
282*a5892487SPeter Brune       if (snes->ls_monitor) {
283*a5892487SPeter Brune         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
284*a5892487SPeter Brune       }
285*a5892487SPeter Brune     }
286*a5892487SPeter Brune     count++;
287*a5892487SPeter Brune   }
288*a5892487SPeter Brune   theend1:
289*a5892487SPeter Brune   /* Optional user-defined check for line search step validity */
290*a5892487SPeter Brune   if (snes->ops->postcheckstep && *flag) {
291*a5892487SPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
292*a5892487SPeter Brune     if (changed_y) {
293*a5892487SPeter Brune       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
294*a5892487SPeter Brune     }
295*a5892487SPeter Brune     if (changed_y || changed_w) { /* recompute the function if the step has changed */
296*a5892487SPeter Brune       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
297*a5892487SPeter Brune       if (snes->domainerror) {
298*a5892487SPeter Brune         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
299*a5892487SPeter Brune         PetscFunctionReturn(0);
300*a5892487SPeter Brune       }
301*a5892487SPeter Brune       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
302*a5892487SPeter Brune       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
303*a5892487SPeter Brune       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
304*a5892487SPeter Brune       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
305*a5892487SPeter Brune       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
306*a5892487SPeter Brune     }
307*a5892487SPeter Brune   }
308*a5892487SPeter Brune   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
309*a5892487SPeter Brune   PetscFunctionReturn(0);
310*a5892487SPeter Brune }
311*a5892487SPeter Brune /* -------------------------------------------------------------------------- */
312*a5892487SPeter Brune 
313*a5892487SPeter Brune #undef __FUNCT__
314*a5892487SPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_LS"
315*a5892487SPeter Brune /*@C
316*a5892487SPeter Brune    SNESLineSearchQuadratic_LS - Performs a quadratic line search.
317*a5892487SPeter Brune 
318*a5892487SPeter Brune    Collective on SNES and Vec
319*a5892487SPeter Brune 
320*a5892487SPeter Brune    Input Parameters:
321*a5892487SPeter Brune +  snes - the SNES context
322*a5892487SPeter Brune .  lsctx - optional context for line search (not used here)
323*a5892487SPeter Brune .  x - current iterate
324*a5892487SPeter Brune .  f - residual evaluated at x
325*a5892487SPeter Brune .  y - search direction
326*a5892487SPeter Brune .  fnorm - 2-norm of f
327*a5892487SPeter Brune -  xnorm - norm of x if known, otherwise 0
328*a5892487SPeter Brune 
329*a5892487SPeter Brune    Output Parameters:
330*a5892487SPeter Brune +  g - residual evaluated at new iterate w
331*a5892487SPeter Brune .  w - new iterate (x + lambda*y)
332*a5892487SPeter Brune .  gnorm - 2-norm of g
333*a5892487SPeter Brune .  ynorm - 2-norm of search length
334*a5892487SPeter Brune -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
335*a5892487SPeter Brune 
336*a5892487SPeter Brune    Options Database Keys:
337*a5892487SPeter Brune +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
338*a5892487SPeter Brune .   -snes_ls_alpha <alpha> - Sets alpha
339*a5892487SPeter Brune .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
340*a5892487SPeter Brune -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
341*a5892487SPeter Brune 
342*a5892487SPeter Brune    Notes:
343*a5892487SPeter Brune    Use SNESLineSearchSet() to set this routine within the SNESLS method.
344*a5892487SPeter Brune 
345*a5892487SPeter Brune    Level: advanced
346*a5892487SPeter Brune 
347*a5892487SPeter Brune .keywords: SNES, nonlinear, quadratic, line search
348*a5892487SPeter Brune 
349*a5892487SPeter Brune .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
350*a5892487SPeter Brune @*/
351*a5892487SPeter Brune PetscErrorCode  SNESLineSearchQuadratic_LS(SNES snes,void *lsctx,Vec x,Vec f,Vec y,PetscReal fnorm,PetscReal xnorm,Vec g,Vec w,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
352*a5892487SPeter Brune {
353*a5892487SPeter Brune   /*
354*a5892487SPeter Brune      Note that for line search purposes we work with with the related
355*a5892487SPeter Brune      minimization problem:
356*a5892487SPeter Brune         min  z(x):  R^n -> R,
357*a5892487SPeter Brune      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
358*a5892487SPeter Brune    */
359*a5892487SPeter Brune   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
360*a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
361*a5892487SPeter Brune   PetscScalar    cinitslope;
362*a5892487SPeter Brune #endif
363*a5892487SPeter Brune   PetscErrorCode ierr;
364*a5892487SPeter Brune   PetscInt       count;
365*a5892487SPeter Brune   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
366*a5892487SPeter Brune 
367*a5892487SPeter Brune   PetscFunctionBegin;
368*a5892487SPeter Brune   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
369*a5892487SPeter Brune   *flag   = PETSC_TRUE;
370*a5892487SPeter Brune 
371*a5892487SPeter Brune   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
372*a5892487SPeter Brune   if (*ynorm == 0.0) {
373*a5892487SPeter Brune     if (snes->ls_monitor) {
374*a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
375*a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
376*a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
377*a5892487SPeter Brune     }
378*a5892487SPeter Brune     *gnorm = fnorm;
379*a5892487SPeter Brune     ierr   = VecCopy(x,w);CHKERRQ(ierr);
380*a5892487SPeter Brune     ierr   = VecCopy(f,g);CHKERRQ(ierr);
381*a5892487SPeter Brune     *flag  = PETSC_FALSE;
382*a5892487SPeter Brune     goto theend2;
383*a5892487SPeter Brune   }
384*a5892487SPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
385*a5892487SPeter Brune     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
386*a5892487SPeter Brune     *ynorm = snes->maxstep;
387*a5892487SPeter Brune   }
388*a5892487SPeter Brune   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
389*a5892487SPeter Brune   minlambda = snes->steptol/rellength;
390*a5892487SPeter Brune   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
391*a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
392*a5892487SPeter Brune   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
393*a5892487SPeter Brune   initslope = PetscRealPart(cinitslope);
394*a5892487SPeter Brune #else
395*a5892487SPeter Brune   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
396*a5892487SPeter Brune #endif
397*a5892487SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
398*a5892487SPeter Brune   if (initslope == 0.0) initslope = -1.0;
399*a5892487SPeter Brune   ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr);
400*a5892487SPeter Brune 
401*a5892487SPeter Brune   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
402*a5892487SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
403*a5892487SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
404*a5892487SPeter Brune     *flag = PETSC_FALSE;
405*a5892487SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
406*a5892487SPeter Brune     goto theend2;
407*a5892487SPeter Brune   }
408*a5892487SPeter Brune   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
409*a5892487SPeter Brune   if (snes->domainerror) {
410*a5892487SPeter Brune     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
411*a5892487SPeter Brune     PetscFunctionReturn(0);
412*a5892487SPeter Brune   }
413*a5892487SPeter Brune   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
414*a5892487SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
415*a5892487SPeter Brune   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
416*a5892487SPeter Brune     if (snes->ls_monitor) {
417*a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
418*a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
419*a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
420*a5892487SPeter Brune     }
421*a5892487SPeter Brune     goto theend2;
422*a5892487SPeter Brune   }
423*a5892487SPeter Brune 
424*a5892487SPeter Brune   /* Fit points with quadratic */
425*a5892487SPeter Brune   lambda = 1.0;
426*a5892487SPeter Brune   count = 1;
427*a5892487SPeter Brune   while (PETSC_TRUE) {
428*a5892487SPeter Brune     if (lambda <= minlambda) { /* bad luck; use full step */
429*a5892487SPeter Brune       if (snes->ls_monitor) {
430*a5892487SPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
431*a5892487SPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
432*a5892487SPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: fnorm=%14.12e, gnorm=%14.12e, ynorm=%14.12e, lambda=%14.12e, initial slope=%14.12e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
433*a5892487SPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
434*a5892487SPeter Brune       }
435*a5892487SPeter Brune       ierr = VecCopy(x,w);CHKERRQ(ierr);
436*a5892487SPeter Brune       *flag = PETSC_FALSE;
437*a5892487SPeter Brune       break;
438*a5892487SPeter Brune     }
439*a5892487SPeter Brune     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
440*a5892487SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
441*a5892487SPeter Brune     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
442*a5892487SPeter Brune     else                         lambda     = lambdatemp;
443*a5892487SPeter Brune 
444*a5892487SPeter Brune     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
445*a5892487SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
446*a5892487SPeter Brune       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
447*a5892487SPeter Brune       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
448*a5892487SPeter Brune       *flag = PETSC_FALSE;
449*a5892487SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
450*a5892487SPeter Brune       break;
451*a5892487SPeter Brune     }
452*a5892487SPeter Brune     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
453*a5892487SPeter Brune     if (snes->domainerror) {
454*a5892487SPeter Brune       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
455*a5892487SPeter Brune       PetscFunctionReturn(0);
456*a5892487SPeter Brune     }
457*a5892487SPeter Brune     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
458*a5892487SPeter Brune     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
459*a5892487SPeter Brune     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
460*a5892487SPeter Brune       if (snes->ls_monitor) {
461*a5892487SPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
462*a5892487SPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr);
463*a5892487SPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
464*a5892487SPeter Brune       }
465*a5892487SPeter Brune       break;
466*a5892487SPeter Brune     }
467*a5892487SPeter Brune     count++;
468*a5892487SPeter Brune   }
469*a5892487SPeter Brune   theend2:
470*a5892487SPeter Brune   /* Optional user-defined check for line search step validity */
471*a5892487SPeter Brune   if (snes->ops->postcheckstep) {
472*a5892487SPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
473*a5892487SPeter Brune     if (changed_y) {
474*a5892487SPeter Brune       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
475*a5892487SPeter Brune     }
476*a5892487SPeter Brune     if (changed_y || changed_w) { /* recompute the function if the step has changed */
477*a5892487SPeter Brune       ierr = SNESComputeFunction(snes,w,g);
478*a5892487SPeter Brune       if (snes->domainerror) {
479*a5892487SPeter Brune         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
480*a5892487SPeter Brune         PetscFunctionReturn(0);
481*a5892487SPeter Brune       }
482*a5892487SPeter Brune       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
483*a5892487SPeter Brune       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
484*a5892487SPeter Brune       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
485*a5892487SPeter Brune       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
486*a5892487SPeter Brune       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
487*a5892487SPeter Brune     }
488*a5892487SPeter Brune   }
489*a5892487SPeter Brune   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
490*a5892487SPeter Brune   PetscFunctionReturn(0);
491*a5892487SPeter Brune }
492*a5892487SPeter Brune 
493*a5892487SPeter Brune 
49404d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
49504d965bbSLois Curfman McInnes 
49604d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
49794b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
49804d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
49904d965bbSLois Curfman McInnes      respectively.
50004d965bbSLois Curfman McInnes 
501fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
50204d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
50304d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
504fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
50504d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
50604d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
5074b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
5084b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
50904d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
51004d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
51104d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
51204d965bbSLois Curfman McInnes      for all nonlinear solvers.
51304d965bbSLois Curfman McInnes 
51404d965bbSLois Curfman McInnes      Another key routine is:
51504d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
51604d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
51704d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
51804d965bbSLois Curfman McInnes      SNESSolve() if necessary.
51904d965bbSLois Curfman McInnes 
52004d965bbSLois Curfman McInnes      Additional basic routines are:
52104d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
52204d965bbSLois Curfman McInnes                                       have actually been used.
52304d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
524186905e3SBarry Smith      SNESView().
52504d965bbSLois Curfman McInnes 
52604d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
52704d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
52804d965bbSLois Curfman McInnes      above description applies to these categories also.
52904d965bbSLois Curfman McInnes 
53004d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
5315e42470aSBarry Smith /*
5324b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
53304d965bbSLois Curfman McInnes    method with a line search.
5345e76c431SBarry Smith 
53504d965bbSLois Curfman McInnes    Input Parameters:
53604d965bbSLois Curfman McInnes .  snes - the SNES context
5375e76c431SBarry Smith 
53804d965bbSLois Curfman McInnes    Output Parameter:
539c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
54004d965bbSLois Curfman McInnes 
54104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
5425e76c431SBarry Smith 
5435e76c431SBarry Smith    Notes:
5445e76c431SBarry Smith    This implements essentially a truncated Newton method with a
5455e76c431SBarry Smith    line search.  By default a cubic backtracking line search
5465e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
5475e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
548393d2d9aSLois Curfman McInnes    and Schnabel.
5495e76c431SBarry Smith */
5504a2ae208SSatish Balay #undef __FUNCT__
5514b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
552dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
5535e76c431SBarry Smith {
5546849ba73SBarry Smith   PetscErrorCode     ierr;
555a7cc72afSBarry Smith   PetscInt           maxits,i,lits;
556ace3abfcSBarry Smith   PetscBool          lssucceed;
557112a2221SBarry Smith   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
55885385478SLisandro Dalcin   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
55985385478SLisandro Dalcin   Vec                Y,X,F,G,W;
5603d4c4710SBarry Smith   KSPConvergedReason kspreason;
5615e76c431SBarry Smith 
5623a40ed3dSBarry Smith   PetscFunctionBegin;
5633d4c4710SBarry Smith   snes->numFailures            = 0;
5643d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
565184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
566184914b5SBarry Smith 
5675e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
5685e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
56939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
5705e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
5715e42470aSBarry Smith   G		= snes->work[1];
5725e42470aSBarry Smith   W		= snes->work[2];
5735e76c431SBarry Smith 
5744c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5754c49b128SBarry Smith   snes->iter = 0;
57675567043SBarry Smith   snes->norm = 0.0;
5774c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
57819717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
5794936397dSBarry Smith   if (snes->domainerror) {
5804936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
5814936397dSBarry Smith     PetscFunctionReturn(0);
5824936397dSBarry Smith   }
5832613ca53SBarry Smith   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
5842613ca53SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
5852613ca53SBarry Smith   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
5862613ca53SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
58762cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5880f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
5895e42470aSBarry Smith   snes->norm = fnorm;
5900f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
59184c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
5927a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
5933f149594SLisandro Dalcin 
5943f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
5953f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
59685385478SLisandro Dalcin   /* test convergence */
59785385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
59806ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
599d034289bSBarry Smith 
6005e76c431SBarry Smith   for (i=0; i<maxits; i++) {
6015e76c431SBarry Smith 
602329e7e40SMatthew Knepley     /* Call general purpose update function */
603e7788613SBarry Smith     if (snes->ops->update) {
604e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
605de8cb200SMatthew Knepley     }
606329e7e40SMatthew Knepley 
607ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
6085334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
60994b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
61071f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
6113d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
6123d4c4710SBarry Smith     if (kspreason < 0) {
6133d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
614ef998cc9SBarry Smith         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
6153d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
616ab81109fSSatish Balay         break;
6173d4c4710SBarry Smith       }
6183d4c4710SBarry Smith     }
6193d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
6203f149594SLisandro Dalcin     snes->linear_its += lits;
6213f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
62274637425SBarry Smith 
623ea630c6eSPeter Brune     if (snes->ops->precheckstep) {
624ace3abfcSBarry Smith       PetscBool  changed_y = PETSC_FALSE;
625ea630c6eSPeter Brune       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
6269c3ca13aSBarry Smith     }
6279c3ca13aSBarry Smith 
628b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
6291e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
63085471664SBarry Smith     }
63174637425SBarry Smith 
632ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
633ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
634e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
635ea4d3ed3SLois Curfman McInnes     */
63685385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
6373f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
638ea630c6eSPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
6398f1a2a5eSBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
6404a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
6414936397dSBarry Smith     if (snes->domainerror) {
6424936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
6434936397dSBarry Smith       PetscFunctionReturn(0);
6444936397dSBarry Smith     }
645a7cc72afSBarry Smith     if (!lssucceed) {
64650ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
647ace3abfcSBarry Smith 	PetscBool  ismin;
648647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6491e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
6508a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
6513505fcc1SBarry Smith         break;
6523505fcc1SBarry Smith       }
65350ffb88aSMatthew Knepley     }
65485385478SLisandro Dalcin     /* Update function and solution vectors */
65585385478SLisandro Dalcin     fnorm = gnorm;
65685385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
65785385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
65885385478SLisandro Dalcin     /* Monitor convergence */
65985385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
66085385478SLisandro Dalcin     snes->iter = i+1;
66185385478SLisandro Dalcin     snes->norm = fnorm;
66285385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
66385385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
6647a03ce2fSLisandro Dalcin     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
66585385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
66685385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
667e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
6683f149594SLisandro Dalcin     if (snes->reason) break;
66916c95cacSBarry Smith   }
67052392280SLois Curfman McInnes   if (i == maxits) {
671ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
67285385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
67352392280SLois Curfman McInnes   }
6743a40ed3dSBarry Smith   PetscFunctionReturn(0);
6755e76c431SBarry Smith }
67604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
67704d965bbSLois Curfman McInnes /*
6784b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
6794b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
68004d965bbSLois Curfman McInnes 
68104d965bbSLois Curfman McInnes    Input Parameter:
68204d965bbSLois Curfman McInnes .  snes - the SNES context
68304d965bbSLois Curfman McInnes .  x - the solution vector
68404d965bbSLois Curfman McInnes 
68504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
68604d965bbSLois Curfman McInnes 
68704d965bbSLois Curfman McInnes    Notes:
6884b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
68904d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
69004d965bbSLois Curfman McInnes    the call to SNESSolve().
69104d965bbSLois Curfman McInnes  */
6924a2ae208SSatish Balay #undef __FUNCT__
6934b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
694dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
6955e76c431SBarry Smith {
696dfbe8321SBarry Smith   PetscErrorCode ierr;
6973a40ed3dSBarry Smith 
6983a40ed3dSBarry Smith   PetscFunctionBegin;
69958c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
7003a40ed3dSBarry Smith   PetscFunctionReturn(0);
7015e76c431SBarry Smith }
70204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7036b8b9a38SLisandro Dalcin 
7046b8b9a38SLisandro Dalcin #undef __FUNCT__
7056b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS"
7066b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes)
7076b8b9a38SLisandro Dalcin {
7086b8b9a38SLisandro Dalcin   PetscErrorCode ierr;
7096b8b9a38SLisandro Dalcin 
7106b8b9a38SLisandro Dalcin   PetscFunctionBegin;
7116b8b9a38SLisandro Dalcin   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
7126b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
7136b8b9a38SLisandro Dalcin }
7146b8b9a38SLisandro Dalcin 
71504d965bbSLois Curfman McInnes /*
7164b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
7174b27c08aSLois Curfman McInnes    with SNESCreate_LS().
71804d965bbSLois Curfman McInnes 
71904d965bbSLois Curfman McInnes    Input Parameter:
72004d965bbSLois Curfman McInnes .  snes - the SNES context
72104d965bbSLois Curfman McInnes 
722de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
72304d965bbSLois Curfman McInnes  */
7244a2ae208SSatish Balay #undef __FUNCT__
7254b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
726dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
7275e76c431SBarry Smith {
728dfbe8321SBarry Smith   PetscErrorCode ierr;
7293a40ed3dSBarry Smith 
7303a40ed3dSBarry Smith   PetscFunctionBegin;
7316b8b9a38SLisandro Dalcin   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
7325c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
733557d3f75SLisandro Dalcin 
7343a40ed3dSBarry Smith   PetscFunctionReturn(0);
7355e76c431SBarry Smith }
73604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
73794298653SBarry Smith 
738329e7e40SMatthew Knepley /*
7394b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
74004d965bbSLois Curfman McInnes 
74104d965bbSLois Curfman McInnes    Input Parameters:
74204d965bbSLois Curfman McInnes .  SNES - the SNES context
74304d965bbSLois Curfman McInnes .  viewer - visualization context
74404d965bbSLois Curfman McInnes 
74504d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
74604d965bbSLois Curfman McInnes */
7474a2ae208SSatish Balay #undef __FUNCT__
7484b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
7496849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
750a935fc98SLois Curfman McInnes {
7512fc52814SBarry Smith   const char     *cstr;
752dfbe8321SBarry Smith   PetscErrorCode ierr;
753ace3abfcSBarry Smith   PetscBool      iascii;
754a935fc98SLois Curfman McInnes 
7553a40ed3dSBarry Smith   PetscFunctionBegin;
7562692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
75732077d6dSBarry Smith   if (iascii) {
75815f5eeeaSPeter Brune     cstr = SNESLineSearchTypeName(snes->ls_type);
759b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
760ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)snes->ls_alpha,(double)snes->maxstep,(double)snes->steptol);CHKERRQ(ierr);
761ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  damping factor=%14.12e\n",(double)snes->damping);CHKERRQ(ierr);
76219bcc07fSBarry Smith   }
7633a40ed3dSBarry Smith   PetscFunctionReturn(0);
764a935fc98SLois Curfman McInnes }
76504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
76604d965bbSLois Curfman McInnes /*
7674b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
76804d965bbSLois Curfman McInnes 
76904d965bbSLois Curfman McInnes    Input Parameter:
77004d965bbSLois Curfman McInnes .  snes - the SNES context
77104d965bbSLois Curfman McInnes 
77204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
77304d965bbSLois Curfman McInnes */
7744a2ae208SSatish Balay #undef __FUNCT__
7754b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
7766849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
7775e42470aSBarry Smith {
778dfbe8321SBarry Smith   PetscErrorCode     ierr;
7795e42470aSBarry Smith 
7803a40ed3dSBarry Smith   PetscFunctionBegin;
781b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
782b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7833a40ed3dSBarry Smith   PetscFunctionReturn(0);
7845e42470aSBarry Smith }
7854827ddcaSBarry Smith 
7864827ddcaSBarry Smith EXTERN_C_BEGIN
7874827ddcaSBarry Smith extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
7884827ddcaSBarry Smith EXTERN_C_END
7894827ddcaSBarry Smith 
79004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
791ebe3b25bSBarry Smith /*MC
792ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
79304d965bbSLois Curfman McInnes 
794ebe3b25bSBarry Smith    Options Database:
795ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
79655d4ea13SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
797e106eecfSBarry Smith .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
798acbee50cSBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
79955d4ea13SBarry Smith .   -snes_ls_monitor - print information about progress of line searches
80055d4ea13SBarry Smith -   -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms
801acbee50cSBarry Smith 
80204d965bbSLois Curfman McInnes 
803ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
80404d965bbSLois Curfman McInnes 
805ee3001cbSBarry Smith    Level: beginner
806ee3001cbSBarry Smith 
80717bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
80817bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
809b3dd4ab5SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
810ebe3b25bSBarry Smith 
811ebe3b25bSBarry Smith M*/
812fb2e594dSBarry Smith EXTERN_C_BEGIN
8134a2ae208SSatish Balay #undef __FUNCT__
8144b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
8157087cfbeSBarry Smith PetscErrorCode  SNESCreate_LS(SNES snes)
8165e42470aSBarry Smith {
817dfbe8321SBarry Smith   PetscErrorCode ierr;
8184b27c08aSLois Curfman McInnes   SNES_LS        *neP;
8195e42470aSBarry Smith 
8203a40ed3dSBarry Smith   PetscFunctionBegin;
821e7788613SBarry Smith   snes->ops->setup           = SNESSetUp_LS;
822e7788613SBarry Smith   snes->ops->solve           = SNESSolve_LS;
823e7788613SBarry Smith   snes->ops->destroy         = SNESDestroy_LS;
824e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
825e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
8266b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_LS;
8275e42470aSBarry Smith 
82842f4f86dSBarry Smith   snes->usesksp                      = PETSC_TRUE;
82942f4f86dSBarry Smith   snes->usespc                       = PETSC_FALSE;
83038f2d2fdSLisandro Dalcin   ierr                               = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
8315e42470aSBarry Smith   snes->data                         = (void*)neP;
83215f5eeeaSPeter Brune   snes->ops->linesearchno            = SNESLineSearchNo;
83315f5eeeaSPeter Brune   snes->ops->linesearchnonorms       = SNESLineSearchNoNorms;
834*a5892487SPeter Brune   snes->ops->linesearchquadratic     = SNESLineSearchQuadratic_LS;
835*a5892487SPeter Brune   snes->ops->linesearchcubic         = SNESLineSearchCubic_LS;
836ea630c6eSPeter Brune   snes->ops->linesearchexact         = PETSC_NULL;
837ea630c6eSPeter Brune   snes->ops->linesearchtest          = PETSC_NULL;
838ea630c6eSPeter Brune   snes->lsP                          = PETSC_NULL;
839ea630c6eSPeter Brune   snes->ops->postcheckstep           = PETSC_NULL;
840ea630c6eSPeter Brune   snes->postcheck                    = PETSC_NULL;
841ea630c6eSPeter Brune   snes->ops->precheckstep            = PETSC_NULL;
842ea630c6eSPeter Brune   snes->precheck                     = PETSC_NULL;
84315f5eeeaSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr);
8443a40ed3dSBarry Smith   PetscFunctionReturn(0);
8455e42470aSBarry Smith }
846fb2e594dSBarry Smith EXTERN_C_END
847