xref: /petsc/src/snes/impls/ls/ls.c (revision 1a3481a4f86b15936ff2dcb29ec4403ee784c865)
15e76c431SBarry Smith #ifndef lint
2*1a3481a4SLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.33 1995/07/20 15:34:54 curfman Exp curfman $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
55e76c431SBarry Smith #include <math.h>
6fbe28522SBarry Smith #include "ls.h"
7a935fc98SLois Curfman McInnes #include "pviewer.h"
8bbb6d6a8SBarry Smith #if defined(HAVE_STRING_H)
9bbb6d6a8SBarry Smith #include <string.h>
10bbb6d6a8SBarry Smith #endif
115e42470aSBarry Smith 
125e42470aSBarry Smith /*
135e42470aSBarry Smith      Implements a line search variant of Newton's Method
145e76c431SBarry Smith     for solving systems of nonlinear equations.
155e76c431SBarry Smith 
165e76c431SBarry Smith     Input parameters:
175e42470aSBarry Smith .   snes - nonlinear context obtained from SNESCreate()
185e76c431SBarry Smith 
195e42470aSBarry Smith     Output Parameters:
20a935fc98SLois Curfman McInnes .   outits  - Number of global iterations until termination.
215e76c431SBarry Smith 
225e76c431SBarry Smith     Notes:
235e76c431SBarry Smith     This implements essentially a truncated Newton method with a
245e76c431SBarry Smith     line search.  By default a cubic backtracking line search
255e76c431SBarry Smith     is employed, as described in the text "Numerical Methods for
265e76c431SBarry Smith     Unconstrained Optimization and Nonlinear Equations" by Dennis
275e42470aSBarry Smith     and Schnabel.  See the examples in src/snes/examples.
285e76c431SBarry Smith */
295e76c431SBarry Smith 
305e42470aSBarry Smith int SNESSolve_LS(SNES snes,int *outits)
315e76c431SBarry Smith {
325e42470aSBarry Smith   SNES_LS      *neP = (SNES_LS *) snes->data;
33761aaf1bSLois Curfman McInnes   int          maxits, i, history_len, ierr, lits, lsfail;
34df60cc22SBarry Smith   MatStructure flg = ALLMAT_DIFFERENT_NONZERO_PATTERN;
356b5873e3SBarry Smith   double       fnorm, gnorm, xnorm, ynorm, *history;
365e42470aSBarry Smith   Vec          Y, X, F, G, W, TMP;
375e76c431SBarry Smith 
385e42470aSBarry Smith   history	= snes->conv_hist;	/* convergence history */
395e42470aSBarry Smith   history_len	= snes->conv_hist_len;	/* 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 
4778b31e54SBarry Smith   ierr = SNESComputeInitialGuess(snes,X); CHKERRQ(ierr); /* X <- X_0 */
485e42470aSBarry Smith   VecNorm(X,&xnorm);		                         /* xnorm = || X || */
4978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);   /* (+/-) F(X) */
505e42470aSBarry Smith   VecNorm(F,&fnorm);	        	                 /* fnorm <- ||F|| */
515e42470aSBarry Smith   snes->norm = fnorm;
525e76c431SBarry Smith   if (history && history_len > 0) history[0] = fnorm;
53bbb6d6a8SBarry Smith   if (snes->monitor)
54bbb6d6a8SBarry Smith     {ierr = (*snes->monitor)(snes,0,fnorm,snes->monP); CHKERRQ(ierr);}
555e76c431SBarry Smith 
565e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
575e42470aSBarry Smith        snes->iter = i+1;
585e76c431SBarry Smith 
595e76c431SBarry Smith        /* Solve J Y = -F, where J is Jacobian matrix */
60df60cc22SBarry Smith        ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,
61bbb6d6a8SBarry Smith                                 &flg); CHKERRQ(ierr);
62a935fc98SLois Curfman McInnes        ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,
63a935fc98SLois Curfman McInnes                                 flg); CHKERRQ(ierr);
6478b31e54SBarry Smith        ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
6581b6cf68SLois Curfman McInnes        ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
66761aaf1bSLois Curfman McInnes        ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);
67761aaf1bSLois Curfman McInnes        if (lsfail) snes->nfailures++;
6878b31e54SBarry Smith        CHKERRQ(ierr);
695e76c431SBarry Smith 
7039e2f89bSBarry Smith        TMP = F; F = G; snes->vec_func_always = F; G = TMP;
7139e2f89bSBarry Smith        TMP = X; X = Y; snes->vec_sol_always = X; Y = TMP;
725e76c431SBarry Smith        fnorm = gnorm;
735e76c431SBarry Smith 
745e42470aSBarry Smith        snes->norm = fnorm;
755e76c431SBarry Smith        if (history && history_len > i+1) history[i+1] = fnorm;
765e42470aSBarry Smith        VecNorm(X,&xnorm);		/* xnorm = || X || */
77bbb6d6a8SBarry Smith        if (snes->monitor)
78bbb6d6a8SBarry Smith          {(*snes->monitor)(snes,i+1,fnorm,snes->monP); CHKERRQ(ierr);}
795e76c431SBarry Smith 
805e76c431SBarry Smith        /* Test for convergence */
81bbb6d6a8SBarry Smith        if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
8239e2f89bSBarry Smith            if (X != snes->vec_sol) {
8339e2f89bSBarry Smith              VecCopy(X,snes->vec_sol);
8439e2f89bSBarry Smith              snes->vec_sol_always = snes->vec_sol;
8539e2f89bSBarry Smith              snes->vec_func_always = snes->vec_func;
8639e2f89bSBarry Smith            }
875e76c431SBarry Smith            break;
885e76c431SBarry Smith        }
895e76c431SBarry Smith   }
905e76c431SBarry Smith   if (i == maxits) i--;
915e42470aSBarry Smith   *outits = i+1;
925e42470aSBarry Smith   return 0;
935e76c431SBarry Smith }
945e76c431SBarry Smith /* ------------------------------------------------------------ */
955e42470aSBarry Smith int SNESSetUp_LS(SNES snes )
965e76c431SBarry Smith {
975e42470aSBarry Smith   int ierr;
9881b6cf68SLois Curfman McInnes   snes->nwork = 4;
9978b31e54SBarry Smith   ierr = VecGetVecs(snes->vec_sol,snes->nwork,&snes->work); CHKERRQ(ierr);
1005e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work );
10181b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1025e42470aSBarry Smith   return 0;
1035e76c431SBarry Smith }
1045e76c431SBarry Smith /* ------------------------------------------------------------ */
1055e42470aSBarry Smith int SNESDestroy_LS(PetscObject obj)
1065e76c431SBarry Smith {
1075e42470aSBarry Smith   SNES snes = (SNES) obj;
1085e42470aSBarry Smith   VecFreeVecs(snes->work, snes->nwork );
10978b31e54SBarry Smith   PETSCFREE(snes->data);
1105e42470aSBarry Smith   return 0;
1115e76c431SBarry Smith }
1125e42470aSBarry Smith /*@
1130de89847SLois Curfman McInnes    SNESDefaultMonitor - Default SNES monitoring routine.  Prints the
1145e42470aSBarry Smith    residual norm at each iteration.
1155e42470aSBarry Smith 
1165e42470aSBarry Smith    Input Parameters:
1170de89847SLois Curfman McInnes .  snes - the SNES context
1180de89847SLois Curfman McInnes .  its - iteration number
1190de89847SLois Curfman McInnes .  fnorm - 2-norm residual value (may be estimated)
1200de89847SLois Curfman McInnes .  dummy - unused context
1215e42470aSBarry Smith 
1225e42470aSBarry Smith    Notes:
1235e42470aSBarry Smith    f is either the residual or its negative, depending on the user's
124a9581db2SBarry Smith    preference, as set with SNESSetFunction().
125a67c8522SLois Curfman McInnes 
126a67c8522SLois Curfman McInnes .keywords: SNES, nonlinear, default, monitor, residual, norm
127a67c8522SLois Curfman McInnes 
128a67c8522SLois Curfman McInnes .seealso: SNESSetMonitor()
1295e42470aSBarry Smith @*/
130eafb4bcbSBarry Smith int SNESDefaultMonitor(SNES snes,int its,double fnorm,void *dummy)
1315e42470aSBarry Smith {
132614e12f8SBarry Smith   MPIU_printf(snes->comm, "iter = %d, Function norm %g \n",its,fnorm);
133614e12f8SBarry Smith   return 0;
134614e12f8SBarry Smith }
135614e12f8SBarry Smith 
136614e12f8SBarry Smith int SNESDefaultSMonitor(SNES snes,int its, double fnorm,void *dummy)
137614e12f8SBarry Smith {
138614e12f8SBarry Smith   if (fnorm > 1.e-9 || fnorm == 0.0) {
1397f813858SBarry Smith     MPIU_printf(snes->comm, "iter = %d, Function norm %g \n",its,fnorm);
14069e4de80SBarry Smith   }
141614e12f8SBarry Smith   else if (fnorm > 1.e-11){
14269e4de80SBarry Smith     MPIU_printf(snes->comm, "iter = %d, Function norm %5.3e \n",its,fnorm);
14369e4de80SBarry Smith   }
144614e12f8SBarry Smith   else {
145614e12f8SBarry Smith     MPIU_printf(snes->comm, "iter = %d, Function norm < 1.e-11\n",its);
146614e12f8SBarry Smith   }
1475e42470aSBarry Smith   return 0;
1485e42470aSBarry Smith }
1495e42470aSBarry Smith 
1505e42470aSBarry Smith /*@
1515e42470aSBarry Smith    SNESDefaultConverged - Default test for monitoring the convergence
1520de89847SLois Curfman McInnes    of the SNES solvers.
1535e42470aSBarry Smith 
1545e42470aSBarry Smith    Input Parameters:
1550de89847SLois Curfman McInnes .  snes - the SNES context
1565e42470aSBarry Smith .  xnorm - 2-norm of current iterate
1575e42470aSBarry Smith .  pnorm - 2-norm of current step
1585e42470aSBarry Smith .  fnorm - 2-norm of residual
1590de89847SLois Curfman McInnes .  dummy - unused context
1605e42470aSBarry Smith 
1615e42470aSBarry Smith    Returns:
1625e42470aSBarry Smith $  2  if  ( fnorm < atol ),
1635e42470aSBarry Smith $  3  if  ( pnorm < xtol*xnorm ),
1645e42470aSBarry Smith $ -2  if  ( nres > max_res ),
1655e42470aSBarry Smith $  0  otherwise,
1665e42470aSBarry Smith 
1675e42470aSBarry Smith    where
1685e42470aSBarry Smith $    atol    - absolute residual norm tolerance,
1690de89847SLois Curfman McInnes $              set with SNESSetAbsoluteTolerance()
1705e42470aSBarry Smith $    max_res - maximum number of residual evaluations,
1710de89847SLois Curfman McInnes $              set with SNESSetMaxResidualEvaluations()
1725e42470aSBarry Smith $    nres    - number of residual evaluations
1735e42470aSBarry Smith $    xtol    - relative residual norm tolerance,
1740de89847SLois Curfman McInnes $              set with SNESSetRelativeTolerance()
1755e42470aSBarry Smith 
1760de89847SLois Curfman McInnes .keywords: SNES, nonlinear, default, converged, convergence
1775e42470aSBarry Smith 
1780de89847SLois Curfman McInnes .seealso: SNESSetConvergenceTest()
1795e42470aSBarry Smith @*/
1805e42470aSBarry Smith int SNESDefaultConverged(SNES snes,double xnorm,double pnorm,double fnorm,
1815e42470aSBarry Smith                          void *dummy)
1825e42470aSBarry Smith {
1835e42470aSBarry Smith   if (fnorm < snes->atol) {
184a67c8522SLois Curfman McInnes     PLogInfo((PetscObject)snes,
185a67c8522SLois Curfman McInnes       "Converged due to absolute residual norm %g < %g\n",fnorm,snes->atol);
1865e42470aSBarry Smith     return 2;
1875e42470aSBarry Smith   }
1885e42470aSBarry Smith   if (pnorm < snes->xtol*(xnorm)) {
189a67c8522SLois Curfman McInnes     PLogInfo((PetscObject)snes,
190a67c8522SLois Curfman McInnes       "Converged due to small update length  %g < %g*%g\n",
1915e42470aSBarry Smith        pnorm,snes->xtol,xnorm);
1925e42470aSBarry Smith     return 3;
1935e42470aSBarry Smith   }
19449e3953dSBarry Smith   if (snes->nfuncs > snes->max_funcs) {
195a67c8522SLois Curfman McInnes     PLogInfo((PetscObject)snes,
196ddbbdb52SLois Curfman McInnes       "Exceeded maximum number of function evaluations: %d > %d\n",
19749e3953dSBarry Smith        snes->nfuncs, snes->max_funcs );
198a67c8522SLois Curfman McInnes     return -2;
199a67c8522SLois Curfman McInnes   }
2005e42470aSBarry Smith   return 0;
2015e42470aSBarry Smith }
2025e42470aSBarry Smith 
2035e76c431SBarry Smith /* ------------------------------------------------------------ */
2045e76c431SBarry Smith /*ARGSUSED*/
2055e76c431SBarry Smith /*@
2065e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2075e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2085e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2095e76c431SBarry Smith 
2105e76c431SBarry Smith    Input Parameters:
2115e42470aSBarry Smith .  snes - nonlinear context
2125e76c431SBarry Smith .  x - current iterate
2135e76c431SBarry Smith .  f - residual evaluated at x
2145e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2155e76c431SBarry Smith .  w - work vector
2165e76c431SBarry Smith .  fnorm - 2-norm of f
2175e76c431SBarry Smith 
218c4a48953SLois Curfman McInnes    Output Parameters:
2195e76c431SBarry Smith .  g - residual evaluated at new iterate y
2205e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2215e76c431SBarry Smith .  gnorm - 2-norm of g
2225e76c431SBarry Smith .  ynorm - 2-norm of search length
223761aaf1bSLois Curfman McInnes .  flag - set to 0, indicating a successful line search
2245e76c431SBarry Smith 
225c4a48953SLois Curfman McInnes    Options Database Key:
226c4a48953SLois Curfman McInnes $  -snes_line_search basic
227c4a48953SLois Curfman McInnes 
22828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
22928ae5a21SLois Curfman McInnes 
230f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
231a935fc98SLois Curfman McInnes           SNESSetLineSearchRoutine()
2325e76c431SBarry Smith @*/
2335e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
234761aaf1bSLois Curfman McInnes               double fnorm, double *ynorm, double *gnorm,int *flag )
2355e76c431SBarry Smith {
2365e42470aSBarry Smith   int    ierr;
2375e42470aSBarry Smith   Scalar one = 1.0;
238761aaf1bSLois Curfman McInnes   *flag = 0;
2397857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
2405e42470aSBarry Smith   VecNorm(y,ynorm);                            /* ynorm = || y || */
2415e42470aSBarry Smith   VecAXPY(&one,x,y);	                       /* y <- x + y      */
24278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr);
2435e42470aSBarry Smith   VecNorm(g,gnorm); 	                       /* gnorm = || g || */
2447857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
245761aaf1bSLois Curfman McInnes   return 0;
2465e76c431SBarry Smith }
2475e76c431SBarry Smith /* ------------------------------------------------------------------ */
2485e76c431SBarry Smith /*@
249f59ffdedSLois Curfman McInnes    SNESCubicLineSearch - This routine performs a cubic line search and
250f59ffdedSLois Curfman McInnes    is the default line search method.
2515e76c431SBarry Smith 
2525e76c431SBarry Smith    Input Parameters:
2535e42470aSBarry Smith .  snes - nonlinear context
2545e76c431SBarry Smith .  x - current iterate
2555e76c431SBarry Smith .  f - residual evaluated at x
2565e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2575e76c431SBarry Smith .  w - work vector
2585e76c431SBarry Smith .  fnorm - 2-norm of f
2595e76c431SBarry Smith 
2605e76c431SBarry Smith    Output parameters:
2615e76c431SBarry Smith .  g - residual evaluated at new iterate y
2625e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2635e76c431SBarry Smith .  gnorm - 2-norm of g
2645e76c431SBarry Smith .  ynorm - 2-norm of search length
265761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
2665e76c431SBarry Smith 
267c4a48953SLois Curfman McInnes    Options Database Key:
268c4a48953SLois Curfman McInnes $  -snes_line_search cubic
269c4a48953SLois Curfman McInnes 
2705e76c431SBarry Smith    Notes:
2715e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
2725e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
2735e76c431SBarry Smith 
27428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
27528ae5a21SLois Curfman McInnes 
276f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine()
2775e76c431SBarry Smith @*/
2785e42470aSBarry Smith int SNESCubicLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
279761aaf1bSLois Curfman McInnes                 double fnorm, double *ynorm, double *gnorm,int *flag)
2805e76c431SBarry Smith {
2815e42470aSBarry Smith   double  steptol, initslope;
2825e42470aSBarry Smith   double  lambdaprev, gnormprev;
2835e76c431SBarry Smith   double  a, b, d, t1, t2;
2846b5873e3SBarry Smith #if defined(PETSC_COMPLEX)
2855e42470aSBarry Smith   Scalar  cinitslope,clambda;
2866b5873e3SBarry Smith #endif
2875e42470aSBarry Smith   int     ierr,count;
2885e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
2895e42470aSBarry Smith   Scalar  one = 1.0,scale;
2905e42470aSBarry Smith   double  maxstep,minlambda,alpha,lambda,lambdatemp;
2915e76c431SBarry Smith 
2927857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
293761aaf1bSLois Curfman McInnes   *flag = 0;
2945e76c431SBarry Smith   alpha   = neP->alpha;
2955e76c431SBarry Smith   maxstep = neP->maxstep;
2965e76c431SBarry Smith   steptol = neP->steptol;
2975e76c431SBarry Smith 
2985e42470aSBarry Smith   VecNorm(y, ynorm );
2995e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3005e42470aSBarry Smith     scale = maxstep/(*ynorm);
3016b5873e3SBarry Smith #if defined(PETSC_COMPLEX)
3026b5873e3SBarry Smith     PLogInfo((PetscObject)snes,"Scaling step by %g\n",real(scale));
3036b5873e3SBarry Smith #else
3045e42470aSBarry Smith     PLogInfo((PetscObject)snes,"Scaling step by %g\n",scale);
3056b5873e3SBarry Smith #endif
3065e42470aSBarry Smith     VecScale(&scale, y );
3075e76c431SBarry Smith     *ynorm = maxstep;
3085e76c431SBarry Smith   }
3095e76c431SBarry Smith   minlambda = steptol/(*ynorm);
3105e42470aSBarry Smith #if defined(PETSC_COMPLEX)
3115e42470aSBarry Smith   VecDot(f, y, &cinitslope );
3125e42470aSBarry Smith   initslope = real(cinitslope);
3135e42470aSBarry Smith #else
3145e42470aSBarry Smith   VecDot(f, y, &initslope );
3155e42470aSBarry Smith #endif
3165e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3175e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3185e76c431SBarry Smith 
3195e42470aSBarry Smith   VecCopy(y, w );
3205e42470aSBarry Smith   VecAXPY(&one, x, w );
32178b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
3225e42470aSBarry Smith   VecNorm(g, gnorm );
3235e76c431SBarry Smith   if (*gnorm <= fnorm + alpha*initslope) {	/* Sufficient reduction */
3245e42470aSBarry Smith       VecCopy(w, y );
3255e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Using full step\n");
326d93a2b8dSBarry Smith       goto theend;
3275e76c431SBarry Smith   }
3285e76c431SBarry Smith 
3295e76c431SBarry Smith   /* Fit points with quadratic */
3305e76c431SBarry Smith   lambda = 1.0; count = 0;
3315e76c431SBarry Smith   lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope));
3325e76c431SBarry Smith   lambdaprev = lambda;
3335e76c431SBarry Smith   gnormprev = *gnorm;
3345e76c431SBarry Smith   if (lambdatemp <= .1*lambda) {
3355e76c431SBarry Smith       lambda = .1*lambda;
3365e76c431SBarry Smith   } else lambda = lambdatemp;
3375e42470aSBarry Smith   VecCopy(x, w );
3385e42470aSBarry Smith #if defined(PETSC_COMPLEX)
3395e42470aSBarry Smith   clambda = lambda; VecAXPY(&clambda, y, w );
3405e42470aSBarry Smith #else
3415e42470aSBarry Smith   VecAXPY(&lambda, y, w );
3425e42470aSBarry Smith #endif
34378b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
3445e42470aSBarry Smith   VecNorm(g, gnorm );
3455e76c431SBarry Smith   if (*gnorm <= fnorm + alpha*initslope) {      /* sufficient reduction */
3465e42470aSBarry Smith       VecCopy(w, y );
3475e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda);
348d93a2b8dSBarry Smith       goto theend;
3495e76c431SBarry Smith   }
3505e76c431SBarry Smith 
3515e76c431SBarry Smith   /* Fit points with cubic */
3525e76c431SBarry Smith   count = 1;
3535e76c431SBarry Smith   while (1) {
3545e76c431SBarry Smith       if (lambda <= minlambda) { /* bad luck; use full step */
3555e42470aSBarry Smith            PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count);
3565e42470aSBarry Smith            PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n",
3575e76c431SBarry Smith                    fnorm,*gnorm, *ynorm,lambda);
3585e42470aSBarry Smith            VecCopy(w, y );
359761aaf1bSLois Curfman McInnes            *flag = -1; break;
3605e76c431SBarry Smith       }
3615e76c431SBarry Smith       t1 = *gnorm - fnorm - lambda*initslope;
3625e76c431SBarry Smith       t2 = gnormprev  - fnorm - lambdaprev*initslope;
3635e76c431SBarry Smith       a = (t1/(lambda*lambda) -
3645e76c431SBarry Smith                 t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
3655e76c431SBarry Smith       b = (-lambdaprev*t1/(lambda*lambda) +
3665e76c431SBarry Smith                 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;
3835e42470aSBarry Smith       VecCopy( x, w );
3845e42470aSBarry Smith #if defined(PETSC_COMPLEX)
3855e42470aSBarry Smith       VecAXPY(&clambda, y, w );
3865e42470aSBarry Smith #else
3875e42470aSBarry Smith       VecAXPY(&lambda, y, w );
3885e42470aSBarry Smith #endif
38978b31e54SBarry Smith       ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
3905e42470aSBarry Smith       VecNorm(g, gnorm );
3915e76c431SBarry Smith       if (*gnorm <= fnorm + alpha*initslope) {      /* is reduction enough */
3925e42470aSBarry Smith          VecCopy(w, y );
3935e42470aSBarry Smith          PLogInfo((PetscObject)snes,"Cubically determined step, lambda %g\n",lambda);
394761aaf1bSLois Curfman McInnes          *flag = -1; break;
3955e76c431SBarry Smith       }
3965e76c431SBarry Smith       count++;
3975e76c431SBarry Smith    }
398d93a2b8dSBarry Smith   theend:
3997857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4005e42470aSBarry Smith   return 0;
4015e76c431SBarry Smith }
4025e42470aSBarry Smith /*@
4035e42470aSBarry Smith    SNESQuadraticLineSearch - This routine performs a cubic line search.
4045e76c431SBarry Smith 
4055e42470aSBarry Smith    Input Parameters:
406f59ffdedSLois Curfman McInnes .  snes - the SNES context
4075e42470aSBarry Smith .  x - current iterate
4085e42470aSBarry Smith .  f - residual evaluated at x
4095e42470aSBarry Smith .  y - search direction (contains new iterate on output)
4105e42470aSBarry Smith .  w - work vector
4115e42470aSBarry Smith .  fnorm - 2-norm of f
4125e42470aSBarry Smith 
413c4a48953SLois Curfman McInnes    Output Parameters:
4145e42470aSBarry Smith .  g - residual evaluated at new iterate y
4155e42470aSBarry Smith .  y - new iterate (contains search direction on input)
4165e42470aSBarry Smith .  gnorm - 2-norm of g
4175e42470aSBarry Smith .  ynorm - 2-norm of search length
418761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
4195e42470aSBarry Smith 
420c4a48953SLois Curfman McInnes    Options Database Key:
421c4a48953SLois Curfman McInnes $  -snes_line_search quadratic
422c4a48953SLois Curfman McInnes 
4235e42470aSBarry Smith    Notes:
424f59ffdedSLois Curfman McInnes    Use SNESSetLineSearchRoutine()
425f59ffdedSLois Curfman McInnes    to set this routine within the SNES_NLS method.
4265e42470aSBarry Smith 
427f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
428f59ffdedSLois Curfman McInnes 
429f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearchRoutine()
4305e42470aSBarry Smith @*/
4315e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
432761aaf1bSLois Curfman McInnes                     double fnorm, double *ynorm, double *gnorm,int *flag)
4335e76c431SBarry Smith {
4345e42470aSBarry Smith   double  steptol, initslope;
4355e42470aSBarry Smith   double  lambdaprev, gnormprev;
4366b5873e3SBarry Smith #if defined(PETSC_COMPLEX)
4375e42470aSBarry Smith   Scalar  cinitslope,clambda;
4386b5873e3SBarry Smith #endif
4395e42470aSBarry Smith   int     ierr,count;
4405e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
4415e42470aSBarry Smith   Scalar  one = 1.0,scale;
4425e42470aSBarry Smith   double  maxstep,minlambda,alpha,lambda,lambdatemp;
4435e76c431SBarry Smith 
4447857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
445761aaf1bSLois Curfman McInnes   *flag = 0;
4465e42470aSBarry Smith   alpha   = neP->alpha;
4475e42470aSBarry Smith   maxstep = neP->maxstep;
4485e42470aSBarry Smith   steptol = neP->steptol;
4495e76c431SBarry Smith 
4505e42470aSBarry Smith   VecNorm(y, ynorm );
4515e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4525e42470aSBarry Smith     scale = maxstep/(*ynorm);
4535e42470aSBarry Smith     VecScale(&scale, y );
4545e42470aSBarry Smith     *ynorm = maxstep;
4555e76c431SBarry Smith   }
4565e42470aSBarry Smith   minlambda = steptol/(*ynorm);
4575e42470aSBarry Smith #if defined(PETSC_COMPLEX)
4585e42470aSBarry Smith   VecDot(f, y, &cinitslope );
4595e42470aSBarry Smith   initslope = real(cinitslope);
4605e42470aSBarry Smith #else
4615e42470aSBarry Smith   VecDot(f, y, &initslope );
4625e42470aSBarry Smith #endif
4635e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
4645e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
4655e42470aSBarry Smith 
4665e42470aSBarry Smith   VecCopy(y, w );
4675e42470aSBarry Smith   VecAXPY(&one, x, w );
46878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
4695e42470aSBarry Smith   VecNorm(g, gnorm );
4705e42470aSBarry Smith   if (*gnorm <= fnorm + alpha*initslope) {	/* Sufficient reduction */
4715e42470aSBarry Smith       VecCopy(w, y );
4725e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Using full step\n");
473d93a2b8dSBarry Smith       goto theend;
4745e42470aSBarry Smith   }
4755e42470aSBarry Smith 
4765e42470aSBarry Smith   /* Fit points with quadratic */
4775e42470aSBarry Smith   lambda = 1.0; count = 0;
4785e42470aSBarry Smith   count = 1;
4795e42470aSBarry Smith   while (1) {
4805e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4815e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Unable to find good step length! %d \n",count);
4825e42470aSBarry Smith       PLogInfo((PetscObject)snes, "f %g fnew %g ynorm %g lambda %g \n",
4835e42470aSBarry Smith                    fnorm,*gnorm, *ynorm,lambda);
4845e42470aSBarry Smith       VecCopy(w, y );
485761aaf1bSLois Curfman McInnes       *flag = -1; break;
4865e42470aSBarry Smith     }
4875e42470aSBarry Smith     lambdatemp = -initslope/(2.0*(*gnorm - fnorm - initslope));
4885e42470aSBarry Smith     lambdaprev = lambda;
4895e42470aSBarry Smith     gnormprev = *gnorm;
4905e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
4915e42470aSBarry Smith       lambda = .1*lambda;
4925e42470aSBarry Smith     } else lambda = lambdatemp;
4935e42470aSBarry Smith     VecCopy(x, w );
4945e42470aSBarry Smith #if defined(PETSC_COMPLEX)
4955e42470aSBarry Smith     clambda = lambda; VecAXPY(&clambda, y, w );
4965e42470aSBarry Smith #else
4975e42470aSBarry Smith     VecAXPY(&lambda, y, w );
4985e42470aSBarry Smith #endif
49978b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
5005e42470aSBarry Smith     VecNorm(g, gnorm );
5015e42470aSBarry Smith     if (*gnorm <= fnorm + alpha*initslope) {      /* sufficient reduction */
5025e42470aSBarry Smith       VecCopy(w, y );
5035e42470aSBarry Smith       PLogInfo((PetscObject)snes,"Quadratically determined step, lambda %g\n",lambda);
50406259719SBarry Smith       break;
5055e42470aSBarry Smith     }
5065e42470aSBarry Smith     count++;
5075e42470aSBarry Smith   }
508d93a2b8dSBarry Smith   theend:
5097857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5105e42470aSBarry Smith   return 0;
5115e76c431SBarry Smith }
5125e76c431SBarry Smith /* ------------------------------------------------------------ */
513b1f0a012SBarry Smith /*@
5145e42470aSBarry Smith    SNESSetLineSearchRoutine - Sets the line search routine to be used
515fbe28522SBarry Smith    by the method SNES_LS.
5165e76c431SBarry Smith 
5175e76c431SBarry Smith    Input Parameters:
518eafb4bcbSBarry Smith .  snes - nonlinear context obtained from SNESCreate()
5195e76c431SBarry Smith .  func - pointer to int function
5205e76c431SBarry Smith 
521c4a48953SLois Curfman McInnes    Available Routines:
522f59ffdedSLois Curfman McInnes .  SNESCubicLineSearch() - default line search
523f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
524f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
5255e76c431SBarry Smith 
526c4a48953SLois Curfman McInnes     Options Database Keys:
527c4a48953SLois Curfman McInnes $   -snes_line_search [basic,quadratic,cubic]
528c4a48953SLois Curfman McInnes 
5295e76c431SBarry Smith    Calling sequence of func:
530f59ffdedSLois Curfman McInnes    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
531761aaf1bSLois Curfman McInnes          Vec w, double fnorm, double *ynorm,
532761aaf1bSLois Curfman McInnes          double *gnorm, *flag)
5335e76c431SBarry Smith 
5345e76c431SBarry Smith     Input parameters for func:
5355e42470aSBarry Smith .   snes - nonlinear context
5365e76c431SBarry Smith .   x - current iterate
5375e76c431SBarry Smith .   f - residual evaluated at x
5385e76c431SBarry Smith .   y - search direction (contains new iterate on output)
5395e76c431SBarry Smith .   w - work vector
5405e76c431SBarry Smith .   fnorm - 2-norm of f
5415e76c431SBarry Smith 
5425e76c431SBarry Smith     Output parameters for func:
5435e76c431SBarry Smith .   g - residual evaluated at new iterate y
5445e76c431SBarry Smith .   y - new iterate (contains search direction on input)
5455e76c431SBarry Smith .   gnorm - 2-norm of g
5465e76c431SBarry Smith .   ynorm - 2-norm of search length
547761aaf1bSLois Curfman McInnes .   flag - set to 0 if the line search succeeds; a nonzero integer
548761aaf1bSLois Curfman McInnes            on failure.
549f59ffdedSLois Curfman McInnes 
550f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
551f59ffdedSLois Curfman McInnes 
552f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
5535e76c431SBarry Smith @*/
5545e42470aSBarry Smith int SNESSetLineSearchRoutine(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
555761aaf1bSLois Curfman McInnes                              double,double *,double*,int*) )
5565e76c431SBarry Smith {
557fbe28522SBarry Smith   if ((snes)->type == SNES_NLS)
5585e42470aSBarry Smith     ((SNES_LS *)(snes->data))->LineSearch = func;
5595e42470aSBarry Smith   return 0;
5605e76c431SBarry Smith }
5615e42470aSBarry Smith 
5625e42470aSBarry Smith static int SNESPrintHelp_LS(SNES snes)
5635e42470aSBarry Smith {
5645e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
56525c7317fSLois Curfman McInnes   MPIU_fprintf(snes->comm,stdout," method ls:\n");
56625c7317fSLois Curfman McInnes   MPIU_fprintf(snes->comm,stdout,"   -snes_line_search [basic,quadratic,cubic]\n");
56725c7317fSLois Curfman McInnes   MPIU_fprintf(snes->comm,stdout,"   -snes_line_search_alpha alpha (default %g)\n",ls->alpha);
56825c7317fSLois Curfman McInnes   MPIU_fprintf(snes->comm,stdout,"   -snes_line_search_maxstep max (default %g)\n",ls->maxstep);
56925c7317fSLois Curfman McInnes   MPIU_fprintf(snes->comm,stdout,"   -snes_line_search_steptol tol (default %g)\n",ls->steptol);
5706b5873e3SBarry Smith   return 0;
5715e42470aSBarry Smith }
5725e42470aSBarry Smith 
573a935fc98SLois Curfman McInnes static int SNESView_LS(PetscObject obj,Viewer viewer)
574a935fc98SLois Curfman McInnes {
575a935fc98SLois Curfman McInnes   SNES    snes = (SNES)obj;
576a935fc98SLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
577a935fc98SLois Curfman McInnes   FILE    *fd = ViewerFileGetPointer_Private(viewer);
578a935fc98SLois Curfman McInnes   char    *cstring;
579a935fc98SLois Curfman McInnes 
580a935fc98SLois Curfman McInnes   if (ls->LineSearch == SNESNoLineSearch)
581a935fc98SLois Curfman McInnes     cstring = "SNESNoLineSearch";
582a935fc98SLois Curfman McInnes   else if (ls->LineSearch == SNESQuadraticLineSearch)
583a935fc98SLois Curfman McInnes     cstring = "SNESQuadraticLineSearch";
584a935fc98SLois Curfman McInnes   else if (ls->LineSearch == SNESCubicLineSearch)
585a935fc98SLois Curfman McInnes     cstring = "SNESCubicLineSearch";
586a935fc98SLois Curfman McInnes   else
587a935fc98SLois Curfman McInnes     cstring = "unknown";
588a935fc98SLois Curfman McInnes   MPIU_fprintf(snes->comm,fd,"    line search variant: %s\n",cstring);
589a935fc98SLois Curfman McInnes   MPIU_fprintf(snes->comm,fd,"    alpha=%g, maxstep=%g, steptol=%g\n",
590a935fc98SLois Curfman McInnes      ls->alpha,ls->maxstep,ls->steptol);
591a935fc98SLois Curfman McInnes   return 0;
592a935fc98SLois Curfman McInnes }
593a935fc98SLois Curfman McInnes 
5945e42470aSBarry Smith static int SNESSetFromOptions_LS(SNES snes)
5955e42470aSBarry Smith {
5965e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
5975e42470aSBarry Smith   char    ver[16];
5985e42470aSBarry Smith   double  tmp;
5995e42470aSBarry Smith 
600df60cc22SBarry Smith   if (OptionsGetDouble(snes->prefix,"-snes_line_search_alpa",&tmp)) {
6015e42470aSBarry Smith     ls->alpha = tmp;
6025e42470aSBarry Smith   }
603df60cc22SBarry Smith   if (OptionsGetDouble(snes->prefix,"-snes_line_search_maxstep",&tmp)) {
6045e42470aSBarry Smith     ls->maxstep = tmp;
6055e42470aSBarry Smith   }
606df60cc22SBarry Smith   if (OptionsGetDouble(snes->prefix,"-snes_line_search_steptol",&tmp)) {
6075e42470aSBarry Smith     ls->steptol = tmp;
6085e42470aSBarry Smith   }
609df60cc22SBarry Smith   if (OptionsGetString(snes->prefix,"-snes_line_search",ver,16)) {
6105e42470aSBarry Smith     if (!strcmp(ver,"basic")) {
6115e42470aSBarry Smith       SNESSetLineSearchRoutine(snes,SNESNoLineSearch);
6125e42470aSBarry Smith     }
6135e42470aSBarry Smith     else if (!strcmp(ver,"quadratic")) {
6145e42470aSBarry Smith       SNESSetLineSearchRoutine(snes,SNESQuadraticLineSearch);
6155e42470aSBarry Smith     }
6165e42470aSBarry Smith     else if (!strcmp(ver,"cubic")) {
6175e42470aSBarry Smith       SNESSetLineSearchRoutine(snes,SNESCubicLineSearch);
6185e42470aSBarry Smith     }
6198c48e012SBarry Smith     else {SETERRQ(1,"SNESSetFromOptions_LS: Unknown line search");}
6205e42470aSBarry Smith   }
6215e42470aSBarry Smith   return 0;
6225e42470aSBarry Smith }
6235e42470aSBarry Smith 
6245e42470aSBarry Smith /* ------------------------------------------------------------ */
6255e42470aSBarry Smith int SNESCreate_LS(SNES  snes )
6265e42470aSBarry Smith {
6275e42470aSBarry Smith   SNES_LS *neP;
6285e42470aSBarry Smith 
629*1a3481a4SLois Curfman McInnes   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,
630*1a3481a4SLois Curfman McInnes     "SNESCreate_LS: Valid for SNES_NONLINEAR_EQUATIONS problems only");
63125c7317fSLois Curfman McInnes   snes->type		= SNES_EQ_NLS;
63249e3953dSBarry Smith   snes->setup		= SNESSetUp_LS;
63349e3953dSBarry Smith   snes->solve		= SNESSolve_LS;
6345e42470aSBarry Smith   snes->destroy		= SNESDestroy_LS;
635bbb6d6a8SBarry Smith   snes->converged	= SNESDefaultConverged;
63649e3953dSBarry Smith   snes->printhelp       = SNESPrintHelp_LS;
63749e3953dSBarry Smith   snes->setfromoptions  = SNESSetFromOptions_LS;
638a935fc98SLois Curfman McInnes   snes->view            = SNESView_LS;
6395e42470aSBarry Smith 
64078b31e54SBarry Smith   neP			= PETSCNEW(SNES_LS);   CHKPTRQ(neP);
6415e42470aSBarry Smith   snes->data    	= (void *) neP;
6425e42470aSBarry Smith   neP->alpha		= 1.e-4;
6435e42470aSBarry Smith   neP->maxstep		= 1.e8;
6445e42470aSBarry Smith   neP->steptol		= 1.e-12;
6455e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
6465e42470aSBarry Smith   return 0;
6475e42470aSBarry Smith }
6485e42470aSBarry Smith 
6495e42470aSBarry Smith 
6505e42470aSBarry Smith 
6515e42470aSBarry Smith 
652