xref: /petsc/src/snes/impls/ls/ls.c (revision a3b388050d4e8a4c38128e174a7572fa24ee1ba5)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*a3b38805SLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.138 1999/06/08 22:54:40 bsmith Exp curfman $";
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);
245*a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
246ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
2478f99978dSLois Curfman McInnes   if (neP->CheckStep) {
2488f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
2498f99978dSLois Curfman McInnes   }
250ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
251*a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
2527857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
2545e76c431SBarry Smith }
25504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2565615d1e5SSatish Balay #undef __FUNC__
25729e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
25882bf6240SBarry Smith 
25929e0b56fSBarry Smith /*@C
26029e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
26129e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
26229e0b56fSBarry Smith    even compute the norm of the function or search direction; this
26329e0b56fSBarry Smith    is intended only when you know the full step is fine and are
26429e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
265c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
266c7afd0dbSLois Curfman McInnes 
267c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
26829e0b56fSBarry Smith 
26929e0b56fSBarry Smith    Input Parameters:
270c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
271af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
27229e0b56fSBarry Smith .  x - current iterate
27329e0b56fSBarry Smith .  f - residual evaluated at x
27429e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
27529e0b56fSBarry Smith .  w - work vector
276c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
27729e0b56fSBarry Smith 
27829e0b56fSBarry Smith    Output Parameters:
279c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
28029e0b56fSBarry Smith .  gnorm - not changed
28129e0b56fSBarry Smith .  ynorm - not changed
282c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
283fee21e36SBarry Smith 
28429e0b56fSBarry Smith    Options Database Key:
285c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
28629e0b56fSBarry Smith 
2878cbba510SBarry Smith    Notes:
288ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
289ea56f5baSLois Curfman McInnes    either the options
290ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
291ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
292ea56f5baSLois Curfman McInnes    SNESSetConvergenceTest(); otherwise, the SNES solver will
293ea56f5baSLois Curfman McInnes    generate an error.
2948cbba510SBarry Smith 
295ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
296ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
297ea56f5baSLois Curfman McInnes    correct, since they are not computed.
298ea56f5baSLois Curfman McInnes 
299ea56f5baSLois Curfman McInnes    Level: advanced
3008cbba510SBarry Smith 
30129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
30229e0b56fSBarry Smith 
30329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
30429e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
30529e0b56fSBarry Smith @*/
306af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
30729e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
30829e0b56fSBarry Smith {
30929e0b56fSBarry Smith   int        ierr;
31029e0b56fSBarry Smith   Scalar     mone = -1.0;
3118f99978dSLois Curfman McInnes   SNES_LS    *neP = (SNES_LS *) snes->data;
3128f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
31329e0b56fSBarry Smith 
3143a40ed3dSBarry Smith   PetscFunctionBegin;
3158cbba510SBarry Smith   *flag = 0;
31629e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
31729e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3188f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3198f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3208f99978dSLois Curfman McInnes   }
32129e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
32229e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
32429e0b56fSBarry Smith }
32504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
32629e0b56fSBarry Smith #undef __FUNC__
3275615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3284b828684SBarry Smith /*@C
329f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3305e76c431SBarry Smith 
331c7afd0dbSLois Curfman McInnes    Collective on SNES
332c7afd0dbSLois Curfman McInnes 
3335e76c431SBarry Smith    Input Parameters:
334c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
335af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3365e76c431SBarry Smith .  x - current iterate
3375e76c431SBarry Smith .  f - residual evaluated at x
3385e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3395e76c431SBarry Smith .  w - work vector
340c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3415e76c431SBarry Smith 
342393d2d9aSLois Curfman McInnes    Output Parameters:
343c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3445e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3455e76c431SBarry Smith .  gnorm - 2-norm of g
3465e76c431SBarry Smith .  ynorm - 2-norm of search length
347c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
348fee21e36SBarry Smith 
349c4a48953SLois Curfman McInnes    Options Database Key:
350c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
351c4a48953SLois Curfman McInnes 
3525e76c431SBarry Smith    Notes:
3535e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3545e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3555e76c431SBarry Smith 
35636851e7fSLois Curfman McInnes    Level: advanced
35736851e7fSLois Curfman McInnes 
35828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
35928ae5a21SLois Curfman McInnes 
360af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3615e76c431SBarry Smith @*/
362af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
363761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3645e76c431SBarry Smith {
365406556e6SLois Curfman McInnes   /*
366406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
367406556e6SLois Curfman McInnes      minimization problem:
368406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
369406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
370406556e6SLois Curfman McInnes    */
371406556e6SLois Curfman McInnes 
372ddd12767SBarry Smith   double     steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
373ea4d3ed3SLois Curfman McInnes   double     maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
374aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
3755e42470aSBarry Smith   Scalar     cinitslope, clambda;
3766b5873e3SBarry Smith #endif
3775e42470aSBarry Smith   int        ierr, count;
3785e42470aSBarry Smith   SNES_LS    *neP = (SNES_LS *) snes->data;
3795334005bSBarry Smith   Scalar     mone = -1.0, scale;
3808f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
3815e76c431SBarry Smith 
3823a40ed3dSBarry Smith   PetscFunctionBegin;
3837857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
384761aaf1bSLois Curfman McInnes   *flag   = 0;
3855e76c431SBarry Smith   alpha   = neP->alpha;
3865e76c431SBarry Smith   maxstep = neP->maxstep;
3875e76c431SBarry Smith   steptol = neP->steptol;
3885e76c431SBarry Smith 
389cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
390a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
391a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
392a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
393a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
394a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
395ad922ac9SBarry Smith     goto theend1;
39694a9d846SBarry Smith   }
3975e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3985e42470aSBarry Smith     scale = maxstep/(*ynorm);
399aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
40092318cfeSSatish Balay     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
4016b5873e3SBarry Smith #else
40294a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
4036b5873e3SBarry Smith #endif
404393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
4055e76c431SBarry Smith     *ynorm = maxstep;
4065e76c431SBarry Smith   }
4075e76c431SBarry Smith   minlambda = steptol/(*ynorm);
408a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
409aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
410a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
41192318cfeSSatish Balay   initslope = PetscReal(cinitslope);
4125e42470aSBarry Smith #else
413a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
4145e42470aSBarry Smith #endif
4155e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
4165e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
4175e76c431SBarry Smith 
418393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
4195334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
42078b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
421313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
422406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
423393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
42494a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
425ad922ac9SBarry Smith     goto theend1;
4265e76c431SBarry Smith   }
4275e76c431SBarry Smith 
4285e76c431SBarry Smith   /* Fit points with quadratic */
429313b5538SBarry Smith   lambda = 1.0;
430a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4315e76c431SBarry Smith   lambdaprev = lambda;
4325e76c431SBarry Smith   gnormprev = *gnorm;
433ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
434ddd12767SBarry Smith   else lambda = lambdatemp;
435393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
436ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
437aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
438ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4395e42470aSBarry Smith #else
440ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4415e42470aSBarry Smith #endif
44278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
443cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
444406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
445393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
446ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
447ad922ac9SBarry Smith     goto theend1;
4485e76c431SBarry Smith   }
4495e76c431SBarry Smith 
4505e76c431SBarry Smith   /* Fit points with cubic */
4515e76c431SBarry Smith   count = 1;
4525e76c431SBarry Smith   while (1) {
4535e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4542b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4552b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
456ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
457393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
458761aaf1bSLois Curfman McInnes       *flag = -1; break;
4595e76c431SBarry Smith     }
460406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
461406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
462ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4632b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4645e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4655e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4665e76c431SBarry Smith     if (a == 0.0) {
4675e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4685e76c431SBarry Smith     } else {
4695e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4705e76c431SBarry Smith     }
4715e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4725e76c431SBarry Smith       lambdatemp = .5*lambda;
4735e76c431SBarry Smith     }
4745e76c431SBarry Smith     lambdaprev = lambda;
4755e76c431SBarry Smith     gnormprev = *gnorm;
4765e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4775e76c431SBarry Smith       lambda = .1*lambda;
4785e76c431SBarry Smith     }
4795e76c431SBarry Smith     else lambda = lambdatemp;
480393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w );CHKERRQ(ierr);
481ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
483ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
484393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4855e42470aSBarry Smith #else
486ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4875e42470aSBarry Smith #endif
48878b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
489cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
490406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
491393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
492ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
4932715565aSLois Curfman McInnes       goto theend1;
4942b022350SLois Curfman McInnes     } else {
4952b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
4965e76c431SBarry Smith     }
4975e76c431SBarry Smith     count++;
4985e76c431SBarry Smith   }
499ad922ac9SBarry Smith   theend1:
5008f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5018f99978dSLois Curfman McInnes   if (neP->CheckStep) {
5028f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
5038f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
5048f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
5058f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
5068f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
5078f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
5088f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
5098f99978dSLois Curfman McInnes     }
5108f99978dSLois Curfman McInnes   }
5117857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5123a40ed3dSBarry Smith   PetscFunctionReturn(0);
5135e76c431SBarry Smith }
51404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
5155615d1e5SSatish Balay #undef __FUNC__
5165615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
5174b828684SBarry Smith /*@C
518f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
5195e76c431SBarry Smith 
520c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
521c7afd0dbSLois Curfman McInnes 
5225e42470aSBarry Smith    Input Parameters:
523c7afd0dbSLois Curfman McInnes +  snes - the SNES context
524af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
5255e42470aSBarry Smith .  x - current iterate
5265e42470aSBarry Smith .  f - residual evaluated at x
5275e42470aSBarry Smith .  y - search direction (contains new iterate on output)
5285e42470aSBarry Smith .  w - work vector
529c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
5305e42470aSBarry Smith 
531c4a48953SLois Curfman McInnes    Output Parameters:
532c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
5335e42470aSBarry Smith .  y - new iterate (contains search direction on input)
5345e42470aSBarry Smith .  gnorm - 2-norm of g
5355e42470aSBarry Smith .  ynorm - 2-norm of search length
536c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
537fee21e36SBarry Smith 
538c4a48953SLois Curfman McInnes    Options Database Key:
539c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
540c4a48953SLois Curfman McInnes 
5415e42470aSBarry Smith    Notes:
542c7afd0dbSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
5435e42470aSBarry Smith 
54436851e7fSLois Curfman McInnes    Level: advanced
54536851e7fSLois Curfman McInnes 
546f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
547f59ffdedSLois Curfman McInnes 
548af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
5495e42470aSBarry Smith @*/
550af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
551761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5525e76c431SBarry Smith {
553406556e6SLois Curfman McInnes   /*
554406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
555406556e6SLois Curfman McInnes      minimization problem:
556406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
557406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
558406556e6SLois Curfman McInnes    */
55940011551SBarry Smith   double     steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
560aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5615e42470aSBarry Smith   Scalar     cinitslope,clambda;
5626b5873e3SBarry Smith #endif
5635e42470aSBarry Smith   int        ierr, count;
5645e42470aSBarry Smith   SNES_LS    *neP = (SNES_LS *) snes->data;
5655334005bSBarry Smith   Scalar     mone = -1.0,scale;
5668f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
5675e76c431SBarry Smith 
5683a40ed3dSBarry Smith   PetscFunctionBegin;
5697857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
570761aaf1bSLois Curfman McInnes   *flag   = 0;
5715e42470aSBarry Smith   alpha   = neP->alpha;
5725e42470aSBarry Smith   maxstep = neP->maxstep;
5735e42470aSBarry Smith   steptol = neP->steptol;
5745e76c431SBarry Smith 
575cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
576b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
57794a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
578b37302c6SLois Curfman McInnes     *gnorm = fnorm;
579b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
580b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
581ad922ac9SBarry Smith     goto theend2;
58294a9d846SBarry Smith   }
5835e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5845e42470aSBarry Smith     scale = maxstep/(*ynorm);
585393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5865e42470aSBarry Smith     *ynorm = maxstep;
5875e76c431SBarry Smith   }
5885e42470aSBarry Smith   minlambda = steptol/(*ynorm);
589a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
590aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
591a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
59292318cfeSSatish Balay   initslope = PetscReal(cinitslope);
5935e42470aSBarry Smith #else
594a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5955e42470aSBarry Smith #endif
5965e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
5975e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
5985e42470aSBarry Smith 
599393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6005334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
60178b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
602cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
603406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
604393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
60594a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
606ad922ac9SBarry Smith     goto theend2;
6075e42470aSBarry Smith   }
6085e42470aSBarry Smith 
6095e42470aSBarry Smith   /* Fit points with quadratic */
610313b5538SBarry Smith   lambda = 1.0;
6115e42470aSBarry Smith   count = 1;
6125e42470aSBarry Smith   while (1) {
6135e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
614981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
615981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
616ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
617393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
618761aaf1bSLois Curfman McInnes       *flag = -1; break;
6195e42470aSBarry Smith     }
620a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
6215e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
6225e42470aSBarry Smith       lambda = .1*lambda;
6235e42470aSBarry Smith     } else lambda = lambdatemp;
624393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
6255334005bSBarry Smith     lambda = -lambda;
626aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
627393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6285e42470aSBarry Smith #else
629393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w);CHKERRQ(ierr);
6305e42470aSBarry Smith #endif
63178b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
632cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
633406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
634393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
635981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
63606259719SBarry Smith       break;
6375e42470aSBarry Smith     }
6385e42470aSBarry Smith     count++;
6395e42470aSBarry Smith   }
640ad922ac9SBarry Smith   theend2:
6418f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6428f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6438f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6448f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6458f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6468f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6478f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6488f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6498f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6508f99978dSLois Curfman McInnes     }
6518f99978dSLois Curfman McInnes   }
6527857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
6533a40ed3dSBarry Smith   PetscFunctionReturn(0);
6545e76c431SBarry Smith }
65504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6565615d1e5SSatish Balay #undef __FUNC__
657d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
658c9e6a524SLois Curfman McInnes /*@C
65977c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
660f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
6615e76c431SBarry Smith 
6625e76c431SBarry Smith    Input Parameters:
663c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
664af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
665c7afd0dbSLois Curfman McInnes -  func - pointer to int function
6665e76c431SBarry Smith 
667fee21e36SBarry Smith    Collective on SNES
668fee21e36SBarry Smith 
669c4a48953SLois Curfman McInnes    Available Routines:
670c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
671f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
672f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
673af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
6745e76c431SBarry Smith 
675c4a48953SLois Curfman McInnes     Options Database Keys:
676af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
677c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
678c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
679c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
680c4a48953SLois Curfman McInnes 
6815e76c431SBarry Smith    Calling sequence of func:
682bd208895SLois Curfman McInnes .vb
683af2d14d2SLois Curfman McInnes    func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
684bd208895SLois Curfman McInnes          double fnorm, double *ynorm, double *gnorm, *flag)
685bd208895SLois Curfman McInnes .ve
6865e76c431SBarry Smith 
6875e76c431SBarry Smith     Input parameters for func:
688c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
689af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
6905e76c431SBarry Smith .   x - current iterate
6915e76c431SBarry Smith .   f - residual evaluated at x
6925e76c431SBarry Smith .   y - search direction (contains new iterate on output)
6935e76c431SBarry Smith .   w - work vector
694c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
6955e76c431SBarry Smith 
6965e76c431SBarry Smith     Output parameters for func:
697c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
6985e76c431SBarry Smith .   y - new iterate (contains search direction on input)
6995e76c431SBarry Smith .   gnorm - 2-norm of g
7005e76c431SBarry Smith .   ynorm - 2-norm of search length
701c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
702761aaf1bSLois Curfman McInnes            on failure.
703f59ffdedSLois Curfman McInnes 
70436851e7fSLois Curfman McInnes     Level: advanced
70536851e7fSLois Curfman McInnes 
706f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
707f59ffdedSLois Curfman McInnes 
7084ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
7095e76c431SBarry Smith @*/
710af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx)
7115e76c431SBarry Smith {
712af2d14d2SLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*);
71382bf6240SBarry Smith 
7143a40ed3dSBarry Smith   PetscFunctionBegin;
71594d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
71682bf6240SBarry Smith   if (f) {
717af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
71882bf6240SBarry Smith   }
7193a40ed3dSBarry Smith   PetscFunctionReturn(0);
7205e76c431SBarry Smith }
72104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
722fb2e594dSBarry Smith EXTERN_C_BEGIN
72382bf6240SBarry Smith #undef __FUNC__
72482bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
725af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
726af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
72782bf6240SBarry Smith {
72882bf6240SBarry Smith   PetscFunctionBegin;
72982bf6240SBarry Smith   ((SNES_LS *)(snes->data))->LineSearch = func;
730af2d14d2SLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP = lsctx;
73182bf6240SBarry Smith   PetscFunctionReturn(0);
73282bf6240SBarry Smith }
733fb2e594dSBarry Smith EXTERN_C_END
73404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
735c8dd0c0dSLois Curfman McInnes #undef __FUNC__
736c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck"
737c8dd0c0dSLois Curfman McInnes /*@C
738530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
7398f99978dSLois Curfman McInnes    by the line search routine in the Newton-based method SNES_EQ_LS.
740c8dd0c0dSLois Curfman McInnes 
741c8dd0c0dSLois Curfman McInnes    Input Parameters:
742c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
743c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
744c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
745c8dd0c0dSLois Curfman McInnes 
746c8dd0c0dSLois Curfman McInnes    Collective on SNES
747c8dd0c0dSLois Curfman McInnes 
748c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
749c8dd0c0dSLois Curfman McInnes .vb
750b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
751c8dd0c0dSLois Curfman McInnes .ve
752b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
753b1ae27eaSLois Curfman McInnes    on failure.
754c8dd0c0dSLois Curfman McInnes 
755c8dd0c0dSLois Curfman McInnes    Input parameters for func:
756c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
757c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
7588f99978dSLois Curfman McInnes -  x - current candidate iterate
759c8dd0c0dSLois Curfman McInnes 
760c8dd0c0dSLois Curfman McInnes    Output parameters for func:
761c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
762c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
763c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
764c8dd0c0dSLois Curfman McInnes 
765c8dd0c0dSLois Curfman McInnes    Level: advanced
766c8dd0c0dSLois Curfman McInnes 
7678f99978dSLois Curfman McInnes    Notes:
768b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
769b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
770b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
771b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
772ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
7738f99978dSLois Curfman McInnes 
774b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
775b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
776b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
777ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
778ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
779ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
780ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
781b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
7828f99978dSLois Curfman McInnes 
783c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
784c8dd0c0dSLois Curfman McInnes 
785c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
786c8dd0c0dSLois Curfman McInnes @*/
7878f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
788c8dd0c0dSLois Curfman McInnes {
7898f99978dSLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
790c8dd0c0dSLois Curfman McInnes 
791c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
792c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
793c8dd0c0dSLois Curfman McInnes   if (f) {
794c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
795c8dd0c0dSLois Curfman McInnes   }
796c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
797c8dd0c0dSLois Curfman McInnes }
798c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
799c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
800c8dd0c0dSLois Curfman McInnes #undef __FUNC__
801c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS"
8028f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
803c8dd0c0dSLois Curfman McInnes {
804c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
805c8dd0c0dSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
806c8dd0c0dSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP = checkctx;
807c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
808c8dd0c0dSLois Curfman McInnes }
809c8dd0c0dSLois Curfman McInnes EXTERN_C_END
810c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
81104d965bbSLois Curfman McInnes /*
81204d965bbSLois Curfman McInnes    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
81382bf6240SBarry Smith 
81404d965bbSLois Curfman McInnes    Input Parameter:
81504d965bbSLois Curfman McInnes .  snes - the SNES context
81604d965bbSLois Curfman McInnes 
81704d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
81804d965bbSLois Curfman McInnes */
8195615d1e5SSatish Balay #undef __FUNC__
820d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
821f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
8225e42470aSBarry Smith {
8235e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
824d132466eSBarry Smith   int     ierr;
8256daaf66cSBarry Smith 
8263a40ed3dSBarry Smith   PetscFunctionBegin;
827d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr);
828d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr);
829d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr);
830d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr);
831d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr);
8323a40ed3dSBarry Smith   PetscFunctionReturn(0);
8335e42470aSBarry Smith }
83404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
83504d965bbSLois Curfman McInnes /*
836de49b36dSLois Curfman McInnes    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
83704d965bbSLois Curfman McInnes 
83804d965bbSLois Curfman McInnes    Input Parameters:
83904d965bbSLois Curfman McInnes .  SNES - the SNES context
84004d965bbSLois Curfman McInnes .  viewer - visualization context
84104d965bbSLois Curfman McInnes 
84204d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
84304d965bbSLois Curfman McInnes */
8445615d1e5SSatish Balay #undef __FUNC__
845d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
846e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
847a935fc98SLois Curfman McInnes {
848a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
84919bcc07fSBarry Smith   char       *cstr;
85051695f54SLois Curfman McInnes   int        ierr;
85119bcc07fSBarry Smith   ViewerType vtype;
852a935fc98SLois Curfman McInnes 
8533a40ed3dSBarry Smith   PetscFunctionBegin;
85419bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr);
8553f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
85619bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
85719bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
85819bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
85919bcc07fSBarry Smith     else                                                cstr = "unknown";
860d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
861d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
8625cd90555SBarry Smith   } else {
8635cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
86419bcc07fSBarry Smith   }
8653a40ed3dSBarry Smith   PetscFunctionReturn(0);
866a935fc98SLois Curfman McInnes }
86704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
86804d965bbSLois Curfman McInnes /*
86904d965bbSLois Curfman McInnes    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
87004d965bbSLois Curfman McInnes 
87104d965bbSLois Curfman McInnes    Input Parameter:
87204d965bbSLois Curfman McInnes .  snes - the SNES context
87304d965bbSLois Curfman McInnes 
87404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
87504d965bbSLois Curfman McInnes */
8765615d1e5SSatish Balay #undef __FUNC__
8775615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
878f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
8795e42470aSBarry Smith {
8805e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
8815e42470aSBarry Smith   char    ver[16];
8825e42470aSBarry Smith   double  tmp;
88325cdf11fSBarry Smith   int     ierr,flg;
8845e42470aSBarry Smith 
8853a40ed3dSBarry Smith   PetscFunctionBegin;
88609d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
88725cdf11fSBarry Smith   if (flg) {
8885e42470aSBarry Smith     ls->alpha = tmp;
8895e42470aSBarry Smith   }
89009d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
89125cdf11fSBarry Smith   if (flg) {
8925e42470aSBarry Smith     ls->maxstep = tmp;
8935e42470aSBarry Smith   }
89409d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
89525cdf11fSBarry Smith   if (flg) {
8965e42470aSBarry Smith     ls->steptol = tmp;
8975e42470aSBarry Smith   }
89809d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg);CHKERRQ(ierr);
89925cdf11fSBarry Smith   if (flg) {
90048d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
901af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
902b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"basicnonorms")) {
903af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
904b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"quadratic")) {
905af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
906b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"cubic")) {
907af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
9085e42470aSBarry Smith     }
909a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
9105e42470aSBarry Smith   }
9113a40ed3dSBarry Smith   PetscFunctionReturn(0);
9125e42470aSBarry Smith }
91304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
91404d965bbSLois Curfman McInnes /*
91504d965bbSLois Curfman McInnes    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
91604d965bbSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
91704d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
91804d965bbSLois Curfman McInnes 
91904d965bbSLois Curfman McInnes 
92004d965bbSLois Curfman McInnes    Input Parameter:
92104d965bbSLois Curfman McInnes .  snes - the SNES context
92204d965bbSLois Curfman McInnes 
92304d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
92404d965bbSLois Curfman McInnes  */
925fb2e594dSBarry Smith EXTERN_C_BEGIN
9265615d1e5SSatish Balay #undef __FUNC__
9275615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
928f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
9295e42470aSBarry Smith {
93082bf6240SBarry Smith   int     ierr;
9315e42470aSBarry Smith   SNES_LS *neP;
9325e42470aSBarry Smith 
9333a40ed3dSBarry Smith   PetscFunctionBegin;
934a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
935a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
936a8c6a408SBarry Smith   }
93782bf6240SBarry Smith 
938f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
939f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
940f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
941f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
942f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
943f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
944f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
9455baf8537SBarry Smith   snes->nwork           = 0;
9465e42470aSBarry Smith 
9470452661fSBarry Smith   neP			= PetscNew(SNES_LS);CHKPTRQ(neP);
948ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
9495e42470aSBarry Smith   snes->data    	= (void *) neP;
9505e42470aSBarry Smith   neP->alpha		= 1.e-4;
9515e42470aSBarry Smith   neP->maxstep		= 1.e8;
9525e42470aSBarry Smith   neP->steptol		= 1.e-12;
9535e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
954c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
955c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
956c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
95782bf6240SBarry Smith 
95894d884c6SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
959e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
960c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
961c8dd0c0dSLois Curfman McInnes                     (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
96282bf6240SBarry Smith 
9633a40ed3dSBarry Smith   PetscFunctionReturn(0);
9645e42470aSBarry Smith }
965fb2e594dSBarry Smith EXTERN_C_END
96682bf6240SBarry Smith 
96782bf6240SBarry Smith 
96882bf6240SBarry Smith 
969