xref: /petsc/src/snes/impls/ls/ls.c (revision abb0e124a6c24cb1ba78480aa3c5e0caea43b62b)
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 
33*abb0e124SMatthew Knepley     ierr = VecSetRandom(W,PETSC_NULL);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);
622dcb1b2aSMatthew Knepley     ierr = VecAXPY(W1,mone,F);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 
1859c3ca13aSBarry Smith     if (neP->precheckstep) {
1869c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
1879c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
1889c3ca13aSBarry Smith     }
1899c3ca13aSBarry 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;
275e74804ceSBarry Smith   if (!snes->work) {
27681b6cf68SLois Curfman McInnes     snes->nwork = 4;
277d7e8b826SBarry Smith     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
278efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
27981b6cf68SLois Curfman McInnes     snes->vec_sol_update_always = snes->work[3];
280e74804ceSBarry Smith   }
2813a40ed3dSBarry Smith   PetscFunctionReturn(0);
2825e76c431SBarry Smith }
28304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28404d965bbSLois Curfman McInnes /*
2854b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2864b27c08aSLois Curfman McInnes    with SNESCreate_LS().
28704d965bbSLois Curfman McInnes 
28804d965bbSLois Curfman McInnes    Input Parameter:
28904d965bbSLois Curfman McInnes .  snes - the SNES context
29004d965bbSLois Curfman McInnes 
291de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29204d965bbSLois Curfman McInnes  */
2934a2ae208SSatish Balay #undef __FUNCT__
2944b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
295dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2965e76c431SBarry Smith {
297dfbe8321SBarry Smith   PetscErrorCode ierr;
2983a40ed3dSBarry Smith 
2993a40ed3dSBarry Smith   PetscFunctionBegin;
3005baf8537SBarry Smith   if (snes->nwork) {
3014b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
30221c89e3eSBarry Smith   }
3035c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3043a40ed3dSBarry Smith   PetscFunctionReturn(0);
3055e76c431SBarry Smith }
30604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3074a2ae208SSatish Balay #undef __FUNCT__
30817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
30982bf6240SBarry Smith 
3104b828684SBarry Smith /*@C
31117bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3125e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3135e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3145e76c431SBarry Smith 
315c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
316c7afd0dbSLois Curfman McInnes 
3175e76c431SBarry Smith    Input Parameters:
318c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
319af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3205e76c431SBarry Smith .  x - current iterate
3215e76c431SBarry Smith .  f - residual evaluated at x
3223c632250SBarry Smith .  y - search direction
3235e76c431SBarry Smith .  w - work vector
324c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3255e76c431SBarry Smith 
326c4a48953SLois Curfman McInnes    Output Parameters:
327c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3283c632250SBarry Smith .  w - new iterate
3295e76c431SBarry Smith .  gnorm - 2-norm of g
3305e76c431SBarry Smith .  ynorm - 2-norm of search length
331a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
332fee21e36SBarry Smith 
333c4a48953SLois Curfman McInnes    Options Database Key:
33417bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
335c4a48953SLois Curfman McInnes 
33636851e7fSLois Curfman McInnes    Level: advanced
33736851e7fSLois Curfman McInnes 
33828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33928ae5a21SLois Curfman McInnes 
34017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
34117bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3425e76c431SBarry Smith @*/
34317bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
3445e76c431SBarry Smith {
345dfbe8321SBarry Smith   PetscErrorCode ierr;
346ea709b57SSatish Balay   PetscScalar    mone = -1.0;
3474b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3483c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3495334005bSBarry Smith 
3503a40ed3dSBarry Smith   PetscFunctionBegin;
351a7cc72afSBarry Smith   *flag = PETSC_TRUE;
352d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
353a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
3542dcb1b2aSMatthew Knepley   ierr = VecWAXPY(w,mone,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3553c632250SBarry Smith   if (neP->postcheckstep) {
3563c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3578f99978dSLois Curfman McInnes   }
3583c632250SBarry Smith   if (changed_y) {
3592dcb1b2aSMatthew Knepley     ierr = VecWAXPY(w,mone,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3603c632250SBarry Smith   }
3613c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
362d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
36319717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
36419717074SBarry Smith   }
365d5e45103SBarry Smith   CHKERRQ(ierr);
366d5e45103SBarry Smith 
367a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
36843e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
369d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3703a40ed3dSBarry Smith   PetscFunctionReturn(0);
3715e76c431SBarry Smith }
37204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3734a2ae208SSatish Balay #undef __FUNCT__
37417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
37582bf6240SBarry Smith 
37629e0b56fSBarry Smith /*@C
37717bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
37829e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
37929e0b56fSBarry Smith    even compute the norm of the function or search direction; this
38029e0b56fSBarry Smith    is intended only when you know the full step is fine and are
38129e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
382c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
383c7afd0dbSLois Curfman McInnes 
384c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
38529e0b56fSBarry Smith 
38629e0b56fSBarry Smith    Input Parameters:
387c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
388af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
38929e0b56fSBarry Smith .  x - current iterate
39029e0b56fSBarry Smith .  f - residual evaluated at x
3913c632250SBarry Smith .  y - search direction
39229e0b56fSBarry Smith .  w - work vector
393c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
39429e0b56fSBarry Smith 
39529e0b56fSBarry Smith    Output Parameters:
396c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3973c632250SBarry Smith .  w - new iterate
39829e0b56fSBarry Smith .  gnorm - not changed
39929e0b56fSBarry Smith .  ynorm - not changed
400a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
401fee21e36SBarry Smith 
40229e0b56fSBarry Smith    Options Database Key:
40317bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
40429e0b56fSBarry Smith 
4058cbba510SBarry Smith    Notes:
40617bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
407ea56f5baSLois Curfman McInnes    either the options
408ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
409ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
410329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
411329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
412329f5518SBarry Smith 
413329f5518SBarry Smith    During the final iteration this will not evaluate the function at
414329f5518SBarry Smith    the solution point. This is to save a function evaluation while
415329f5518SBarry Smith    using pseudo-timestepping.
4168cbba510SBarry Smith 
417ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
418ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
419ea56f5baSLois Curfman McInnes    correct, since they are not computed.
420ea56f5baSLois Curfman McInnes 
421ea56f5baSLois Curfman McInnes    Level: advanced
4228cbba510SBarry Smith 
42329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
42429e0b56fSBarry Smith 
42517bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
42617bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
42729e0b56fSBarry Smith @*/
42817bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
42929e0b56fSBarry Smith {
430dfbe8321SBarry Smith   PetscErrorCode ierr;
431ea709b57SSatish Balay   PetscScalar    mone = -1.0;
4324b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4333c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
43429e0b56fSBarry Smith 
4353a40ed3dSBarry Smith   PetscFunctionBegin;
436a7cc72afSBarry Smith   *flag = PETSC_TRUE;
437d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4382dcb1b2aSMatthew Knepley   ierr = VecWAXPY(w,mone,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4393c632250SBarry Smith   if (neP->postcheckstep) {
4403c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4413c632250SBarry Smith   }
4423c632250SBarry Smith   if (changed_y) {
4432dcb1b2aSMatthew Knepley     ierr = VecWAXPY(w,mone,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4448f99978dSLois Curfman McInnes   }
445329f5518SBarry Smith 
446329f5518SBarry Smith   /* don't evaluate function the last time through */
447329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4483c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
449d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
45019717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
45119717074SBarry Smith     }
452d5e45103SBarry Smith     CHKERRQ(ierr);
453329f5518SBarry Smith   }
454d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4553a40ed3dSBarry Smith   PetscFunctionReturn(0);
45629e0b56fSBarry Smith }
45704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4584a2ae208SSatish Balay #undef __FUNCT__
45917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4604b828684SBarry Smith /*@C
46117bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4625e76c431SBarry Smith 
463c7afd0dbSLois Curfman McInnes    Collective on SNES
464c7afd0dbSLois Curfman McInnes 
4655e76c431SBarry Smith    Input Parameters:
466c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
467af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4685e76c431SBarry Smith .  x - current iterate
4695e76c431SBarry Smith .  f - residual evaluated at x
4703c632250SBarry Smith .  y - search direction
4715e76c431SBarry Smith .  w - work vector
472c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4735e76c431SBarry Smith 
474393d2d9aSLois Curfman McInnes    Output Parameters:
475c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4763c632250SBarry Smith .  w - new iterate
4775e76c431SBarry Smith .  gnorm - 2-norm of g
4785e76c431SBarry Smith .  ynorm - 2-norm of search length
479a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
480fee21e36SBarry Smith 
481c4a48953SLois Curfman McInnes    Options Database Key:
48217bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
483c4a48953SLois Curfman McInnes 
4845e76c431SBarry Smith    Notes:
4855e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4865e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4875e76c431SBarry Smith 
48836851e7fSLois Curfman McInnes    Level: advanced
48936851e7fSLois Curfman McInnes 
49028ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49128ae5a21SLois Curfman McInnes 
49217bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
4935e76c431SBarry Smith @*/
49417bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
4955e76c431SBarry Smith {
496406556e6SLois Curfman McInnes   /*
497406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
498406556e6SLois Curfman McInnes      minimization problem:
499406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
500406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
501406556e6SLois Curfman McInnes    */
502406556e6SLois Curfman McInnes 
5035ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
504329f5518SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
505aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
506ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5076b5873e3SBarry Smith #endif
508dfbe8321SBarry Smith   PetscErrorCode ierr;
509a7cc72afSBarry Smith   PetscInt       count;
5104b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
511ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
5123c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5135e76c431SBarry Smith 
5143a40ed3dSBarry Smith   PetscFunctionBegin;
515d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
516a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5175e76c431SBarry Smith   alpha   = neP->alpha;
5185e76c431SBarry Smith   maxstep = neP->maxstep;
5195e76c431SBarry Smith   steptol = neP->steptol;
5205e76c431SBarry Smith 
521cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5229e247f21SBarry Smith   if (!*ynorm) {
52317bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Search direction and size is 0\n"));CHKERRQ(ierr);
524a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5253c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
526a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
527a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
528ad922ac9SBarry Smith     goto theend1;
52994a9d846SBarry Smith   }
5305e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5315e42470aSBarry Smith     scale = maxstep/(*ynorm);
532aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
53317bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr);
5346b5873e3SBarry Smith #else
53517bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr);
5366b5873e3SBarry Smith #endif
5372dcb1b2aSMatthew Knepley     ierr = VecScale(y,scale);CHKERRQ(ierr);
5385e76c431SBarry Smith     *ynorm = maxstep;
5395e76c431SBarry Smith   }
540ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5415ca10a37SBarry Smith   minlambda = steptol/rellength;
542a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
543aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
544a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
545329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5465e42470aSBarry Smith #else
547a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5485e42470aSBarry Smith #endif
5495e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5505e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5515e76c431SBarry Smith 
552393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5532dcb1b2aSMatthew Knepley   ierr = VecAYPX(w,mone,x);CHKERRQ(ierr);
55443e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
55517bae607SBarry Smith     ierr  = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
55643e71028SBarry Smith     *flag = PETSC_FALSE;
55743e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
55843e71028SBarry Smith     goto theend1;
55943e71028SBarry Smith   }
56019717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
561d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
56219717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
56319717074SBarry Smith   }
564d5e45103SBarry Smith   CHKERRQ(ierr);
565313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56643e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5675d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
56817bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Using full step\n"));CHKERRQ(ierr);
569ad922ac9SBarry Smith     goto theend1;
5705e76c431SBarry Smith   }
5715e76c431SBarry Smith 
5725e76c431SBarry Smith   /* Fit points with quadratic */
573313b5538SBarry Smith   lambda     = 1.0;
574a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5755e76c431SBarry Smith   lambdaprev = lambda;
5765e76c431SBarry Smith   gnormprev  = *gnorm;
577329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
578ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
579ddd12767SBarry Smith   else                         lambda = lambdatemp;
580393d2d9aSLois Curfman McInnes   ierr      = VecCopy(x,w);CHKERRQ(ierr);
581ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
582aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5832dcb1b2aSMatthew Knepley   clambda   = lambdaneg; ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr);
5845e42470aSBarry Smith #else
5852dcb1b2aSMatthew Knepley   ierr      = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr);
5865e42470aSBarry Smith #endif
58743e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
58817bae607SBarry Smith     ierr  = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr);
58943e71028SBarry Smith     *flag = PETSC_FALSE;
59043e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
59143e71028SBarry Smith     goto theend1;
59243e71028SBarry Smith   }
59319717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
594d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59519717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59619717074SBarry Smith   }
597d5e45103SBarry Smith   CHKERRQ(ierr);
598cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
59943e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6005ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
60117bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
602ad922ac9SBarry Smith     goto theend1;
6035e76c431SBarry Smith   }
6045e76c431SBarry Smith 
6055e76c431SBarry Smith   /* Fit points with cubic */
6065e76c431SBarry Smith   count = 1;
6078229bfc2SMatthew Knepley   while (count < 10000) {
6085e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
60917bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
61017bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
61143e71028SBarry Smith       *flag = PETSC_FALSE;
61243e71028SBarry Smith       break;
6135e76c431SBarry Smith     }
6145d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6155d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
616ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6172b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6185e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6195e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6205e76c431SBarry Smith     if (a == 0.0) {
6215e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6225e76c431SBarry Smith     } else {
6235e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6245e76c431SBarry Smith     }
6255e76c431SBarry Smith     lambdaprev = lambda;
6265e76c431SBarry Smith     gnormprev  = *gnorm;
627329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
628329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6295e76c431SBarry Smith     else                         lambda     = lambdatemp;
630393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
631ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
632aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
633ea4d3ed3SLois Curfman McInnes     clambda   = lambdaneg;
6342dcb1b2aSMatthew Knepley     ierr      = VecAXPY(w,clambda,y);CHKERRQ(ierr);
6355e42470aSBarry Smith #else
6362dcb1b2aSMatthew Knepley     ierr      = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr);
6375e42470aSBarry Smith #endif
63843e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
63917bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
64017bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
641ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
64243e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
643ed98166eSMatthew Knepley       break;
644ed98166eSMatthew Knepley     }
64519717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
646d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64719717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64819717074SBarry Smith     }
649d5e45103SBarry Smith     CHKERRQ(ierr);
650cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
65143e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6525ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
65317bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
65443e71028SBarry Smith       break;
6552b022350SLois Curfman McInnes     } else {
65617bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda));CHKERRQ(ierr);
6575e76c431SBarry Smith     }
6585e76c431SBarry Smith     count++;
6595e76c431SBarry Smith   }
6608229bfc2SMatthew Knepley   if (count >= 10000) {
661cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6628229bfc2SMatthew Knepley   }
663ad922ac9SBarry Smith   theend1:
6648f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6653c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6663c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6673c632250SBarry Smith     if (changed_y) {
6682dcb1b2aSMatthew Knepley       ierr = VecWAXPY(w,mone,y,x);CHKERRQ(ierr);
6693c632250SBarry Smith     }
6703c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6713c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
672d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
67319717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
67419717074SBarry Smith       }
675d5e45103SBarry Smith       CHKERRQ(ierr);
6768f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67743e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6783c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6798f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6803c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6818f99978dSLois Curfman McInnes     }
6828f99978dSLois Curfman McInnes   }
683d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6843a40ed3dSBarry Smith   PetscFunctionReturn(0);
6855e76c431SBarry Smith }
68604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6874a2ae208SSatish Balay #undef __FUNCT__
68817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6894b828684SBarry Smith /*@C
69017bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6915e76c431SBarry Smith 
692c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
693c7afd0dbSLois Curfman McInnes 
6945e42470aSBarry Smith    Input Parameters:
695c7afd0dbSLois Curfman McInnes +  snes - the SNES context
696af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6975e42470aSBarry Smith .  x - current iterate
6985e42470aSBarry Smith .  f - residual evaluated at x
6993c632250SBarry Smith .  y - search direction
7005e42470aSBarry Smith .  w - work vector
701c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
7025e42470aSBarry Smith 
703c4a48953SLois Curfman McInnes    Output Parameters:
704c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
7053c632250SBarry Smith .  w - new iterate
7065e42470aSBarry Smith .  gnorm - 2-norm of g
7075e42470aSBarry Smith .  ynorm - 2-norm of search length
708a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
709fee21e36SBarry Smith 
710c4a48953SLois Curfman McInnes    Options Database Key:
71117bae607SBarry Smith .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
712c4a48953SLois Curfman McInnes 
7135e42470aSBarry Smith    Notes:
71417bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7155e42470aSBarry Smith 
71636851e7fSLois Curfman McInnes    Level: advanced
71736851e7fSLois Curfman McInnes 
718f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
719f59ffdedSLois Curfman McInnes 
72017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7215e42470aSBarry Smith @*/
72217bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
7235e76c431SBarry Smith {
724406556e6SLois Curfman McInnes   /*
725406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
726406556e6SLois Curfman McInnes      minimization problem:
727406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
728406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
729406556e6SLois Curfman McInnes    */
7305ca10a37SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
731aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
732ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7336b5873e3SBarry Smith #endif
734dfbe8321SBarry Smith   PetscErrorCode ierr;
735a7cc72afSBarry Smith   PetscInt       count;
7364b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
737ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
7383c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7395e76c431SBarry Smith 
7403a40ed3dSBarry Smith   PetscFunctionBegin;
741d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
742a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7435e42470aSBarry Smith   alpha   = neP->alpha;
7445e42470aSBarry Smith   maxstep = neP->maxstep;
7455e42470aSBarry Smith   steptol = neP->steptol;
7465e76c431SBarry Smith 
7473505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74863b9fa5eSBarry Smith   if (*ynorm == 0.0) {
74917bae607SBarry Smith     ierr   = PetscLogInfo((snes,"SNESLineSearchQuadratic: Search direction and size is 0\n"));CHKERRQ(ierr);
750b37302c6SLois Curfman McInnes     *gnorm = fnorm;
751b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
752b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
753a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
754ad922ac9SBarry Smith     goto theend2;
75594a9d846SBarry Smith   }
7565e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7575e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7582dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7595e42470aSBarry Smith     *ynorm = maxstep;
7605e76c431SBarry Smith   }
761ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7625ca10a37SBarry Smith   minlambda = steptol/rellength;
763a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
764aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
765a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
766329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7675e42470aSBarry Smith #else
768a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7695e42470aSBarry Smith #endif
7705e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7715e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7725e42470aSBarry Smith 
773393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7742dcb1b2aSMatthew Knepley   ierr = VecAYPX(w,mone,x);CHKERRQ(ierr);
77543e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
77617bae607SBarry Smith     ierr  = PetscLogInfo((snes,"SNESLineSearchQuadratic:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
77743e71028SBarry Smith     ierr  = VecCopy(w,y);CHKERRQ(ierr);
77843e71028SBarry Smith     *flag = PETSC_FALSE;
77943e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
78043e71028SBarry Smith     goto theend2;
78143e71028SBarry Smith   }
78219717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
783d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
78419717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
78519717074SBarry Smith   }
786d5e45103SBarry Smith   CHKERRQ(ierr);
787cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
78843e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7895d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
790393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
79117bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic: Using full step\n"));CHKERRQ(ierr);
792ad922ac9SBarry Smith     goto theend2;
7935e42470aSBarry Smith   }
7945e42470aSBarry Smith 
7955e42470aSBarry Smith   /* Fit points with quadratic */
796313b5538SBarry Smith   lambda = 1.0;
7975e42470aSBarry Smith   count = 1;
7985ca10a37SBarry Smith   while (PETSC_TRUE) {
7995e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
80017bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
80117bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
802f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
80343e71028SBarry Smith       *flag = PETSC_FALSE;
80443e71028SBarry Smith       break;
8055e42470aSBarry Smith     }
806a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
807329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
808329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
809329f5518SBarry Smith     else                         lambda     = lambdatemp;
810393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
8113505fcc1SBarry Smith     lambdaneg = -lambda;
812aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
8132dcb1b2aSMatthew Knepley     clambda   = lambdaneg; ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr);
8145e42470aSBarry Smith #else
8152dcb1b2aSMatthew Knepley     ierr      = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr);
8165e42470aSBarry Smith #endif
81743e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
81817bae607SBarry Smith       ierr  = PetscLogInfo((snes,"SNESLineSearchQuadratic:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
81917bae607SBarry Smith       ierr  = PetscLogInfo((snes,"SNESLineSearchQuadratic:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
82043e71028SBarry Smith       ierr  = VecCopy(w,y);CHKERRQ(ierr);
821ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
82243e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
823ed98166eSMatthew Knepley       break;
824ed98166eSMatthew Knepley     }
82519717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
826d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
82719717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
82819717074SBarry Smith     }
829d5e45103SBarry Smith     CHKERRQ(ierr);
830cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
83143e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8325ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
833393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
83417bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr);
83506259719SBarry Smith       break;
8365e42470aSBarry Smith     }
8375e42470aSBarry Smith     count++;
8385e42470aSBarry Smith   }
839ad922ac9SBarry Smith   theend2:
8408f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8413c632250SBarry Smith   if (neP->postcheckstep) {
8423c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8433c632250SBarry Smith     if (changed_y) {
8442dcb1b2aSMatthew Knepley       ierr = VecWAXPY(w,mone,y,x);CHKERRQ(ierr);
8453c632250SBarry Smith     }
8463c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8473c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
848d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
84919717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
85019717074SBarry Smith       }
851d5e45103SBarry Smith       CHKERRQ(ierr);
8528f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8533c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8548f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8553c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8563c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8578f99978dSLois Curfman McInnes     }
8588f99978dSLois Curfman McInnes   }
859d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8603a40ed3dSBarry Smith   PetscFunctionReturn(0);
8615e76c431SBarry Smith }
8622343ba6eSBarry Smith 
86304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8644a2ae208SSatish Balay #undef __FUNCT__
86517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
866c9e6a524SLois Curfman McInnes /*@C
86717bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8684b27c08aSLois Curfman McInnes    by the method SNESLS.
8695e76c431SBarry Smith 
8705e76c431SBarry Smith    Input Parameters:
871c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
872af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
873c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8745e76c431SBarry Smith 
875fee21e36SBarry Smith    Collective on SNES
876fee21e36SBarry Smith 
877c4a48953SLois Curfman McInnes    Available Routines:
87817bae607SBarry Smith +  SNESLineSearchCubic() - default line search
87917bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
88017bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
88117bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8825e76c431SBarry Smith 
883c4a48953SLois Curfman McInnes     Options Database Keys:
8844b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8854b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8864b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8874b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8883304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8893304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
890c4a48953SLois Curfman McInnes 
8915e76c431SBarry Smith    Calling sequence of func:
892bd208895SLois Curfman McInnes .vb
893af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
894a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
895bd208895SLois Curfman McInnes .ve
8965e76c431SBarry Smith 
8975e76c431SBarry Smith     Input parameters for func:
898c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
899af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
9005e76c431SBarry Smith .   x - current iterate
9015e76c431SBarry Smith .   f - residual evaluated at x
9023c632250SBarry Smith .   y - search direction
9035e76c431SBarry Smith .   w - work vector
904c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9055e76c431SBarry Smith 
9065e76c431SBarry Smith     Output parameters for func:
907c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9083c632250SBarry Smith .   w - new iterate
9095e76c431SBarry Smith .   gnorm - 2-norm of g
9105e76c431SBarry Smith .   ynorm - 2-norm of search length
911a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
912f59ffdedSLois Curfman McInnes 
91336851e7fSLois Curfman McInnes     Level: advanced
91436851e7fSLois Curfman McInnes 
915f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
916f59ffdedSLois Curfman McInnes 
91717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
91817bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9195e76c431SBarry Smith @*/
92017bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9215e76c431SBarry Smith {
922a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
92382bf6240SBarry Smith 
9243a40ed3dSBarry Smith   PetscFunctionBegin;
92517bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
92682bf6240SBarry Smith   if (f) {
927af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92882bf6240SBarry Smith   }
9293a40ed3dSBarry Smith   PetscFunctionReturn(0);
9305e76c431SBarry Smith }
9318e019c35SBarry Smith 
932a7cc72afSBarry 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*/
93304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
934fb2e594dSBarry Smith EXTERN_C_BEGIN
9354a2ae208SSatish Balay #undef __FUNCT__
93617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
93717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
93882bf6240SBarry Smith {
93982bf6240SBarry Smith   PetscFunctionBegin;
9404b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9414b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
94282bf6240SBarry Smith   PetscFunctionReturn(0);
94382bf6240SBarry Smith }
944fb2e594dSBarry Smith EXTERN_C_END
94504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9464a2ae208SSatish Balay #undef __FUNCT__
9473c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
948c8dd0c0dSLois Curfman McInnes /*@C
9493c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9504b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
951c8dd0c0dSLois Curfman McInnes 
952c8dd0c0dSLois Curfman McInnes    Input Parameters:
953c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9543c632250SBarry Smith .  func - pointer to function
955c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
956c8dd0c0dSLois Curfman McInnes 
957c8dd0c0dSLois Curfman McInnes    Collective on SNES
958c8dd0c0dSLois Curfman McInnes 
959c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
960c8dd0c0dSLois Curfman McInnes .vb
9613c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
962c8dd0c0dSLois Curfman McInnes .ve
963b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
964b1ae27eaSLois Curfman McInnes    on failure.
965c8dd0c0dSLois Curfman McInnes 
966c8dd0c0dSLois Curfman McInnes    Input parameters for func:
967c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
968c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9693c632250SBarry Smith .  x - previous iterate
9703c632250SBarry Smith .  y - new search direction and length
9713c632250SBarry Smith -  w - current candidate iterate
972c8dd0c0dSLois Curfman McInnes 
973c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9743c632250SBarry Smith +  y - search direction (possibly changed)
9753c632250SBarry Smith .  w - current iterate (possibly modified)
9763c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9773c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
978c8dd0c0dSLois Curfman McInnes 
979c8dd0c0dSLois Curfman McInnes    Level: advanced
980c8dd0c0dSLois Curfman McInnes 
9819e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9829e247f21SBarry Smith 
9833c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9843c632250SBarry Smith 
9853c632250SBarry Smith    On input w = x + y
9863c632250SBarry Smith 
98717bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
988b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
989ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9908f99978dSLois Curfman McInnes 
99117bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
992ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
993ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
994ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9959e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
996b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9978f99978dSLois Curfman McInnes 
998c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
999c8dd0c0dSLois Curfman McInnes 
100017bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1001c8dd0c0dSLois Curfman McInnes @*/
10023c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
1003c8dd0c0dSLois Curfman McInnes {
10043c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
1005c8dd0c0dSLois Curfman McInnes 
1006c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10073c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1008c8dd0c0dSLois Curfman McInnes   if (f) {
1009c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1010c8dd0c0dSLois Curfman McInnes   }
1011c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1012c8dd0c0dSLois Curfman McInnes }
10139c3ca13aSBarry Smith #undef __FUNCT__
10149c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10159c3ca13aSBarry Smith /*@C
10169c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10179c3ca13aSBarry Smith 
10189c3ca13aSBarry Smith    Input Parameters:
10199c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10209c3ca13aSBarry Smith .  func - pointer to function
10219c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10229c3ca13aSBarry Smith 
10239c3ca13aSBarry Smith    Collective on SNES
10249c3ca13aSBarry Smith 
10259c3ca13aSBarry Smith    Calling sequence of func:
10269c3ca13aSBarry Smith .vb
10279c3ca13aSBarry Smith    int func (SNES snes, Vec x,Vec y,,void *checkctx, PetscTruth *changed_y)
10289c3ca13aSBarry Smith .ve
10299c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10309c3ca13aSBarry Smith    on failure.
10319c3ca13aSBarry Smith 
10329c3ca13aSBarry Smith    Input parameters for func:
10339c3ca13aSBarry Smith +  snes - nonlinear context
10349c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10359c3ca13aSBarry Smith .  x - previous iterate
10369c3ca13aSBarry Smith -  y - new search direction and length
10379c3ca13aSBarry Smith 
10389c3ca13aSBarry Smith    Output parameters for func:
10399c3ca13aSBarry Smith +  y - search direction (possibly changed)
10409c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10419c3ca13aSBarry Smith 
10429c3ca13aSBarry Smith    Level: advanced
10439c3ca13aSBarry Smith 
10449c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10459c3ca13aSBarry Smith 
10469c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10479c3ca13aSBarry Smith 
104817bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck()
10499c3ca13aSBarry Smith @*/
10509c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10519c3ca13aSBarry Smith {
10529c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10539c3ca13aSBarry Smith 
10549c3ca13aSBarry Smith   PetscFunctionBegin;
10559c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10569c3ca13aSBarry Smith   if (f) {
10579c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10589c3ca13aSBarry Smith   }
10599c3ca13aSBarry Smith   PetscFunctionReturn(0);
10609c3ca13aSBarry Smith }
10619c3ca13aSBarry Smith 
1062c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10639c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10649c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1065c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10664a2ae208SSatish Balay #undef __FUNCT__
10673c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10689c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1069c8dd0c0dSLois Curfman McInnes {
1070c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10713c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10723c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1073c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1074c8dd0c0dSLois Curfman McInnes }
1075c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10769c3ca13aSBarry Smith 
10779c3ca13aSBarry Smith EXTERN_C_BEGIN
10789c3ca13aSBarry Smith #undef __FUNCT__
10799c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10809c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10819c3ca13aSBarry Smith {
10829c3ca13aSBarry Smith   PetscFunctionBegin;
10839c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10849c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10859c3ca13aSBarry Smith   PetscFunctionReturn(0);
10869c3ca13aSBarry Smith }
10879c3ca13aSBarry Smith EXTERN_C_END
1088c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
108904d965bbSLois Curfman McInnes /*
10904b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
1091329e7e40SMatthew Knepley 
1092329e7e40SMatthew Knepley    Input Parameter:
1093329e7e40SMatthew Knepley .  snes - the SNES context
1094329e7e40SMatthew Knepley 
1095329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
1096329e7e40SMatthew Knepley */
1097329e7e40SMatthew Knepley #undef __FUNCT__
10984b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
10996849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
1100329e7e40SMatthew Knepley {
11014b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
1102329e7e40SMatthew Knepley 
1103329e7e40SMatthew Knepley   PetscFunctionBegin;
11044b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
11054b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
11064b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
11074b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
11084b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
1109329e7e40SMatthew Knepley   PetscFunctionReturn(0);
1110329e7e40SMatthew Knepley }
1111329e7e40SMatthew Knepley 
1112329e7e40SMatthew Knepley /*
11134b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
111404d965bbSLois Curfman McInnes 
111504d965bbSLois Curfman McInnes    Input Parameters:
111604d965bbSLois Curfman McInnes .  SNES - the SNES context
111704d965bbSLois Curfman McInnes .  viewer - visualization context
111804d965bbSLois Curfman McInnes 
111904d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
112004d965bbSLois Curfman McInnes */
11214a2ae208SSatish Balay #undef __FUNCT__
11224b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11236849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1124a935fc98SLois Curfman McInnes {
11254b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11262fc52814SBarry Smith   const char     *cstr;
1127dfbe8321SBarry Smith   PetscErrorCode ierr;
112832077d6dSBarry Smith   PetscTruth     iascii;
1129a935fc98SLois Curfman McInnes 
11303a40ed3dSBarry Smith   PetscFunctionBegin;
113132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
113232077d6dSBarry Smith   if (iascii) {
113317bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
113417bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
113517bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
113619bcc07fSBarry Smith     else                                                cstr = "unknown";
1137b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1138b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11395cd90555SBarry Smith   } else {
114079a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
114119bcc07fSBarry Smith   }
11423a40ed3dSBarry Smith   PetscFunctionReturn(0);
1143a935fc98SLois Curfman McInnes }
114404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
114504d965bbSLois Curfman McInnes /*
11464b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
114704d965bbSLois Curfman McInnes 
114804d965bbSLois Curfman McInnes    Input Parameter:
114904d965bbSLois Curfman McInnes .  snes - the SNES context
115004d965bbSLois Curfman McInnes 
115104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
115204d965bbSLois Curfman McInnes */
11534a2ae208SSatish Balay #undef __FUNCT__
11544b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11556849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11565e42470aSBarry Smith {
11574b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1158e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1159dfbe8321SBarry Smith   PetscErrorCode ierr;
1160a7cc72afSBarry Smith   PetscInt       indx;
1161f1af5d2fSBarry Smith   PetscTruth     flg;
11625e42470aSBarry Smith 
11633a40ed3dSBarry Smith   PetscFunctionBegin;
1164b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11654b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11664b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11674b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1168186905e3SBarry Smith 
116917bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
117025cdf11fSBarry Smith     if (flg) {
117122e36657SBarry Smith       switch (indx) {
1172b49fd9e1SBarry Smith       case 0:
117317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1174b49fd9e1SBarry Smith         break;
1175b49fd9e1SBarry Smith       case 1:
117617bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1177b49fd9e1SBarry Smith         break;
1178b49fd9e1SBarry Smith       case 2:
117917bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1180b49fd9e1SBarry Smith         break;
1181b49fd9e1SBarry Smith       case 3:
118217bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1183b49fd9e1SBarry Smith         break;
11845e42470aSBarry Smith       }
11855e42470aSBarry Smith     }
1186b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11873a40ed3dSBarry Smith   PetscFunctionReturn(0);
11885e42470aSBarry Smith }
118904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1190ebe3b25bSBarry Smith /*MC
1191ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
119204d965bbSLois Curfman McInnes 
1193ebe3b25bSBarry Smith    Options Database:
1194ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1195ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1196ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1197ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1198ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1199ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
120004d965bbSLois Curfman McInnes 
1201ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
120204d965bbSLois Curfman McInnes 
1203ee3001cbSBarry Smith    Level: beginner
1204ee3001cbSBarry Smith 
120517bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
120617bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
120717bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1208ebe3b25bSBarry Smith 
1209ebe3b25bSBarry Smith M*/
1210fb2e594dSBarry Smith EXTERN_C_BEGIN
12114a2ae208SSatish Balay #undef __FUNCT__
12124b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
121363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
12145e42470aSBarry Smith {
1215dfbe8321SBarry Smith   PetscErrorCode ierr;
12164b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12175e42470aSBarry Smith 
12183a40ed3dSBarry Smith   PetscFunctionBegin;
12194b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
12204b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
12214b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
12224b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
12234b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
12244b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
12254b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
12265baf8537SBarry Smith   snes->nwork           = 0;
12275e42470aSBarry Smith 
12284b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
122952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
12305e42470aSBarry Smith   snes->data    	= (void*)neP;
12315e42470aSBarry Smith   neP->alpha		= 1.e-4;
12325e42470aSBarry Smith   neP->maxstep		= 1.e8;
12335e42470aSBarry Smith   neP->steptol		= 1.e-12;
123417bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1235c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12363c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12373c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12383c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12393c632250SBarry Smith   neP->precheck         = PETSC_NULL;
124082bf6240SBarry Smith 
124117bae607SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
12423c632250SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
12439c3ca13aSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
124482bf6240SBarry Smith 
12453a40ed3dSBarry Smith   PetscFunctionReturn(0);
12465e42470aSBarry Smith }
1247fb2e594dSBarry Smith EXTERN_C_END
124882bf6240SBarry Smith 
124982bf6240SBarry Smith 
125082bf6240SBarry Smith 
1251