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