xref: /petsc/src/snes/impls/ls/ls.c (revision e106eecf280d2eadc5a817a2dc6d6d2b72b33833)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
25e76c431SBarry Smith 
37c4f633dSBarry Smith #include "../src/snes/impls/ls/ls.h"
45e42470aSBarry Smith 
58a5d9ee4SBarry Smith /*
620f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
720f52e01SBarry Smith     || F(u) ||_2 but not a zero, F(u) = 0. 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
920f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
108a5d9ee4SBarry Smith */
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
131e2582c4SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,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);
261e2582c4SBarry Smith     ierr = PetscInfo1(snes,"|| 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);
401e2582c4SBarry Smith     ierr = PetscInfo1(snes,"(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"
511e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,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) {
681e2582c4SBarry Smith       ierr = PetscInfo1(snes,"||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;
13985385478SLisandro Dalcin   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
14085385478SLisandro Dalcin   Vec                Y,X,F,G,W;
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;
15785385478SLisandro Dalcin   snes->norm = 0;
1584c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15919717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1604936397dSBarry Smith   if (snes->domainerror) {
1614936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1624936397dSBarry Smith     PetscFunctionReturn(0);
1634936397dSBarry Smith   }
1642613ca53SBarry Smith   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1652613ca53SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
1662613ca53SBarry Smith   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
1672613ca53SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
168f66fdb6dSSatish Balay   if PetscIsInfOrNanReal(fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1690f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1705e42470aSBarry Smith   snes->norm = fnorm;
1710f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17284c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
17394a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1743f149594SLisandro Dalcin 
1753f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
1763f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
17785385478SLisandro Dalcin   /* test convergence */
17885385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
17906ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
180d034289bSBarry Smith 
1815e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1825e76c431SBarry Smith 
183329e7e40SMatthew Knepley     /* Call general purpose update function */
184e7788613SBarry Smith     if (snes->ops->update) {
185e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
186de8cb200SMatthew Knepley     }
187329e7e40SMatthew Knepley 
188ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1895334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
19094b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
19171f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
1923d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1933d4c4710SBarry Smith     if (kspreason < 0) {
1943d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
195ef998cc9SBarry Smith         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
1963d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
197ab81109fSSatish Balay         break;
1983d4c4710SBarry Smith       }
1993d4c4710SBarry Smith     }
2003d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2013f149594SLisandro Dalcin     snes->linear_its += lits;
2023f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
20374637425SBarry Smith 
2049c3ca13aSBarry Smith     if (neP->precheckstep) {
2059c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
2069c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
2079c3ca13aSBarry Smith     }
2089c3ca13aSBarry Smith 
209b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
2101e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
21185471664SBarry Smith     }
21274637425SBarry Smith 
213ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
214ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
215e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
216ea4d3ed3SLois Curfman McInnes     */
21785385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
2183f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
219dc357ad2SBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
220ae15b995SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
2214a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2224936397dSBarry Smith     if (snes->domainerror) {
2234936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2244936397dSBarry Smith       PetscFunctionReturn(0);
2254936397dSBarry Smith     }
226a7cc72afSBarry Smith     if (!lssucceed) {
22750ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
22885385478SLisandro Dalcin 	PetscTruth ismin;
2293505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2301e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
2318a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2323505fcc1SBarry Smith         break;
2333505fcc1SBarry Smith       }
23450ffb88aSMatthew Knepley     }
23585385478SLisandro Dalcin     /* Update function and solution vectors */
23685385478SLisandro Dalcin     fnorm = gnorm;
23785385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
23885385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
23985385478SLisandro Dalcin     /* Monitor convergence */
24085385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
24185385478SLisandro Dalcin     snes->iter = i+1;
24285385478SLisandro Dalcin     snes->norm = fnorm;
24385385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
24485385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
24585385478SLisandro Dalcin     SNESMonitor(snes,snes->iter,snes->norm);
24685385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
24785385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
248e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2493f149594SLisandro Dalcin     if (snes->reason) break;
25016c95cacSBarry Smith   }
25152392280SLois Curfman McInnes   if (i == maxits) {
252ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
25385385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
25452392280SLois Curfman McInnes   }
2553a40ed3dSBarry Smith   PetscFunctionReturn(0);
2565e76c431SBarry Smith }
25704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25804d965bbSLois Curfman McInnes /*
2594b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2604b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
26104d965bbSLois Curfman McInnes 
26204d965bbSLois Curfman McInnes    Input Parameter:
26304d965bbSLois Curfman McInnes .  snes - the SNES context
26404d965bbSLois Curfman McInnes .  x - the solution vector
26504d965bbSLois Curfman McInnes 
26604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26704d965bbSLois Curfman McInnes 
26804d965bbSLois Curfman McInnes    Notes:
2694b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
27004d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
27104d965bbSLois Curfman McInnes    the call to SNESSolve().
27204d965bbSLois Curfman McInnes  */
2734a2ae208SSatish Balay #undef __FUNCT__
2744b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
275dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2765e76c431SBarry Smith {
277dfbe8321SBarry Smith   PetscErrorCode ierr;
2783a40ed3dSBarry Smith 
2793a40ed3dSBarry Smith   PetscFunctionBegin;
28085385478SLisandro Dalcin   if (!snes->vec_sol_update) {
28185385478SLisandro Dalcin     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
28285385478SLisandro Dalcin     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
28385385478SLisandro Dalcin   }
284e74804ceSBarry Smith   if (!snes->work) {
28585385478SLisandro Dalcin     snes->nwork = 3;
28685385478SLisandro Dalcin     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
287efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
288e74804ceSBarry Smith   }
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
2905e76c431SBarry Smith }
29104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
29204d965bbSLois Curfman McInnes /*
2934b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2944b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29504d965bbSLois Curfman McInnes 
29604d965bbSLois Curfman McInnes    Input Parameter:
29704d965bbSLois Curfman McInnes .  snes - the SNES context
29804d965bbSLois Curfman McInnes 
299de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30004d965bbSLois Curfman McInnes  */
3014a2ae208SSatish Balay #undef __FUNCT__
3024b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
303dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3045e76c431SBarry Smith {
305dfbe8321SBarry Smith   PetscErrorCode ierr;
3063a40ed3dSBarry Smith 
3073a40ed3dSBarry Smith   PetscFunctionBegin;
30885385478SLisandro Dalcin   if (snes->vec_sol_update) {
30985385478SLisandro Dalcin     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
31085385478SLisandro Dalcin     snes->vec_sol_update = PETSC_NULL;
31185385478SLisandro Dalcin   }
3125baf8537SBarry Smith   if (snes->nwork) {
3134b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
31485385478SLisandro Dalcin     snes->nwork = 0;
31585385478SLisandro Dalcin     snes->work  = PETSC_NULL;
31621c89e3eSBarry Smith   }
3175c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
318557d3f75SLisandro Dalcin 
319557d3f75SLisandro Dalcin   /* clear composed functions */
320557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
321557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
322557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
323557d3f75SLisandro Dalcin 
3243a40ed3dSBarry Smith   PetscFunctionReturn(0);
3255e76c431SBarry Smith }
32604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3274a2ae208SSatish Balay #undef __FUNCT__
32817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
32982bf6240SBarry Smith 
3304b828684SBarry Smith /*@C
33117bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3325e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3335e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3345e76c431SBarry Smith 
335c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
336c7afd0dbSLois Curfman McInnes 
3375e76c431SBarry Smith    Input Parameters:
338c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
339af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3405e76c431SBarry Smith .  x - current iterate
3415e76c431SBarry Smith .  f - residual evaluated at x
3423c632250SBarry Smith .  y - search direction
3435e76c431SBarry Smith .  w - work vector
344dc357ad2SBarry Smith .  fnorm - 2-norm of f
345dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
3465e76c431SBarry Smith 
347c4a48953SLois Curfman McInnes    Output Parameters:
348c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3493c632250SBarry Smith .  w - new iterate
3505e76c431SBarry Smith .  gnorm - 2-norm of g
3515e76c431SBarry Smith .  ynorm - 2-norm of search length
352a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
353fee21e36SBarry Smith 
354c4a48953SLois Curfman McInnes    Options Database Key:
35517bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
356c4a48953SLois Curfman McInnes 
35736851e7fSLois Curfman McInnes    Level: advanced
35836851e7fSLois Curfman McInnes 
35928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
36028ae5a21SLois Curfman McInnes 
36117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
36217bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3635e76c431SBarry Smith @*/
364dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
3655e76c431SBarry Smith {
366dfbe8321SBarry Smith   PetscErrorCode ierr;
3674b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3683c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3695334005bSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
371a7cc72afSBarry Smith   *flag = PETSC_TRUE;
372d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
373a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
37479f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3753c632250SBarry Smith   if (neP->postcheckstep) {
3763c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3778f99978dSLois Curfman McInnes   }
3783c632250SBarry Smith   if (changed_y) {
37979f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3803c632250SBarry Smith   }
3814936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
3824936397dSBarry Smith   if (!snes->domainerror) {
383a3b38805SLois Curfman McInnes     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
384f66fdb6dSSatish Balay     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3854936397dSBarry Smith   }
386d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
3885e76c431SBarry Smith }
38904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3904a2ae208SSatish Balay #undef __FUNCT__
39117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
39282bf6240SBarry Smith 
39329e0b56fSBarry Smith /*@C
39417bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
39529e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
39629e0b56fSBarry Smith    even compute the norm of the function or search direction; this
39729e0b56fSBarry Smith    is intended only when you know the full step is fine and are
39829e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
399c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
400c7afd0dbSLois Curfman McInnes 
401c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
40229e0b56fSBarry Smith 
40329e0b56fSBarry Smith    Input Parameters:
404c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
405af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
40629e0b56fSBarry Smith .  x - current iterate
40729e0b56fSBarry Smith .  f - residual evaluated at x
4083c632250SBarry Smith .  y - search direction
40929e0b56fSBarry Smith .  w - work vector
410dc357ad2SBarry Smith .  fnorm - 2-norm of f
411dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
41229e0b56fSBarry Smith 
41329e0b56fSBarry Smith    Output Parameters:
414c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4153c632250SBarry Smith .  w - new iterate
41629e0b56fSBarry Smith .  gnorm - not changed
41729e0b56fSBarry Smith .  ynorm - not changed
418a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
419fee21e36SBarry Smith 
42029e0b56fSBarry Smith    Options Database Key:
42117bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
42229e0b56fSBarry Smith 
4238cbba510SBarry Smith    Notes:
42417bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
425ea56f5baSLois Curfman McInnes    either the options
426ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
427ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
428329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
429329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
430329f5518SBarry Smith 
431329f5518SBarry Smith    During the final iteration this will not evaluate the function at
432329f5518SBarry Smith    the solution point. This is to save a function evaluation while
433329f5518SBarry Smith    using pseudo-timestepping.
4348cbba510SBarry Smith 
435ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
436a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
437ea56f5baSLois Curfman McInnes    correct, since they are not computed.
438ea56f5baSLois Curfman McInnes 
439ea56f5baSLois Curfman McInnes    Level: advanced
4408cbba510SBarry Smith 
44129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
44229e0b56fSBarry Smith 
44317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
44417bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
44529e0b56fSBarry Smith @*/
446dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
44729e0b56fSBarry Smith {
448dfbe8321SBarry Smith   PetscErrorCode ierr;
4494b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4503c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
45129e0b56fSBarry Smith 
4523a40ed3dSBarry Smith   PetscFunctionBegin;
453a7cc72afSBarry Smith   *flag = PETSC_TRUE;
454d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
45579f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4563c632250SBarry Smith   if (neP->postcheckstep) {
4573c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4583c632250SBarry Smith   }
4593c632250SBarry Smith   if (changed_y) {
46079f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4618f99978dSLois Curfman McInnes   }
462329f5518SBarry Smith 
463329f5518SBarry Smith   /* don't evaluate function the last time through */
464329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4654936397dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
466329f5518SBarry Smith   }
467d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
46929e0b56fSBarry Smith }
47004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4714a2ae208SSatish Balay #undef __FUNCT__
47217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4734b828684SBarry Smith /*@C
47417bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4755e76c431SBarry Smith 
476c7afd0dbSLois Curfman McInnes    Collective on SNES
477c7afd0dbSLois Curfman McInnes 
4785e76c431SBarry Smith    Input Parameters:
479c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
480af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4815e76c431SBarry Smith .  x - current iterate
4825e76c431SBarry Smith .  f - residual evaluated at x
4833c632250SBarry Smith .  y - search direction
4845e76c431SBarry Smith .  w - work vector
485dc357ad2SBarry Smith .  fnorm - 2-norm of f
486dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
4875e76c431SBarry Smith 
488393d2d9aSLois Curfman McInnes    Output Parameters:
489c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4903c632250SBarry Smith .  w - new iterate
4915e76c431SBarry Smith .  gnorm - 2-norm of g
4925e76c431SBarry Smith .  ynorm - 2-norm of search length
493a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
494fee21e36SBarry Smith 
495c4a48953SLois Curfman McInnes    Options Database Key:
496*e106eecfSBarry Smith +  -snes_ls cubic - Activates SNESLineSearchCubic()
497*e106eecfSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
498*e106eecfSBarry Smith .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
499*e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
500*e106eecfSBarry Smith 
501c4a48953SLois Curfman McInnes 
5025e76c431SBarry Smith    Notes:
5035e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
5045e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
5055e76c431SBarry Smith 
50636851e7fSLois Curfman McInnes    Level: advanced
50736851e7fSLois Curfman McInnes 
50828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
50928ae5a21SLois Curfman McInnes 
51017bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5115e76c431SBarry Smith @*/
512dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
5135e76c431SBarry Smith {
514406556e6SLois Curfman McInnes   /*
515406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
516406556e6SLois Curfman McInnes      minimization problem:
517406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
518406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
519406556e6SLois Curfman McInnes    */
520406556e6SLois Curfman McInnes 
521*e106eecfSBarry Smith   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
522*e106eecfSBarry Smith   PetscReal      minlambda,lambda,lambdatemp;
523aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
524f5b6597dSBarry Smith   PetscScalar    cinitslope;
5256b5873e3SBarry Smith #endif
526dfbe8321SBarry Smith   PetscErrorCode ierr;
527a7cc72afSBarry Smith   PetscInt       count;
5284b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
5293c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5305e76c431SBarry Smith 
5313a40ed3dSBarry Smith   PetscFunctionBegin;
532d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
533a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5345e76c431SBarry Smith 
535cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5369e247f21SBarry Smith   if (!*ynorm) {
537ae15b995SBarry Smith     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
538a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5393c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
540a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
541a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
542ad922ac9SBarry Smith     goto theend1;
54394a9d846SBarry Smith   }
544*e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
545*e106eecfSBarry Smith     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(neP->maxstep/(*ynorm)),*ynorm);CHKERRQ(ierr);
546*e106eecfSBarry Smith     ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
547*e106eecfSBarry Smith     *ynorm = neP->maxstep;
5485e76c431SBarry Smith   }
549ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
550*e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
551a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
552aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
553a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
554329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5555e42470aSBarry Smith #else
556a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5575e42470aSBarry Smith #endif
5585e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5595e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5605e76c431SBarry Smith 
561e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
56243e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
563ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
56443e71028SBarry Smith     *flag = PETSC_FALSE;
56543e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
56643e71028SBarry Smith     goto theend1;
56743e71028SBarry Smith   }
5684936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5694936397dSBarry Smith   if (snes->domainerror) {
5704936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5714936397dSBarry Smith     PetscFunctionReturn(0);
57219717074SBarry Smith   }
573313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
574f66fdb6dSSatish Balay   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
575f5b6597dSBarry Smith   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
576*e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
577f5b6597dSBarry Smith     ierr = PetscInfo2(snes,"Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
578ad922ac9SBarry Smith     goto theend1;
5795e76c431SBarry Smith   }
5805e76c431SBarry Smith 
5815e76c431SBarry Smith   /* Fit points with quadratic */
582313b5538SBarry Smith   lambda     = 1.0;
583a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5845e76c431SBarry Smith   lambdaprev = lambda;
5855e76c431SBarry Smith   gnormprev  = *gnorm;
586329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
587ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
588ddd12767SBarry Smith   else                         lambda = lambdatemp;
589275d25c3SBarry Smith 
590e68848bdSBarry Smith   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
59143e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
592ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
59343e71028SBarry Smith     *flag = PETSC_FALSE;
59443e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
59543e71028SBarry Smith     goto theend1;
59643e71028SBarry Smith   }
597f5b6597dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5984936397dSBarry Smith   if (snes->domainerror) {
5994936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6004936397dSBarry Smith     PetscFunctionReturn(0);
60119717074SBarry Smith   }
602cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
603f66fdb6dSSatish Balay   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
604f5b6597dSBarry Smith   ierr = PetscInfo1(snes,"gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
605*e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
606ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
607ad922ac9SBarry Smith     goto theend1;
6085e76c431SBarry Smith   }
6095e76c431SBarry Smith 
6105e76c431SBarry Smith   /* Fit points with cubic */
6115e76c431SBarry Smith   count = 1;
612bb9195b6SBarry Smith   while (PETSC_TRUE) {
613dc357ad2SBarry Smith     if (lambda <= minlambda) {
614dc357ad2SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
615*e106eecfSBarry Smith       ierr = PetscInfo6(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,minlambda,lambda,initslope);CHKERRQ(ierr);
61643e71028SBarry Smith       *flag = PETSC_FALSE;
61743e71028SBarry Smith       break;
6185e76c431SBarry Smith     }
6195d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6205d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
621ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6222b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6235e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6245e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6255e76c431SBarry Smith     if (a == 0.0) {
6265e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6275e76c431SBarry Smith     } else {
6285e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6295e76c431SBarry Smith     }
6305e76c431SBarry Smith     lambdaprev = lambda;
6315e76c431SBarry Smith     gnormprev  = *gnorm;
632329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
633329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6345e76c431SBarry Smith     else                         lambda     = lambdatemp;
635e68848bdSBarry Smith     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
63643e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
637ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
638ae15b995SBarry 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);
639ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
64043e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
641ed98166eSMatthew Knepley       break;
642ed98166eSMatthew Knepley     }
643f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6444936397dSBarry Smith     if (snes->domainerror) {
6454936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6464936397dSBarry Smith       PetscFunctionReturn(0);
64719717074SBarry Smith     }
648cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
649f66fdb6dSSatish Balay     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
650*e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */
651f5b6597dSBarry Smith       ierr = PetscInfo2(snes,"Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
65243e71028SBarry Smith       break;
6532b022350SLois Curfman McInnes     } else {
654f5b6597dSBarry Smith       ierr = PetscInfo2(snes,"Cubic step no good, shrinking lambda, current gnorem %G lambda=%18.16e\n",*gnorm,lambda);CHKERRQ(ierr);
6555e76c431SBarry Smith     }
6565e76c431SBarry Smith     count++;
6575e76c431SBarry Smith   }
658ad922ac9SBarry Smith   theend1:
6598f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6603c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6613c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6623c632250SBarry Smith     if (changed_y) {
66379f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6643c632250SBarry Smith     }
6653c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
666f5b6597dSBarry Smith       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6674936397dSBarry Smith       if (snes->domainerror) {
6684936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6694936397dSBarry Smith         PetscFunctionReturn(0);
67019717074SBarry Smith       }
6718f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
672f66fdb6dSSatish Balay       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
673a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6748f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
675a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6768f99978dSLois Curfman McInnes     }
6778f99978dSLois Curfman McInnes   }
678d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6793a40ed3dSBarry Smith   PetscFunctionReturn(0);
6805e76c431SBarry Smith }
68104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6824a2ae208SSatish Balay #undef __FUNCT__
68317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6844b828684SBarry Smith /*@C
68517bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6865e76c431SBarry Smith 
687c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
688c7afd0dbSLois Curfman McInnes 
6895e42470aSBarry Smith    Input Parameters:
690c7afd0dbSLois Curfman McInnes +  snes - the SNES context
691af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6925e42470aSBarry Smith .  x - current iterate
6935e42470aSBarry Smith .  f - residual evaluated at x
6943c632250SBarry Smith .  y - search direction
6955e42470aSBarry Smith .  w - work vector
696dc357ad2SBarry Smith .  fnorm - 2-norm of f
697dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
6985e42470aSBarry Smith 
699c4a48953SLois Curfman McInnes    Output Parameters:
7007f3332b4SBarry Smith +  g - residual evaluated at new iterate w
701*e106eecfSBarry Smith .  w - new iterate (x + lambda*y)
7025e42470aSBarry Smith .  gnorm - 2-norm of g
7035e42470aSBarry Smith .  ynorm - 2-norm of search length
704a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
705fee21e36SBarry Smith 
706ce9499c7SBarry Smith    Options Database Keys:
707ce9499c7SBarry Smith +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
708ce9499c7SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
709*e106eecfSBarry Smith .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
710*e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
711*e106eecfSBarry Smith 
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 @*/
721dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,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    */
729*e106eecfSBarry Smith   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
730aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
731f5b6597dSBarry Smith   PetscScalar    cinitslope;
7326b5873e3SBarry Smith #endif
733dfbe8321SBarry Smith   PetscErrorCode ierr;
734a7cc72afSBarry Smith   PetscInt       count;
7354b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
7363c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7375e76c431SBarry Smith 
7383a40ed3dSBarry Smith   PetscFunctionBegin;
739d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
740a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7415e76c431SBarry Smith 
7423505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74363b9fa5eSBarry Smith   if (*ynorm == 0.0) {
744ae15b995SBarry Smith     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
745b37302c6SLois Curfman McInnes     *gnorm = fnorm;
746e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
747b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
748a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
749ad922ac9SBarry Smith     goto theend2;
75094a9d846SBarry Smith   }
751*e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
752*e106eecfSBarry Smith     ierr   = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
753*e106eecfSBarry Smith     *ynorm = neP->maxstep;
7545e76c431SBarry Smith   }
755ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
756*e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
757a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
758aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
759a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
760329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7615e42470aSBarry Smith #else
762a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7635e42470aSBarry Smith #endif
7645e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7655e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
766f8833479SBarry Smith   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
7675e42470aSBarry Smith 
768e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
76943e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
770ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
77143e71028SBarry Smith     *flag = PETSC_FALSE;
77243e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77343e71028SBarry Smith     goto theend2;
77443e71028SBarry Smith   }
7754936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
7764936397dSBarry Smith   if (snes->domainerror) {
7774936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7784936397dSBarry Smith     PetscFunctionReturn(0);
77919717074SBarry Smith   }
780cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
781f66fdb6dSSatish Balay   if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
782*e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
783ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
784ad922ac9SBarry Smith     goto theend2;
7855e42470aSBarry Smith   }
7865e42470aSBarry Smith 
7875e42470aSBarry Smith   /* Fit points with quadratic */
788313b5538SBarry Smith   lambda = 1.0;
7895e42470aSBarry Smith   count = 1;
7905ca10a37SBarry Smith   while (PETSC_TRUE) {
7915e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
792ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
793ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
794e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
79543e71028SBarry Smith       *flag = PETSC_FALSE;
79643e71028SBarry Smith       break;
7975e42470aSBarry Smith     }
798a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
799329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
800329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
801329f5518SBarry Smith     else                         lambda     = lambdatemp;
802275d25c3SBarry Smith 
803e68848bdSBarry Smith     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
80443e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
805ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
806ae15b995SBarry 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);
807ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
80843e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
809ed98166eSMatthew Knepley       break;
810ed98166eSMatthew Knepley     }
811f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8124936397dSBarry Smith     if (snes->domainerror) {
8134936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8144936397dSBarry Smith       PetscFunctionReturn(0);
81519717074SBarry Smith     }
816cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
817f66fdb6dSSatish Balay     if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
818*e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
819ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
82006259719SBarry Smith       break;
8215e42470aSBarry Smith     }
8225e42470aSBarry Smith     count++;
8235e42470aSBarry Smith   }
824ad922ac9SBarry Smith   theend2:
8258f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8263c632250SBarry Smith   if (neP->postcheckstep) {
8273c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8283c632250SBarry Smith     if (changed_y) {
82979f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8303c632250SBarry Smith     }
8313c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8323c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
8334936397dSBarry Smith       if (snes->domainerror) {
8344936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8354936397dSBarry Smith         PetscFunctionReturn(0);
83619717074SBarry Smith       }
8378f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
838a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
8398f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
840a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
841f66fdb6dSSatish Balay       if PetscIsInfOrNanReal(*gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8428f99978dSLois Curfman McInnes     }
8438f99978dSLois Curfman McInnes   }
844d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
8465e76c431SBarry Smith }
8472343ba6eSBarry Smith 
84804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8494a2ae208SSatish Balay #undef __FUNCT__
85017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
851c9e6a524SLois Curfman McInnes /*@C
85217bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8534b27c08aSLois Curfman McInnes    by the method SNESLS.
8545e76c431SBarry Smith 
8555e76c431SBarry Smith    Input Parameters:
856c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
857af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
858c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8595e76c431SBarry Smith 
860fee21e36SBarry Smith    Collective on SNES
861fee21e36SBarry Smith 
862c4a48953SLois Curfman McInnes    Available Routines:
86317bae607SBarry Smith +  SNESLineSearchCubic() - default line search
86417bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
86517bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
86617bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8675e76c431SBarry Smith 
868c4a48953SLois Curfman McInnes     Options Database Keys:
8694b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8704b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
871*e106eecfSBarry Smith .   -snes_ls_maxstep <maxstep> - Sets maximum step the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
872*e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
873c4a48953SLois Curfman McInnes 
8745e76c431SBarry Smith    Calling sequence of func:
875bd208895SLois Curfman McInnes .vb
876dc357ad2SBarry Smith    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
877bd208895SLois Curfman McInnes .ve
8785e76c431SBarry Smith 
8795e76c431SBarry Smith     Input parameters for func:
880c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
881af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8825e76c431SBarry Smith .   x - current iterate
8835e76c431SBarry Smith .   f - residual evaluated at x
8843c632250SBarry Smith .   y - search direction
8855e76c431SBarry Smith .   w - work vector
886c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8875e76c431SBarry Smith 
8885e76c431SBarry Smith     Output parameters for func:
889c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8903c632250SBarry Smith .   w - new iterate
8915e76c431SBarry Smith .   gnorm - 2-norm of g
8925e76c431SBarry Smith .   ynorm - 2-norm of search length
893a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
894f59ffdedSLois Curfman McInnes 
89536851e7fSLois Curfman McInnes     Level: advanced
89636851e7fSLois Curfman McInnes 
897f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
898f59ffdedSLois Curfman McInnes 
89917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
90017bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9015e76c431SBarry Smith @*/
902dc357ad2SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9035e76c431SBarry Smith {
904dc357ad2SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
90582bf6240SBarry Smith 
9063a40ed3dSBarry Smith   PetscFunctionBegin;
90717bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
90882bf6240SBarry Smith   if (f) {
909af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
91082bf6240SBarry Smith   }
9113a40ed3dSBarry Smith   PetscFunctionReturn(0);
9125e76c431SBarry Smith }
9138e019c35SBarry Smith 
914dc357ad2SBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
91504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
916fb2e594dSBarry Smith EXTERN_C_BEGIN
9174a2ae208SSatish Balay #undef __FUNCT__
91817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
91917bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
92082bf6240SBarry Smith {
92182bf6240SBarry Smith   PetscFunctionBegin;
9224b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9234b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
92482bf6240SBarry Smith   PetscFunctionReturn(0);
92582bf6240SBarry Smith }
926fb2e594dSBarry Smith EXTERN_C_END
92704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9284a2ae208SSatish Balay #undef __FUNCT__
9293c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
930c8dd0c0dSLois Curfman McInnes /*@C
9313c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9324b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
933c8dd0c0dSLois Curfman McInnes 
934c8dd0c0dSLois Curfman McInnes    Input Parameters:
935c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9363c632250SBarry Smith .  func - pointer to function
937c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
938c8dd0c0dSLois Curfman McInnes 
939c8dd0c0dSLois Curfman McInnes    Collective on SNES
940c8dd0c0dSLois Curfman McInnes 
941c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
942c8dd0c0dSLois Curfman McInnes .vb
9433c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
944c8dd0c0dSLois Curfman McInnes .ve
945b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
946b1ae27eaSLois Curfman McInnes    on failure.
947c8dd0c0dSLois Curfman McInnes 
948c8dd0c0dSLois Curfman McInnes    Input parameters for func:
949c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
950c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9513c632250SBarry Smith .  x - previous iterate
9523c632250SBarry Smith .  y - new search direction and length
9533c632250SBarry Smith -  w - current candidate iterate
954c8dd0c0dSLois Curfman McInnes 
955c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9563c632250SBarry Smith +  y - search direction (possibly changed)
9573c632250SBarry Smith .  w - current iterate (possibly modified)
9583c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9593c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
960c8dd0c0dSLois Curfman McInnes 
961c8dd0c0dSLois Curfman McInnes    Level: advanced
962c8dd0c0dSLois Curfman McInnes 
9639e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9649e247f21SBarry Smith 
9653c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9663c632250SBarry Smith 
9673c632250SBarry Smith    On input w = x + y
9683c632250SBarry Smith 
96917bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
970b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
971ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9728f99978dSLois Curfman McInnes 
97317bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
974ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
975ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
976ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9779e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
978b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9798f99978dSLois Curfman McInnes 
980c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
981c8dd0c0dSLois Curfman McInnes 
98217bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
983c8dd0c0dSLois Curfman McInnes @*/
9843c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
985c8dd0c0dSLois Curfman McInnes {
9863c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
987c8dd0c0dSLois Curfman McInnes 
988c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9893c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
990c8dd0c0dSLois Curfman McInnes   if (f) {
991c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
992c8dd0c0dSLois Curfman McInnes   }
993c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
994c8dd0c0dSLois Curfman McInnes }
9959c3ca13aSBarry Smith #undef __FUNCT__
9969c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
9979c3ca13aSBarry Smith /*@C
9989c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
9997e4bb74cSBarry Smith          before the line search is called.
10009c3ca13aSBarry Smith 
10019c3ca13aSBarry Smith    Input Parameters:
10029c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10039c3ca13aSBarry Smith .  func - pointer to function
10049c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10059c3ca13aSBarry Smith 
10069c3ca13aSBarry Smith    Collective on SNES
10079c3ca13aSBarry Smith 
10089c3ca13aSBarry Smith    Calling sequence of func:
10099c3ca13aSBarry Smith .vb
10101e3347e8SBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
10119c3ca13aSBarry Smith .ve
10129c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10139c3ca13aSBarry Smith    on failure.
10149c3ca13aSBarry Smith 
10159c3ca13aSBarry Smith    Input parameters for func:
10169c3ca13aSBarry Smith +  snes - nonlinear context
10179c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10189c3ca13aSBarry Smith .  x - previous iterate
10199c3ca13aSBarry Smith -  y - new search direction and length
10209c3ca13aSBarry Smith 
10219c3ca13aSBarry Smith    Output parameters for func:
10229c3ca13aSBarry Smith +  y - search direction (possibly changed)
10239c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10249c3ca13aSBarry Smith 
10259c3ca13aSBarry Smith    Level: advanced
10269c3ca13aSBarry Smith 
10279c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10289c3ca13aSBarry Smith 
10299c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10309c3ca13aSBarry Smith 
10317e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10329c3ca13aSBarry Smith @*/
10339c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10349c3ca13aSBarry Smith {
10359c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10369c3ca13aSBarry Smith 
10379c3ca13aSBarry Smith   PetscFunctionBegin;
10389c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10399c3ca13aSBarry Smith   if (f) {
10409c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10419c3ca13aSBarry Smith   }
10429c3ca13aSBarry Smith   PetscFunctionReturn(0);
10439c3ca13aSBarry Smith }
10449c3ca13aSBarry Smith 
1045c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10469c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10479c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1048c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10494a2ae208SSatish Balay #undef __FUNCT__
10503c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10519c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1052c8dd0c0dSLois Curfman McInnes {
1053c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10543c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10553c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1056c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1057c8dd0c0dSLois Curfman McInnes }
1058c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10599c3ca13aSBarry Smith 
10609c3ca13aSBarry Smith EXTERN_C_BEGIN
10619c3ca13aSBarry Smith #undef __FUNCT__
10629c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10639c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10649c3ca13aSBarry Smith {
10659c3ca13aSBarry Smith   PetscFunctionBegin;
10669c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10679c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10689c3ca13aSBarry Smith   PetscFunctionReturn(0);
10699c3ca13aSBarry Smith }
10709c3ca13aSBarry Smith EXTERN_C_END
1071329e7e40SMatthew Knepley 
1072329e7e40SMatthew Knepley /*
10734b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
107404d965bbSLois Curfman McInnes 
107504d965bbSLois Curfman McInnes    Input Parameters:
107604d965bbSLois Curfman McInnes .  SNES - the SNES context
107704d965bbSLois Curfman McInnes .  viewer - visualization context
107804d965bbSLois Curfman McInnes 
107904d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
108004d965bbSLois Curfman McInnes */
10814a2ae208SSatish Balay #undef __FUNCT__
10824b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10836849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1084a935fc98SLois Curfman McInnes {
10854b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
10862fc52814SBarry Smith   const char     *cstr;
1087dfbe8321SBarry Smith   PetscErrorCode ierr;
108832077d6dSBarry Smith   PetscTruth     iascii;
1089a935fc98SLois Curfman McInnes 
10903a40ed3dSBarry Smith   PetscFunctionBegin;
109132077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
109232077d6dSBarry Smith   if (iascii) {
109317bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
109417bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
109517bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
109619bcc07fSBarry Smith     else                                                cstr = "unknown";
1097b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1098*e106eecfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr);
10995cd90555SBarry Smith   } else {
110079a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
110119bcc07fSBarry Smith   }
11023a40ed3dSBarry Smith   PetscFunctionReturn(0);
1103a935fc98SLois Curfman McInnes }
110404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
110504d965bbSLois Curfman McInnes /*
11064b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
110704d965bbSLois Curfman McInnes 
110804d965bbSLois Curfman McInnes    Input Parameter:
110904d965bbSLois Curfman McInnes .  snes - the SNES context
111004d965bbSLois Curfman McInnes 
111104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
111204d965bbSLois Curfman McInnes */
11134a2ae208SSatish Balay #undef __FUNCT__
11144b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11156849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11165e42470aSBarry Smith {
11174b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1118e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1119dfbe8321SBarry Smith   PetscErrorCode ierr;
1120a7cc72afSBarry Smith   PetscInt       indx;
1121f1af5d2fSBarry Smith   PetscTruth     flg;
11225e42470aSBarry Smith 
11233a40ed3dSBarry Smith   PetscFunctionBegin;
1124b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11254b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11264b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1127*e106eecfSBarry Smith     ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr);
1128186905e3SBarry Smith 
112917bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
113025cdf11fSBarry Smith     if (flg) {
113122e36657SBarry Smith       switch (indx) {
1132b49fd9e1SBarry Smith       case 0:
113317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1134b49fd9e1SBarry Smith         break;
1135b49fd9e1SBarry Smith       case 1:
113617bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1137b49fd9e1SBarry Smith         break;
1138b49fd9e1SBarry Smith       case 2:
113917bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1140b49fd9e1SBarry Smith         break;
1141b49fd9e1SBarry Smith       case 3:
114217bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1143b49fd9e1SBarry Smith         break;
11445e42470aSBarry Smith       }
11455e42470aSBarry Smith     }
1146b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11473a40ed3dSBarry Smith   PetscFunctionReturn(0);
11485e42470aSBarry Smith }
114904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1150ebe3b25bSBarry Smith /*MC
1151ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
115204d965bbSLois Curfman McInnes 
1153ebe3b25bSBarry Smith    Options Database:
1154ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1155ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1156*e106eecfSBarry Smith .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
1157*e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
115804d965bbSLois Curfman McInnes 
1159ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
116004d965bbSLois Curfman McInnes 
1161ee3001cbSBarry Smith    Level: beginner
1162ee3001cbSBarry Smith 
116317bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
116417bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
116517bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1166ebe3b25bSBarry Smith 
1167ebe3b25bSBarry Smith M*/
1168fb2e594dSBarry Smith EXTERN_C_BEGIN
11694a2ae208SSatish Balay #undef __FUNCT__
11704b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
117163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
11725e42470aSBarry Smith {
1173dfbe8321SBarry Smith   PetscErrorCode ierr;
11744b27c08aSLois Curfman McInnes   SNES_LS        *neP;
11755e42470aSBarry Smith 
11763a40ed3dSBarry Smith   PetscFunctionBegin;
1177e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1178e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1179e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1180e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1181e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
11825e42470aSBarry Smith 
118338f2d2fdSLisandro Dalcin   ierr                  = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
11845e42470aSBarry Smith   snes->data    	= (void*)neP;
11855e42470aSBarry Smith   neP->alpha		= 1.e-4;
11865e42470aSBarry Smith   neP->maxstep		= 1.e8;
1187*e106eecfSBarry Smith   neP->minlambda        = 1.e-12;
118817bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1189c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
11903c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
11913c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
11923c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
11933c632250SBarry Smith   neP->precheck         = PETSC_NULL;
119482bf6240SBarry Smith 
1195557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1196557d3f75SLisandro Dalcin 					   "SNESLineSearchSet_LS",
1197557d3f75SLisandro Dalcin 					   SNESLineSearchSet_LS);CHKERRQ(ierr);
1198557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1199557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPostCheck_LS",
1200557d3f75SLisandro Dalcin 					   SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
1201557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1202557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPreCheck_LS",
1203557d3f75SLisandro Dalcin 					   SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
120482bf6240SBarry Smith 
12053a40ed3dSBarry Smith   PetscFunctionReturn(0);
12065e42470aSBarry Smith }
1207fb2e594dSBarry Smith EXTERN_C_END
120882bf6240SBarry Smith 
120982bf6240SBarry Smith 
121082bf6240SBarry Smith 
1211