xref: /petsc/src/snes/impls/ls/ls.c (revision 79f364605d8610b93e7f13639ad7dac12d01ffca)
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 
33abb0e124SMatthew 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;
5674637425SBarry Smith 
5774637425SBarry Smith   PetscFunctionBegin;
5874637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5974637425SBarry Smith   if (hastranspose) {
6074637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
61*79f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6274637425SBarry Smith 
6374637425SBarry Smith     /* Compute || J^T W|| */
6474637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6574637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6674637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
679e247f21SBarry Smith     if (a1) {
6863ba0a88SBarry Smith       ierr = PetscLogInfo((0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1));CHKERRQ(ierr);
6985471664SBarry Smith     }
7074637425SBarry Smith   }
7174637425SBarry Smith   PetscFunctionReturn(0);
7274637425SBarry Smith }
7374637425SBarry Smith 
7404d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7504d965bbSLois Curfman McInnes 
7604d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7794b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7804d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
7904d965bbSLois Curfman McInnes      respectively.
8004d965bbSLois Curfman McInnes 
81fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8204d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8304d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
84fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8504d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8604d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
874b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
884b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
8904d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9004d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9104d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9204d965bbSLois Curfman McInnes      for all nonlinear solvers.
9304d965bbSLois Curfman McInnes 
9404d965bbSLois Curfman McInnes      Another key routine is:
9504d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9604d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9704d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9804d965bbSLois Curfman McInnes      SNESSolve() if necessary.
9904d965bbSLois Curfman McInnes 
10004d965bbSLois Curfman McInnes      Additional basic routines are:
10104d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10204d965bbSLois Curfman McInnes                                       have actually been used.
10304d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
104186905e3SBarry Smith      SNESView().
10504d965bbSLois Curfman McInnes 
10604d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10704d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10804d965bbSLois Curfman McInnes      above description applies to these categories also.
10904d965bbSLois Curfman McInnes 
11004d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1115e42470aSBarry Smith /*
1124b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11304d965bbSLois Curfman McInnes    method with a line search.
1145e76c431SBarry Smith 
11504d965bbSLois Curfman McInnes    Input Parameters:
11604d965bbSLois Curfman McInnes .  snes - the SNES context
1175e76c431SBarry Smith 
11804d965bbSLois Curfman McInnes    Output Parameter:
119c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12004d965bbSLois Curfman McInnes 
12104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1225e76c431SBarry Smith 
1235e76c431SBarry Smith    Notes:
1245e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1255e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1265e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1275e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
128393d2d9aSLois Curfman McInnes    and Schnabel.
1295e76c431SBarry Smith */
1304a2ae208SSatish Balay #undef __FUNCT__
1314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
132dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1335e76c431SBarry Smith {
1344b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
1356849ba73SBarry Smith   PetscErrorCode ierr;
136a7cc72afSBarry Smith   PetscInt       maxits,i,lits;
137a7cc72afSBarry Smith   PetscTruth     lssucceed;
138112a2221SBarry Smith   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
139329f5518SBarry Smith   PetscReal      fnorm,gnorm,xnorm,ynorm;
1405e42470aSBarry Smith   Vec            Y,X,F,G,W,TMP;
141c293cc10SBarry Smith   KSP            ksp;
1425e76c431SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
14494b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
145184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
146184914b5SBarry Smith 
1475e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1485e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1505e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1515e42470aSBarry Smith   G		= snes->work[1];
1525e42470aSBarry Smith   W		= snes->work[2];
1535e76c431SBarry Smith 
1544c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1554c49b128SBarry Smith   snes->iter = 0;
1564c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15719717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
158cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
15943e71028SBarry Smith   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1600f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1615e42470aSBarry Smith   snes->norm = fnorm;
1620f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16384c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16494a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1655e76c431SBarry Smith 
16670441072SBarry Smith   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16794a9d846SBarry Smith 
168d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
169d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
170d034289bSBarry Smith 
1715e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1725e76c431SBarry Smith 
173329e7e40SMatthew Knepley     /* Call general purpose update function */
174abc0a331SBarry Smith     if (snes->update) {
175329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
176de8cb200SMatthew Knepley     }
177329e7e40SMatthew Knepley 
178ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1795334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18094b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18123ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
182c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
18374637425SBarry Smith 
1849c3ca13aSBarry Smith     if (neP->precheckstep) {
1859c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
1869c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
1879c3ca13aSBarry Smith     }
1889c3ca13aSBarry Smith 
189b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
19074637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
19185471664SBarry Smith     }
19274637425SBarry Smith 
19390cbd584SBarry Smith     /* should check what happened to the linear solve? */
1943505fcc1SBarry Smith     snes->linear_its += lits;
19563ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr);
196ea4d3ed3SLois Curfman McInnes 
197ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
198ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
199ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
200ea4d3ed3SLois Curfman McInnes     */
20181b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
202a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
20363ba0a88SBarry 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);
20443e71028SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
205a3b891d8SBarry Smith 
206a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
2073c632250SBarry Smith     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
208a3b891d8SBarry Smith     fnorm = gnorm;
209a3b891d8SBarry Smith 
2105ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2115ed2d596SBarry Smith     snes->iter = i+1;
2125ed2d596SBarry Smith     snes->norm = fnorm;
2135ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2145ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2155ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2165ed2d596SBarry Smith 
217a7cc72afSBarry Smith     if (!lssucceed) {
2188a5d9ee4SBarry Smith       PetscTruth ismin;
21950ffb88aSMatthew Knepley 
22050ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2213505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2228a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2238a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2243505fcc1SBarry Smith         break;
2253505fcc1SBarry Smith       }
22650ffb88aSMatthew Knepley     }
2275e76c431SBarry Smith 
2285e76c431SBarry Smith     /* Test for convergence */
22929e0b56fSBarry Smith     if (snes->converged) {
23029e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2313505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2323505fcc1SBarry Smith       if (snes->reason) {
23316c95cacSBarry Smith         break;
23416c95cacSBarry Smith       }
23516c95cacSBarry Smith     }
23629e0b56fSBarry Smith   }
23739e2f89bSBarry Smith   if (X != snes->vec_sol) {
238393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
239e15f7bb5SBarry Smith   }
240e15f7bb5SBarry Smith   if (F != snes->vec_func) {
241e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
242e15f7bb5SBarry Smith   }
24339e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
24439e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
24552392280SLois Curfman McInnes   if (i == maxits) {
24663ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr);
2473505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
24852392280SLois Curfman McInnes   }
2493a40ed3dSBarry Smith   PetscFunctionReturn(0);
2505e76c431SBarry Smith }
25104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25204d965bbSLois Curfman McInnes /*
2534b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2544b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25504d965bbSLois Curfman McInnes 
25604d965bbSLois Curfman McInnes    Input Parameter:
25704d965bbSLois Curfman McInnes .  snes - the SNES context
25804d965bbSLois Curfman McInnes .  x - the solution vector
25904d965bbSLois Curfman McInnes 
26004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26104d965bbSLois Curfman McInnes 
26204d965bbSLois Curfman McInnes    Notes:
2634b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26404d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26504d965bbSLois Curfman McInnes    the call to SNESSolve().
26604d965bbSLois Curfman McInnes  */
2674a2ae208SSatish Balay #undef __FUNCT__
2684b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
269dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2705e76c431SBarry Smith {
271dfbe8321SBarry Smith   PetscErrorCode ierr;
2723a40ed3dSBarry Smith 
2733a40ed3dSBarry Smith   PetscFunctionBegin;
274e74804ceSBarry Smith   if (!snes->work) {
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];
279e74804ceSBarry Smith   }
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
2815e76c431SBarry Smith }
28204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28304d965bbSLois Curfman McInnes /*
2844b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2854b27c08aSLois Curfman McInnes    with SNESCreate_LS().
28604d965bbSLois Curfman McInnes 
28704d965bbSLois Curfman McInnes    Input Parameter:
28804d965bbSLois Curfman McInnes .  snes - the SNES context
28904d965bbSLois Curfman McInnes 
290de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29104d965bbSLois Curfman McInnes  */
2924a2ae208SSatish Balay #undef __FUNCT__
2934b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
294dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2955e76c431SBarry Smith {
296dfbe8321SBarry Smith   PetscErrorCode ierr;
2973a40ed3dSBarry Smith 
2983a40ed3dSBarry Smith   PetscFunctionBegin;
2995baf8537SBarry Smith   if (snes->nwork) {
3004b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
30121c89e3eSBarry Smith   }
3025c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3033a40ed3dSBarry Smith   PetscFunctionReturn(0);
3045e76c431SBarry Smith }
30504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3064a2ae208SSatish Balay #undef __FUNCT__
30717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
30882bf6240SBarry Smith 
3094b828684SBarry Smith /*@C
31017bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3115e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3125e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3135e76c431SBarry Smith 
314c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
315c7afd0dbSLois Curfman McInnes 
3165e76c431SBarry Smith    Input Parameters:
317c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
318af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3195e76c431SBarry Smith .  x - current iterate
3205e76c431SBarry Smith .  f - residual evaluated at x
3213c632250SBarry Smith .  y - search direction
3225e76c431SBarry Smith .  w - work vector
323c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3245e76c431SBarry Smith 
325c4a48953SLois Curfman McInnes    Output Parameters:
326c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3273c632250SBarry Smith .  w - new iterate
3285e76c431SBarry Smith .  gnorm - 2-norm of g
3295e76c431SBarry Smith .  ynorm - 2-norm of search length
330a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
331fee21e36SBarry Smith 
332c4a48953SLois Curfman McInnes    Options Database Key:
33317bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
334c4a48953SLois Curfman McInnes 
33536851e7fSLois Curfman McInnes    Level: advanced
33636851e7fSLois Curfman McInnes 
33728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33828ae5a21SLois Curfman McInnes 
33917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
34017bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3415e76c431SBarry Smith @*/
34217bae607SBarry 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)
3435e76c431SBarry Smith {
344dfbe8321SBarry Smith   PetscErrorCode ierr;
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 || */
352*79f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);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) {
357*79f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);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__
37217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
37382bf6240SBarry Smith 
37429e0b56fSBarry Smith /*@C
37517bae607SBarry Smith    SNESLineSearchNoNorms - 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:
40117bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
40229e0b56fSBarry Smith 
4038cbba510SBarry Smith    Notes:
40417bae607SBarry Smith    SNESLineSearchNoNorms() 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 
42317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
42417bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
42529e0b56fSBarry Smith @*/
42617bae607SBarry 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)
42729e0b56fSBarry Smith {
428dfbe8321SBarry Smith   PetscErrorCode ierr;
4294b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4303c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
43129e0b56fSBarry Smith 
4323a40ed3dSBarry Smith   PetscFunctionBegin;
433a7cc72afSBarry Smith   *flag = PETSC_TRUE;
434d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
435*79f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4363c632250SBarry Smith   if (neP->postcheckstep) {
4373c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4383c632250SBarry Smith   }
4393c632250SBarry Smith   if (changed_y) {
440*79f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4418f99978dSLois Curfman McInnes   }
442329f5518SBarry Smith 
443329f5518SBarry Smith   /* don't evaluate function the last time through */
444329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4453c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
446d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
44719717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
44819717074SBarry Smith     }
449d5e45103SBarry Smith     CHKERRQ(ierr);
450329f5518SBarry Smith   }
451d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4523a40ed3dSBarry Smith   PetscFunctionReturn(0);
45329e0b56fSBarry Smith }
45404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4554a2ae208SSatish Balay #undef __FUNCT__
45617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4574b828684SBarry Smith /*@C
45817bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4595e76c431SBarry Smith 
460c7afd0dbSLois Curfman McInnes    Collective on SNES
461c7afd0dbSLois Curfman McInnes 
4625e76c431SBarry Smith    Input Parameters:
463c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
464af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4655e76c431SBarry Smith .  x - current iterate
4665e76c431SBarry Smith .  f - residual evaluated at x
4673c632250SBarry Smith .  y - search direction
4685e76c431SBarry Smith .  w - work vector
469c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4705e76c431SBarry Smith 
471393d2d9aSLois Curfman McInnes    Output Parameters:
472c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4733c632250SBarry Smith .  w - new iterate
4745e76c431SBarry Smith .  gnorm - 2-norm of g
4755e76c431SBarry Smith .  ynorm - 2-norm of search length
476a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
477fee21e36SBarry Smith 
478c4a48953SLois Curfman McInnes    Options Database Key:
47917bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
480c4a48953SLois Curfman McInnes 
4815e76c431SBarry Smith    Notes:
4825e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4835e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4845e76c431SBarry Smith 
48536851e7fSLois Curfman McInnes    Level: advanced
48636851e7fSLois Curfman McInnes 
48728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
48828ae5a21SLois Curfman McInnes 
48917bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
4905e76c431SBarry Smith @*/
49117bae607SBarry 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)
4925e76c431SBarry Smith {
493406556e6SLois Curfman McInnes   /*
494406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
495406556e6SLois Curfman McInnes      minimization problem:
496406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
497406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
498406556e6SLois Curfman McInnes    */
499406556e6SLois Curfman McInnes 
5005ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
501329f5518SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
502aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
503ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5046b5873e3SBarry Smith #endif
505dfbe8321SBarry Smith   PetscErrorCode ierr;
506a7cc72afSBarry Smith   PetscInt       count;
5074b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
508*79f36460SBarry Smith   PetscScalar    scale;
5093c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5105e76c431SBarry Smith 
5113a40ed3dSBarry Smith   PetscFunctionBegin;
512d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
513a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5145e76c431SBarry Smith   alpha   = neP->alpha;
5155e76c431SBarry Smith   maxstep = neP->maxstep;
5165e76c431SBarry Smith   steptol = neP->steptol;
5175e76c431SBarry Smith 
518cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5199e247f21SBarry Smith   if (!*ynorm) {
52017bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Search direction and size is 0\n"));CHKERRQ(ierr);
521a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5223c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
523a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
524a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
525ad922ac9SBarry Smith     goto theend1;
52694a9d846SBarry Smith   }
5275e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5285e42470aSBarry Smith     scale = maxstep/(*ynorm);
529aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
53017bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr);
5316b5873e3SBarry Smith #else
53217bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr);
5336b5873e3SBarry Smith #endif
5342dcb1b2aSMatthew Knepley     ierr = VecScale(y,scale);CHKERRQ(ierr);
5355e76c431SBarry Smith     *ynorm = maxstep;
5365e76c431SBarry Smith   }
537ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5385ca10a37SBarry Smith   minlambda = steptol/rellength;
539a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
540aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
541a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
542329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5435e42470aSBarry Smith #else
544a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5455e42470aSBarry Smith #endif
5465e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5475e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5485e76c431SBarry Smith 
549393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
550*79f36460SBarry Smith   ierr = VecAYPX(w,-1.0,x);CHKERRQ(ierr);
55143e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
55217bae607SBarry Smith     ierr  = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
55343e71028SBarry Smith     *flag = PETSC_FALSE;
55443e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
55543e71028SBarry Smith     goto theend1;
55643e71028SBarry Smith   }
55719717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
558d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
55919717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
56019717074SBarry Smith   }
561d5e45103SBarry Smith   CHKERRQ(ierr);
562313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56343e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5645d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
56517bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Using full step\n"));CHKERRQ(ierr);
566ad922ac9SBarry Smith     goto theend1;
5675e76c431SBarry Smith   }
5685e76c431SBarry Smith 
5695e76c431SBarry Smith   /* Fit points with quadratic */
570313b5538SBarry Smith   lambda     = 1.0;
571a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5725e76c431SBarry Smith   lambdaprev = lambda;
5735e76c431SBarry Smith   gnormprev  = *gnorm;
574329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
575ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
576ddd12767SBarry Smith   else                         lambda = lambdatemp;
577393d2d9aSLois Curfman McInnes   ierr      = VecCopy(x,w);CHKERRQ(ierr);
578ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
579aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
5802dcb1b2aSMatthew Knepley   clambda   = lambdaneg; ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr);
5815e42470aSBarry Smith #else
5822dcb1b2aSMatthew Knepley   ierr      = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr);
5835e42470aSBarry Smith #endif
58443e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
58517bae607SBarry Smith     ierr  = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr);
58643e71028SBarry Smith     *flag = PETSC_FALSE;
58743e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
58843e71028SBarry Smith     goto theend1;
58943e71028SBarry Smith   }
59019717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
591d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59219717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59319717074SBarry Smith   }
594d5e45103SBarry Smith   CHKERRQ(ierr);
595cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
59643e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5975ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
59817bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
599ad922ac9SBarry Smith     goto theend1;
6005e76c431SBarry Smith   }
6015e76c431SBarry Smith 
6025e76c431SBarry Smith   /* Fit points with cubic */
6035e76c431SBarry Smith   count = 1;
6048229bfc2SMatthew Knepley   while (count < 10000) {
6055e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
60617bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
60717bae607SBarry 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);
60843e71028SBarry Smith       *flag = PETSC_FALSE;
60943e71028SBarry Smith       break;
6105e76c431SBarry Smith     }
6115d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6125d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
613ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6142b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6155e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6165e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6175e76c431SBarry Smith     if (a == 0.0) {
6185e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6195e76c431SBarry Smith     } else {
6205e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6215e76c431SBarry Smith     }
6225e76c431SBarry Smith     lambdaprev = lambda;
6235e76c431SBarry Smith     gnormprev  = *gnorm;
624329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
625329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6265e76c431SBarry Smith     else                         lambda     = lambdatemp;
627393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
628ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
629aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
630ea4d3ed3SLois Curfman McInnes     clambda   = lambdaneg;
6312dcb1b2aSMatthew Knepley     ierr      = VecAXPY(w,clambda,y);CHKERRQ(ierr);
6325e42470aSBarry Smith #else
6332dcb1b2aSMatthew Knepley     ierr      = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr);
6345e42470aSBarry Smith #endif
63543e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
63617bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
63717bae607SBarry 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);
638ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
63943e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
640ed98166eSMatthew Knepley       break;
641ed98166eSMatthew Knepley     }
64219717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
643d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64419717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64519717074SBarry Smith     }
646d5e45103SBarry Smith     CHKERRQ(ierr);
647cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
64843e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6495ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
65017bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
65143e71028SBarry Smith       break;
6522b022350SLois Curfman McInnes     } else {
65317bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchCubic: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda));CHKERRQ(ierr);
6545e76c431SBarry Smith     }
6555e76c431SBarry Smith     count++;
6565e76c431SBarry Smith   }
6578229bfc2SMatthew Knepley   if (count >= 10000) {
658cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6598229bfc2SMatthew Knepley   }
660ad922ac9SBarry Smith   theend1:
6618f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6623c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6633c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6643c632250SBarry Smith     if (changed_y) {
665*79f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6663c632250SBarry Smith     }
6673c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6683c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
669d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
67019717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
67119717074SBarry Smith       }
672d5e45103SBarry Smith       CHKERRQ(ierr);
6738f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67443e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6753c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6768f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6773c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6788f99978dSLois Curfman McInnes     }
6798f99978dSLois Curfman McInnes   }
680d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6813a40ed3dSBarry Smith   PetscFunctionReturn(0);
6825e76c431SBarry Smith }
68304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6844a2ae208SSatish Balay #undef __FUNCT__
68517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6864b828684SBarry Smith /*@C
68717bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6885e76c431SBarry Smith 
689c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
690c7afd0dbSLois Curfman McInnes 
6915e42470aSBarry Smith    Input Parameters:
692c7afd0dbSLois Curfman McInnes +  snes - the SNES context
693af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6945e42470aSBarry Smith .  x - current iterate
6955e42470aSBarry Smith .  f - residual evaluated at x
6963c632250SBarry Smith .  y - search direction
6975e42470aSBarry Smith .  w - work vector
698c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6995e42470aSBarry Smith 
700c4a48953SLois Curfman McInnes    Output Parameters:
701c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
7023c632250SBarry Smith .  w - new iterate
7035e42470aSBarry Smith .  gnorm - 2-norm of g
7045e42470aSBarry Smith .  ynorm - 2-norm of search length
705a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
706fee21e36SBarry Smith 
707c4a48953SLois Curfman McInnes    Options Database Key:
70817bae607SBarry Smith .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
709c4a48953SLois Curfman McInnes 
7105e42470aSBarry Smith    Notes:
71117bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7125e42470aSBarry Smith 
71336851e7fSLois Curfman McInnes    Level: advanced
71436851e7fSLois Curfman McInnes 
715f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
716f59ffdedSLois Curfman McInnes 
71717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7185e42470aSBarry Smith @*/
71917bae607SBarry 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)
7205e76c431SBarry Smith {
721406556e6SLois Curfman McInnes   /*
722406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
723406556e6SLois Curfman McInnes      minimization problem:
724406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
725406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
726406556e6SLois Curfman McInnes    */
7275ca10a37SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
728aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
729ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7306b5873e3SBarry Smith #endif
731dfbe8321SBarry Smith   PetscErrorCode ierr;
732a7cc72afSBarry Smith   PetscInt       count;
7334b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
734*79f36460SBarry Smith   PetscScalar    scale;
7353c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7365e76c431SBarry Smith 
7373a40ed3dSBarry Smith   PetscFunctionBegin;
738d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
739a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7405e42470aSBarry Smith   alpha   = neP->alpha;
7415e42470aSBarry Smith   maxstep = neP->maxstep;
7425e42470aSBarry Smith   steptol = neP->steptol;
7435e76c431SBarry Smith 
7443505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74563b9fa5eSBarry Smith   if (*ynorm == 0.0) {
74617bae607SBarry Smith     ierr   = PetscLogInfo((snes,"SNESLineSearchQuadratic: Search direction and size is 0\n"));CHKERRQ(ierr);
747b37302c6SLois Curfman McInnes     *gnorm = fnorm;
748b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
749b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
750a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
751ad922ac9SBarry Smith     goto theend2;
75294a9d846SBarry Smith   }
7535e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7545e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7552dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7565e42470aSBarry Smith     *ynorm = maxstep;
7575e76c431SBarry Smith   }
758ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7595ca10a37SBarry Smith   minlambda = steptol/rellength;
760a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
761aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
762a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
763329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7645e42470aSBarry Smith #else
765a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7665e42470aSBarry Smith #endif
7675e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7685e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7695e42470aSBarry Smith 
770393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
771*79f36460SBarry Smith   ierr = VecAYPX(w,-1.0,x);CHKERRQ(ierr);
77243e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
77317bae607SBarry Smith     ierr  = PetscLogInfo((snes,"SNESLineSearchQuadratic:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
77443e71028SBarry Smith     ierr  = VecCopy(w,y);CHKERRQ(ierr);
77543e71028SBarry Smith     *flag = PETSC_FALSE;
77643e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77743e71028SBarry Smith     goto theend2;
77843e71028SBarry Smith   }
77919717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
780d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
78119717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
78219717074SBarry Smith   }
783d5e45103SBarry Smith   CHKERRQ(ierr);
784cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
78543e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7865d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
787393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
78817bae607SBarry Smith     ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic: Using full step\n"));CHKERRQ(ierr);
789ad922ac9SBarry Smith     goto theend2;
7905e42470aSBarry Smith   }
7915e42470aSBarry Smith 
7925e42470aSBarry Smith   /* Fit points with quadratic */
793313b5538SBarry Smith   lambda = 1.0;
7945e42470aSBarry Smith   count = 1;
7955ca10a37SBarry Smith   while (PETSC_TRUE) {
7965e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
79717bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
79817bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
799f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
80043e71028SBarry Smith       *flag = PETSC_FALSE;
80143e71028SBarry Smith       break;
8025e42470aSBarry Smith     }
803a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
804329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
805329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
806329f5518SBarry Smith     else                         lambda     = lambdatemp;
807393d2d9aSLois Curfman McInnes     ierr      = VecCopy(x,w);CHKERRQ(ierr);
8083505fcc1SBarry Smith     lambdaneg = -lambda;
809aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
8102dcb1b2aSMatthew Knepley     clambda   = lambdaneg; ierr = VecAXPY(w,clambda,y);CHKERRQ(ierr);
8115e42470aSBarry Smith #else
8122dcb1b2aSMatthew Knepley     ierr      = VecAXPY(w,lambdaneg,y);CHKERRQ(ierr);
8135e42470aSBarry Smith #endif
81443e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
81517bae607SBarry Smith       ierr  = PetscLogInfo((snes,"SNESLineSearchQuadratic:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
81617bae607SBarry 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);
81743e71028SBarry Smith       ierr  = VecCopy(w,y);CHKERRQ(ierr);
818ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
81943e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
820ed98166eSMatthew Knepley       break;
821ed98166eSMatthew Knepley     }
82219717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
823d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
82419717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
82519717074SBarry Smith     }
826d5e45103SBarry Smith     CHKERRQ(ierr);
827cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82843e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8295ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
830393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
83117bae607SBarry Smith       ierr = PetscLogInfo((snes,"SNESLineSearchQuadratic:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr);
83206259719SBarry Smith       break;
8335e42470aSBarry Smith     }
8345e42470aSBarry Smith     count++;
8355e42470aSBarry Smith   }
836ad922ac9SBarry Smith   theend2:
8378f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8383c632250SBarry Smith   if (neP->postcheckstep) {
8393c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8403c632250SBarry Smith     if (changed_y) {
841*79f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8423c632250SBarry Smith     }
8433c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8443c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
845d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
84619717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
84719717074SBarry Smith       }
848d5e45103SBarry Smith       CHKERRQ(ierr);
8498f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8503c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8518f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8523c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8533c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8548f99978dSLois Curfman McInnes     }
8558f99978dSLois Curfman McInnes   }
856d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8573a40ed3dSBarry Smith   PetscFunctionReturn(0);
8585e76c431SBarry Smith }
8592343ba6eSBarry Smith 
86004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8614a2ae208SSatish Balay #undef __FUNCT__
86217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
863c9e6a524SLois Curfman McInnes /*@C
86417bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8654b27c08aSLois Curfman McInnes    by the method SNESLS.
8665e76c431SBarry Smith 
8675e76c431SBarry Smith    Input Parameters:
868c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
869af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
870c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8715e76c431SBarry Smith 
872fee21e36SBarry Smith    Collective on SNES
873fee21e36SBarry Smith 
874c4a48953SLois Curfman McInnes    Available Routines:
87517bae607SBarry Smith +  SNESLineSearchCubic() - default line search
87617bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
87717bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
87817bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8795e76c431SBarry Smith 
880c4a48953SLois Curfman McInnes     Options Database Keys:
8814b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8824b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8834b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8844b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8853304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8863304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
887c4a48953SLois Curfman McInnes 
8885e76c431SBarry Smith    Calling sequence of func:
889bd208895SLois Curfman McInnes .vb
890af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
891a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
892bd208895SLois Curfman McInnes .ve
8935e76c431SBarry Smith 
8945e76c431SBarry Smith     Input parameters for func:
895c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
896af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8975e76c431SBarry Smith .   x - current iterate
8985e76c431SBarry Smith .   f - residual evaluated at x
8993c632250SBarry Smith .   y - search direction
9005e76c431SBarry Smith .   w - work vector
901c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9025e76c431SBarry Smith 
9035e76c431SBarry Smith     Output parameters for func:
904c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9053c632250SBarry Smith .   w - new iterate
9065e76c431SBarry Smith .   gnorm - 2-norm of g
9075e76c431SBarry Smith .   ynorm - 2-norm of search length
908a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
909f59ffdedSLois Curfman McInnes 
91036851e7fSLois Curfman McInnes     Level: advanced
91136851e7fSLois Curfman McInnes 
912f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
913f59ffdedSLois Curfman McInnes 
91417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
91517bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9165e76c431SBarry Smith @*/
91717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9185e76c431SBarry Smith {
919a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
92082bf6240SBarry Smith 
9213a40ed3dSBarry Smith   PetscFunctionBegin;
92217bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
92382bf6240SBarry Smith   if (f) {
924af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92582bf6240SBarry Smith   }
9263a40ed3dSBarry Smith   PetscFunctionReturn(0);
9275e76c431SBarry Smith }
9288e019c35SBarry Smith 
929a7cc72afSBarry 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*/
93004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
931fb2e594dSBarry Smith EXTERN_C_BEGIN
9324a2ae208SSatish Balay #undef __FUNCT__
93317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
93417bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
93582bf6240SBarry Smith {
93682bf6240SBarry Smith   PetscFunctionBegin;
9374b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9384b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
93982bf6240SBarry Smith   PetscFunctionReturn(0);
94082bf6240SBarry Smith }
941fb2e594dSBarry Smith EXTERN_C_END
94204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9434a2ae208SSatish Balay #undef __FUNCT__
9443c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
945c8dd0c0dSLois Curfman McInnes /*@C
9463c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9474b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
948c8dd0c0dSLois Curfman McInnes 
949c8dd0c0dSLois Curfman McInnes    Input Parameters:
950c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9513c632250SBarry Smith .  func - pointer to function
952c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
953c8dd0c0dSLois Curfman McInnes 
954c8dd0c0dSLois Curfman McInnes    Collective on SNES
955c8dd0c0dSLois Curfman McInnes 
956c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
957c8dd0c0dSLois Curfman McInnes .vb
9583c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
959c8dd0c0dSLois Curfman McInnes .ve
960b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
961b1ae27eaSLois Curfman McInnes    on failure.
962c8dd0c0dSLois Curfman McInnes 
963c8dd0c0dSLois Curfman McInnes    Input parameters for func:
964c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
965c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9663c632250SBarry Smith .  x - previous iterate
9673c632250SBarry Smith .  y - new search direction and length
9683c632250SBarry Smith -  w - current candidate iterate
969c8dd0c0dSLois Curfman McInnes 
970c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9713c632250SBarry Smith +  y - search direction (possibly changed)
9723c632250SBarry Smith .  w - current iterate (possibly modified)
9733c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9743c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
975c8dd0c0dSLois Curfman McInnes 
976c8dd0c0dSLois Curfman McInnes    Level: advanced
977c8dd0c0dSLois Curfman McInnes 
9789e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9799e247f21SBarry Smith 
9803c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9813c632250SBarry Smith 
9823c632250SBarry Smith    On input w = x + y
9833c632250SBarry Smith 
98417bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
985b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
986ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9878f99978dSLois Curfman McInnes 
98817bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
989ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
990ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
991ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9929e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
993b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9948f99978dSLois Curfman McInnes 
995c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
996c8dd0c0dSLois Curfman McInnes 
99717bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
998c8dd0c0dSLois Curfman McInnes @*/
9993c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
1000c8dd0c0dSLois Curfman McInnes {
10013c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
1002c8dd0c0dSLois Curfman McInnes 
1003c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10043c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1005c8dd0c0dSLois Curfman McInnes   if (f) {
1006c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1007c8dd0c0dSLois Curfman McInnes   }
1008c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1009c8dd0c0dSLois Curfman McInnes }
10109c3ca13aSBarry Smith #undef __FUNCT__
10119c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10129c3ca13aSBarry Smith /*@C
10139c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10149c3ca13aSBarry Smith 
10159c3ca13aSBarry Smith    Input Parameters:
10169c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10179c3ca13aSBarry Smith .  func - pointer to function
10189c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10199c3ca13aSBarry Smith 
10209c3ca13aSBarry Smith    Collective on SNES
10219c3ca13aSBarry Smith 
10229c3ca13aSBarry Smith    Calling sequence of func:
10239c3ca13aSBarry Smith .vb
10249c3ca13aSBarry Smith    int func (SNES snes, Vec x,Vec y,,void *checkctx, PetscTruth *changed_y)
10259c3ca13aSBarry Smith .ve
10269c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10279c3ca13aSBarry Smith    on failure.
10289c3ca13aSBarry Smith 
10299c3ca13aSBarry Smith    Input parameters for func:
10309c3ca13aSBarry Smith +  snes - nonlinear context
10319c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10329c3ca13aSBarry Smith .  x - previous iterate
10339c3ca13aSBarry Smith -  y - new search direction and length
10349c3ca13aSBarry Smith 
10359c3ca13aSBarry Smith    Output parameters for func:
10369c3ca13aSBarry Smith +  y - search direction (possibly changed)
10379c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10389c3ca13aSBarry Smith 
10399c3ca13aSBarry Smith    Level: advanced
10409c3ca13aSBarry Smith 
10419c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10429c3ca13aSBarry Smith 
10439c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10449c3ca13aSBarry Smith 
104517bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck()
10469c3ca13aSBarry Smith @*/
10479c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10489c3ca13aSBarry Smith {
10499c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10509c3ca13aSBarry Smith 
10519c3ca13aSBarry Smith   PetscFunctionBegin;
10529c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10539c3ca13aSBarry Smith   if (f) {
10549c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10559c3ca13aSBarry Smith   }
10569c3ca13aSBarry Smith   PetscFunctionReturn(0);
10579c3ca13aSBarry Smith }
10589c3ca13aSBarry Smith 
1059c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10609c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10619c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1062c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10634a2ae208SSatish Balay #undef __FUNCT__
10643c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10659c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1066c8dd0c0dSLois Curfman McInnes {
1067c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10683c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10693c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1070c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1071c8dd0c0dSLois Curfman McInnes }
1072c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10739c3ca13aSBarry Smith 
10749c3ca13aSBarry Smith EXTERN_C_BEGIN
10759c3ca13aSBarry Smith #undef __FUNCT__
10769c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10779c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10789c3ca13aSBarry Smith {
10799c3ca13aSBarry Smith   PetscFunctionBegin;
10809c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10819c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10829c3ca13aSBarry Smith   PetscFunctionReturn(0);
10839c3ca13aSBarry Smith }
10849c3ca13aSBarry Smith EXTERN_C_END
1085c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
108604d965bbSLois Curfman McInnes /*
10874b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
1088329e7e40SMatthew Knepley 
1089329e7e40SMatthew Knepley    Input Parameter:
1090329e7e40SMatthew Knepley .  snes - the SNES context
1091329e7e40SMatthew Knepley 
1092329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
1093329e7e40SMatthew Knepley */
1094329e7e40SMatthew Knepley #undef __FUNCT__
10954b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
10966849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
1097329e7e40SMatthew Knepley {
10984b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
1099329e7e40SMatthew Knepley 
1100329e7e40SMatthew Knepley   PetscFunctionBegin;
11014b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
11024b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
11034b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
11044b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
11054b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
1106329e7e40SMatthew Knepley   PetscFunctionReturn(0);
1107329e7e40SMatthew Knepley }
1108329e7e40SMatthew Knepley 
1109329e7e40SMatthew Knepley /*
11104b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
111104d965bbSLois Curfman McInnes 
111204d965bbSLois Curfman McInnes    Input Parameters:
111304d965bbSLois Curfman McInnes .  SNES - the SNES context
111404d965bbSLois Curfman McInnes .  viewer - visualization context
111504d965bbSLois Curfman McInnes 
111604d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
111704d965bbSLois Curfman McInnes */
11184a2ae208SSatish Balay #undef __FUNCT__
11194b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11206849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1121a935fc98SLois Curfman McInnes {
11224b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11232fc52814SBarry Smith   const char     *cstr;
1124dfbe8321SBarry Smith   PetscErrorCode ierr;
112532077d6dSBarry Smith   PetscTruth     iascii;
1126a935fc98SLois Curfman McInnes 
11273a40ed3dSBarry Smith   PetscFunctionBegin;
112832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
112932077d6dSBarry Smith   if (iascii) {
113017bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
113117bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
113217bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
113319bcc07fSBarry Smith     else                                                cstr = "unknown";
1134b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1135b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11365cd90555SBarry Smith   } else {
113779a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
113819bcc07fSBarry Smith   }
11393a40ed3dSBarry Smith   PetscFunctionReturn(0);
1140a935fc98SLois Curfman McInnes }
114104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
114204d965bbSLois Curfman McInnes /*
11434b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
114404d965bbSLois Curfman McInnes 
114504d965bbSLois Curfman McInnes    Input Parameter:
114604d965bbSLois Curfman McInnes .  snes - the SNES context
114704d965bbSLois Curfman McInnes 
114804d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
114904d965bbSLois Curfman McInnes */
11504a2ae208SSatish Balay #undef __FUNCT__
11514b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11526849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11535e42470aSBarry Smith {
11544b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1155e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1156dfbe8321SBarry Smith   PetscErrorCode ierr;
1157a7cc72afSBarry Smith   PetscInt       indx;
1158f1af5d2fSBarry Smith   PetscTruth     flg;
11595e42470aSBarry Smith 
11603a40ed3dSBarry Smith   PetscFunctionBegin;
1161b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11624b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11634b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11644b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1165186905e3SBarry Smith 
116617bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
116725cdf11fSBarry Smith     if (flg) {
116822e36657SBarry Smith       switch (indx) {
1169b49fd9e1SBarry Smith       case 0:
117017bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1171b49fd9e1SBarry Smith         break;
1172b49fd9e1SBarry Smith       case 1:
117317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1174b49fd9e1SBarry Smith         break;
1175b49fd9e1SBarry Smith       case 2:
117617bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1177b49fd9e1SBarry Smith         break;
1178b49fd9e1SBarry Smith       case 3:
117917bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1180b49fd9e1SBarry Smith         break;
11815e42470aSBarry Smith       }
11825e42470aSBarry Smith     }
1183b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11843a40ed3dSBarry Smith   PetscFunctionReturn(0);
11855e42470aSBarry Smith }
118604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1187ebe3b25bSBarry Smith /*MC
1188ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
118904d965bbSLois Curfman McInnes 
1190ebe3b25bSBarry Smith    Options Database:
1191ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1192ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1193ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1194ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1195ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1196ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
119704d965bbSLois Curfman McInnes 
1198ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
119904d965bbSLois Curfman McInnes 
1200ee3001cbSBarry Smith    Level: beginner
1201ee3001cbSBarry Smith 
120217bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
120317bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
120417bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1205ebe3b25bSBarry Smith 
1206ebe3b25bSBarry Smith M*/
1207fb2e594dSBarry Smith EXTERN_C_BEGIN
12084a2ae208SSatish Balay #undef __FUNCT__
12094b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
121063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
12115e42470aSBarry Smith {
1212dfbe8321SBarry Smith   PetscErrorCode ierr;
12134b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12145e42470aSBarry Smith 
12153a40ed3dSBarry Smith   PetscFunctionBegin;
12164b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
12174b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
12184b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
12194b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
12204b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
12214b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
12224b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
12235baf8537SBarry Smith   snes->nwork           = 0;
12245e42470aSBarry Smith 
12254b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
122652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
12275e42470aSBarry Smith   snes->data    	= (void*)neP;
12285e42470aSBarry Smith   neP->alpha		= 1.e-4;
12295e42470aSBarry Smith   neP->maxstep		= 1.e8;
12305e42470aSBarry Smith   neP->steptol		= 1.e-12;
123117bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1232c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12333c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12343c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12353c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12363c632250SBarry Smith   neP->precheck         = PETSC_NULL;
123782bf6240SBarry Smith 
123817bae607SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
12393c632250SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
12409c3ca13aSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
124182bf6240SBarry Smith 
12423a40ed3dSBarry Smith   PetscFunctionReturn(0);
12435e42470aSBarry Smith }
1244fb2e594dSBarry Smith EXTERN_C_END
124582bf6240SBarry Smith 
124682bf6240SBarry Smith 
124782bf6240SBarry Smith 
1248