xref: /petsc/src/snes/impls/ls/ls.c (revision 8a5d9ee4d611ecbd142307055b81a4eaf96cbeec)
1*8a5d9ee4SBarry Smith /*$Id: ls.c,v 1.156 2000/07/23 16:42:19 bsmith Exp bsmith $*/
25e76c431SBarry Smith 
370f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
45e42470aSBarry Smith 
5*8a5d9ee4SBarry Smith /*
6*8a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
7*8a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
8*8a5d9ee4SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J
9*8a5d9ee4SBarry Smith */
10*8a5d9ee4SBarry Smith #undef __FUNC__
11*8a5d9ee4SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckLocalMin_Private"></a>*/"SNESLSCheckLocalMin_Private"
12*8a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
13*8a5d9ee4SBarry Smith {
14*8a5d9ee4SBarry Smith   PetscReal a1;
15*8a5d9ee4SBarry Smith   int       ierr;
16*8a5d9ee4SBarry Smith 
17*8a5d9ee4SBarry Smith   PetscFunctionBegin;
18*8a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
19*8a5d9ee4SBarry Smith     /* Compute || J^T F|| */
20*8a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
21*8a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
22*8a5d9ee4SBarry Smith     PLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g || F || %g ||J^T F||/|| F || %g\n",a1,fnorm,a1/fnorm);
23*8a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
24*8a5d9ee4SBarry Smith   PetscFunctionReturn(0);
25*8a5d9ee4SBarry Smith }
26*8a5d9ee4SBarry Smith 
2704d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
2804d965bbSLois Curfman McInnes 
2904d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
3004d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
3104d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
3204d965bbSLois Curfman McInnes      respectively.
3304d965bbSLois Curfman McInnes 
34fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
3504d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
3604d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
37fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
3804d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
3904d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
4004d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
4104d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
4204d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
4304d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
4404d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
4504d965bbSLois Curfman McInnes      for all nonlinear solvers.
4604d965bbSLois Curfman McInnes 
4704d965bbSLois Curfman McInnes      Another key routine is:
4804d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
4904d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
5004d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
5104d965bbSLois Curfman McInnes      SNESSolve() if necessary.
5204d965bbSLois Curfman McInnes 
5304d965bbSLois Curfman McInnes      Additional basic routines are:
5404d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
5504d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
5604d965bbSLois Curfman McInnes                                       have actually been used.
5704d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
5804d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
5904d965bbSLois Curfman McInnes 
6004d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
6104d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
6204d965bbSLois Curfman McInnes      above description applies to these categories also.
6304d965bbSLois Curfman McInnes 
6404d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
655e42470aSBarry Smith /*
6604d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
6704d965bbSLois Curfman McInnes    method with a line search.
685e76c431SBarry Smith 
6904d965bbSLois Curfman McInnes    Input Parameters:
7004d965bbSLois Curfman McInnes .  snes - the SNES context
715e76c431SBarry Smith 
7204d965bbSLois Curfman McInnes    Output Parameter:
73c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
7404d965bbSLois Curfman McInnes 
7504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
765e76c431SBarry Smith 
775e76c431SBarry Smith    Notes:
785e76c431SBarry Smith    This implements essentially a truncated Newton method with a
795e76c431SBarry Smith    line search.  By default a cubic backtracking line search
805e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
815e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
82393d2d9aSLois Curfman McInnes    and Schnabel.
835e76c431SBarry Smith */
845615d1e5SSatish Balay #undef __FUNC__
85b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS"
86f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
875e76c431SBarry Smith {
886831982aSBarry Smith   SNES_EQ_LS          *neP = (SNES_EQ_LS*)snes->data;
8984c569c9SBarry Smith   int                 maxits,i,ierr,lits,lsfail;
90112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
91329f5518SBarry Smith   PetscReal           fnorm,gnorm,xnorm,ynorm;
925e42470aSBarry Smith   Vec                 Y,X,F,G,W,TMP;
935e76c431SBarry Smith 
943a40ed3dSBarry Smith   PetscFunctionBegin;
95184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
96184914b5SBarry Smith 
975e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
985e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
9939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1005e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1015e42470aSBarry Smith   G		= snes->work[1];
1025e42470aSBarry Smith   W		= snes->work[2];
1035e76c431SBarry Smith 
1044c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1054c49b128SBarry Smith   snes->iter = 0;
1064c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1075334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
108cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1090f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1105e42470aSBarry Smith   snes->norm = fnorm;
1110f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
11284c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
11394a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1145e76c431SBarry Smith 
115184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
11694a9d846SBarry Smith 
117d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
118d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
119d034289bSBarry Smith 
1205e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1215e76c431SBarry Smith 
122ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1235334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1245334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
12578b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
12690cbd584SBarry Smith     /* should check what happened to the linear solve? */
1273505fcc1SBarry Smith     snes->linear_its += lits;
128af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
129ea4d3ed3SLois Curfman McInnes 
130ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
131ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
132ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
133ea4d3ed3SLois Curfman McInnes     */
13481b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
135af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
136af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
1373505fcc1SBarry Smith     if (lsfail) {
138*8a5d9ee4SBarry Smith       PetscTruth ismin;
1393505fcc1SBarry Smith       snes->nfailures++;
1403505fcc1SBarry Smith       snes->reason = SNES_DIVERGED_LS_FAILURE;
141*8a5d9ee4SBarry Smith       ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
142*8a5d9ee4SBarry Smith       if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
1433505fcc1SBarry Smith       break;
1443505fcc1SBarry Smith     }
1455e76c431SBarry Smith 
14639e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
14739e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1485e76c431SBarry Smith     fnorm = gnorm;
1495e76c431SBarry Smith 
1500f5bd95cSBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
15184c569c9SBarry Smith     snes->iter = i+1;
1525e42470aSBarry Smith     snes->norm = fnorm;
1530f5bd95cSBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15484c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
15594a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1565e76c431SBarry Smith 
1575e76c431SBarry Smith     /* Test for convergence */
15829e0b56fSBarry Smith     if (snes->converged) {
15929e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
1603505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1613505fcc1SBarry Smith       if (snes->reason) {
16216c95cacSBarry Smith         break;
16316c95cacSBarry Smith       }
16416c95cacSBarry Smith     }
16529e0b56fSBarry Smith   }
16639e2f89bSBarry Smith   if (X != snes->vec_sol) {
167393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
168e15f7bb5SBarry Smith   }
169e15f7bb5SBarry Smith   if (F != snes->vec_func) {
170e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
171e15f7bb5SBarry Smith   }
17239e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
17339e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
17452392280SLois Curfman McInnes   if (i == maxits) {
175981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
17652392280SLois Curfman McInnes     i--;
1773505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
17852392280SLois Curfman McInnes   }
1790f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1800f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1815e42470aSBarry Smith   *outits = i+1;
1823a40ed3dSBarry Smith   PetscFunctionReturn(0);
1835e76c431SBarry Smith }
18404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
18504d965bbSLois Curfman McInnes /*
18604d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
1876831982aSBarry Smith    of the SNESEQLS nonlinear solver.
18804d965bbSLois Curfman McInnes 
18904d965bbSLois Curfman McInnes    Input Parameter:
19004d965bbSLois Curfman McInnes .  snes - the SNES context
19104d965bbSLois Curfman McInnes .  x - the solution vector
19204d965bbSLois Curfman McInnes 
19304d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
19404d965bbSLois Curfman McInnes 
19504d965bbSLois Curfman McInnes    Notes:
19604d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
19704d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
19804d965bbSLois Curfman McInnes    the call to SNESSolve().
19904d965bbSLois Curfman McInnes  */
2005615d1e5SSatish Balay #undef __FUNC__
201b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS"
202f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
2035e76c431SBarry Smith {
2045e42470aSBarry Smith   int ierr;
2053a40ed3dSBarry Smith 
2063a40ed3dSBarry Smith   PetscFunctionBegin;
20781b6cf68SLois Curfman McInnes   snes->nwork = 4;
208d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
2095e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
21081b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
2125e76c431SBarry Smith }
21304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
21404d965bbSLois Curfman McInnes /*
2156831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
21604d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
21704d965bbSLois Curfman McInnes 
21804d965bbSLois Curfman McInnes    Input Parameter:
21904d965bbSLois Curfman McInnes .  snes - the SNES context
22004d965bbSLois Curfman McInnes 
221de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
22204d965bbSLois Curfman McInnes  */
2235615d1e5SSatish Balay #undef __FUNC__
224b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS"
225e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
2265e76c431SBarry Smith {
227393d2d9aSLois Curfman McInnes   int  ierr;
2283a40ed3dSBarry Smith 
2293a40ed3dSBarry Smith   PetscFunctionBegin;
2305baf8537SBarry Smith   if (snes->nwork) {
2314b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
23221c89e3eSBarry Smith   }
2335c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2343a40ed3dSBarry Smith   PetscFunctionReturn(0);
2355e76c431SBarry Smith }
23604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2375615d1e5SSatish Balay #undef __FUNC__
238b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch"
23982bf6240SBarry Smith 
2404b828684SBarry Smith /*@C
2415e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2425e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2435e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2445e76c431SBarry Smith 
245c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
246c7afd0dbSLois Curfman McInnes 
2475e76c431SBarry Smith    Input Parameters:
248c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
249af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
2505e76c431SBarry Smith .  x - current iterate
2515e76c431SBarry Smith .  f - residual evaluated at x
2525e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2535e76c431SBarry Smith .  w - work vector
254c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
2555e76c431SBarry Smith 
256c4a48953SLois Curfman McInnes    Output Parameters:
257c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
2585e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2595e76c431SBarry Smith .  gnorm - 2-norm of g
2605e76c431SBarry Smith .  ynorm - 2-norm of search length
261c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
262fee21e36SBarry Smith 
263c4a48953SLois Curfman McInnes    Options Database Key:
264c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
265c4a48953SLois Curfman McInnes 
26636851e7fSLois Curfman McInnes    Level: advanced
26736851e7fSLois Curfman McInnes 
26828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
26928ae5a21SLois Curfman McInnes 
270f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
271af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
2725e76c431SBarry Smith @*/
273af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
274329f5518SBarry Smith                      PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
2755e76c431SBarry Smith {
2765e42470aSBarry Smith   int        ierr;
2775334005bSBarry Smith   Scalar     mone = -1.0;
2786831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
2798f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
2805334005bSBarry Smith 
2813a40ed3dSBarry Smith   PetscFunctionBegin;
282761aaf1bSLois Curfman McInnes   *flag = 0;
2837857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
284a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
285ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
2868f99978dSLois Curfman McInnes   if (neP->CheckStep) {
2878f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
2888f99978dSLois Curfman McInnes   }
289ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
290a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
2917857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2923a40ed3dSBarry Smith   PetscFunctionReturn(0);
2935e76c431SBarry Smith }
29404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2955615d1e5SSatish Balay #undef __FUNC__
296b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms"
29782bf6240SBarry Smith 
29829e0b56fSBarry Smith /*@C
29929e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
30029e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
30129e0b56fSBarry Smith    even compute the norm of the function or search direction; this
30229e0b56fSBarry Smith    is intended only when you know the full step is fine and are
30329e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
304c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
305c7afd0dbSLois Curfman McInnes 
306c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
30729e0b56fSBarry Smith 
30829e0b56fSBarry Smith    Input Parameters:
309c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
310af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
31129e0b56fSBarry Smith .  x - current iterate
31229e0b56fSBarry Smith .  f - residual evaluated at x
31329e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
31429e0b56fSBarry Smith .  w - work vector
315c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
31629e0b56fSBarry Smith 
31729e0b56fSBarry Smith    Output Parameters:
318c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
31929e0b56fSBarry Smith .  gnorm - not changed
32029e0b56fSBarry Smith .  ynorm - not changed
321c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
322fee21e36SBarry Smith 
32329e0b56fSBarry Smith    Options Database Key:
324c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
32529e0b56fSBarry Smith 
3268cbba510SBarry Smith    Notes:
327ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
328ea56f5baSLois Curfman McInnes    either the options
329ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
330ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
331329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
332329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
333329f5518SBarry Smith 
334329f5518SBarry Smith    During the final iteration this will not evaluate the function at
335329f5518SBarry Smith    the solution point. This is to save a function evaluation while
336329f5518SBarry Smith    using pseudo-timestepping.
3378cbba510SBarry Smith 
338ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
339ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
340ea56f5baSLois Curfman McInnes    correct, since they are not computed.
341ea56f5baSLois Curfman McInnes 
342ea56f5baSLois Curfman McInnes    Level: advanced
3438cbba510SBarry Smith 
34429e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
34529e0b56fSBarry Smith 
34629e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
34729e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
34829e0b56fSBarry Smith @*/
349af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
350329f5518SBarry Smith                             PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
35129e0b56fSBarry Smith {
35229e0b56fSBarry Smith   int        ierr;
35329e0b56fSBarry Smith   Scalar     mone = -1.0;
3546831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
3558f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
35629e0b56fSBarry Smith 
3573a40ed3dSBarry Smith   PetscFunctionBegin;
3588cbba510SBarry Smith   *flag = 0;
35929e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
36029e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3618f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3628f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3638f99978dSLois Curfman McInnes   }
364329f5518SBarry Smith 
365329f5518SBarry Smith   /* don't evaluate function the last time through */
366329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
36729e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
368329f5518SBarry Smith   }
36929e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3703a40ed3dSBarry Smith   PetscFunctionReturn(0);
37129e0b56fSBarry Smith }
37204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
37329e0b56fSBarry Smith #undef __FUNC__
374b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch"
3754b828684SBarry Smith /*@C
376f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3775e76c431SBarry Smith 
378c7afd0dbSLois Curfman McInnes    Collective on SNES
379c7afd0dbSLois Curfman McInnes 
3805e76c431SBarry Smith    Input Parameters:
381c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
382af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3835e76c431SBarry Smith .  x - current iterate
3845e76c431SBarry Smith .  f - residual evaluated at x
3855e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3865e76c431SBarry Smith .  w - work vector
387c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3885e76c431SBarry Smith 
389393d2d9aSLois Curfman McInnes    Output Parameters:
390c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3915e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3925e76c431SBarry Smith .  gnorm - 2-norm of g
3935e76c431SBarry Smith .  ynorm - 2-norm of search length
394c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
395fee21e36SBarry Smith 
396c4a48953SLois Curfman McInnes    Options Database Key:
397c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
398c4a48953SLois Curfman McInnes 
3995e76c431SBarry Smith    Notes:
4005e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4015e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4025e76c431SBarry Smith 
40336851e7fSLois Curfman McInnes    Level: advanced
40436851e7fSLois Curfman McInnes 
40528ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
40628ae5a21SLois Curfman McInnes 
407af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4085e76c431SBarry Smith @*/
409af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
410329f5518SBarry Smith                         PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
4115e76c431SBarry Smith {
412406556e6SLois Curfman McInnes   /*
413406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
414406556e6SLois Curfman McInnes      minimization problem:
415406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
416406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
417406556e6SLois Curfman McInnes    */
418406556e6SLois Curfman McInnes 
419329f5518SBarry Smith   PetscReal  steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
420329f5518SBarry Smith   PetscReal  maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
421aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4225e42470aSBarry Smith   Scalar     cinitslope,clambda;
4236b5873e3SBarry Smith #endif
4245e42470aSBarry Smith   int        ierr,count;
4256831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
4265334005bSBarry Smith   Scalar     mone = -1.0,scale;
4278f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
4285e76c431SBarry Smith 
4293a40ed3dSBarry Smith   PetscFunctionBegin;
4307857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
431761aaf1bSLois Curfman McInnes   *flag   = 0;
4325e76c431SBarry Smith   alpha   = neP->alpha;
4335e76c431SBarry Smith   maxstep = neP->maxstep;
4345e76c431SBarry Smith   steptol = neP->steptol;
4355e76c431SBarry Smith 
436cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
437a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
438a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
439a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
440a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
441a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
442ad922ac9SBarry Smith     goto theend1;
44394a9d846SBarry Smith   }
4445e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4455e42470aSBarry Smith     scale = maxstep/(*ynorm);
446aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
44790cbd584SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
4486b5873e3SBarry Smith #else
44990cbd584SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
4506b5873e3SBarry Smith #endif
451393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
4525e76c431SBarry Smith     *ynorm = maxstep;
4535e76c431SBarry Smith   }
4545e76c431SBarry Smith   minlambda = steptol/(*ynorm);
455a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
456aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
457a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
458329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
4595e42470aSBarry Smith #else
460a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
4615e42470aSBarry Smith #endif
4625e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
4635e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
4645e76c431SBarry Smith 
465393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
4665334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
46778b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
468313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
469406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
470393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
47194a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
472ad922ac9SBarry Smith     goto theend1;
4735e76c431SBarry Smith   }
4745e76c431SBarry Smith 
4755e76c431SBarry Smith   /* Fit points with quadratic */
476313b5538SBarry Smith   lambda = 1.0;
477a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4785e76c431SBarry Smith   lambdaprev = lambda;
4795e76c431SBarry Smith   gnormprev = *gnorm;
480329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
481ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
482ddd12767SBarry Smith   else                         lambda = lambdatemp;
483393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
484ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
485aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
486ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
4875e42470aSBarry Smith #else
488ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
4895e42470aSBarry Smith #endif
49078b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
491cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
492406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
493393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
494ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
495ad922ac9SBarry Smith     goto theend1;
4965e76c431SBarry Smith   }
4975e76c431SBarry Smith 
4985e76c431SBarry Smith   /* Fit points with cubic */
4995e76c431SBarry Smith   count = 1;
5005e76c431SBarry Smith   while (1) {
5015e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
5022b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
50390cbd584SBarry Smith       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
504393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
505761aaf1bSLois Curfman McInnes       *flag = -1; break;
5065e76c431SBarry Smith     }
507406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
508406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
509ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5102b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5115e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5125e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5135e76c431SBarry Smith     if (a == 0.0) {
5145e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5155e76c431SBarry Smith     } else {
5165e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5175e76c431SBarry Smith     }
5185e76c431SBarry Smith     lambdaprev = lambda;
5195e76c431SBarry Smith     gnormprev  = *gnorm;
520329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
521329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
5225e76c431SBarry Smith     else                         lambda     = lambdatemp;
523393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
524ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
525aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
526ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
527393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5285e42470aSBarry Smith #else
529ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5305e42470aSBarry Smith #endif
53178b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
532cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
533406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
534393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
535ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
5362715565aSLois Curfman McInnes       goto theend1;
5372b022350SLois Curfman McInnes     } else {
5382b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
5395e76c431SBarry Smith     }
5405e76c431SBarry Smith     count++;
5415e76c431SBarry Smith   }
542ad922ac9SBarry Smith   theend1:
5438f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5448f99978dSLois Curfman McInnes   if (neP->CheckStep) {
5458f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
5468f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
5478f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
5488f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
5498f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
5508f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
5518f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
5528f99978dSLois Curfman McInnes     }
5538f99978dSLois Curfman McInnes   }
5547857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5553a40ed3dSBarry Smith   PetscFunctionReturn(0);
5565e76c431SBarry Smith }
55704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
5585615d1e5SSatish Balay #undef __FUNC__
559b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch"
5604b828684SBarry Smith /*@C
561f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
5625e76c431SBarry Smith 
563c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
564c7afd0dbSLois Curfman McInnes 
5655e42470aSBarry Smith    Input Parameters:
566c7afd0dbSLois Curfman McInnes +  snes - the SNES context
567af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
5685e42470aSBarry Smith .  x - current iterate
5695e42470aSBarry Smith .  f - residual evaluated at x
5705e42470aSBarry Smith .  y - search direction (contains new iterate on output)
5715e42470aSBarry Smith .  w - work vector
572c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
5735e42470aSBarry Smith 
574c4a48953SLois Curfman McInnes    Output Parameters:
575c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
5765e42470aSBarry Smith .  y - new iterate (contains search direction on input)
5775e42470aSBarry Smith .  gnorm - 2-norm of g
5785e42470aSBarry Smith .  ynorm - 2-norm of search length
579c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
580fee21e36SBarry Smith 
581c4a48953SLois Curfman McInnes    Options Database Key:
582c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
583c4a48953SLois Curfman McInnes 
5845e42470aSBarry Smith    Notes:
5856831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
5865e42470aSBarry Smith 
58736851e7fSLois Curfman McInnes    Level: advanced
58836851e7fSLois Curfman McInnes 
589f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
590f59ffdedSLois Curfman McInnes 
591af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
5925e42470aSBarry Smith @*/
593af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
594329f5518SBarry Smith                            PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
5955e76c431SBarry Smith {
596406556e6SLois Curfman McInnes   /*
597406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
598406556e6SLois Curfman McInnes      minimization problem:
599406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
600406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
601406556e6SLois Curfman McInnes    */
602329f5518SBarry Smith   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
603aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6045e42470aSBarry Smith   Scalar     cinitslope,clambda;
6056b5873e3SBarry Smith #endif
6065e42470aSBarry Smith   int        ierr,count;
6076831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
6085334005bSBarry Smith   Scalar     mone = -1.0,scale;
6098f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
6105e76c431SBarry Smith 
6113a40ed3dSBarry Smith   PetscFunctionBegin;
6127857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
613761aaf1bSLois Curfman McInnes   *flag   = 0;
6145e42470aSBarry Smith   alpha   = neP->alpha;
6155e42470aSBarry Smith   maxstep = neP->maxstep;
6165e42470aSBarry Smith   steptol = neP->steptol;
6175e76c431SBarry Smith 
6183505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
619b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
62094a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
621b37302c6SLois Curfman McInnes     *gnorm = fnorm;
622b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
623b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
624ad922ac9SBarry Smith     goto theend2;
62594a9d846SBarry Smith   }
6265e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6275e42470aSBarry Smith     scale = maxstep/(*ynorm);
628393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6295e42470aSBarry Smith     *ynorm = maxstep;
6305e76c431SBarry Smith   }
6315e42470aSBarry Smith   minlambda = steptol/(*ynorm);
632a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
633aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
634a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
635329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
6365e42470aSBarry Smith #else
637a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6385e42470aSBarry Smith #endif
6395e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
6405e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6415e42470aSBarry Smith 
642393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6435334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
64478b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
645cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
646406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
647393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
64894a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
649ad922ac9SBarry Smith     goto theend2;
6505e42470aSBarry Smith   }
6515e42470aSBarry Smith 
6525e42470aSBarry Smith   /* Fit points with quadratic */
653313b5538SBarry Smith   lambda = 1.0;
6545e42470aSBarry Smith   count = 1;
6555e42470aSBarry Smith   while (1) {
6565e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
657981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
658981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
659ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
660393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
661761aaf1bSLois Curfman McInnes       *flag = -1; break;
6625e42470aSBarry Smith     }
663a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
664329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
665329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
666329f5518SBarry Smith     else                         lambda     = lambdatemp;
667393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
6683505fcc1SBarry Smith     lambdaneg = -lambda;
669aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6703505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6715e42470aSBarry Smith #else
6723505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6735e42470aSBarry Smith #endif
67478b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
675cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
676406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
677393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
678981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
67906259719SBarry Smith       break;
6805e42470aSBarry Smith     }
6815e42470aSBarry Smith     count++;
6825e42470aSBarry Smith   }
683ad922ac9SBarry Smith   theend2:
6848f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6858f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6868f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6878f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6888f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6898f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6908f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6918f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6928f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6938f99978dSLois Curfman McInnes     }
6948f99978dSLois Curfman McInnes   }
6957857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
6963a40ed3dSBarry Smith   PetscFunctionReturn(0);
6975e76c431SBarry Smith }
69804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6995615d1e5SSatish Balay #undef __FUNC__
700b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch"
701c9e6a524SLois Curfman McInnes /*@C
70277c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7036831982aSBarry Smith    by the method SNESEQLS.
7045e76c431SBarry Smith 
7055e76c431SBarry Smith    Input Parameters:
706c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
707af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
708c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7095e76c431SBarry Smith 
710fee21e36SBarry Smith    Collective on SNES
711fee21e36SBarry Smith 
712c4a48953SLois Curfman McInnes    Available Routines:
713c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
714f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
715f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
716af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7175e76c431SBarry Smith 
718c4a48953SLois Curfman McInnes     Options Database Keys:
719af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
720c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
721c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
722c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
723c4a48953SLois Curfman McInnes 
7245e76c431SBarry Smith    Calling sequence of func:
725bd208895SLois Curfman McInnes .vb
726af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
727329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
728bd208895SLois Curfman McInnes .ve
7295e76c431SBarry Smith 
7305e76c431SBarry Smith     Input parameters for func:
731c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
732af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
7335e76c431SBarry Smith .   x - current iterate
7345e76c431SBarry Smith .   f - residual evaluated at x
7355e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7365e76c431SBarry Smith .   w - work vector
737c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
7385e76c431SBarry Smith 
7395e76c431SBarry Smith     Output parameters for func:
740c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7415e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7425e76c431SBarry Smith .   gnorm - 2-norm of g
7435e76c431SBarry Smith .   ynorm - 2-norm of search length
744c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
745761aaf1bSLois Curfman McInnes            on failure.
746f59ffdedSLois Curfman McInnes 
74736851e7fSLois Curfman McInnes     Level: advanced
74836851e7fSLois Curfman McInnes 
749f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
750f59ffdedSLois Curfman McInnes 
7514ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
7525e76c431SBarry Smith @*/
753329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
7545e76c431SBarry Smith {
755329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
75682bf6240SBarry Smith 
7573a40ed3dSBarry Smith   PetscFunctionBegin;
75894d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
75982bf6240SBarry Smith   if (f) {
760af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
76182bf6240SBarry Smith   }
7623a40ed3dSBarry Smith   PetscFunctionReturn(0);
7635e76c431SBarry Smith }
76404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
765fb2e594dSBarry Smith EXTERN_C_BEGIN
76682bf6240SBarry Smith #undef __FUNC__
767b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS"
768af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
769af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
77082bf6240SBarry Smith {
77182bf6240SBarry Smith   PetscFunctionBegin;
7726831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
7736831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
77482bf6240SBarry Smith   PetscFunctionReturn(0);
77582bf6240SBarry Smith }
776fb2e594dSBarry Smith EXTERN_C_END
77704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
778c8dd0c0dSLois Curfman McInnes #undef __FUNC__
779b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck"
780c8dd0c0dSLois Curfman McInnes /*@C
781530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
7826831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
783c8dd0c0dSLois Curfman McInnes 
784c8dd0c0dSLois Curfman McInnes    Input Parameters:
785c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
786c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
787c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
788c8dd0c0dSLois Curfman McInnes 
789c8dd0c0dSLois Curfman McInnes    Collective on SNES
790c8dd0c0dSLois Curfman McInnes 
791c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
792c8dd0c0dSLois Curfman McInnes .vb
793b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
794c8dd0c0dSLois Curfman McInnes .ve
795b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
796b1ae27eaSLois Curfman McInnes    on failure.
797c8dd0c0dSLois Curfman McInnes 
798c8dd0c0dSLois Curfman McInnes    Input parameters for func:
799c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
800c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8018f99978dSLois Curfman McInnes -  x - current candidate iterate
802c8dd0c0dSLois Curfman McInnes 
803c8dd0c0dSLois Curfman McInnes    Output parameters for func:
804c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
805c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
806c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
807c8dd0c0dSLois Curfman McInnes 
808c8dd0c0dSLois Curfman McInnes    Level: advanced
809c8dd0c0dSLois Curfman McInnes 
8108f99978dSLois Curfman McInnes    Notes:
811b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
812b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
813b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
814b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
815ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8168f99978dSLois Curfman McInnes 
817b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
818b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
819b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
820ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
821ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
822ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
823ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
824b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8258f99978dSLois Curfman McInnes 
826c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
827c8dd0c0dSLois Curfman McInnes 
828c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
829c8dd0c0dSLois Curfman McInnes @*/
8308f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
831c8dd0c0dSLois Curfman McInnes {
8328f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
833c8dd0c0dSLois Curfman McInnes 
834c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
835c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
836c8dd0c0dSLois Curfman McInnes   if (f) {
837c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
838c8dd0c0dSLois Curfman McInnes   }
839c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
840c8dd0c0dSLois Curfman McInnes }
841c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
842c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
843c8dd0c0dSLois Curfman McInnes #undef __FUNC__
844b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS"
8458f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
846c8dd0c0dSLois Curfman McInnes {
847c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
8486831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
8496831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
850c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
851c8dd0c0dSLois Curfman McInnes }
852c8dd0c0dSLois Curfman McInnes EXTERN_C_END
853c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
85404d965bbSLois Curfman McInnes /*
8556831982aSBarry Smith    SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method.
85682bf6240SBarry Smith 
85704d965bbSLois Curfman McInnes    Input Parameter:
85804d965bbSLois Curfman McInnes .  snes - the SNES context
85904d965bbSLois Curfman McInnes 
86004d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
86104d965bbSLois Curfman McInnes */
8625615d1e5SSatish Balay #undef __FUNC__
863b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS"
864f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
8655e42470aSBarry Smith {
8666831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
867d132466eSBarry Smith   int     ierr;
8686daaf66cSBarry Smith 
8693a40ed3dSBarry Smith   PetscFunctionBegin;
8706831982aSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr);
871d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr);
872d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr);
873d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr);
874d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr);
8753a40ed3dSBarry Smith   PetscFunctionReturn(0);
8765e42470aSBarry Smith }
87704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
87804d965bbSLois Curfman McInnes /*
8796831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
88004d965bbSLois Curfman McInnes 
88104d965bbSLois Curfman McInnes    Input Parameters:
88204d965bbSLois Curfman McInnes .  SNES - the SNES context
88304d965bbSLois Curfman McInnes .  viewer - visualization context
88404d965bbSLois Curfman McInnes 
88504d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
88604d965bbSLois Curfman McInnes */
8875615d1e5SSatish Balay #undef __FUNC__
888b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS"
889e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
890a935fc98SLois Curfman McInnes {
8916831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
89219bcc07fSBarry Smith   char       *cstr;
89351695f54SLois Curfman McInnes   int        ierr;
8946831982aSBarry Smith   PetscTruth isascii;
895a935fc98SLois Curfman McInnes 
8963a40ed3dSBarry Smith   PetscFunctionBegin;
8976831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
8980f5bd95cSBarry Smith   if (isascii) {
89919bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
90019bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
90119bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
90219bcc07fSBarry Smith     else                                                cstr = "unknown";
903d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
904d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9055cd90555SBarry Smith   } else {
9060f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
90719bcc07fSBarry Smith   }
9083a40ed3dSBarry Smith   PetscFunctionReturn(0);
909a935fc98SLois Curfman McInnes }
91004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
91104d965bbSLois Curfman McInnes /*
9126831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
91304d965bbSLois Curfman McInnes 
91404d965bbSLois Curfman McInnes    Input Parameter:
91504d965bbSLois Curfman McInnes .  snes - the SNES context
91604d965bbSLois Curfman McInnes 
91704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
91804d965bbSLois Curfman McInnes */
9195615d1e5SSatish Balay #undef __FUNC__
920b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS"
921f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
9225e42470aSBarry Smith {
9236831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
9245e42470aSBarry Smith   char       ver[16];
9255e42470aSBarry Smith   double     tmp;
926f1af5d2fSBarry Smith   int        ierr;
927f1af5d2fSBarry Smith   PetscTruth flg;
9285e42470aSBarry Smith 
9293a40ed3dSBarry Smith   PetscFunctionBegin;
93009d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr);
93125cdf11fSBarry Smith   if (flg) {
9325e42470aSBarry Smith     ls->alpha = tmp;
9335e42470aSBarry Smith   }
93409d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr);
93525cdf11fSBarry Smith   if (flg) {
9365e42470aSBarry Smith     ls->maxstep = tmp;
9375e42470aSBarry Smith   }
93809d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr);
93925cdf11fSBarry Smith   if (flg) {
9405e42470aSBarry Smith     ls->steptol = tmp;
9415e42470aSBarry Smith   }
94209d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr);
94325cdf11fSBarry Smith   if (flg) {
944c38d4ed2SBarry Smith     PetscTruth isbasic,isnonorms,isquad,iscubic;
9450f5bd95cSBarry Smith 
946c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr);
947c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr);
948c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr);
949c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr);
9500f5bd95cSBarry Smith 
9510f5bd95cSBarry Smith     if (isbasic) {
952af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
9530f5bd95cSBarry Smith     } else if (isnonorms) {
954af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
9550f5bd95cSBarry Smith     } else if (isquad) {
956af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
9570f5bd95cSBarry Smith     } else if (iscubic) {
958af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
9595e42470aSBarry Smith     }
960a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
9615e42470aSBarry Smith   }
9623a40ed3dSBarry Smith   PetscFunctionReturn(0);
9635e42470aSBarry Smith }
96404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
96504d965bbSLois Curfman McInnes /*
9666831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
9676831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
96804d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
96904d965bbSLois Curfman McInnes 
97004d965bbSLois Curfman McInnes 
97104d965bbSLois Curfman McInnes    Input Parameter:
97204d965bbSLois Curfman McInnes .  snes - the SNES context
97304d965bbSLois Curfman McInnes 
97404d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
97504d965bbSLois Curfman McInnes  */
976fb2e594dSBarry Smith EXTERN_C_BEGIN
9775615d1e5SSatish Balay #undef __FUNC__
978b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS"
979f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
9805e42470aSBarry Smith {
98182bf6240SBarry Smith   int        ierr;
9826831982aSBarry Smith   SNES_EQ_LS *neP;
9835e42470aSBarry Smith 
9843a40ed3dSBarry Smith   PetscFunctionBegin;
985a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
986a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
987a8c6a408SBarry Smith   }
98882bf6240SBarry Smith 
989f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
990f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
991f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
992f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
993f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
994f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
995f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
9965baf8537SBarry Smith   snes->nwork           = 0;
9975e42470aSBarry Smith 
9986831982aSBarry Smith   neP			= PetscNew(SNES_EQ_LS);CHKPTRQ(neP);
9996831982aSBarry Smith   PLogObjectMemory(snes,sizeof(SNES_EQ_LS));
10005e42470aSBarry Smith   snes->data    	= (void*)neP;
10015e42470aSBarry Smith   neP->alpha		= 1.e-4;
10025e42470aSBarry Smith   neP->maxstep		= 1.e8;
10035e42470aSBarry Smith   neP->steptol		= 1.e-12;
10045e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1005c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1006c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1007c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
100882bf6240SBarry Smith 
1009f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
10100c97f09dSSatish Balay                     SNESSetLineSearch_LS);CHKERRQ(ierr);
1011f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
10120c97f09dSSatish Balay                     SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
101382bf6240SBarry Smith 
10143a40ed3dSBarry Smith   PetscFunctionReturn(0);
10155e42470aSBarry Smith }
1016fb2e594dSBarry Smith EXTERN_C_END
101782bf6240SBarry Smith 
101882bf6240SBarry Smith 
101982bf6240SBarry Smith 
1020