xref: /petsc/src/snes/impls/ls/ls.c (revision 6831982abb6453c2f3e25efecb5d0951d333371e)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*6831982aSBarry Smith static char vcid[] = "$Id: ls.c,v 1.144 1999/10/13 20:38:29 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 {
68*6831982aSBarry Smith   SNES_EQ_LS          *neP = (SNES_EQ_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;
73184914b5SBarry Smith   SNESConvergedReason reason;
745e76c431SBarry Smith 
753a40ed3dSBarry Smith   PetscFunctionBegin;
76184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
77184914b5SBarry Smith 
785e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
795e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
8039e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
815e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
825e42470aSBarry Smith   G		= snes->work[1];
835e42470aSBarry Smith   W		= snes->work[2];
845e76c431SBarry Smith 
855334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
86cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
870f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
8884c569c9SBarry Smith   snes->iter = 0;
895e42470aSBarry Smith   snes->norm = fnorm;
900f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
9184c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
9294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
935e76c431SBarry Smith 
94184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
9594a9d846SBarry Smith 
96d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
97d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
98d034289bSBarry Smith 
995e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
1005e76c431SBarry Smith 
101ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1025334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1035334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
10478b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
1057a00f4a9SLois Curfman McInnes     snes->linear_its += PetscAbsInt(lits);
106af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
107ea4d3ed3SLois Curfman McInnes 
108ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
109ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
110ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
111ea4d3ed3SLois Curfman McInnes     */
11281b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
113af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
114af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
115761aaf1bSLois Curfman McInnes     if (lsfail) snes->nfailures++;
1165e76c431SBarry Smith 
11739e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
11839e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1195e76c431SBarry Smith     fnorm = gnorm;
1205e76c431SBarry Smith 
1210f5bd95cSBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
12284c569c9SBarry Smith     snes->iter = i+1;
1235e42470aSBarry Smith     snes->norm = fnorm;
1240f5bd95cSBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
12584c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
12694a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1275e76c431SBarry Smith 
1285e76c431SBarry Smith     /* Test for convergence */
12929e0b56fSBarry Smith     if (snes->converged) {
13029e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
131184914b5SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&reason,snes->cnvP);CHKERRQ(ierr);
132184914b5SBarry Smith       if (reason) {
13316c95cacSBarry Smith         break;
13416c95cacSBarry Smith       }
13516c95cacSBarry Smith     }
13629e0b56fSBarry Smith   }
13739e2f89bSBarry Smith   if (X != snes->vec_sol) {
138393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
13939e2f89bSBarry Smith     snes->vec_sol_always  = snes->vec_sol;
14039e2f89bSBarry Smith     snes->vec_func_always = snes->vec_func;
14139e2f89bSBarry Smith   }
14252392280SLois Curfman McInnes   if (i == maxits) {
143981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
14452392280SLois Curfman McInnes     i--;
145184914b5SBarry Smith     reason = SNES_DIVERGED_MAX_IT;
14652392280SLois Curfman McInnes   }
1470f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
148184914b5SBarry Smith   snes->reason = reason;
1490f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1505e42470aSBarry Smith   *outits = i+1;
1513a40ed3dSBarry Smith   PetscFunctionReturn(0);
1525e76c431SBarry Smith }
15304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
15404d965bbSLois Curfman McInnes /*
15504d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
156*6831982aSBarry Smith    of the SNESEQLS nonlinear solver.
15704d965bbSLois Curfman McInnes 
15804d965bbSLois Curfman McInnes    Input Parameter:
15904d965bbSLois Curfman McInnes .  snes - the SNES context
16004d965bbSLois Curfman McInnes .  x - the solution vector
16104d965bbSLois Curfman McInnes 
16204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
16304d965bbSLois Curfman McInnes 
16404d965bbSLois Curfman McInnes    Notes:
16504d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
16604d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
16704d965bbSLois Curfman McInnes    the call to SNESSolve().
16804d965bbSLois Curfman McInnes  */
1695615d1e5SSatish Balay #undef __FUNC__
1705615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
171f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
1725e76c431SBarry Smith {
1735e42470aSBarry Smith   int ierr;
1743a40ed3dSBarry Smith 
1753a40ed3dSBarry Smith   PetscFunctionBegin;
17681b6cf68SLois Curfman McInnes   snes->nwork = 4;
177d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1785e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
17981b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1803a40ed3dSBarry Smith   PetscFunctionReturn(0);
1815e76c431SBarry Smith }
18204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
18304d965bbSLois Curfman McInnes /*
184*6831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
18504d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
18604d965bbSLois Curfman McInnes 
18704d965bbSLois Curfman McInnes    Input Parameter:
18804d965bbSLois Curfman McInnes .  snes - the SNES context
18904d965bbSLois Curfman McInnes 
190de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
19104d965bbSLois Curfman McInnes  */
1925615d1e5SSatish Balay #undef __FUNC__
193d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
194e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
1955e76c431SBarry Smith {
196393d2d9aSLois Curfman McInnes   int  ierr;
1973a40ed3dSBarry Smith 
1983a40ed3dSBarry Smith   PetscFunctionBegin;
1995baf8537SBarry Smith   if (snes->nwork) {
2004b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
20121c89e3eSBarry Smith   }
2025c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2033a40ed3dSBarry Smith   PetscFunctionReturn(0);
2045e76c431SBarry Smith }
20504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2065615d1e5SSatish Balay #undef __FUNC__
2075615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
20882bf6240SBarry Smith 
2094b828684SBarry Smith /*@C
2105e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2115e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2125e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2135e76c431SBarry Smith 
214c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
215c7afd0dbSLois Curfman McInnes 
2165e76c431SBarry Smith    Input Parameters:
217c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
218af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
2195e76c431SBarry Smith .  x - current iterate
2205e76c431SBarry Smith .  f - residual evaluated at x
2215e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2225e76c431SBarry Smith .  w - work vector
223c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
2245e76c431SBarry Smith 
225c4a48953SLois Curfman McInnes    Output Parameters:
226c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
2275e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2285e76c431SBarry Smith .  gnorm - 2-norm of g
2295e76c431SBarry Smith .  ynorm - 2-norm of search length
230c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
231fee21e36SBarry Smith 
232c4a48953SLois Curfman McInnes    Options Database Key:
233c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
234c4a48953SLois Curfman McInnes 
23536851e7fSLois Curfman McInnes    Level: advanced
23636851e7fSLois Curfman McInnes 
23728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
23828ae5a21SLois Curfman McInnes 
239f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
240af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
2415e76c431SBarry Smith @*/
242af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
243761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
2445e76c431SBarry Smith {
2455e42470aSBarry Smith   int        ierr;
2465334005bSBarry Smith   Scalar     mone = -1.0;
247*6831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data;
2488f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
2495334005bSBarry Smith 
2503a40ed3dSBarry Smith   PetscFunctionBegin;
251761aaf1bSLois Curfman McInnes   *flag = 0;
2527857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
253a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
254ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
2558f99978dSLois Curfman McInnes   if (neP->CheckStep) {
2568f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
2578f99978dSLois Curfman McInnes   }
258ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
259a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
2607857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
2625e76c431SBarry Smith }
26304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2645615d1e5SSatish Balay #undef __FUNC__
26529e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
26682bf6240SBarry Smith 
26729e0b56fSBarry Smith /*@C
26829e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
26929e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
27029e0b56fSBarry Smith    even compute the norm of the function or search direction; this
27129e0b56fSBarry Smith    is intended only when you know the full step is fine and are
27229e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
273c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
274c7afd0dbSLois Curfman McInnes 
275c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
27629e0b56fSBarry Smith 
27729e0b56fSBarry Smith    Input Parameters:
278c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
279af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
28029e0b56fSBarry Smith .  x - current iterate
28129e0b56fSBarry Smith .  f - residual evaluated at x
28229e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
28329e0b56fSBarry Smith .  w - work vector
284c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
28529e0b56fSBarry Smith 
28629e0b56fSBarry Smith    Output Parameters:
287c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
28829e0b56fSBarry Smith .  gnorm - not changed
28929e0b56fSBarry Smith .  ynorm - not changed
290c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
291fee21e36SBarry Smith 
29229e0b56fSBarry Smith    Options Database Key:
293c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
29429e0b56fSBarry Smith 
2958cbba510SBarry Smith    Notes:
296ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
297ea56f5baSLois Curfman McInnes    either the options
298ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
299ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
300ea56f5baSLois Curfman McInnes    SNESSetConvergenceTest(); otherwise, the SNES solver will
301ea56f5baSLois Curfman McInnes    generate an error.
3028cbba510SBarry Smith 
303ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
304ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
305ea56f5baSLois Curfman McInnes    correct, since they are not computed.
306ea56f5baSLois Curfman McInnes 
307ea56f5baSLois Curfman McInnes    Level: advanced
3088cbba510SBarry Smith 
30929e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
31029e0b56fSBarry Smith 
31129e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
31229e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
31329e0b56fSBarry Smith @*/
314af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
31529e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
31629e0b56fSBarry Smith {
31729e0b56fSBarry Smith   int        ierr;
31829e0b56fSBarry Smith   Scalar     mone = -1.0;
319*6831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data;
3208f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
32129e0b56fSBarry Smith 
3223a40ed3dSBarry Smith   PetscFunctionBegin;
3238cbba510SBarry Smith   *flag = 0;
32429e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
32529e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3268f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3278f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3288f99978dSLois Curfman McInnes   }
32929e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
33029e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3313a40ed3dSBarry Smith   PetscFunctionReturn(0);
33229e0b56fSBarry Smith }
33304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
33429e0b56fSBarry Smith #undef __FUNC__
3355615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3364b828684SBarry Smith /*@C
337f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3385e76c431SBarry Smith 
339c7afd0dbSLois Curfman McInnes    Collective on SNES
340c7afd0dbSLois Curfman McInnes 
3415e76c431SBarry Smith    Input Parameters:
342c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
343af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3445e76c431SBarry Smith .  x - current iterate
3455e76c431SBarry Smith .  f - residual evaluated at x
3465e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3475e76c431SBarry Smith .  w - work vector
348c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3495e76c431SBarry Smith 
350393d2d9aSLois Curfman McInnes    Output Parameters:
351c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3525e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3535e76c431SBarry Smith .  gnorm - 2-norm of g
3545e76c431SBarry Smith .  ynorm - 2-norm of search length
355c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
356fee21e36SBarry Smith 
357c4a48953SLois Curfman McInnes    Options Database Key:
358c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
359c4a48953SLois Curfman McInnes 
3605e76c431SBarry Smith    Notes:
3615e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3625e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3635e76c431SBarry Smith 
36436851e7fSLois Curfman McInnes    Level: advanced
36536851e7fSLois Curfman McInnes 
36628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
36728ae5a21SLois Curfman McInnes 
368af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3695e76c431SBarry Smith @*/
370af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
371761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3725e76c431SBarry Smith {
373406556e6SLois Curfman McInnes   /*
374406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
375406556e6SLois Curfman McInnes      minimization problem:
376406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
377406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
378406556e6SLois Curfman McInnes    */
379406556e6SLois Curfman McInnes 
380ddd12767SBarry Smith   double     steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
381ea4d3ed3SLois Curfman McInnes   double     maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
382aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
3835e42470aSBarry Smith   Scalar     cinitslope, clambda;
3846b5873e3SBarry Smith #endif
3855e42470aSBarry Smith   int        ierr, count;
386*6831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data;
3875334005bSBarry Smith   Scalar     mone = -1.0, scale;
3888f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
3895e76c431SBarry Smith 
3903a40ed3dSBarry Smith   PetscFunctionBegin;
3917857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
392761aaf1bSLois Curfman McInnes   *flag   = 0;
3935e76c431SBarry Smith   alpha   = neP->alpha;
3945e76c431SBarry Smith   maxstep = neP->maxstep;
3955e76c431SBarry Smith   steptol = neP->steptol;
3965e76c431SBarry Smith 
397cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
398a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
399a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
400a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
401a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
402a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
403ad922ac9SBarry Smith     goto theend1;
40494a9d846SBarry Smith   }
4055e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4065e42470aSBarry Smith     scale = maxstep/(*ynorm);
407aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
40892318cfeSSatish Balay     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
4096b5873e3SBarry Smith #else
41094a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
4116b5873e3SBarry Smith #endif
412393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
4135e76c431SBarry Smith     *ynorm = maxstep;
4145e76c431SBarry Smith   }
4155e76c431SBarry Smith   minlambda = steptol/(*ynorm);
416a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
417aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
418a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
41992318cfeSSatish Balay   initslope = PetscReal(cinitslope);
4205e42470aSBarry Smith #else
421a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
4225e42470aSBarry Smith #endif
4235e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
4245e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
4255e76c431SBarry Smith 
426393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
4275334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
42878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
429313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
430406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
431393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
43294a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
433ad922ac9SBarry Smith     goto theend1;
4345e76c431SBarry Smith   }
4355e76c431SBarry Smith 
4365e76c431SBarry Smith   /* Fit points with quadratic */
437313b5538SBarry Smith   lambda = 1.0;
438a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4395e76c431SBarry Smith   lambdaprev = lambda;
4405e76c431SBarry Smith   gnormprev = *gnorm;
441ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
442ddd12767SBarry Smith   else lambda = lambdatemp;
443393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
444ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
445aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
446ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4475e42470aSBarry Smith #else
448ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4495e42470aSBarry Smith #endif
45078b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
451cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
452406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
453393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
454ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
455ad922ac9SBarry Smith     goto theend1;
4565e76c431SBarry Smith   }
4575e76c431SBarry Smith 
4585e76c431SBarry Smith   /* Fit points with cubic */
4595e76c431SBarry Smith   count = 1;
4605e76c431SBarry Smith   while (1) {
4615e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4622b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4632b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
464ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
465393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
466761aaf1bSLois Curfman McInnes       *flag = -1; break;
4675e76c431SBarry Smith     }
468406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
469406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
470ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4712b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4725e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4735e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4745e76c431SBarry Smith     if (a == 0.0) {
4755e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4765e76c431SBarry Smith     } else {
4775e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4785e76c431SBarry Smith     }
4795e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4805e76c431SBarry Smith       lambdatemp = .5*lambda;
4815e76c431SBarry Smith     }
4825e76c431SBarry Smith     lambdaprev = lambda;
4835e76c431SBarry Smith     gnormprev = *gnorm;
4845e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4855e76c431SBarry Smith       lambda = .1*lambda;
4865e76c431SBarry Smith     }
4875e76c431SBarry Smith     else lambda = lambdatemp;
488393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w );CHKERRQ(ierr);
489ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
490aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
491ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
492393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4935e42470aSBarry Smith #else
494ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4955e42470aSBarry Smith #endif
49678b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
497cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
498406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
499393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
500ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
5012715565aSLois Curfman McInnes       goto theend1;
5022b022350SLois Curfman McInnes     } else {
5032b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
5045e76c431SBarry Smith     }
5055e76c431SBarry Smith     count++;
5065e76c431SBarry Smith   }
507ad922ac9SBarry Smith   theend1:
5088f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5098f99978dSLois Curfman McInnes   if (neP->CheckStep) {
5108f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
5118f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
5128f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
5138f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
5148f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
5158f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
5168f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
5178f99978dSLois Curfman McInnes     }
5188f99978dSLois Curfman McInnes   }
5197857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5203a40ed3dSBarry Smith   PetscFunctionReturn(0);
5215e76c431SBarry Smith }
52204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
5235615d1e5SSatish Balay #undef __FUNC__
5245615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
5254b828684SBarry Smith /*@C
526f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
5275e76c431SBarry Smith 
528c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
529c7afd0dbSLois Curfman McInnes 
5305e42470aSBarry Smith    Input Parameters:
531c7afd0dbSLois Curfman McInnes +  snes - the SNES context
532af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
5335e42470aSBarry Smith .  x - current iterate
5345e42470aSBarry Smith .  f - residual evaluated at x
5355e42470aSBarry Smith .  y - search direction (contains new iterate on output)
5365e42470aSBarry Smith .  w - work vector
537c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
5385e42470aSBarry Smith 
539c4a48953SLois Curfman McInnes    Output Parameters:
540c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
5415e42470aSBarry Smith .  y - new iterate (contains search direction on input)
5425e42470aSBarry Smith .  gnorm - 2-norm of g
5435e42470aSBarry Smith .  ynorm - 2-norm of search length
544c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
545fee21e36SBarry Smith 
546c4a48953SLois Curfman McInnes    Options Database Key:
547c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
548c4a48953SLois Curfman McInnes 
5495e42470aSBarry Smith    Notes:
550*6831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
5515e42470aSBarry Smith 
55236851e7fSLois Curfman McInnes    Level: advanced
55336851e7fSLois Curfman McInnes 
554f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
555f59ffdedSLois Curfman McInnes 
556af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
5575e42470aSBarry Smith @*/
558af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
559761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5605e76c431SBarry Smith {
561406556e6SLois Curfman McInnes   /*
562406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
563406556e6SLois Curfman McInnes      minimization problem:
564406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
565406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
566406556e6SLois Curfman McInnes    */
56740011551SBarry Smith   double     steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
568aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5695e42470aSBarry Smith   Scalar     cinitslope,clambda;
5706b5873e3SBarry Smith #endif
5715e42470aSBarry Smith   int        ierr, count;
572*6831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS *) snes->data;
5735334005bSBarry Smith   Scalar     mone = -1.0,scale;
5748f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
5755e76c431SBarry Smith 
5763a40ed3dSBarry Smith   PetscFunctionBegin;
5777857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
578761aaf1bSLois Curfman McInnes   *flag   = 0;
5795e42470aSBarry Smith   alpha   = neP->alpha;
5805e42470aSBarry Smith   maxstep = neP->maxstep;
5815e42470aSBarry Smith   steptol = neP->steptol;
5825e76c431SBarry Smith 
583cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
584b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
58594a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
586b37302c6SLois Curfman McInnes     *gnorm = fnorm;
587b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
588b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
589ad922ac9SBarry Smith     goto theend2;
59094a9d846SBarry Smith   }
5915e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5925e42470aSBarry Smith     scale = maxstep/(*ynorm);
593393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5945e42470aSBarry Smith     *ynorm = maxstep;
5955e76c431SBarry Smith   }
5965e42470aSBarry Smith   minlambda = steptol/(*ynorm);
597a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
598aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
599a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
60092318cfeSSatish Balay   initslope = PetscReal(cinitslope);
6015e42470aSBarry Smith #else
602a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6035e42470aSBarry Smith #endif
6045e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
6055e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6065e42470aSBarry Smith 
607393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6085334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
60978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
610cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
611406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
612393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
61394a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
614ad922ac9SBarry Smith     goto theend2;
6155e42470aSBarry Smith   }
6165e42470aSBarry Smith 
6175e42470aSBarry Smith   /* Fit points with quadratic */
618313b5538SBarry Smith   lambda = 1.0;
6195e42470aSBarry Smith   count = 1;
6205e42470aSBarry Smith   while (1) {
6215e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
622981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
623981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
624ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
625393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
626761aaf1bSLois Curfman McInnes       *flag = -1; break;
6275e42470aSBarry Smith     }
628a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
6295e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
6305e42470aSBarry Smith       lambda = .1*lambda;
6315e42470aSBarry Smith     } else lambda = lambdatemp;
632393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
6335334005bSBarry Smith     lambda = -lambda;
634aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
635393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6365e42470aSBarry Smith #else
637393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w);CHKERRQ(ierr);
6385e42470aSBarry Smith #endif
63978b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
640cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
641406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
642393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
643981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
64406259719SBarry Smith       break;
6455e42470aSBarry Smith     }
6465e42470aSBarry Smith     count++;
6475e42470aSBarry Smith   }
648ad922ac9SBarry Smith   theend2:
6498f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6508f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6518f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6528f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6538f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6548f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6558f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6568f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6578f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6588f99978dSLois Curfman McInnes     }
6598f99978dSLois Curfman McInnes   }
6607857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
6613a40ed3dSBarry Smith   PetscFunctionReturn(0);
6625e76c431SBarry Smith }
66304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6645615d1e5SSatish Balay #undef __FUNC__
665d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
666c9e6a524SLois Curfman McInnes /*@C
66777c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
668*6831982aSBarry Smith    by the method SNESEQLS.
6695e76c431SBarry Smith 
6705e76c431SBarry Smith    Input Parameters:
671c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
672af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
673c7afd0dbSLois Curfman McInnes -  func - pointer to int function
6745e76c431SBarry Smith 
675fee21e36SBarry Smith    Collective on SNES
676fee21e36SBarry Smith 
677c4a48953SLois Curfman McInnes    Available Routines:
678c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
679f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
680f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
681af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
6825e76c431SBarry Smith 
683c4a48953SLois Curfman McInnes     Options Database Keys:
684af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
685c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
686c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
687c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
688c4a48953SLois Curfman McInnes 
6895e76c431SBarry Smith    Calling sequence of func:
690bd208895SLois Curfman McInnes .vb
691af2d14d2SLois Curfman McInnes    func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
692bd208895SLois Curfman McInnes          double fnorm, double *ynorm, double *gnorm, *flag)
693bd208895SLois Curfman McInnes .ve
6945e76c431SBarry Smith 
6955e76c431SBarry Smith     Input parameters for func:
696c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
697af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
6985e76c431SBarry Smith .   x - current iterate
6995e76c431SBarry Smith .   f - residual evaluated at x
7005e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7015e76c431SBarry Smith .   w - work vector
702c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
7035e76c431SBarry Smith 
7045e76c431SBarry Smith     Output parameters for func:
705c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7065e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7075e76c431SBarry Smith .   gnorm - 2-norm of g
7085e76c431SBarry Smith .   ynorm - 2-norm of search length
709c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
710761aaf1bSLois Curfman McInnes            on failure.
711f59ffdedSLois Curfman McInnes 
71236851e7fSLois Curfman McInnes     Level: advanced
71336851e7fSLois Curfman McInnes 
714f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
715f59ffdedSLois Curfman McInnes 
7164ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
7175e76c431SBarry Smith @*/
718af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx)
7195e76c431SBarry Smith {
720af2d14d2SLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*);
72182bf6240SBarry Smith 
7223a40ed3dSBarry Smith   PetscFunctionBegin;
72394d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
72482bf6240SBarry Smith   if (f) {
725af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
72682bf6240SBarry Smith   }
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
7285e76c431SBarry Smith }
72904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
730fb2e594dSBarry Smith EXTERN_C_BEGIN
73182bf6240SBarry Smith #undef __FUNC__
73282bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
733af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
734af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
73582bf6240SBarry Smith {
73682bf6240SBarry Smith   PetscFunctionBegin;
737*6831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
738*6831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
73982bf6240SBarry Smith   PetscFunctionReturn(0);
74082bf6240SBarry Smith }
741fb2e594dSBarry Smith EXTERN_C_END
74204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
743c8dd0c0dSLois Curfman McInnes #undef __FUNC__
744c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck"
745c8dd0c0dSLois Curfman McInnes /*@C
746530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
747*6831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
748c8dd0c0dSLois Curfman McInnes 
749c8dd0c0dSLois Curfman McInnes    Input Parameters:
750c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
751c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
752c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
753c8dd0c0dSLois Curfman McInnes 
754c8dd0c0dSLois Curfman McInnes    Collective on SNES
755c8dd0c0dSLois Curfman McInnes 
756c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
757c8dd0c0dSLois Curfman McInnes .vb
758b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
759c8dd0c0dSLois Curfman McInnes .ve
760b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
761b1ae27eaSLois Curfman McInnes    on failure.
762c8dd0c0dSLois Curfman McInnes 
763c8dd0c0dSLois Curfman McInnes    Input parameters for func:
764c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
765c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
7668f99978dSLois Curfman McInnes -  x - current candidate iterate
767c8dd0c0dSLois Curfman McInnes 
768c8dd0c0dSLois Curfman McInnes    Output parameters for func:
769c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
770c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
771c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
772c8dd0c0dSLois Curfman McInnes 
773c8dd0c0dSLois Curfman McInnes    Level: advanced
774c8dd0c0dSLois Curfman McInnes 
7758f99978dSLois Curfman McInnes    Notes:
776b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
777b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
778b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
779b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
780ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
7818f99978dSLois Curfman McInnes 
782b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
783b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
784b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
785ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
786ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
787ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
788ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
789b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
7908f99978dSLois Curfman McInnes 
791c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
792c8dd0c0dSLois Curfman McInnes 
793c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
794c8dd0c0dSLois Curfman McInnes @*/
7958f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
796c8dd0c0dSLois Curfman McInnes {
7978f99978dSLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
798c8dd0c0dSLois Curfman McInnes 
799c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
800c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
801c8dd0c0dSLois Curfman McInnes   if (f) {
802c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
803c8dd0c0dSLois Curfman McInnes   }
804c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
805c8dd0c0dSLois Curfman McInnes }
806c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
807c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
808c8dd0c0dSLois Curfman McInnes #undef __FUNC__
809c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS"
8108f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
811c8dd0c0dSLois Curfman McInnes {
812c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
813*6831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
814*6831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
815c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
816c8dd0c0dSLois Curfman McInnes }
817c8dd0c0dSLois Curfman McInnes EXTERN_C_END
818c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
81904d965bbSLois Curfman McInnes /*
820*6831982aSBarry Smith    SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method.
82182bf6240SBarry Smith 
82204d965bbSLois Curfman McInnes    Input Parameter:
82304d965bbSLois Curfman McInnes .  snes - the SNES context
82404d965bbSLois Curfman McInnes 
82504d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
82604d965bbSLois Curfman McInnes */
8275615d1e5SSatish Balay #undef __FUNC__
828d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
829f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
8305e42470aSBarry Smith {
831*6831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
832d132466eSBarry Smith   int     ierr;
8336daaf66cSBarry Smith 
8343a40ed3dSBarry Smith   PetscFunctionBegin;
835*6831982aSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr);
836d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr);
837d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr);
838d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr);
839d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr);
8403a40ed3dSBarry Smith   PetscFunctionReturn(0);
8415e42470aSBarry Smith }
84204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
84304d965bbSLois Curfman McInnes /*
844*6831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
84504d965bbSLois Curfman McInnes 
84604d965bbSLois Curfman McInnes    Input Parameters:
84704d965bbSLois Curfman McInnes .  SNES - the SNES context
84804d965bbSLois Curfman McInnes .  viewer - visualization context
84904d965bbSLois Curfman McInnes 
85004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
85104d965bbSLois Curfman McInnes */
8525615d1e5SSatish Balay #undef __FUNC__
853d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
854e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
855a935fc98SLois Curfman McInnes {
856*6831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
85719bcc07fSBarry Smith   char       *cstr;
85851695f54SLois Curfman McInnes   int        ierr;
859*6831982aSBarry Smith   PetscTruth isascii;
860a935fc98SLois Curfman McInnes 
8613a40ed3dSBarry Smith   PetscFunctionBegin;
862*6831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
8630f5bd95cSBarry Smith   if (isascii) {
86419bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
86519bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
86619bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
86719bcc07fSBarry Smith     else                                                cstr = "unknown";
868d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
869d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
8705cd90555SBarry Smith   } else {
8710f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
87219bcc07fSBarry Smith   }
8733a40ed3dSBarry Smith   PetscFunctionReturn(0);
874a935fc98SLois Curfman McInnes }
87504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
87604d965bbSLois Curfman McInnes /*
877*6831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
87804d965bbSLois Curfman McInnes 
87904d965bbSLois Curfman McInnes    Input Parameter:
88004d965bbSLois Curfman McInnes .  snes - the SNES context
88104d965bbSLois Curfman McInnes 
88204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
88304d965bbSLois Curfman McInnes */
8845615d1e5SSatish Balay #undef __FUNC__
8855615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
886f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
8875e42470aSBarry Smith {
888*6831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
8895e42470aSBarry Smith   char       ver[16];
8905e42470aSBarry Smith   double     tmp;
89125cdf11fSBarry Smith   int        ierr,flg;
8925e42470aSBarry Smith 
8933a40ed3dSBarry Smith   PetscFunctionBegin;
89409d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
89525cdf11fSBarry Smith   if (flg) {
8965e42470aSBarry Smith     ls->alpha = tmp;
8975e42470aSBarry Smith   }
89809d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
89925cdf11fSBarry Smith   if (flg) {
9005e42470aSBarry Smith     ls->maxstep = tmp;
9015e42470aSBarry Smith   }
90209d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
90325cdf11fSBarry Smith   if (flg) {
9045e42470aSBarry Smith     ls->steptol = tmp;
9055e42470aSBarry Smith   }
90609d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg);CHKERRQ(ierr);
90725cdf11fSBarry Smith   if (flg) {
9080f5bd95cSBarry Smith     int isbasic,isnonorms,isquad,iscubic;
9090f5bd95cSBarry Smith 
9100f5bd95cSBarry Smith     isbasic   = !PetscStrcmp(ver,"basic");
9110f5bd95cSBarry Smith     isnonorms = !PetscStrcmp(ver,"basicnonorms");
9120f5bd95cSBarry Smith     isquad    = !PetscStrcmp(ver,"quadratic");
9130f5bd95cSBarry Smith     iscubic   = !PetscStrcmp(ver,"cubic");
9140f5bd95cSBarry Smith 
9150f5bd95cSBarry Smith     if (isbasic) {
916af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
9170f5bd95cSBarry Smith     } else if (isnonorms) {
918af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
9190f5bd95cSBarry Smith     } else if (isquad) {
920af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
9210f5bd95cSBarry Smith     } else if (iscubic) {
922af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
9235e42470aSBarry Smith     }
924a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
9255e42470aSBarry Smith   }
9263a40ed3dSBarry Smith   PetscFunctionReturn(0);
9275e42470aSBarry Smith }
92804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
92904d965bbSLois Curfman McInnes /*
930*6831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
931*6831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
93204d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
93304d965bbSLois Curfman McInnes 
93404d965bbSLois Curfman McInnes 
93504d965bbSLois Curfman McInnes    Input Parameter:
93604d965bbSLois Curfman McInnes .  snes - the SNES context
93704d965bbSLois Curfman McInnes 
93804d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
93904d965bbSLois Curfman McInnes  */
940fb2e594dSBarry Smith EXTERN_C_BEGIN
9415615d1e5SSatish Balay #undef __FUNC__
9425615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
943f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
9445e42470aSBarry Smith {
94582bf6240SBarry Smith   int        ierr;
946*6831982aSBarry Smith   SNES_EQ_LS *neP;
9475e42470aSBarry Smith 
9483a40ed3dSBarry Smith   PetscFunctionBegin;
949a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
950a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
951a8c6a408SBarry Smith   }
95282bf6240SBarry Smith 
953f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
954f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
955f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
956f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
957f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
958f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
959f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
9605baf8537SBarry Smith   snes->nwork           = 0;
9615e42470aSBarry Smith 
962*6831982aSBarry Smith   neP			= PetscNew(SNES_EQ_LS);CHKPTRQ(neP);
963*6831982aSBarry Smith   PLogObjectMemory(snes,sizeof(SNES_EQ_LS));
9645e42470aSBarry Smith   snes->data    	= (void *) neP;
9655e42470aSBarry Smith   neP->alpha		= 1.e-4;
9665e42470aSBarry Smith   neP->maxstep		= 1.e8;
9675e42470aSBarry Smith   neP->steptol		= 1.e-12;
9685e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
969c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
970c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
971c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
97282bf6240SBarry Smith 
97394d884c6SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
974e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
975c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
976c8dd0c0dSLois Curfman McInnes                     (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
97782bf6240SBarry Smith 
9783a40ed3dSBarry Smith   PetscFunctionReturn(0);
9795e42470aSBarry Smith }
980fb2e594dSBarry Smith EXTERN_C_END
98182bf6240SBarry Smith 
98282bf6240SBarry Smith 
98382bf6240SBarry Smith 
984