xref: /petsc/src/snes/impls/ls/ls.c (revision 3f1db9ec2fd39765c6c3a00831044586630c4cca)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: ls.c,v 1.117 1998/12/03 04:05:32 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, 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 .keywords: SNES, nonlinear, line search, cubic
226 
227 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
228           SNESSetLineSearch()
229 @*/
230 int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
231                      double fnorm, double *ynorm, double *gnorm,int *flag )
232 {
233   int    ierr;
234   Scalar mone = -1.0;
235 
236   PetscFunctionBegin;
237   *flag = 0;
238   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
239   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
240   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
241   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
242   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
243   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
244   PetscFunctionReturn(0);
245 }
246 /* -------------------------------------------------------------------------- */
247 #undef __FUNC__
248 #define __FUNC__ "SNESNoLineSearchNoNorms"
249 
250 /*@C
251    SNESNoLineSearchNoNorms - This routine is not a line search at
252    all; it simply uses the full Newton step. This version does not
253    even compute the norm of the function or search direction; this
254    is intended only when you know the full step is fine and are
255    not checking for convergence of the nonlinear iteration (for
256    example, you are running always for a fixed number of Newton steps).
257 
258    Collective on SNES and Vec
259 
260    Input Parameters:
261 +  snes - nonlinear context
262 .  x - current iterate
263 .  f - residual evaluated at x
264 .  y - search direction (contains new iterate on output)
265 .  w - work vector
266 -  fnorm - 2-norm of f
267 
268    Output Parameters:
269 +  g - residual evaluated at new iterate y
270 .  gnorm - not changed
271 .  ynorm - not changed
272 -  flag - set to 0, indicating a successful line search
273 
274    Options Database Key:
275 .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
276 
277 .keywords: SNES, nonlinear, line search, cubic
278 
279 .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
280           SNESSetLineSearch(), SNESNoLineSearch()
281 @*/
282 int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
283                      double fnorm, double *ynorm, double *gnorm,int *flag )
284 {
285   int    ierr;
286   Scalar mone = -1.0;
287 
288   PetscFunctionBegin;
289   *flag = 0;
290   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
291   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
292   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
293   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
294   PetscFunctionReturn(0);
295 }
296 /* -------------------------------------------------------------------------- */
297 #undef __FUNC__
298 #define __FUNC__ "SNESCubicLineSearch"
299 /*@C
300    SNESCubicLineSearch - Performs a cubic line search (default line search method).
301 
302    Collective on SNES
303 
304    Input Parameters:
305 +  snes - nonlinear context
306 .  x - current iterate
307 .  f - residual evaluated at x
308 .  y - search direction (contains new iterate on output)
309 .  w - work vector
310 -  fnorm - 2-norm of f
311 
312    Output Parameters:
313 +  g - residual evaluated at new iterate y
314 .  y - new iterate (contains search direction on input)
315 .  gnorm - 2-norm of g
316 .  ynorm - 2-norm of search length
317 -  flag - 0 if line search succeeds; -1 on failure.
318 
319    Options Database Key:
320 $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
321 
322    Notes:
323    This line search is taken from "Numerical Methods for Unconstrained
324    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
325 
326 .keywords: SNES, nonlinear, line search, cubic
327 
328 .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
329 @*/
330 int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
331                         double fnorm,double *ynorm,double *gnorm,int *flag)
332 {
333   /*
334      Note that for line search purposes we work with with the related
335      minimization problem:
336         min  z(x):  R^n -> R,
337      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
338    */
339 
340   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
341   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
342 #if defined(USE_PETSC_COMPLEX)
343   Scalar  cinitslope, clambda;
344 #endif
345   int     ierr, count;
346   SNES_LS *neP = (SNES_LS *) snes->data;
347   Scalar  mone = -1.0,scale;
348 
349   PetscFunctionBegin;
350   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
351   *flag   = 0;
352   alpha   = neP->alpha;
353   maxstep = neP->maxstep;
354   steptol = neP->steptol;
355 
356   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
357   if (*ynorm < snes->atol) {
358     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
359     *gnorm = fnorm;
360     ierr = VecCopy(x,y); CHKERRQ(ierr);
361     ierr = VecCopy(f,g); CHKERRQ(ierr);
362     goto theend1;
363   }
364   if (*ynorm > maxstep) {	/* Step too big, so scale back */
365     scale = maxstep/(*ynorm);
366 #if defined(USE_PETSC_COMPLEX)
367     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
368 #else
369     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
370 #endif
371     ierr = VecScale(&scale,y); CHKERRQ(ierr);
372     *ynorm = maxstep;
373   }
374   minlambda = steptol/(*ynorm);
375   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
376 #if defined(USE_PETSC_COMPLEX)
377   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
378   initslope = PetscReal(cinitslope);
379 #else
380   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
381 #endif
382   if (initslope > 0.0) initslope = -initslope;
383   if (initslope == 0.0) initslope = -1.0;
384 
385   ierr = VecCopy(y,w); CHKERRQ(ierr);
386   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
387   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
388   ierr = VecNorm(g,NORM_2,gnorm);
389   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
390     ierr = VecCopy(w,y); CHKERRQ(ierr);
391     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
392     goto theend1;
393   }
394 
395   /* Fit points with quadratic */
396   lambda = 1.0; count = 0;
397   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
398   lambdaprev = lambda;
399   gnormprev = *gnorm;
400   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
401   else lambda = lambdatemp;
402   ierr   = VecCopy(x,w); CHKERRQ(ierr);
403   lambdaneg = -lambda;
404 #if defined(USE_PETSC_COMPLEX)
405   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
406 #else
407   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
408 #endif
409   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
410   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
411   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
412     ierr = VecCopy(w,y); CHKERRQ(ierr);
413     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
414     goto theend1;
415   }
416 
417   /* Fit points with cubic */
418   count = 1;
419   while (1) {
420     if (lambda <= minlambda) { /* bad luck; use full step */
421       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
422       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
423                fnorm,*gnorm,*ynorm,lambda,initslope);
424       ierr = VecCopy(w,y); CHKERRQ(ierr);
425       *flag = -1; break;
426     }
427     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
428     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
429     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
430     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
431     d  = b*b - 3*a*initslope;
432     if (d < 0.0) d = 0.0;
433     if (a == 0.0) {
434       lambdatemp = -initslope/(2.0*b);
435     } else {
436       lambdatemp = (-b + sqrt(d))/(3.0*a);
437     }
438     if (lambdatemp > .5*lambda) {
439       lambdatemp = .5*lambda;
440     }
441     lambdaprev = lambda;
442     gnormprev = *gnorm;
443     if (lambdatemp <= .1*lambda) {
444       lambda = .1*lambda;
445     }
446     else lambda = lambdatemp;
447     ierr = VecCopy( x, w ); CHKERRQ(ierr);
448     lambdaneg = -lambda;
449 #if defined(USE_PETSC_COMPLEX)
450     clambda = lambdaneg;
451     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
452 #else
453     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
454 #endif
455     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
456     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
457     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
458       ierr = VecCopy(w,y); CHKERRQ(ierr);
459       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
460       goto theend1;
461     } else {
462       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
463     }
464     count++;
465   }
466   theend1:
467   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
468   PetscFunctionReturn(0);
469 }
470 /* -------------------------------------------------------------------------- */
471 #undef __FUNC__
472 #define __FUNC__ "SNESQuadraticLineSearch"
473 /*@C
474    SNESQuadraticLineSearch - Performs a quadratic line search.
475 
476    Collective on SNES and Vec
477 
478    Input Parameters:
479 +  snes - the SNES context
480 .  x - current iterate
481 .  f - residual evaluated at x
482 .  y - search direction (contains new iterate on output)
483 .  w - work vector
484 -  fnorm - 2-norm of f
485 
486    Output Parameters:
487 +  g - residual evaluated at new iterate y
488 .  y - new iterate (contains search direction on input)
489 .  gnorm - 2-norm of g
490 .  ynorm - 2-norm of search length
491 -  flag - 0 if line search succeeds; -1 on failure.
492 
493    Options Database Key:
494 .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
495 
496    Notes:
497    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
498 
499 .keywords: SNES, nonlinear, quadratic, line search
500 
501 .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
502 @*/
503 int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
504                            double fnorm, double *ynorm, double *gnorm,int *flag)
505 {
506   /*
507      Note that for line search purposes we work with with the related
508      minimization problem:
509         min  z(x):  R^n -> R,
510      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
511    */
512   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
513 #if defined(USE_PETSC_COMPLEX)
514   Scalar  cinitslope,clambda;
515 #endif
516   int     ierr,count;
517   SNES_LS *neP = (SNES_LS *) snes->data;
518   Scalar  mone = -1.0,scale;
519 
520   PetscFunctionBegin;
521   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
522   *flag   = 0;
523   alpha   = neP->alpha;
524   maxstep = neP->maxstep;
525   steptol = neP->steptol;
526 
527   VecNorm(y, NORM_2,ynorm );
528   if (*ynorm < snes->atol) {
529     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
530     *gnorm = fnorm;
531     ierr = VecCopy(x,y); CHKERRQ(ierr);
532     ierr = VecCopy(f,g); CHKERRQ(ierr);
533     goto theend2;
534   }
535   if (*ynorm > maxstep) {	/* Step too big, so scale back */
536     scale = maxstep/(*ynorm);
537     ierr = VecScale(&scale,y); CHKERRQ(ierr);
538     *ynorm = maxstep;
539   }
540   minlambda = steptol/(*ynorm);
541   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
542 #if defined(USE_PETSC_COMPLEX)
543   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
544   initslope = PetscReal(cinitslope);
545 #else
546   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
547 #endif
548   if (initslope > 0.0) initslope = -initslope;
549   if (initslope == 0.0) initslope = -1.0;
550 
551   ierr = VecCopy(y,w); CHKERRQ(ierr);
552   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
553   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
554   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
555   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
556     ierr = VecCopy(w,y); CHKERRQ(ierr);
557     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
558     goto theend2;
559   }
560 
561   /* Fit points with quadratic */
562   lambda = 1.0; count = 0;
563   count = 1;
564   while (1) {
565     if (lambda <= minlambda) { /* bad luck; use full step */
566       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
567       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
568                fnorm,*gnorm,*ynorm,lambda,initslope);
569       ierr = VecCopy(w,y); CHKERRQ(ierr);
570       *flag = -1; break;
571     }
572     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
573     if (lambdatemp <= .1*lambda) {
574       lambda = .1*lambda;
575     } else lambda = lambdatemp;
576     ierr = VecCopy(x,w); CHKERRQ(ierr);
577     lambda = -lambda;
578 #if defined(USE_PETSC_COMPLEX)
579     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
580 #else
581     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
582 #endif
583     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
584     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
585     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
586       ierr = VecCopy(w,y); CHKERRQ(ierr);
587       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
588       break;
589     }
590     count++;
591   }
592   theend2:
593   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
594   PetscFunctionReturn(0);
595 }
596 /* -------------------------------------------------------------------------- */
597 #undef __FUNC__
598 #define __FUNC__ "SNESSetLineSearch"
599 /*@C
600    SNESSetLineSearch - Sets the line search routine to be used
601    by the method SNES_EQ_LS.
602 
603    Input Parameters:
604 +  snes - nonlinear context obtained from SNESCreate()
605 -  func - pointer to int function
606 
607    Collective on SNES
608 
609    Available Routines:
610 +  SNESCubicLineSearch() - default line search
611 .  SNESQuadraticLineSearch() - quadratic line search
612 .  SNESNoLineSearch() - the full Newton step (actually not a line search)
613 -  SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel)
614 
615     Options Database Keys:
616 +   -snes_eq_ls [basic,quadratic,cubic] - Selects line search
617 .   -snes_eq_ls_alpha <alpha> - Sets alpha
618 .   -snes_eq_ls_maxstep <max> - Sets maxstep
619 -   -snes_eq_ls_steptol <steptol> - Sets steptol
620 
621    Calling sequence of func:
622 .vb
623    func (SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
624          double fnorm, double *ynorm, double *gnorm, *flag)
625 .ve
626 
627     Input parameters for func:
628 +   snes - nonlinear context
629 .   x - current iterate
630 .   f - residual evaluated at x
631 .   y - search direction (contains new iterate on output)
632 .   w - work vector
633 -   fnorm - 2-norm of f
634 
635     Output parameters for func:
636 +   g - residual evaluated at new iterate y
637 .   y - new iterate (contains search direction on input)
638 .   gnorm - 2-norm of g
639 .   ynorm - 2-norm of search length
640 -   flag - set to 0 if the line search succeeds; a nonzero integer
641            on failure.
642 
643 .keywords: SNES, nonlinear, set, line search, routine
644 
645 .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
646 @*/
647 int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
648                              double,double*,double*,int*))
649 {
650   int ierr, (*f)(SNES,int (*)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*));
651 
652   PetscFunctionBegin;
653   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
654   if (f) {
655     ierr = (*f)(snes,func);CHKERRQ(ierr);
656   }
657   PetscFunctionReturn(0);
658 }
659 /* -------------------------------------------------------------------------- */
660 EXTERN_C_BEGIN
661 #undef __FUNC__
662 #define __FUNC__ "SNESSetLineSearch_LS"
663 int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
664                          double,double*,double*,int*))
665 {
666   PetscFunctionBegin;
667   ((SNES_LS *)(snes->data))->LineSearch = func;
668   PetscFunctionReturn(0);
669 }
670 EXTERN_C_END
671 /* -------------------------------------------------------------------------- */
672 /*
673    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
674 
675    Input Parameter:
676 .  snes - the SNES context
677 
678    Application Interface Routine: SNESPrintHelp()
679 */
680 #undef __FUNC__
681 #define __FUNC__ "SNESPrintHelp_EQ_LS"
682 static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
683 {
684   SNES_LS *ls = (SNES_LS *)snes->data;
685 
686   PetscFunctionBegin;
687   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
688   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
689   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
690   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
691   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
692   PetscFunctionReturn(0);
693 }
694 /* -------------------------------------------------------------------------- */
695 /*
696    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
697 
698    Input Parameters:
699 .  SNES - the SNES context
700 .  viewer - visualization context
701 
702    Application Interface Routine: SNESView()
703 */
704 #undef __FUNC__
705 #define __FUNC__ "SNESView_EQ_LS"
706 static int SNESView_EQ_LS(SNES snes,Viewer viewer)
707 {
708   SNES_LS    *ls = (SNES_LS *)snes->data;
709   FILE       *fd;
710   char       *cstr;
711   int        ierr;
712   ViewerType vtype;
713 
714   PetscFunctionBegin;
715   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
716   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
717     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
718     if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
719     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
720     else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
721     else cstr = "unknown";
722     PetscFPrintf(snes->comm,fd,"    line search variant: %s\n",cstr);
723     PetscFPrintf(snes->comm,fd,"    alpha=%g, maxstep=%g, steptol=%g\n",
724                  ls->alpha,ls->maxstep,ls->steptol);
725   } else {
726     SETERRQ(1,1,"Viewer type not supported for this object");
727   }
728   PetscFunctionReturn(0);
729 }
730 /* -------------------------------------------------------------------------- */
731 /*
732    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
733 
734    Input Parameter:
735 .  snes - the SNES context
736 
737    Application Interface Routine: SNESSetFromOptions()
738 */
739 #undef __FUNC__
740 #define __FUNC__ "SNESSetFromOptions_EQ_LS"
741 static int SNESSetFromOptions_EQ_LS(SNES snes)
742 {
743   SNES_LS *ls = (SNES_LS *)snes->data;
744   char    ver[16];
745   double  tmp;
746   int     ierr,flg;
747 
748   PetscFunctionBegin;
749   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
750   if (flg) {
751     ls->alpha = tmp;
752   }
753   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
754   if (flg) {
755     ls->maxstep = tmp;
756   }
757   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
758   if (flg) {
759     ls->steptol = tmp;
760   }
761   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
762   if (flg) {
763     if (!PetscStrcmp(ver,"basic")) {
764       ierr = SNESSetLineSearch(snes,SNESNoLineSearch);CHKERRQ(ierr);
765     } else if (!PetscStrcmp(ver,"basicnonorms")) {
766       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);CHKERRQ(ierr);
767     } else if (!PetscStrcmp(ver,"quadratic")) {
768       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch);CHKERRQ(ierr);
769     } else if (!PetscStrcmp(ver,"cubic")) {
770       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch);CHKERRQ(ierr);
771     }
772     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
773   }
774   PetscFunctionReturn(0);
775 }
776 /* -------------------------------------------------------------------------- */
777 /*
778    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
779    SNES_LS, and sets this as the private data within the generic nonlinear solver
780    context, SNES, that was created within SNESCreate().
781 
782 
783    Input Parameter:
784 .  snes - the SNES context
785 
786    Application Interface Routine: SNESCreate()
787  */
788 EXTERN_C_BEGIN
789 #undef __FUNC__
790 #define __FUNC__ "SNESCreate_EQ_LS"
791 int SNESCreate_EQ_LS(SNES snes)
792 {
793   int     ierr;
794   SNES_LS *neP;
795 
796   PetscFunctionBegin;
797   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
798     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
799   }
800 
801   snes->setup		= SNESSetUp_EQ_LS;
802   snes->solve		= SNESSolve_EQ_LS;
803   snes->destroy		= SNESDestroy_EQ_LS;
804   snes->converged	= SNESConverged_EQ_LS;
805   snes->printhelp       = SNESPrintHelp_EQ_LS;
806   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
807   snes->view            = SNESView_EQ_LS;
808   snes->nwork           = 0;
809 
810   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
811   PLogObjectMemory(snes,sizeof(SNES_LS));
812   snes->data    	= (void *) neP;
813   neP->alpha		= 1.e-4;
814   neP->maxstep		= 1.e8;
815   neP->steptol		= 1.e-12;
816   neP->LineSearch       = SNESCubicLineSearch;
817 
818   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
819                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
820 
821   PetscFunctionReturn(0);
822 }
823 EXTERN_C_END
824 
825 
826 
827