xref: /petsc/src/snes/impls/ls/ls.c (revision 69e4de806daeac5cd3538e86d1ec7d6cbbb057ff)
15e76c431SBarry Smith #ifndef lint
2*69e4de80SBarry Smith static char vcid[] = "$Id: ls.c,v 1.19 1995/05/25 22:48:47 bsmith Exp bsmith $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
55e76c431SBarry Smith #include <math.h>
6fbe28522SBarry Smith #include "ls.h"
75e42470aSBarry Smith 
85e42470aSBarry Smith /*
95e42470aSBarry Smith      Implements a line search variant of Newton's Method
105e76c431SBarry Smith     for solving systems of nonlinear equations.
115e76c431SBarry Smith 
125e76c431SBarry Smith     Input parameters:
135e42470aSBarry Smith .   snes - nonlinear context obtained from SNESCreate()
145e76c431SBarry Smith 
155e42470aSBarry Smith     Output Parameters:
165e42470aSBarry Smith .   its  - Number of global iterations until termination.
175e76c431SBarry Smith 
185e76c431SBarry Smith     Notes:
195e42470aSBarry Smith     See SNESCreate() and SNESSetUp() for information on the definition and
205e76c431SBarry Smith     initialization of the nonlinear solver context.
215e76c431SBarry Smith 
225e76c431SBarry Smith     This implements essentially a truncated Newton method with a
235e76c431SBarry Smith     line search.  By default a cubic backtracking line search
245e76c431SBarry Smith     is employed, as described in the text "Numerical Methods for
255e76c431SBarry Smith     Unconstrained Optimization and Nonlinear Equations" by Dennis
265e42470aSBarry Smith     and Schnabel.  See the examples in src/snes/examples.
275e76c431SBarry Smith */
285e76c431SBarry Smith 
295e42470aSBarry Smith int SNESSolve_LS( SNES snes, int *outits )
305e76c431SBarry Smith {
315e42470aSBarry Smith   SNES_LS      *neP = (SNES_LS *) snes->data;
32edd2f0e1SBarry Smith   int          maxits, i, history_len,ierr,lits;
33df60cc22SBarry Smith   MatStructure flg = ALLMAT_DIFFERENT_NONZERO_PATTERN;
346b5873e3SBarry Smith   double       fnorm, gnorm, xnorm, ynorm, *history;
355e42470aSBarry Smith   Vec          Y, X, F, G, W, TMP;
365e76c431SBarry Smith 
375e42470aSBarry Smith   history	= snes->conv_hist;	/* convergence history */
385e42470aSBarry Smith   history_len	= snes->conv_hist_len;	/* convergence history length */
395e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
405e42470aSBarry Smith   X		= snes->vec_sol;		/* solution vector */
4139e2f89bSBarry Smith   F		= snes->vec_func;		/* residual vector */
425e42470aSBarry Smith   Y		= snes->work[0];		/* work vectors */
435e42470aSBarry Smith   G		= snes->work[1];
445e42470aSBarry Smith   W		= snes->work[2];
455e76c431SBarry Smith 
465e42470aSBarry Smith   ierr = SNESComputeInitialGuess(snes,X); CHKERR(ierr);  /* X <- X_0 */
475e42470aSBarry Smith   VecNorm( X, &xnorm );		       /* xnorm = || X || */
48a9581db2SBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERR(ierr); /* (+/-) F(X) */
495e42470aSBarry Smith   VecNorm(F, &fnorm );	        	/* fnorm <- || F || */
505e42470aSBarry Smith   snes->norm = fnorm;
515e76c431SBarry Smith   if (history && history_len > 0) history[0] = fnorm;
52eafb4bcbSBarry Smith   if (snes->Monitor) (*snes->Monitor)(snes,0,fnorm,snes->monP);
535e76c431SBarry Smith 
545e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
555e42470aSBarry Smith        snes->iter = i+1;
565e76c431SBarry Smith 
575e76c431SBarry Smith        /* Solve J Y = -F, where J is Jacobian matrix */
58df60cc22SBarry Smith        ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,
59df60cc22SBarry Smith                                 &flg,snes->jacP); CHKERR(ierr);
6023242f5aSBarry Smith        ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);
615e42470aSBarry Smith        ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERR(ierr);
625e42470aSBarry Smith        ierr = (*neP->LineSearch)(snes, X, F, G, Y, W, fnorm, &ynorm, &gnorm );
6323242f5aSBarry Smith        CHKERR(ierr);
645e76c431SBarry Smith 
6539e2f89bSBarry Smith        TMP = F; F = G; snes->vec_func_always = F; G = TMP;
6639e2f89bSBarry Smith        TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
675e76c431SBarry Smith        fnorm = gnorm;
685e76c431SBarry Smith 
695e42470aSBarry Smith        snes->norm = fnorm;
705e76c431SBarry Smith        if (history && history_len > i+1) history[i+1] = fnorm;
715e42470aSBarry Smith        VecNorm( X, &xnorm );		/* xnorm = || X || */
72eafb4bcbSBarry Smith        if (snes->Monitor) (*snes->Monitor)(snes,i+1,fnorm,snes->monP);
735e76c431SBarry Smith 
745e76c431SBarry Smith        /* Test for convergence */
755e42470aSBarry Smith        if ((*snes->Converged)( snes, xnorm, ynorm, fnorm,snes->cnvP )) {
7639e2f89bSBarry Smith            if (X != snes->vec_sol) {
7739e2f89bSBarry Smith              VecCopy( X, snes->vec_sol );
7839e2f89bSBarry Smith              snes->vec_sol_always = snes->vec_sol;
7939e2f89bSBarry Smith              snes->vec_func_always = snes->vec_func;
8039e2f89bSBarry Smith            }
815e76c431SBarry Smith            break;
825e76c431SBarry Smith        }
835e76c431SBarry Smith   }
845e76c431SBarry Smith   if (i == maxits) i--;
855e42470aSBarry Smith   *outits = i+1;
865e42470aSBarry Smith   return 0;
875e76c431SBarry Smith }
885e76c431SBarry Smith /* ------------------------------------------------------------ */
895e42470aSBarry Smith int SNESSetUp_LS(SNES snes )
905e76c431SBarry Smith {
915e42470aSBarry Smith   int ierr;
925e42470aSBarry Smith   snes->nwork = 3;
935e42470aSBarry Smith   ierr = VecGetVecs( snes->vec_sol, snes->nwork,&snes->work ); CHKERR(ierr);
945e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work );
955e42470aSBarry Smith   return 0;
965e76c431SBarry Smith }
975e76c431SBarry Smith /* ------------------------------------------------------------ */
985e42470aSBarry Smith int SNESDestroy_LS(PetscObject obj)
995e76c431SBarry Smith {
1005e42470aSBarry Smith   SNES snes = (SNES) obj;
1015e42470aSBarry Smith   VecFreeVecs(snes->work, snes->nwork );
102df60cc22SBarry Smith   FREE(snes->data);
1035e42470aSBarry Smith   return 0;
1045e76c431SBarry Smith }
1055e42470aSBarry Smith /*@
1060de89847SLois Curfman McInnes    SNESDefaultMonitor - Default SNES monitoring routine.  Prints the
1075e42470aSBarry Smith    residual norm at each iteration.
1085e42470aSBarry Smith 
1095e42470aSBarry Smith    Input Parameters:
1100de89847SLois Curfman McInnes .  snes - the SNES context
1110de89847SLois Curfman McInnes .  its - iteration number
1120de89847SLois Curfman McInnes .  fnorm - 2-norm residual value (may be estimated)
1130de89847SLois Curfman McInnes .  dummy - unused context
1145e42470aSBarry Smith 
1155e42470aSBarry Smith    Notes:
1165e42470aSBarry Smith    f is either the residual or its negative, depending on the user's
117a9581db2SBarry Smith    preference, as set with SNESSetFunction().
118a67c8522SLois Curfman McInnes 
119a67c8522SLois Curfman McInnes .keywords: SNES, nonlinear, default, monitor, residual, norm
120a67c8522SLois Curfman McInnes 
121a67c8522SLois Curfman McInnes .seealso: SNESSetMonitor()
1225e42470aSBarry Smith @*/
123eafb4bcbSBarry Smith int SNESDefaultMonitor(SNES snes,int its, double fnorm,void *dummy)
1245e42470aSBarry Smith {
125*69e4de80SBarry Smith   if (fnorm > 1.e-9 || fnorm < -1.e-9 || fnorm == 0.0) {
1267f813858SBarry Smith     MPIU_printf(snes->comm, "iter = %d, Function norm %g \n",its,fnorm);
127*69e4de80SBarry Smith   }
128*69e4de80SBarry Smith   else {
129*69e4de80SBarry Smith     MPIU_printf(snes->comm, "iter = %d, Function norm %5.3e \n",its,fnorm);
130*69e4de80SBarry Smith   }
1315e42470aSBarry Smith   return 0;
1325e42470aSBarry Smith }
1335e42470aSBarry Smith 
1345e42470aSBarry Smith /*@
1355e42470aSBarry Smith    SNESDefaultConverged - Default test for monitoring the convergence
1360de89847SLois Curfman McInnes    of the SNES solvers.
1375e42470aSBarry Smith 
1385e42470aSBarry Smith    Input Parameters:
1390de89847SLois Curfman McInnes .  snes - the SNES context
1405e42470aSBarry Smith .  xnorm - 2-norm of current iterate
1415e42470aSBarry Smith .  pnorm - 2-norm of current step
1425e42470aSBarry Smith .  fnorm - 2-norm of residual
1430de89847SLois Curfman McInnes .  dummy - unused context
1445e42470aSBarry Smith 
1455e42470aSBarry Smith    Returns:
1465e42470aSBarry Smith $  2  if  ( fnorm < atol ),
1475e42470aSBarry Smith $  3  if  ( pnorm < xtol*xnorm ),
1485e42470aSBarry Smith $ -2  if  ( nres > max_res ),
1495e42470aSBarry Smith $  0  otherwise,
1505e42470aSBarry Smith 
1515e42470aSBarry Smith    where
1525e42470aSBarry Smith $    atol    - absolute residual norm tolerance,
1530de89847SLois Curfman McInnes $              set with SNESSetAbsoluteTolerance()
1545e42470aSBarry Smith $    max_res - maximum number of residual evaluations,
1550de89847SLois Curfman McInnes $              set with SNESSetMaxResidualEvaluations()
1565e42470aSBarry Smith $    nres    - number of residual evaluations
1575e42470aSBarry Smith $    xtol    - relative residual norm tolerance,
1580de89847SLois Curfman McInnes $              set with SNESSetRelativeTolerance()
1595e42470aSBarry Smith 
1600de89847SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
1615e42470aSBarry Smith 
1620de89847SLois Curfman McInnes .seealso: SNESSetConvergenceTest()
1635e42470aSBarry Smith @*/
1645e42470aSBarry Smith int SNESDefaultConverged(SNES snes,double xnorm,double pnorm,double fnorm,
1655e42470aSBarry Smith                          void *dummy)
1665e42470aSBarry Smith {
1675e42470aSBarry Smith   if (fnorm < snes->atol) {
168a67c8522SLois Curfman McInnes     PLogInfo((PetscObject)snes,
169a67c8522SLois Curfman McInnes       "Converged due to absolute residual norm %g < %g\n",fnorm,snes->atol);
1705e42470aSBarry Smith     return 2;
1715e42470aSBarry Smith   }
1725e42470aSBarry Smith   if (pnorm < snes->xtol*(xnorm)) {
173a67c8522SLois Curfman McInnes     PLogInfo((PetscObject)snes,
174a67c8522SLois Curfman McInnes       "Converged due to small update length  %g < %g*%g\n",
1755e42470aSBarry Smith        pnorm,snes->xtol,xnorm);
1765e42470aSBarry Smith     return 3;
1775e42470aSBarry Smith   }
17849e3953dSBarry Smith   if (snes->nfuncs > snes->max_funcs) {
179a67c8522SLois Curfman McInnes     PLogInfo((PetscObject)snes,
180a67c8522SLois Curfman McInnes       "Exceeded maximum number of residual evaluations: %d > %d\n",
18149e3953dSBarry Smith        snes->nfuncs, snes->max_funcs );
182a67c8522SLois Curfman McInnes     return -2;
183a67c8522SLois Curfman McInnes   }
1845e42470aSBarry Smith   return 0;
1855e42470aSBarry Smith }
1865e42470aSBarry Smith 
1875e76c431SBarry Smith /* ------------------------------------------------------------ */
1885e76c431SBarry Smith /*ARGSUSED*/
1895e76c431SBarry Smith /*@
1905e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
1915e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
1925e76c431SBarry Smith    to serve as a template and is not recommended for general use.
1935e76c431SBarry Smith 
1945e76c431SBarry Smith    Input Parameters:
1955e42470aSBarry Smith .  snes - nonlinear context
1965e76c431SBarry Smith .  x - current iterate
1975e76c431SBarry Smith .  f - residual evaluated at x
1985e76c431SBarry Smith .  y - search direction (contains new iterate on output)
1995e76c431SBarry Smith .  w - work vector
2005e76c431SBarry Smith .  fnorm - 2-norm of f
2015e76c431SBarry Smith 
202c4a48953SLois Curfman McInnes    Output Parameters:
2035e76c431SBarry Smith .  g - residual evaluated at new iterate y
2045e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2055e76c431SBarry Smith .  gnorm - 2-norm of g
2065e76c431SBarry Smith .  ynorm - 2-norm of search length
2075e76c431SBarry Smith 
208c4a48953SLois Curfman McInnes    Options Database Key:
209c4a48953SLois Curfman McInnes $  -snes_line_search basic
210c4a48953SLois Curfman McInnes 
2115e76c431SBarry Smith    Returns:
2125e76c431SBarry Smith    1, indicating success of the step.
2135e76c431SBarry Smith 
21428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
21528ae5a21SLois Curfman McInnes 
216f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
217f59ffdedSLois Curfman McInnes .seealso: SNESSetLineSearchRoutine()
2185e76c431SBarry Smith @*/
2195e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
2205e42470aSBarry Smith                              double fnorm, double *ynorm, double *gnorm )
2215e76c431SBarry Smith {
2225e42470aSBarry Smith   int    ierr;
2235e42470aSBarry Smith   Scalar one = 1.0;
2247857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
2255e42470aSBarry Smith   VecNorm(y, ynorm );	/* ynorm = || y ||    */
2265e42470aSBarry Smith   VecAXPY(&one, x, y );	/* y <- x + y         */
227a9581db2SBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERR(ierr);
2285e42470aSBarry Smith   VecNorm( g, gnorm ); 	/* gnorm = || g ||    */
2297857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2305e76c431SBarry Smith   return 1;
2315e76c431SBarry Smith }
2325e76c431SBarry Smith /* ------------------------------------------------------------------ */
2335e76c431SBarry Smith /*@
234f59ffdedSLois Curfman McInnes    SNESCubicLineSearch - This routine performs a cubic line search and
235f59ffdedSLois Curfman McInnes    is the default line search method.
2365e76c431SBarry Smith 
2375e76c431SBarry Smith    Input Parameters:
2385e42470aSBarry Smith .  snes - nonlinear context
2395e76c431SBarry Smith .  x - current iterate
2405e76c431SBarry Smith .  f - residual evaluated at x
2415e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2425e76c431SBarry Smith .  w - work vector
2435e76c431SBarry Smith .  fnorm - 2-norm of f
2445e76c431SBarry Smith 
2455e76c431SBarry Smith    Output parameters:
2465e76c431SBarry Smith .  g - residual evaluated at new iterate y
2475e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2485e76c431SBarry Smith .  gnorm - 2-norm of g
2495e76c431SBarry Smith .  ynorm - 2-norm of search length
2505e76c431SBarry Smith 
2515e76c431SBarry Smith    Returns:
2525e76c431SBarry Smith    1 if the line search succeeds; 0 if the line search fails.
2535e76c431SBarry Smith 
254c4a48953SLois Curfman McInnes    Options Database Key:
255c4a48953SLois Curfman McInnes $  -snes_line_search cubic
256c4a48953SLois Curfman McInnes 
2575e76c431SBarry Smith    Notes:
2585e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
2595e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
2605e76c431SBarry Smith 
26128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
26228ae5a21SLois Curfman McInnes 
263f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine()
2645e76c431SBarry Smith @*/
2655e42470aSBarry Smith int SNESCubicLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
2665e42470aSBarry Smith                               double fnorm, double *ynorm, double *gnorm )
2675e76c431SBarry Smith {
2685e42470aSBarry Smith   double  steptol, initslope;
2695e42470aSBarry Smith   double  lambdaprev, gnormprev;
2705e76c431SBarry Smith   double  a, b, d, t1, t2;
2716b5873e3SBarry Smith #if defined(PETSC_COMPLEX)
2725e42470aSBarry Smith   Scalar  cinitslope,clambda;
2736b5873e3SBarry Smith #endif
2745e42470aSBarry Smith   int     ierr,count;
2755e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
2765e42470aSBarry Smith   Scalar  one = 1.0,scale;
2775e42470aSBarry Smith   double  maxstep,minlambda,alpha,lambda,lambdatemp;
2785e76c431SBarry Smith 
2797857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
2805e76c431SBarry Smith   alpha   = neP->alpha;
2815e76c431SBarry Smith   maxstep = neP->maxstep;
2825e76c431SBarry Smith   steptol = neP->steptol;
2835e76c431SBarry Smith 
2845e42470aSBarry Smith   VecNorm(y, ynorm );
2855e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
2865e42470aSBarry Smith     scale = maxstep/(*ynorm);
2876b5873e3SBarry Smith #if defined(PETSC_COMPLEX)
2886b5873e3SBarry Smith     PLogInfo((PetscObject)snes,"Scaling step by %g\n",real(scale));
2896b5873e3SBarry Smith #else
2905e42470aSBarry Smith     PLogInfo((PetscObject)snes,"Scaling step by %g\n",scale);
2916b5873e3SBarry Smith #endif
2925e42470aSBarry Smith     VecScale(&scale, y );
2935e76c431SBarry Smith     *ynorm = maxstep;
2945e76c431SBarry Smith   }
2955e76c431SBarry Smith   minlambda = steptol/(*ynorm);
2965e42470aSBarry Smith #if defined(PETSC_COMPLEX)
2975e42470aSBarry Smith   VecDot(f, y, &cinitslope );
2985e42470aSBarry Smith   initslope = real(cinitslope);
2995e42470aSBarry Smith #else
3005e42470aSBarry Smith   VecDot(f, y, &initslope );
3015e42470aSBarry Smith #endif
3025e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3035e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3045e76c431SBarry Smith 
3055e42470aSBarry Smith   VecCopy(y, w );
3065e42470aSBarry Smith   VecAXPY(&one, x, w );
307a9581db2SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
3085e42470aSBarry Smith   VecNorm(g, gnorm );
3095e76c431SBarry Smith   if (*gnorm <= fnorm + alpha*initslope) {	/* Sufficient reduction */
3105e42470aSBarry Smith       VecCopy(w, y );
3115e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Using full step\n");
3127857610eSBarry Smith       PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3135e42470aSBarry Smith       return 0;
3145e76c431SBarry Smith   }
3155e76c431SBarry Smith 
3165e76c431SBarry Smith   /* Fit points with quadratic */
3175e76c431SBarry Smith   lambda = 1.0; count = 0;
3185e76c431SBarry Smith   lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope));
3195e76c431SBarry Smith   lambdaprev = lambda;
3205e76c431SBarry Smith   gnormprev = *gnorm;
3215e76c431SBarry Smith   if (lambdatemp <= .1*lambda) {
3225e76c431SBarry Smith       lambda = .1*lambda;
3235e76c431SBarry Smith   } else lambda = lambdatemp;
3245e42470aSBarry Smith   VecCopy(x, w );
3255e42470aSBarry Smith #if defined(PETSC_COMPLEX)
3265e42470aSBarry Smith   clambda = lambda; VecAXPY(&clambda, y, w );
3275e42470aSBarry Smith #else
3285e42470aSBarry Smith   VecAXPY(&lambda, y, w );
3295e42470aSBarry Smith #endif
330a9581db2SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
3315e42470aSBarry Smith   VecNorm(g, gnorm );
3325e76c431SBarry Smith   if (*gnorm <= fnorm + alpha*initslope) {      /* sufficient reduction */
3335e42470aSBarry Smith       VecCopy(w, y );
3345e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda);
3357857610eSBarry Smith       PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3365e42470aSBarry Smith       return 0;
3375e76c431SBarry Smith   }
3385e76c431SBarry Smith 
3395e76c431SBarry Smith   /* Fit points with cubic */
3405e76c431SBarry Smith   count = 1;
3415e76c431SBarry Smith   while (1) {
3425e76c431SBarry Smith       if (lambda <= minlambda) { /* bad luck; use full step */
3435e42470aSBarry Smith            PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count);
3445e42470aSBarry Smith            PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n",
3455e76c431SBarry Smith                    fnorm,*gnorm, *ynorm,lambda);
3465e42470aSBarry Smith            VecCopy(w, y );
3477857610eSBarry Smith            PLogEventEnd(SNES_LineSearch,snes,x,f,g);
34823242f5aSBarry Smith            return -1;
3495e76c431SBarry Smith       }
3505e76c431SBarry Smith       t1 = *gnorm - fnorm - lambda*initslope;
3515e76c431SBarry Smith       t2 = gnormprev  - fnorm - lambdaprev*initslope;
3525e76c431SBarry Smith       a = (t1/(lambda*lambda) -
3535e76c431SBarry Smith                 t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
3545e76c431SBarry Smith       b = (-lambdaprev*t1/(lambda*lambda) +
3555e76c431SBarry Smith                 lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
3565e76c431SBarry Smith       d = b*b - 3*a*initslope;
3575e76c431SBarry Smith       if (d < 0.0) d = 0.0;
3585e76c431SBarry Smith       if (a == 0.0) {
3595e76c431SBarry Smith          lambdatemp = -initslope/(2.0*b);
3605e76c431SBarry Smith       } else {
3615e76c431SBarry Smith          lambdatemp = (-b + sqrt(d))/(3.0*a);
3625e76c431SBarry Smith       }
3635e76c431SBarry Smith       if (lambdatemp > .5*lambda) {
3645e76c431SBarry Smith          lambdatemp = .5*lambda;
3655e76c431SBarry Smith       }
3665e76c431SBarry Smith       lambdaprev = lambda;
3675e76c431SBarry Smith       gnormprev = *gnorm;
3685e76c431SBarry Smith       if (lambdatemp <= .1*lambda) {
3695e76c431SBarry Smith          lambda = .1*lambda;
3705e76c431SBarry Smith       }
3715e76c431SBarry Smith       else lambda = lambdatemp;
3725e42470aSBarry Smith       VecCopy( x, w );
3735e42470aSBarry Smith #if defined(PETSC_COMPLEX)
3745e42470aSBarry Smith       VecAXPY(&clambda, y, w );
3755e42470aSBarry Smith #else
3765e42470aSBarry Smith       VecAXPY(&lambda, y, w );
3775e42470aSBarry Smith #endif
378a9581db2SBarry Smith       ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
3795e42470aSBarry Smith       VecNorm(g, gnorm );
3805e76c431SBarry Smith       if (*gnorm <= fnorm + alpha*initslope) {      /* is reduction enough */
3815e42470aSBarry Smith          VecCopy(w, y );
3825e42470aSBarry Smith          PLogInfo((PetscObject)snes,"Cubically determined step, lambda %g\n",lambda);
3837857610eSBarry Smith          PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3845e42470aSBarry Smith          return 0;
3855e76c431SBarry Smith       }
3865e76c431SBarry Smith       count++;
3875e76c431SBarry Smith    }
3887857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3895e42470aSBarry Smith   return 0;
3905e76c431SBarry Smith }
3915e42470aSBarry Smith /*@
3925e42470aSBarry Smith    SNESQuadraticLineSearch - This routine performs a cubic line search.
3935e76c431SBarry Smith 
3945e42470aSBarry Smith    Input Parameters:
395f59ffdedSLois Curfman McInnes .  snes - the SNES context
3965e42470aSBarry Smith .  x - current iterate
3975e42470aSBarry Smith .  f - residual evaluated at x
3985e42470aSBarry Smith .  y - search direction (contains new iterate on output)
3995e42470aSBarry Smith .  w - work vector
4005e42470aSBarry Smith .  fnorm - 2-norm of f
4015e42470aSBarry Smith 
402c4a48953SLois Curfman McInnes    Output Parameters:
4035e42470aSBarry Smith .  g - residual evaluated at new iterate y
4045e42470aSBarry Smith .  y - new iterate (contains search direction on input)
4055e42470aSBarry Smith .  gnorm - 2-norm of g
4065e42470aSBarry Smith .  ynorm - 2-norm of search length
4075e42470aSBarry Smith 
4085e42470aSBarry Smith    Returns:
4095e42470aSBarry Smith    1 if the line search succeeds; 0 if the line search fails.
4105e42470aSBarry Smith 
411c4a48953SLois Curfman McInnes    Options Database Key:
412c4a48953SLois Curfman McInnes $  -snes_line_search quadratic
413c4a48953SLois Curfman McInnes 
4145e42470aSBarry Smith    Notes:
415f59ffdedSLois Curfman McInnes    Use SNESSetLineSearchRoutine()
416f59ffdedSLois Curfman McInnes    to set this routine within the SNES_NLS method.
4175e42470aSBarry Smith 
418f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
419f59ffdedSLois Curfman McInnes 
420f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine()
4215e42470aSBarry Smith @*/
4225e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
4235e42470aSBarry Smith                               double fnorm, double *ynorm, double *gnorm )
4245e76c431SBarry Smith {
4255e42470aSBarry Smith   double  steptol, initslope;
4265e42470aSBarry Smith   double  lambdaprev, gnormprev;
4276b5873e3SBarry Smith #if defined(PETSC_COMPLEX)
4285e42470aSBarry Smith   Scalar  cinitslope,clambda;
4296b5873e3SBarry Smith #endif
4305e42470aSBarry Smith   int     ierr,count;
4315e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
4325e42470aSBarry Smith   Scalar  one = 1.0,scale;
4335e42470aSBarry Smith   double  maxstep,minlambda,alpha,lambda,lambdatemp;
4345e76c431SBarry Smith 
4357857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
4365e42470aSBarry Smith   alpha   = neP->alpha;
4375e42470aSBarry Smith   maxstep = neP->maxstep;
4385e42470aSBarry Smith   steptol = neP->steptol;
4395e76c431SBarry Smith 
4405e42470aSBarry Smith   VecNorm(y, ynorm );
4415e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4425e42470aSBarry Smith     scale = maxstep/(*ynorm);
4435e42470aSBarry Smith     VecScale(&scale, y );
4445e42470aSBarry Smith     *ynorm = maxstep;
4455e76c431SBarry Smith   }
4465e42470aSBarry Smith   minlambda = steptol/(*ynorm);
4475e42470aSBarry Smith #if defined(PETSC_COMPLEX)
4485e42470aSBarry Smith   VecDot(f, y, &cinitslope );
4495e42470aSBarry Smith   initslope = real(cinitslope);
4505e42470aSBarry Smith #else
4515e42470aSBarry Smith   VecDot(f, y, &initslope );
4525e42470aSBarry Smith #endif
4535e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
4545e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
4555e42470aSBarry Smith 
4565e42470aSBarry Smith   VecCopy(y, w );
4575e42470aSBarry Smith   VecAXPY(&one, x, w );
458a9581db2SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
4595e42470aSBarry Smith   VecNorm(g, gnorm );
4605e42470aSBarry Smith   if (*gnorm <= fnorm + alpha*initslope) {	/* Sufficient reduction */
4615e42470aSBarry Smith       VecCopy(w, y );
4625e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Using full step\n");
4637857610eSBarry Smith       PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4645e42470aSBarry Smith       return 0;
4655e42470aSBarry Smith   }
4665e42470aSBarry Smith 
4675e42470aSBarry Smith   /* Fit points with quadratic */
4685e42470aSBarry Smith   lambda = 1.0; count = 0;
4695e42470aSBarry Smith   count = 1;
4705e42470aSBarry Smith   while (1) {
4715e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4725e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count);
4735e42470aSBarry Smith       PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n",
4745e42470aSBarry Smith                    fnorm,*gnorm, *ynorm,lambda);
4755e42470aSBarry Smith       VecCopy(w, y );
4767857610eSBarry Smith       PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4775e42470aSBarry Smith       return 0;
4785e42470aSBarry Smith     }
4795e42470aSBarry Smith     lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope));
4805e42470aSBarry Smith     lambdaprev = lambda;
4815e42470aSBarry Smith     gnormprev = *gnorm;
4825e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
4835e42470aSBarry Smith       lambda = .1*lambda;
4845e42470aSBarry Smith     } else lambda = lambdatemp;
4855e42470aSBarry Smith     VecCopy(x, w );
4865e42470aSBarry Smith #if defined(PETSC_COMPLEX)
4875e42470aSBarry Smith     clambda = lambda; VecAXPY(&clambda, y, w );
4885e42470aSBarry Smith #else
4895e42470aSBarry Smith     VecAXPY(&lambda, y, w );
4905e42470aSBarry Smith #endif
491a9581db2SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERR(ierr);
4925e42470aSBarry Smith     VecNorm(g, gnorm );
4935e42470aSBarry Smith     if (*gnorm <= fnorm + alpha*initslope) {      /* sufficient reduction */
4945e42470aSBarry Smith       VecCopy(w, y );
4955e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda);
4967857610eSBarry Smith       PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4975e42470aSBarry Smith       return 0;
4985e42470aSBarry Smith     }
4995e42470aSBarry Smith     count++;
5005e42470aSBarry Smith   }
5015e42470aSBarry Smith 
5027857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5035e42470aSBarry Smith   return 0;
5045e76c431SBarry Smith }
5055e76c431SBarry Smith /* ------------------------------------------------------------ */
5065e76c431SBarry Smith /*@C
5075e42470aSBarry Smith    SNESSetLineSearchRoutine - Sets the line search routine to be used
508fbe28522SBarry Smith    by the method SNES_LS.
5095e76c431SBarry Smith 
5105e76c431SBarry Smith    Input Parameters:
511eafb4bcbSBarry Smith .  snes - nonlinear context obtained from SNESCreate()
5125e76c431SBarry Smith .  func - pointer to int function
5135e76c431SBarry Smith 
514c4a48953SLois Curfman McInnes    Available Routines:
515f59ffdedSLois Curfman McInnes .  SNESCubicLineSearch() - default line search
516f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
517f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
5185e76c431SBarry Smith 
519c4a48953SLois Curfman McInnes     Options Database Keys:
520c4a48953SLois Curfman McInnes $   -snes_line_search [basic,quadratic,cubic]
521c4a48953SLois Curfman McInnes 
5225e76c431SBarry Smith    Calling sequence of func:
523f59ffdedSLois Curfman McInnes    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
5245e42470aSBarry Smith          Vec w, double fnorm, double *ynorm, double *gnorm)
5255e76c431SBarry Smith 
5265e76c431SBarry Smith     Input parameters for func:
5275e42470aSBarry Smith .   snes - nonlinear context
5285e76c431SBarry Smith .   x - current iterate
5295e76c431SBarry Smith .   f - residual evaluated at x
5305e76c431SBarry Smith .   y - search direction (contains new iterate on output)
5315e76c431SBarry Smith .   w - work vector
5325e76c431SBarry Smith .   fnorm - 2-norm of f
5335e76c431SBarry Smith 
5345e76c431SBarry Smith     Output parameters for func:
5355e76c431SBarry Smith .   g - residual evaluated at new iterate y
5365e76c431SBarry Smith .   y - new iterate (contains search direction on input)
5375e76c431SBarry Smith .   gnorm - 2-norm of g
5385e76c431SBarry Smith .   ynorm - 2-norm of search length
5395e76c431SBarry Smith 
5405e76c431SBarry Smith     Returned by func:
5415e76c431SBarry Smith     1 if the line search succeeds; 0 if the line search fails.
542f59ffdedSLois Curfman McInnes 
543f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
544f59ffdedSLois Curfman McInnes 
545f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
5465e76c431SBarry Smith @*/
5475e42470aSBarry Smith int SNESSetLineSearchRoutine(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
5485e42470aSBarry Smith                              double,double *,double*) )
5495e76c431SBarry Smith {
550fbe28522SBarry Smith   if ((snes)->type == SNES_NLS)
5515e42470aSBarry Smith     ((SNES_LS *)(snes->data))->LineSearch = func;
5525e42470aSBarry Smith   return 0;
5535e76c431SBarry Smith }
5545e42470aSBarry Smith 
5555e42470aSBarry Smith static int SNESPrintHelp_LS(SNES snes)
5565e42470aSBarry Smith {
5575e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
5585e42470aSBarry Smith   fprintf(stderr,"-snes_line_search [basic,quadratic,cubic]\n");
5595e42470aSBarry Smith   fprintf(stderr,"-snes_line_search_alpha alpha (default %g)\n",ls->alpha);
5605e42470aSBarry Smith   fprintf(stderr,"-snes_line_search_maxstep max (default %g)\n",ls->maxstep);
5615e42470aSBarry Smith   fprintf(stderr,"-snes_line_search_steptol tol (default %g)\n",ls->steptol);
5626b5873e3SBarry Smith   return 0;
5635e42470aSBarry Smith }
5645e42470aSBarry Smith 
5655e42470aSBarry Smith static int SNESSetFromOptions_LS(SNES snes)
5665e42470aSBarry Smith {
5675e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
5685e42470aSBarry Smith   char    ver[16];
5695e42470aSBarry Smith   double  tmp;
5705e42470aSBarry Smith 
571df60cc22SBarry Smith   if (OptionsGetDouble(snes->prefix,"-snes_line_search_alpa",&tmp)) {
5725e42470aSBarry Smith     ls->alpha = tmp;
5735e42470aSBarry Smith   }
574df60cc22SBarry Smith   if (OptionsGetDouble(snes->prefix,"-snes_line_search_maxstep",&tmp)) {
5755e42470aSBarry Smith     ls->maxstep = tmp;
5765e42470aSBarry Smith   }
577df60cc22SBarry Smith   if (OptionsGetDouble(snes->prefix,"-snes_line_search_steptol",&tmp)) {
5785e42470aSBarry Smith     ls->steptol = tmp;
5795e42470aSBarry Smith   }
580df60cc22SBarry Smith   if (OptionsGetString(snes->prefix,"-snes_line_search",ver,16)) {
5815e42470aSBarry Smith     if (!strcmp(ver,"basic")) {
5825e42470aSBarry Smith       SNESSetLineSearchRoutine(snes,SNESNoLineSearch);
5835e42470aSBarry Smith     }
5845e42470aSBarry Smith     else if (!strcmp(ver,"quadratic")) {
5855e42470aSBarry Smith       SNESSetLineSearchRoutine(snes,SNESQuadraticLineSearch);
5865e42470aSBarry Smith     }
5875e42470aSBarry Smith     else if (!strcmp(ver,"cubic")) {
5885e42470aSBarry Smith       SNESSetLineSearchRoutine(snes,SNESCubicLineSearch);
5895e42470aSBarry Smith     }
5905e42470aSBarry Smith     else {SETERR(1,"Unknown line search?");}
5915e42470aSBarry Smith   }
5925e42470aSBarry Smith   return 0;
5935e42470aSBarry Smith }
5945e42470aSBarry Smith 
5955e42470aSBarry Smith /* ------------------------------------------------------------ */
5965e42470aSBarry Smith int SNESCreate_LS(SNES  snes )
5975e42470aSBarry Smith {
5985e42470aSBarry Smith   SNES_LS *neP;
5995e42470aSBarry Smith 
600fbe28522SBarry Smith   snes->type		= SNES_NLS;
60149e3953dSBarry Smith   snes->setup		= SNESSetUp_LS;
60249e3953dSBarry Smith   snes->solve		= SNESSolve_LS;
6035e42470aSBarry Smith   snes->destroy		= SNESDestroy_LS;
6045e42470aSBarry Smith   snes->Converged	= SNESDefaultConverged;
60549e3953dSBarry Smith   snes->printhelp       = SNESPrintHelp_LS;
60649e3953dSBarry Smith   snes->setfromoptions  = SNESSetFromOptions_LS;
6075e42470aSBarry Smith 
6085e42470aSBarry Smith   neP			= NEW(SNES_LS);   CHKPTR(neP);
6095e42470aSBarry Smith   snes->data    	= (void *) neP;
6105e42470aSBarry Smith   neP->alpha		= 1.e-4;
6115e42470aSBarry Smith   neP->maxstep		= 1.e8;
6125e42470aSBarry Smith   neP->steptol		= 1.e-12;
6135e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
6145e42470aSBarry Smith   return 0;
6155e42470aSBarry Smith }
6165e42470aSBarry Smith 
6175e42470aSBarry Smith 
6185e42470aSBarry Smith 
6195e42470aSBarry Smith 
620