xref: /petsc/src/snes/impls/ls/ls.c (revision 3a40ed3dce77c081171d005ae1a6ff4bb9d13b6f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*3a40ed3dSBarry Smith static char vcid[] = "$Id: ls.c,v 1.95 1997/08/22 15:17:59 bsmith Exp bsmith $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
55e76c431SBarry Smith #include <math.h>
670f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
7b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
85e42470aSBarry Smith 
95e42470aSBarry Smith /*
105e42470aSBarry Smith      Implements a line search variant of Newton's Method
115e76c431SBarry Smith     for solving systems of nonlinear equations.
125e76c431SBarry Smith 
135e76c431SBarry Smith     Input parameters:
145e42470aSBarry Smith .   snes - nonlinear context obtained from SNESCreate()
155e76c431SBarry Smith 
165e42470aSBarry Smith     Output Parameters:
17a935fc98SLois Curfman McInnes .   outits  - Number of global iterations until termination.
185e76c431SBarry Smith 
195e76c431SBarry Smith     Notes:
205e76c431SBarry Smith     This implements essentially a truncated Newton method with a
215e76c431SBarry Smith     line search.  By default a cubic backtracking line search
225e76c431SBarry Smith     is employed, as described in the text "Numerical Methods for
235e76c431SBarry Smith     Unconstrained Optimization and Nonlinear Equations" by Dennis
24393d2d9aSLois Curfman McInnes     and Schnabel.
255e76c431SBarry Smith */
265e76c431SBarry Smith 
275615d1e5SSatish Balay #undef __FUNC__
285615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS"
29f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
305e76c431SBarry Smith {
315e42470aSBarry Smith   SNES_LS       *neP = (SNES_LS *) snes->data;
32761aaf1bSLois Curfman McInnes   int           maxits, i, history_len, ierr, lits, lsfail;
33112a2221SBarry Smith   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
346b5873e3SBarry Smith   double        fnorm, gnorm, xnorm, ynorm, *history;
355e42470aSBarry Smith   Vec           Y, X, F, G, W, TMP;
365e76c431SBarry Smith 
37*3a40ed3dSBarry Smith   PetscFunctionBegin;
385e42470aSBarry Smith   history	= snes->conv_hist;	/* convergence history */
3951979daaSLois Curfman McInnes   history_len	= snes->conv_hist_size;	/* convergence history length */
405e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
415e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
4239e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
435e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
445e42470aSBarry Smith   G		= snes->work[1];
455e42470aSBarry Smith   W		= snes->work[2];
465e76c431SBarry Smith 
4725ed9cd1SBarry Smith   snes->iter = 0;
485334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
49cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
505e42470aSBarry Smith   snes->norm = fnorm;
5151979daaSLois Curfman McInnes   if (history) history[0] = fnorm;
5294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
535e76c431SBarry Smith 
54*3a40ed3dSBarry Smith   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
5594a9d846SBarry Smith 
56d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
57d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
58d034289bSBarry Smith 
595e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
605e42470aSBarry Smith     snes->iter = i+1;
615e76c431SBarry Smith 
62ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
635334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
645334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
6578b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
667a00f4a9SLois Curfman McInnes     snes->linear_its += PetscAbsInt(lits);
67af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
68ea4d3ed3SLois Curfman McInnes 
69ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
70ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
71ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
72ea4d3ed3SLois Curfman McInnes     */
7381b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
74ddd12767SBarry Smith     ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
75af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
76761aaf1bSLois Curfman McInnes     if (lsfail) snes->nfailures++;
775e76c431SBarry Smith 
7839e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
7939e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
805e76c431SBarry Smith     fnorm = gnorm;
815e76c431SBarry Smith 
825e42470aSBarry Smith     snes->norm = fnorm;
835e76c431SBarry Smith     if (history && history_len > i+1) history[i+1] = fnorm;
8494a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
855e76c431SBarry Smith 
865e76c431SBarry Smith     /* Test for convergence */
8729e0b56fSBarry Smith     if (snes->converged) {
8829e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
89bbb6d6a8SBarry Smith       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
9016c95cacSBarry Smith         break;
9116c95cacSBarry Smith       }
9216c95cacSBarry Smith     }
9329e0b56fSBarry Smith   }
9439e2f89bSBarry Smith   if (X != snes->vec_sol) {
95393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
9639e2f89bSBarry Smith     snes->vec_sol_always  = snes->vec_sol;
9739e2f89bSBarry Smith     snes->vec_func_always = snes->vec_func;
9839e2f89bSBarry Smith   }
9952392280SLois Curfman McInnes   if (i == maxits) {
10094a424c1SBarry Smith     PLogInfo(snes,
101f63b844aSLois Curfman McInnes       "SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
10252392280SLois Curfman McInnes     i--;
10352392280SLois Curfman McInnes   }
10451979daaSLois Curfman McInnes   if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1;
1055e42470aSBarry Smith   *outits = i+1;
106*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
1075e76c431SBarry Smith }
1085e76c431SBarry Smith /* ------------------------------------------------------------ */
1095615d1e5SSatish Balay #undef __FUNC__
1105615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
111f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes )
1125e76c431SBarry Smith {
1135e42470aSBarry Smith   int ierr;
114*3a40ed3dSBarry Smith 
115*3a40ed3dSBarry Smith   PetscFunctionBegin;
11681b6cf68SLois Curfman McInnes   snes->nwork = 4;
117d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1185e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
11981b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
120*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
1215e76c431SBarry Smith }
1225e76c431SBarry Smith /* ------------------------------------------------------------ */
1235615d1e5SSatish Balay #undef __FUNC__
124d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
125f63b844aSLois Curfman McInnes int SNESDestroy_EQ_LS(PetscObject obj)
1265e76c431SBarry Smith {
1275e42470aSBarry Smith   SNES snes = (SNES) obj;
128393d2d9aSLois Curfman McInnes   int  ierr;
129*3a40ed3dSBarry Smith 
130*3a40ed3dSBarry Smith   PetscFunctionBegin;
1315baf8537SBarry Smith   if (snes->nwork) {
1324b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
13321c89e3eSBarry Smith   }
1340452661fSBarry Smith   PetscFree(snes->data);
135*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
1365e76c431SBarry Smith }
1375e76c431SBarry Smith /* ------------------------------------------------------------ */
1385615d1e5SSatish Balay #undef __FUNC__
1395615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
1405e76c431SBarry Smith /*ARGSUSED*/
1414b828684SBarry Smith /*@C
1425e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
1435e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
1445e76c431SBarry Smith    to serve as a template and is not recommended for general use.
1455e76c431SBarry Smith 
1465e76c431SBarry Smith    Input Parameters:
1475e42470aSBarry Smith .  snes - nonlinear context
1485e76c431SBarry Smith .  x - current iterate
1495e76c431SBarry Smith .  f - residual evaluated at x
1505e76c431SBarry Smith .  y - search direction (contains new iterate on output)
1515e76c431SBarry Smith .  w - work vector
1525e76c431SBarry Smith .  fnorm - 2-norm of f
1535e76c431SBarry Smith 
154c4a48953SLois Curfman McInnes    Output Parameters:
1555e76c431SBarry Smith .  g - residual evaluated at new iterate y
1565e76c431SBarry Smith .  y - new iterate (contains search direction on input)
1575e76c431SBarry Smith .  gnorm - 2-norm of g
1585e76c431SBarry Smith .  ynorm - 2-norm of search length
159761aaf1bSLois Curfman McInnes .  flag - set to 0, indicating a successful line search
1605e76c431SBarry Smith 
161c4a48953SLois Curfman McInnes    Options Database Key:
16209d61ba7SLois Curfman McInnes $  -snes_eq_ls basic
163c4a48953SLois Curfman McInnes 
16428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
16528ae5a21SLois Curfman McInnes 
166f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
16777c4ece6SBarry Smith           SNESSetLineSearch()
1685e76c431SBarry Smith @*/
1695e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
170761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
1715e76c431SBarry Smith {
1725e42470aSBarry Smith   int    ierr;
1735334005bSBarry Smith   Scalar mone = -1.0;
1745334005bSBarry Smith 
175*3a40ed3dSBarry Smith   PetscFunctionBegin;
176761aaf1bSLois Curfman McInnes   *flag = 0;
1777857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
178cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
179ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
180ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
181cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
1827857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
183*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
1845e76c431SBarry Smith }
1855e76c431SBarry Smith /* ------------------------------------------------------------------ */
1865615d1e5SSatish Balay #undef __FUNC__
18729e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
18829e0b56fSBarry Smith /*ARGSUSED*/
18929e0b56fSBarry Smith /*@C
19029e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
19129e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
19229e0b56fSBarry Smith    even compute the norm of the function or search direction; this
19329e0b56fSBarry Smith    is intended only when you know the full step is fine and are
19429e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
19529e0b56fSBarry Smith    example, you are running always for a fixed number of Newton
19629e0b56fSBarry Smith    steps).
19729e0b56fSBarry Smith 
19829e0b56fSBarry Smith    Input Parameters:
19929e0b56fSBarry Smith .  snes - nonlinear context
20029e0b56fSBarry Smith .  x - current iterate
20129e0b56fSBarry Smith .  f - residual evaluated at x
20229e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
20329e0b56fSBarry Smith .  w - work vector
20429e0b56fSBarry Smith .  fnorm - 2-norm of f
20529e0b56fSBarry Smith 
20629e0b56fSBarry Smith    Output Parameters:
20729e0b56fSBarry Smith .  g - residual evaluated at new iterate y
20829e0b56fSBarry Smith .  gnorm - not changed
20929e0b56fSBarry Smith .  ynorm - not changed
21029e0b56fSBarry Smith .  flag - set to 0, indicating a successful line search
21129e0b56fSBarry Smith 
21229e0b56fSBarry Smith    Options Database Key:
21329e0b56fSBarry Smith $  -snes_eq_ls basicnonorms
21429e0b56fSBarry Smith 
21529e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
21629e0b56fSBarry Smith 
21729e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
21829e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
21929e0b56fSBarry Smith @*/
22029e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
22129e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
22229e0b56fSBarry Smith {
22329e0b56fSBarry Smith   int    ierr;
22429e0b56fSBarry Smith   Scalar mone = -1.0;
22529e0b56fSBarry Smith 
226*3a40ed3dSBarry Smith   PetscFunctionBegin;
22729e0b56fSBarry Smith   *flag = 0;
22829e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
22929e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
23029e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
23129e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
232*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
23329e0b56fSBarry Smith }
23429e0b56fSBarry Smith /* ------------------------------------------------------------------ */
23529e0b56fSBarry Smith #undef __FUNC__
2365615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
2374b828684SBarry Smith /*@C
238f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
2395e76c431SBarry Smith 
2405e76c431SBarry Smith    Input Parameters:
2415e42470aSBarry Smith .  snes - nonlinear context
2425e76c431SBarry Smith .  x - current iterate
2435e76c431SBarry Smith .  f - residual evaluated at x
2445e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2455e76c431SBarry Smith .  w - work vector
2465e76c431SBarry Smith .  fnorm - 2-norm of f
2475e76c431SBarry Smith 
248393d2d9aSLois Curfman McInnes    Output Parameters:
2495e76c431SBarry Smith .  g - residual evaluated at new iterate y
2505e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2515e76c431SBarry Smith .  gnorm - 2-norm of g
2525e76c431SBarry Smith .  ynorm - 2-norm of search length
253761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
2545e76c431SBarry Smith 
255c4a48953SLois Curfman McInnes    Options Database Key:
25609d61ba7SLois Curfman McInnes $  -snes_eq_ls cubic
257c4a48953SLois Curfman McInnes 
2585e76c431SBarry Smith    Notes:
2595e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
2605e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
2615e76c431SBarry Smith 
26228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
26328ae5a21SLois Curfman McInnes 
26477c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
2655e76c431SBarry Smith @*/
2665e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
267761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
2685e76c431SBarry Smith {
269406556e6SLois Curfman McInnes   /*
270406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
271406556e6SLois Curfman McInnes      minimization problem:
272406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
273406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
274406556e6SLois Curfman McInnes    */
275406556e6SLois Curfman McInnes 
276ddd12767SBarry Smith   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
277ea4d3ed3SLois Curfman McInnes   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
278*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
2795e42470aSBarry Smith   Scalar  cinitslope, clambda;
2806b5873e3SBarry Smith #endif
2815e42470aSBarry Smith   int     ierr, count;
2825e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
2835334005bSBarry Smith   Scalar  mone = -1.0,scale;
2845e76c431SBarry Smith 
285*3a40ed3dSBarry Smith   PetscFunctionBegin;
2867857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
287761aaf1bSLois Curfman McInnes   *flag   = 0;
2885e76c431SBarry Smith   alpha   = neP->alpha;
2895e76c431SBarry Smith   maxstep = neP->maxstep;
2905e76c431SBarry Smith   steptol = neP->steptol;
2915e76c431SBarry Smith 
292cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
293a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
294a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
295a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
296a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
297a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
298ad922ac9SBarry Smith     goto theend1;
29994a9d846SBarry Smith   }
3005e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3015e42470aSBarry Smith     scale = maxstep/(*ynorm);
302*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
30394a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",real(scale));
3046b5873e3SBarry Smith #else
30594a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
3066b5873e3SBarry Smith #endif
307393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
3085e76c431SBarry Smith     *ynorm = maxstep;
3095e76c431SBarry Smith   }
3105e76c431SBarry Smith   minlambda = steptol/(*ynorm);
311a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
312*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
313a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
3145e42470aSBarry Smith   initslope = real(cinitslope);
3155e42470aSBarry Smith #else
316a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
3175e42470aSBarry Smith #endif
3185e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3195e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3205e76c431SBarry Smith 
321393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
3225334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
32378b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
324cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);
325406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
326393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
32794a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
328ad922ac9SBarry Smith     goto theend1;
3295e76c431SBarry Smith   }
3305e76c431SBarry Smith 
3315e76c431SBarry Smith   /* Fit points with quadratic */
3325e76c431SBarry Smith   lambda = 1.0; count = 0;
333a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
3345e76c431SBarry Smith   lambdaprev = lambda;
3355e76c431SBarry Smith   gnormprev = *gnorm;
336ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
337ddd12767SBarry Smith   else lambda = lambdatemp;
338393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w); CHKERRQ(ierr);
339ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
340*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
341ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
3425e42470aSBarry Smith #else
343ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
3445e42470aSBarry Smith #endif
34578b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
346cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
347406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
348393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
349ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
350ad922ac9SBarry Smith     goto theend1;
3515e76c431SBarry Smith   }
3525e76c431SBarry Smith 
3535e76c431SBarry Smith   /* Fit points with cubic */
3545e76c431SBarry Smith   count = 1;
3555e76c431SBarry Smith   while (1) {
3565e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
3572b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
3582b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
359ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
360393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
361761aaf1bSLois Curfman McInnes       *flag = -1; break;
3625e76c431SBarry Smith     }
363406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
364406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
365ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
3662b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
3675e76c431SBarry Smith     d  = b*b - 3*a*initslope;
3685e76c431SBarry Smith     if (d < 0.0) d = 0.0;
3695e76c431SBarry Smith     if (a == 0.0) {
3705e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
3715e76c431SBarry Smith     } else {
3725e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
3735e76c431SBarry Smith     }
3745e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
3755e76c431SBarry Smith       lambdatemp = .5*lambda;
3765e76c431SBarry Smith     }
3775e76c431SBarry Smith     lambdaprev = lambda;
3785e76c431SBarry Smith     gnormprev = *gnorm;
3795e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
3805e76c431SBarry Smith       lambda = .1*lambda;
3815e76c431SBarry Smith     }
3825e76c431SBarry Smith     else lambda = lambdatemp;
383393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w ); CHKERRQ(ierr);
384ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
385*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
386ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
387393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
3885e42470aSBarry Smith #else
389ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
3905e42470aSBarry Smith #endif
39178b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
392cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
393406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
394393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
395ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
3962715565aSLois Curfman McInnes       goto theend1;
3972b022350SLois Curfman McInnes     } else {
3982b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
3995e76c431SBarry Smith     }
4005e76c431SBarry Smith     count++;
4015e76c431SBarry Smith   }
402ad922ac9SBarry Smith   theend1:
4037857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
404*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
4055e76c431SBarry Smith }
40652392280SLois Curfman McInnes /* ------------------------------------------------------------------ */
4075615d1e5SSatish Balay #undef __FUNC__
4085615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
4094b828684SBarry Smith /*@C
410f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
4115e76c431SBarry Smith 
4125e42470aSBarry Smith    Input Parameters:
413f59ffdedSLois Curfman McInnes .  snes - the SNES context
4145e42470aSBarry Smith .  x - current iterate
4155e42470aSBarry Smith .  f - residual evaluated at x
4165e42470aSBarry Smith .  y - search direction (contains new iterate on output)
4175e42470aSBarry Smith .  w - work vector
4185e42470aSBarry Smith .  fnorm - 2-norm of f
4195e42470aSBarry Smith 
420c4a48953SLois Curfman McInnes    Output Parameters:
4215e42470aSBarry Smith .  g - residual evaluated at new iterate y
4225e42470aSBarry Smith .  y - new iterate (contains search direction on input)
4235e42470aSBarry Smith .  gnorm - 2-norm of g
4245e42470aSBarry Smith .  ynorm - 2-norm of search length
425761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
4265e42470aSBarry Smith 
427c4a48953SLois Curfman McInnes    Options Database Key:
42809d61ba7SLois Curfman McInnes $  -snes_eq_ls quadratic
429c4a48953SLois Curfman McInnes 
4305e42470aSBarry Smith    Notes:
43177c4ece6SBarry Smith    Use SNESSetLineSearch()
432f63b844aSLois Curfman McInnes    to set this routine within the SNES_EQ_LS method.
4335e42470aSBarry Smith 
434f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
435f59ffdedSLois Curfman McInnes 
43677c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
4375e42470aSBarry Smith @*/
4385e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
439761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
4405e76c431SBarry Smith {
441406556e6SLois Curfman McInnes   /*
442406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
443406556e6SLois Curfman McInnes      minimization problem:
444406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
445406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
446406556e6SLois Curfman McInnes    */
447ddd12767SBarry Smith   double  steptol,initslope,lambdaprev,gnormprev,maxstep,minlambda,alpha,lambda,lambdatemp;
448*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
4495e42470aSBarry Smith   Scalar  cinitslope,clambda;
4506b5873e3SBarry Smith #endif
4515e42470aSBarry Smith   int     ierr,count;
4525e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
4535334005bSBarry Smith   Scalar  mone = -1.0,scale;
4545e76c431SBarry Smith 
455*3a40ed3dSBarry Smith   PetscFunctionBegin;
4567857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
457761aaf1bSLois Curfman McInnes   *flag   = 0;
4585e42470aSBarry Smith   alpha   = neP->alpha;
4595e42470aSBarry Smith   maxstep = neP->maxstep;
4605e42470aSBarry Smith   steptol = neP->steptol;
4615e76c431SBarry Smith 
462cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
463b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
46494a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
465b37302c6SLois Curfman McInnes     *gnorm = fnorm;
466b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
467b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
468ad922ac9SBarry Smith     goto theend2;
46994a9d846SBarry Smith   }
4705e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4715e42470aSBarry Smith     scale = maxstep/(*ynorm);
472393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
4735e42470aSBarry Smith     *ynorm = maxstep;
4745e76c431SBarry Smith   }
4755e42470aSBarry Smith   minlambda = steptol/(*ynorm);
476a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
477*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
478a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
4795e42470aSBarry Smith   initslope = real(cinitslope);
4805e42470aSBarry Smith #else
481a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
4825e42470aSBarry Smith #endif
4835e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
4845e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
4855e42470aSBarry Smith 
486393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
4875334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
48878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
489cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
490406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
491393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
49294a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
493ad922ac9SBarry Smith     goto theend2;
4945e42470aSBarry Smith   }
4955e42470aSBarry Smith 
4965e42470aSBarry Smith   /* Fit points with quadratic */
4975e42470aSBarry Smith   lambda = 1.0; count = 0;
4985e42470aSBarry Smith   count = 1;
4995e42470aSBarry Smith   while (1) {
5005e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
50194a424c1SBarry Smith       PLogInfo(snes,
5024b828684SBarry Smith           "SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
50394a424c1SBarry Smith       PLogInfo(snes,
504ea4d3ed3SLois Curfman McInnes       "SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
505ea4d3ed3SLois Curfman McInnes           fnorm,*gnorm,*ynorm,lambda,initslope);
506393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
507761aaf1bSLois Curfman McInnes       *flag = -1; break;
5085e42470aSBarry Smith     }
509a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5105e42470aSBarry Smith     lambdaprev = lambda;
5115e42470aSBarry Smith     gnormprev = *gnorm;
5125e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
5135e42470aSBarry Smith       lambda = .1*lambda;
5145e42470aSBarry Smith     } else lambda = lambdatemp;
515393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w); CHKERRQ(ierr);
5165334005bSBarry Smith     lambda = -lambda;
517*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
518393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
5195e42470aSBarry Smith #else
520393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
5215e42470aSBarry Smith #endif
52278b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
523cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
524406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
525393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
52694a424c1SBarry Smith       PLogInfo(snes,
527ea4d3ed3SLois Curfman McInnes         "SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
52806259719SBarry Smith       break;
5295e42470aSBarry Smith     }
5305e42470aSBarry Smith     count++;
5315e42470aSBarry Smith   }
532ad922ac9SBarry Smith   theend2:
5337857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
534*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
5355e76c431SBarry Smith }
5365e76c431SBarry Smith /* ------------------------------------------------------------ */
5375615d1e5SSatish Balay #undef __FUNC__
538d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
539c9e6a524SLois Curfman McInnes /*@C
54077c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
541f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
5425e76c431SBarry Smith 
5435e76c431SBarry Smith    Input Parameters:
544eafb4bcbSBarry Smith .  snes - nonlinear context obtained from SNESCreate()
5455e76c431SBarry Smith .  func - pointer to int function
5465e76c431SBarry Smith 
547c4a48953SLois Curfman McInnes    Available Routines:
548f59ffdedSLois Curfman McInnes .  SNESCubicLineSearch() - default line search
549f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
550f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
5515e76c431SBarry Smith 
552c4a48953SLois Curfman McInnes     Options Database Keys:
55309d61ba7SLois Curfman McInnes $   -snes_eq_ls [basic,quadratic,cubic]
554912ebf9aSLois Curfman McInnes $   -snes_eq_ls_alpha <alpha>
555912ebf9aSLois Curfman McInnes $   -snes_eq_ls_maxstep <max>
556912ebf9aSLois Curfman McInnes $   -snes_eq_ls_steptol <steptol>
557c4a48953SLois Curfman McInnes 
5585e76c431SBarry Smith    Calling sequence of func:
559f59ffdedSLois Curfman McInnes    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
560761aaf1bSLois Curfman McInnes          Vec w, double fnorm, double *ynorm,
561761aaf1bSLois Curfman McInnes          double *gnorm, *flag)
5625e76c431SBarry Smith 
5635e76c431SBarry Smith     Input parameters for func:
5645e42470aSBarry Smith .   snes - nonlinear context
5655e76c431SBarry Smith .   x - current iterate
5665e76c431SBarry Smith .   f - residual evaluated at x
5675e76c431SBarry Smith .   y - search direction (contains new iterate on output)
5685e76c431SBarry Smith .   w - work vector
5695e76c431SBarry Smith .   fnorm - 2-norm of f
5705e76c431SBarry Smith 
5715e76c431SBarry Smith     Output parameters for func:
5725e76c431SBarry Smith .   g - residual evaluated at new iterate y
5735e76c431SBarry Smith .   y - new iterate (contains search direction on input)
5745e76c431SBarry Smith .   gnorm - 2-norm of g
5755e76c431SBarry Smith .   ynorm - 2-norm of search length
576761aaf1bSLois Curfman McInnes .   flag - set to 0 if the line search succeeds; a nonzero integer
577761aaf1bSLois Curfman McInnes            on failure.
578f59ffdedSLois Curfman McInnes 
579f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
580f59ffdedSLois Curfman McInnes 
581f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
5825e76c431SBarry Smith @*/
58377c4ece6SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
584761aaf1bSLois Curfman McInnes                              double,double*,double*,int*))
5855e76c431SBarry Smith {
586*3a40ed3dSBarry Smith   PetscFunctionBegin;
587f63b844aSLois Curfman McInnes   if ((snes)->type == SNES_EQ_LS) ((SNES_LS *)(snes->data))->LineSearch = func;
588*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
5895e76c431SBarry Smith }
59052392280SLois Curfman McInnes /* ------------------------------------------------------------------ */
5915615d1e5SSatish Balay #undef __FUNC__
592d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
593f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
5945e42470aSBarry Smith {
5955e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
5966daaf66cSBarry Smith 
597*3a40ed3dSBarry Smith   PetscFunctionBegin;
598f63b844aSLois Curfman McInnes   PetscPrintf(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
59909d61ba7SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
60009d61ba7SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
60109d61ba7SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
60209d61ba7SLois Curfman McInnes   PetscPrintf(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
603*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
6045e42470aSBarry Smith }
60552392280SLois Curfman McInnes /* ------------------------------------------------------------------ */
6065615d1e5SSatish Balay #undef __FUNC__
607d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
608f63b844aSLois Curfman McInnes static int SNESView_EQ_LS(PetscObject obj,Viewer viewer)
609a935fc98SLois Curfman McInnes {
610a935fc98SLois Curfman McInnes   SNES       snes = (SNES)obj;
611a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
612d636dbe3SBarry Smith   FILE       *fd;
61319bcc07fSBarry Smith   char       *cstr;
61451695f54SLois Curfman McInnes   int        ierr;
61519bcc07fSBarry Smith   ViewerType vtype;
616a935fc98SLois Curfman McInnes 
617*3a40ed3dSBarry Smith   PetscFunctionBegin;
61819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
61919bcc07fSBarry Smith   if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
62090ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
62119bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
62219bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
62319bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
62419bcc07fSBarry Smith     else cstr = "unknown";
62577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"    line search variant: %s\n",cstr);
62677c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"    alpha=%g, maxstep=%g, steptol=%g\n",
627a935fc98SLois Curfman McInnes                  ls->alpha,ls->maxstep,ls->steptol);
62819bcc07fSBarry Smith   }
629*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
630a935fc98SLois Curfman McInnes }
63152392280SLois Curfman McInnes /* ------------------------------------------------------------------ */
6325615d1e5SSatish Balay #undef __FUNC__
6335615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
634f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
6355e42470aSBarry Smith {
6365e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
6375e42470aSBarry Smith   char    ver[16];
6385e42470aSBarry Smith   double  tmp;
63925cdf11fSBarry Smith   int     ierr,flg;
6405e42470aSBarry Smith 
641*3a40ed3dSBarry Smith   PetscFunctionBegin;
64209d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
64325cdf11fSBarry Smith   if (flg) {
6445e42470aSBarry Smith     ls->alpha = tmp;
6455e42470aSBarry Smith   }
64609d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
64725cdf11fSBarry Smith   if (flg) {
6485e42470aSBarry Smith     ls->maxstep = tmp;
6495e42470aSBarry Smith   }
65009d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
65125cdf11fSBarry Smith   if (flg) {
6525e42470aSBarry Smith     ls->steptol = tmp;
6535e42470aSBarry Smith   }
65409d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
65525cdf11fSBarry Smith   if (flg) {
65648d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
65777c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESNoLineSearch);
6585e42470aSBarry Smith     }
65929e0b56fSBarry Smith     else if (!PetscStrcmp(ver,"basicnonorms")) {
66029e0b56fSBarry Smith       SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);
66129e0b56fSBarry Smith     }
66248d91487SBarry Smith     else if (!PetscStrcmp(ver,"quadratic")) {
66377c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESQuadraticLineSearch);
6645e42470aSBarry Smith     }
66548d91487SBarry Smith     else if (!PetscStrcmp(ver,"cubic")) {
66677c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESCubicLineSearch);
6675e42470aSBarry Smith     }
668e3372554SBarry Smith     else {SETERRQ(1,0,"Unknown line search");}
6695e42470aSBarry Smith   }
670*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
6715e42470aSBarry Smith }
6725e42470aSBarry Smith /* ------------------------------------------------------------ */
6735615d1e5SSatish Balay #undef __FUNC__
6745615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
675f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES  snes )
6765e42470aSBarry Smith {
6775e42470aSBarry Smith   SNES_LS *neP;
6785e42470aSBarry Smith 
679*3a40ed3dSBarry Smith   PetscFunctionBegin;
68029e0b56fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only");
681f63b844aSLois Curfman McInnes   snes->type		= SNES_EQ_LS;
682f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
683f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
684f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
685f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
686f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
687f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
688f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
6895baf8537SBarry Smith   snes->nwork           = 0;
6905e42470aSBarry Smith 
6910452661fSBarry Smith   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
692ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
6935e42470aSBarry Smith   snes->data    	= (void *) neP;
6945e42470aSBarry Smith   neP->alpha		= 1.e-4;
6955e42470aSBarry Smith   neP->maxstep		= 1.e8;
6965e42470aSBarry Smith   neP->steptol		= 1.e-12;
6975e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
698*3a40ed3dSBarry Smith   PetscFunctionReturn(0);
6995e42470aSBarry Smith }
7005e42470aSBarry Smith 
701