xref: /petsc/src/snes/impls/ls/ls.c (revision 4ccb76dd2a06ade17c2800eec40bc3f6b9abd66f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*4ccb76ddSBarry Smith static char vcid[] = "$Id: ls.c,v 1.137 1999/06/08 22:53:15 bsmith Exp bsmith $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
570f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
65e42470aSBarry Smith 
704d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
804d965bbSLois Curfman McInnes 
904d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
1004d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
1104d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
1204d965bbSLois Curfman McInnes      respectively.
1304d965bbSLois Curfman McInnes 
14fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
1504d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
1604d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
17fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
1804d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
1904d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
2004d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
2104d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
2204d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
2304d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
2404d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
2504d965bbSLois Curfman McInnes      for all nonlinear solvers.
2604d965bbSLois Curfman McInnes 
2704d965bbSLois Curfman McInnes      Another key routine is:
2804d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
2904d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
3004d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
3104d965bbSLois Curfman McInnes      SNESSolve() if necessary.
3204d965bbSLois Curfman McInnes 
3304d965bbSLois Curfman McInnes      Additional basic routines are:
3404d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
3504d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
3604d965bbSLois Curfman McInnes                                       have actually been used.
3704d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
3804d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
3904d965bbSLois Curfman McInnes 
4004d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
4104d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
4204d965bbSLois Curfman McInnes      above description applies to these categories also.
4304d965bbSLois Curfman McInnes 
4404d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
455e42470aSBarry Smith /*
4604d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
4704d965bbSLois Curfman McInnes    method with a line search.
485e76c431SBarry Smith 
4904d965bbSLois Curfman McInnes    Input Parameters:
5004d965bbSLois Curfman McInnes .  snes - the SNES context
515e76c431SBarry Smith 
5204d965bbSLois Curfman McInnes    Output Parameter:
53c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
5404d965bbSLois Curfman McInnes 
5504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
565e76c431SBarry Smith 
575e76c431SBarry Smith    Notes:
585e76c431SBarry Smith    This implements essentially a truncated Newton method with a
595e76c431SBarry Smith    line search.  By default a cubic backtracking line search
605e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
615e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
62393d2d9aSLois Curfman McInnes    and Schnabel.
635e76c431SBarry Smith */
645615d1e5SSatish Balay #undef __FUNC__
655615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS"
66f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
675e76c431SBarry Smith {
685e42470aSBarry Smith   SNES_LS       *neP = (SNES_LS *) snes->data;
6984c569c9SBarry Smith   int           maxits, i, ierr, lits, lsfail;
70112a2221SBarry Smith   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
7184c569c9SBarry Smith   double        fnorm, gnorm, xnorm, ynorm;
725e42470aSBarry Smith   Vec           Y, X, F, G, W, TMP;
735e76c431SBarry Smith 
743a40ed3dSBarry Smith   PetscFunctionBegin;
755e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
765e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
7739e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
785e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
795e42470aSBarry Smith   G		= snes->work[1];
805e42470aSBarry Smith   W		= snes->work[2];
815e76c431SBarry Smith 
825334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
83cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
8406d0569aSBarry Smith   PetscAMSTakeAccess(snes);
8584c569c9SBarry Smith   snes->iter = 0;
865e42470aSBarry Smith   snes->norm = fnorm;
8706d0569aSBarry Smith   PetscAMSGrantAccess(snes);
8884c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
8994a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
905e76c431SBarry Smith 
913a40ed3dSBarry Smith   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
9294a9d846SBarry Smith 
93d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
94d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
95d034289bSBarry Smith 
965e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
975e76c431SBarry Smith 
98ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
995334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1005334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
10178b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
1027a00f4a9SLois Curfman McInnes     snes->linear_its += PetscAbsInt(lits);
103af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
104ea4d3ed3SLois Curfman McInnes 
105ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
106ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
107ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
108ea4d3ed3SLois Curfman McInnes     */
10981b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
110af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
111af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
112761aaf1bSLois Curfman McInnes     if (lsfail) snes->nfailures++;
1135e76c431SBarry Smith 
11439e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
11539e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1165e76c431SBarry Smith     fnorm = gnorm;
1175e76c431SBarry Smith 
11806d0569aSBarry Smith     PetscAMSTakeAccess(snes);
11984c569c9SBarry Smith     snes->iter = i+1;
1205e42470aSBarry Smith     snes->norm = fnorm;
12106d0569aSBarry Smith     PetscAMSGrantAccess(snes);
12284c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
12394a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1245e76c431SBarry Smith 
1255e76c431SBarry Smith     /* Test for convergence */
12629e0b56fSBarry Smith     if (snes->converged) {
12729e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
128bbb6d6a8SBarry Smith       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
12916c95cacSBarry Smith         break;
13016c95cacSBarry Smith       }
13116c95cacSBarry Smith     }
13229e0b56fSBarry Smith   }
13339e2f89bSBarry Smith   if (X != snes->vec_sol) {
134393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
13539e2f89bSBarry Smith     snes->vec_sol_always  = snes->vec_sol;
13639e2f89bSBarry Smith     snes->vec_func_always = snes->vec_func;
13739e2f89bSBarry Smith   }
13852392280SLois Curfman McInnes   if (i == maxits) {
139981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
14052392280SLois Curfman McInnes     i--;
14152392280SLois Curfman McInnes   }
1425e42470aSBarry Smith   *outits = i+1;
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1445e76c431SBarry Smith }
14504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
14604d965bbSLois Curfman McInnes /*
14704d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
14804d965bbSLois Curfman McInnes    of the SNES_EQ_LS nonlinear solver.
14904d965bbSLois Curfman McInnes 
15004d965bbSLois Curfman McInnes    Input Parameter:
15104d965bbSLois Curfman McInnes .  snes - the SNES context
15204d965bbSLois Curfman McInnes .  x - the solution vector
15304d965bbSLois Curfman McInnes 
15404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
15504d965bbSLois Curfman McInnes 
15604d965bbSLois Curfman McInnes    Notes:
15704d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
15804d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
15904d965bbSLois Curfman McInnes    the call to SNESSolve().
16004d965bbSLois Curfman McInnes  */
1615615d1e5SSatish Balay #undef __FUNC__
1625615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
163f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
1645e76c431SBarry Smith {
1655e42470aSBarry Smith   int ierr;
1663a40ed3dSBarry Smith 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
16881b6cf68SLois Curfman McInnes   snes->nwork = 4;
169d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1705e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
17181b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1735e76c431SBarry Smith }
17404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
17504d965bbSLois Curfman McInnes /*
17604d965bbSLois Curfman McInnes    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
17704d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
17804d965bbSLois Curfman McInnes 
17904d965bbSLois Curfman McInnes    Input Parameter:
18004d965bbSLois Curfman McInnes .  snes - the SNES context
18104d965bbSLois Curfman McInnes 
182de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
18304d965bbSLois Curfman McInnes  */
1845615d1e5SSatish Balay #undef __FUNC__
185d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
186e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
1875e76c431SBarry Smith {
188393d2d9aSLois Curfman McInnes   int  ierr;
1893a40ed3dSBarry Smith 
1903a40ed3dSBarry Smith   PetscFunctionBegin;
1915baf8537SBarry Smith   if (snes->nwork) {
1924b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
19321c89e3eSBarry Smith   }
1940452661fSBarry Smith   PetscFree(snes->data);
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1965e76c431SBarry Smith }
19704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1985615d1e5SSatish Balay #undef __FUNC__
1995615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
20082bf6240SBarry Smith 
2014b828684SBarry Smith /*@C
2025e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2035e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2045e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2055e76c431SBarry Smith 
206c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
207c7afd0dbSLois Curfman McInnes 
2085e76c431SBarry Smith    Input Parameters:
209c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
210af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
2115e76c431SBarry Smith .  x - current iterate
2125e76c431SBarry Smith .  f - residual evaluated at x
2135e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2145e76c431SBarry Smith .  w - work vector
215c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
2165e76c431SBarry Smith 
217c4a48953SLois Curfman McInnes    Output Parameters:
218c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
2195e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2205e76c431SBarry Smith .  gnorm - 2-norm of g
2215e76c431SBarry Smith .  ynorm - 2-norm of search length
222c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
223fee21e36SBarry Smith 
224c4a48953SLois Curfman McInnes    Options Database Key:
225c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
226c4a48953SLois Curfman McInnes 
22736851e7fSLois Curfman McInnes    Level: advanced
22836851e7fSLois Curfman McInnes 
22928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
23028ae5a21SLois Curfman McInnes 
231f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
232af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
2335e76c431SBarry Smith @*/
234af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
235761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
2365e76c431SBarry Smith {
2375e42470aSBarry Smith   int        ierr;
2385334005bSBarry Smith   Scalar     mone = -1.0;
2398f99978dSLois Curfman McInnes   SNES_LS    *neP = (SNES_LS *) snes->data;
2408f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
2415334005bSBarry Smith 
2423a40ed3dSBarry Smith   PetscFunctionBegin;
243761aaf1bSLois Curfman McInnes   *flag = 0;
2447857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
245ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
2468f99978dSLois Curfman McInnes   if (neP->CheckStep) {
2478f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
2488f99978dSLois Curfman McInnes   }
249ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
2508f99978dSLois Curfman McInnes   ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
2518f99978dSLois Curfman McInnes   ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
2528f99978dSLois Curfman McInnes   ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
2538f99978dSLois Curfman McInnes   ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
2547857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2553a40ed3dSBarry Smith   PetscFunctionReturn(0);
2565e76c431SBarry Smith }
25704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2585615d1e5SSatish Balay #undef __FUNC__
25929e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
26082bf6240SBarry Smith 
26129e0b56fSBarry Smith /*@C
26229e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
26329e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
26429e0b56fSBarry Smith    even compute the norm of the function or search direction; this
26529e0b56fSBarry Smith    is intended only when you know the full step is fine and are
26629e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
267c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
268c7afd0dbSLois Curfman McInnes 
269c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
27029e0b56fSBarry Smith 
27129e0b56fSBarry Smith    Input Parameters:
272c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
273af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
27429e0b56fSBarry Smith .  x - current iterate
27529e0b56fSBarry Smith .  f - residual evaluated at x
27629e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
27729e0b56fSBarry Smith .  w - work vector
278c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
27929e0b56fSBarry Smith 
28029e0b56fSBarry Smith    Output Parameters:
281c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
28229e0b56fSBarry Smith .  gnorm - not changed
28329e0b56fSBarry Smith .  ynorm - not changed
284c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
285fee21e36SBarry Smith 
28629e0b56fSBarry Smith    Options Database Key:
287c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
28829e0b56fSBarry Smith 
2898cbba510SBarry Smith    Notes:
290ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
291ea56f5baSLois Curfman McInnes    either the options
292ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
293ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
294ea56f5baSLois Curfman McInnes    SNESSetConvergenceTest(); otherwise, the SNES solver will
295ea56f5baSLois Curfman McInnes    generate an error.
2968cbba510SBarry Smith 
297ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
298ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
299ea56f5baSLois Curfman McInnes    correct, since they are not computed.
300ea56f5baSLois Curfman McInnes 
301ea56f5baSLois Curfman McInnes    Level: advanced
3028cbba510SBarry Smith 
30329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
30429e0b56fSBarry Smith 
30529e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
30629e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
30729e0b56fSBarry Smith @*/
308af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
30929e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
31029e0b56fSBarry Smith {
31129e0b56fSBarry Smith   int        ierr;
31229e0b56fSBarry Smith   Scalar     mone = -1.0;
3138f99978dSLois Curfman McInnes   SNES_LS    *neP = (SNES_LS *) snes->data;
3148f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
31529e0b56fSBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
3178cbba510SBarry Smith   *flag = 0;
31829e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
31929e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3208f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3218f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3228f99978dSLois Curfman McInnes   }
32329e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
32429e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3253a40ed3dSBarry Smith   PetscFunctionReturn(0);
32629e0b56fSBarry Smith }
32704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
32829e0b56fSBarry Smith #undef __FUNC__
3295615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3304b828684SBarry Smith /*@C
331f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3325e76c431SBarry Smith 
333c7afd0dbSLois Curfman McInnes    Collective on SNES
334c7afd0dbSLois Curfman McInnes 
3355e76c431SBarry Smith    Input Parameters:
336c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
337af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3385e76c431SBarry Smith .  x - current iterate
3395e76c431SBarry Smith .  f - residual evaluated at x
3405e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3415e76c431SBarry Smith .  w - work vector
342c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3435e76c431SBarry Smith 
344393d2d9aSLois Curfman McInnes    Output Parameters:
345c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3465e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3475e76c431SBarry Smith .  gnorm - 2-norm of g
3485e76c431SBarry Smith .  ynorm - 2-norm of search length
349c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
350fee21e36SBarry Smith 
351c4a48953SLois Curfman McInnes    Options Database Key:
352c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
353c4a48953SLois Curfman McInnes 
3545e76c431SBarry Smith    Notes:
3555e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3565e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3575e76c431SBarry Smith 
35836851e7fSLois Curfman McInnes    Level: advanced
35936851e7fSLois Curfman McInnes 
36028ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
36128ae5a21SLois Curfman McInnes 
362af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3635e76c431SBarry Smith @*/
364af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
365761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3665e76c431SBarry Smith {
367406556e6SLois Curfman McInnes   /*
368406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
369406556e6SLois Curfman McInnes      minimization problem:
370406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
371406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
372406556e6SLois Curfman McInnes    */
373406556e6SLois Curfman McInnes 
374ddd12767SBarry Smith   double     steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
375ea4d3ed3SLois Curfman McInnes   double     maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
376aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
3775e42470aSBarry Smith   Scalar     cinitslope, clambda;
3786b5873e3SBarry Smith #endif
3795e42470aSBarry Smith   int        ierr, count;
3805e42470aSBarry Smith   SNES_LS    *neP = (SNES_LS *) snes->data;
3815334005bSBarry Smith   Scalar     mone = -1.0, scale;
3828f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
3835e76c431SBarry Smith 
3843a40ed3dSBarry Smith   PetscFunctionBegin;
3857857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
386761aaf1bSLois Curfman McInnes   *flag   = 0;
3875e76c431SBarry Smith   alpha   = neP->alpha;
3885e76c431SBarry Smith   maxstep = neP->maxstep;
3895e76c431SBarry Smith   steptol = neP->steptol;
3905e76c431SBarry Smith 
391cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
392a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
393a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
394a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
395a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
396a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
397ad922ac9SBarry Smith     goto theend1;
39894a9d846SBarry Smith   }
3995e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4005e42470aSBarry Smith     scale = maxstep/(*ynorm);
401aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
40292318cfeSSatish Balay     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
4036b5873e3SBarry Smith #else
40494a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
4056b5873e3SBarry Smith #endif
406393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
4075e76c431SBarry Smith     *ynorm = maxstep;
4085e76c431SBarry Smith   }
4095e76c431SBarry Smith   minlambda = steptol/(*ynorm);
410a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
411aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
412a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
41392318cfeSSatish Balay   initslope = PetscReal(cinitslope);
4145e42470aSBarry Smith #else
415a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
4165e42470aSBarry Smith #endif
4175e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
4185e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
4195e76c431SBarry Smith 
420393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
4215334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
42278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
423313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
424406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
425393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
42694a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
427ad922ac9SBarry Smith     goto theend1;
4285e76c431SBarry Smith   }
4295e76c431SBarry Smith 
4305e76c431SBarry Smith   /* Fit points with quadratic */
431313b5538SBarry Smith   lambda = 1.0;
432a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4335e76c431SBarry Smith   lambdaprev = lambda;
4345e76c431SBarry Smith   gnormprev = *gnorm;
435ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
436ddd12767SBarry Smith   else lambda = lambdatemp;
437393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
438ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
439aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
440ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4415e42470aSBarry Smith #else
442ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4435e42470aSBarry Smith #endif
44478b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
445cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
446406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
447393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
448ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
449ad922ac9SBarry Smith     goto theend1;
4505e76c431SBarry Smith   }
4515e76c431SBarry Smith 
4525e76c431SBarry Smith   /* Fit points with cubic */
4535e76c431SBarry Smith   count = 1;
4545e76c431SBarry Smith   while (1) {
4555e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4562b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4572b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
458ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
459393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
460761aaf1bSLois Curfman McInnes       *flag = -1; break;
4615e76c431SBarry Smith     }
462406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
463406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
464ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4652b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4665e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4675e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4685e76c431SBarry Smith     if (a == 0.0) {
4695e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4705e76c431SBarry Smith     } else {
4715e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4725e76c431SBarry Smith     }
4735e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4745e76c431SBarry Smith       lambdatemp = .5*lambda;
4755e76c431SBarry Smith     }
4765e76c431SBarry Smith     lambdaprev = lambda;
4775e76c431SBarry Smith     gnormprev = *gnorm;
4785e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4795e76c431SBarry Smith       lambda = .1*lambda;
4805e76c431SBarry Smith     }
4815e76c431SBarry Smith     else lambda = lambdatemp;
482393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w );CHKERRQ(ierr);
483ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
484aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
485ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
486393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4875e42470aSBarry Smith #else
488ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4895e42470aSBarry Smith #endif
49078b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
491cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
492406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
493393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
494ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
4952715565aSLois Curfman McInnes       goto theend1;
4962b022350SLois Curfman McInnes     } else {
4972b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
4985e76c431SBarry Smith     }
4995e76c431SBarry Smith     count++;
5005e76c431SBarry Smith   }
501ad922ac9SBarry Smith   theend1:
5028f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5038f99978dSLois Curfman McInnes   if (neP->CheckStep) {
5048f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
5058f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
5068f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
5078f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
5088f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
5098f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
5108f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
5118f99978dSLois Curfman McInnes     }
5128f99978dSLois Curfman McInnes   }
5137857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5143a40ed3dSBarry Smith   PetscFunctionReturn(0);
5155e76c431SBarry Smith }
51604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
5175615d1e5SSatish Balay #undef __FUNC__
5185615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
5194b828684SBarry Smith /*@C
520f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
5215e76c431SBarry Smith 
522c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
523c7afd0dbSLois Curfman McInnes 
5245e42470aSBarry Smith    Input Parameters:
525c7afd0dbSLois Curfman McInnes +  snes - the SNES context
526af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
5275e42470aSBarry Smith .  x - current iterate
5285e42470aSBarry Smith .  f - residual evaluated at x
5295e42470aSBarry Smith .  y - search direction (contains new iterate on output)
5305e42470aSBarry Smith .  w - work vector
531c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
5325e42470aSBarry Smith 
533c4a48953SLois Curfman McInnes    Output Parameters:
534c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
5355e42470aSBarry Smith .  y - new iterate (contains search direction on input)
5365e42470aSBarry Smith .  gnorm - 2-norm of g
5375e42470aSBarry Smith .  ynorm - 2-norm of search length
538c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
539fee21e36SBarry Smith 
540c4a48953SLois Curfman McInnes    Options Database Key:
541c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
542c4a48953SLois Curfman McInnes 
5435e42470aSBarry Smith    Notes:
544c7afd0dbSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
5455e42470aSBarry Smith 
54636851e7fSLois Curfman McInnes    Level: advanced
54736851e7fSLois Curfman McInnes 
548f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
549f59ffdedSLois Curfman McInnes 
550af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
5515e42470aSBarry Smith @*/
552af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
553761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5545e76c431SBarry Smith {
555406556e6SLois Curfman McInnes   /*
556406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
557406556e6SLois Curfman McInnes      minimization problem:
558406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
559406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
560406556e6SLois Curfman McInnes    */
56140011551SBarry Smith   double     steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
562aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5635e42470aSBarry Smith   Scalar     cinitslope,clambda;
5646b5873e3SBarry Smith #endif
5655e42470aSBarry Smith   int        ierr, count;
5665e42470aSBarry Smith   SNES_LS    *neP = (SNES_LS *) snes->data;
5675334005bSBarry Smith   Scalar     mone = -1.0,scale;
5688f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
5695e76c431SBarry Smith 
5703a40ed3dSBarry Smith   PetscFunctionBegin;
5717857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
572761aaf1bSLois Curfman McInnes   *flag   = 0;
5735e42470aSBarry Smith   alpha   = neP->alpha;
5745e42470aSBarry Smith   maxstep = neP->maxstep;
5755e42470aSBarry Smith   steptol = neP->steptol;
5765e76c431SBarry Smith 
577cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
578b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
57994a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
580b37302c6SLois Curfman McInnes     *gnorm = fnorm;
581b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
582b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
583ad922ac9SBarry Smith     goto theend2;
58494a9d846SBarry Smith   }
5855e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5865e42470aSBarry Smith     scale = maxstep/(*ynorm);
587393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5885e42470aSBarry Smith     *ynorm = maxstep;
5895e76c431SBarry Smith   }
5905e42470aSBarry Smith   minlambda = steptol/(*ynorm);
591a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
592aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
593a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
59492318cfeSSatish Balay   initslope = PetscReal(cinitslope);
5955e42470aSBarry Smith #else
596a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5975e42470aSBarry Smith #endif
5985e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
5995e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6005e42470aSBarry Smith 
601393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6025334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
60378b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
604cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
605406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
606393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
60794a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
608ad922ac9SBarry Smith     goto theend2;
6095e42470aSBarry Smith   }
6105e42470aSBarry Smith 
6115e42470aSBarry Smith   /* Fit points with quadratic */
612313b5538SBarry Smith   lambda = 1.0;
6135e42470aSBarry Smith   count = 1;
6145e42470aSBarry Smith   while (1) {
6155e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
616981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
617981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
618ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
619393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
620761aaf1bSLois Curfman McInnes       *flag = -1; break;
6215e42470aSBarry Smith     }
622a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
6235e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
6245e42470aSBarry Smith       lambda = .1*lambda;
6255e42470aSBarry Smith     } else lambda = lambdatemp;
626393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
6275334005bSBarry Smith     lambda = -lambda;
628aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
629393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6305e42470aSBarry Smith #else
631393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w);CHKERRQ(ierr);
6325e42470aSBarry Smith #endif
63378b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
634cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
635406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
636393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
637981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
63806259719SBarry Smith       break;
6395e42470aSBarry Smith     }
6405e42470aSBarry Smith     count++;
6415e42470aSBarry Smith   }
642ad922ac9SBarry Smith   theend2:
6438f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6448f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6458f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6468f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6478f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6488f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6498f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6508f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6518f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6528f99978dSLois Curfman McInnes     }
6538f99978dSLois Curfman McInnes   }
6547857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
6553a40ed3dSBarry Smith   PetscFunctionReturn(0);
6565e76c431SBarry Smith }
65704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6585615d1e5SSatish Balay #undef __FUNC__
659d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
660c9e6a524SLois Curfman McInnes /*@C
66177c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
662f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
6635e76c431SBarry Smith 
6645e76c431SBarry Smith    Input Parameters:
665c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
666af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
667c7afd0dbSLois Curfman McInnes -  func - pointer to int function
6685e76c431SBarry Smith 
669fee21e36SBarry Smith    Collective on SNES
670fee21e36SBarry Smith 
671c4a48953SLois Curfman McInnes    Available Routines:
672c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
673f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
674f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
675af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
6765e76c431SBarry Smith 
677c4a48953SLois Curfman McInnes     Options Database Keys:
678af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
679c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
680c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
681c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
682c4a48953SLois Curfman McInnes 
6835e76c431SBarry Smith    Calling sequence of func:
684bd208895SLois Curfman McInnes .vb
685af2d14d2SLois Curfman McInnes    func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
686bd208895SLois Curfman McInnes          double fnorm, double *ynorm, double *gnorm, *flag)
687bd208895SLois Curfman McInnes .ve
6885e76c431SBarry Smith 
6895e76c431SBarry Smith     Input parameters for func:
690c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
691af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
6925e76c431SBarry Smith .   x - current iterate
6935e76c431SBarry Smith .   f - residual evaluated at x
6945e76c431SBarry Smith .   y - search direction (contains new iterate on output)
6955e76c431SBarry Smith .   w - work vector
696c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
6975e76c431SBarry Smith 
6985e76c431SBarry Smith     Output parameters for func:
699c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7005e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7015e76c431SBarry Smith .   gnorm - 2-norm of g
7025e76c431SBarry Smith .   ynorm - 2-norm of search length
703c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
704761aaf1bSLois Curfman McInnes            on failure.
705f59ffdedSLois Curfman McInnes 
70636851e7fSLois Curfman McInnes     Level: advanced
70736851e7fSLois Curfman McInnes 
708f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
709f59ffdedSLois Curfman McInnes 
710*4ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
7115e76c431SBarry Smith @*/
712af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx)
7135e76c431SBarry Smith {
714af2d14d2SLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*);
71582bf6240SBarry Smith 
7163a40ed3dSBarry Smith   PetscFunctionBegin;
71794d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
71882bf6240SBarry Smith   if (f) {
719af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
72082bf6240SBarry Smith   }
7213a40ed3dSBarry Smith   PetscFunctionReturn(0);
7225e76c431SBarry Smith }
72304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
724fb2e594dSBarry Smith EXTERN_C_BEGIN
72582bf6240SBarry Smith #undef __FUNC__
72682bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
727af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
728af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
72982bf6240SBarry Smith {
73082bf6240SBarry Smith   PetscFunctionBegin;
73182bf6240SBarry Smith   ((SNES_LS *)(snes->data))->LineSearch = func;
732af2d14d2SLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP = lsctx;
73382bf6240SBarry Smith   PetscFunctionReturn(0);
73482bf6240SBarry Smith }
735fb2e594dSBarry Smith EXTERN_C_END
73604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
737c8dd0c0dSLois Curfman McInnes #undef __FUNC__
738c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck"
739c8dd0c0dSLois Curfman McInnes /*@C
740530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
7418f99978dSLois Curfman McInnes    by the line search routine in the Newton-based method SNES_EQ_LS.
742c8dd0c0dSLois Curfman McInnes 
743c8dd0c0dSLois Curfman McInnes    Input Parameters:
744c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
745c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
746c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
747c8dd0c0dSLois Curfman McInnes 
748c8dd0c0dSLois Curfman McInnes    Collective on SNES
749c8dd0c0dSLois Curfman McInnes 
750c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
751c8dd0c0dSLois Curfman McInnes .vb
752b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
753c8dd0c0dSLois Curfman McInnes .ve
754b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
755b1ae27eaSLois Curfman McInnes    on failure.
756c8dd0c0dSLois Curfman McInnes 
757c8dd0c0dSLois Curfman McInnes    Input parameters for func:
758c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
759c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
7608f99978dSLois Curfman McInnes -  x - current candidate iterate
761c8dd0c0dSLois Curfman McInnes 
762c8dd0c0dSLois Curfman McInnes    Output parameters for func:
763c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
764c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
765c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
766c8dd0c0dSLois Curfman McInnes 
767c8dd0c0dSLois Curfman McInnes    Level: advanced
768c8dd0c0dSLois Curfman McInnes 
7698f99978dSLois Curfman McInnes    Notes:
770b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
771b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
772b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
773b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
774ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
7758f99978dSLois Curfman McInnes 
776b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
777b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
778b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
779ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
780ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
781ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
782ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
783b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
7848f99978dSLois Curfman McInnes 
785c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
786c8dd0c0dSLois Curfman McInnes 
787c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
788c8dd0c0dSLois Curfman McInnes @*/
7898f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
790c8dd0c0dSLois Curfman McInnes {
7918f99978dSLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
792c8dd0c0dSLois Curfman McInnes 
793c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
794c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
795c8dd0c0dSLois Curfman McInnes   if (f) {
796c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
797c8dd0c0dSLois Curfman McInnes   }
798c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
799c8dd0c0dSLois Curfman McInnes }
800c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
801c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
802c8dd0c0dSLois Curfman McInnes #undef __FUNC__
803c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS"
8048f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
805c8dd0c0dSLois Curfman McInnes {
806c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
807c8dd0c0dSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
808c8dd0c0dSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP = checkctx;
809c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
810c8dd0c0dSLois Curfman McInnes }
811c8dd0c0dSLois Curfman McInnes EXTERN_C_END
812c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
81304d965bbSLois Curfman McInnes /*
81404d965bbSLois Curfman McInnes    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
81582bf6240SBarry Smith 
81604d965bbSLois Curfman McInnes    Input Parameter:
81704d965bbSLois Curfman McInnes .  snes - the SNES context
81804d965bbSLois Curfman McInnes 
81904d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
82004d965bbSLois Curfman McInnes */
8215615d1e5SSatish Balay #undef __FUNC__
822d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
823f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
8245e42470aSBarry Smith {
8255e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
826d132466eSBarry Smith   int     ierr;
8276daaf66cSBarry Smith 
8283a40ed3dSBarry Smith   PetscFunctionBegin;
829d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr);
830d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr);
831d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr);
832d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr);
833d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr);
8343a40ed3dSBarry Smith   PetscFunctionReturn(0);
8355e42470aSBarry Smith }
83604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
83704d965bbSLois Curfman McInnes /*
838de49b36dSLois Curfman McInnes    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
83904d965bbSLois Curfman McInnes 
84004d965bbSLois Curfman McInnes    Input Parameters:
84104d965bbSLois Curfman McInnes .  SNES - the SNES context
84204d965bbSLois Curfman McInnes .  viewer - visualization context
84304d965bbSLois Curfman McInnes 
84404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
84504d965bbSLois Curfman McInnes */
8465615d1e5SSatish Balay #undef __FUNC__
847d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
848e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
849a935fc98SLois Curfman McInnes {
850a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
85119bcc07fSBarry Smith   char       *cstr;
85251695f54SLois Curfman McInnes   int        ierr;
85319bcc07fSBarry Smith   ViewerType vtype;
854a935fc98SLois Curfman McInnes 
8553a40ed3dSBarry Smith   PetscFunctionBegin;
85619bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
8573f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
85819bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
85919bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
86019bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
86119bcc07fSBarry Smith     else                                                cstr = "unknown";
862d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
863d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
8645cd90555SBarry Smith   } else {
8655cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
86619bcc07fSBarry Smith   }
8673a40ed3dSBarry Smith   PetscFunctionReturn(0);
868a935fc98SLois Curfman McInnes }
86904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
87004d965bbSLois Curfman McInnes /*
87104d965bbSLois Curfman McInnes    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
87204d965bbSLois Curfman McInnes 
87304d965bbSLois Curfman McInnes    Input Parameter:
87404d965bbSLois Curfman McInnes .  snes - the SNES context
87504d965bbSLois Curfman McInnes 
87604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
87704d965bbSLois Curfman McInnes */
8785615d1e5SSatish Balay #undef __FUNC__
8795615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
880f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
8815e42470aSBarry Smith {
8825e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
8835e42470aSBarry Smith   char    ver[16];
8845e42470aSBarry Smith   double  tmp;
88525cdf11fSBarry Smith   int     ierr,flg;
8865e42470aSBarry Smith 
8873a40ed3dSBarry Smith   PetscFunctionBegin;
88809d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
88925cdf11fSBarry Smith   if (flg) {
8905e42470aSBarry Smith     ls->alpha = tmp;
8915e42470aSBarry Smith   }
89209d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
89325cdf11fSBarry Smith   if (flg) {
8945e42470aSBarry Smith     ls->maxstep = tmp;
8955e42470aSBarry Smith   }
89609d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
89725cdf11fSBarry Smith   if (flg) {
8985e42470aSBarry Smith     ls->steptol = tmp;
8995e42470aSBarry Smith   }
90009d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg);CHKERRQ(ierr);
90125cdf11fSBarry Smith   if (flg) {
90248d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
903af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
904b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"basicnonorms")) {
905af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
906b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"quadratic")) {
907af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
908b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"cubic")) {
909af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
9105e42470aSBarry Smith     }
911a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
9125e42470aSBarry Smith   }
9133a40ed3dSBarry Smith   PetscFunctionReturn(0);
9145e42470aSBarry Smith }
91504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
91604d965bbSLois Curfman McInnes /*
91704d965bbSLois Curfman McInnes    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
91804d965bbSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
91904d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
92004d965bbSLois Curfman McInnes 
92104d965bbSLois Curfman McInnes 
92204d965bbSLois Curfman McInnes    Input Parameter:
92304d965bbSLois Curfman McInnes .  snes - the SNES context
92404d965bbSLois Curfman McInnes 
92504d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
92604d965bbSLois Curfman McInnes  */
927fb2e594dSBarry Smith EXTERN_C_BEGIN
9285615d1e5SSatish Balay #undef __FUNC__
9295615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
930f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
9315e42470aSBarry Smith {
93282bf6240SBarry Smith   int     ierr;
9335e42470aSBarry Smith   SNES_LS *neP;
9345e42470aSBarry Smith 
9353a40ed3dSBarry Smith   PetscFunctionBegin;
936a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
937a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
938a8c6a408SBarry Smith   }
93982bf6240SBarry Smith 
940f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
941f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
942f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
943f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
944f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
945f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
946f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
9475baf8537SBarry Smith   snes->nwork           = 0;
9485e42470aSBarry Smith 
9490452661fSBarry Smith   neP			= PetscNew(SNES_LS);CHKPTRQ(neP);
950ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
9515e42470aSBarry Smith   snes->data    	= (void *) neP;
9525e42470aSBarry Smith   neP->alpha		= 1.e-4;
9535e42470aSBarry Smith   neP->maxstep		= 1.e8;
9545e42470aSBarry Smith   neP->steptol		= 1.e-12;
9555e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
956c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
957c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
958c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
95982bf6240SBarry Smith 
96094d884c6SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
961e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
962c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
963c8dd0c0dSLois Curfman McInnes                     (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
96482bf6240SBarry Smith 
9653a40ed3dSBarry Smith   PetscFunctionReturn(0);
9665e42470aSBarry Smith }
967fb2e594dSBarry Smith EXTERN_C_END
96882bf6240SBarry Smith 
96982bf6240SBarry Smith 
97082bf6240SBarry Smith 
971