xref: /petsc/src/snes/impls/ls/ls.c (revision 3c632250c5290da0280719623b052aa271290270)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
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 */
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
13dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
148a5d9ee4SBarry Smith {
158a5d9ee4SBarry Smith   PetscReal      a1;
16dfbe8321SBarry Smith   PetscErrorCode 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);
2663ba0a88SBarry Smith     ierr = PetscLogInfo((0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm));CHKERRQ(ierr);
278a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2836669109SBarry Smith   } else {
2936669109SBarry Smith     Vec         work;
30ea709b57SSatish Balay     PetscScalar 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);
4063ba0a88SBarry Smith     ierr = PetscLogInfo((0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1));CHKERRQ(ierr);
4136669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4236669109SBarry Smith   }
438a5d9ee4SBarry Smith   PetscFunctionReturn(0);
448a5d9ee4SBarry Smith }
458a5d9ee4SBarry Smith 
4674637425SBarry Smith /*
475ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4874637425SBarry Smith */
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
51dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5274637425SBarry Smith {
5374637425SBarry Smith   PetscReal      a1,a2;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
5574637425SBarry Smith   PetscTruth     hastranspose;
56ea709b57SSatish Balay   PetscScalar    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);
689e247f21SBarry Smith     if (a1) {
6963ba0a88SBarry Smith       ierr = PetscLogInfo((0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1));CHKERRQ(ierr);
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,
7894b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, 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
884b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
894b27c08aSLois Curfman McInnes      systems of nonlinear equations 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 /*
1134b27c08aSLois Curfman McInnes    SNESSolve_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 */
1314a2ae208SSatish Balay #undef __FUNCT__
1324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
133dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1345e76c431SBarry Smith {
1354b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
1366849ba73SBarry Smith   PetscErrorCode ierr;
137a7cc72afSBarry Smith   PetscInt       maxits,i,lits;
138a7cc72afSBarry Smith   PetscTruth     lssucceed;
139112a2221SBarry Smith   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
140329f5518SBarry Smith   PetscReal      fnorm,gnorm,xnorm,ynorm;
1415e42470aSBarry Smith   Vec            Y,X,F,G,W,TMP;
142c293cc10SBarry Smith   KSP            ksp;
1435e76c431SBarry Smith 
1443a40ed3dSBarry Smith   PetscFunctionBegin;
14594b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
146184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
147184914b5SBarry Smith 
1485e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1495e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
15039e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1515e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1525e42470aSBarry Smith   G		= snes->work[1];
1535e42470aSBarry Smith   W		= snes->work[2];
1545e76c431SBarry Smith 
1554c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1564c49b128SBarry Smith   snes->iter = 0;
1574c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15819717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
159cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
16043e71028SBarry Smith   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1610f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1625e42470aSBarry Smith   snes->norm = fnorm;
1630f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16484c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16594a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1665e76c431SBarry Smith 
16770441072SBarry Smith   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16894a9d846SBarry Smith 
169d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
170d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
171d034289bSBarry Smith 
1725e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1735e76c431SBarry Smith 
174329e7e40SMatthew Knepley     /* Call general purpose update function */
175abc0a331SBarry Smith     if (snes->update) {
176329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
177de8cb200SMatthew Knepley     }
178329e7e40SMatthew Knepley 
179ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1805334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18194b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18223ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
183c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
18474637425SBarry Smith 
185b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
18674637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
18785471664SBarry Smith     }
18874637425SBarry Smith 
18990cbd584SBarry Smith     /* should check what happened to the linear solve? */
1903505fcc1SBarry Smith     snes->linear_its += lits;
19163ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr);
192ea4d3ed3SLois Curfman McInnes 
193ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
194ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
195ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
196ea4d3ed3SLois Curfman McInnes     */
19781b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
198a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
19963ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed));CHKERRQ(ierr);
20043e71028SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
201a3b891d8SBarry Smith 
202a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
203*3c632250SBarry Smith     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
204a3b891d8SBarry Smith     fnorm = gnorm;
205a3b891d8SBarry Smith 
2065ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2075ed2d596SBarry Smith     snes->iter = i+1;
2085ed2d596SBarry Smith     snes->norm = fnorm;
2095ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2105ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2115ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2125ed2d596SBarry Smith 
213a7cc72afSBarry Smith     if (!lssucceed) {
2148a5d9ee4SBarry Smith       PetscTruth ismin;
21550ffb88aSMatthew Knepley 
21650ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2173505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2188a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2198a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2203505fcc1SBarry Smith         break;
2213505fcc1SBarry Smith       }
22250ffb88aSMatthew Knepley     }
2235e76c431SBarry Smith 
2245e76c431SBarry Smith     /* Test for convergence */
22529e0b56fSBarry Smith     if (snes->converged) {
22629e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2273505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2283505fcc1SBarry Smith       if (snes->reason) {
22916c95cacSBarry Smith         break;
23016c95cacSBarry Smith       }
23116c95cacSBarry Smith     }
23229e0b56fSBarry Smith   }
23339e2f89bSBarry Smith   if (X != snes->vec_sol) {
234393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
235e15f7bb5SBarry Smith   }
236e15f7bb5SBarry Smith   if (F != snes->vec_func) {
237e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
238e15f7bb5SBarry Smith   }
23939e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
24039e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
24152392280SLois Curfman McInnes   if (i == maxits) {
24263ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr);
2433505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
24452392280SLois Curfman McInnes   }
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
2465e76c431SBarry Smith }
24704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
24804d965bbSLois Curfman McInnes /*
2494b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2504b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25104d965bbSLois Curfman McInnes 
25204d965bbSLois Curfman McInnes    Input Parameter:
25304d965bbSLois Curfman McInnes .  snes - the SNES context
25404d965bbSLois Curfman McInnes .  x - the solution vector
25504d965bbSLois Curfman McInnes 
25604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25704d965bbSLois Curfman McInnes 
25804d965bbSLois Curfman McInnes    Notes:
2594b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26004d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26104d965bbSLois Curfman McInnes    the call to SNESSolve().
26204d965bbSLois Curfman McInnes  */
2634a2ae208SSatish Balay #undef __FUNCT__
2644b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
265dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2665e76c431SBarry Smith {
267dfbe8321SBarry Smith   PetscErrorCode ierr;
2683a40ed3dSBarry Smith 
2693a40ed3dSBarry Smith   PetscFunctionBegin;
27081b6cf68SLois Curfman McInnes   snes->nwork = 4;
271d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
272efee365bSSatish Balay   ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
27381b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2743a40ed3dSBarry Smith   PetscFunctionReturn(0);
2755e76c431SBarry Smith }
27604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27704d965bbSLois Curfman McInnes /*
2784b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2794b27c08aSLois Curfman McInnes    with SNESCreate_LS().
28004d965bbSLois Curfman McInnes 
28104d965bbSLois Curfman McInnes    Input Parameter:
28204d965bbSLois Curfman McInnes .  snes - the SNES context
28304d965bbSLois Curfman McInnes 
284de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28504d965bbSLois Curfman McInnes  */
2864a2ae208SSatish Balay #undef __FUNCT__
2874b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
288dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2895e76c431SBarry Smith {
290dfbe8321SBarry Smith   PetscErrorCode ierr;
2913a40ed3dSBarry Smith 
2923a40ed3dSBarry Smith   PetscFunctionBegin;
2935baf8537SBarry Smith   if (snes->nwork) {
2944b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29521c89e3eSBarry Smith   }
2965c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
2985e76c431SBarry Smith }
29904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3004a2ae208SSatish Balay #undef __FUNCT__
3014a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
30282bf6240SBarry Smith 
3034b828684SBarry Smith /*@C
3045e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3055e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3065e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3075e76c431SBarry Smith 
308c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
309c7afd0dbSLois Curfman McInnes 
3105e76c431SBarry Smith    Input Parameters:
311c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
312af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3135e76c431SBarry Smith .  x - current iterate
3145e76c431SBarry Smith .  f - residual evaluated at x
315*3c632250SBarry Smith .  y - search direction
3165e76c431SBarry Smith .  w - work vector
317c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3185e76c431SBarry Smith 
319c4a48953SLois Curfman McInnes    Output Parameters:
320c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
321*3c632250SBarry Smith .  w - new iterate
3225e76c431SBarry Smith .  gnorm - 2-norm of g
3235e76c431SBarry Smith .  ynorm - 2-norm of search length
324a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
325fee21e36SBarry Smith 
326c4a48953SLois Curfman McInnes    Options Database Key:
3274b27c08aSLois Curfman McInnes .  -snes_ls basic - Activates SNESNoLineSearch()
328c4a48953SLois Curfman McInnes 
32936851e7fSLois Curfman McInnes    Level: advanced
33036851e7fSLois Curfman McInnes 
33128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33228ae5a21SLois Curfman McInnes 
333f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
334af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3355e76c431SBarry Smith @*/
33663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
3375e76c431SBarry Smith {
338dfbe8321SBarry Smith   PetscErrorCode ierr;
339ea709b57SSatish Balay   PetscScalar    mone = -1.0;
3404b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
341*3c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3425334005bSBarry Smith 
3433a40ed3dSBarry Smith   PetscFunctionBegin;
344a7cc72afSBarry Smith   *flag = PETSC_TRUE;
345d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
346a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
347*3c632250SBarry Smith   ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y   */
348*3c632250SBarry Smith   if (neP->postcheckstep) {
349*3c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3508f99978dSLois Curfman McInnes   }
351*3c632250SBarry Smith   if (changed_y) {
352*3c632250SBarry Smith     ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y   */
353*3c632250SBarry Smith   }
354*3c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
355d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
35619717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
35719717074SBarry Smith   }
358d5e45103SBarry Smith   CHKERRQ(ierr);
359d5e45103SBarry Smith 
360a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
36143e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
362d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3633a40ed3dSBarry Smith   PetscFunctionReturn(0);
3645e76c431SBarry Smith }
36504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3664a2ae208SSatish Balay #undef __FUNCT__
3674a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
36882bf6240SBarry Smith 
36929e0b56fSBarry Smith /*@C
37029e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
37129e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
37229e0b56fSBarry Smith    even compute the norm of the function or search direction; this
37329e0b56fSBarry Smith    is intended only when you know the full step is fine and are
37429e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
375c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
376c7afd0dbSLois Curfman McInnes 
377c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
37829e0b56fSBarry Smith 
37929e0b56fSBarry Smith    Input Parameters:
380c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
381af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
38229e0b56fSBarry Smith .  x - current iterate
38329e0b56fSBarry Smith .  f - residual evaluated at x
384*3c632250SBarry Smith .  y - search direction
38529e0b56fSBarry Smith .  w - work vector
386c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
38729e0b56fSBarry Smith 
38829e0b56fSBarry Smith    Output Parameters:
389c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
390*3c632250SBarry Smith .  w - new iterate
39129e0b56fSBarry Smith .  gnorm - not changed
39229e0b56fSBarry Smith .  ynorm - not changed
393a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
394fee21e36SBarry Smith 
39529e0b56fSBarry Smith    Options Database Key:
3964b27c08aSLois Curfman McInnes .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
39729e0b56fSBarry Smith 
3988cbba510SBarry Smith    Notes:
399ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
400ea56f5baSLois Curfman McInnes    either the options
401ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
402ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
403329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
404329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
405329f5518SBarry Smith 
406329f5518SBarry Smith    During the final iteration this will not evaluate the function at
407329f5518SBarry Smith    the solution point. This is to save a function evaluation while
408329f5518SBarry Smith    using pseudo-timestepping.
4098cbba510SBarry Smith 
410ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
411ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
412ea56f5baSLois Curfman McInnes    correct, since they are not computed.
413ea56f5baSLois Curfman McInnes 
414ea56f5baSLois Curfman McInnes    Level: advanced
4158cbba510SBarry Smith 
41629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
41729e0b56fSBarry Smith 
41829e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
41929e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
42029e0b56fSBarry Smith @*/
42163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
42229e0b56fSBarry Smith {
423dfbe8321SBarry Smith   PetscErrorCode ierr;
424ea709b57SSatish Balay   PetscScalar    mone = -1.0;
4254b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
426*3c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
42729e0b56fSBarry Smith 
4283a40ed3dSBarry Smith   PetscFunctionBegin;
429a7cc72afSBarry Smith   *flag = PETSC_TRUE;
430d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
431*3c632250SBarry Smith   ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y      */
432*3c632250SBarry Smith   if (neP->postcheckstep) {
433*3c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
434*3c632250SBarry Smith   }
435*3c632250SBarry Smith   if (changed_y) {
436*3c632250SBarry Smith     ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y   */
4378f99978dSLois Curfman McInnes   }
438329f5518SBarry Smith 
439329f5518SBarry Smith   /* don't evaluate function the last time through */
440329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
441*3c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
442d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
44319717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
44419717074SBarry Smith     }
445d5e45103SBarry Smith     CHKERRQ(ierr);
446329f5518SBarry Smith   }
447d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4483a40ed3dSBarry Smith   PetscFunctionReturn(0);
44929e0b56fSBarry Smith }
45004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4514a2ae208SSatish Balay #undef __FUNCT__
4524a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4534b828684SBarry Smith /*@C
454f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4555e76c431SBarry Smith 
456c7afd0dbSLois Curfman McInnes    Collective on SNES
457c7afd0dbSLois Curfman McInnes 
4585e76c431SBarry Smith    Input Parameters:
459c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
460af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4615e76c431SBarry Smith .  x - current iterate
4625e76c431SBarry Smith .  f - residual evaluated at x
463*3c632250SBarry Smith .  y - search direction
4645e76c431SBarry Smith .  w - work vector
465c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4665e76c431SBarry Smith 
467393d2d9aSLois Curfman McInnes    Output Parameters:
468c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
469*3c632250SBarry Smith .  w - new iterate
4705e76c431SBarry Smith .  gnorm - 2-norm of g
4715e76c431SBarry Smith .  ynorm - 2-norm of search length
472a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
473fee21e36SBarry Smith 
474c4a48953SLois Curfman McInnes    Options Database Key:
4754b27c08aSLois Curfman McInnes $  -snes_ls cubic - Activates SNESCubicLineSearch()
476c4a48953SLois Curfman McInnes 
4775e76c431SBarry Smith    Notes:
4785e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4795e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4805e76c431SBarry Smith 
48136851e7fSLois Curfman McInnes    Level: advanced
48236851e7fSLois Curfman McInnes 
48328ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
48428ae5a21SLois Curfman McInnes 
485af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4865e76c431SBarry Smith @*/
48763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
4885e76c431SBarry Smith {
489406556e6SLois Curfman McInnes   /*
490406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
491406556e6SLois Curfman McInnes      minimization problem:
492406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
493406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
494406556e6SLois Curfman McInnes    */
495406556e6SLois Curfman McInnes 
4965ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
497329f5518SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
498aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
499ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5006b5873e3SBarry Smith #endif
501dfbe8321SBarry Smith   PetscErrorCode ierr;
502a7cc72afSBarry Smith   PetscInt       count;
5034b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
504ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
505*3c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5065e76c431SBarry Smith 
5073a40ed3dSBarry Smith   PetscFunctionBegin;
508d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
509a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5105e76c431SBarry Smith   alpha   = neP->alpha;
5115e76c431SBarry Smith   maxstep = neP->maxstep;
5125e76c431SBarry Smith   steptol = neP->steptol;
5135e76c431SBarry Smith 
514cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5159e247f21SBarry Smith   if (!*ynorm) {
51663ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
517a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
518*3c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
519a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
520a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
521ad922ac9SBarry Smith     goto theend1;
52294a9d846SBarry Smith   }
5235e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5245e42470aSBarry Smith     scale = maxstep/(*ynorm);
525aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
52663ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr);
5276b5873e3SBarry Smith #else
52863ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr);
5296b5873e3SBarry Smith #endif
530393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5315e76c431SBarry Smith     *ynorm = maxstep;
5325e76c431SBarry Smith   }
533ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5345ca10a37SBarry Smith   minlambda = steptol/rellength;
535a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
536aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
537a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
538329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5395e42470aSBarry Smith #else
540a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5415e42470aSBarry Smith #endif
5425e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5435e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5445e76c431SBarry Smith 
545393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5465334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
54743e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
54863ba0a88SBarry Smith     ierr  = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
54943e71028SBarry Smith     *flag = PETSC_FALSE;
55043e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
55143e71028SBarry Smith     goto theend1;
55243e71028SBarry Smith   }
55319717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
554d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
55519717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
55619717074SBarry Smith   }
557d5e45103SBarry Smith   CHKERRQ(ierr);
558313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
55943e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5605d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
56163ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Using full step\n"));CHKERRQ(ierr);
562ad922ac9SBarry Smith     goto theend1;
5635e76c431SBarry Smith   }
5645e76c431SBarry Smith 
5655e76c431SBarry Smith   /* Fit points with quadratic */
566313b5538SBarry Smith   lambda     = 1.0;
567a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5685e76c431SBarry Smith   lambdaprev = lambda;
5695e76c431SBarry Smith   gnormprev  = *gnorm;
570329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
571ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
572ddd12767SBarry Smith   else                         lambda = lambdatemp;
573393d2d9aSLois Curfman McInnes   ierr      = VecCopy(x,w);CHKERRQ(ierr);
574ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
575aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
576ea4d3ed3SLois Curfman McInnes   clambda   = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5775e42470aSBarry Smith #else
578ea4d3ed3SLois Curfman McInnes   ierr      = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5795e42470aSBarry Smith #endif
58043e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
58163ba0a88SBarry Smith     ierr  = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr);
58243e71028SBarry Smith     *flag = PETSC_FALSE;
58343e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
58443e71028SBarry Smith     goto theend1;
58543e71028SBarry Smith   }
58619717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
587d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
58819717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
58919717074SBarry Smith   }
590d5e45103SBarry Smith   CHKERRQ(ierr);
591cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
59243e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5935ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
59463ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
595ad922ac9SBarry Smith     goto theend1;
5965e76c431SBarry Smith   }
5975e76c431SBarry Smith 
5985e76c431SBarry Smith   /* Fit points with cubic */
5995e76c431SBarry Smith   count = 1;
6008229bfc2SMatthew Knepley   while (count < 10000) {
6015e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
60263ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
60363ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
60443e71028SBarry Smith       *flag = PETSC_FALSE;
60543e71028SBarry Smith       break;
6065e76c431SBarry Smith     }
6075d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6085d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
609ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6102b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6115e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6125e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6135e76c431SBarry Smith     if (a == 0.0) {
6145e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6155e76c431SBarry Smith     } else {
6165e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6175e76c431SBarry Smith     }
6185e76c431SBarry Smith     lambdaprev = lambda;
6195e76c431SBarry Smith     gnormprev  = *gnorm;
620329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
621329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6225e76c431SBarry Smith     else                         lambda     = lambdatemp;
623393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
624ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
625aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
626ea4d3ed3SLois Curfman McInnes     clambda   = lambdaneg;
627393d2d9aSLois Curfman McInnes     ierr      = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6285e42470aSBarry Smith #else
629ea4d3ed3SLois Curfman McInnes     ierr      = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6305e42470aSBarry Smith #endif
63143e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
63263ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
63363ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
634ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
63543e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
636ed98166eSMatthew Knepley       break;
637ed98166eSMatthew Knepley     }
63819717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
639d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64019717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64119717074SBarry Smith     }
642d5e45103SBarry Smith     CHKERRQ(ierr);
643cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
64443e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6455ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
64663ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
64743e71028SBarry Smith       break;
6482b022350SLois Curfman McInnes     } else {
64963ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda));CHKERRQ(ierr);
6505e76c431SBarry Smith     }
6515e76c431SBarry Smith     count++;
6525e76c431SBarry Smith   }
6538229bfc2SMatthew Knepley   if (count >= 10000) {
654cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6558229bfc2SMatthew Knepley   }
656ad922ac9SBarry Smith   theend1:
6578f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
658*3c632250SBarry Smith   if (neP->postcheckstep && *flag) {
659*3c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
660*3c632250SBarry Smith     if (changed_y) {
661*3c632250SBarry Smith       ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);
662*3c632250SBarry Smith     }
663*3c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
664*3c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
665d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
66619717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
66719717074SBarry Smith       }
668d5e45103SBarry Smith       CHKERRQ(ierr);
6698f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67043e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
671*3c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6728f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
673*3c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6748f99978dSLois Curfman McInnes     }
6758f99978dSLois Curfman McInnes   }
676d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6773a40ed3dSBarry Smith   PetscFunctionReturn(0);
6785e76c431SBarry Smith }
67904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6804a2ae208SSatish Balay #undef __FUNCT__
6814a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6824b828684SBarry Smith /*@C
683f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6845e76c431SBarry Smith 
685c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
686c7afd0dbSLois Curfman McInnes 
6875e42470aSBarry Smith    Input Parameters:
688c7afd0dbSLois Curfman McInnes +  snes - the SNES context
689af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6905e42470aSBarry Smith .  x - current iterate
6915e42470aSBarry Smith .  f - residual evaluated at x
692*3c632250SBarry Smith .  y - search direction
6935e42470aSBarry Smith .  w - work vector
694c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6955e42470aSBarry Smith 
696c4a48953SLois Curfman McInnes    Output Parameters:
697c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
698*3c632250SBarry Smith .  w - new iterate
6995e42470aSBarry Smith .  gnorm - 2-norm of g
7005e42470aSBarry Smith .  ynorm - 2-norm of search length
701a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
702fee21e36SBarry Smith 
703c4a48953SLois Curfman McInnes    Options Database Key:
7044b27c08aSLois Curfman McInnes .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
705c4a48953SLois Curfman McInnes 
7065e42470aSBarry Smith    Notes:
7074b27c08aSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNESLS method.
7085e42470aSBarry Smith 
70936851e7fSLois Curfman McInnes    Level: advanced
71036851e7fSLois Curfman McInnes 
711f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
712f59ffdedSLois Curfman McInnes 
713af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
7145e42470aSBarry Smith @*/
71563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
7165e76c431SBarry Smith {
717406556e6SLois Curfman McInnes   /*
718406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
719406556e6SLois Curfman McInnes      minimization problem:
720406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
721406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
722406556e6SLois Curfman McInnes    */
7235ca10a37SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
724aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
725ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7266b5873e3SBarry Smith #endif
727dfbe8321SBarry Smith   PetscErrorCode ierr;
728a7cc72afSBarry Smith   PetscInt       count;
7294b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
730ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
731*3c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7325e76c431SBarry Smith 
7333a40ed3dSBarry Smith   PetscFunctionBegin;
734d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
735a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7365e42470aSBarry Smith   alpha   = neP->alpha;
7375e42470aSBarry Smith   maxstep = neP->maxstep;
7385e42470aSBarry Smith   steptol = neP->steptol;
7395e76c431SBarry Smith 
7403505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74163b9fa5eSBarry Smith   if (*ynorm == 0.0) {
74263ba0a88SBarry Smith     ierr   = PetscLogInfo((snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
743b37302c6SLois Curfman McInnes     *gnorm = fnorm;
744b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
745b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
746a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
747ad922ac9SBarry Smith     goto theend2;
74894a9d846SBarry Smith   }
7495e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7505e42470aSBarry Smith     scale  = maxstep/(*ynorm);
751393d2d9aSLois Curfman McInnes     ierr   = VecScale(&scale,y);CHKERRQ(ierr);
7525e42470aSBarry Smith     *ynorm = maxstep;
7535e76c431SBarry Smith   }
754ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7555ca10a37SBarry Smith   minlambda = steptol/rellength;
756a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
757aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
758a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
759329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7605e42470aSBarry Smith #else
761a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7625e42470aSBarry Smith #endif
7635e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7645e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7655e42470aSBarry Smith 
766393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7675334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
76843e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
76963ba0a88SBarry Smith     ierr  = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
77043e71028SBarry Smith     ierr  = VecCopy(w,y);CHKERRQ(ierr);
77143e71028SBarry Smith     *flag = PETSC_FALSE;
77243e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77343e71028SBarry Smith     goto theend2;
77443e71028SBarry Smith   }
77519717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
776d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
77719717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
77819717074SBarry Smith   }
779d5e45103SBarry Smith   CHKERRQ(ierr);
780cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
78143e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7825d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
783393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
78463ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Using full step\n"));CHKERRQ(ierr);
785ad922ac9SBarry Smith     goto theend2;
7865e42470aSBarry Smith   }
7875e42470aSBarry Smith 
7885e42470aSBarry Smith   /* Fit points with quadratic */
789313b5538SBarry Smith   lambda = 1.0;
7905e42470aSBarry Smith   count = 1;
7915ca10a37SBarry Smith   while (PETSC_TRUE) {
7925e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
79363ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
79463ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
795f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
79643e71028SBarry Smith       *flag = PETSC_FALSE;
79743e71028SBarry Smith       break;
7985e42470aSBarry Smith     }
799a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
800329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
801329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
802329f5518SBarry Smith     else                         lambda     = lambdatemp;
803393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
8043505fcc1SBarry Smith     lambdaneg = -lambda;
805aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
8063505fcc1SBarry Smith     clambda   = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
8075e42470aSBarry Smith #else
8083505fcc1SBarry Smith     ierr      = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
8095e42470aSBarry Smith #endif
81043e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
81163ba0a88SBarry Smith       ierr  = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
81263ba0a88SBarry Smith       ierr  = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
81343e71028SBarry Smith       ierr  = VecCopy(w,y);CHKERRQ(ierr);
814ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
81543e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
816ed98166eSMatthew Knepley       break;
817ed98166eSMatthew Knepley     }
81819717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
819d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
82019717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
82119717074SBarry Smith     }
822d5e45103SBarry Smith     CHKERRQ(ierr);
823cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82443e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8255ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
826393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
82763ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr);
82806259719SBarry Smith       break;
8295e42470aSBarry Smith     }
8305e42470aSBarry Smith     count++;
8315e42470aSBarry Smith   }
832ad922ac9SBarry Smith   theend2:
8338f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
834*3c632250SBarry Smith   if (neP->postcheckstep) {
835*3c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
836*3c632250SBarry Smith     if (changed_y) {
837*3c632250SBarry Smith       ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);
838*3c632250SBarry Smith     }
839*3c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
840*3c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
841d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
84219717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
84319717074SBarry Smith       }
844d5e45103SBarry Smith       CHKERRQ(ierr);
8458f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
846*3c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8478f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
848*3c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
849*3c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8508f99978dSLois Curfman McInnes     }
8518f99978dSLois Curfman McInnes   }
852d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8533a40ed3dSBarry Smith   PetscFunctionReturn(0);
8545e76c431SBarry Smith }
8552343ba6eSBarry Smith 
85604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8574a2ae208SSatish Balay #undef __FUNCT__
8584a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
859c9e6a524SLois Curfman McInnes /*@C
86077c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
8614b27c08aSLois Curfman McInnes    by the method SNESLS.
8625e76c431SBarry Smith 
8635e76c431SBarry Smith    Input Parameters:
864c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
865af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
866c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8675e76c431SBarry Smith 
868fee21e36SBarry Smith    Collective on SNES
869fee21e36SBarry Smith 
870c4a48953SLois Curfman McInnes    Available Routines:
871c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
872f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
873f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
874af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8755e76c431SBarry Smith 
876c4a48953SLois Curfman McInnes     Options Database Keys:
8774b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8784b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8794b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8804b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8813304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8823304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
883c4a48953SLois Curfman McInnes 
8845e76c431SBarry Smith    Calling sequence of func:
885bd208895SLois Curfman McInnes .vb
886af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
887a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
888bd208895SLois Curfman McInnes .ve
8895e76c431SBarry Smith 
8905e76c431SBarry Smith     Input parameters for func:
891c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
892af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8935e76c431SBarry Smith .   x - current iterate
8945e76c431SBarry Smith .   f - residual evaluated at x
895*3c632250SBarry Smith .   y - search direction
8965e76c431SBarry Smith .   w - work vector
897c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8985e76c431SBarry Smith 
8995e76c431SBarry Smith     Output parameters for func:
900c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
901*3c632250SBarry Smith .   w - new iterate
9025e76c431SBarry Smith .   gnorm - 2-norm of g
9035e76c431SBarry Smith .   ynorm - 2-norm of search length
904a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
905f59ffdedSLois Curfman McInnes 
90636851e7fSLois Curfman McInnes     Level: advanced
90736851e7fSLois Curfman McInnes 
908f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
909f59ffdedSLois Curfman McInnes 
9105d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
911*3c632250SBarry Smith           SNESLineSearchSetPostCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
9125e76c431SBarry Smith @*/
91363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9145e76c431SBarry Smith {
915a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
91682bf6240SBarry Smith 
9173a40ed3dSBarry Smith   PetscFunctionBegin;
918c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
91982bf6240SBarry Smith   if (f) {
920af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92182bf6240SBarry Smith   }
9223a40ed3dSBarry Smith   PetscFunctionReturn(0);
9235e76c431SBarry Smith }
9248e019c35SBarry Smith 
925a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
92604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
927fb2e594dSBarry Smith EXTERN_C_BEGIN
9284a2ae208SSatish Balay #undef __FUNCT__
9294a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
93063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
93182bf6240SBarry Smith {
93282bf6240SBarry Smith   PetscFunctionBegin;
9334b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9344b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
93582bf6240SBarry Smith   PetscFunctionReturn(0);
93682bf6240SBarry Smith }
937fb2e594dSBarry Smith EXTERN_C_END
93804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9394a2ae208SSatish Balay #undef __FUNCT__
940*3c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
941c8dd0c0dSLois Curfman McInnes /*@C
942*3c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9434b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
944c8dd0c0dSLois Curfman McInnes 
945c8dd0c0dSLois Curfman McInnes    Input Parameters:
946c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
947*3c632250SBarry Smith .  func - pointer to function
948c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
949c8dd0c0dSLois Curfman McInnes 
950c8dd0c0dSLois Curfman McInnes    Collective on SNES
951c8dd0c0dSLois Curfman McInnes 
952c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
953c8dd0c0dSLois Curfman McInnes .vb
954*3c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
955c8dd0c0dSLois Curfman McInnes .ve
956b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
957b1ae27eaSLois Curfman McInnes    on failure.
958c8dd0c0dSLois Curfman McInnes 
959c8dd0c0dSLois Curfman McInnes    Input parameters for func:
960c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
961c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
962*3c632250SBarry Smith .  x - previous iterate
963*3c632250SBarry Smith .  y - new search direction and length
964*3c632250SBarry Smith -  w - current candidate iterate
965c8dd0c0dSLois Curfman McInnes 
966c8dd0c0dSLois Curfman McInnes    Output parameters for func:
967*3c632250SBarry Smith +  y - search direction (possibly changed)
968*3c632250SBarry Smith .  w - current iterate (possibly modified)
969*3c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
970*3c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
971c8dd0c0dSLois Curfman McInnes 
972c8dd0c0dSLois Curfman McInnes    Level: advanced
973c8dd0c0dSLois Curfman McInnes 
9749e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9759e247f21SBarry Smith 
976*3c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
977*3c632250SBarry Smith 
978*3c632250SBarry Smith    On input w = x + y
979*3c632250SBarry Smith 
9809e247f21SBarry Smith    SNESNoLineSearch() and SNESNoLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
981b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
982ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9838f99978dSLois Curfman McInnes 
9849e247f21SBarry Smith    SNESQuadraticLineSearch() and SNESCubicLineSearch() (1) compute a candidate iterate u_{i+1} as well as a
985ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
986ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
987ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9889e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
989b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9908f99978dSLois Curfman McInnes 
991c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
992c8dd0c0dSLois Curfman McInnes 
993c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
994c8dd0c0dSLois Curfman McInnes @*/
995*3c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
996c8dd0c0dSLois Curfman McInnes {
997*3c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
998c8dd0c0dSLois Curfman McInnes 
999c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
1000*3c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1001c8dd0c0dSLois Curfman McInnes   if (f) {
1002c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1003c8dd0c0dSLois Curfman McInnes   }
1004c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1005c8dd0c0dSLois Curfman McInnes }
1006c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1007*3c632250SBarry Smith typedef PetscErrorCode (*FCN)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
1008c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10094a2ae208SSatish Balay #undef __FUNCT__
1010*3c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
1011*3c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN func,void *checkctx)
1012c8dd0c0dSLois Curfman McInnes {
1013c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
1014*3c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
1015*3c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1016c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1017c8dd0c0dSLois Curfman McInnes }
1018c8dd0c0dSLois Curfman McInnes EXTERN_C_END
1019c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
102004d965bbSLois Curfman McInnes /*
10214b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
1022329e7e40SMatthew Knepley 
1023329e7e40SMatthew Knepley    Input Parameter:
1024329e7e40SMatthew Knepley .  snes - the SNES context
1025329e7e40SMatthew Knepley 
1026329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
1027329e7e40SMatthew Knepley */
1028329e7e40SMatthew Knepley #undef __FUNCT__
10294b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
10306849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
1031329e7e40SMatthew Knepley {
10324b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
1033329e7e40SMatthew Knepley 
1034329e7e40SMatthew Knepley   PetscFunctionBegin;
10354b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
10364b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
10374b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
10384b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
10394b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
1040329e7e40SMatthew Knepley   PetscFunctionReturn(0);
1041329e7e40SMatthew Knepley }
1042329e7e40SMatthew Knepley 
1043329e7e40SMatthew Knepley /*
10444b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
104504d965bbSLois Curfman McInnes 
104604d965bbSLois Curfman McInnes    Input Parameters:
104704d965bbSLois Curfman McInnes .  SNES - the SNES context
104804d965bbSLois Curfman McInnes .  viewer - visualization context
104904d965bbSLois Curfman McInnes 
105004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
105104d965bbSLois Curfman McInnes */
10524a2ae208SSatish Balay #undef __FUNCT__
10534b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10546849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1055a935fc98SLois Curfman McInnes {
10564b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
10572fc52814SBarry Smith   const char     *cstr;
1058dfbe8321SBarry Smith   PetscErrorCode ierr;
105932077d6dSBarry Smith   PetscTruth     iascii;
1060a935fc98SLois Curfman McInnes 
10613a40ed3dSBarry Smith   PetscFunctionBegin;
106232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
106332077d6dSBarry Smith   if (iascii) {
106419bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
106519bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
106619bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
106719bcc07fSBarry Smith     else                                                cstr = "unknown";
1068b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1069b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
10705cd90555SBarry Smith   } else {
107179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
107219bcc07fSBarry Smith   }
10733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1074a935fc98SLois Curfman McInnes }
107504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
107604d965bbSLois Curfman McInnes /*
10774b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
107804d965bbSLois Curfman McInnes 
107904d965bbSLois Curfman McInnes    Input Parameter:
108004d965bbSLois Curfman McInnes .  snes - the SNES context
108104d965bbSLois Curfman McInnes 
108204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
108304d965bbSLois Curfman McInnes */
10844a2ae208SSatish Balay #undef __FUNCT__
10854b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
10866849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
10875e42470aSBarry Smith {
10884b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1089e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1090dfbe8321SBarry Smith   PetscErrorCode ierr;
1091a7cc72afSBarry Smith   PetscInt       indx;
1092f1af5d2fSBarry Smith   PetscTruth     flg;
10935e42470aSBarry Smith 
10943a40ed3dSBarry Smith   PetscFunctionBegin;
1095b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
10964b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
10974b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
10984b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1099186905e3SBarry Smith 
110022e36657SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
110125cdf11fSBarry Smith     if (flg) {
110222e36657SBarry Smith       switch (indx) {
1103b49fd9e1SBarry Smith       case 0:
1104af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1105b49fd9e1SBarry Smith         break;
1106b49fd9e1SBarry Smith       case 1:
1107af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1108b49fd9e1SBarry Smith         break;
1109b49fd9e1SBarry Smith       case 2:
1110af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1111b49fd9e1SBarry Smith         break;
1112b49fd9e1SBarry Smith       case 3:
1113af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1114b49fd9e1SBarry Smith         break;
11155e42470aSBarry Smith       }
11165e42470aSBarry Smith     }
1117b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11183a40ed3dSBarry Smith   PetscFunctionReturn(0);
11195e42470aSBarry Smith }
112004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1121ebe3b25bSBarry Smith /*MC
1122ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
112304d965bbSLois Curfman McInnes 
1124ebe3b25bSBarry Smith    Options Database:
1125ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1126ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1127ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1128ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1129ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1130ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
113104d965bbSLois Curfman McInnes 
1132ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
113304d965bbSLois Curfman McInnes 
1134ee3001cbSBarry Smith    Level: beginner
1135ee3001cbSBarry Smith 
1136ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1137*3c632250SBarry Smith            SNESLineSearchSetPostCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1138ebe3b25bSBarry Smith           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1139ebe3b25bSBarry Smith 
1140ebe3b25bSBarry Smith M*/
1141fb2e594dSBarry Smith EXTERN_C_BEGIN
11424a2ae208SSatish Balay #undef __FUNCT__
11434b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
114463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
11455e42470aSBarry Smith {
1146dfbe8321SBarry Smith   PetscErrorCode ierr;
11474b27c08aSLois Curfman McInnes   SNES_LS        *neP;
11485e42470aSBarry Smith 
11493a40ed3dSBarry Smith   PetscFunctionBegin;
11504b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
11514b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
11524b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
11534b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
11544b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
11554b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
11564b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
11575baf8537SBarry Smith   snes->nwork           = 0;
11585e42470aSBarry Smith 
11594b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
116052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
11615e42470aSBarry Smith   snes->data    	= (void*)neP;
11625e42470aSBarry Smith   neP->alpha		= 1.e-4;
11635e42470aSBarry Smith   neP->maxstep		= 1.e8;
11645e42470aSBarry Smith   neP->steptol		= 1.e-12;
11655e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1166c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1167*3c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
1168*3c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
1169*3c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
1170*3c632250SBarry Smith   neP->precheck         = PETSC_NULL;
117182bf6240SBarry Smith 
11725d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
1173*3c632250SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
117482bf6240SBarry Smith 
11753a40ed3dSBarry Smith   PetscFunctionReturn(0);
11765e42470aSBarry Smith }
1177fb2e594dSBarry Smith EXTERN_C_END
117882bf6240SBarry Smith 
117982bf6240SBarry Smith 
118082bf6240SBarry Smith 
1181