xref: /petsc/src/snes/impls/ls/ls.c (revision 71f87433c2d5ca8a2055e7bac4e3ff8159d57028)
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);
179*71f87433Sdalcinl     ierr = SNES_KSPSolve(snes,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);
307557d3f75SLisandro Dalcin 
308557d3f75SLisandro Dalcin   /* clear composed functions */
309557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
310557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
311557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
312557d3f75SLisandro Dalcin 
3133a40ed3dSBarry Smith   PetscFunctionReturn(0);
3145e76c431SBarry Smith }
31504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3164a2ae208SSatish Balay #undef __FUNCT__
31717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
31882bf6240SBarry Smith 
3194b828684SBarry Smith /*@C
32017bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3215e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3225e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3235e76c431SBarry Smith 
324c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
325c7afd0dbSLois Curfman McInnes 
3265e76c431SBarry Smith    Input Parameters:
327c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
328af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3295e76c431SBarry Smith .  x - current iterate
3305e76c431SBarry Smith .  f - residual evaluated at x
3313c632250SBarry Smith .  y - search direction
3325e76c431SBarry Smith .  w - work vector
333c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3345e76c431SBarry Smith 
335c4a48953SLois Curfman McInnes    Output Parameters:
336c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3373c632250SBarry Smith .  w - new iterate
3385e76c431SBarry Smith .  gnorm - 2-norm of g
3395e76c431SBarry Smith .  ynorm - 2-norm of search length
340a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
341fee21e36SBarry Smith 
342c4a48953SLois Curfman McInnes    Options Database Key:
34317bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
344c4a48953SLois Curfman McInnes 
34536851e7fSLois Curfman McInnes    Level: advanced
34636851e7fSLois Curfman McInnes 
34728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
34828ae5a21SLois Curfman McInnes 
34917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
35017bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3515e76c431SBarry Smith @*/
35217bae607SBarry 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)
3535e76c431SBarry Smith {
354dfbe8321SBarry Smith   PetscErrorCode ierr;
3554b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3563c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3575334005bSBarry Smith 
3583a40ed3dSBarry Smith   PetscFunctionBegin;
359a7cc72afSBarry Smith   *flag = PETSC_TRUE;
360d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
361a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
36279f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3633c632250SBarry Smith   if (neP->postcheckstep) {
3643c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3658f99978dSLois Curfman McInnes   }
3663c632250SBarry Smith   if (changed_y) {
36779f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3683c632250SBarry Smith   }
3693c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
370d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
37119717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
37219717074SBarry Smith   }
373d5e45103SBarry Smith   CHKERRQ(ierr);
374d5e45103SBarry Smith 
375a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
37643e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
377d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3783a40ed3dSBarry Smith   PetscFunctionReturn(0);
3795e76c431SBarry Smith }
38004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3814a2ae208SSatish Balay #undef __FUNCT__
38217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
38382bf6240SBarry Smith 
38429e0b56fSBarry Smith /*@C
38517bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
38629e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
38729e0b56fSBarry Smith    even compute the norm of the function or search direction; this
38829e0b56fSBarry Smith    is intended only when you know the full step is fine and are
38929e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
390c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
391c7afd0dbSLois Curfman McInnes 
392c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
39329e0b56fSBarry Smith 
39429e0b56fSBarry Smith    Input Parameters:
395c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
396af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
39729e0b56fSBarry Smith .  x - current iterate
39829e0b56fSBarry Smith .  f - residual evaluated at x
3993c632250SBarry Smith .  y - search direction
40029e0b56fSBarry Smith .  w - work vector
401c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
40229e0b56fSBarry Smith 
40329e0b56fSBarry Smith    Output Parameters:
404c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4053c632250SBarry Smith .  w - new iterate
40629e0b56fSBarry Smith .  gnorm - not changed
40729e0b56fSBarry Smith .  ynorm - not changed
408a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
409fee21e36SBarry Smith 
41029e0b56fSBarry Smith    Options Database Key:
41117bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
41229e0b56fSBarry Smith 
4138cbba510SBarry Smith    Notes:
41417bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
415ea56f5baSLois Curfman McInnes    either the options
416ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
417ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
418329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
419329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
420329f5518SBarry Smith 
421329f5518SBarry Smith    During the final iteration this will not evaluate the function at
422329f5518SBarry Smith    the solution point. This is to save a function evaluation while
423329f5518SBarry Smith    using pseudo-timestepping.
4248cbba510SBarry Smith 
425ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
426a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
427ea56f5baSLois Curfman McInnes    correct, since they are not computed.
428ea56f5baSLois Curfman McInnes 
429ea56f5baSLois Curfman McInnes    Level: advanced
4308cbba510SBarry Smith 
43129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
43229e0b56fSBarry Smith 
43317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
43417bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
43529e0b56fSBarry Smith @*/
43617bae607SBarry 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)
43729e0b56fSBarry Smith {
438dfbe8321SBarry Smith   PetscErrorCode ierr;
4394b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4403c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
44129e0b56fSBarry Smith 
4423a40ed3dSBarry Smith   PetscFunctionBegin;
443a7cc72afSBarry Smith   *flag = PETSC_TRUE;
444d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
44579f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4463c632250SBarry Smith   if (neP->postcheckstep) {
4473c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4483c632250SBarry Smith   }
4493c632250SBarry Smith   if (changed_y) {
45079f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4518f99978dSLois Curfman McInnes   }
452329f5518SBarry Smith 
453329f5518SBarry Smith   /* don't evaluate function the last time through */
454329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4553c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
456d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
45719717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
45819717074SBarry Smith     }
459d5e45103SBarry Smith     CHKERRQ(ierr);
460329f5518SBarry Smith   }
461d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4623a40ed3dSBarry Smith   PetscFunctionReturn(0);
46329e0b56fSBarry Smith }
46404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4654a2ae208SSatish Balay #undef __FUNCT__
46617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4674b828684SBarry Smith /*@C
46817bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4695e76c431SBarry Smith 
470c7afd0dbSLois Curfman McInnes    Collective on SNES
471c7afd0dbSLois Curfman McInnes 
4725e76c431SBarry Smith    Input Parameters:
473c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
474af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4755e76c431SBarry Smith .  x - current iterate
4765e76c431SBarry Smith .  f - residual evaluated at x
4773c632250SBarry Smith .  y - search direction
4785e76c431SBarry Smith .  w - work vector
479c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4805e76c431SBarry Smith 
481393d2d9aSLois Curfman McInnes    Output Parameters:
482c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4833c632250SBarry Smith .  w - new iterate
4845e76c431SBarry Smith .  gnorm - 2-norm of g
4855e76c431SBarry Smith .  ynorm - 2-norm of search length
486a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
487fee21e36SBarry Smith 
488c4a48953SLois Curfman McInnes    Options Database Key:
48917bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
490c4a48953SLois Curfman McInnes 
4915e76c431SBarry Smith    Notes:
4925e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4935e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4945e76c431SBarry Smith 
49536851e7fSLois Curfman McInnes    Level: advanced
49636851e7fSLois Curfman McInnes 
49728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49828ae5a21SLois Curfman McInnes 
49917bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5005e76c431SBarry Smith @*/
50117bae607SBarry 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)
5025e76c431SBarry Smith {
503406556e6SLois Curfman McInnes   /*
504406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
505406556e6SLois Curfman McInnes      minimization problem:
506406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
507406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
508406556e6SLois Curfman McInnes    */
509406556e6SLois Curfman McInnes 
5105ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
511275d25c3SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
512aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
513ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5146b5873e3SBarry Smith #endif
515dfbe8321SBarry Smith   PetscErrorCode ierr;
516a7cc72afSBarry Smith   PetscInt       count;
5174b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
51879f36460SBarry Smith   PetscScalar    scale;
5193c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5205e76c431SBarry Smith 
5213a40ed3dSBarry Smith   PetscFunctionBegin;
522d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
523a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5245e76c431SBarry Smith   alpha   = neP->alpha;
5255e76c431SBarry Smith   maxstep = neP->maxstep;
5265e76c431SBarry Smith   steptol = neP->steptol;
5275e76c431SBarry Smith 
528cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5299e247f21SBarry Smith   if (!*ynorm) {
530ae15b995SBarry Smith     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
531a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5323c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
533a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
534a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
535ad922ac9SBarry Smith     goto theend1;
53694a9d846SBarry Smith   }
5375e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5385e42470aSBarry Smith     scale = maxstep/(*ynorm);
539ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr);
5402dcb1b2aSMatthew Knepley     ierr = VecScale(y,scale);CHKERRQ(ierr);
5415e76c431SBarry Smith     *ynorm = maxstep;
5425e76c431SBarry Smith   }
543ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5445ca10a37SBarry Smith   minlambda = steptol/rellength;
545a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
546aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
547a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
548329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5495e42470aSBarry Smith #else
550a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5515e42470aSBarry Smith #endif
5525e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5535e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5545e76c431SBarry Smith 
555e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
55643e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
557ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
55843e71028SBarry Smith     *flag = PETSC_FALSE;
55943e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
56043e71028SBarry Smith     goto theend1;
56143e71028SBarry Smith   }
56219717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
563d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
56419717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
56519717074SBarry Smith   }
566d5e45103SBarry Smith   CHKERRQ(ierr);
567313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56843e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5695d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
570ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
571ad922ac9SBarry Smith     goto theend1;
5725e76c431SBarry Smith   }
5735e76c431SBarry Smith 
5745e76c431SBarry Smith   /* Fit points with quadratic */
575313b5538SBarry Smith   lambda     = 1.0;
576a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5775e76c431SBarry Smith   lambdaprev = lambda;
5785e76c431SBarry Smith   gnormprev  = *gnorm;
579329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
580ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
581ddd12767SBarry Smith   else                         lambda = lambdatemp;
582275d25c3SBarry Smith 
583aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
584a23f76dfSSatish Balay   clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
5855e42470aSBarry Smith #else
586e68848bdSBarry Smith   ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
5875e42470aSBarry Smith #endif
58843e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
589ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
59043e71028SBarry Smith     *flag = PETSC_FALSE;
59143e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
59243e71028SBarry Smith     goto theend1;
59343e71028SBarry Smith   }
59419717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
595d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59619717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59719717074SBarry Smith   }
598d5e45103SBarry Smith   CHKERRQ(ierr);
599cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
60043e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6015ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
602ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
603ad922ac9SBarry Smith     goto theend1;
6045e76c431SBarry Smith   }
6055e76c431SBarry Smith 
6065e76c431SBarry Smith   /* Fit points with cubic */
6075e76c431SBarry Smith   count = 1;
6088229bfc2SMatthew Knepley   while (count < 10000) {
6095e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
610ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
611ae15b995SBarry 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);
61243e71028SBarry Smith       *flag = PETSC_FALSE;
61343e71028SBarry Smith       break;
6145e76c431SBarry Smith     }
6155d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6165d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
617ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6182b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6195e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6205e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6215e76c431SBarry Smith     if (a == 0.0) {
6225e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6235e76c431SBarry Smith     } else {
6245e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6255e76c431SBarry Smith     }
6265e76c431SBarry Smith     lambdaprev = lambda;
6275e76c431SBarry Smith     gnormprev  = *gnorm;
628329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
629329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6305e76c431SBarry Smith     else                         lambda     = lambdatemp;
631aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
632275d25c3SBarry Smith     clambda   = lambda;
633e68848bdSBarry Smith     ierr      = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
6345e42470aSBarry Smith #else
635e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
6365e42470aSBarry Smith #endif
63743e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
638ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
639ae15b995SBarry 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);
640ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
64143e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
642ed98166eSMatthew Knepley       break;
643ed98166eSMatthew Knepley     }
64419717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
645d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64619717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64719717074SBarry Smith     }
648d5e45103SBarry Smith     CHKERRQ(ierr);
649cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
65043e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6515ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
652ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
65343e71028SBarry Smith       break;
6542b022350SLois Curfman McInnes     } else {
655ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);CHKERRQ(ierr);
6565e76c431SBarry Smith     }
6575e76c431SBarry Smith     count++;
6585e76c431SBarry Smith   }
6598229bfc2SMatthew Knepley   if (count >= 10000) {
660cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6618229bfc2SMatthew Knepley   }
662ad922ac9SBarry Smith   theend1:
6638f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6643c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6653c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6663c632250SBarry Smith     if (changed_y) {
66779f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6683c632250SBarry Smith     }
6693c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6703c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
671d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
67219717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
67319717074SBarry Smith       }
674d5e45103SBarry Smith       CHKERRQ(ierr);
6758f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67643e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6773c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6788f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6793c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6808f99978dSLois Curfman McInnes     }
6818f99978dSLois Curfman McInnes   }
682d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6833a40ed3dSBarry Smith   PetscFunctionReturn(0);
6845e76c431SBarry Smith }
68504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6864a2ae208SSatish Balay #undef __FUNCT__
68717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6884b828684SBarry Smith /*@C
68917bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6905e76c431SBarry Smith 
691c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
692c7afd0dbSLois Curfman McInnes 
6935e42470aSBarry Smith    Input Parameters:
694c7afd0dbSLois Curfman McInnes +  snes - the SNES context
695af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6965e42470aSBarry Smith .  x - current iterate
6975e42470aSBarry Smith .  f - residual evaluated at x
6983c632250SBarry Smith .  y - search direction
6995e42470aSBarry Smith .  w - work vector
700c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
7015e42470aSBarry Smith 
702c4a48953SLois Curfman McInnes    Output Parameters:
7037f3332b4SBarry Smith +  g - residual evaluated at new iterate w
7047f3332b4SBarry Smith .  w - new iterate (x + alpha*y)
7055e42470aSBarry Smith .  gnorm - 2-norm of g
7065e42470aSBarry Smith .  ynorm - 2-norm of search length
707a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
708fee21e36SBarry Smith 
709c4a48953SLois Curfman McInnes    Options Database Key:
71017bae607SBarry Smith .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
711c4a48953SLois Curfman McInnes 
7125e42470aSBarry Smith    Notes:
71317bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7145e42470aSBarry Smith 
71536851e7fSLois Curfman McInnes    Level: advanced
71636851e7fSLois Curfman McInnes 
717f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
718f59ffdedSLois Curfman McInnes 
71917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7205e42470aSBarry Smith @*/
72117bae607SBarry 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)
7225e76c431SBarry Smith {
723406556e6SLois Curfman McInnes   /*
724406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
725406556e6SLois Curfman McInnes      minimization problem:
726406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
727406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
728406556e6SLois Curfman McInnes    */
729275d25c3SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
730aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
731ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7326b5873e3SBarry Smith #endif
733dfbe8321SBarry Smith   PetscErrorCode ierr;
734a7cc72afSBarry Smith   PetscInt       count;
7354b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
73679f36460SBarry Smith   PetscScalar    scale;
7373c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7385e76c431SBarry Smith 
7393a40ed3dSBarry Smith   PetscFunctionBegin;
740d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
741a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7425e42470aSBarry Smith   alpha   = neP->alpha;
7435e42470aSBarry Smith   maxstep = neP->maxstep;
7445e42470aSBarry Smith   steptol = neP->steptol;
7455e76c431SBarry Smith 
7463505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74763b9fa5eSBarry Smith   if (*ynorm == 0.0) {
748ae15b995SBarry Smith     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
749b37302c6SLois Curfman McInnes     *gnorm = fnorm;
750e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
751b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
752a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
753ad922ac9SBarry Smith     goto theend2;
75494a9d846SBarry Smith   }
7555e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7565e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7572dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7585e42470aSBarry Smith     *ynorm = maxstep;
7595e76c431SBarry Smith   }
760ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7615ca10a37SBarry Smith   minlambda = steptol/rellength;
762a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
763aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
764a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
765329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7665e42470aSBarry Smith #else
767a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7685e42470aSBarry Smith #endif
7695e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7705e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7715e42470aSBarry Smith 
772e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
77343e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
774ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");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 */
787ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
788ad922ac9SBarry Smith     goto theend2;
7895e42470aSBarry Smith   }
7905e42470aSBarry Smith 
7915e42470aSBarry Smith   /* Fit points with quadratic */
792313b5538SBarry Smith   lambda = 1.0;
7935e42470aSBarry Smith   count = 1;
7945ca10a37SBarry Smith   while (PETSC_TRUE) {
7955e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
796ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
797ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
798e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
79943e71028SBarry Smith       *flag = PETSC_FALSE;
80043e71028SBarry Smith       break;
8015e42470aSBarry Smith     }
802a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
803329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
804329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
805329f5518SBarry Smith     else                         lambda     = lambdatemp;
806275d25c3SBarry Smith 
807aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
808e68848bdSBarry Smith     clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
8095e42470aSBarry Smith #else
810e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8115e42470aSBarry Smith #endif
81243e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
813ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
814ae15b995SBarry 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);
815ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
81643e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
817ed98166eSMatthew Knepley       break;
818ed98166eSMatthew Knepley     }
81919717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
820d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
82119717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
82219717074SBarry Smith     }
823d5e45103SBarry Smith     CHKERRQ(ierr);
824cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82543e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8265ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
827ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
82806259719SBarry Smith       break;
8295e42470aSBarry Smith     }
8305e42470aSBarry Smith     count++;
8315e42470aSBarry Smith   }
832ad922ac9SBarry Smith   theend2:
8338f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8343c632250SBarry Smith   if (neP->postcheckstep) {
8353c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8363c632250SBarry Smith     if (changed_y) {
83779f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8383c632250SBarry Smith     }
8393c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8403c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
841d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
84219717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
84319717074SBarry Smith       }
844d5e45103SBarry Smith       CHKERRQ(ierr);
8458f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8463c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8478f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8483c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8493c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8508f99978dSLois Curfman McInnes     }
8518f99978dSLois Curfman McInnes   }
852d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8533a40ed3dSBarry Smith   PetscFunctionReturn(0);
8545e76c431SBarry Smith }
8552343ba6eSBarry Smith 
85604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8574a2ae208SSatish Balay #undef __FUNCT__
85817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
859c9e6a524SLois Curfman McInnes /*@C
86017bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8614b27c08aSLois Curfman McInnes    by the method SNESLS.
8625e76c431SBarry Smith 
8635e76c431SBarry Smith    Input Parameters:
864c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
865af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
866c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8675e76c431SBarry Smith 
868fee21e36SBarry Smith    Collective on SNES
869fee21e36SBarry Smith 
870c4a48953SLois Curfman McInnes    Available Routines:
87117bae607SBarry Smith +  SNESLineSearchCubic() - default line search
87217bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
87317bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
87417bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8755e76c431SBarry Smith 
876c4a48953SLois Curfman McInnes     Options Database Keys:
8774b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8784b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8794b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8804b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8813304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8823304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
883c4a48953SLois Curfman McInnes 
8845e76c431SBarry Smith    Calling sequence of func:
885bd208895SLois Curfman McInnes .vb
886af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
887a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
888bd208895SLois Curfman McInnes .ve
8895e76c431SBarry Smith 
8905e76c431SBarry Smith     Input parameters for func:
891c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
892af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8935e76c431SBarry Smith .   x - current iterate
8945e76c431SBarry Smith .   f - residual evaluated at x
8953c632250SBarry Smith .   y - search direction
8965e76c431SBarry Smith .   w - work vector
897c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8985e76c431SBarry Smith 
8995e76c431SBarry Smith     Output parameters for func:
900c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9013c632250SBarry Smith .   w - new iterate
9025e76c431SBarry Smith .   gnorm - 2-norm of g
9035e76c431SBarry Smith .   ynorm - 2-norm of search length
904a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
905f59ffdedSLois Curfman McInnes 
90636851e7fSLois Curfman McInnes     Level: advanced
90736851e7fSLois Curfman McInnes 
908f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
909f59ffdedSLois Curfman McInnes 
91017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
91117bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9125e76c431SBarry Smith @*/
91317bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9145e76c431SBarry Smith {
915a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
91682bf6240SBarry Smith 
9173a40ed3dSBarry Smith   PetscFunctionBegin;
91817bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
91982bf6240SBarry Smith   if (f) {
920af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92182bf6240SBarry Smith   }
9223a40ed3dSBarry Smith   PetscFunctionReturn(0);
9235e76c431SBarry Smith }
9248e019c35SBarry Smith 
925a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
92604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
927fb2e594dSBarry Smith EXTERN_C_BEGIN
9284a2ae208SSatish Balay #undef __FUNCT__
92917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
93017bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
93182bf6240SBarry Smith {
93282bf6240SBarry Smith   PetscFunctionBegin;
9334b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9344b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
93582bf6240SBarry Smith   PetscFunctionReturn(0);
93682bf6240SBarry Smith }
937fb2e594dSBarry Smith EXTERN_C_END
93804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9394a2ae208SSatish Balay #undef __FUNCT__
9403c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
941c8dd0c0dSLois Curfman McInnes /*@C
9423c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9434b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
944c8dd0c0dSLois Curfman McInnes 
945c8dd0c0dSLois Curfman McInnes    Input Parameters:
946c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9473c632250SBarry Smith .  func - pointer to function
948c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
949c8dd0c0dSLois Curfman McInnes 
950c8dd0c0dSLois Curfman McInnes    Collective on SNES
951c8dd0c0dSLois Curfman McInnes 
952c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
953c8dd0c0dSLois Curfman McInnes .vb
9543c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
955c8dd0c0dSLois Curfman McInnes .ve
956b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
957b1ae27eaSLois Curfman McInnes    on failure.
958c8dd0c0dSLois Curfman McInnes 
959c8dd0c0dSLois Curfman McInnes    Input parameters for func:
960c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
961c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9623c632250SBarry Smith .  x - previous iterate
9633c632250SBarry Smith .  y - new search direction and length
9643c632250SBarry Smith -  w - current candidate iterate
965c8dd0c0dSLois Curfman McInnes 
966c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9673c632250SBarry Smith +  y - search direction (possibly changed)
9683c632250SBarry Smith .  w - current iterate (possibly modified)
9693c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9703c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
971c8dd0c0dSLois Curfman McInnes 
972c8dd0c0dSLois Curfman McInnes    Level: advanced
973c8dd0c0dSLois Curfman McInnes 
9749e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9759e247f21SBarry Smith 
9763c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9773c632250SBarry Smith 
9783c632250SBarry Smith    On input w = x + y
9793c632250SBarry Smith 
98017bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
981b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
982ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9838f99978dSLois Curfman McInnes 
98417bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
985ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
986ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
987ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9889e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
989b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9908f99978dSLois Curfman McInnes 
991c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
992c8dd0c0dSLois Curfman McInnes 
99317bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
994c8dd0c0dSLois Curfman McInnes @*/
9953c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
996c8dd0c0dSLois Curfman McInnes {
9973c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
998c8dd0c0dSLois Curfman McInnes 
999c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10003c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1001c8dd0c0dSLois Curfman McInnes   if (f) {
1002c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1003c8dd0c0dSLois Curfman McInnes   }
1004c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1005c8dd0c0dSLois Curfman McInnes }
10069c3ca13aSBarry Smith #undef __FUNCT__
10079c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10089c3ca13aSBarry Smith /*@C
10099c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10107e4bb74cSBarry Smith          before the line search is called.
10119c3ca13aSBarry Smith 
10129c3ca13aSBarry Smith    Input Parameters:
10139c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10149c3ca13aSBarry Smith .  func - pointer to function
10159c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10169c3ca13aSBarry Smith 
10179c3ca13aSBarry Smith    Collective on SNES
10189c3ca13aSBarry Smith 
10199c3ca13aSBarry Smith    Calling sequence of func:
10209c3ca13aSBarry Smith .vb
10211e3347e8SBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
10229c3ca13aSBarry Smith .ve
10239c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10249c3ca13aSBarry Smith    on failure.
10259c3ca13aSBarry Smith 
10269c3ca13aSBarry Smith    Input parameters for func:
10279c3ca13aSBarry Smith +  snes - nonlinear context
10289c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10299c3ca13aSBarry Smith .  x - previous iterate
10309c3ca13aSBarry Smith -  y - new search direction and length
10319c3ca13aSBarry Smith 
10329c3ca13aSBarry Smith    Output parameters for func:
10339c3ca13aSBarry Smith +  y - search direction (possibly changed)
10349c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10359c3ca13aSBarry Smith 
10369c3ca13aSBarry Smith    Level: advanced
10379c3ca13aSBarry Smith 
10389c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10399c3ca13aSBarry Smith 
10409c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10419c3ca13aSBarry Smith 
10427e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10439c3ca13aSBarry Smith @*/
10449c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10459c3ca13aSBarry Smith {
10469c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10479c3ca13aSBarry Smith 
10489c3ca13aSBarry Smith   PetscFunctionBegin;
10499c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10509c3ca13aSBarry Smith   if (f) {
10519c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10529c3ca13aSBarry Smith   }
10539c3ca13aSBarry Smith   PetscFunctionReturn(0);
10549c3ca13aSBarry Smith }
10559c3ca13aSBarry Smith 
1056c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10579c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10589c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1059c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10604a2ae208SSatish Balay #undef __FUNCT__
10613c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10629c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1063c8dd0c0dSLois Curfman McInnes {
1064c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10653c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10663c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1067c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1068c8dd0c0dSLois Curfman McInnes }
1069c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10709c3ca13aSBarry Smith 
10719c3ca13aSBarry Smith EXTERN_C_BEGIN
10729c3ca13aSBarry Smith #undef __FUNCT__
10739c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10749c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10759c3ca13aSBarry Smith {
10769c3ca13aSBarry Smith   PetscFunctionBegin;
10779c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10789c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10799c3ca13aSBarry Smith   PetscFunctionReturn(0);
10809c3ca13aSBarry Smith }
10819c3ca13aSBarry Smith EXTERN_C_END
1082329e7e40SMatthew Knepley 
1083329e7e40SMatthew Knepley /*
10844b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
108504d965bbSLois Curfman McInnes 
108604d965bbSLois Curfman McInnes    Input Parameters:
108704d965bbSLois Curfman McInnes .  SNES - the SNES context
108804d965bbSLois Curfman McInnes .  viewer - visualization context
108904d965bbSLois Curfman McInnes 
109004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
109104d965bbSLois Curfman McInnes */
10924a2ae208SSatish Balay #undef __FUNCT__
10934b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10946849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1095a935fc98SLois Curfman McInnes {
10964b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
10972fc52814SBarry Smith   const char     *cstr;
1098dfbe8321SBarry Smith   PetscErrorCode ierr;
109932077d6dSBarry Smith   PetscTruth     iascii;
1100a935fc98SLois Curfman McInnes 
11013a40ed3dSBarry Smith   PetscFunctionBegin;
110232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
110332077d6dSBarry Smith   if (iascii) {
110417bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
110517bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
110617bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
110719bcc07fSBarry Smith     else                                                cstr = "unknown";
1108b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1109a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11105cd90555SBarry Smith   } else {
111179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
111219bcc07fSBarry Smith   }
11133a40ed3dSBarry Smith   PetscFunctionReturn(0);
1114a935fc98SLois Curfman McInnes }
111504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
111604d965bbSLois Curfman McInnes /*
11174b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
111804d965bbSLois Curfman McInnes 
111904d965bbSLois Curfman McInnes    Input Parameter:
112004d965bbSLois Curfman McInnes .  snes - the SNES context
112104d965bbSLois Curfman McInnes 
112204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
112304d965bbSLois Curfman McInnes */
11244a2ae208SSatish Balay #undef __FUNCT__
11254b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11266849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11275e42470aSBarry Smith {
11284b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1129e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1130dfbe8321SBarry Smith   PetscErrorCode ierr;
1131a7cc72afSBarry Smith   PetscInt       indx;
1132f1af5d2fSBarry Smith   PetscTruth     flg;
11335e42470aSBarry Smith 
11343a40ed3dSBarry Smith   PetscFunctionBegin;
1135b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11364b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11374b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11384b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1139186905e3SBarry Smith 
114017bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
114125cdf11fSBarry Smith     if (flg) {
114222e36657SBarry Smith       switch (indx) {
1143b49fd9e1SBarry Smith       case 0:
114417bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1145b49fd9e1SBarry Smith         break;
1146b49fd9e1SBarry Smith       case 1:
114717bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1148b49fd9e1SBarry Smith         break;
1149b49fd9e1SBarry Smith       case 2:
115017bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1151b49fd9e1SBarry Smith         break;
1152b49fd9e1SBarry Smith       case 3:
115317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1154b49fd9e1SBarry Smith         break;
11555e42470aSBarry Smith       }
11565e42470aSBarry Smith     }
1157b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11583a40ed3dSBarry Smith   PetscFunctionReturn(0);
11595e42470aSBarry Smith }
116004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1161ebe3b25bSBarry Smith /*MC
1162ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
116304d965bbSLois Curfman McInnes 
1164ebe3b25bSBarry Smith    Options Database:
1165ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1166ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1167ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1168ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1169ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1170ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
117104d965bbSLois Curfman McInnes 
1172ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
117304d965bbSLois Curfman McInnes 
1174ee3001cbSBarry Smith    Level: beginner
1175ee3001cbSBarry Smith 
117617bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
117717bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
117817bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1179ebe3b25bSBarry Smith 
1180ebe3b25bSBarry Smith M*/
1181fb2e594dSBarry Smith EXTERN_C_BEGIN
11824a2ae208SSatish Balay #undef __FUNCT__
11834b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
118463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
11855e42470aSBarry Smith {
1186dfbe8321SBarry Smith   PetscErrorCode ierr;
11874b27c08aSLois Curfman McInnes   SNES_LS        *neP;
11885e42470aSBarry Smith 
11893a40ed3dSBarry Smith   PetscFunctionBegin;
1190e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1191e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1192e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1193e7788613SBarry Smith   snes->ops->converged	     = SNESConverged_LS;
1194e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1195e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
11965baf8537SBarry Smith   snes->nwork                = 0;
11975e42470aSBarry Smith 
11984b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
119952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
12005e42470aSBarry Smith   snes->data    	= (void*)neP;
12015e42470aSBarry Smith   neP->alpha		= 1.e-4;
12025e42470aSBarry Smith   neP->maxstep		= 1.e8;
12035e42470aSBarry Smith   neP->steptol		= 1.e-12;
120417bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1205c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12063c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12073c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12083c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12093c632250SBarry Smith   neP->precheck         = PETSC_NULL;
121082bf6240SBarry Smith 
1211557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1212557d3f75SLisandro Dalcin 					   "SNESLineSearchSet_LS",
1213557d3f75SLisandro Dalcin 					   SNESLineSearchSet_LS);CHKERRQ(ierr);
1214557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1215557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPostCheck_LS",
1216557d3f75SLisandro Dalcin 					   SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
1217557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1218557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPreCheck_LS",
1219557d3f75SLisandro Dalcin 					   SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
122082bf6240SBarry Smith 
12213a40ed3dSBarry Smith   PetscFunctionReturn(0);
12225e42470aSBarry Smith }
1223fb2e594dSBarry Smith EXTERN_C_END
122482bf6240SBarry Smith 
122582bf6240SBarry Smith 
122682bf6240SBarry Smith 
1227