xref: /petsc/src/snes/impls/ls/ls.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1*29bbc08cSBarry Smith /*$Id: ls.c,v 1.165 2000/09/07 15:17:50 balay 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
836669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
936669109SBarry 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;
1736669109SBarry Smith   PetscTruth hastranspose;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2136669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2236669109SBarry 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);
265d1a10b1SBarry Smith     PLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm);
278a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2836669109SBarry Smith   } else {
2936669109SBarry Smith     Vec       work;
3036669109SBarry Smith     Scalar    result;
3136669109SBarry Smith     PetscReal wnorm;
3236669109SBarry Smith 
3336669109SBarry Smith     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3736669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
3836669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
3936669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
405d1a10b1SBarry Smith     PLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1);
4136669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4236669109SBarry Smith   }
438a5d9ee4SBarry Smith   PetscFunctionReturn(0);
448a5d9ee4SBarry Smith }
458a5d9ee4SBarry Smith 
4674637425SBarry Smith /*
4774637425SBarry Smith      Checks if J^T(F - AX) = 0
4874637425SBarry Smith */
4974637425SBarry Smith #undef __FUNC__
5074637425SBarry Smith #define __FUNC__ /*<a name="SNESLSCheckResidual_Private"></a>*/"SNESLSCheckResidual_Private"
5174637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5274637425SBarry Smith {
5374637425SBarry Smith   PetscReal  a1,a2;
5474637425SBarry Smith   int        ierr;
5574637425SBarry Smith   PetscTruth hastranspose;
5674637425SBarry Smith   Scalar     mone = -1.0;
5774637425SBarry Smith 
5874637425SBarry Smith   PetscFunctionBegin;
5974637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
6074637425SBarry Smith   if (hastranspose) {
6174637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6274637425SBarry Smith     ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr);
6374637425SBarry Smith 
6474637425SBarry Smith     /* Compute || J^T W|| */
6574637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6674637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6774637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
689994e62eSSatish Balay     if (a1 != 0) {
6985471664SBarry Smith       PLogInfo(0,"SNESSolve_EQ_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1);
7085471664SBarry Smith     }
7174637425SBarry Smith   }
7274637425SBarry Smith   PetscFunctionReturn(0);
7374637425SBarry Smith }
7474637425SBarry Smith 
7504d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7604d965bbSLois Curfman McInnes 
7704d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7804d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
7904d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8004d965bbSLois Curfman McInnes      respectively.
8104d965bbSLois Curfman McInnes 
82fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8304d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8404d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
85fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8604d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8704d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
8804d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
8904d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
9004d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9104d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9204d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9304d965bbSLois Curfman McInnes      for all nonlinear solvers.
9404d965bbSLois Curfman McInnes 
9504d965bbSLois Curfman McInnes      Another key routine is:
9604d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9704d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9804d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9904d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10004d965bbSLois Curfman McInnes 
10104d965bbSLois Curfman McInnes      Additional basic routines are:
10204d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10304d965bbSLois Curfman McInnes                                       have actually been used.
10404d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
105186905e3SBarry Smith      SNESView().
10604d965bbSLois Curfman McInnes 
10704d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10804d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10904d965bbSLois Curfman McInnes      above description applies to these categories also.
11004d965bbSLois Curfman McInnes 
11104d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1125e42470aSBarry Smith /*
11304d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
11404d965bbSLois Curfman McInnes    method with a line search.
1155e76c431SBarry Smith 
11604d965bbSLois Curfman McInnes    Input Parameters:
11704d965bbSLois Curfman McInnes .  snes - the SNES context
1185e76c431SBarry Smith 
11904d965bbSLois Curfman McInnes    Output Parameter:
120c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12104d965bbSLois Curfman McInnes 
12204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1235e76c431SBarry Smith 
1245e76c431SBarry Smith    Notes:
1255e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1265e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1275e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1285e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
129393d2d9aSLois Curfman McInnes    and Schnabel.
1305e76c431SBarry Smith */
1315615d1e5SSatish Balay #undef __FUNC__
132b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve_EQ_LS"
133f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
1345e76c431SBarry Smith {
1356831982aSBarry Smith   SNES_EQ_LS          *neP = (SNES_EQ_LS*)snes->data;
13684c569c9SBarry Smith   int                 maxits,i,ierr,lits,lsfail;
137112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
138329f5518SBarry Smith   PetscReal           fnorm,gnorm,xnorm,ynorm;
1395e42470aSBarry Smith   Vec                 Y,X,F,G,W,TMP;
1405e76c431SBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
142184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
143184914b5SBarry Smith 
1445e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1455e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14639e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1475e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1485e42470aSBarry Smith   G		= snes->work[1];
1495e42470aSBarry Smith   W		= snes->work[2];
1505e76c431SBarry Smith 
1514c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1524c49b128SBarry Smith   snes->iter = 0;
1534c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1545334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
155cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1560f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1575e42470aSBarry Smith   snes->norm = fnorm;
1580f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15984c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16094a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1615e76c431SBarry Smith 
162184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16394a9d846SBarry Smith 
164d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
165d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
166d034289bSBarry Smith 
1675e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1685e76c431SBarry Smith 
169ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1705334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1715334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
17278b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
17374637425SBarry Smith 
17485471664SBarry Smith     if (PLogPrintInfo){
17574637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
17685471664SBarry Smith     }
17774637425SBarry Smith 
17890cbd584SBarry Smith     /* should check what happened to the linear solve? */
1793505fcc1SBarry Smith     snes->linear_its += lits;
180af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
181ea4d3ed3SLois Curfman McInnes 
182ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
183ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
184ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
185ea4d3ed3SLois Curfman McInnes     */
18681b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
187af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
188af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
1893505fcc1SBarry Smith     if (lsfail) {
1908a5d9ee4SBarry Smith       PetscTruth ismin;
1913505fcc1SBarry Smith       snes->nfailures++;
1923505fcc1SBarry Smith       snes->reason = SNES_DIVERGED_LS_FAILURE;
1938a5d9ee4SBarry Smith       ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
1948a5d9ee4SBarry Smith       if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
1953505fcc1SBarry Smith       break;
1963505fcc1SBarry Smith     }
1975e76c431SBarry Smith 
19839e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
19939e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
2005e76c431SBarry Smith     fnorm = gnorm;
2015e76c431SBarry Smith 
2020f5bd95cSBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
20384c569c9SBarry Smith     snes->iter = i+1;
2045e42470aSBarry Smith     snes->norm = fnorm;
2050f5bd95cSBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
20684c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
20794a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2085e76c431SBarry Smith 
2095e76c431SBarry Smith     /* Test for convergence */
21029e0b56fSBarry Smith     if (snes->converged) {
21129e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2123505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2133505fcc1SBarry Smith       if (snes->reason) {
21416c95cacSBarry Smith         break;
21516c95cacSBarry Smith       }
21616c95cacSBarry Smith     }
21729e0b56fSBarry Smith   }
21839e2f89bSBarry Smith   if (X != snes->vec_sol) {
219393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
220e15f7bb5SBarry Smith   }
221e15f7bb5SBarry Smith   if (F != snes->vec_func) {
222e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
223e15f7bb5SBarry Smith   }
22439e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
22539e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
22652392280SLois Curfman McInnes   if (i == maxits) {
227981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
22852392280SLois Curfman McInnes     i--;
2293505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
23052392280SLois Curfman McInnes   }
2310f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2320f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2335e42470aSBarry Smith   *outits = i+1;
2343a40ed3dSBarry Smith   PetscFunctionReturn(0);
2355e76c431SBarry Smith }
23604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
23704d965bbSLois Curfman McInnes /*
23804d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
2396831982aSBarry Smith    of the SNESEQLS nonlinear solver.
24004d965bbSLois Curfman McInnes 
24104d965bbSLois Curfman McInnes    Input Parameter:
24204d965bbSLois Curfman McInnes .  snes - the SNES context
24304d965bbSLois Curfman McInnes .  x - the solution vector
24404d965bbSLois Curfman McInnes 
24504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
24604d965bbSLois Curfman McInnes 
24704d965bbSLois Curfman McInnes    Notes:
24804d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
24904d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
25004d965bbSLois Curfman McInnes    the call to SNESSolve().
25104d965bbSLois Curfman McInnes  */
2525615d1e5SSatish Balay #undef __FUNC__
253b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp_EQ_LS"
254f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
2555e76c431SBarry Smith {
2565e42470aSBarry Smith   int ierr;
2573a40ed3dSBarry Smith 
2583a40ed3dSBarry Smith   PetscFunctionBegin;
25981b6cf68SLois Curfman McInnes   snes->nwork = 4;
260d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
2615e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
26281b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2633a40ed3dSBarry Smith   PetscFunctionReturn(0);
2645e76c431SBarry Smith }
26504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
26604d965bbSLois Curfman McInnes /*
2676831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
26804d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
26904d965bbSLois Curfman McInnes 
27004d965bbSLois Curfman McInnes    Input Parameter:
27104d965bbSLois Curfman McInnes .  snes - the SNES context
27204d965bbSLois Curfman McInnes 
273de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
27404d965bbSLois Curfman McInnes  */
2755615d1e5SSatish Balay #undef __FUNC__
276b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy_EQ_LS"
277e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
2785e76c431SBarry Smith {
279393d2d9aSLois Curfman McInnes   int  ierr;
2803a40ed3dSBarry Smith 
2813a40ed3dSBarry Smith   PetscFunctionBegin;
2825baf8537SBarry Smith   if (snes->nwork) {
2834b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
28421c89e3eSBarry Smith   }
2855c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2863a40ed3dSBarry Smith   PetscFunctionReturn(0);
2875e76c431SBarry Smith }
28804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2895615d1e5SSatish Balay #undef __FUNC__
290b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearch"
29182bf6240SBarry Smith 
2924b828684SBarry Smith /*@C
2935e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2945e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2955e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2965e76c431SBarry Smith 
297c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
298c7afd0dbSLois Curfman McInnes 
2995e76c431SBarry Smith    Input Parameters:
300c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
301af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3025e76c431SBarry Smith .  x - current iterate
3035e76c431SBarry Smith .  f - residual evaluated at x
3045e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3055e76c431SBarry Smith .  w - work vector
306c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3075e76c431SBarry Smith 
308c4a48953SLois Curfman McInnes    Output Parameters:
309c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3105e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3115e76c431SBarry Smith .  gnorm - 2-norm of g
3125e76c431SBarry Smith .  ynorm - 2-norm of search length
313c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
314fee21e36SBarry Smith 
315c4a48953SLois Curfman McInnes    Options Database Key:
316c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
317c4a48953SLois Curfman McInnes 
31836851e7fSLois Curfman McInnes    Level: advanced
31936851e7fSLois Curfman McInnes 
32028ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
32128ae5a21SLois Curfman McInnes 
322f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
323af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3245e76c431SBarry Smith @*/
3255d1a10b1SBarry Smith int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
3265e76c431SBarry Smith {
3275e42470aSBarry Smith   int        ierr;
3285334005bSBarry Smith   Scalar     mone = -1.0;
3296831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
3308f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
3315334005bSBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
333761aaf1bSLois Curfman McInnes   *flag = 0;
33481bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
335a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
336ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3378f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3388f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3398f99978dSLois Curfman McInnes   }
340ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
341a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
34281bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3433a40ed3dSBarry Smith   PetscFunctionReturn(0);
3445e76c431SBarry Smith }
34504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3465615d1e5SSatish Balay #undef __FUNC__
347b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESNoLineSearchNoNorms"
34882bf6240SBarry Smith 
34929e0b56fSBarry Smith /*@C
35029e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
35129e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
35229e0b56fSBarry Smith    even compute the norm of the function or search direction; this
35329e0b56fSBarry Smith    is intended only when you know the full step is fine and are
35429e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
355c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
356c7afd0dbSLois Curfman McInnes 
357c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
35829e0b56fSBarry Smith 
35929e0b56fSBarry Smith    Input Parameters:
360c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
361af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
36229e0b56fSBarry Smith .  x - current iterate
36329e0b56fSBarry Smith .  f - residual evaluated at x
36429e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
36529e0b56fSBarry Smith .  w - work vector
366c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
36729e0b56fSBarry Smith 
36829e0b56fSBarry Smith    Output Parameters:
369c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
37029e0b56fSBarry Smith .  gnorm - not changed
37129e0b56fSBarry Smith .  ynorm - not changed
372c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
373fee21e36SBarry Smith 
37429e0b56fSBarry Smith    Options Database Key:
375c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
37629e0b56fSBarry Smith 
3778cbba510SBarry Smith    Notes:
378ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
379ea56f5baSLois Curfman McInnes    either the options
380ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
381ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
382329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
383329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
384329f5518SBarry Smith 
385329f5518SBarry Smith    During the final iteration this will not evaluate the function at
386329f5518SBarry Smith    the solution point. This is to save a function evaluation while
387329f5518SBarry Smith    using pseudo-timestepping.
3888cbba510SBarry Smith 
389ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
390ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
391ea56f5baSLois Curfman McInnes    correct, since they are not computed.
392ea56f5baSLois Curfman McInnes 
393ea56f5baSLois Curfman McInnes    Level: advanced
3948cbba510SBarry Smith 
39529e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
39629e0b56fSBarry Smith 
39729e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
39829e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
39929e0b56fSBarry Smith @*/
4005d1a10b1SBarry Smith int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
40129e0b56fSBarry Smith {
40229e0b56fSBarry Smith   int        ierr;
40329e0b56fSBarry Smith   Scalar     mone = -1.0;
4046831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
4058f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
40629e0b56fSBarry Smith 
4073a40ed3dSBarry Smith   PetscFunctionBegin;
4088cbba510SBarry Smith   *flag = 0;
40981bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
41029e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4118f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4128f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4138f99978dSLois Curfman McInnes   }
414329f5518SBarry Smith 
415329f5518SBarry Smith   /* don't evaluate function the last time through */
416329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
41729e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
418329f5518SBarry Smith   }
41981bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4203a40ed3dSBarry Smith   PetscFunctionReturn(0);
42129e0b56fSBarry Smith }
42204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
42329e0b56fSBarry Smith #undef __FUNC__
424b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCubicLineSearch"
4254b828684SBarry Smith /*@C
426f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4275e76c431SBarry Smith 
428c7afd0dbSLois Curfman McInnes    Collective on SNES
429c7afd0dbSLois Curfman McInnes 
4305e76c431SBarry Smith    Input Parameters:
431c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
432af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4335e76c431SBarry Smith .  x - current iterate
4345e76c431SBarry Smith .  f - residual evaluated at x
4355e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4365e76c431SBarry Smith .  w - work vector
437c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4385e76c431SBarry Smith 
439393d2d9aSLois Curfman McInnes    Output Parameters:
440c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4415e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4425e76c431SBarry Smith .  gnorm - 2-norm of g
4435e76c431SBarry Smith .  ynorm - 2-norm of search length
444c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
445fee21e36SBarry Smith 
446c4a48953SLois Curfman McInnes    Options Database Key:
447c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
448c4a48953SLois Curfman McInnes 
4495e76c431SBarry Smith    Notes:
4505e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4515e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4525e76c431SBarry Smith 
45336851e7fSLois Curfman McInnes    Level: advanced
45436851e7fSLois Curfman McInnes 
45528ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
45628ae5a21SLois Curfman McInnes 
457af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4585e76c431SBarry Smith @*/
4595d1a10b1SBarry Smith int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
4605e76c431SBarry Smith {
461406556e6SLois Curfman McInnes   /*
462406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
463406556e6SLois Curfman McInnes      minimization problem:
464406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
465406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
466406556e6SLois Curfman McInnes    */
467406556e6SLois Curfman McInnes 
468329f5518SBarry Smith   PetscReal  steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
469329f5518SBarry Smith   PetscReal  maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
470aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
4715e42470aSBarry Smith   Scalar     cinitslope,clambda;
4726b5873e3SBarry Smith #endif
4735e42470aSBarry Smith   int        ierr,count;
4746831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
4755334005bSBarry Smith   Scalar     mone = -1.0,scale;
4768f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
4775e76c431SBarry Smith 
4783a40ed3dSBarry Smith   PetscFunctionBegin;
47981bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
480761aaf1bSLois Curfman McInnes   *flag   = 0;
4815e76c431SBarry Smith   alpha   = neP->alpha;
4825e76c431SBarry Smith   maxstep = neP->maxstep;
4835e76c431SBarry Smith   steptol = neP->steptol;
4845e76c431SBarry Smith 
485cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
486a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
487a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
488a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
489a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
490a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
491ad922ac9SBarry Smith     goto theend1;
49294a9d846SBarry Smith   }
4935e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
4945e42470aSBarry Smith     scale = maxstep/(*ynorm);
495aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
49690cbd584SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
4976b5873e3SBarry Smith #else
49890cbd584SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
4996b5873e3SBarry Smith #endif
500393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5015e76c431SBarry Smith     *ynorm = maxstep;
5025e76c431SBarry Smith   }
5035e76c431SBarry Smith   minlambda = steptol/(*ynorm);
504a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
505aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
506a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
507329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5085e42470aSBarry Smith #else
509a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5105e42470aSBarry Smith #endif
5115e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5125e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5135e76c431SBarry Smith 
514393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5155334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
51678b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
517313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5185d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
519393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
52094a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
521ad922ac9SBarry Smith     goto theend1;
5225e76c431SBarry Smith   }
5235e76c431SBarry Smith 
5245e76c431SBarry Smith   /* Fit points with quadratic */
525313b5538SBarry Smith   lambda = 1.0;
526a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5275e76c431SBarry Smith   lambdaprev = lambda;
5285e76c431SBarry Smith   gnormprev = *gnorm;
529329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
530ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
531ddd12767SBarry Smith   else                         lambda = lambdatemp;
532393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
533ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
534aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
535ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5365e42470aSBarry Smith #else
537ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5385e42470aSBarry Smith #endif
53978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
540cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5415d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
542393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
543ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
544ad922ac9SBarry Smith     goto theend1;
5455e76c431SBarry Smith   }
5465e76c431SBarry Smith 
5475e76c431SBarry Smith   /* Fit points with cubic */
5485e76c431SBarry Smith   count = 1;
5495e76c431SBarry Smith   while (1) {
5505e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
5512b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
55290cbd584SBarry Smith       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
553393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
554761aaf1bSLois Curfman McInnes       *flag = -1; break;
5555e76c431SBarry Smith     }
5565d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
5575d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
558ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5592b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5605e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5615e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5625e76c431SBarry Smith     if (a == 0.0) {
5635e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5645e76c431SBarry Smith     } else {
5655e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5665e76c431SBarry Smith     }
5675e76c431SBarry Smith     lambdaprev = lambda;
5685e76c431SBarry Smith     gnormprev  = *gnorm;
569329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
570329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
5715e76c431SBarry Smith     else                         lambda     = lambdatemp;
572393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
573ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
574aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
575ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
576393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5775e42470aSBarry Smith #else
578ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5795e42470aSBarry Smith #endif
58078b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
581cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5825d1a10b1SBarry Smith     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
583393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
584ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
5852715565aSLois Curfman McInnes       goto theend1;
5862b022350SLois Curfman McInnes     } else {
5872b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
5885e76c431SBarry Smith     }
5895e76c431SBarry Smith     count++;
5905e76c431SBarry Smith   }
591ad922ac9SBarry Smith   theend1:
5928f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5938f99978dSLois Curfman McInnes   if (neP->CheckStep) {
5948f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
5958f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
5968f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
5978f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
5988f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
5998f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6008f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6018f99978dSLois Curfman McInnes     }
6028f99978dSLois Curfman McInnes   }
60381bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6043a40ed3dSBarry Smith   PetscFunctionReturn(0);
6055e76c431SBarry Smith }
60604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6075615d1e5SSatish Balay #undef __FUNC__
608b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESQuadraticLineSearch"
6094b828684SBarry Smith /*@C
610f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6115e76c431SBarry Smith 
612c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
613c7afd0dbSLois Curfman McInnes 
6145e42470aSBarry Smith    Input Parameters:
615c7afd0dbSLois Curfman McInnes +  snes - the SNES context
616af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6175e42470aSBarry Smith .  x - current iterate
6185e42470aSBarry Smith .  f - residual evaluated at x
6195e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6205e42470aSBarry Smith .  w - work vector
621c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6225e42470aSBarry Smith 
623c4a48953SLois Curfman McInnes    Output Parameters:
624c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6255e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6265e42470aSBarry Smith .  gnorm - 2-norm of g
6275e42470aSBarry Smith .  ynorm - 2-norm of search length
628c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
629fee21e36SBarry Smith 
630c4a48953SLois Curfman McInnes    Options Database Key:
631c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
632c4a48953SLois Curfman McInnes 
6335e42470aSBarry Smith    Notes:
6346831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
6355e42470aSBarry Smith 
63636851e7fSLois Curfman McInnes    Level: advanced
63736851e7fSLois Curfman McInnes 
638f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
639f59ffdedSLois Curfman McInnes 
640af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6415e42470aSBarry Smith @*/
6425d1a10b1SBarry Smith int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
6435e76c431SBarry Smith {
644406556e6SLois Curfman McInnes   /*
645406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
646406556e6SLois Curfman McInnes      minimization problem:
647406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
648406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
649406556e6SLois Curfman McInnes    */
650329f5518SBarry Smith   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
651aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
6525e42470aSBarry Smith   Scalar     cinitslope,clambda;
6536b5873e3SBarry Smith #endif
6545e42470aSBarry Smith   int        ierr,count;
6556831982aSBarry Smith   SNES_EQ_LS *neP = (SNES_EQ_LS*)snes->data;
6565334005bSBarry Smith   Scalar     mone = -1.0,scale;
6578f99978dSLois Curfman McInnes   PetscTruth change_y = PETSC_FALSE;
6585e76c431SBarry Smith 
6593a40ed3dSBarry Smith   PetscFunctionBegin;
66081bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
661761aaf1bSLois Curfman McInnes   *flag   = 0;
6625e42470aSBarry Smith   alpha   = neP->alpha;
6635e42470aSBarry Smith   maxstep = neP->maxstep;
6645e42470aSBarry Smith   steptol = neP->steptol;
6655e76c431SBarry Smith 
6663505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
667b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
66894a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
669b37302c6SLois Curfman McInnes     *gnorm = fnorm;
670b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
671b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
672ad922ac9SBarry Smith     goto theend2;
67394a9d846SBarry Smith   }
6745e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6755e42470aSBarry Smith     scale = maxstep/(*ynorm);
676393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6775e42470aSBarry Smith     *ynorm = maxstep;
6785e76c431SBarry Smith   }
6795e42470aSBarry Smith   minlambda = steptol/(*ynorm);
680a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
681aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
682a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
683329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
6845e42470aSBarry Smith #else
685a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6865e42470aSBarry Smith #endif
6875e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
6885e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6895e42470aSBarry Smith 
690393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6915334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
69278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
693cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
6945d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
695393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
69694a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
697ad922ac9SBarry Smith     goto theend2;
6985e42470aSBarry Smith   }
6995e42470aSBarry Smith 
7005e42470aSBarry Smith   /* Fit points with quadratic */
701313b5538SBarry Smith   lambda = 1.0;
7025e42470aSBarry Smith   count = 1;
7035e42470aSBarry Smith   while (1) {
7045e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
705981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
7065d1a10b1SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
707393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
708761aaf1bSLois Curfman McInnes       *flag = -1; break;
7095e42470aSBarry Smith     }
710a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
711329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
712329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
713329f5518SBarry Smith     else                         lambda     = lambdatemp;
714393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7153505fcc1SBarry Smith     lambdaneg = -lambda;
716aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7173505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7185e42470aSBarry Smith #else
7193505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7205e42470aSBarry Smith #endif
72178b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
722cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7235d1a10b1SBarry Smith     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
724393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
725981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
72606259719SBarry Smith       break;
7275e42470aSBarry Smith     }
7285e42470aSBarry Smith     count++;
7295e42470aSBarry Smith   }
730ad922ac9SBarry Smith   theend2:
7318f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7328f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7338f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7348f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7358f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7368f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7378f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7388f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7398f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7408f99978dSLois Curfman McInnes     }
7418f99978dSLois Curfman McInnes   }
74281bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7433a40ed3dSBarry Smith   PetscFunctionReturn(0);
7445e76c431SBarry Smith }
74504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7465615d1e5SSatish Balay #undef __FUNC__
747b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch"
748c9e6a524SLois Curfman McInnes /*@C
74977c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7506831982aSBarry Smith    by the method SNESEQLS.
7515e76c431SBarry Smith 
7525e76c431SBarry Smith    Input Parameters:
753c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
754af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
755c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7565e76c431SBarry Smith 
757fee21e36SBarry Smith    Collective on SNES
758fee21e36SBarry Smith 
759c4a48953SLois Curfman McInnes    Available Routines:
760c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
761f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
762f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
763af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7645e76c431SBarry Smith 
765c4a48953SLois Curfman McInnes     Options Database Keys:
766af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
767c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
768c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
769c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
770c4a48953SLois Curfman McInnes 
7715e76c431SBarry Smith    Calling sequence of func:
772bd208895SLois Curfman McInnes .vb
773af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
774329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
775bd208895SLois Curfman McInnes .ve
7765e76c431SBarry Smith 
7775e76c431SBarry Smith     Input parameters for func:
778c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
779af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
7805e76c431SBarry Smith .   x - current iterate
7815e76c431SBarry Smith .   f - residual evaluated at x
7825e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7835e76c431SBarry Smith .   w - work vector
784c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
7855e76c431SBarry Smith 
7865e76c431SBarry Smith     Output parameters for func:
787c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7885e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7895e76c431SBarry Smith .   gnorm - 2-norm of g
7905e76c431SBarry Smith .   ynorm - 2-norm of search length
791c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
792761aaf1bSLois Curfman McInnes            on failure.
793f59ffdedSLois Curfman McInnes 
79436851e7fSLois Curfman McInnes     Level: advanced
79536851e7fSLois Curfman McInnes 
796f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
797f59ffdedSLois Curfman McInnes 
7985d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
7995d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8005e76c431SBarry Smith @*/
801329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
8025e76c431SBarry Smith {
803329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
80482bf6240SBarry Smith 
8053a40ed3dSBarry Smith   PetscFunctionBegin;
80694d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
80782bf6240SBarry Smith   if (f) {
808af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
80982bf6240SBarry Smith   }
8103a40ed3dSBarry Smith   PetscFunctionReturn(0);
8115e76c431SBarry Smith }
81204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
813fb2e594dSBarry Smith EXTERN_C_BEGIN
81482bf6240SBarry Smith #undef __FUNC__
815b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearch_LS"
816af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
817af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
81882bf6240SBarry Smith {
81982bf6240SBarry Smith   PetscFunctionBegin;
8206831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
8216831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
82282bf6240SBarry Smith   PetscFunctionReturn(0);
82382bf6240SBarry Smith }
824fb2e594dSBarry Smith EXTERN_C_END
82504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
826c8dd0c0dSLois Curfman McInnes #undef __FUNC__
827b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck"
828c8dd0c0dSLois Curfman McInnes /*@C
829530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8306831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
831c8dd0c0dSLois Curfman McInnes 
832c8dd0c0dSLois Curfman McInnes    Input Parameters:
833c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
834c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
835c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
836c8dd0c0dSLois Curfman McInnes 
837c8dd0c0dSLois Curfman McInnes    Collective on SNES
838c8dd0c0dSLois Curfman McInnes 
839c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
840c8dd0c0dSLois Curfman McInnes .vb
841b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
842c8dd0c0dSLois Curfman McInnes .ve
843b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
844b1ae27eaSLois Curfman McInnes    on failure.
845c8dd0c0dSLois Curfman McInnes 
846c8dd0c0dSLois Curfman McInnes    Input parameters for func:
847c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
848c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8498f99978dSLois Curfman McInnes -  x - current candidate iterate
850c8dd0c0dSLois Curfman McInnes 
851c8dd0c0dSLois Curfman McInnes    Output parameters for func:
852c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
853c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
854c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
855c8dd0c0dSLois Curfman McInnes 
856c8dd0c0dSLois Curfman McInnes    Level: advanced
857c8dd0c0dSLois Curfman McInnes 
8588f99978dSLois Curfman McInnes    Notes:
859b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
860b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
861b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
862b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
863ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8648f99978dSLois Curfman McInnes 
865b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
866b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
867b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
868ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
869ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
870ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
871ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
872b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8738f99978dSLois Curfman McInnes 
874c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
875c8dd0c0dSLois Curfman McInnes 
876c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
877c8dd0c0dSLois Curfman McInnes @*/
8788f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
879c8dd0c0dSLois Curfman McInnes {
8808f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
881c8dd0c0dSLois Curfman McInnes 
882c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
883c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
884c8dd0c0dSLois Curfman McInnes   if (f) {
885c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
886c8dd0c0dSLois Curfman McInnes   }
887c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
888c8dd0c0dSLois Curfman McInnes }
889c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
890c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
891c8dd0c0dSLois Curfman McInnes #undef __FUNC__
892b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetLineSearchCheck_LS"
8938f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
894c8dd0c0dSLois Curfman McInnes {
895c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
8966831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
8976831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
898c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
899c8dd0c0dSLois Curfman McInnes }
900c8dd0c0dSLois Curfman McInnes EXTERN_C_END
901c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
90204d965bbSLois Curfman McInnes /*
9036831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
90404d965bbSLois Curfman McInnes 
90504d965bbSLois Curfman McInnes    Input Parameters:
90604d965bbSLois Curfman McInnes .  SNES - the SNES context
90704d965bbSLois Curfman McInnes .  viewer - visualization context
90804d965bbSLois Curfman McInnes 
90904d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
91004d965bbSLois Curfman McInnes */
9115615d1e5SSatish Balay #undef __FUNC__
912b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView_EQ_LS"
913e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
914a935fc98SLois Curfman McInnes {
9156831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
91619bcc07fSBarry Smith   char       *cstr;
91751695f54SLois Curfman McInnes   int        ierr;
9186831982aSBarry Smith   PetscTruth isascii;
919a935fc98SLois Curfman McInnes 
9203a40ed3dSBarry Smith   PetscFunctionBegin;
9216831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
9220f5bd95cSBarry Smith   if (isascii) {
92319bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
92419bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
92519bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
92619bcc07fSBarry Smith     else                                                cstr = "unknown";
927d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
928d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9295cd90555SBarry Smith   } else {
930*29bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
93119bcc07fSBarry Smith   }
9323a40ed3dSBarry Smith   PetscFunctionReturn(0);
933a935fc98SLois Curfman McInnes }
93404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
93504d965bbSLois Curfman McInnes /*
9366831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
93704d965bbSLois Curfman McInnes 
93804d965bbSLois Curfman McInnes    Input Parameter:
93904d965bbSLois Curfman McInnes .  snes - the SNES context
94004d965bbSLois Curfman McInnes 
94104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
94204d965bbSLois Curfman McInnes */
9435615d1e5SSatish Balay #undef __FUNC__
944b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions_EQ_LS"
945f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
9465e42470aSBarry Smith {
9476831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
948186905e3SBarry Smith   char       ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
949f1af5d2fSBarry Smith   int        ierr;
950f1af5d2fSBarry Smith   PetscTruth flg;
9515e42470aSBarry Smith 
9523a40ed3dSBarry Smith   PetscFunctionBegin;
9534bbc92c1SBarry Smith   ierr = OptionsHead("SNES Line search options");CHKERRQ(ierr);
954186905e3SBarry Smith     ierr = OptionsDouble("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
955186905e3SBarry Smith     ierr = OptionsDouble("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
956186905e3SBarry Smith     ierr = OptionsDouble("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
957186905e3SBarry Smith 
958186905e3SBarry Smith     ierr = OptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);CHKERRQ(ierr);
95925cdf11fSBarry Smith     if (flg) {
960c38d4ed2SBarry Smith       PetscTruth isbasic,isnonorms,isquad,iscubic;
9610f5bd95cSBarry Smith 
962186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[0],&isbasic);CHKERRQ(ierr);
963186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[1],&isnonorms);CHKERRQ(ierr);
964186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[2],&isquad);CHKERRQ(ierr);
965186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[3],&iscubic);CHKERRQ(ierr);
9660f5bd95cSBarry Smith 
9670f5bd95cSBarry Smith       if (isbasic) {
968af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
9690f5bd95cSBarry Smith       } else if (isnonorms) {
970af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
9710f5bd95cSBarry Smith       } else if (isquad) {
972af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
9730f5bd95cSBarry Smith       } else if (iscubic) {
974af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
9755e42470aSBarry Smith       }
976*29bbc08cSBarry Smith       else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
9775e42470aSBarry Smith     }
9784bbc92c1SBarry Smith   ierr = OptionsTail();CHKERRQ(ierr);
9793a40ed3dSBarry Smith   PetscFunctionReturn(0);
9805e42470aSBarry Smith }
98104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
98204d965bbSLois Curfman McInnes /*
9836831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
9846831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
98504d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
98604d965bbSLois Curfman McInnes 
98704d965bbSLois Curfman McInnes 
98804d965bbSLois Curfman McInnes    Input Parameter:
98904d965bbSLois Curfman McInnes .  snes - the SNES context
99004d965bbSLois Curfman McInnes 
99104d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
99204d965bbSLois Curfman McInnes  */
993fb2e594dSBarry Smith EXTERN_C_BEGIN
9945615d1e5SSatish Balay #undef __FUNC__
995b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate_EQ_LS"
996f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
9975e42470aSBarry Smith {
99882bf6240SBarry Smith   int        ierr;
9996831982aSBarry Smith   SNES_EQ_LS *neP;
10005e42470aSBarry Smith 
10013a40ed3dSBarry Smith   PetscFunctionBegin;
1002a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1003*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1004a8c6a408SBarry Smith   }
100582bf6240SBarry Smith 
1006f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
1007f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
1008f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
1009f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
1010f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
1011f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
10125baf8537SBarry Smith   snes->nwork           = 0;
10135e42470aSBarry Smith 
10146831982aSBarry Smith   neP			= PetscNew(SNES_EQ_LS);CHKPTRQ(neP);
10156831982aSBarry Smith   PLogObjectMemory(snes,sizeof(SNES_EQ_LS));
10165e42470aSBarry Smith   snes->data    	= (void*)neP;
10175e42470aSBarry Smith   neP->alpha		= 1.e-4;
10185e42470aSBarry Smith   neP->maxstep		= 1.e8;
10195e42470aSBarry Smith   neP->steptol		= 1.e-12;
10205e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1021c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1022c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1023c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
102482bf6240SBarry Smith 
10255d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
10265d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
102782bf6240SBarry Smith 
10283a40ed3dSBarry Smith   PetscFunctionReturn(0);
10295e42470aSBarry Smith }
1030fb2e594dSBarry Smith EXTERN_C_END
103182bf6240SBarry Smith 
103282bf6240SBarry Smith 
103382bf6240SBarry Smith 
1034