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