xref: /petsc/src/snes/impls/ls/ls.c (revision 91e9ee9fc0200d5446579c7185d9e655300e8525)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: ls.c,v 1.123 1999/02/09 23:32:19 bsmith Exp bsmith $";
3 #endif
4 
5 #include "src/snes/impls/ls/ls.h"
6 
7 /*  --------------------------------------------------------------------
8 
9      This file implements a truncated Newton method with a line search,
10      for solving a system of nonlinear equations, using the SLES, Vec,
11      and Mat interfaces for linear solvers, vectors, and matrices,
12      respectively.
13 
14      The following basic routines are required for each nonlinear solver:
15           SNESCreate_XXX()          - Creates a nonlinear solver context
16           SNESSetFromOptions_XXX()  - Sets runtime options
17           SNESSolve_XXX()           - Solves the nonlinear system
18           SNESDestroy_XXX()         - Destroys the nonlinear solver context
19      The suffix "_XXX" denotes a particular implementation, in this case
20      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
21      systems of nonlinear equations (EQ) with a line search (LS) method.
22      These routines are actually called via the common user interface
23      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
24      SNESDestroy(), so the application code interface remains identical
25      for all nonlinear solvers.
26 
27      Another key routine is:
28           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
29      by setting data structures and options.   The interface routine SNESSetUp()
30      is not usually called directly by the user, but instead is called by
31      SNESSolve() if necessary.
32 
33      Additional basic routines are:
34           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
35           SNESView_XXX()            - Prints details of runtime options that
36                                       have actually been used.
37      These are called by application codes via the interface routines
38      SNESPrintHelp() and SNESView().
39 
40      The various types of solvers (preconditioners, Krylov subspace methods,
41      nonlinear solvers, timesteppers) are all organized similarly, so the
42      above description applies to these categories also.
43 
44     -------------------------------------------------------------------- */
45 /*
46    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
47    method with a line search.
48 
49    Input Parameters:
50 .  snes - the SNES context
51 
52    Output Parameter:
53 .  outits - number of iterations until termination
54 
55    Application Interface Routine: SNESSolve()
56 
57    Notes:
58    This implements essentially a truncated Newton method with a
59    line search.  By default a cubic backtracking line search
60    is employed, as described in the text "Numerical Methods for
61    Unconstrained Optimization and Nonlinear Equations" by Dennis
62    and Schnabel.
63 */
64 #undef __FUNC__
65 #define __FUNC__ "SNESSolve_EQ_LS"
66 int SNESSolve_EQ_LS(SNES snes,int *outits)
67 {
68   SNES_LS       *neP = (SNES_LS *) snes->data;
69   int           maxits, i, ierr, lits, lsfail;
70   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
71   double        fnorm, gnorm, xnorm, ynorm;
72   Vec           Y, X, F, G, W, TMP;
73 
74   PetscFunctionBegin;
75   maxits	= snes->max_its;	/* maximum number of iterations */
76   X		= snes->vec_sol;	/* solution vector */
77   F		= snes->vec_func;	/* residual vector */
78   Y		= snes->work[0];	/* work vectors */
79   G		= snes->work[1];
80   W		= snes->work[2];
81 
82   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
83   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
84   PetscAMSTakeAccess(snes);
85   snes->iter = 0;
86   snes->norm = fnorm;
87   PetscAMSGrantAccess(snes);
88   SNESLogConvHistory(snes,fnorm,0);
89   SNESMonitor(snes,0,fnorm);
90 
91   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
92 
93   /* set parameter for default relative tolerance convergence test */
94   snes->ttol = fnorm*snes->rtol;
95 
96   for ( i=0; i<maxits; i++ ) {
97 
98     /* Solve J Y = F, where J is Jacobian matrix */
99     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
100     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
101     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
102     snes->linear_its += PetscAbsInt(lits);
103     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
104 
105     /* Compute a (scaled) negative update in the line search routine:
106          Y <- X - lambda*Y
107        and evaluate G(Y) = function(Y))
108     */
109     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
110     ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
111     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
112     if (lsfail) snes->nfailures++;
113 
114     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
115     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
116     fnorm = gnorm;
117 
118     PetscAMSTakeAccess(snes);
119     snes->iter = i+1;
120     snes->norm = fnorm;
121     PetscAMSGrantAccess(snes);
122     SNESLogConvHistory(snes,fnorm,lits);
123     SNESMonitor(snes,i+1,fnorm);
124 
125     /* Test for convergence */
126     if (snes->converged) {
127       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
128       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
129         break;
130       }
131     }
132   }
133   if (X != snes->vec_sol) {
134     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
135     snes->vec_sol_always  = snes->vec_sol;
136     snes->vec_func_always = snes->vec_func;
137   }
138   if (i == maxits) {
139     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
140     i--;
141   }
142   *outits = i+1;
143   PetscFunctionReturn(0);
144 }
145 /* -------------------------------------------------------------------------- */
146 /*
147    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
148    of the SNES_EQ_LS nonlinear solver.
149 
150    Input Parameter:
151 .  snes - the SNES context
152 .  x - the solution vector
153 
154    Application Interface Routine: SNESSetUp()
155 
156    Notes:
157    For basic use of the SNES solvers the user need not explicitly call
158    SNESSetUp(), since these actions will automatically occur during
159    the call to SNESSolve().
160  */
161 #undef __FUNC__
162 #define __FUNC__ "SNESSetUp_EQ_LS"
163 int SNESSetUp_EQ_LS(SNES snes)
164 {
165   int ierr;
166 
167   PetscFunctionBegin;
168   snes->nwork = 4;
169   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
170   PLogObjectParents(snes,snes->nwork,snes->work);
171   snes->vec_sol_update_always = snes->work[3];
172   PetscFunctionReturn(0);
173 }
174 /* -------------------------------------------------------------------------- */
175 /*
176    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
177    with SNESCreate_EQ_LS().
178 
179    Input Parameter:
180 .  snes - the SNES context
181 
182    Application Interface Routine: SNESDestroy()
183  */
184 #undef __FUNC__
185 #define __FUNC__ "SNESDestroy_EQ_LS"
186 int SNESDestroy_EQ_LS(SNES snes)
187 {
188   int  ierr;
189 
190   PetscFunctionBegin;
191   if (snes->nwork) {
192     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
193   }
194   PetscFree(snes->data);
195   PetscFunctionReturn(0);
196 }
197 /* -------------------------------------------------------------------------- */
198 #undef __FUNC__
199 #define __FUNC__ "SNESNoLineSearch"
200 
201 /*@C
202    SNESNoLineSearch - This routine is not a line search at all;
203    it simply uses the full Newton step.  Thus, this routine is intended
204    to serve as a template and is not recommended for general use.
205 
206    Collective on SNES and Vec
207 
208    Input Parameters:
209 +  snes - nonlinear context
210 .  x - current iterate
211 .  f - residual evaluated at x
212 .  y - search direction (contains new iterate on output)
213 .  w - work vector
214 -  fnorm - 2-norm of f
215 
216    Output Parameters:
217 +  g - residual evaluated at new iterate y
218 .  y - new iterate (contains search direction on input)
219 .  gnorm - 2-norm of g
220 .  ynorm - 2-norm of search length
221 -  flag - set to 0, indicating a successful line search
222 
223    Options Database Key:
224 .  -snes_eq_ls basic - Activates SNESNoLineSearch()
225 
226    Level: advanced
227 
228 .keywords: SNES, nonlinear, line search, cubic
229 
230 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
231           SNESSetLineSearch()
232 @*/
233 int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
234                      double fnorm, double *ynorm, double *gnorm,int *flag )
235 {
236   int    ierr;
237   Scalar mone = -1.0;
238 
239   PetscFunctionBegin;
240   *flag = 0;
241   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
242   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
243   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
244   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
245   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
246   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
247   PetscFunctionReturn(0);
248 }
249 /* -------------------------------------------------------------------------- */
250 #undef __FUNC__
251 #define __FUNC__ "SNESNoLineSearchNoNorms"
252 
253 /*@C
254    SNESNoLineSearchNoNorms - This routine is not a line search at
255    all; it simply uses the full Newton step. This version does not
256    even compute the norm of the function or search direction; this
257    is intended only when you know the full step is fine and are
258    not checking for convergence of the nonlinear iteration (for
259    example, you are running always for a fixed number of Newton steps).
260 
261    Collective on SNES and Vec
262 
263    Input Parameters:
264 +  snes - nonlinear context
265 .  x - current iterate
266 .  f - residual evaluated at x
267 .  y - search direction (contains new iterate on output)
268 .  w - work vector
269 -  fnorm - 2-norm of f
270 
271    Output Parameters:
272 +  g - residual evaluated at new iterate y
273 .  gnorm - not changed
274 .  ynorm - not changed
275 -  flag - set to 0, indicating a successful line search
276 
277    Options Database Key:
278 .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
279 
280    Level: advanced
281 
282 .keywords: SNES, nonlinear, line search, cubic
283 
284 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
285           SNESSetLineSearch(), SNESNoLineSearch()
286 @*/
287 int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
288                      double fnorm, double *ynorm, double *gnorm,int *flag )
289 {
290   int    ierr;
291   Scalar mone = -1.0;
292 
293   PetscFunctionBegin;
294   *flag = 0;
295   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
296   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
297   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
298   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
299   PetscFunctionReturn(0);
300 }
301 /* -------------------------------------------------------------------------- */
302 #undef __FUNC__
303 #define __FUNC__ "SNESCubicLineSearch"
304 /*@C
305    SNESCubicLineSearch - Performs a cubic line search (default line search method).
306 
307    Collective on SNES
308 
309    Input Parameters:
310 +  snes - nonlinear context
311 .  x - current iterate
312 .  f - residual evaluated at x
313 .  y - search direction (contains new iterate on output)
314 .  w - work vector
315 -  fnorm - 2-norm of f
316 
317    Output Parameters:
318 +  g - residual evaluated at new iterate y
319 .  y - new iterate (contains search direction on input)
320 .  gnorm - 2-norm of g
321 .  ynorm - 2-norm of search length
322 -  flag - 0 if line search succeeds; -1 on failure.
323 
324    Options Database Key:
325 $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
326 
327    Notes:
328    This line search is taken from "Numerical Methods for Unconstrained
329    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
330 
331    Level: advanced
332 
333 .keywords: SNES, nonlinear, line search, cubic
334 
335 .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
336 @*/
337 int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
338                         double fnorm,double *ynorm,double *gnorm,int *flag)
339 {
340   /*
341      Note that for line search purposes we work with with the related
342      minimization problem:
343         min  z(x):  R^n -> R,
344      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
345    */
346 
347   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
348   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
349 #if defined(USE_PETSC_COMPLEX)
350   Scalar  cinitslope, clambda;
351 #endif
352   int     ierr, count;
353   SNES_LS *neP = (SNES_LS *) snes->data;
354   Scalar  mone = -1.0,scale;
355 
356   PetscFunctionBegin;
357   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
358   *flag   = 0;
359   alpha   = neP->alpha;
360   maxstep = neP->maxstep;
361   steptol = neP->steptol;
362 
363   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
364   if (*ynorm < snes->atol) {
365     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
366     *gnorm = fnorm;
367     ierr = VecCopy(x,y); CHKERRQ(ierr);
368     ierr = VecCopy(f,g); CHKERRQ(ierr);
369     goto theend1;
370   }
371   if (*ynorm > maxstep) {	/* Step too big, so scale back */
372     scale = maxstep/(*ynorm);
373 #if defined(USE_PETSC_COMPLEX)
374     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
375 #else
376     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
377 #endif
378     ierr = VecScale(&scale,y); CHKERRQ(ierr);
379     *ynorm = maxstep;
380   }
381   minlambda = steptol/(*ynorm);
382   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
383 #if defined(USE_PETSC_COMPLEX)
384   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
385   initslope = PetscReal(cinitslope);
386 #else
387   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
388 #endif
389   if (initslope > 0.0) initslope = -initslope;
390   if (initslope == 0.0) initslope = -1.0;
391 
392   ierr = VecCopy(y,w); CHKERRQ(ierr);
393   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
394   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
395   ierr = VecNorm(g,NORM_2,gnorm);
396   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
397     ierr = VecCopy(w,y); CHKERRQ(ierr);
398     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
399     goto theend1;
400   }
401 
402   /* Fit points with quadratic */
403   lambda = 1.0; count = 0;
404   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
405   lambdaprev = lambda;
406   gnormprev = *gnorm;
407   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
408   else lambda = lambdatemp;
409   ierr   = VecCopy(x,w); CHKERRQ(ierr);
410   lambdaneg = -lambda;
411 #if defined(USE_PETSC_COMPLEX)
412   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
413 #else
414   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
415 #endif
416   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
417   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
418   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
419     ierr = VecCopy(w,y); CHKERRQ(ierr);
420     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
421     goto theend1;
422   }
423 
424   /* Fit points with cubic */
425   count = 1;
426   while (1) {
427     if (lambda <= minlambda) { /* bad luck; use full step */
428       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
429       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
430                fnorm,*gnorm,*ynorm,lambda,initslope);
431       ierr = VecCopy(w,y); CHKERRQ(ierr);
432       *flag = -1; break;
433     }
434     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
435     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
436     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
437     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
438     d  = b*b - 3*a*initslope;
439     if (d < 0.0) d = 0.0;
440     if (a == 0.0) {
441       lambdatemp = -initslope/(2.0*b);
442     } else {
443       lambdatemp = (-b + sqrt(d))/(3.0*a);
444     }
445     if (lambdatemp > .5*lambda) {
446       lambdatemp = .5*lambda;
447     }
448     lambdaprev = lambda;
449     gnormprev = *gnorm;
450     if (lambdatemp <= .1*lambda) {
451       lambda = .1*lambda;
452     }
453     else lambda = lambdatemp;
454     ierr = VecCopy( x, w ); CHKERRQ(ierr);
455     lambdaneg = -lambda;
456 #if defined(USE_PETSC_COMPLEX)
457     clambda = lambdaneg;
458     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
459 #else
460     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
461 #endif
462     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
463     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
464     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
465       ierr = VecCopy(w,y); CHKERRQ(ierr);
466       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
467       goto theend1;
468     } else {
469       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
470     }
471     count++;
472   }
473   theend1:
474   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
475   PetscFunctionReturn(0);
476 }
477 /* -------------------------------------------------------------------------- */
478 #undef __FUNC__
479 #define __FUNC__ "SNESQuadraticLineSearch"
480 /*@C
481    SNESQuadraticLineSearch - Performs a quadratic line search.
482 
483    Collective on SNES and Vec
484 
485    Input Parameters:
486 +  snes - the SNES context
487 .  x - current iterate
488 .  f - residual evaluated at x
489 .  y - search direction (contains new iterate on output)
490 .  w - work vector
491 -  fnorm - 2-norm of f
492 
493    Output Parameters:
494 +  g - residual evaluated at new iterate y
495 .  y - new iterate (contains search direction on input)
496 .  gnorm - 2-norm of g
497 .  ynorm - 2-norm of search length
498 -  flag - 0 if line search succeeds; -1 on failure.
499 
500    Options Database Key:
501 .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
502 
503    Notes:
504    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
505 
506    Level: advanced
507 
508 .keywords: SNES, nonlinear, quadratic, line search
509 
510 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
511 @*/
512 int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
513                            double fnorm, double *ynorm, double *gnorm,int *flag)
514 {
515   /*
516      Note that for line search purposes we work with with the related
517      minimization problem:
518         min  z(x):  R^n -> R,
519      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
520    */
521   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
522 #if defined(USE_PETSC_COMPLEX)
523   Scalar  cinitslope,clambda;
524 #endif
525   int     ierr,count;
526   SNES_LS *neP = (SNES_LS *) snes->data;
527   Scalar  mone = -1.0,scale;
528 
529   PetscFunctionBegin;
530   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
531   *flag   = 0;
532   alpha   = neP->alpha;
533   maxstep = neP->maxstep;
534   steptol = neP->steptol;
535 
536   VecNorm(y, NORM_2,ynorm );
537   if (*ynorm < snes->atol) {
538     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
539     *gnorm = fnorm;
540     ierr = VecCopy(x,y); CHKERRQ(ierr);
541     ierr = VecCopy(f,g); CHKERRQ(ierr);
542     goto theend2;
543   }
544   if (*ynorm > maxstep) {	/* Step too big, so scale back */
545     scale = maxstep/(*ynorm);
546     ierr = VecScale(&scale,y); CHKERRQ(ierr);
547     *ynorm = maxstep;
548   }
549   minlambda = steptol/(*ynorm);
550   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
551 #if defined(USE_PETSC_COMPLEX)
552   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
553   initslope = PetscReal(cinitslope);
554 #else
555   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
556 #endif
557   if (initslope > 0.0) initslope = -initslope;
558   if (initslope == 0.0) initslope = -1.0;
559 
560   ierr = VecCopy(y,w); CHKERRQ(ierr);
561   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
562   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
563   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
564   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
565     ierr = VecCopy(w,y); CHKERRQ(ierr);
566     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
567     goto theend2;
568   }
569 
570   /* Fit points with quadratic */
571   lambda = 1.0; count = 0;
572   count = 1;
573   while (1) {
574     if (lambda <= minlambda) { /* bad luck; use full step */
575       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
576       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
577                fnorm,*gnorm,*ynorm,lambda,initslope);
578       ierr = VecCopy(w,y); CHKERRQ(ierr);
579       *flag = -1; break;
580     }
581     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
582     if (lambdatemp <= .1*lambda) {
583       lambda = .1*lambda;
584     } else lambda = lambdatemp;
585     ierr = VecCopy(x,w); CHKERRQ(ierr);
586     lambda = -lambda;
587 #if defined(USE_PETSC_COMPLEX)
588     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
589 #else
590     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
591 #endif
592     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
593     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
594     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
595       ierr = VecCopy(w,y); CHKERRQ(ierr);
596       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
597       break;
598     }
599     count++;
600   }
601   theend2:
602   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
603   PetscFunctionReturn(0);
604 }
605 /* -------------------------------------------------------------------------- */
606 #undef __FUNC__
607 #define __FUNC__ "SNESSetLineSearch"
608 /*@C
609    SNESSetLineSearch - Sets the line search routine to be used
610    by the method SNES_EQ_LS.
611 
612    Input Parameters:
613 +  snes - nonlinear context obtained from SNESCreate()
614 -  func - pointer to int function
615 
616    Collective on SNES
617 
618    Available Routines:
619 +  SNESCubicLineSearch() - default line search
620 .  SNESQuadraticLineSearch() - quadratic line search
621 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
622 -  SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel)
623 
624     Options Database Keys:
625 +   -snes_eq_ls [basic,quadratic,cubic] - Selects line search
626 .   -snes_eq_ls_alpha <alpha> - Sets alpha
627 .   -snes_eq_ls_maxstep <max> - Sets maxstep
628 -   -snes_eq_ls_steptol <steptol> - Sets steptol
629 
630    Calling sequence of func:
631 .vb
632    func (SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
633          double fnorm, double *ynorm, double *gnorm, *flag)
634 .ve
635 
636     Input parameters for func:
637 +   snes - nonlinear context
638 .   x - current iterate
639 .   f - residual evaluated at x
640 .   y - search direction (contains new iterate on output)
641 .   w - work vector
642 -   fnorm - 2-norm of f
643 
644     Output parameters for func:
645 +   g - residual evaluated at new iterate y
646 .   y - new iterate (contains search direction on input)
647 .   gnorm - 2-norm of g
648 .   ynorm - 2-norm of search length
649 -   flag - set to 0 if the line search succeeds; a nonzero integer
650            on failure.
651 
652     Level: advanced
653 
654 .keywords: SNES, nonlinear, set, line search, routine
655 
656 .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
657 @*/
658 int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*))
659 {
660   int ierr, (*f)(SNES,int (*)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*));
661 
662   PetscFunctionBegin;
663   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
664   if (f) {
665     ierr = (*f)(snes,func);CHKERRQ(ierr);
666   }
667   PetscFunctionReturn(0);
668 }
669 /* -------------------------------------------------------------------------- */
670 EXTERN_C_BEGIN
671 #undef __FUNC__
672 #define __FUNC__ "SNESSetLineSearch_LS"
673 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
674                          double,double*,double*,int*))
675 {
676   PetscFunctionBegin;
677   ((SNES_LS *)(snes->data))->LineSearch = func;
678   PetscFunctionReturn(0);
679 }
680 EXTERN_C_END
681 /* -------------------------------------------------------------------------- */
682 /*
683    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
684 
685    Input Parameter:
686 .  snes - the SNES context
687 
688    Application Interface Routine: SNESPrintHelp()
689 */
690 #undef __FUNC__
691 #define __FUNC__ "SNESPrintHelp_EQ_LS"
692 static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
693 {
694   SNES_LS *ls = (SNES_LS *)snes->data;
695 
696   PetscFunctionBegin;
697   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
698   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
699   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
700   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
701   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
702   PetscFunctionReturn(0);
703 }
704 /* -------------------------------------------------------------------------- */
705 /*
706    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
707 
708    Input Parameters:
709 .  SNES - the SNES context
710 .  viewer - visualization context
711 
712    Application Interface Routine: SNESView()
713 */
714 #undef __FUNC__
715 #define __FUNC__ "SNESView_EQ_LS"
716 static int SNESView_EQ_LS(SNES snes,Viewer viewer)
717 {
718   SNES_LS    *ls = (SNES_LS *)snes->data;
719   char       *cstr;
720   int        ierr;
721   ViewerType vtype;
722 
723   PetscFunctionBegin;
724   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
725   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
726     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
727     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
728     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
729     else                                                cstr = "unknown";
730     ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);
731     ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);
732   } else {
733     SETERRQ(1,1,"Viewer type not supported for this object");
734   }
735   PetscFunctionReturn(0);
736 }
737 /* -------------------------------------------------------------------------- */
738 /*
739    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
740 
741    Input Parameter:
742 .  snes - the SNES context
743 
744    Application Interface Routine: SNESSetFromOptions()
745 */
746 #undef __FUNC__
747 #define __FUNC__ "SNESSetFromOptions_EQ_LS"
748 static int SNESSetFromOptions_EQ_LS(SNES snes)
749 {
750   SNES_LS *ls = (SNES_LS *)snes->data;
751   char    ver[16];
752   double  tmp;
753   int     ierr,flg;
754 
755   PetscFunctionBegin;
756   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
757   if (flg) {
758     ls->alpha = tmp;
759   }
760   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
761   if (flg) {
762     ls->maxstep = tmp;
763   }
764   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
765   if (flg) {
766     ls->steptol = tmp;
767   }
768   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
769   if (flg) {
770     if (!PetscStrcmp(ver,"basic")) {
771       ierr = SNESSetLineSearch(snes,SNESNoLineSearch);CHKERRQ(ierr);
772     } else if (!PetscStrcmp(ver,"basicnonorms")) {
773       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);CHKERRQ(ierr);
774     } else if (!PetscStrcmp(ver,"quadratic")) {
775       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch);CHKERRQ(ierr);
776     } else if (!PetscStrcmp(ver,"cubic")) {
777       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch);CHKERRQ(ierr);
778     }
779     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
780   }
781   PetscFunctionReturn(0);
782 }
783 /* -------------------------------------------------------------------------- */
784 /*
785    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
786    SNES_LS, and sets this as the private data within the generic nonlinear solver
787    context, SNES, that was created within SNESCreate().
788 
789 
790    Input Parameter:
791 .  snes - the SNES context
792 
793    Application Interface Routine: SNESCreate()
794  */
795 EXTERN_C_BEGIN
796 #undef __FUNC__
797 #define __FUNC__ "SNESCreate_EQ_LS"
798 int SNESCreate_EQ_LS(SNES snes)
799 {
800   int     ierr;
801   SNES_LS *neP;
802 
803   PetscFunctionBegin;
804   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
805     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
806   }
807 
808   snes->setup		= SNESSetUp_EQ_LS;
809   snes->solve		= SNESSolve_EQ_LS;
810   snes->destroy		= SNESDestroy_EQ_LS;
811   snes->converged	= SNESConverged_EQ_LS;
812   snes->printhelp       = SNESPrintHelp_EQ_LS;
813   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
814   snes->view            = SNESView_EQ_LS;
815   snes->nwork           = 0;
816 
817   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
818   PLogObjectMemory(snes,sizeof(SNES_LS));
819   snes->data    	= (void *) neP;
820   neP->alpha		= 1.e-4;
821   neP->maxstep		= 1.e8;
822   neP->steptol		= 1.e-12;
823   neP->LineSearch       = SNESCubicLineSearch;
824 
825   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
826                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
827 
828   PetscFunctionReturn(0);
829 }
830 EXTERN_C_END
831 
832 
833 
834