xref: /petsc/src/snes/impls/ls/ls.c (revision 9c3ca13afd1a6531fa13996d6426b76a28ba2985)
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 
185*9c3ca13aSBarry Smith     if (neP->precheckstep) {
186*9c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
187*9c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
188*9c3ca13aSBarry Smith     }
189*9c3ca13aSBarry Smith 
190b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
19174637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
19285471664SBarry Smith     }
19374637425SBarry Smith 
19490cbd584SBarry Smith     /* should check what happened to the linear solve? */
1953505fcc1SBarry Smith     snes->linear_its += lits;
19663ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr);
197ea4d3ed3SLois Curfman McInnes 
198ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
199ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
200ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
201ea4d3ed3SLois Curfman McInnes     */
20281b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
203a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
20463ba0a88SBarry 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);
20543e71028SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
206a3b891d8SBarry Smith 
207a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
2083c632250SBarry Smith     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
209a3b891d8SBarry Smith     fnorm = gnorm;
210a3b891d8SBarry Smith 
2115ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2125ed2d596SBarry Smith     snes->iter = i+1;
2135ed2d596SBarry Smith     snes->norm = fnorm;
2145ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2155ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2165ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2175ed2d596SBarry Smith 
218a7cc72afSBarry Smith     if (!lssucceed) {
2198a5d9ee4SBarry Smith       PetscTruth ismin;
22050ffb88aSMatthew Knepley 
22150ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2223505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2238a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2248a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2253505fcc1SBarry Smith         break;
2263505fcc1SBarry Smith       }
22750ffb88aSMatthew Knepley     }
2285e76c431SBarry Smith 
2295e76c431SBarry Smith     /* Test for convergence */
23029e0b56fSBarry Smith     if (snes->converged) {
23129e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2323505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2333505fcc1SBarry Smith       if (snes->reason) {
23416c95cacSBarry Smith         break;
23516c95cacSBarry Smith       }
23616c95cacSBarry Smith     }
23729e0b56fSBarry Smith   }
23839e2f89bSBarry Smith   if (X != snes->vec_sol) {
239393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
240e15f7bb5SBarry Smith   }
241e15f7bb5SBarry Smith   if (F != snes->vec_func) {
242e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
243e15f7bb5SBarry Smith   }
24439e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
24539e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
24652392280SLois Curfman McInnes   if (i == maxits) {
24763ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr);
2483505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
24952392280SLois Curfman McInnes   }
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
2515e76c431SBarry Smith }
25204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25304d965bbSLois Curfman McInnes /*
2544b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2554b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25604d965bbSLois Curfman McInnes 
25704d965bbSLois Curfman McInnes    Input Parameter:
25804d965bbSLois Curfman McInnes .  snes - the SNES context
25904d965bbSLois Curfman McInnes .  x - the solution vector
26004d965bbSLois Curfman McInnes 
26104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26204d965bbSLois Curfman McInnes 
26304d965bbSLois Curfman McInnes    Notes:
2644b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26504d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26604d965bbSLois Curfman McInnes    the call to SNESSolve().
26704d965bbSLois Curfman McInnes  */
2684a2ae208SSatish Balay #undef __FUNCT__
2694b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
270dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2715e76c431SBarry Smith {
272dfbe8321SBarry Smith   PetscErrorCode ierr;
2733a40ed3dSBarry Smith 
2743a40ed3dSBarry Smith   PetscFunctionBegin;
27581b6cf68SLois Curfman McInnes   snes->nwork = 4;
276d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
277efee365bSSatish Balay   ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
27881b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
2805e76c431SBarry Smith }
28104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28204d965bbSLois Curfman McInnes /*
2834b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2844b27c08aSLois Curfman McInnes    with SNESCreate_LS().
28504d965bbSLois Curfman McInnes 
28604d965bbSLois Curfman McInnes    Input Parameter:
28704d965bbSLois Curfman McInnes .  snes - the SNES context
28804d965bbSLois Curfman McInnes 
289de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29004d965bbSLois Curfman McInnes  */
2914a2ae208SSatish Balay #undef __FUNCT__
2924b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
293dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2945e76c431SBarry Smith {
295dfbe8321SBarry Smith   PetscErrorCode ierr;
2963a40ed3dSBarry Smith 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2985baf8537SBarry Smith   if (snes->nwork) {
2994b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
30021c89e3eSBarry Smith   }
3015c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3023a40ed3dSBarry Smith   PetscFunctionReturn(0);
3035e76c431SBarry Smith }
30404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3054a2ae208SSatish Balay #undef __FUNCT__
3064a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
30782bf6240SBarry Smith 
3084b828684SBarry Smith /*@C
3095e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3105e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3115e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3125e76c431SBarry Smith 
313c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
314c7afd0dbSLois Curfman McInnes 
3155e76c431SBarry Smith    Input Parameters:
316c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
317af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3185e76c431SBarry Smith .  x - current iterate
3195e76c431SBarry Smith .  f - residual evaluated at x
3203c632250SBarry Smith .  y - search direction
3215e76c431SBarry Smith .  w - work vector
322c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3235e76c431SBarry Smith 
324c4a48953SLois Curfman McInnes    Output Parameters:
325c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3263c632250SBarry Smith .  w - new iterate
3275e76c431SBarry Smith .  gnorm - 2-norm of g
3285e76c431SBarry Smith .  ynorm - 2-norm of search length
329a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
330fee21e36SBarry Smith 
331c4a48953SLois Curfman McInnes    Options Database Key:
3324b27c08aSLois Curfman McInnes .  -snes_ls basic - Activates SNESNoLineSearch()
333c4a48953SLois Curfman McInnes 
33436851e7fSLois Curfman McInnes    Level: advanced
33536851e7fSLois Curfman McInnes 
33628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33728ae5a21SLois Curfman McInnes 
338f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
339af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3405e76c431SBarry Smith @*/
34163dd3a1aSKris 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)
3425e76c431SBarry Smith {
343dfbe8321SBarry Smith   PetscErrorCode ierr;
344ea709b57SSatish Balay   PetscScalar    mone = -1.0;
3454b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3463c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3475334005bSBarry Smith 
3483a40ed3dSBarry Smith   PetscFunctionBegin;
349a7cc72afSBarry Smith   *flag = PETSC_TRUE;
350d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
351a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
3523c632250SBarry Smith   ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y   */
3533c632250SBarry Smith   if (neP->postcheckstep) {
3543c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3558f99978dSLois Curfman McInnes   }
3563c632250SBarry Smith   if (changed_y) {
3573c632250SBarry Smith     ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y   */
3583c632250SBarry Smith   }
3593c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
360d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
36119717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
36219717074SBarry Smith   }
363d5e45103SBarry Smith   CHKERRQ(ierr);
364d5e45103SBarry Smith 
365a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
36643e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
367d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3683a40ed3dSBarry Smith   PetscFunctionReturn(0);
3695e76c431SBarry Smith }
37004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3714a2ae208SSatish Balay #undef __FUNCT__
3724a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
37382bf6240SBarry Smith 
37429e0b56fSBarry Smith /*@C
37529e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
37629e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
37729e0b56fSBarry Smith    even compute the norm of the function or search direction; this
37829e0b56fSBarry Smith    is intended only when you know the full step is fine and are
37929e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
380c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
381c7afd0dbSLois Curfman McInnes 
382c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
38329e0b56fSBarry Smith 
38429e0b56fSBarry Smith    Input Parameters:
385c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
386af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
38729e0b56fSBarry Smith .  x - current iterate
38829e0b56fSBarry Smith .  f - residual evaluated at x
3893c632250SBarry Smith .  y - search direction
39029e0b56fSBarry Smith .  w - work vector
391c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
39229e0b56fSBarry Smith 
39329e0b56fSBarry Smith    Output Parameters:
394c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3953c632250SBarry Smith .  w - new iterate
39629e0b56fSBarry Smith .  gnorm - not changed
39729e0b56fSBarry Smith .  ynorm - not changed
398a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
399fee21e36SBarry Smith 
40029e0b56fSBarry Smith    Options Database Key:
4014b27c08aSLois Curfman McInnes .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
40229e0b56fSBarry Smith 
4038cbba510SBarry Smith    Notes:
404ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
405ea56f5baSLois Curfman McInnes    either the options
406ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
407ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
408329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
409329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
410329f5518SBarry Smith 
411329f5518SBarry Smith    During the final iteration this will not evaluate the function at
412329f5518SBarry Smith    the solution point. This is to save a function evaluation while
413329f5518SBarry Smith    using pseudo-timestepping.
4148cbba510SBarry Smith 
415ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
416ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
417ea56f5baSLois Curfman McInnes    correct, since they are not computed.
418ea56f5baSLois Curfman McInnes 
419ea56f5baSLois Curfman McInnes    Level: advanced
4208cbba510SBarry Smith 
42129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
42229e0b56fSBarry Smith 
42329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
42429e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
42529e0b56fSBarry Smith @*/
42663dd3a1aSKris 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)
42729e0b56fSBarry Smith {
428dfbe8321SBarry Smith   PetscErrorCode ierr;
429ea709b57SSatish Balay   PetscScalar    mone = -1.0;
4304b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4313c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
43229e0b56fSBarry Smith 
4333a40ed3dSBarry Smith   PetscFunctionBegin;
434a7cc72afSBarry Smith   *flag = PETSC_TRUE;
435d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4363c632250SBarry Smith   ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y      */
4373c632250SBarry Smith   if (neP->postcheckstep) {
4383c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4393c632250SBarry Smith   }
4403c632250SBarry Smith   if (changed_y) {
4413c632250SBarry Smith     ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);            /* w <- x - y   */
4428f99978dSLois Curfman McInnes   }
443329f5518SBarry Smith 
444329f5518SBarry Smith   /* don't evaluate function the last time through */
445329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4463c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
447d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
44819717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
44919717074SBarry Smith     }
450d5e45103SBarry Smith     CHKERRQ(ierr);
451329f5518SBarry Smith   }
452d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4533a40ed3dSBarry Smith   PetscFunctionReturn(0);
45429e0b56fSBarry Smith }
45504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4564a2ae208SSatish Balay #undef __FUNCT__
4574a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4584b828684SBarry Smith /*@C
459f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4605e76c431SBarry Smith 
461c7afd0dbSLois Curfman McInnes    Collective on SNES
462c7afd0dbSLois Curfman McInnes 
4635e76c431SBarry Smith    Input Parameters:
464c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
465af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4665e76c431SBarry Smith .  x - current iterate
4675e76c431SBarry Smith .  f - residual evaluated at x
4683c632250SBarry Smith .  y - search direction
4695e76c431SBarry Smith .  w - work vector
470c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4715e76c431SBarry Smith 
472393d2d9aSLois Curfman McInnes    Output Parameters:
473c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4743c632250SBarry Smith .  w - new iterate
4755e76c431SBarry Smith .  gnorm - 2-norm of g
4765e76c431SBarry Smith .  ynorm - 2-norm of search length
477a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
478fee21e36SBarry Smith 
479c4a48953SLois Curfman McInnes    Options Database Key:
4804b27c08aSLois Curfman McInnes $  -snes_ls cubic - Activates SNESCubicLineSearch()
481c4a48953SLois Curfman McInnes 
4825e76c431SBarry Smith    Notes:
4835e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4845e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4855e76c431SBarry Smith 
48636851e7fSLois Curfman McInnes    Level: advanced
48736851e7fSLois Curfman McInnes 
48828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
48928ae5a21SLois Curfman McInnes 
490af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4915e76c431SBarry Smith @*/
49263dd3a1aSKris 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)
4935e76c431SBarry Smith {
494406556e6SLois Curfman McInnes   /*
495406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
496406556e6SLois Curfman McInnes      minimization problem:
497406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
498406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
499406556e6SLois Curfman McInnes    */
500406556e6SLois Curfman McInnes 
5015ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
502329f5518SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
503aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
504ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5056b5873e3SBarry Smith #endif
506dfbe8321SBarry Smith   PetscErrorCode ierr;
507a7cc72afSBarry Smith   PetscInt       count;
5084b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
509ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
5103c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5115e76c431SBarry Smith 
5123a40ed3dSBarry Smith   PetscFunctionBegin;
513d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
514a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5155e76c431SBarry Smith   alpha   = neP->alpha;
5165e76c431SBarry Smith   maxstep = neP->maxstep;
5175e76c431SBarry Smith   steptol = neP->steptol;
5185e76c431SBarry Smith 
519cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5209e247f21SBarry Smith   if (!*ynorm) {
52163ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
522a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5233c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
524a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
525a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
526ad922ac9SBarry Smith     goto theend1;
52794a9d846SBarry Smith   }
5285e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5295e42470aSBarry Smith     scale = maxstep/(*ynorm);
530aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
53163ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr);
5326b5873e3SBarry Smith #else
53363ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr);
5346b5873e3SBarry Smith #endif
535393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5365e76c431SBarry Smith     *ynorm = maxstep;
5375e76c431SBarry Smith   }
538ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5395ca10a37SBarry Smith   minlambda = steptol/rellength;
540a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
541aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
542a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
543329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5445e42470aSBarry Smith #else
545a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5465e42470aSBarry Smith #endif
5475e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5485e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5495e76c431SBarry Smith 
550393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5515334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
55243e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
55363ba0a88SBarry Smith     ierr  = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
55443e71028SBarry Smith     *flag = PETSC_FALSE;
55543e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
55643e71028SBarry Smith     goto theend1;
55743e71028SBarry Smith   }
55819717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
559d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
56019717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
56119717074SBarry Smith   }
562d5e45103SBarry Smith   CHKERRQ(ierr);
563313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56443e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5655d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
56663ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Using full step\n"));CHKERRQ(ierr);
567ad922ac9SBarry Smith     goto theend1;
5685e76c431SBarry Smith   }
5695e76c431SBarry Smith 
5705e76c431SBarry Smith   /* Fit points with quadratic */
571313b5538SBarry Smith   lambda     = 1.0;
572a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5735e76c431SBarry Smith   lambdaprev = lambda;
5745e76c431SBarry Smith   gnormprev  = *gnorm;
575329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
576ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
577ddd12767SBarry Smith   else                         lambda = lambdatemp;
578393d2d9aSLois Curfman McInnes   ierr      = VecCopy(x,w);CHKERRQ(ierr);
579ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
580aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
581ea4d3ed3SLois Curfman McInnes   clambda   = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5825e42470aSBarry Smith #else
583ea4d3ed3SLois Curfman McInnes   ierr      = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5845e42470aSBarry Smith #endif
58543e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
58663ba0a88SBarry Smith     ierr  = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr);
58743e71028SBarry Smith     *flag = PETSC_FALSE;
58843e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
58943e71028SBarry Smith     goto theend1;
59043e71028SBarry Smith   }
59119717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
592d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59319717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59419717074SBarry Smith   }
595d5e45103SBarry Smith   CHKERRQ(ierr);
596cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
59743e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5985ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
59963ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
600ad922ac9SBarry Smith     goto theend1;
6015e76c431SBarry Smith   }
6025e76c431SBarry Smith 
6035e76c431SBarry Smith   /* Fit points with cubic */
6045e76c431SBarry Smith   count = 1;
6058229bfc2SMatthew Knepley   while (count < 10000) {
6065e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
60763ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
60863ba0a88SBarry 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);
60943e71028SBarry Smith       *flag = PETSC_FALSE;
61043e71028SBarry Smith       break;
6115e76c431SBarry Smith     }
6125d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6135d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
614ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6152b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6165e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6175e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6185e76c431SBarry Smith     if (a == 0.0) {
6195e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6205e76c431SBarry Smith     } else {
6215e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6225e76c431SBarry Smith     }
6235e76c431SBarry Smith     lambdaprev = lambda;
6245e76c431SBarry Smith     gnormprev  = *gnorm;
625329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
626329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6275e76c431SBarry Smith     else                         lambda     = lambdatemp;
628393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
629ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
630aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
631ea4d3ed3SLois Curfman McInnes     clambda   = lambdaneg;
632393d2d9aSLois Curfman McInnes     ierr      = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6335e42470aSBarry Smith #else
634ea4d3ed3SLois Curfman McInnes     ierr      = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6355e42470aSBarry Smith #endif
63643e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
63763ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
63863ba0a88SBarry 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);
639ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
64043e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
641ed98166eSMatthew Knepley       break;
642ed98166eSMatthew Knepley     }
64319717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
644d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64519717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64619717074SBarry Smith     }
647d5e45103SBarry Smith     CHKERRQ(ierr);
648cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
64943e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6505ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
65163ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
65243e71028SBarry Smith       break;
6532b022350SLois Curfman McInnes     } else {
65463ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda));CHKERRQ(ierr);
6555e76c431SBarry Smith     }
6565e76c431SBarry Smith     count++;
6575e76c431SBarry Smith   }
6588229bfc2SMatthew Knepley   if (count >= 10000) {
659cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6608229bfc2SMatthew Knepley   }
661ad922ac9SBarry Smith   theend1:
6628f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6633c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6643c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6653c632250SBarry Smith     if (changed_y) {
6663c632250SBarry Smith       ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);
6673c632250SBarry Smith     }
6683c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6693c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
670d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
67119717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
67219717074SBarry Smith       }
673d5e45103SBarry Smith       CHKERRQ(ierr);
6748f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67543e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6763c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6778f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6783c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6798f99978dSLois Curfman McInnes     }
6808f99978dSLois Curfman McInnes   }
681d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6823a40ed3dSBarry Smith   PetscFunctionReturn(0);
6835e76c431SBarry Smith }
68404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6854a2ae208SSatish Balay #undef __FUNCT__
6864a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6874b828684SBarry Smith /*@C
688f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6895e76c431SBarry Smith 
690c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
691c7afd0dbSLois Curfman McInnes 
6925e42470aSBarry Smith    Input Parameters:
693c7afd0dbSLois Curfman McInnes +  snes - the SNES context
694af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6955e42470aSBarry Smith .  x - current iterate
6965e42470aSBarry Smith .  f - residual evaluated at x
6973c632250SBarry Smith .  y - search direction
6985e42470aSBarry Smith .  w - work vector
699c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
7005e42470aSBarry Smith 
701c4a48953SLois Curfman McInnes    Output Parameters:
702c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
7033c632250SBarry Smith .  w - new iterate
7045e42470aSBarry Smith .  gnorm - 2-norm of g
7055e42470aSBarry Smith .  ynorm - 2-norm of search length
706a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
707fee21e36SBarry Smith 
708c4a48953SLois Curfman McInnes    Options Database Key:
7094b27c08aSLois Curfman McInnes .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
710c4a48953SLois Curfman McInnes 
7115e42470aSBarry Smith    Notes:
7124b27c08aSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNESLS method.
7135e42470aSBarry Smith 
71436851e7fSLois Curfman McInnes    Level: advanced
71536851e7fSLois Curfman McInnes 
716f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
717f59ffdedSLois Curfman McInnes 
718af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
7195e42470aSBarry Smith @*/
72063dd3a1aSKris 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)
7215e76c431SBarry Smith {
722406556e6SLois Curfman McInnes   /*
723406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
724406556e6SLois Curfman McInnes      minimization problem:
725406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
726406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
727406556e6SLois Curfman McInnes    */
7285ca10a37SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
729aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
730ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7316b5873e3SBarry Smith #endif
732dfbe8321SBarry Smith   PetscErrorCode ierr;
733a7cc72afSBarry Smith   PetscInt       count;
7344b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
735ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
7363c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7375e76c431SBarry Smith 
7383a40ed3dSBarry Smith   PetscFunctionBegin;
739d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
740a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7415e42470aSBarry Smith   alpha   = neP->alpha;
7425e42470aSBarry Smith   maxstep = neP->maxstep;
7435e42470aSBarry Smith   steptol = neP->steptol;
7445e76c431SBarry Smith 
7453505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74663b9fa5eSBarry Smith   if (*ynorm == 0.0) {
74763ba0a88SBarry Smith     ierr   = PetscLogInfo((snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
748b37302c6SLois Curfman McInnes     *gnorm = fnorm;
749b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
750b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
751a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
752ad922ac9SBarry Smith     goto theend2;
75394a9d846SBarry Smith   }
7545e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7555e42470aSBarry Smith     scale  = maxstep/(*ynorm);
756393d2d9aSLois Curfman McInnes     ierr   = VecScale(&scale,y);CHKERRQ(ierr);
7575e42470aSBarry Smith     *ynorm = maxstep;
7585e76c431SBarry Smith   }
759ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7605ca10a37SBarry Smith   minlambda = steptol/rellength;
761a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
762aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
763a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
764329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7655e42470aSBarry Smith #else
766a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7675e42470aSBarry Smith #endif
7685e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7695e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7705e42470aSBarry Smith 
771393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7725334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
77343e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
77463ba0a88SBarry Smith     ierr  = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
77543e71028SBarry Smith     ierr  = VecCopy(w,y);CHKERRQ(ierr);
77643e71028SBarry Smith     *flag = PETSC_FALSE;
77743e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77843e71028SBarry Smith     goto theend2;
77943e71028SBarry Smith   }
78019717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
781d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
78219717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
78319717074SBarry Smith   }
784d5e45103SBarry Smith   CHKERRQ(ierr);
785cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
78643e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7875d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
788393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
78963ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Using full step\n"));CHKERRQ(ierr);
790ad922ac9SBarry Smith     goto theend2;
7915e42470aSBarry Smith   }
7925e42470aSBarry Smith 
7935e42470aSBarry Smith   /* Fit points with quadratic */
794313b5538SBarry Smith   lambda = 1.0;
7955e42470aSBarry Smith   count = 1;
7965ca10a37SBarry Smith   while (PETSC_TRUE) {
7975e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
79863ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
79963ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
800f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
80143e71028SBarry Smith       *flag = PETSC_FALSE;
80243e71028SBarry Smith       break;
8035e42470aSBarry Smith     }
804a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
805329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
806329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
807329f5518SBarry Smith     else                         lambda     = lambdatemp;
808393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
8093505fcc1SBarry Smith     lambdaneg = -lambda;
810aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
8113505fcc1SBarry Smith     clambda   = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
8125e42470aSBarry Smith #else
8133505fcc1SBarry Smith     ierr      = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
8145e42470aSBarry Smith #endif
81543e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
81663ba0a88SBarry Smith       ierr  = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
81763ba0a88SBarry 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);
81843e71028SBarry Smith       ierr  = VecCopy(w,y);CHKERRQ(ierr);
819ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
82043e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
821ed98166eSMatthew Knepley       break;
822ed98166eSMatthew Knepley     }
82319717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
824d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
82519717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
82619717074SBarry Smith     }
827d5e45103SBarry Smith     CHKERRQ(ierr);
828cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82943e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8305ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
831393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
83263ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr);
83306259719SBarry Smith       break;
8345e42470aSBarry Smith     }
8355e42470aSBarry Smith     count++;
8365e42470aSBarry Smith   }
837ad922ac9SBarry Smith   theend2:
8388f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8393c632250SBarry Smith   if (neP->postcheckstep) {
8403c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8413c632250SBarry Smith     if (changed_y) {
8423c632250SBarry Smith       ierr = VecWAXPY(&mone,y,x,w);CHKERRQ(ierr);
8433c632250SBarry Smith     }
8443c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8453c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
846d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
84719717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
84819717074SBarry Smith       }
849d5e45103SBarry Smith       CHKERRQ(ierr);
8508f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8513c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8528f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8533c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8543c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8558f99978dSLois Curfman McInnes     }
8568f99978dSLois Curfman McInnes   }
857d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8583a40ed3dSBarry Smith   PetscFunctionReturn(0);
8595e76c431SBarry Smith }
8602343ba6eSBarry Smith 
86104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8624a2ae208SSatish Balay #undef __FUNCT__
8634a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
864c9e6a524SLois Curfman McInnes /*@C
86577c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
8664b27c08aSLois Curfman McInnes    by the method SNESLS.
8675e76c431SBarry Smith 
8685e76c431SBarry Smith    Input Parameters:
869c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
870af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
871c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8725e76c431SBarry Smith 
873fee21e36SBarry Smith    Collective on SNES
874fee21e36SBarry Smith 
875c4a48953SLois Curfman McInnes    Available Routines:
876c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
877f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
878f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
879af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8805e76c431SBarry Smith 
881c4a48953SLois Curfman McInnes     Options Database Keys:
8824b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8834b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8844b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8854b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8863304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8873304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
888c4a48953SLois Curfman McInnes 
8895e76c431SBarry Smith    Calling sequence of func:
890bd208895SLois Curfman McInnes .vb
891af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
892a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
893bd208895SLois Curfman McInnes .ve
8945e76c431SBarry Smith 
8955e76c431SBarry Smith     Input parameters for func:
896c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
897af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8985e76c431SBarry Smith .   x - current iterate
8995e76c431SBarry Smith .   f - residual evaluated at x
9003c632250SBarry Smith .   y - search direction
9015e76c431SBarry Smith .   w - work vector
902c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9035e76c431SBarry Smith 
9045e76c431SBarry Smith     Output parameters for func:
905c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9063c632250SBarry Smith .   w - new iterate
9075e76c431SBarry Smith .   gnorm - 2-norm of g
9085e76c431SBarry Smith .   ynorm - 2-norm of search length
909a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
910f59ffdedSLois Curfman McInnes 
91136851e7fSLois Curfman McInnes     Level: advanced
91236851e7fSLois Curfman McInnes 
913f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
914f59ffdedSLois Curfman McInnes 
9155d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
916*9c3ca13aSBarry Smith           SNESLineSearchSetPostCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams(), SNESLineSearchSetPreCheck()
9175e76c431SBarry Smith @*/
91863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9195e76c431SBarry Smith {
920a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
92182bf6240SBarry Smith 
9223a40ed3dSBarry Smith   PetscFunctionBegin;
923c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
92482bf6240SBarry Smith   if (f) {
925af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92682bf6240SBarry Smith   }
9273a40ed3dSBarry Smith   PetscFunctionReturn(0);
9285e76c431SBarry Smith }
9298e019c35SBarry Smith 
930a7cc72afSBarry 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*/
93104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
932fb2e594dSBarry Smith EXTERN_C_BEGIN
9334a2ae208SSatish Balay #undef __FUNCT__
9344a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
93563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
93682bf6240SBarry Smith {
93782bf6240SBarry Smith   PetscFunctionBegin;
9384b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9394b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
94082bf6240SBarry Smith   PetscFunctionReturn(0);
94182bf6240SBarry Smith }
942fb2e594dSBarry Smith EXTERN_C_END
94304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9444a2ae208SSatish Balay #undef __FUNCT__
9453c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
946c8dd0c0dSLois Curfman McInnes /*@C
9473c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9484b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
949c8dd0c0dSLois Curfman McInnes 
950c8dd0c0dSLois Curfman McInnes    Input Parameters:
951c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9523c632250SBarry Smith .  func - pointer to function
953c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
954c8dd0c0dSLois Curfman McInnes 
955c8dd0c0dSLois Curfman McInnes    Collective on SNES
956c8dd0c0dSLois Curfman McInnes 
957c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
958c8dd0c0dSLois Curfman McInnes .vb
9593c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
960c8dd0c0dSLois Curfman McInnes .ve
961b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
962b1ae27eaSLois Curfman McInnes    on failure.
963c8dd0c0dSLois Curfman McInnes 
964c8dd0c0dSLois Curfman McInnes    Input parameters for func:
965c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
966c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9673c632250SBarry Smith .  x - previous iterate
9683c632250SBarry Smith .  y - new search direction and length
9693c632250SBarry Smith -  w - current candidate iterate
970c8dd0c0dSLois Curfman McInnes 
971c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9723c632250SBarry Smith +  y - search direction (possibly changed)
9733c632250SBarry Smith .  w - current iterate (possibly modified)
9743c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9753c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
976c8dd0c0dSLois Curfman McInnes 
977c8dd0c0dSLois Curfman McInnes    Level: advanced
978c8dd0c0dSLois Curfman McInnes 
9799e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9809e247f21SBarry Smith 
9813c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9823c632250SBarry Smith 
9833c632250SBarry Smith    On input w = x + y
9843c632250SBarry Smith 
9859e247f21SBarry Smith    SNESNoLineSearch() and SNESNoLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
986b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
987ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9888f99978dSLois Curfman McInnes 
9899e247f21SBarry Smith    SNESQuadraticLineSearch() and SNESCubicLineSearch() (1) compute a candidate iterate u_{i+1} as well as a
990ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
991ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
992ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9939e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
994b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9958f99978dSLois Curfman McInnes 
996c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
997c8dd0c0dSLois Curfman McInnes 
998*9c3ca13aSBarry Smith .seealso: SNESSetLineSearch(), SNESLineSearchSetPreCheck()
999c8dd0c0dSLois Curfman McInnes @*/
10003c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
1001c8dd0c0dSLois Curfman McInnes {
10023c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
1003c8dd0c0dSLois Curfman McInnes 
1004c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10053c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1006c8dd0c0dSLois Curfman McInnes   if (f) {
1007c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1008c8dd0c0dSLois Curfman McInnes   }
1009c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1010c8dd0c0dSLois Curfman McInnes }
1011*9c3ca13aSBarry Smith #undef __FUNCT__
1012*9c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
1013*9c3ca13aSBarry Smith /*@C
1014*9c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
1015*9c3ca13aSBarry Smith 
1016*9c3ca13aSBarry Smith    Input Parameters:
1017*9c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
1018*9c3ca13aSBarry Smith .  func - pointer to function
1019*9c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
1020*9c3ca13aSBarry Smith 
1021*9c3ca13aSBarry Smith    Collective on SNES
1022*9c3ca13aSBarry Smith 
1023*9c3ca13aSBarry Smith    Calling sequence of func:
1024*9c3ca13aSBarry Smith .vb
1025*9c3ca13aSBarry Smith    int func (SNES snes, Vec x,Vec y,,void *checkctx, PetscTruth *changed_y)
1026*9c3ca13aSBarry Smith .ve
1027*9c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
1028*9c3ca13aSBarry Smith    on failure.
1029*9c3ca13aSBarry Smith 
1030*9c3ca13aSBarry Smith    Input parameters for func:
1031*9c3ca13aSBarry Smith +  snes - nonlinear context
1032*9c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
1033*9c3ca13aSBarry Smith .  x - previous iterate
1034*9c3ca13aSBarry Smith -  y - new search direction and length
1035*9c3ca13aSBarry Smith 
1036*9c3ca13aSBarry Smith    Output parameters for func:
1037*9c3ca13aSBarry Smith +  y - search direction (possibly changed)
1038*9c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
1039*9c3ca13aSBarry Smith 
1040*9c3ca13aSBarry Smith    Level: advanced
1041*9c3ca13aSBarry Smith 
1042*9c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
1043*9c3ca13aSBarry Smith 
1044*9c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
1045*9c3ca13aSBarry Smith 
1046*9c3ca13aSBarry Smith .seealso: SNESSetLineSearch(), SNESLineSearchSetPostCheck()
1047*9c3ca13aSBarry Smith @*/
1048*9c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
1049*9c3ca13aSBarry Smith {
1050*9c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
1051*9c3ca13aSBarry Smith 
1052*9c3ca13aSBarry Smith   PetscFunctionBegin;
1053*9c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1054*9c3ca13aSBarry Smith   if (f) {
1055*9c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1056*9c3ca13aSBarry Smith   }
1057*9c3ca13aSBarry Smith   PetscFunctionReturn(0);
1058*9c3ca13aSBarry Smith }
1059*9c3ca13aSBarry Smith 
1060c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1061*9c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
1062*9c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1063c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10644a2ae208SSatish Balay #undef __FUNCT__
10653c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
1066*9c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1067c8dd0c0dSLois Curfman McInnes {
1068c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10693c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10703c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1071c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1072c8dd0c0dSLois Curfman McInnes }
1073c8dd0c0dSLois Curfman McInnes EXTERN_C_END
1074*9c3ca13aSBarry Smith 
1075*9c3ca13aSBarry Smith EXTERN_C_BEGIN
1076*9c3ca13aSBarry Smith #undef __FUNCT__
1077*9c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
1078*9c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
1079*9c3ca13aSBarry Smith {
1080*9c3ca13aSBarry Smith   PetscFunctionBegin;
1081*9c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
1082*9c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
1083*9c3ca13aSBarry Smith   PetscFunctionReturn(0);
1084*9c3ca13aSBarry Smith }
1085*9c3ca13aSBarry Smith EXTERN_C_END
1086c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
108704d965bbSLois Curfman McInnes /*
10884b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
1089329e7e40SMatthew Knepley 
1090329e7e40SMatthew Knepley    Input Parameter:
1091329e7e40SMatthew Knepley .  snes - the SNES context
1092329e7e40SMatthew Knepley 
1093329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
1094329e7e40SMatthew Knepley */
1095329e7e40SMatthew Knepley #undef __FUNCT__
10964b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
10976849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
1098329e7e40SMatthew Knepley {
10994b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
1100329e7e40SMatthew Knepley 
1101329e7e40SMatthew Knepley   PetscFunctionBegin;
11024b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
11034b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
11044b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
11054b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
11064b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
1107329e7e40SMatthew Knepley   PetscFunctionReturn(0);
1108329e7e40SMatthew Knepley }
1109329e7e40SMatthew Knepley 
1110329e7e40SMatthew Knepley /*
11114b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
111204d965bbSLois Curfman McInnes 
111304d965bbSLois Curfman McInnes    Input Parameters:
111404d965bbSLois Curfman McInnes .  SNES - the SNES context
111504d965bbSLois Curfman McInnes .  viewer - visualization context
111604d965bbSLois Curfman McInnes 
111704d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
111804d965bbSLois Curfman McInnes */
11194a2ae208SSatish Balay #undef __FUNCT__
11204b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11216849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1122a935fc98SLois Curfman McInnes {
11234b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11242fc52814SBarry Smith   const char     *cstr;
1125dfbe8321SBarry Smith   PetscErrorCode ierr;
112632077d6dSBarry Smith   PetscTruth     iascii;
1127a935fc98SLois Curfman McInnes 
11283a40ed3dSBarry Smith   PetscFunctionBegin;
112932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
113032077d6dSBarry Smith   if (iascii) {
113119bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
113219bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
113319bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
113419bcc07fSBarry Smith     else                                                cstr = "unknown";
1135b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1136b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11375cd90555SBarry Smith   } else {
113879a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
113919bcc07fSBarry Smith   }
11403a40ed3dSBarry Smith   PetscFunctionReturn(0);
1141a935fc98SLois Curfman McInnes }
114204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
114304d965bbSLois Curfman McInnes /*
11444b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
114504d965bbSLois Curfman McInnes 
114604d965bbSLois Curfman McInnes    Input Parameter:
114704d965bbSLois Curfman McInnes .  snes - the SNES context
114804d965bbSLois Curfman McInnes 
114904d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
115004d965bbSLois Curfman McInnes */
11514a2ae208SSatish Balay #undef __FUNCT__
11524b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11536849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11545e42470aSBarry Smith {
11554b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1156e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1157dfbe8321SBarry Smith   PetscErrorCode ierr;
1158a7cc72afSBarry Smith   PetscInt       indx;
1159f1af5d2fSBarry Smith   PetscTruth     flg;
11605e42470aSBarry Smith 
11613a40ed3dSBarry Smith   PetscFunctionBegin;
1162b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11634b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11644b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11654b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1166186905e3SBarry Smith 
116722e36657SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
116825cdf11fSBarry Smith     if (flg) {
116922e36657SBarry Smith       switch (indx) {
1170b49fd9e1SBarry Smith       case 0:
1171af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1172b49fd9e1SBarry Smith         break;
1173b49fd9e1SBarry Smith       case 1:
1174af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1175b49fd9e1SBarry Smith         break;
1176b49fd9e1SBarry Smith       case 2:
1177af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1178b49fd9e1SBarry Smith         break;
1179b49fd9e1SBarry Smith       case 3:
1180af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1181b49fd9e1SBarry Smith         break;
11825e42470aSBarry Smith       }
11835e42470aSBarry Smith     }
1184b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11853a40ed3dSBarry Smith   PetscFunctionReturn(0);
11865e42470aSBarry Smith }
118704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1188ebe3b25bSBarry Smith /*MC
1189ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
119004d965bbSLois Curfman McInnes 
1191ebe3b25bSBarry Smith    Options Database:
1192ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1193ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1194ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1195ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1196ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1197ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
119804d965bbSLois Curfman McInnes 
1199ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
120004d965bbSLois Curfman McInnes 
1201ee3001cbSBarry Smith    Level: beginner
1202ee3001cbSBarry Smith 
1203ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
12043c632250SBarry Smith            SNESLineSearchSetPostCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1205*9c3ca13aSBarry Smith           SNESSetLineSearch(), SNESNoLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1206ebe3b25bSBarry Smith 
1207ebe3b25bSBarry Smith M*/
1208fb2e594dSBarry Smith EXTERN_C_BEGIN
12094a2ae208SSatish Balay #undef __FUNCT__
12104b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
121163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
12125e42470aSBarry Smith {
1213dfbe8321SBarry Smith   PetscErrorCode ierr;
12144b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12155e42470aSBarry Smith 
12163a40ed3dSBarry Smith   PetscFunctionBegin;
12174b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
12184b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
12194b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
12204b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
12214b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
12224b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
12234b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
12245baf8537SBarry Smith   snes->nwork           = 0;
12255e42470aSBarry Smith 
12264b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
122752e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
12285e42470aSBarry Smith   snes->data    	= (void*)neP;
12295e42470aSBarry Smith   neP->alpha		= 1.e-4;
12305e42470aSBarry Smith   neP->maxstep		= 1.e8;
12315e42470aSBarry Smith   neP->steptol		= 1.e-12;
12325e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1233c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12343c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12353c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12363c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12373c632250SBarry Smith   neP->precheck         = PETSC_NULL;
123882bf6240SBarry Smith 
12395d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
12403c632250SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
1241*9c3ca13aSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
124282bf6240SBarry Smith 
12433a40ed3dSBarry Smith   PetscFunctionReturn(0);
12445e42470aSBarry Smith }
1245fb2e594dSBarry Smith EXTERN_C_END
124682bf6240SBarry Smith 
124782bf6240SBarry Smith 
124882bf6240SBarry Smith 
1249