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