xref: /petsc/src/snes/impls/ls/ls.c (revision 70d3b23b377b60e976d343365f7c92b8b2f0f13f)
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 
73a5892487SPeter Brune #undef __FUNCT__
74a5892487SPeter Brune #define __FUNCT__ "SNESLineSearchCubic_LS"
75a5892487SPeter Brune /*@C
76a5892487SPeter Brune    SNESLineSearchCubic - Performs a cubic line search (default line search method).
77a5892487SPeter Brune 
78a5892487SPeter Brune    Collective on SNES
79a5892487SPeter Brune 
80a5892487SPeter Brune    Input Parameters:
81a5892487SPeter Brune +  snes - nonlinear context
82a5892487SPeter Brune .  lsctx - optional context for line search (not used here)
83a5892487SPeter Brune .  x - current iterate
84a5892487SPeter Brune .  f - residual evaluated at x
85a5892487SPeter Brune .  y - search direction
86a5892487SPeter Brune .  fnorm - 2-norm of f
87a5892487SPeter Brune -  xnorm - norm of x if known, otherwise 0
88a5892487SPeter Brune 
89a5892487SPeter Brune    Output Parameters:
90a5892487SPeter Brune +  g - residual evaluated at new iterate y
91a5892487SPeter Brune .  w - new iterate
92a5892487SPeter Brune .  gnorm - 2-norm of g
93a5892487SPeter Brune .  ynorm - 2-norm of search length
94a5892487SPeter Brune -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
95a5892487SPeter Brune 
96a5892487SPeter Brune    Options Database Key:
97a5892487SPeter Brune +  -snes_ls cubic - Activates SNESLineSearchCubic()
98a5892487SPeter Brune .   -snes_ls_alpha <alpha> - Sets alpha
99a5892487SPeter 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)
100a5892487SPeter Brune -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
101a5892487SPeter Brune 
102a5892487SPeter Brune 
103a5892487SPeter Brune    Notes:
104a5892487SPeter Brune    This line search is taken from "Numerical Methods for Unconstrained
105a5892487SPeter Brune    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
106a5892487SPeter Brune 
107a5892487SPeter Brune    Level: advanced
108a5892487SPeter Brune 
109a5892487SPeter Brune .keywords: SNES, nonlinear, line search, cubic
110a5892487SPeter Brune 
111a5892487SPeter Brune .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
112a5892487SPeter Brune @*/
113a5892487SPeter 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)
114a5892487SPeter Brune {
115a5892487SPeter Brune   /*
116a5892487SPeter Brune      Note that for line search purposes we work with with the related
117a5892487SPeter Brune      minimization problem:
118a5892487SPeter Brune         min  z(x):  R^n -> R,
119a5892487SPeter Brune      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
120a5892487SPeter Brune    */
121a5892487SPeter Brune 
122a5892487SPeter Brune   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
123a5892487SPeter Brune   PetscReal      minlambda,lambda,lambdatemp;
124a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
125a5892487SPeter Brune   PetscScalar    cinitslope;
126a5892487SPeter Brune #endif
127a5892487SPeter Brune   PetscErrorCode ierr;
128a5892487SPeter Brune   PetscInt       count;
129a5892487SPeter Brune   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
130a5892487SPeter Brune   MPI_Comm       comm;
131a5892487SPeter Brune 
132a5892487SPeter Brune   PetscFunctionBegin;
133a5892487SPeter Brune   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
134a5892487SPeter Brune   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
135a5892487SPeter Brune   *flag   = PETSC_TRUE;
136a5892487SPeter Brune 
137a5892487SPeter Brune   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
138a5892487SPeter Brune   if (*ynorm == 0.0) {
139a5892487SPeter Brune     if (snes->ls_monitor) {
140a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
141a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
142a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
143a5892487SPeter Brune     }
144a5892487SPeter Brune     *gnorm = fnorm;
145a5892487SPeter Brune     ierr   = VecCopy(x,w);CHKERRQ(ierr);
146a5892487SPeter Brune     ierr   = VecCopy(f,g);CHKERRQ(ierr);
147a5892487SPeter Brune     *flag  = PETSC_FALSE;
148a5892487SPeter Brune     goto theend1;
149a5892487SPeter Brune   }
150a5892487SPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
151a5892487SPeter Brune     if (snes->ls_monitor) {
152a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
153a5892487SPeter 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);
154a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
155a5892487SPeter Brune     }
156a5892487SPeter Brune     ierr = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
157a5892487SPeter Brune     *ynorm = snes->maxstep;
158a5892487SPeter Brune   }
159a5892487SPeter Brune   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
160a5892487SPeter Brune   minlambda = snes->steptol/rellength;
161a5892487SPeter Brune   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
162a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
163a5892487SPeter Brune   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
164a5892487SPeter Brune   initslope = PetscRealPart(cinitslope);
165a5892487SPeter Brune #else
166a5892487SPeter Brune   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
167a5892487SPeter Brune #endif
168a5892487SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
169a5892487SPeter Brune   if (initslope == 0.0) initslope = -1.0;
170a5892487SPeter Brune 
171a5892487SPeter Brune   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
172a5892487SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
173a5892487SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
174a5892487SPeter Brune     *flag = PETSC_FALSE;
175a5892487SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
176a5892487SPeter Brune     goto theend1;
177a5892487SPeter Brune   }
178a5892487SPeter Brune   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
179a5892487SPeter Brune   if (snes->domainerror) {
180a5892487SPeter Brune     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
181a5892487SPeter Brune     PetscFunctionReturn(0);
182a5892487SPeter Brune   }
183a5892487SPeter Brune   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
184a5892487SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
185a5892487SPeter Brune   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
186a5892487SPeter Brune   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
187a5892487SPeter Brune     if (snes->ls_monitor) {
188a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
189a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
190a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
191a5892487SPeter Brune     }
192a5892487SPeter Brune     goto theend1;
193a5892487SPeter Brune   }
194a5892487SPeter Brune 
195a5892487SPeter Brune   /* Fit points with quadratic */
196a5892487SPeter Brune   lambda     = 1.0;
197a5892487SPeter Brune   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
198a5892487SPeter Brune   lambdaprev = lambda;
199a5892487SPeter Brune   gnormprev  = *gnorm;
200a5892487SPeter Brune   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
201a5892487SPeter Brune   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
202a5892487SPeter Brune   else                         lambda = lambdatemp;
203a5892487SPeter Brune 
204a5892487SPeter Brune   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
205a5892487SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
206a5892487SPeter Brune     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
207a5892487SPeter Brune     *flag = PETSC_FALSE;
208a5892487SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
209a5892487SPeter Brune     goto theend1;
210a5892487SPeter Brune   }
211a5892487SPeter Brune   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
212a5892487SPeter Brune   if (snes->domainerror) {
213a5892487SPeter Brune     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
214a5892487SPeter Brune     PetscFunctionReturn(0);
215a5892487SPeter Brune   }
216a5892487SPeter Brune   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
217a5892487SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
218a5892487SPeter Brune   if (snes->ls_monitor) {
219a5892487SPeter Brune     ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
220a5892487SPeter Brune     ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr);
221a5892487SPeter Brune     ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
222a5892487SPeter Brune   }
223a5892487SPeter Brune   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
224a5892487SPeter Brune     if (snes->ls_monitor) {
225a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
226a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
227a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
228a5892487SPeter Brune     }
229a5892487SPeter Brune     goto theend1;
230a5892487SPeter Brune   }
231a5892487SPeter Brune 
232a5892487SPeter Brune   /* Fit points with cubic */
233a5892487SPeter Brune   count = 1;
234a5892487SPeter Brune   while (PETSC_TRUE) {
235a5892487SPeter Brune     if (lambda <= minlambda) {
236a5892487SPeter Brune       if (snes->ls_monitor) {
237a5892487SPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
238a5892487SPeter Brune 	ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
239a5892487SPeter 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);
240a5892487SPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
241a5892487SPeter Brune       }
242a5892487SPeter Brune       *flag = PETSC_FALSE;
243a5892487SPeter Brune       break;
244a5892487SPeter Brune     }
245a5892487SPeter Brune     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
246a5892487SPeter Brune     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
247a5892487SPeter Brune     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
248a5892487SPeter Brune     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
249a5892487SPeter Brune     d  = b*b - 3*a*initslope;
250a5892487SPeter Brune     if (d < 0.0) d = 0.0;
251a5892487SPeter Brune     if (a == 0.0) {
252a5892487SPeter Brune       lambdatemp = -initslope/(2.0*b);
253a5892487SPeter Brune     } else {
254a5892487SPeter Brune       lambdatemp = (-b + sqrt(d))/(3.0*a);
255a5892487SPeter Brune     }
256a5892487SPeter Brune     lambdaprev = lambda;
257a5892487SPeter Brune     gnormprev  = *gnorm;
258a5892487SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
259a5892487SPeter Brune     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
260a5892487SPeter Brune     else                         lambda     = lambdatemp;
261a5892487SPeter Brune     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
262a5892487SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
263a5892487SPeter Brune       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
264a5892487SPeter 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);
265a5892487SPeter Brune       *flag = PETSC_FALSE;
266a5892487SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
267a5892487SPeter Brune       break;
268a5892487SPeter Brune     }
269a5892487SPeter Brune     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
270a5892487SPeter Brune     if (snes->domainerror) {
271a5892487SPeter Brune       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
272a5892487SPeter Brune       PetscFunctionReturn(0);
273a5892487SPeter Brune     }
274a5892487SPeter Brune     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
275a5892487SPeter Brune     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
276a5892487SPeter Brune     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* is reduction enough? */
277a5892487SPeter Brune       if (snes->ls_monitor) {
278a5892487SPeter Brune 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
279a5892487SPeter Brune       }
280a5892487SPeter Brune       break;
281a5892487SPeter Brune     } else {
282a5892487SPeter Brune       if (snes->ls_monitor) {
283a5892487SPeter 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);
284a5892487SPeter Brune       }
285a5892487SPeter Brune     }
286a5892487SPeter Brune     count++;
287a5892487SPeter Brune   }
288a5892487SPeter Brune   theend1:
289a5892487SPeter Brune   /* Optional user-defined check for line search step validity */
290a5892487SPeter Brune   if (snes->ops->postcheckstep && *flag) {
291a5892487SPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
292a5892487SPeter Brune     if (changed_y) {
293a5892487SPeter Brune       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
294a5892487SPeter Brune     }
295a5892487SPeter Brune     if (changed_y || changed_w) { /* recompute the function if the step has changed */
296a5892487SPeter Brune       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
297a5892487SPeter Brune       if (snes->domainerror) {
298a5892487SPeter Brune         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
299a5892487SPeter Brune         PetscFunctionReturn(0);
300a5892487SPeter Brune       }
301a5892487SPeter Brune       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
302a5892487SPeter Brune       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
303a5892487SPeter Brune       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
304a5892487SPeter Brune       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
305a5892487SPeter Brune       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
306a5892487SPeter Brune     }
307a5892487SPeter Brune   }
308a5892487SPeter Brune   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
309a5892487SPeter Brune   PetscFunctionReturn(0);
310a5892487SPeter Brune }
311a5892487SPeter Brune /* -------------------------------------------------------------------------- */
312a5892487SPeter Brune 
313a5892487SPeter Brune #undef __FUNCT__
314a5892487SPeter Brune #define __FUNCT__ "SNESLineSearchQuadratic_LS"
315a5892487SPeter Brune /*@C
316a5892487SPeter Brune    SNESLineSearchQuadratic_LS - Performs a quadratic line search.
317a5892487SPeter Brune 
318a5892487SPeter Brune    Collective on SNES and Vec
319a5892487SPeter Brune 
320a5892487SPeter Brune    Input Parameters:
321a5892487SPeter Brune +  snes - the SNES context
322a5892487SPeter Brune .  lsctx - optional context for line search (not used here)
323a5892487SPeter Brune .  x - current iterate
324a5892487SPeter Brune .  f - residual evaluated at x
325a5892487SPeter Brune .  y - search direction
326a5892487SPeter Brune .  fnorm - 2-norm of f
327a5892487SPeter Brune -  xnorm - norm of x if known, otherwise 0
328a5892487SPeter Brune 
329a5892487SPeter Brune    Output Parameters:
330a5892487SPeter Brune +  g - residual evaluated at new iterate w
331a5892487SPeter Brune .  w - new iterate (x + lambda*y)
332a5892487SPeter Brune .  gnorm - 2-norm of g
333a5892487SPeter Brune .  ynorm - 2-norm of search length
334a5892487SPeter Brune -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
335a5892487SPeter Brune 
336a5892487SPeter Brune    Options Database Keys:
337a5892487SPeter Brune +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
338a5892487SPeter Brune .   -snes_ls_alpha <alpha> - Sets alpha
339a5892487SPeter 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)
340a5892487SPeter Brune -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
341a5892487SPeter Brune 
342a5892487SPeter Brune    Notes:
343a5892487SPeter Brune    Use SNESLineSearchSet() to set this routine within the SNESLS method.
344a5892487SPeter Brune 
345a5892487SPeter Brune    Level: advanced
346a5892487SPeter Brune 
347a5892487SPeter Brune .keywords: SNES, nonlinear, quadratic, line search
348a5892487SPeter Brune 
349a5892487SPeter Brune .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
350a5892487SPeter Brune @*/
351a5892487SPeter 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)
352a5892487SPeter Brune {
353a5892487SPeter Brune   /*
354a5892487SPeter Brune      Note that for line search purposes we work with with the related
355a5892487SPeter Brune      minimization problem:
356a5892487SPeter Brune         min  z(x):  R^n -> R,
357a5892487SPeter Brune      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
358a5892487SPeter Brune    */
359a5892487SPeter Brune   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
360a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
361a5892487SPeter Brune   PetscScalar    cinitslope;
362a5892487SPeter Brune #endif
363a5892487SPeter Brune   PetscErrorCode ierr;
364a5892487SPeter Brune   PetscInt       count;
365a5892487SPeter Brune   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
366a5892487SPeter Brune 
367a5892487SPeter Brune   PetscFunctionBegin;
368a5892487SPeter Brune   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
369a5892487SPeter Brune   *flag   = PETSC_TRUE;
370a5892487SPeter Brune 
371a5892487SPeter Brune   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
372a5892487SPeter Brune   if (*ynorm == 0.0) {
373a5892487SPeter Brune     if (snes->ls_monitor) {
374a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
375a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
376a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
377a5892487SPeter Brune     }
378a5892487SPeter Brune     *gnorm = fnorm;
379a5892487SPeter Brune     ierr   = VecCopy(x,w);CHKERRQ(ierr);
380a5892487SPeter Brune     ierr   = VecCopy(f,g);CHKERRQ(ierr);
381a5892487SPeter Brune     *flag  = PETSC_FALSE;
382a5892487SPeter Brune     goto theend2;
383a5892487SPeter Brune   }
384a5892487SPeter Brune   if (*ynorm > snes->maxstep) {	/* Step too big, so scale back */
385a5892487SPeter Brune     ierr   = VecScale(y,snes->maxstep/(*ynorm));CHKERRQ(ierr);
386a5892487SPeter Brune     *ynorm = snes->maxstep;
387a5892487SPeter Brune   }
388a5892487SPeter Brune   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
389a5892487SPeter Brune   minlambda = snes->steptol/rellength;
390a5892487SPeter Brune   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
391a5892487SPeter Brune #if defined(PETSC_USE_COMPLEX)
392a5892487SPeter Brune   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
393a5892487SPeter Brune   initslope = PetscRealPart(cinitslope);
394a5892487SPeter Brune #else
395a5892487SPeter Brune   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
396a5892487SPeter Brune #endif
397a5892487SPeter Brune   if (initslope > 0.0)  initslope = -initslope;
398a5892487SPeter Brune   if (initslope == 0.0) initslope = -1.0;
399a5892487SPeter Brune   ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr);
400a5892487SPeter Brune 
401a5892487SPeter Brune   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
402a5892487SPeter Brune   if (snes->nfuncs >= snes->max_funcs) {
403a5892487SPeter Brune     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
404a5892487SPeter Brune     *flag = PETSC_FALSE;
405a5892487SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
406a5892487SPeter Brune     goto theend2;
407a5892487SPeter Brune   }
408a5892487SPeter Brune   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
409a5892487SPeter Brune   if (snes->domainerror) {
410a5892487SPeter Brune     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
411a5892487SPeter Brune     PetscFunctionReturn(0);
412a5892487SPeter Brune   }
413a5892487SPeter Brune   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
414a5892487SPeter Brune   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
415a5892487SPeter Brune   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + snes->ls_alpha*initslope) { /* Sufficient reduction */
416a5892487SPeter Brune     if (snes->ls_monitor) {
417a5892487SPeter Brune       ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
418a5892487SPeter Brune       ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
419a5892487SPeter Brune       ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
420a5892487SPeter Brune     }
421a5892487SPeter Brune     goto theend2;
422a5892487SPeter Brune   }
423a5892487SPeter Brune 
424a5892487SPeter Brune   /* Fit points with quadratic */
425a5892487SPeter Brune   lambda = 1.0;
426a5892487SPeter Brune   count = 1;
427a5892487SPeter Brune   while (PETSC_TRUE) {
428a5892487SPeter Brune     if (lambda <= minlambda) { /* bad luck; use full step */
429a5892487SPeter Brune       if (snes->ls_monitor) {
430a5892487SPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
431a5892487SPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
432a5892487SPeter 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);
433a5892487SPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
434a5892487SPeter Brune       }
435a5892487SPeter Brune       ierr = VecCopy(x,w);CHKERRQ(ierr);
436a5892487SPeter Brune       *flag = PETSC_FALSE;
437a5892487SPeter Brune       break;
438a5892487SPeter Brune     }
439a5892487SPeter Brune     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
440a5892487SPeter Brune     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
441a5892487SPeter Brune     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
442a5892487SPeter Brune     else                         lambda     = lambdatemp;
443a5892487SPeter Brune 
444a5892487SPeter Brune     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
445a5892487SPeter Brune     if (snes->nfuncs >= snes->max_funcs) {
446a5892487SPeter Brune       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
447a5892487SPeter 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);
448a5892487SPeter Brune       *flag = PETSC_FALSE;
449a5892487SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
450a5892487SPeter Brune       break;
451a5892487SPeter Brune     }
452a5892487SPeter Brune     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
453a5892487SPeter Brune     if (snes->domainerror) {
454a5892487SPeter Brune       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
455a5892487SPeter Brune       PetscFunctionReturn(0);
456a5892487SPeter Brune     }
457a5892487SPeter Brune     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
458a5892487SPeter Brune     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
459a5892487SPeter Brune     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*snes->ls_alpha*initslope) { /* sufficient reduction */
460a5892487SPeter Brune       if (snes->ls_monitor) {
461a5892487SPeter Brune         ierr = PetscViewerASCIIAddTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
462a5892487SPeter Brune         ierr = PetscViewerASCIIPrintf(snes->ls_monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr);
463a5892487SPeter Brune         ierr = PetscViewerASCIISubtractTab(snes->ls_monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
464a5892487SPeter Brune       }
465a5892487SPeter Brune       break;
466a5892487SPeter Brune     }
467a5892487SPeter Brune     count++;
468a5892487SPeter Brune   }
469a5892487SPeter Brune   theend2:
470a5892487SPeter Brune   /* Optional user-defined check for line search step validity */
471a5892487SPeter Brune   if (snes->ops->postcheckstep) {
472a5892487SPeter Brune     ierr = (*snes->ops->postcheckstep)(snes,x,y,w,snes->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
473a5892487SPeter Brune     if (changed_y) {
474a5892487SPeter Brune       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
475a5892487SPeter Brune     }
476a5892487SPeter Brune     if (changed_y || changed_w) { /* recompute the function if the step has changed */
477a5892487SPeter Brune       ierr = SNESComputeFunction(snes,w,g);
478a5892487SPeter Brune       if (snes->domainerror) {
479a5892487SPeter Brune         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
480a5892487SPeter Brune         PetscFunctionReturn(0);
481a5892487SPeter Brune       }
482a5892487SPeter Brune       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
483a5892487SPeter Brune       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
484a5892487SPeter Brune       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
485a5892487SPeter Brune       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
486a5892487SPeter Brune       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
487a5892487SPeter Brune     }
488a5892487SPeter Brune   }
489a5892487SPeter Brune   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
490a5892487SPeter Brune   PetscFunctionReturn(0);
491a5892487SPeter Brune }
492a5892487SPeter Brune 
493a5892487SPeter Brune 
494*70d3b23bSPeter Brune EXTERN_C_BEGIN
495*70d3b23bSPeter Brune #undef __FUNCT__
496*70d3b23bSPeter Brune #define __FUNCT__ "SNESLineSearchSetType_LS"
497*70d3b23bSPeter Brune PetscErrorCode  SNESLineSearchSetType_LS(SNES snes, SNESLineSearchType type)
498*70d3b23bSPeter Brune {
499*70d3b23bSPeter Brune   PetscErrorCode ierr;
500*70d3b23bSPeter Brune   PetscFunctionBegin;
501*70d3b23bSPeter Brune 
502*70d3b23bSPeter Brune   switch (type) {
503*70d3b23bSPeter Brune   case SNES_LS_BASIC:
504*70d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
505*70d3b23bSPeter Brune     break;
506*70d3b23bSPeter Brune   case SNES_LS_BASIC_NONORMS:
507*70d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
508*70d3b23bSPeter Brune     break;
509*70d3b23bSPeter Brune   case SNES_LS_QUADRATIC:
510*70d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic_LS,PETSC_NULL);CHKERRQ(ierr);
511*70d3b23bSPeter Brune     break;
512*70d3b23bSPeter Brune   case SNES_LS_CUBIC:
513*70d3b23bSPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchCubic_LS,PETSC_NULL);CHKERRQ(ierr);
514*70d3b23bSPeter Brune     break;
515*70d3b23bSPeter Brune   default:
516*70d3b23bSPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type.");
517*70d3b23bSPeter Brune     break;
518*70d3b23bSPeter Brune   }
519*70d3b23bSPeter Brune   snes->ls_type = type;
520*70d3b23bSPeter Brune   PetscFunctionReturn(0);
521*70d3b23bSPeter Brune }
522*70d3b23bSPeter Brune EXTERN_C_END
523*70d3b23bSPeter Brune 
524*70d3b23bSPeter Brune 
52504d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
52604d965bbSLois Curfman McInnes 
52704d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
52894b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
52904d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
53004d965bbSLois Curfman McInnes      respectively.
53104d965bbSLois Curfman McInnes 
532fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
53304d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
53404d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
535fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
53604d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
53704d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
5384b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
5394b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
54004d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
54104d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
54204d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
54304d965bbSLois Curfman McInnes      for all nonlinear solvers.
54404d965bbSLois Curfman McInnes 
54504d965bbSLois Curfman McInnes      Another key routine is:
54604d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
54704d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
54804d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
54904d965bbSLois Curfman McInnes      SNESSolve() if necessary.
55004d965bbSLois Curfman McInnes 
55104d965bbSLois Curfman McInnes      Additional basic routines are:
55204d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
55304d965bbSLois Curfman McInnes                                       have actually been used.
55404d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
555186905e3SBarry Smith      SNESView().
55604d965bbSLois Curfman McInnes 
55704d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
55804d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
55904d965bbSLois Curfman McInnes      above description applies to these categories also.
56004d965bbSLois Curfman McInnes 
56104d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
5625e42470aSBarry Smith /*
5634b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
56404d965bbSLois Curfman McInnes    method with a line search.
5655e76c431SBarry Smith 
56604d965bbSLois Curfman McInnes    Input Parameters:
56704d965bbSLois Curfman McInnes .  snes - the SNES context
5685e76c431SBarry Smith 
56904d965bbSLois Curfman McInnes    Output Parameter:
570c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
57104d965bbSLois Curfman McInnes 
57204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
5735e76c431SBarry Smith 
5745e76c431SBarry Smith    Notes:
5755e76c431SBarry Smith    This implements essentially a truncated Newton method with a
5765e76c431SBarry Smith    line search.  By default a cubic backtracking line search
5775e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
5785e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
579393d2d9aSLois Curfman McInnes    and Schnabel.
5805e76c431SBarry Smith */
5814a2ae208SSatish Balay #undef __FUNCT__
5824b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
583dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
5845e76c431SBarry Smith {
5856849ba73SBarry Smith   PetscErrorCode     ierr;
586a7cc72afSBarry Smith   PetscInt           maxits,i,lits;
587ace3abfcSBarry Smith   PetscBool          lssucceed;
588112a2221SBarry Smith   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
58985385478SLisandro Dalcin   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
59085385478SLisandro Dalcin   Vec                Y,X,F,G,W;
5913d4c4710SBarry Smith   KSPConvergedReason kspreason;
5925e76c431SBarry Smith 
5933a40ed3dSBarry Smith   PetscFunctionBegin;
5943d4c4710SBarry Smith   snes->numFailures            = 0;
5953d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
596184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
597184914b5SBarry Smith 
5985e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
5995e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
60039e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
6015e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
6025e42470aSBarry Smith   G		= snes->work[1];
6035e42470aSBarry Smith   W		= snes->work[2];
6045e76c431SBarry Smith 
6054c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
6064c49b128SBarry Smith   snes->iter = 0;
60775567043SBarry Smith   snes->norm = 0.0;
6084c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
60919717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
6104936397dSBarry Smith   if (snes->domainerror) {
6114936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
6124936397dSBarry Smith     PetscFunctionReturn(0);
6134936397dSBarry Smith   }
6142613ca53SBarry Smith   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
6152613ca53SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
6162613ca53SBarry Smith   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
6172613ca53SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
61862cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6190f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
6205e42470aSBarry Smith   snes->norm = fnorm;
6210f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
62284c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
6237a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
6243f149594SLisandro Dalcin 
6253f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
6263f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
62785385478SLisandro Dalcin   /* test convergence */
62885385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
62906ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
630d034289bSBarry Smith 
6315e76c431SBarry Smith   for (i=0; i<maxits; i++) {
6325e76c431SBarry Smith 
633329e7e40SMatthew Knepley     /* Call general purpose update function */
634e7788613SBarry Smith     if (snes->ops->update) {
635e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
636de8cb200SMatthew Knepley     }
637329e7e40SMatthew Knepley 
638ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
6395334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
64094b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
64171f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
6423d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
6433d4c4710SBarry Smith     if (kspreason < 0) {
6443d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
645ef998cc9SBarry 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);
6463d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
647ab81109fSSatish Balay         break;
6483d4c4710SBarry Smith       }
6493d4c4710SBarry Smith     }
6503d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
6513f149594SLisandro Dalcin     snes->linear_its += lits;
6523f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
65374637425SBarry Smith 
654ea630c6eSPeter Brune     if (snes->ops->precheckstep) {
655ace3abfcSBarry Smith       PetscBool  changed_y = PETSC_FALSE;
656ea630c6eSPeter Brune       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
6579c3ca13aSBarry Smith     }
6589c3ca13aSBarry Smith 
659b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
6601e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
66185471664SBarry Smith     }
66274637425SBarry Smith 
663ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
664ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
665e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
666ea4d3ed3SLois Curfman McInnes     */
66785385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
6683f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
669ea630c6eSPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
6708f1a2a5eSBarry 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);
6714a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
6724936397dSBarry Smith     if (snes->domainerror) {
6734936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
6744936397dSBarry Smith       PetscFunctionReturn(0);
6754936397dSBarry Smith     }
676a7cc72afSBarry Smith     if (!lssucceed) {
67750ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
678ace3abfcSBarry Smith 	PetscBool  ismin;
679647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6801e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
6818a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
6823505fcc1SBarry Smith         break;
6833505fcc1SBarry Smith       }
68450ffb88aSMatthew Knepley     }
68585385478SLisandro Dalcin     /* Update function and solution vectors */
68685385478SLisandro Dalcin     fnorm = gnorm;
68785385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
68885385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
68985385478SLisandro Dalcin     /* Monitor convergence */
69085385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
69185385478SLisandro Dalcin     snes->iter = i+1;
69285385478SLisandro Dalcin     snes->norm = fnorm;
69385385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
69485385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
6957a03ce2fSLisandro Dalcin     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
69685385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
69785385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
698e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
6993f149594SLisandro Dalcin     if (snes->reason) break;
70016c95cacSBarry Smith   }
70152392280SLois Curfman McInnes   if (i == maxits) {
702ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
70385385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
70452392280SLois Curfman McInnes   }
7053a40ed3dSBarry Smith   PetscFunctionReturn(0);
7065e76c431SBarry Smith }
70704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
70804d965bbSLois Curfman McInnes /*
7094b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
7104b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
71104d965bbSLois Curfman McInnes 
71204d965bbSLois Curfman McInnes    Input Parameter:
71304d965bbSLois Curfman McInnes .  snes - the SNES context
71404d965bbSLois Curfman McInnes .  x - the solution vector
71504d965bbSLois Curfman McInnes 
71604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
71704d965bbSLois Curfman McInnes 
71804d965bbSLois Curfman McInnes    Notes:
7194b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
72004d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
72104d965bbSLois Curfman McInnes    the call to SNESSolve().
72204d965bbSLois Curfman McInnes  */
7234a2ae208SSatish Balay #undef __FUNCT__
7244b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
725dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
7265e76c431SBarry Smith {
727dfbe8321SBarry Smith   PetscErrorCode ierr;
7283a40ed3dSBarry Smith 
7293a40ed3dSBarry Smith   PetscFunctionBegin;
73058c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
7313a40ed3dSBarry Smith   PetscFunctionReturn(0);
7325e76c431SBarry Smith }
73304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7346b8b9a38SLisandro Dalcin 
7356b8b9a38SLisandro Dalcin #undef __FUNCT__
7366b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS"
7376b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes)
7386b8b9a38SLisandro Dalcin {
7396b8b9a38SLisandro Dalcin   PetscErrorCode ierr;
7406b8b9a38SLisandro Dalcin 
7416b8b9a38SLisandro Dalcin   PetscFunctionBegin;
7426b8b9a38SLisandro Dalcin   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
7436b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
7446b8b9a38SLisandro Dalcin }
7456b8b9a38SLisandro Dalcin 
74604d965bbSLois Curfman McInnes /*
7474b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
7484b27c08aSLois Curfman McInnes    with SNESCreate_LS().
74904d965bbSLois Curfman McInnes 
75004d965bbSLois Curfman McInnes    Input Parameter:
75104d965bbSLois Curfman McInnes .  snes - the SNES context
75204d965bbSLois Curfman McInnes 
753de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
75404d965bbSLois Curfman McInnes  */
7554a2ae208SSatish Balay #undef __FUNCT__
7564b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
757dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
7585e76c431SBarry Smith {
759dfbe8321SBarry Smith   PetscErrorCode ierr;
7603a40ed3dSBarry Smith 
7613a40ed3dSBarry Smith   PetscFunctionBegin;
7626b8b9a38SLisandro Dalcin   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
7635c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
764557d3f75SLisandro Dalcin 
7653a40ed3dSBarry Smith   PetscFunctionReturn(0);
7665e76c431SBarry Smith }
76704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
76894298653SBarry Smith 
769329e7e40SMatthew Knepley /*
7704b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
77104d965bbSLois Curfman McInnes 
77204d965bbSLois Curfman McInnes    Input Parameters:
77304d965bbSLois Curfman McInnes .  SNES - the SNES context
77404d965bbSLois Curfman McInnes .  viewer - visualization context
77504d965bbSLois Curfman McInnes 
77604d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
77704d965bbSLois Curfman McInnes */
7784a2ae208SSatish Balay #undef __FUNCT__
7794b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
7806849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
781a935fc98SLois Curfman McInnes {
7822fc52814SBarry Smith   const char     *cstr;
783dfbe8321SBarry Smith   PetscErrorCode ierr;
784ace3abfcSBarry Smith   PetscBool      iascii;
785a935fc98SLois Curfman McInnes 
7863a40ed3dSBarry Smith   PetscFunctionBegin;
7872692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
78832077d6dSBarry Smith   if (iascii) {
78915f5eeeaSPeter Brune     cstr = SNESLineSearchTypeName(snes->ls_type);
790b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
791ea630c6eSPeter 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);
792ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  damping factor=%14.12e\n",(double)snes->damping);CHKERRQ(ierr);
79350f6de3fSJed Brown   }
79450f6de3fSJed Brown   PetscFunctionReturn(0);
79550f6de3fSJed Brown }
79650f6de3fSJed Brown 
79704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
79804d965bbSLois Curfman McInnes /*
7994b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
80004d965bbSLois Curfman McInnes 
80104d965bbSLois Curfman McInnes    Input Parameter:
80204d965bbSLois Curfman McInnes .  snes - the SNES context
80304d965bbSLois Curfman McInnes 
80404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
80504d965bbSLois Curfman McInnes */
8064a2ae208SSatish Balay #undef __FUNCT__
8074b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
8086849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
8095e42470aSBarry Smith {
810dfbe8321SBarry Smith   PetscErrorCode ierr;
8115e42470aSBarry Smith 
8123a40ed3dSBarry Smith   PetscFunctionBegin;
813b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
814b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
8153a40ed3dSBarry Smith   PetscFunctionReturn(0);
8165e42470aSBarry Smith }
8174827ddcaSBarry Smith 
81804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
819ebe3b25bSBarry Smith /*MC
820ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
82104d965bbSLois Curfman McInnes 
822ebe3b25bSBarry Smith    Options Database:
823ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
82455d4ea13SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
825e106eecfSBarry 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)
826acbee50cSBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
82755d4ea13SBarry Smith .   -snes_ls_monitor - print information about progress of line searches
82855d4ea13SBarry Smith -   -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms
829acbee50cSBarry Smith 
83004d965bbSLois Curfman McInnes 
831ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
83204d965bbSLois Curfman McInnes 
833ee3001cbSBarry Smith    Level: beginner
834ee3001cbSBarry Smith 
83517bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
83617bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
837b3dd4ab5SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
838ebe3b25bSBarry Smith 
839ebe3b25bSBarry Smith M*/
840fb2e594dSBarry Smith EXTERN_C_BEGIN
8414a2ae208SSatish Balay #undef __FUNCT__
8424b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
8437087cfbeSBarry Smith PetscErrorCode  SNESCreate_LS(SNES snes)
8445e42470aSBarry Smith {
845dfbe8321SBarry Smith   PetscErrorCode ierr;
8464b27c08aSLois Curfman McInnes   SNES_LS        *neP;
8475e42470aSBarry Smith 
8483a40ed3dSBarry Smith   PetscFunctionBegin;
849e7788613SBarry Smith   snes->ops->setup           = SNESSetUp_LS;
850e7788613SBarry Smith   snes->ops->solve           = SNESSolve_LS;
851e7788613SBarry Smith   snes->ops->destroy         = SNESDestroy_LS;
852e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
853e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
8546b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_LS;
8555e42470aSBarry Smith 
85642f4f86dSBarry Smith   snes->usesksp                      = PETSC_TRUE;
85742f4f86dSBarry Smith   snes->usespc                       = PETSC_FALSE;
85838f2d2fdSLisandro Dalcin   ierr                               = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
8595e42470aSBarry Smith   snes->data                         = (void*)neP;
86015f5eeeaSPeter Brune   snes->ops->linesearchno            = SNESLineSearchNo;
861ea630c6eSPeter Brune   snes->lsP                          = PETSC_NULL;
862ea630c6eSPeter Brune   snes->ops->postcheckstep           = PETSC_NULL;
863ea630c6eSPeter Brune   snes->postcheck                    = PETSC_NULL;
864ea630c6eSPeter Brune   snes->ops->precheckstep            = PETSC_NULL;
865ea630c6eSPeter Brune   snes->precheck                     = PETSC_NULL;
866*70d3b23bSPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_LS",SNESLineSearchSetType_LS);CHKERRQ(ierr);
86715f5eeeaSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr);
8683a40ed3dSBarry Smith   PetscFunctionReturn(0);
8695e42470aSBarry Smith }
870fb2e594dSBarry Smith EXTERN_C_END
871