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