xref: /petsc/src/snes/impls/ls/ls.c (revision 36669109cfbc5b99b88f1bf1ba0194db2654cc8e)
1*36669109SBarry Smith /*$Id: ls.c,v 1.157 2000/07/24 03:32:47 bsmith Exp bsmith $*/
25e76c431SBarry Smith 
370f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
45e42470aSBarry Smith 
58a5d9ee4SBarry Smith /*
68a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
78a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
8*36669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
9*36669109SBarry Smith     for this trick.
108a5d9ee4SBarry Smith */
118a5d9ee4SBarry Smith #undef __FUNC__
128a5d9ee4SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckLocalMin_Private"></a>*/"SNESLSCheckLocalMin_Private"
138a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
148a5d9ee4SBarry Smith {
158a5d9ee4SBarry Smith   PetscReal  a1;
168a5d9ee4SBarry Smith   int        ierr;
17*36669109SBarry Smith   PetscTruth hastranspose;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
21*36669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
22*36669109SBarry Smith   if (hastranspose) {
238a5d9ee4SBarry Smith     /* Compute || J^T F|| */
248a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
258a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
268a5d9ee4SBarry Smith     PLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g || F || %g ||J^T F||/|| F || %g\n",a1,fnorm,a1/fnorm);
278a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
28*36669109SBarry Smith   } else {
29*36669109SBarry Smith     Vec       work;
30*36669109SBarry Smith     Scalar    result;
31*36669109SBarry Smith     PetscReal wnorm;
32*36669109SBarry Smith 
33*36669109SBarry Smith     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
34*36669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
35*36669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
36*36669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
37*36669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
38*36669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
39*36669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
40*36669109SBarry Smith     PLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g\n",a1);
41*36669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
42*36669109SBarry Smith   }
438a5d9ee4SBarry Smith   PetscFunctionReturn(0);
448a5d9ee4SBarry Smith }
458a5d9ee4SBarry Smith 
4604d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
4704d965bbSLois Curfman McInnes 
4804d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
4904d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
5004d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
5104d965bbSLois Curfman McInnes      respectively.
5204d965bbSLois Curfman McInnes 
53fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
5404d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
5504d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
56fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
5704d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
5804d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
5904d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
6004d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
6104d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
6204d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
6304d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
6404d965bbSLois Curfman McInnes      for all nonlinear solvers.
6504d965bbSLois Curfman McInnes 
6604d965bbSLois Curfman McInnes      Another key routine is:
6704d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
6804d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
6904d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
7004d965bbSLois Curfman McInnes      SNESSolve() if necessary.
7104d965bbSLois Curfman McInnes 
7204d965bbSLois Curfman McInnes      Additional basic routines are:
7304d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
7404d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
7504d965bbSLois Curfman McInnes                                       have actually been used.
7604d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
7704d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
7804d965bbSLois Curfman McInnes 
7904d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
8004d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
8104d965bbSLois Curfman McInnes      above description applies to these categories also.
8204d965bbSLois Curfman McInnes 
8304d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
845e42470aSBarry Smith /*
8504d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
8604d965bbSLois Curfman McInnes    method with a line search.
875e76c431SBarry Smith 
8804d965bbSLois Curfman McInnes    Input Parameters:
8904d965bbSLois Curfman McInnes .  snes - the SNES context
905e76c431SBarry Smith 
9104d965bbSLois Curfman McInnes    Output Parameter:
92c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
9304d965bbSLois Curfman McInnes 
9404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
955e76c431SBarry Smith 
965e76c431SBarry Smith    Notes:
975e76c431SBarry Smith    This implements essentially a truncated Newton method with a
985e76c431SBarry Smith    line search.  By default a cubic backtracking line search
995e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1005e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
101393d2d9aSLois Curfman McInnes    and Schnabel.
1025e76c431SBarry Smith */
1035615d1e5SSatish Balay #undef __FUNC__
104b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS"
105f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
1065e76c431SBarry Smith {
1076831982aSBarry Smith   SNES_EQ_LS          *neP = (SNES_EQ_LS*)snes->data;
10884c569c9SBarry Smith   int                 maxits,i,ierr,lits,lsfail;
109112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
110329f5518SBarry Smith   PetscReal           fnorm,gnorm,xnorm,ynorm;
1115e42470aSBarry Smith   Vec                 Y,X,F,G,W,TMP;
1125e76c431SBarry Smith 
1133a40ed3dSBarry Smith   PetscFunctionBegin;
114184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
115184914b5SBarry Smith 
1165e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1175e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
11839e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1195e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1205e42470aSBarry Smith   G		= snes->work[1];
1215e42470aSBarry Smith   W		= snes->work[2];
1225e76c431SBarry Smith 
1234c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1244c49b128SBarry Smith   snes->iter = 0;
1254c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1265334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
127cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1280f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1295e42470aSBarry Smith   snes->norm = fnorm;
1300f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
13184c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
13294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1335e76c431SBarry Smith 
134184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
13594a9d846SBarry Smith 
136d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
137d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
138d034289bSBarry Smith 
1395e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1405e76c431SBarry Smith 
141ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1425334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1435334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
14478b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
14590cbd584SBarry Smith     /* should check what happened to the linear solve? */
1463505fcc1SBarry Smith     snes->linear_its += lits;
147af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
148ea4d3ed3SLois Curfman McInnes 
149ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
150ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
151ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
152ea4d3ed3SLois Curfman McInnes     */
15381b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
154af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
155af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
1563505fcc1SBarry Smith     if (lsfail) {
1578a5d9ee4SBarry Smith       PetscTruth ismin;
1583505fcc1SBarry Smith       snes->nfailures++;
1593505fcc1SBarry Smith       snes->reason = SNES_DIVERGED_LS_FAILURE;
1608a5d9ee4SBarry Smith       ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
1618a5d9ee4SBarry Smith       if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
1623505fcc1SBarry Smith       break;
1633505fcc1SBarry Smith     }
1645e76c431SBarry Smith 
16539e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
16639e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1675e76c431SBarry Smith     fnorm = gnorm;
1685e76c431SBarry Smith 
1690f5bd95cSBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
17084c569c9SBarry Smith     snes->iter = i+1;
1715e42470aSBarry Smith     snes->norm = fnorm;
1720f5bd95cSBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17384c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
17494a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1755e76c431SBarry Smith 
1765e76c431SBarry Smith     /* Test for convergence */
17729e0b56fSBarry Smith     if (snes->converged) {
17829e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
1793505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1803505fcc1SBarry Smith       if (snes->reason) {
18116c95cacSBarry Smith         break;
18216c95cacSBarry Smith       }
18316c95cacSBarry Smith     }
18429e0b56fSBarry Smith   }
18539e2f89bSBarry Smith   if (X != snes->vec_sol) {
186393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
187e15f7bb5SBarry Smith   }
188e15f7bb5SBarry Smith   if (F != snes->vec_func) {
189e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
190e15f7bb5SBarry Smith   }
19139e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
19239e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
19352392280SLois Curfman McInnes   if (i == maxits) {
194981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
19552392280SLois Curfman McInnes     i--;
1963505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
19752392280SLois Curfman McInnes   }
1980f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1990f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2005e42470aSBarry Smith   *outits = i+1;
2013a40ed3dSBarry Smith   PetscFunctionReturn(0);
2025e76c431SBarry Smith }
20304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
20404d965bbSLois Curfman McInnes /*
20504d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
2066831982aSBarry Smith    of the SNESEQLS nonlinear solver.
20704d965bbSLois Curfman McInnes 
20804d965bbSLois Curfman McInnes    Input Parameter:
20904d965bbSLois Curfman McInnes .  snes - the SNES context
21004d965bbSLois Curfman McInnes .  x - the solution vector
21104d965bbSLois Curfman McInnes 
21204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
21304d965bbSLois Curfman McInnes 
21404d965bbSLois Curfman McInnes    Notes:
21504d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
21604d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
21704d965bbSLois Curfman McInnes    the call to SNESSolve().
21804d965bbSLois Curfman McInnes  */
2195615d1e5SSatish Balay #undef __FUNC__
220b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS"
221f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
2225e76c431SBarry Smith {
2235e42470aSBarry Smith   int ierr;
2243a40ed3dSBarry Smith 
2253a40ed3dSBarry Smith   PetscFunctionBegin;
22681b6cf68SLois Curfman McInnes   snes->nwork = 4;
227d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
2285e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
22981b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2303a40ed3dSBarry Smith   PetscFunctionReturn(0);
2315e76c431SBarry Smith }
23204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
23304d965bbSLois Curfman McInnes /*
2346831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
23504d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
23604d965bbSLois Curfman McInnes 
23704d965bbSLois Curfman McInnes    Input Parameter:
23804d965bbSLois Curfman McInnes .  snes - the SNES context
23904d965bbSLois Curfman McInnes 
240de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
24104d965bbSLois Curfman McInnes  */
2425615d1e5SSatish Balay #undef __FUNC__
243b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS"
244e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
2455e76c431SBarry Smith {
246393d2d9aSLois Curfman McInnes   int  ierr;
2473a40ed3dSBarry Smith 
2483a40ed3dSBarry Smith   PetscFunctionBegin;
2495baf8537SBarry Smith   if (snes->nwork) {
2504b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
25121c89e3eSBarry Smith   }
2525c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
2545e76c431SBarry Smith }
25504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2565615d1e5SSatish Balay #undef __FUNC__
257b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch"
25882bf6240SBarry Smith 
2594b828684SBarry Smith /*@C
2605e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2615e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2625e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2635e76c431SBarry Smith 
264c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
265c7afd0dbSLois Curfman McInnes 
2665e76c431SBarry Smith    Input Parameters:
267c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
268af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
2695e76c431SBarry Smith .  x - current iterate
2705e76c431SBarry Smith .  f - residual evaluated at x
2715e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2725e76c431SBarry Smith .  w - work vector
273c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
2745e76c431SBarry Smith 
275c4a48953SLois Curfman McInnes    Output Parameters:
276c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
2775e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2785e76c431SBarry Smith .  gnorm - 2-norm of g
2795e76c431SBarry Smith .  ynorm - 2-norm of search length
280c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
281fee21e36SBarry Smith 
282c4a48953SLois Curfman McInnes    Options Database Key:
283c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
284c4a48953SLois Curfman McInnes 
28536851e7fSLois Curfman McInnes    Level: advanced
28636851e7fSLois Curfman McInnes 
28728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
28828ae5a21SLois Curfman McInnes 
289f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
290af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
2915e76c431SBarry Smith @*/
292af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
293329f5518SBarry Smith                      PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
2945e76c431SBarry Smith {
2955e42470aSBarry Smith   int        ierr;
2965334005bSBarry Smith   Scalar     mone = -1.0;
2976831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
2988f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
2995334005bSBarry Smith 
3003a40ed3dSBarry Smith   PetscFunctionBegin;
301761aaf1bSLois Curfman McInnes   *flag = 0;
3027857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
303a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
304ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3058f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3068f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3078f99978dSLois Curfman McInnes   }
308ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
309a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
3107857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3113a40ed3dSBarry Smith   PetscFunctionReturn(0);
3125e76c431SBarry Smith }
31304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3145615d1e5SSatish Balay #undef __FUNC__
315b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms"
31682bf6240SBarry Smith 
31729e0b56fSBarry Smith /*@C
31829e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
31929e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
32029e0b56fSBarry Smith    even compute the norm of the function or search direction; this
32129e0b56fSBarry Smith    is intended only when you know the full step is fine and are
32229e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
323c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
324c7afd0dbSLois Curfman McInnes 
325c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
32629e0b56fSBarry Smith 
32729e0b56fSBarry Smith    Input Parameters:
328c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
329af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
33029e0b56fSBarry Smith .  x - current iterate
33129e0b56fSBarry Smith .  f - residual evaluated at x
33229e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
33329e0b56fSBarry Smith .  w - work vector
334c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
33529e0b56fSBarry Smith 
33629e0b56fSBarry Smith    Output Parameters:
337c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
33829e0b56fSBarry Smith .  gnorm - not changed
33929e0b56fSBarry Smith .  ynorm - not changed
340c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
341fee21e36SBarry Smith 
34229e0b56fSBarry Smith    Options Database Key:
343c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
34429e0b56fSBarry Smith 
3458cbba510SBarry Smith    Notes:
346ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
347ea56f5baSLois Curfman McInnes    either the options
348ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
349ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
350329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
351329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
352329f5518SBarry Smith 
353329f5518SBarry Smith    During the final iteration this will not evaluate the function at
354329f5518SBarry Smith    the solution point. This is to save a function evaluation while
355329f5518SBarry Smith    using pseudo-timestepping.
3568cbba510SBarry Smith 
357ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
358ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
359ea56f5baSLois Curfman McInnes    correct, since they are not computed.
360ea56f5baSLois Curfman McInnes 
361ea56f5baSLois Curfman McInnes    Level: advanced
3628cbba510SBarry Smith 
36329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
36429e0b56fSBarry Smith 
36529e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
36629e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
36729e0b56fSBarry Smith @*/
368af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
369329f5518SBarry Smith                             PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
37029e0b56fSBarry Smith {
37129e0b56fSBarry Smith   int        ierr;
37229e0b56fSBarry Smith   Scalar     mone = -1.0;
3736831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
3748f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
37529e0b56fSBarry Smith 
3763a40ed3dSBarry Smith   PetscFunctionBegin;
3778cbba510SBarry Smith   *flag = 0;
37829e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
37929e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3808f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3818f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3828f99978dSLois Curfman McInnes   }
383329f5518SBarry Smith 
384329f5518SBarry Smith   /* don't evaluate function the last time through */
385329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
38629e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
387329f5518SBarry Smith   }
38829e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3893a40ed3dSBarry Smith   PetscFunctionReturn(0);
39029e0b56fSBarry Smith }
39104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
39229e0b56fSBarry Smith #undef __FUNC__
393b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch"
3944b828684SBarry Smith /*@C
395f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3965e76c431SBarry Smith 
397c7afd0dbSLois Curfman McInnes    Collective on SNES
398c7afd0dbSLois Curfman McInnes 
3995e76c431SBarry Smith    Input Parameters:
400c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
401af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4025e76c431SBarry Smith .  x - current iterate
4035e76c431SBarry Smith .  f - residual evaluated at x
4045e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4055e76c431SBarry Smith .  w - work vector
406c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4075e76c431SBarry Smith 
408393d2d9aSLois Curfman McInnes    Output Parameters:
409c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4105e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4115e76c431SBarry Smith .  gnorm - 2-norm of g
4125e76c431SBarry Smith .  ynorm - 2-norm of search length
413c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
414fee21e36SBarry Smith 
415c4a48953SLois Curfman McInnes    Options Database Key:
416c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
417c4a48953SLois Curfman McInnes 
4185e76c431SBarry Smith    Notes:
4195e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4205e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4215e76c431SBarry Smith 
42236851e7fSLois Curfman McInnes    Level: advanced
42336851e7fSLois Curfman McInnes 
42428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
42528ae5a21SLois Curfman McInnes 
426af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4275e76c431SBarry Smith @*/
428af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
429329f5518SBarry Smith                         PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
4305e76c431SBarry Smith {
431406556e6SLois Curfman McInnes   /*
432406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
433406556e6SLois Curfman McInnes      minimization problem:
434406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
435406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
436406556e6SLois Curfman McInnes    */
437406556e6SLois Curfman McInnes 
438329f5518SBarry Smith   PetscReal  steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
439329f5518SBarry Smith   PetscReal  maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
440aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4415e42470aSBarry Smith   Scalar     cinitslope,clambda;
4426b5873e3SBarry Smith #endif
4435e42470aSBarry Smith   int        ierr,count;
4446831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
4455334005bSBarry Smith   Scalar     mone = -1.0,scale;
4468f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
4475e76c431SBarry Smith 
4483a40ed3dSBarry Smith   PetscFunctionBegin;
4497857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
450761aaf1bSLois Curfman McInnes   *flag   = 0;
4515e76c431SBarry Smith   alpha   = neP->alpha;
4525e76c431SBarry Smith   maxstep = neP->maxstep;
4535e76c431SBarry Smith   steptol = neP->steptol;
4545e76c431SBarry Smith 
455cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
456a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
457a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
458a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
459a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
460a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
461ad922ac9SBarry Smith     goto theend1;
46294a9d846SBarry Smith   }
4635e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4645e42470aSBarry Smith     scale = maxstep/(*ynorm);
465aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
46690cbd584SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
4676b5873e3SBarry Smith #else
46890cbd584SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
4696b5873e3SBarry Smith #endif
470393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
4715e76c431SBarry Smith     *ynorm = maxstep;
4725e76c431SBarry Smith   }
4735e76c431SBarry Smith   minlambda = steptol/(*ynorm);
474a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
475aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
476a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
477329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
4785e42470aSBarry Smith #else
479a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
4805e42470aSBarry Smith #endif
4815e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
4825e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
4835e76c431SBarry Smith 
484393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
4855334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
48678b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
487313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
488406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
489393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
49094a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
491ad922ac9SBarry Smith     goto theend1;
4925e76c431SBarry Smith   }
4935e76c431SBarry Smith 
4945e76c431SBarry Smith   /* Fit points with quadratic */
495313b5538SBarry Smith   lambda = 1.0;
496a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4975e76c431SBarry Smith   lambdaprev = lambda;
4985e76c431SBarry Smith   gnormprev = *gnorm;
499329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
500ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
501ddd12767SBarry Smith   else                         lambda = lambdatemp;
502393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
503ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
504aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
505ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5065e42470aSBarry Smith #else
507ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5085e42470aSBarry Smith #endif
50978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
510cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
511406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
512393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
513ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
514ad922ac9SBarry Smith     goto theend1;
5155e76c431SBarry Smith   }
5165e76c431SBarry Smith 
5175e76c431SBarry Smith   /* Fit points with cubic */
5185e76c431SBarry Smith   count = 1;
5195e76c431SBarry Smith   while (1) {
5205e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
5212b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
52290cbd584SBarry Smith       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
523393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
524761aaf1bSLois Curfman McInnes       *flag = -1; break;
5255e76c431SBarry Smith     }
526406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
527406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
528ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5292b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5305e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5315e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5325e76c431SBarry Smith     if (a == 0.0) {
5335e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5345e76c431SBarry Smith     } else {
5355e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5365e76c431SBarry Smith     }
5375e76c431SBarry Smith     lambdaprev = lambda;
5385e76c431SBarry Smith     gnormprev  = *gnorm;
539329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
540329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
5415e76c431SBarry Smith     else                         lambda     = lambdatemp;
542393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
543ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
544aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
545ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
546393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5475e42470aSBarry Smith #else
548ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5495e42470aSBarry Smith #endif
55078b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
551cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
552406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
553393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
554ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
5552715565aSLois Curfman McInnes       goto theend1;
5562b022350SLois Curfman McInnes     } else {
5572b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
5585e76c431SBarry Smith     }
5595e76c431SBarry Smith     count++;
5605e76c431SBarry Smith   }
561ad922ac9SBarry Smith   theend1:
5628f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5638f99978dSLois Curfman McInnes   if (neP->CheckStep) {
5648f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
5658f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
5668f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
5678f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
5688f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
5698f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
5708f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
5718f99978dSLois Curfman McInnes     }
5728f99978dSLois Curfman McInnes   }
5737857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5743a40ed3dSBarry Smith   PetscFunctionReturn(0);
5755e76c431SBarry Smith }
57604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
5775615d1e5SSatish Balay #undef __FUNC__
578b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch"
5794b828684SBarry Smith /*@C
580f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
5815e76c431SBarry Smith 
582c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
583c7afd0dbSLois Curfman McInnes 
5845e42470aSBarry Smith    Input Parameters:
585c7afd0dbSLois Curfman McInnes +  snes - the SNES context
586af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
5875e42470aSBarry Smith .  x - current iterate
5885e42470aSBarry Smith .  f - residual evaluated at x
5895e42470aSBarry Smith .  y - search direction (contains new iterate on output)
5905e42470aSBarry Smith .  w - work vector
591c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
5925e42470aSBarry Smith 
593c4a48953SLois Curfman McInnes    Output Parameters:
594c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
5955e42470aSBarry Smith .  y - new iterate (contains search direction on input)
5965e42470aSBarry Smith .  gnorm - 2-norm of g
5975e42470aSBarry Smith .  ynorm - 2-norm of search length
598c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
599fee21e36SBarry Smith 
600c4a48953SLois Curfman McInnes    Options Database Key:
601c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
602c4a48953SLois Curfman McInnes 
6035e42470aSBarry Smith    Notes:
6046831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
6055e42470aSBarry Smith 
60636851e7fSLois Curfman McInnes    Level: advanced
60736851e7fSLois Curfman McInnes 
608f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
609f59ffdedSLois Curfman McInnes 
610af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6115e42470aSBarry Smith @*/
612af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
613329f5518SBarry Smith                            PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
6145e76c431SBarry Smith {
615406556e6SLois Curfman McInnes   /*
616406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
617406556e6SLois Curfman McInnes      minimization problem:
618406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
619406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
620406556e6SLois Curfman McInnes    */
621329f5518SBarry Smith   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
622aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6235e42470aSBarry Smith   Scalar     cinitslope,clambda;
6246b5873e3SBarry Smith #endif
6255e42470aSBarry Smith   int        ierr,count;
6266831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
6275334005bSBarry Smith   Scalar     mone = -1.0,scale;
6288f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
6295e76c431SBarry Smith 
6303a40ed3dSBarry Smith   PetscFunctionBegin;
6317857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
632761aaf1bSLois Curfman McInnes   *flag   = 0;
6335e42470aSBarry Smith   alpha   = neP->alpha;
6345e42470aSBarry Smith   maxstep = neP->maxstep;
6355e42470aSBarry Smith   steptol = neP->steptol;
6365e76c431SBarry Smith 
6373505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
638b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
63994a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
640b37302c6SLois Curfman McInnes     *gnorm = fnorm;
641b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
642b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
643ad922ac9SBarry Smith     goto theend2;
64494a9d846SBarry Smith   }
6455e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6465e42470aSBarry Smith     scale = maxstep/(*ynorm);
647393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6485e42470aSBarry Smith     *ynorm = maxstep;
6495e76c431SBarry Smith   }
6505e42470aSBarry Smith   minlambda = steptol/(*ynorm);
651a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
652aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
653a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
654329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
6555e42470aSBarry Smith #else
656a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6575e42470aSBarry Smith #endif
6585e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
6595e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6605e42470aSBarry Smith 
661393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6625334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
66378b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
664cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
665406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
666393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
66794a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
668ad922ac9SBarry Smith     goto theend2;
6695e42470aSBarry Smith   }
6705e42470aSBarry Smith 
6715e42470aSBarry Smith   /* Fit points with quadratic */
672313b5538SBarry Smith   lambda = 1.0;
6735e42470aSBarry Smith   count = 1;
6745e42470aSBarry Smith   while (1) {
6755e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
676981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
677981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
678ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
679393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
680761aaf1bSLois Curfman McInnes       *flag = -1; break;
6815e42470aSBarry Smith     }
682a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
683329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
684329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
685329f5518SBarry Smith     else                         lambda     = lambdatemp;
686393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
6873505fcc1SBarry Smith     lambdaneg = -lambda;
688aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6893505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6905e42470aSBarry Smith #else
6913505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6925e42470aSBarry Smith #endif
69378b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
694cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
695406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
696393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
697981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
69806259719SBarry Smith       break;
6995e42470aSBarry Smith     }
7005e42470aSBarry Smith     count++;
7015e42470aSBarry Smith   }
702ad922ac9SBarry Smith   theend2:
7038f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7048f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7058f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7068f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7078f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7088f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7098f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7108f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7118f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7128f99978dSLois Curfman McInnes     }
7138f99978dSLois Curfman McInnes   }
7147857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
7153a40ed3dSBarry Smith   PetscFunctionReturn(0);
7165e76c431SBarry Smith }
71704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7185615d1e5SSatish Balay #undef __FUNC__
719b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch"
720c9e6a524SLois Curfman McInnes /*@C
72177c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7226831982aSBarry Smith    by the method SNESEQLS.
7235e76c431SBarry Smith 
7245e76c431SBarry Smith    Input Parameters:
725c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
726af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
727c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7285e76c431SBarry Smith 
729fee21e36SBarry Smith    Collective on SNES
730fee21e36SBarry Smith 
731c4a48953SLois Curfman McInnes    Available Routines:
732c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
733f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
734f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
735af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7365e76c431SBarry Smith 
737c4a48953SLois Curfman McInnes     Options Database Keys:
738af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
739c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
740c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
741c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
742c4a48953SLois Curfman McInnes 
7435e76c431SBarry Smith    Calling sequence of func:
744bd208895SLois Curfman McInnes .vb
745af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
746329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
747bd208895SLois Curfman McInnes .ve
7485e76c431SBarry Smith 
7495e76c431SBarry Smith     Input parameters for func:
750c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
751af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
7525e76c431SBarry Smith .   x - current iterate
7535e76c431SBarry Smith .   f - residual evaluated at x
7545e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7555e76c431SBarry Smith .   w - work vector
756c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
7575e76c431SBarry Smith 
7585e76c431SBarry Smith     Output parameters for func:
759c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7605e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7615e76c431SBarry Smith .   gnorm - 2-norm of g
7625e76c431SBarry Smith .   ynorm - 2-norm of search length
763c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
764761aaf1bSLois Curfman McInnes            on failure.
765f59ffdedSLois Curfman McInnes 
76636851e7fSLois Curfman McInnes     Level: advanced
76736851e7fSLois Curfman McInnes 
768f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
769f59ffdedSLois Curfman McInnes 
7704ccb76ddSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(), SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
7715e76c431SBarry Smith @*/
772329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
7735e76c431SBarry Smith {
774329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
77582bf6240SBarry Smith 
7763a40ed3dSBarry Smith   PetscFunctionBegin;
77794d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
77882bf6240SBarry Smith   if (f) {
779af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
78082bf6240SBarry Smith   }
7813a40ed3dSBarry Smith   PetscFunctionReturn(0);
7825e76c431SBarry Smith }
78304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
784fb2e594dSBarry Smith EXTERN_C_BEGIN
78582bf6240SBarry Smith #undef __FUNC__
786b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS"
787af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
788af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
78982bf6240SBarry Smith {
79082bf6240SBarry Smith   PetscFunctionBegin;
7916831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
7926831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
79382bf6240SBarry Smith   PetscFunctionReturn(0);
79482bf6240SBarry Smith }
795fb2e594dSBarry Smith EXTERN_C_END
79604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
797c8dd0c0dSLois Curfman McInnes #undef __FUNC__
798b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck"
799c8dd0c0dSLois Curfman McInnes /*@C
800530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8016831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
802c8dd0c0dSLois Curfman McInnes 
803c8dd0c0dSLois Curfman McInnes    Input Parameters:
804c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
805c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
806c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
807c8dd0c0dSLois Curfman McInnes 
808c8dd0c0dSLois Curfman McInnes    Collective on SNES
809c8dd0c0dSLois Curfman McInnes 
810c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
811c8dd0c0dSLois Curfman McInnes .vb
812b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
813c8dd0c0dSLois Curfman McInnes .ve
814b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
815b1ae27eaSLois Curfman McInnes    on failure.
816c8dd0c0dSLois Curfman McInnes 
817c8dd0c0dSLois Curfman McInnes    Input parameters for func:
818c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
819c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8208f99978dSLois Curfman McInnes -  x - current candidate iterate
821c8dd0c0dSLois Curfman McInnes 
822c8dd0c0dSLois Curfman McInnes    Output parameters for func:
823c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
824c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
825c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
826c8dd0c0dSLois Curfman McInnes 
827c8dd0c0dSLois Curfman McInnes    Level: advanced
828c8dd0c0dSLois Curfman McInnes 
8298f99978dSLois Curfman McInnes    Notes:
830b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
831b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
832b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
833b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
834ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8358f99978dSLois Curfman McInnes 
836b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
837b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
838b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
839ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
840ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
841ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
842ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
843b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8448f99978dSLois Curfman McInnes 
845c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
846c8dd0c0dSLois Curfman McInnes 
847c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
848c8dd0c0dSLois Curfman McInnes @*/
8498f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
850c8dd0c0dSLois Curfman McInnes {
8518f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
852c8dd0c0dSLois Curfman McInnes 
853c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
854c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
855c8dd0c0dSLois Curfman McInnes   if (f) {
856c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
857c8dd0c0dSLois Curfman McInnes   }
858c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
859c8dd0c0dSLois Curfman McInnes }
860c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
861c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
862c8dd0c0dSLois Curfman McInnes #undef __FUNC__
863b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS"
8648f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
865c8dd0c0dSLois Curfman McInnes {
866c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
8676831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
8686831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
869c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
870c8dd0c0dSLois Curfman McInnes }
871c8dd0c0dSLois Curfman McInnes EXTERN_C_END
872c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
87304d965bbSLois Curfman McInnes /*
8746831982aSBarry Smith    SNESPrintHelp_EQ_LS - Prints all options for the SNESEQLS method.
87582bf6240SBarry Smith 
87604d965bbSLois Curfman McInnes    Input Parameter:
87704d965bbSLois Curfman McInnes .  snes - the SNES context
87804d965bbSLois Curfman McInnes 
87904d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
88004d965bbSLois Curfman McInnes */
8815615d1e5SSatish Balay #undef __FUNC__
882b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPrintHelp_EQ_LS"
883f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
8845e42470aSBarry Smith {
8856831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
886d132466eSBarry Smith   int     ierr;
8876daaf66cSBarry Smith 
8883a40ed3dSBarry Smith   PetscFunctionBegin;
8896831982aSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm," method SNESEQLS (ls) for systems of nonlinear equations:\n");CHKERRQ(ierr);
890d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);CHKERRQ(ierr);
891d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);CHKERRQ(ierr);
892d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);CHKERRQ(ierr);
893d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);CHKERRQ(ierr);
8943a40ed3dSBarry Smith   PetscFunctionReturn(0);
8955e42470aSBarry Smith }
89604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
89704d965bbSLois Curfman McInnes /*
8986831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
89904d965bbSLois Curfman McInnes 
90004d965bbSLois Curfman McInnes    Input Parameters:
90104d965bbSLois Curfman McInnes .  SNES - the SNES context
90204d965bbSLois Curfman McInnes .  viewer - visualization context
90304d965bbSLois Curfman McInnes 
90404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
90504d965bbSLois Curfman McInnes */
9065615d1e5SSatish Balay #undef __FUNC__
907b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS"
908e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
909a935fc98SLois Curfman McInnes {
9106831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
91119bcc07fSBarry Smith   char       *cstr;
91251695f54SLois Curfman McInnes   int        ierr;
9136831982aSBarry Smith   PetscTruth isascii;
914a935fc98SLois Curfman McInnes 
9153a40ed3dSBarry Smith   PetscFunctionBegin;
9166831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
9170f5bd95cSBarry Smith   if (isascii) {
91819bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
91919bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
92019bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
92119bcc07fSBarry Smith     else                                                cstr = "unknown";
922d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
923d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9245cd90555SBarry Smith   } else {
9250f5bd95cSBarry Smith     SETERRQ1(1,1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
92619bcc07fSBarry Smith   }
9273a40ed3dSBarry Smith   PetscFunctionReturn(0);
928a935fc98SLois Curfman McInnes }
92904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
93004d965bbSLois Curfman McInnes /*
9316831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
93204d965bbSLois Curfman McInnes 
93304d965bbSLois Curfman McInnes    Input Parameter:
93404d965bbSLois Curfman McInnes .  snes - the SNES context
93504d965bbSLois Curfman McInnes 
93604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
93704d965bbSLois Curfman McInnes */
9385615d1e5SSatish Balay #undef __FUNC__
939b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS"
940f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
9415e42470aSBarry Smith {
9426831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
9435e42470aSBarry Smith   char       ver[16];
9445e42470aSBarry Smith   double     tmp;
945f1af5d2fSBarry Smith   int        ierr;
946f1af5d2fSBarry Smith   PetscTruth flg;
9475e42470aSBarry Smith 
9483a40ed3dSBarry Smith   PetscFunctionBegin;
94909d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp,&flg);CHKERRQ(ierr);
95025cdf11fSBarry Smith   if (flg) {
9515e42470aSBarry Smith     ls->alpha = tmp;
9525e42470aSBarry Smith   }
95309d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp,&flg);CHKERRQ(ierr);
95425cdf11fSBarry Smith   if (flg) {
9555e42470aSBarry Smith     ls->maxstep = tmp;
9565e42470aSBarry Smith   }
95709d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp,&flg);CHKERRQ(ierr);
95825cdf11fSBarry Smith   if (flg) {
9595e42470aSBarry Smith     ls->steptol = tmp;
9605e42470aSBarry Smith   }
96109d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16,&flg);CHKERRQ(ierr);
96225cdf11fSBarry Smith   if (flg) {
963c38d4ed2SBarry Smith     PetscTruth isbasic,isnonorms,isquad,iscubic;
9640f5bd95cSBarry Smith 
965c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"basic",&isbasic);CHKERRQ(ierr);
966c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"basicnonorms",&isnonorms);CHKERRQ(ierr);
967c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"quadratic",&isquad);CHKERRQ(ierr);
968c38d4ed2SBarry Smith     ierr = PetscStrcmp(ver,"cubic",&iscubic);CHKERRQ(ierr);
9690f5bd95cSBarry Smith 
9700f5bd95cSBarry Smith     if (isbasic) {
971af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
9720f5bd95cSBarry Smith     } else if (isnonorms) {
973af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
9740f5bd95cSBarry Smith     } else if (isquad) {
975af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
9760f5bd95cSBarry Smith     } else if (iscubic) {
977af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
9785e42470aSBarry Smith     }
979a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
9805e42470aSBarry Smith   }
9813a40ed3dSBarry Smith   PetscFunctionReturn(0);
9825e42470aSBarry Smith }
98304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
98404d965bbSLois Curfman McInnes /*
9856831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
9866831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
98704d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
98804d965bbSLois Curfman McInnes 
98904d965bbSLois Curfman McInnes 
99004d965bbSLois Curfman McInnes    Input Parameter:
99104d965bbSLois Curfman McInnes .  snes - the SNES context
99204d965bbSLois Curfman McInnes 
99304d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
99404d965bbSLois Curfman McInnes  */
995fb2e594dSBarry Smith EXTERN_C_BEGIN
9965615d1e5SSatish Balay #undef __FUNC__
997b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS"
998f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
9995e42470aSBarry Smith {
100082bf6240SBarry Smith   int        ierr;
10016831982aSBarry Smith   SNES_EQ_LS *neP;
10025e42470aSBarry Smith 
10033a40ed3dSBarry Smith   PetscFunctionBegin;
1004a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1005a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
1006a8c6a408SBarry Smith   }
100782bf6240SBarry Smith 
1008f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
1009f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
1010f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
1011f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
1012f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
1013f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
1014f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
10155baf8537SBarry Smith   snes->nwork           = 0;
10165e42470aSBarry Smith 
10176831982aSBarry Smith   neP			= PetscNew(SNES_EQ_LS);CHKPTRQ(neP);
10186831982aSBarry Smith   PLogObjectMemory(snes,sizeof(SNES_EQ_LS));
10195e42470aSBarry Smith   snes->data    	= (void*)neP;
10205e42470aSBarry Smith   neP->alpha		= 1.e-4;
10215e42470aSBarry Smith   neP->maxstep		= 1.e8;
10225e42470aSBarry Smith   neP->steptol		= 1.e-12;
10235e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1024c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1025c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1026c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
102782bf6240SBarry Smith 
1028f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
10290c97f09dSSatish Balay                     SNESSetLineSearch_LS);CHKERRQ(ierr);
1030f1af5d2fSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
10310c97f09dSSatish Balay                     SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
103282bf6240SBarry Smith 
10333a40ed3dSBarry Smith   PetscFunctionReturn(0);
10345e42470aSBarry Smith }
1035fb2e594dSBarry Smith EXTERN_C_END
103682bf6240SBarry Smith 
103782bf6240SBarry Smith 
103882bf6240SBarry Smith 
1039