xref: /petsc/src/snes/impls/ls/ls.c (revision d5d37b614d6ea02bb947345c08c45a4da2fd754a)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: ls.c,v 1.120 1998/12/23 22:53:13 bsmith Exp curfman $";
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, history_len, ierr, lits, lsfail;
70   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
71   double        fnorm, gnorm, xnorm, ynorm, *history;
72   Vec           Y, X, F, G, W, TMP;
73 
74   PetscFunctionBegin;
75   history	= snes->conv_hist;	/* convergence history */
76   history_len	= snes->conv_hist_size;	/* convergence history length */
77   maxits	= snes->max_its;	/* maximum number of iterations */
78   X		= snes->vec_sol;	/* solution vector */
79   F		= snes->vec_func;	/* residual vector */
80   Y		= snes->work[0];	/* work vectors */
81   G		= snes->work[1];
82   W		= snes->work[2];
83 
84   snes->iter = 0;
85   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
86   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
87   snes->norm = fnorm;
88   if (history) history[0] = fnorm;
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     snes->iter = i+1;
98 
99     /* Solve J Y = F, where J is Jacobian matrix */
100     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
101     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
102     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
103     snes->linear_its += PetscAbsInt(lits);
104     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
105 
106     /* Compute a (scaled) negative update in the line search routine:
107          Y <- X - lambda*Y
108        and evaluate G(Y) = function(Y))
109     */
110     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
111     ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
112     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
113     if (lsfail) snes->nfailures++;
114 
115     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
116     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
117     fnorm = gnorm;
118 
119     snes->norm = fnorm;
120     if (history && history_len > i+1) history[i+1] = fnorm;
121     SNESMonitor(snes,i+1,fnorm);
122 
123     /* Test for convergence */
124     if (snes->converged) {
125       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
126       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
127         break;
128       }
129     }
130   }
131   if (X != snes->vec_sol) {
132     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
133     snes->vec_sol_always  = snes->vec_sol;
134     snes->vec_func_always = snes->vec_func;
135   }
136   if (i == maxits) {
137     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
138     i--;
139   }
140   if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1;
141   *outits = i+1;
142   PetscFunctionReturn(0);
143 }
144 /* -------------------------------------------------------------------------- */
145 /*
146    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
147    of the SNES_EQ_LS nonlinear solver.
148 
149    Input Parameter:
150 .  snes - the SNES context
151 .  x - the solution vector
152 
153    Application Interface Routine: SNESSetUp()
154 
155    Notes:
156    For basic use of the SNES solvers the user need not explicitly call
157    SNESSetUp(), since these actions will automatically occur during
158    the call to SNESSolve().
159  */
160 #undef __FUNC__
161 #define __FUNC__ "SNESSetUp_EQ_LS"
162 int SNESSetUp_EQ_LS(SNES snes)
163 {
164   int ierr;
165 
166   PetscFunctionBegin;
167   snes->nwork = 4;
168   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
169   PLogObjectParents(snes,snes->nwork,snes->work);
170   snes->vec_sol_update_always = snes->work[3];
171   PetscFunctionReturn(0);
172 }
173 /* -------------------------------------------------------------------------- */
174 /*
175    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
176    with SNESCreate_EQ_LS().
177 
178    Input Parameter:
179 .  snes - the SNES context
180 
181    Application Interface Routine: SNESDestroy()
182  */
183 #undef __FUNC__
184 #define __FUNC__ "SNESDestroy_EQ_LS"
185 int SNESDestroy_EQ_LS(SNES snes)
186 {
187   int  ierr;
188 
189   PetscFunctionBegin;
190   if (snes->nwork) {
191     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
192   }
193   PetscFree(snes->data);
194   PetscFunctionReturn(0);
195 }
196 /* -------------------------------------------------------------------------- */
197 #undef __FUNC__
198 #define __FUNC__ "SNESNoLineSearch"
199 
200 /*@C
201    SNESNoLineSearch - This routine is not a line search at all;
202    it simply uses the full Newton step.  Thus, this routine is intended
203    to serve as a template and is not recommended for general use.
204 
205    Collective on SNES and Vec
206 
207    Input Parameters:
208 +  snes - nonlinear context
209 .  x - current iterate
210 .  f - residual evaluated at x
211 .  y - search direction (contains new iterate on output)
212 .  w - work vector
213 -  fnorm - 2-norm of f
214 
215    Output Parameters:
216 +  g - residual evaluated at new iterate y
217 .  y - new iterate (contains search direction on input)
218 .  gnorm - 2-norm of g
219 .  ynorm - 2-norm of search length
220 -  flag - set to 0, indicating a successful line search
221 
222    Options Database Key:
223 .  -snes_eq_ls basic - Activates SNESNoLineSearch()
224 
225    Level: advanced
226 
227 .keywords: SNES, nonlinear, line search, cubic
228 
229 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
230           SNESSetLineSearch()
231 @*/
232 int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
233                      double fnorm, double *ynorm, double *gnorm,int *flag )
234 {
235   int    ierr;
236   Scalar mone = -1.0;
237 
238   PetscFunctionBegin;
239   *flag = 0;
240   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
241   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
242   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
243   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
244   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
245   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
246   PetscFunctionReturn(0);
247 }
248 /* -------------------------------------------------------------------------- */
249 #undef __FUNC__
250 #define __FUNC__ "SNESNoLineSearchNoNorms"
251 
252 /*@C
253    SNESNoLineSearchNoNorms - This routine is not a line search at
254    all; it simply uses the full Newton step. This version does not
255    even compute the norm of the function or search direction; this
256    is intended only when you know the full step is fine and are
257    not checking for convergence of the nonlinear iteration (for
258    example, you are running always for a fixed number of Newton steps).
259 
260    Collective on SNES and Vec
261 
262    Input Parameters:
263 +  snes - nonlinear context
264 .  x - current iterate
265 .  f - residual evaluated at x
266 .  y - search direction (contains new iterate on output)
267 .  w - work vector
268 -  fnorm - 2-norm of f
269 
270    Output Parameters:
271 +  g - residual evaluated at new iterate y
272 .  gnorm - not changed
273 .  ynorm - not changed
274 -  flag - set to 0, indicating a successful line search
275 
276    Options Database Key:
277 .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
278 
279    Level: advanced
280 
281 .keywords: SNES, nonlinear, line search, cubic
282 
283 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
284           SNESSetLineSearch(), SNESNoLineSearch()
285 @*/
286 int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
287                      double fnorm, double *ynorm, double *gnorm,int *flag )
288 {
289   int    ierr;
290   Scalar mone = -1.0;
291 
292   PetscFunctionBegin;
293   *flag = 0;
294   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
295   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
296   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
297   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
298   PetscFunctionReturn(0);
299 }
300 /* -------------------------------------------------------------------------- */
301 #undef __FUNC__
302 #define __FUNC__ "SNESCubicLineSearch"
303 /*@C
304    SNESCubicLineSearch - Performs a cubic line search (default line search method).
305 
306    Collective on SNES
307 
308    Input Parameters:
309 +  snes - nonlinear context
310 .  x - current iterate
311 .  f - residual evaluated at x
312 .  y - search direction (contains new iterate on output)
313 .  w - work vector
314 -  fnorm - 2-norm of f
315 
316    Output Parameters:
317 +  g - residual evaluated at new iterate y
318 .  y - new iterate (contains search direction on input)
319 .  gnorm - 2-norm of g
320 .  ynorm - 2-norm of search length
321 -  flag - 0 if line search succeeds; -1 on failure.
322 
323    Options Database Key:
324 $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
325 
326    Notes:
327    This line search is taken from "Numerical Methods for Unconstrained
328    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
329 
330    Level: advanced
331 
332 .keywords: SNES, nonlinear, line search, cubic
333 
334 .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
335 @*/
336 int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
337                         double fnorm,double *ynorm,double *gnorm,int *flag)
338 {
339   /*
340      Note that for line search purposes we work with with the related
341      minimization problem:
342         min  z(x):  R^n -> R,
343      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
344    */
345 
346   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
347   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
348 #if defined(USE_PETSC_COMPLEX)
349   Scalar  cinitslope, clambda;
350 #endif
351   int     ierr, count;
352   SNES_LS *neP = (SNES_LS *) snes->data;
353   Scalar  mone = -1.0,scale;
354 
355   PetscFunctionBegin;
356   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
357   *flag   = 0;
358   alpha   = neP->alpha;
359   maxstep = neP->maxstep;
360   steptol = neP->steptol;
361 
362   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
363   if (*ynorm < snes->atol) {
364     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
365     *gnorm = fnorm;
366     ierr = VecCopy(x,y); CHKERRQ(ierr);
367     ierr = VecCopy(f,g); CHKERRQ(ierr);
368     goto theend1;
369   }
370   if (*ynorm > maxstep) {	/* Step too big, so scale back */
371     scale = maxstep/(*ynorm);
372 #if defined(USE_PETSC_COMPLEX)
373     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
374 #else
375     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
376 #endif
377     ierr = VecScale(&scale,y); CHKERRQ(ierr);
378     *ynorm = maxstep;
379   }
380   minlambda = steptol/(*ynorm);
381   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
382 #if defined(USE_PETSC_COMPLEX)
383   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
384   initslope = PetscReal(cinitslope);
385 #else
386   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
387 #endif
388   if (initslope > 0.0) initslope = -initslope;
389   if (initslope == 0.0) initslope = -1.0;
390 
391   ierr = VecCopy(y,w); CHKERRQ(ierr);
392   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
393   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
394   ierr = VecNorm(g,NORM_2,gnorm);
395   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
396     ierr = VecCopy(w,y); CHKERRQ(ierr);
397     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
398     goto theend1;
399   }
400 
401   /* Fit points with quadratic */
402   lambda = 1.0; count = 0;
403   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
404   lambdaprev = lambda;
405   gnormprev = *gnorm;
406   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
407   else lambda = lambdatemp;
408   ierr   = VecCopy(x,w); CHKERRQ(ierr);
409   lambdaneg = -lambda;
410 #if defined(USE_PETSC_COMPLEX)
411   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
412 #else
413   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
414 #endif
415   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
416   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
417   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
418     ierr = VecCopy(w,y); CHKERRQ(ierr);
419     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
420     goto theend1;
421   }
422 
423   /* Fit points with cubic */
424   count = 1;
425   while (1) {
426     if (lambda <= minlambda) { /* bad luck; use full step */
427       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
428       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
429                fnorm,*gnorm,*ynorm,lambda,initslope);
430       ierr = VecCopy(w,y); CHKERRQ(ierr);
431       *flag = -1; break;
432     }
433     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
434     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
435     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
436     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
437     d  = b*b - 3*a*initslope;
438     if (d < 0.0) d = 0.0;
439     if (a == 0.0) {
440       lambdatemp = -initslope/(2.0*b);
441     } else {
442       lambdatemp = (-b + sqrt(d))/(3.0*a);
443     }
444     if (lambdatemp > .5*lambda) {
445       lambdatemp = .5*lambda;
446     }
447     lambdaprev = lambda;
448     gnormprev = *gnorm;
449     if (lambdatemp <= .1*lambda) {
450       lambda = .1*lambda;
451     }
452     else lambda = lambdatemp;
453     ierr = VecCopy( x, w ); CHKERRQ(ierr);
454     lambdaneg = -lambda;
455 #if defined(USE_PETSC_COMPLEX)
456     clambda = lambdaneg;
457     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
458 #else
459     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
460 #endif
461     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
462     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
463     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
464       ierr = VecCopy(w,y); CHKERRQ(ierr);
465       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
466       goto theend1;
467     } else {
468       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
469     }
470     count++;
471   }
472   theend1:
473   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
474   PetscFunctionReturn(0);
475 }
476 /* -------------------------------------------------------------------------- */
477 #undef __FUNC__
478 #define __FUNC__ "SNESQuadraticLineSearch"
479 /*@C
480    SNESQuadraticLineSearch - Performs a quadratic line search.
481 
482    Collective on SNES and Vec
483 
484    Input Parameters:
485 +  snes - the SNES context
486 .  x - current iterate
487 .  f - residual evaluated at x
488 .  y - search direction (contains new iterate on output)
489 .  w - work vector
490 -  fnorm - 2-norm of f
491 
492    Output Parameters:
493 +  g - residual evaluated at new iterate y
494 .  y - new iterate (contains search direction on input)
495 .  gnorm - 2-norm of g
496 .  ynorm - 2-norm of search length
497 -  flag - 0 if line search succeeds; -1 on failure.
498 
499    Options Database Key:
500 .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
501 
502    Notes:
503    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
504 
505    Level: advanced
506 
507 .keywords: SNES, nonlinear, quadratic, line search
508 
509 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
510 @*/
511 int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
512                            double fnorm, double *ynorm, double *gnorm,int *flag)
513 {
514   /*
515      Note that for line search purposes we work with with the related
516      minimization problem:
517         min  z(x):  R^n -> R,
518      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
519    */
520   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
521 #if defined(USE_PETSC_COMPLEX)
522   Scalar  cinitslope,clambda;
523 #endif
524   int     ierr,count;
525   SNES_LS *neP = (SNES_LS *) snes->data;
526   Scalar  mone = -1.0,scale;
527 
528   PetscFunctionBegin;
529   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
530   *flag   = 0;
531   alpha   = neP->alpha;
532   maxstep = neP->maxstep;
533   steptol = neP->steptol;
534 
535   VecNorm(y, NORM_2,ynorm );
536   if (*ynorm < snes->atol) {
537     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
538     *gnorm = fnorm;
539     ierr = VecCopy(x,y); CHKERRQ(ierr);
540     ierr = VecCopy(f,g); CHKERRQ(ierr);
541     goto theend2;
542   }
543   if (*ynorm > maxstep) {	/* Step too big, so scale back */
544     scale = maxstep/(*ynorm);
545     ierr = VecScale(&scale,y); CHKERRQ(ierr);
546     *ynorm = maxstep;
547   }
548   minlambda = steptol/(*ynorm);
549   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
550 #if defined(USE_PETSC_COMPLEX)
551   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
552   initslope = PetscReal(cinitslope);
553 #else
554   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
555 #endif
556   if (initslope > 0.0) initslope = -initslope;
557   if (initslope == 0.0) initslope = -1.0;
558 
559   ierr = VecCopy(y,w); CHKERRQ(ierr);
560   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
561   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
562   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
563   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
564     ierr = VecCopy(w,y); CHKERRQ(ierr);
565     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
566     goto theend2;
567   }
568 
569   /* Fit points with quadratic */
570   lambda = 1.0; count = 0;
571   count = 1;
572   while (1) {
573     if (lambda <= minlambda) { /* bad luck; use full step */
574       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
575       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
576                fnorm,*gnorm,*ynorm,lambda,initslope);
577       ierr = VecCopy(w,y); CHKERRQ(ierr);
578       *flag = -1; break;
579     }
580     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
581     if (lambdatemp <= .1*lambda) {
582       lambda = .1*lambda;
583     } else lambda = lambdatemp;
584     ierr = VecCopy(x,w); CHKERRQ(ierr);
585     lambda = -lambda;
586 #if defined(USE_PETSC_COMPLEX)
587     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
588 #else
589     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
590 #endif
591     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
592     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
593     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
594       ierr = VecCopy(w,y); CHKERRQ(ierr);
595       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
596       break;
597     }
598     count++;
599   }
600   theend2:
601   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
602   PetscFunctionReturn(0);
603 }
604 /* -------------------------------------------------------------------------- */
605 #undef __FUNC__
606 #define __FUNC__ "SNESSetLineSearch"
607 /*@C
608    SNESSetLineSearch - Sets the line search routine to be used
609    by the method SNES_EQ_LS.
610 
611    Input Parameters:
612 +  snes - nonlinear context obtained from SNESCreate()
613 -  func - pointer to int function
614 
615    Collective on SNES
616 
617    Available Routines:
618 +  SNESCubicLineSearch() - default line search
619 .  SNESQuadraticLineSearch() - quadratic line search
620 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
621 -  SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel)
622 
623     Options Database Keys:
624 +   -snes_eq_ls [basic,quadratic,cubic] - Selects line search
625 .   -snes_eq_ls_alpha <alpha> - Sets alpha
626 .   -snes_eq_ls_maxstep <max> - Sets maxstep
627 -   -snes_eq_ls_steptol <steptol> - Sets steptol
628 
629    Calling sequence of func:
630 .vb
631    func (SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
632          double fnorm, double *ynorm, double *gnorm, *flag)
633 .ve
634 
635     Input parameters for func:
636 +   snes - nonlinear context
637 .   x - current iterate
638 .   f - residual evaluated at x
639 .   y - search direction (contains new iterate on output)
640 .   w - work vector
641 -   fnorm - 2-norm of f
642 
643     Output parameters for func:
644 +   g - residual evaluated at new iterate y
645 .   y - new iterate (contains search direction on input)
646 .   gnorm - 2-norm of g
647 .   ynorm - 2-norm of search length
648 -   flag - set to 0 if the line search succeeds; a nonzero integer
649            on failure.
650 
651     Level: advanced
652 
653 .keywords: SNES, nonlinear, set, line search, routine
654 
655 .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
656 @*/
657 int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
658                              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