xref: /petsc/src/snes/impls/ls/ls.c (revision a6570f20d5bf426f46356aeda2a42fd24720f376)
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);
26ae15b995SBarry Smith     ierr = PetscInfo1(0,"|| 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);
40ae15b995SBarry Smith     ierr = PetscInfo1(0,"(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);
6179f36460SBarry 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) {
68ae15b995SBarry Smith       ierr = PetscInfo1(0,"||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;
1413d4c4710SBarry Smith   KSPConvergedReason kspreason;
1425e76c431SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
1443d4c4710SBarry Smith   snes->numFailures            = 0;
1453d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
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);
166e7788613SBarry Smith   ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
16706ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
168d034289bSBarry Smith 
1695e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1705e76c431SBarry Smith 
171329e7e40SMatthew Knepley     /* Call general purpose update function */
172e7788613SBarry Smith     if (snes->ops->update) {
173e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
174de8cb200SMatthew Knepley     }
175329e7e40SMatthew Knepley 
176ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1775334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
17894b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
17923ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
1803d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1813d4c4710SBarry Smith     if (kspreason < 0) {
1823d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
1833d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
1843d4c4710SBarry Smith         PetscFunctionReturn(0);
1853d4c4710SBarry Smith       }
1863d4c4710SBarry Smith     }
1873d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
18874637425SBarry Smith 
1899c3ca13aSBarry Smith     if (neP->precheckstep) {
1909c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
1919c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
1929c3ca13aSBarry Smith     }
1939c3ca13aSBarry Smith 
194b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
19574637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
19685471664SBarry Smith     }
19774637425SBarry Smith 
19890cbd584SBarry Smith     /* should check what happened to the linear solve? */
1993505fcc1SBarry Smith     snes->linear_its += lits;
200ae15b995SBarry Smith     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
201ea4d3ed3SLois Curfman McInnes 
202ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
203ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
204e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
205ea4d3ed3SLois Curfman McInnes     */
20681b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
207a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
208ae15b995SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
209a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
2103c632250SBarry Smith     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
211a3b891d8SBarry Smith     fnorm = gnorm;
2124a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
213a3b891d8SBarry Smith 
2145ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2155ed2d596SBarry Smith     snes->iter = i+1;
2165ed2d596SBarry Smith     snes->norm = fnorm;
2175ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2185ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2195ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2205ed2d596SBarry Smith 
221a7cc72afSBarry Smith     if (!lssucceed) {
2228a5d9ee4SBarry Smith       PetscTruth ismin;
22350ffb88aSMatthew Knepley 
22450ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2253505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2268a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2278a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2283505fcc1SBarry Smith         break;
2293505fcc1SBarry Smith       }
23050ffb88aSMatthew Knepley     }
2315e76c431SBarry Smith 
2325e76c431SBarry Smith     /* Test for convergence */
233e7788613SBarry Smith     if (snes->ops->converged) {
23429e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
235e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2363505fcc1SBarry Smith       if (snes->reason) {
23716c95cacSBarry Smith         break;
23816c95cacSBarry Smith       }
23916c95cacSBarry Smith     }
24029e0b56fSBarry Smith   }
24139e2f89bSBarry Smith   if (X != snes->vec_sol) {
242393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
243e15f7bb5SBarry Smith   }
244e15f7bb5SBarry Smith   if (F != snes->vec_func) {
245e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
246e15f7bb5SBarry Smith   }
24739e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
24839e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
24952392280SLois Curfman McInnes   if (i == maxits) {
250ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
2513505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
25252392280SLois Curfman McInnes   }
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
2545e76c431SBarry Smith }
25504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25604d965bbSLois Curfman McInnes /*
2574b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2584b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25904d965bbSLois Curfman McInnes 
26004d965bbSLois Curfman McInnes    Input Parameter:
26104d965bbSLois Curfman McInnes .  snes - the SNES context
26204d965bbSLois Curfman McInnes .  x - the solution vector
26304d965bbSLois Curfman McInnes 
26404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26504d965bbSLois Curfman McInnes 
26604d965bbSLois Curfman McInnes    Notes:
2674b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26804d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26904d965bbSLois Curfman McInnes    the call to SNESSolve().
27004d965bbSLois Curfman McInnes  */
2714a2ae208SSatish Balay #undef __FUNCT__
2724b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
273dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2745e76c431SBarry Smith {
275dfbe8321SBarry Smith   PetscErrorCode ierr;
2763a40ed3dSBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
278e74804ceSBarry Smith   if (!snes->work) {
27981b6cf68SLois Curfman McInnes     snes->nwork = 4;
280d7e8b826SBarry Smith     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
281efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
28281b6cf68SLois Curfman McInnes     snes->vec_sol_update_always = snes->work[3];
283e74804ceSBarry Smith   }
2843a40ed3dSBarry Smith   PetscFunctionReturn(0);
2855e76c431SBarry Smith }
28604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28704d965bbSLois Curfman McInnes /*
2884b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2894b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29004d965bbSLois Curfman McInnes 
29104d965bbSLois Curfman McInnes    Input Parameter:
29204d965bbSLois Curfman McInnes .  snes - the SNES context
29304d965bbSLois Curfman McInnes 
294de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29504d965bbSLois Curfman McInnes  */
2964a2ae208SSatish Balay #undef __FUNCT__
2974b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
298dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2995e76c431SBarry Smith {
300dfbe8321SBarry Smith   PetscErrorCode ierr;
3013a40ed3dSBarry Smith 
3023a40ed3dSBarry Smith   PetscFunctionBegin;
3035baf8537SBarry Smith   if (snes->nwork) {
3044b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
30521c89e3eSBarry Smith   }
3065c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3073a40ed3dSBarry Smith   PetscFunctionReturn(0);
3085e76c431SBarry Smith }
30904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3104a2ae208SSatish Balay #undef __FUNCT__
31117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
31282bf6240SBarry Smith 
3134b828684SBarry Smith /*@C
31417bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3155e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3165e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3175e76c431SBarry Smith 
318c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
319c7afd0dbSLois Curfman McInnes 
3205e76c431SBarry Smith    Input Parameters:
321c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
322af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3235e76c431SBarry Smith .  x - current iterate
3245e76c431SBarry Smith .  f - residual evaluated at x
3253c632250SBarry Smith .  y - search direction
3265e76c431SBarry Smith .  w - work vector
327c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3285e76c431SBarry Smith 
329c4a48953SLois Curfman McInnes    Output Parameters:
330c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3313c632250SBarry Smith .  w - new iterate
3325e76c431SBarry Smith .  gnorm - 2-norm of g
3335e76c431SBarry Smith .  ynorm - 2-norm of search length
334a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
335fee21e36SBarry Smith 
336c4a48953SLois Curfman McInnes    Options Database Key:
33717bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
338c4a48953SLois Curfman McInnes 
33936851e7fSLois Curfman McInnes    Level: advanced
34036851e7fSLois Curfman McInnes 
34128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
34228ae5a21SLois Curfman McInnes 
34317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
34417bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3455e76c431SBarry Smith @*/
34617bae607SBarry 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)
3475e76c431SBarry Smith {
348dfbe8321SBarry Smith   PetscErrorCode ierr;
3494b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3503c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3515334005bSBarry Smith 
3523a40ed3dSBarry Smith   PetscFunctionBegin;
353a7cc72afSBarry Smith   *flag = PETSC_TRUE;
354d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
355a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
35679f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3573c632250SBarry Smith   if (neP->postcheckstep) {
3583c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3598f99978dSLois Curfman McInnes   }
3603c632250SBarry Smith   if (changed_y) {
36179f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3623c632250SBarry Smith   }
3633c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
364d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
36519717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
36619717074SBarry Smith   }
367d5e45103SBarry Smith   CHKERRQ(ierr);
368d5e45103SBarry Smith 
369a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
37043e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
371d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3723a40ed3dSBarry Smith   PetscFunctionReturn(0);
3735e76c431SBarry Smith }
37404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3754a2ae208SSatish Balay #undef __FUNCT__
37617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
37782bf6240SBarry Smith 
37829e0b56fSBarry Smith /*@C
37917bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
38029e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
38129e0b56fSBarry Smith    even compute the norm of the function or search direction; this
38229e0b56fSBarry Smith    is intended only when you know the full step is fine and are
38329e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
384c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
385c7afd0dbSLois Curfman McInnes 
386c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
38729e0b56fSBarry Smith 
38829e0b56fSBarry Smith    Input Parameters:
389c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
390af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
39129e0b56fSBarry Smith .  x - current iterate
39229e0b56fSBarry Smith .  f - residual evaluated at x
3933c632250SBarry Smith .  y - search direction
39429e0b56fSBarry Smith .  w - work vector
395c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
39629e0b56fSBarry Smith 
39729e0b56fSBarry Smith    Output Parameters:
398c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3993c632250SBarry Smith .  w - new iterate
40029e0b56fSBarry Smith .  gnorm - not changed
40129e0b56fSBarry Smith .  ynorm - not changed
402a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
403fee21e36SBarry Smith 
40429e0b56fSBarry Smith    Options Database Key:
40517bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
40629e0b56fSBarry Smith 
4078cbba510SBarry Smith    Notes:
40817bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
409ea56f5baSLois Curfman McInnes    either the options
410ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
411ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
412329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
413329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
414329f5518SBarry Smith 
415329f5518SBarry Smith    During the final iteration this will not evaluate the function at
416329f5518SBarry Smith    the solution point. This is to save a function evaluation while
417329f5518SBarry Smith    using pseudo-timestepping.
4188cbba510SBarry Smith 
419ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
420*a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
421ea56f5baSLois Curfman McInnes    correct, since they are not computed.
422ea56f5baSLois Curfman McInnes 
423ea56f5baSLois Curfman McInnes    Level: advanced
4248cbba510SBarry Smith 
42529e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
42629e0b56fSBarry Smith 
42717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
42817bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
42929e0b56fSBarry Smith @*/
43017bae607SBarry 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)
43129e0b56fSBarry Smith {
432dfbe8321SBarry Smith   PetscErrorCode ierr;
4334b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4343c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
43529e0b56fSBarry Smith 
4363a40ed3dSBarry Smith   PetscFunctionBegin;
437a7cc72afSBarry Smith   *flag = PETSC_TRUE;
438d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
43979f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4403c632250SBarry Smith   if (neP->postcheckstep) {
4413c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4423c632250SBarry Smith   }
4433c632250SBarry Smith   if (changed_y) {
44479f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4458f99978dSLois Curfman McInnes   }
446329f5518SBarry Smith 
447329f5518SBarry Smith   /* don't evaluate function the last time through */
448329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4493c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
450d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
45119717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
45219717074SBarry Smith     }
453d5e45103SBarry Smith     CHKERRQ(ierr);
454329f5518SBarry Smith   }
455d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4563a40ed3dSBarry Smith   PetscFunctionReturn(0);
45729e0b56fSBarry Smith }
45804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4594a2ae208SSatish Balay #undef __FUNCT__
46017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4614b828684SBarry Smith /*@C
46217bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4635e76c431SBarry Smith 
464c7afd0dbSLois Curfman McInnes    Collective on SNES
465c7afd0dbSLois Curfman McInnes 
4665e76c431SBarry Smith    Input Parameters:
467c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
468af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4695e76c431SBarry Smith .  x - current iterate
4705e76c431SBarry Smith .  f - residual evaluated at x
4713c632250SBarry Smith .  y - search direction
4725e76c431SBarry Smith .  w - work vector
473c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4745e76c431SBarry Smith 
475393d2d9aSLois Curfman McInnes    Output Parameters:
476c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4773c632250SBarry Smith .  w - new iterate
4785e76c431SBarry Smith .  gnorm - 2-norm of g
4795e76c431SBarry Smith .  ynorm - 2-norm of search length
480a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
481fee21e36SBarry Smith 
482c4a48953SLois Curfman McInnes    Options Database Key:
48317bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
484c4a48953SLois Curfman McInnes 
4855e76c431SBarry Smith    Notes:
4865e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4875e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4885e76c431SBarry Smith 
48936851e7fSLois Curfman McInnes    Level: advanced
49036851e7fSLois Curfman McInnes 
49128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49228ae5a21SLois Curfman McInnes 
49317bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
4945e76c431SBarry Smith @*/
49517bae607SBarry 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)
4965e76c431SBarry Smith {
497406556e6SLois Curfman McInnes   /*
498406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
499406556e6SLois Curfman McInnes      minimization problem:
500406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
501406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
502406556e6SLois Curfman McInnes    */
503406556e6SLois Curfman McInnes 
5045ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
505275d25c3SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
506aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
507ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5086b5873e3SBarry Smith #endif
509dfbe8321SBarry Smith   PetscErrorCode ierr;
510a7cc72afSBarry Smith   PetscInt       count;
5114b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
51279f36460SBarry Smith   PetscScalar    scale;
5133c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5145e76c431SBarry Smith 
5153a40ed3dSBarry Smith   PetscFunctionBegin;
516d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
517a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5185e76c431SBarry Smith   alpha   = neP->alpha;
5195e76c431SBarry Smith   maxstep = neP->maxstep;
5205e76c431SBarry Smith   steptol = neP->steptol;
5215e76c431SBarry Smith 
522cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5239e247f21SBarry Smith   if (!*ynorm) {
524ae15b995SBarry Smith     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
525a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5263c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
527a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
528a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
529ad922ac9SBarry Smith     goto theend1;
53094a9d846SBarry Smith   }
5315e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5325e42470aSBarry Smith     scale = maxstep/(*ynorm);
533ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr);
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 
549e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
55043e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
551ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
55243e71028SBarry Smith     *flag = PETSC_FALSE;
55343e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
55443e71028SBarry Smith     goto theend1;
55543e71028SBarry Smith   }
55619717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
557d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
55819717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
55919717074SBarry Smith   }
560d5e45103SBarry Smith   CHKERRQ(ierr);
561313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56243e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5635d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
564ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
565ad922ac9SBarry Smith     goto theend1;
5665e76c431SBarry Smith   }
5675e76c431SBarry Smith 
5685e76c431SBarry Smith   /* Fit points with quadratic */
569313b5538SBarry Smith   lambda     = 1.0;
570a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5715e76c431SBarry Smith   lambdaprev = lambda;
5725e76c431SBarry Smith   gnormprev  = *gnorm;
573329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
574ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
575ddd12767SBarry Smith   else                         lambda = lambdatemp;
576275d25c3SBarry Smith 
577aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
578a23f76dfSSatish Balay   clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
5795e42470aSBarry Smith #else
580e68848bdSBarry Smith   ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
5815e42470aSBarry Smith #endif
58243e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
583ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
58443e71028SBarry Smith     *flag = PETSC_FALSE;
58543e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
58643e71028SBarry Smith     goto theend1;
58743e71028SBarry Smith   }
58819717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
589d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59019717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59119717074SBarry Smith   }
592d5e45103SBarry Smith   CHKERRQ(ierr);
593cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
59443e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5955ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
596ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
597ad922ac9SBarry Smith     goto theend1;
5985e76c431SBarry Smith   }
5995e76c431SBarry Smith 
6005e76c431SBarry Smith   /* Fit points with cubic */
6015e76c431SBarry Smith   count = 1;
6028229bfc2SMatthew Knepley   while (count < 10000) {
6035e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
604ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
605ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
60643e71028SBarry Smith       *flag = PETSC_FALSE;
60743e71028SBarry Smith       break;
6085e76c431SBarry Smith     }
6095d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6105d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
611ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6122b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6135e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6145e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6155e76c431SBarry Smith     if (a == 0.0) {
6165e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6175e76c431SBarry Smith     } else {
6185e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6195e76c431SBarry Smith     }
6205e76c431SBarry Smith     lambdaprev = lambda;
6215e76c431SBarry Smith     gnormprev  = *gnorm;
622329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
623329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6245e76c431SBarry Smith     else                         lambda     = lambdatemp;
625aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
626275d25c3SBarry Smith     clambda   = lambda;
627e68848bdSBarry Smith     ierr      = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
6285e42470aSBarry Smith #else
629e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
6305e42470aSBarry Smith #endif
63143e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
632ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
633ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
634ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
63543e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
636ed98166eSMatthew Knepley       break;
637ed98166eSMatthew Knepley     }
63819717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
639d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64019717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64119717074SBarry Smith     }
642d5e45103SBarry Smith     CHKERRQ(ierr);
643cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
64443e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6455ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
646ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
64743e71028SBarry Smith       break;
6482b022350SLois Curfman McInnes     } else {
649ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);CHKERRQ(ierr);
6505e76c431SBarry Smith     }
6515e76c431SBarry Smith     count++;
6525e76c431SBarry Smith   }
6538229bfc2SMatthew Knepley   if (count >= 10000) {
654cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6558229bfc2SMatthew Knepley   }
656ad922ac9SBarry Smith   theend1:
6578f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6583c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6593c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6603c632250SBarry Smith     if (changed_y) {
66179f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6623c632250SBarry Smith     }
6633c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6643c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
665d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
66619717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
66719717074SBarry Smith       }
668d5e45103SBarry Smith       CHKERRQ(ierr);
6698f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67043e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6713c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6728f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6733c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6748f99978dSLois Curfman McInnes     }
6758f99978dSLois Curfman McInnes   }
676d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6773a40ed3dSBarry Smith   PetscFunctionReturn(0);
6785e76c431SBarry Smith }
67904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6804a2ae208SSatish Balay #undef __FUNCT__
68117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6824b828684SBarry Smith /*@C
68317bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6845e76c431SBarry Smith 
685c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
686c7afd0dbSLois Curfman McInnes 
6875e42470aSBarry Smith    Input Parameters:
688c7afd0dbSLois Curfman McInnes +  snes - the SNES context
689af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6905e42470aSBarry Smith .  x - current iterate
6915e42470aSBarry Smith .  f - residual evaluated at x
6923c632250SBarry Smith .  y - search direction
6935e42470aSBarry Smith .  w - work vector
694c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6955e42470aSBarry Smith 
696c4a48953SLois Curfman McInnes    Output Parameters:
6977f3332b4SBarry Smith +  g - residual evaluated at new iterate w
6987f3332b4SBarry Smith .  w - new iterate (x + alpha*y)
6995e42470aSBarry Smith .  gnorm - 2-norm of g
7005e42470aSBarry Smith .  ynorm - 2-norm of search length
701a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
702fee21e36SBarry Smith 
703c4a48953SLois Curfman McInnes    Options Database Key:
70417bae607SBarry Smith .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
705c4a48953SLois Curfman McInnes 
7065e42470aSBarry Smith    Notes:
70717bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7085e42470aSBarry Smith 
70936851e7fSLois Curfman McInnes    Level: advanced
71036851e7fSLois Curfman McInnes 
711f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
712f59ffdedSLois Curfman McInnes 
71317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7145e42470aSBarry Smith @*/
71517bae607SBarry 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)
7165e76c431SBarry Smith {
717406556e6SLois Curfman McInnes   /*
718406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
719406556e6SLois Curfman McInnes      minimization problem:
720406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
721406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
722406556e6SLois Curfman McInnes    */
723275d25c3SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
724aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
725ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7266b5873e3SBarry Smith #endif
727dfbe8321SBarry Smith   PetscErrorCode ierr;
728a7cc72afSBarry Smith   PetscInt       count;
7294b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
73079f36460SBarry Smith   PetscScalar    scale;
7313c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7325e76c431SBarry Smith 
7333a40ed3dSBarry Smith   PetscFunctionBegin;
734d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
735a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7365e42470aSBarry Smith   alpha   = neP->alpha;
7375e42470aSBarry Smith   maxstep = neP->maxstep;
7385e42470aSBarry Smith   steptol = neP->steptol;
7395e76c431SBarry Smith 
7403505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74163b9fa5eSBarry Smith   if (*ynorm == 0.0) {
742ae15b995SBarry Smith     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
743b37302c6SLois Curfman McInnes     *gnorm = fnorm;
744e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
745b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
746a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
747ad922ac9SBarry Smith     goto theend2;
74894a9d846SBarry Smith   }
7495e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7505e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7512dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7525e42470aSBarry Smith     *ynorm = maxstep;
7535e76c431SBarry Smith   }
754ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7555ca10a37SBarry Smith   minlambda = steptol/rellength;
756a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
757aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
758a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
759329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7605e42470aSBarry Smith #else
761a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7625e42470aSBarry Smith #endif
7635e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7645e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7655e42470aSBarry Smith 
766e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
76743e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
768ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
76943e71028SBarry Smith     *flag = PETSC_FALSE;
77043e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77143e71028SBarry Smith     goto theend2;
77243e71028SBarry Smith   }
77319717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
774d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
77519717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
77619717074SBarry Smith   }
777d5e45103SBarry Smith   CHKERRQ(ierr);
778cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
77943e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7805d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
781ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
782ad922ac9SBarry Smith     goto theend2;
7835e42470aSBarry Smith   }
7845e42470aSBarry Smith 
7855e42470aSBarry Smith   /* Fit points with quadratic */
786313b5538SBarry Smith   lambda = 1.0;
7875e42470aSBarry Smith   count = 1;
7885ca10a37SBarry Smith   while (PETSC_TRUE) {
7895e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
790ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
791ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
792e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
79343e71028SBarry Smith       *flag = PETSC_FALSE;
79443e71028SBarry Smith       break;
7955e42470aSBarry Smith     }
796a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
797329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
798329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
799329f5518SBarry Smith     else                         lambda     = lambdatemp;
800275d25c3SBarry Smith 
801aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
802e68848bdSBarry Smith     clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
8035e42470aSBarry Smith #else
804e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8055e42470aSBarry Smith #endif
80643e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
807ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
808ae15b995SBarry Smith       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
809ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
81043e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
811ed98166eSMatthew Knepley       break;
812ed98166eSMatthew Knepley     }
81319717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
814d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
81519717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
81619717074SBarry Smith     }
817d5e45103SBarry Smith     CHKERRQ(ierr);
818cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
81943e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8205ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
821ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
82206259719SBarry Smith       break;
8235e42470aSBarry Smith     }
8245e42470aSBarry Smith     count++;
8255e42470aSBarry Smith   }
826ad922ac9SBarry Smith   theend2:
8278f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8283c632250SBarry Smith   if (neP->postcheckstep) {
8293c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8303c632250SBarry Smith     if (changed_y) {
83179f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8323c632250SBarry Smith     }
8333c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8343c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
835d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
83619717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
83719717074SBarry Smith       }
838d5e45103SBarry Smith       CHKERRQ(ierr);
8398f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8403c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8418f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8423c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8433c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8448f99978dSLois Curfman McInnes     }
8458f99978dSLois Curfman McInnes   }
846d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8473a40ed3dSBarry Smith   PetscFunctionReturn(0);
8485e76c431SBarry Smith }
8492343ba6eSBarry Smith 
85004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8514a2ae208SSatish Balay #undef __FUNCT__
85217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
853c9e6a524SLois Curfman McInnes /*@C
85417bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8554b27c08aSLois Curfman McInnes    by the method SNESLS.
8565e76c431SBarry Smith 
8575e76c431SBarry Smith    Input Parameters:
858c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
859af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
860c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8615e76c431SBarry Smith 
862fee21e36SBarry Smith    Collective on SNES
863fee21e36SBarry Smith 
864c4a48953SLois Curfman McInnes    Available Routines:
86517bae607SBarry Smith +  SNESLineSearchCubic() - default line search
86617bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
86717bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
86817bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8695e76c431SBarry Smith 
870c4a48953SLois Curfman McInnes     Options Database Keys:
8714b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8724b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8734b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8744b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8753304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8763304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
877c4a48953SLois Curfman McInnes 
8785e76c431SBarry Smith    Calling sequence of func:
879bd208895SLois Curfman McInnes .vb
880af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
881a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
882bd208895SLois Curfman McInnes .ve
8835e76c431SBarry Smith 
8845e76c431SBarry Smith     Input parameters for func:
885c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
886af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8875e76c431SBarry Smith .   x - current iterate
8885e76c431SBarry Smith .   f - residual evaluated at x
8893c632250SBarry Smith .   y - search direction
8905e76c431SBarry Smith .   w - work vector
891c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8925e76c431SBarry Smith 
8935e76c431SBarry Smith     Output parameters for func:
894c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8953c632250SBarry Smith .   w - new iterate
8965e76c431SBarry Smith .   gnorm - 2-norm of g
8975e76c431SBarry Smith .   ynorm - 2-norm of search length
898a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
899f59ffdedSLois Curfman McInnes 
90036851e7fSLois Curfman McInnes     Level: advanced
90136851e7fSLois Curfman McInnes 
902f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
903f59ffdedSLois Curfman McInnes 
90417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
90517bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9065e76c431SBarry Smith @*/
90717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9085e76c431SBarry Smith {
909a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
91082bf6240SBarry Smith 
9113a40ed3dSBarry Smith   PetscFunctionBegin;
91217bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
91382bf6240SBarry Smith   if (f) {
914af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
91582bf6240SBarry Smith   }
9163a40ed3dSBarry Smith   PetscFunctionReturn(0);
9175e76c431SBarry Smith }
9188e019c35SBarry Smith 
919a7cc72afSBarry 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*/
92004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
921fb2e594dSBarry Smith EXTERN_C_BEGIN
9224a2ae208SSatish Balay #undef __FUNCT__
92317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
92417bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
92582bf6240SBarry Smith {
92682bf6240SBarry Smith   PetscFunctionBegin;
9274b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9284b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
92982bf6240SBarry Smith   PetscFunctionReturn(0);
93082bf6240SBarry Smith }
931fb2e594dSBarry Smith EXTERN_C_END
93204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9334a2ae208SSatish Balay #undef __FUNCT__
9343c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
935c8dd0c0dSLois Curfman McInnes /*@C
9363c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9374b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
938c8dd0c0dSLois Curfman McInnes 
939c8dd0c0dSLois Curfman McInnes    Input Parameters:
940c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9413c632250SBarry Smith .  func - pointer to function
942c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
943c8dd0c0dSLois Curfman McInnes 
944c8dd0c0dSLois Curfman McInnes    Collective on SNES
945c8dd0c0dSLois Curfman McInnes 
946c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
947c8dd0c0dSLois Curfman McInnes .vb
9483c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
949c8dd0c0dSLois Curfman McInnes .ve
950b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
951b1ae27eaSLois Curfman McInnes    on failure.
952c8dd0c0dSLois Curfman McInnes 
953c8dd0c0dSLois Curfman McInnes    Input parameters for func:
954c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
955c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9563c632250SBarry Smith .  x - previous iterate
9573c632250SBarry Smith .  y - new search direction and length
9583c632250SBarry Smith -  w - current candidate iterate
959c8dd0c0dSLois Curfman McInnes 
960c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9613c632250SBarry Smith +  y - search direction (possibly changed)
9623c632250SBarry Smith .  w - current iterate (possibly modified)
9633c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9643c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
965c8dd0c0dSLois Curfman McInnes 
966c8dd0c0dSLois Curfman McInnes    Level: advanced
967c8dd0c0dSLois Curfman McInnes 
9689e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9699e247f21SBarry Smith 
9703c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9713c632250SBarry Smith 
9723c632250SBarry Smith    On input w = x + y
9733c632250SBarry Smith 
97417bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
975b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
976ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9778f99978dSLois Curfman McInnes 
97817bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
979ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
980ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
981ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9829e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
983b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9848f99978dSLois Curfman McInnes 
985c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
986c8dd0c0dSLois Curfman McInnes 
98717bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
988c8dd0c0dSLois Curfman McInnes @*/
9893c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
990c8dd0c0dSLois Curfman McInnes {
9913c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
992c8dd0c0dSLois Curfman McInnes 
993c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9943c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
995c8dd0c0dSLois Curfman McInnes   if (f) {
996c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
997c8dd0c0dSLois Curfman McInnes   }
998c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
999c8dd0c0dSLois Curfman McInnes }
10009c3ca13aSBarry Smith #undef __FUNCT__
10019c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10029c3ca13aSBarry Smith /*@C
10039c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10047e4bb74cSBarry Smith          before the line search is called.
10059c3ca13aSBarry Smith 
10069c3ca13aSBarry Smith    Input Parameters:
10079c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10089c3ca13aSBarry Smith .  func - pointer to function
10099c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10109c3ca13aSBarry Smith 
10119c3ca13aSBarry Smith    Collective on SNES
10129c3ca13aSBarry Smith 
10139c3ca13aSBarry Smith    Calling sequence of func:
10149c3ca13aSBarry Smith .vb
10151e3347e8SBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
10169c3ca13aSBarry Smith .ve
10179c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10189c3ca13aSBarry Smith    on failure.
10199c3ca13aSBarry Smith 
10209c3ca13aSBarry Smith    Input parameters for func:
10219c3ca13aSBarry Smith +  snes - nonlinear context
10229c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10239c3ca13aSBarry Smith .  x - previous iterate
10249c3ca13aSBarry Smith -  y - new search direction and length
10259c3ca13aSBarry Smith 
10269c3ca13aSBarry Smith    Output parameters for func:
10279c3ca13aSBarry Smith +  y - search direction (possibly changed)
10289c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10299c3ca13aSBarry Smith 
10309c3ca13aSBarry Smith    Level: advanced
10319c3ca13aSBarry Smith 
10329c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10339c3ca13aSBarry Smith 
10349c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10359c3ca13aSBarry Smith 
10367e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10379c3ca13aSBarry Smith @*/
10389c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10399c3ca13aSBarry Smith {
10409c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10419c3ca13aSBarry Smith 
10429c3ca13aSBarry Smith   PetscFunctionBegin;
10439c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10449c3ca13aSBarry Smith   if (f) {
10459c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10469c3ca13aSBarry Smith   }
10479c3ca13aSBarry Smith   PetscFunctionReturn(0);
10489c3ca13aSBarry Smith }
10499c3ca13aSBarry Smith 
1050c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10519c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10529c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1053c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10544a2ae208SSatish Balay #undef __FUNCT__
10553c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10569c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1057c8dd0c0dSLois Curfman McInnes {
1058c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10593c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10603c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1061c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1062c8dd0c0dSLois Curfman McInnes }
1063c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10649c3ca13aSBarry Smith 
10659c3ca13aSBarry Smith EXTERN_C_BEGIN
10669c3ca13aSBarry Smith #undef __FUNCT__
10679c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10689c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10699c3ca13aSBarry Smith {
10709c3ca13aSBarry Smith   PetscFunctionBegin;
10719c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10729c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10739c3ca13aSBarry Smith   PetscFunctionReturn(0);
10749c3ca13aSBarry Smith }
10759c3ca13aSBarry Smith EXTERN_C_END
1076329e7e40SMatthew Knepley 
1077329e7e40SMatthew Knepley /*
10784b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
107904d965bbSLois Curfman McInnes 
108004d965bbSLois Curfman McInnes    Input Parameters:
108104d965bbSLois Curfman McInnes .  SNES - the SNES context
108204d965bbSLois Curfman McInnes .  viewer - visualization context
108304d965bbSLois Curfman McInnes 
108404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
108504d965bbSLois Curfman McInnes */
10864a2ae208SSatish Balay #undef __FUNCT__
10874b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10886849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1089a935fc98SLois Curfman McInnes {
10904b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
10912fc52814SBarry Smith   const char     *cstr;
1092dfbe8321SBarry Smith   PetscErrorCode ierr;
109332077d6dSBarry Smith   PetscTruth     iascii;
1094a935fc98SLois Curfman McInnes 
10953a40ed3dSBarry Smith   PetscFunctionBegin;
109632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
109732077d6dSBarry Smith   if (iascii) {
109817bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
109917bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
110017bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
110119bcc07fSBarry Smith     else                                                cstr = "unknown";
1102b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1103a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11045cd90555SBarry Smith   } else {
110579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
110619bcc07fSBarry Smith   }
11073a40ed3dSBarry Smith   PetscFunctionReturn(0);
1108a935fc98SLois Curfman McInnes }
110904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
111004d965bbSLois Curfman McInnes /*
11114b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
111204d965bbSLois Curfman McInnes 
111304d965bbSLois Curfman McInnes    Input Parameter:
111404d965bbSLois Curfman McInnes .  snes - the SNES context
111504d965bbSLois Curfman McInnes 
111604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
111704d965bbSLois Curfman McInnes */
11184a2ae208SSatish Balay #undef __FUNCT__
11194b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11206849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11215e42470aSBarry Smith {
11224b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1123e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1124dfbe8321SBarry Smith   PetscErrorCode ierr;
1125a7cc72afSBarry Smith   PetscInt       indx;
1126f1af5d2fSBarry Smith   PetscTruth     flg;
11275e42470aSBarry Smith 
11283a40ed3dSBarry Smith   PetscFunctionBegin;
1129b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11304b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11314b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11324b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1133186905e3SBarry Smith 
113417bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
113525cdf11fSBarry Smith     if (flg) {
113622e36657SBarry Smith       switch (indx) {
1137b49fd9e1SBarry Smith       case 0:
113817bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1139b49fd9e1SBarry Smith         break;
1140b49fd9e1SBarry Smith       case 1:
114117bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1142b49fd9e1SBarry Smith         break;
1143b49fd9e1SBarry Smith       case 2:
114417bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1145b49fd9e1SBarry Smith         break;
1146b49fd9e1SBarry Smith       case 3:
114717bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1148b49fd9e1SBarry Smith         break;
11495e42470aSBarry Smith       }
11505e42470aSBarry Smith     }
1151b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11523a40ed3dSBarry Smith   PetscFunctionReturn(0);
11535e42470aSBarry Smith }
115404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1155ebe3b25bSBarry Smith /*MC
1156ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
115704d965bbSLois Curfman McInnes 
1158ebe3b25bSBarry Smith    Options Database:
1159ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1160ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1161ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1162ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1163ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1164ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
116504d965bbSLois Curfman McInnes 
1166ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
116704d965bbSLois Curfman McInnes 
1168ee3001cbSBarry Smith    Level: beginner
1169ee3001cbSBarry Smith 
117017bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
117117bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
117217bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1173ebe3b25bSBarry Smith 
1174ebe3b25bSBarry Smith M*/
1175fb2e594dSBarry Smith EXTERN_C_BEGIN
11764a2ae208SSatish Balay #undef __FUNCT__
11774b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
117863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
11795e42470aSBarry Smith {
1180dfbe8321SBarry Smith   PetscErrorCode ierr;
11814b27c08aSLois Curfman McInnes   SNES_LS        *neP;
11825e42470aSBarry Smith 
11833a40ed3dSBarry Smith   PetscFunctionBegin;
1184e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1185e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1186e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1187e7788613SBarry Smith   snes->ops->converged	     = SNESConverged_LS;
1188e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1189e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
11905baf8537SBarry Smith   snes->nwork                = 0;
11915e42470aSBarry Smith 
11924b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
119352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
11945e42470aSBarry Smith   snes->data    	= (void*)neP;
11955e42470aSBarry Smith   neP->alpha		= 1.e-4;
11965e42470aSBarry Smith   neP->maxstep		= 1.e8;
11975e42470aSBarry Smith   neP->steptol		= 1.e-12;
119817bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1199c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12003c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12013c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12023c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12033c632250SBarry Smith   neP->precheck         = PETSC_NULL;
120482bf6240SBarry Smith 
120517bae607SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
12063c632250SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
12079c3ca13aSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
120882bf6240SBarry Smith 
12093a40ed3dSBarry Smith   PetscFunctionReturn(0);
12105e42470aSBarry Smith }
1211fb2e594dSBarry Smith EXTERN_C_END
121282bf6240SBarry Smith 
121382bf6240SBarry Smith 
121482bf6240SBarry Smith 
1215