xref: /petsc/src/snes/impls/ls/ls.c (revision 4936397d005f346ba90e3a30e7b35f4332998b76)
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"
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);
160*4936397dSBarry Smith   if (snes->domainerror) {
161*4936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
162*4936397dSBarry Smith     PetscFunctionReturn(0);
163*4936397dSBarry Smith   }
164cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
16543e71028SBarry Smith   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1660f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1675e42470aSBarry Smith   snes->norm = fnorm;
1680f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16984c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
17094a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1713f149594SLisandro Dalcin 
1723f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
1733f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
17485385478SLisandro Dalcin   /* test convergence */
17585385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
17606ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
177d034289bSBarry Smith 
1785e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1795e76c431SBarry Smith 
180329e7e40SMatthew Knepley     /* Call general purpose update function */
181e7788613SBarry Smith     if (snes->ops->update) {
182e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
183de8cb200SMatthew Knepley     }
184329e7e40SMatthew Knepley 
185ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1865334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18794b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18871f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
1893d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1903d4c4710SBarry Smith     if (kspreason < 0) {
1913d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
192ef998cc9SBarry 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);
1933d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
194ab81109fSSatish Balay         break;
1953d4c4710SBarry Smith       }
1963d4c4710SBarry Smith     }
1973d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1983f149594SLisandro Dalcin     snes->linear_its += lits;
1993f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
20074637425SBarry Smith 
2019c3ca13aSBarry Smith     if (neP->precheckstep) {
2029c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
2039c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
2049c3ca13aSBarry Smith     }
2059c3ca13aSBarry Smith 
206b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
2071e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
20885471664SBarry Smith     }
20974637425SBarry Smith 
210ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
211ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
212e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
213ea4d3ed3SLois Curfman McInnes     */
21485385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
2153f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
216a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
217ae15b995SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
2184a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
219*4936397dSBarry Smith     if (snes->domainerror) {
220*4936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
221*4936397dSBarry Smith       PetscFunctionReturn(0);
222*4936397dSBarry Smith     }
223a7cc72afSBarry Smith     if (!lssucceed) {
22450ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
22585385478SLisandro Dalcin 	PetscTruth ismin;
2263505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2271e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
2288a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2293505fcc1SBarry Smith         break;
2303505fcc1SBarry Smith       }
23150ffb88aSMatthew Knepley     }
23285385478SLisandro Dalcin     /* Update function and solution vectors */
23385385478SLisandro Dalcin     fnorm = gnorm;
23485385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
23585385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
23685385478SLisandro Dalcin     /* Monitor convergence */
23785385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
23885385478SLisandro Dalcin     snes->iter = i+1;
23985385478SLisandro Dalcin     snes->norm = fnorm;
24085385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
24185385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
24285385478SLisandro Dalcin     SNESMonitor(snes,snes->iter,snes->norm);
24385385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
24485385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
245e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2463f149594SLisandro Dalcin     if (snes->reason) break;
24716c95cacSBarry Smith   }
24852392280SLois Curfman McInnes   if (i == maxits) {
249ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
25085385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
25152392280SLois Curfman McInnes   }
2523a40ed3dSBarry Smith   PetscFunctionReturn(0);
2535e76c431SBarry Smith }
25404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25504d965bbSLois Curfman McInnes /*
2564b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2574b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25804d965bbSLois Curfman McInnes 
25904d965bbSLois Curfman McInnes    Input Parameter:
26004d965bbSLois Curfman McInnes .  snes - the SNES context
26104d965bbSLois Curfman McInnes .  x - the solution vector
26204d965bbSLois Curfman McInnes 
26304d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26404d965bbSLois Curfman McInnes 
26504d965bbSLois Curfman McInnes    Notes:
2664b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26704d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26804d965bbSLois Curfman McInnes    the call to SNESSolve().
26904d965bbSLois Curfman McInnes  */
2704a2ae208SSatish Balay #undef __FUNCT__
2714b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
272dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2735e76c431SBarry Smith {
274dfbe8321SBarry Smith   PetscErrorCode ierr;
2753a40ed3dSBarry Smith 
2763a40ed3dSBarry Smith   PetscFunctionBegin;
27785385478SLisandro Dalcin   if (!snes->vec_sol_update) {
27885385478SLisandro Dalcin     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
27985385478SLisandro Dalcin     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
28085385478SLisandro Dalcin   }
281e74804ceSBarry Smith   if (!snes->work) {
28285385478SLisandro Dalcin     snes->nwork = 3;
28385385478SLisandro Dalcin     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
284efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
285e74804ceSBarry Smith   }
2863a40ed3dSBarry Smith   PetscFunctionReturn(0);
2875e76c431SBarry Smith }
28804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28904d965bbSLois Curfman McInnes /*
2904b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2914b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29204d965bbSLois Curfman McInnes 
29304d965bbSLois Curfman McInnes    Input Parameter:
29404d965bbSLois Curfman McInnes .  snes - the SNES context
29504d965bbSLois Curfman McInnes 
296de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29704d965bbSLois Curfman McInnes  */
2984a2ae208SSatish Balay #undef __FUNCT__
2994b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
300dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3015e76c431SBarry Smith {
302dfbe8321SBarry Smith   PetscErrorCode ierr;
3033a40ed3dSBarry Smith 
3043a40ed3dSBarry Smith   PetscFunctionBegin;
30585385478SLisandro Dalcin   if (snes->vec_sol_update) {
30685385478SLisandro Dalcin     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
30785385478SLisandro Dalcin     snes->vec_sol_update = PETSC_NULL;
30885385478SLisandro Dalcin   }
3095baf8537SBarry Smith   if (snes->nwork) {
3104b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
31185385478SLisandro Dalcin     snes->nwork = 0;
31285385478SLisandro Dalcin     snes->work  = PETSC_NULL;
31321c89e3eSBarry Smith   }
3145c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
315557d3f75SLisandro Dalcin 
316557d3f75SLisandro Dalcin   /* clear composed functions */
317557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
318557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
319557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
320557d3f75SLisandro Dalcin 
3213a40ed3dSBarry Smith   PetscFunctionReturn(0);
3225e76c431SBarry Smith }
32304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3244a2ae208SSatish Balay #undef __FUNCT__
32517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
32682bf6240SBarry Smith 
3274b828684SBarry Smith /*@C
32817bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3295e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3305e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3315e76c431SBarry Smith 
332c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
333c7afd0dbSLois Curfman McInnes 
3345e76c431SBarry Smith    Input Parameters:
335c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
336af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3375e76c431SBarry Smith .  x - current iterate
3385e76c431SBarry Smith .  f - residual evaluated at x
3393c632250SBarry Smith .  y - search direction
3405e76c431SBarry Smith .  w - work vector
341c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3425e76c431SBarry Smith 
343c4a48953SLois Curfman McInnes    Output Parameters:
344c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3453c632250SBarry Smith .  w - new iterate
3465e76c431SBarry Smith .  gnorm - 2-norm of g
3475e76c431SBarry Smith .  ynorm - 2-norm of search length
348a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
349fee21e36SBarry Smith 
350c4a48953SLois Curfman McInnes    Options Database Key:
35117bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
352c4a48953SLois Curfman McInnes 
35336851e7fSLois Curfman McInnes    Level: advanced
35436851e7fSLois Curfman McInnes 
35528ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
35628ae5a21SLois Curfman McInnes 
35717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
35817bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3595e76c431SBarry Smith @*/
36017bae607SBarry 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)
3615e76c431SBarry Smith {
362dfbe8321SBarry Smith   PetscErrorCode ierr;
3634b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3643c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3655334005bSBarry Smith 
3663a40ed3dSBarry Smith   PetscFunctionBegin;
367a7cc72afSBarry Smith   *flag = PETSC_TRUE;
368d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
369a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
37079f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3713c632250SBarry Smith   if (neP->postcheckstep) {
3723c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3738f99978dSLois Curfman McInnes   }
3743c632250SBarry Smith   if (changed_y) {
37579f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3763c632250SBarry Smith   }
377*4936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
378*4936397dSBarry Smith   if (!snes->domainerror) {
379a3b38805SLois Curfman McInnes     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
38043e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
381*4936397dSBarry Smith   }
382d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3833a40ed3dSBarry Smith   PetscFunctionReturn(0);
3845e76c431SBarry Smith }
38504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3864a2ae208SSatish Balay #undef __FUNCT__
38717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
38882bf6240SBarry Smith 
38929e0b56fSBarry Smith /*@C
39017bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
39129e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
39229e0b56fSBarry Smith    even compute the norm of the function or search direction; this
39329e0b56fSBarry Smith    is intended only when you know the full step is fine and are
39429e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
395c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
396c7afd0dbSLois Curfman McInnes 
397c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
39829e0b56fSBarry Smith 
39929e0b56fSBarry Smith    Input Parameters:
400c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
401af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
40229e0b56fSBarry Smith .  x - current iterate
40329e0b56fSBarry Smith .  f - residual evaluated at x
4043c632250SBarry Smith .  y - search direction
40529e0b56fSBarry Smith .  w - work vector
406c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
40729e0b56fSBarry Smith 
40829e0b56fSBarry Smith    Output Parameters:
409c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4103c632250SBarry Smith .  w - new iterate
41129e0b56fSBarry Smith .  gnorm - not changed
41229e0b56fSBarry Smith .  ynorm - not changed
413a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
414fee21e36SBarry Smith 
41529e0b56fSBarry Smith    Options Database Key:
41617bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
41729e0b56fSBarry Smith 
4188cbba510SBarry Smith    Notes:
41917bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
420ea56f5baSLois Curfman McInnes    either the options
421ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
422ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
423329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
424329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
425329f5518SBarry Smith 
426329f5518SBarry Smith    During the final iteration this will not evaluate the function at
427329f5518SBarry Smith    the solution point. This is to save a function evaluation while
428329f5518SBarry Smith    using pseudo-timestepping.
4298cbba510SBarry Smith 
430ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
431a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
432ea56f5baSLois Curfman McInnes    correct, since they are not computed.
433ea56f5baSLois Curfman McInnes 
434ea56f5baSLois Curfman McInnes    Level: advanced
4358cbba510SBarry Smith 
43629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
43729e0b56fSBarry Smith 
43817bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
43917bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
44029e0b56fSBarry Smith @*/
44117bae607SBarry 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)
44229e0b56fSBarry Smith {
443dfbe8321SBarry Smith   PetscErrorCode ierr;
4444b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4453c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
44629e0b56fSBarry Smith 
4473a40ed3dSBarry Smith   PetscFunctionBegin;
448a7cc72afSBarry Smith   *flag = PETSC_TRUE;
449d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
45079f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4513c632250SBarry Smith   if (neP->postcheckstep) {
4523c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4533c632250SBarry Smith   }
4543c632250SBarry Smith   if (changed_y) {
45579f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4568f99978dSLois Curfman McInnes   }
457329f5518SBarry Smith 
458329f5518SBarry Smith   /* don't evaluate function the last time through */
459329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
460*4936397dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
461329f5518SBarry Smith   }
462d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4633a40ed3dSBarry Smith   PetscFunctionReturn(0);
46429e0b56fSBarry Smith }
46504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4664a2ae208SSatish Balay #undef __FUNCT__
46717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4684b828684SBarry Smith /*@C
46917bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4705e76c431SBarry Smith 
471c7afd0dbSLois Curfman McInnes    Collective on SNES
472c7afd0dbSLois Curfman McInnes 
4735e76c431SBarry Smith    Input Parameters:
474c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
475af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4765e76c431SBarry Smith .  x - current iterate
4775e76c431SBarry Smith .  f - residual evaluated at x
4783c632250SBarry Smith .  y - search direction
4795e76c431SBarry Smith .  w - work vector
480c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4815e76c431SBarry Smith 
482393d2d9aSLois Curfman McInnes    Output Parameters:
483c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4843c632250SBarry Smith .  w - new iterate
4855e76c431SBarry Smith .  gnorm - 2-norm of g
4865e76c431SBarry Smith .  ynorm - 2-norm of search length
487a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
488fee21e36SBarry Smith 
489c4a48953SLois Curfman McInnes    Options Database Key:
49017bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
491c4a48953SLois Curfman McInnes 
4925e76c431SBarry Smith    Notes:
4935e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4945e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4955e76c431SBarry Smith 
49636851e7fSLois Curfman McInnes    Level: advanced
49736851e7fSLois Curfman McInnes 
49828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49928ae5a21SLois Curfman McInnes 
50017bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5015e76c431SBarry Smith @*/
50217bae607SBarry 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)
5035e76c431SBarry Smith {
504406556e6SLois Curfman McInnes   /*
505406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
506406556e6SLois Curfman McInnes      minimization problem:
507406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
508406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
509406556e6SLois Curfman McInnes    */
510406556e6SLois Curfman McInnes 
5115ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
512275d25c3SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
513aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
514ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5156b5873e3SBarry Smith #endif
516dfbe8321SBarry Smith   PetscErrorCode ierr;
517a7cc72afSBarry Smith   PetscInt       count;
5184b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
51979f36460SBarry Smith   PetscScalar    scale;
5203c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5215e76c431SBarry Smith 
5223a40ed3dSBarry Smith   PetscFunctionBegin;
523d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
524a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5255e76c431SBarry Smith   alpha   = neP->alpha;
5265e76c431SBarry Smith   maxstep = neP->maxstep;
5275e76c431SBarry Smith   steptol = neP->steptol;
5285e76c431SBarry Smith 
529cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5309e247f21SBarry Smith   if (!*ynorm) {
531ae15b995SBarry Smith     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
532a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5333c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
534a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
535a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
536ad922ac9SBarry Smith     goto theend1;
53794a9d846SBarry Smith   }
5385e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5395e42470aSBarry Smith     scale = maxstep/(*ynorm);
540ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr);
5412dcb1b2aSMatthew Knepley     ierr = VecScale(y,scale);CHKERRQ(ierr);
5425e76c431SBarry Smith     *ynorm = maxstep;
5435e76c431SBarry Smith   }
544ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5455ca10a37SBarry Smith   minlambda = steptol/rellength;
546a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
547aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
548a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
549329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5505e42470aSBarry Smith #else
551a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5525e42470aSBarry Smith #endif
5535e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5545e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5555e76c431SBarry Smith 
556e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
55743e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
558ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
55943e71028SBarry Smith     *flag = PETSC_FALSE;
56043e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
56143e71028SBarry Smith     goto theend1;
56243e71028SBarry Smith   }
563*4936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
564*4936397dSBarry Smith   if (snes->domainerror) {
565*4936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
566*4936397dSBarry Smith     PetscFunctionReturn(0);
56719717074SBarry Smith   }
568313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56943e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5705d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
571ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
572ad922ac9SBarry Smith     goto theend1;
5735e76c431SBarry Smith   }
5745e76c431SBarry Smith 
5755e76c431SBarry Smith   /* Fit points with quadratic */
576313b5538SBarry Smith   lambda     = 1.0;
577a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5785e76c431SBarry Smith   lambdaprev = lambda;
5795e76c431SBarry Smith   gnormprev  = *gnorm;
580329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
581ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
582ddd12767SBarry Smith   else                         lambda = lambdatemp;
583275d25c3SBarry Smith 
584aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
585a23f76dfSSatish Balay   clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
5865e42470aSBarry Smith #else
587e68848bdSBarry Smith   ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
5885e42470aSBarry Smith #endif
58943e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
590ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
59143e71028SBarry Smith     *flag = PETSC_FALSE;
59243e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
59343e71028SBarry Smith     goto theend1;
59443e71028SBarry Smith   }
59519717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
596*4936397dSBarry Smith   if (snes->domainerror) {
597*4936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
598*4936397dSBarry Smith     PetscFunctionReturn(0);
59919717074SBarry Smith   }
600cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
60143e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6025ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
603ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
604ad922ac9SBarry Smith     goto theend1;
6055e76c431SBarry Smith   }
6065e76c431SBarry Smith 
6075e76c431SBarry Smith   /* Fit points with cubic */
6085e76c431SBarry Smith   count = 1;
6098229bfc2SMatthew Knepley   while (count < 10000) {
6105e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
611ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
612ae15b995SBarry 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);
61343e71028SBarry Smith       *flag = PETSC_FALSE;
61443e71028SBarry Smith       break;
6155e76c431SBarry Smith     }
6165d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6175d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
618ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6192b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6205e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6215e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6225e76c431SBarry Smith     if (a == 0.0) {
6235e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6245e76c431SBarry Smith     } else {
6255e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6265e76c431SBarry Smith     }
6275e76c431SBarry Smith     lambdaprev = lambda;
6285e76c431SBarry Smith     gnormprev  = *gnorm;
629329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
630329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6315e76c431SBarry Smith     else                         lambda     = lambdatemp;
632aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
633275d25c3SBarry Smith     clambda   = lambda;
634e68848bdSBarry Smith     ierr      = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
6355e42470aSBarry Smith #else
636e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
6375e42470aSBarry Smith #endif
63843e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
639ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
640ae15b995SBarry 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);
641ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
64243e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
643ed98166eSMatthew Knepley       break;
644ed98166eSMatthew Knepley     }
64519717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
646*4936397dSBarry Smith     if (snes->domainerror) {
647*4936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
648*4936397dSBarry Smith       PetscFunctionReturn(0);
64919717074SBarry Smith     }
650cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
65143e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6525ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
653ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
65443e71028SBarry Smith       break;
6552b022350SLois Curfman McInnes     } else {
656ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);CHKERRQ(ierr);
6575e76c431SBarry Smith     }
6585e76c431SBarry Smith     count++;
6595e76c431SBarry Smith   }
6608229bfc2SMatthew Knepley   if (count >= 10000) {
661cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6628229bfc2SMatthew Knepley   }
663ad922ac9SBarry Smith   theend1:
6648f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6653c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6663c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6673c632250SBarry Smith     if (changed_y) {
66879f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6693c632250SBarry Smith     }
6703c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6713c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
672*4936397dSBarry Smith       if (snes->domainerror) {
673*4936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
674*4936397dSBarry Smith         PetscFunctionReturn(0);
67519717074SBarry Smith       }
6768f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67743e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6783c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6798f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6803c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6818f99978dSLois Curfman McInnes     }
6828f99978dSLois Curfman McInnes   }
683d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6843a40ed3dSBarry Smith   PetscFunctionReturn(0);
6855e76c431SBarry Smith }
68604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6874a2ae208SSatish Balay #undef __FUNCT__
68817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6894b828684SBarry Smith /*@C
69017bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6915e76c431SBarry Smith 
692c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
693c7afd0dbSLois Curfman McInnes 
6945e42470aSBarry Smith    Input Parameters:
695c7afd0dbSLois Curfman McInnes +  snes - the SNES context
696af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6975e42470aSBarry Smith .  x - current iterate
6985e42470aSBarry Smith .  f - residual evaluated at x
6993c632250SBarry Smith .  y - search direction
7005e42470aSBarry Smith .  w - work vector
701c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
7025e42470aSBarry Smith 
703c4a48953SLois Curfman McInnes    Output Parameters:
7047f3332b4SBarry Smith +  g - residual evaluated at new iterate w
7057f3332b4SBarry Smith .  w - new iterate (x + alpha*y)
7065e42470aSBarry Smith .  gnorm - 2-norm of g
7075e42470aSBarry Smith .  ynorm - 2-norm of search length
708a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
709fee21e36SBarry Smith 
710ce9499c7SBarry Smith    Options Database Keys:
711ce9499c7SBarry Smith +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
712ce9499c7SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
713ce9499c7SBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
714ce9499c7SBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
715ce9499c7SBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
716ce9499c7SBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
7175e42470aSBarry Smith    Notes:
71817bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7195e42470aSBarry Smith 
72036851e7fSLois Curfman McInnes    Level: advanced
72136851e7fSLois Curfman McInnes 
722f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
723f59ffdedSLois Curfman McInnes 
72417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7255e42470aSBarry Smith @*/
72617bae607SBarry 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)
7275e76c431SBarry Smith {
728406556e6SLois Curfman McInnes   /*
729406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
730406556e6SLois Curfman McInnes      minimization problem:
731406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
732406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
733406556e6SLois Curfman McInnes    */
734275d25c3SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
735aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
736ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7376b5873e3SBarry Smith #endif
738dfbe8321SBarry Smith   PetscErrorCode ierr;
739a7cc72afSBarry Smith   PetscInt       count;
7404b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
74179f36460SBarry Smith   PetscScalar    scale;
7423c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7435e76c431SBarry Smith 
7443a40ed3dSBarry Smith   PetscFunctionBegin;
745d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
746a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7475e42470aSBarry Smith   alpha   = neP->alpha;
7485e42470aSBarry Smith   maxstep = neP->maxstep;
7495e42470aSBarry Smith   steptol = neP->steptol;
7505e76c431SBarry Smith 
7513505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
75263b9fa5eSBarry Smith   if (*ynorm == 0.0) {
753ae15b995SBarry Smith     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
754b37302c6SLois Curfman McInnes     *gnorm = fnorm;
755e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
756b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
757a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
758ad922ac9SBarry Smith     goto theend2;
75994a9d846SBarry Smith   }
7605e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7615e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7622dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7635e42470aSBarry Smith     *ynorm = maxstep;
7645e76c431SBarry Smith   }
765ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7665ca10a37SBarry Smith   minlambda = steptol/rellength;
767a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
768aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
769a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
770329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7715e42470aSBarry Smith #else
772a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7735e42470aSBarry Smith #endif
7745e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7755e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7765e42470aSBarry Smith 
777e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
77843e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
779ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
78043e71028SBarry Smith     *flag = PETSC_FALSE;
78143e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
78243e71028SBarry Smith     goto theend2;
78343e71028SBarry Smith   }
784*4936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
785*4936397dSBarry Smith   if (snes->domainerror) {
786*4936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
787*4936397dSBarry Smith     PetscFunctionReturn(0);
78819717074SBarry Smith   }
789cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
79043e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7915d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
792ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
793ad922ac9SBarry Smith     goto theend2;
7945e42470aSBarry Smith   }
7955e42470aSBarry Smith 
7965e42470aSBarry Smith   /* Fit points with quadratic */
797313b5538SBarry Smith   lambda = 1.0;
7985e42470aSBarry Smith   count = 1;
7995ca10a37SBarry Smith   while (PETSC_TRUE) {
8005e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
801ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
802ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
803e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
80443e71028SBarry Smith       *flag = PETSC_FALSE;
80543e71028SBarry Smith       break;
8065e42470aSBarry Smith     }
807a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
808329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
809329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
810329f5518SBarry Smith     else                         lambda     = lambdatemp;
811275d25c3SBarry Smith 
812aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
813e68848bdSBarry Smith     clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
8145e42470aSBarry Smith #else
815e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8165e42470aSBarry Smith #endif
81743e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
818ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
819ae15b995SBarry 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);
820ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
82143e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
822ed98166eSMatthew Knepley       break;
823ed98166eSMatthew Knepley     }
82419717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
825*4936397dSBarry Smith     if (snes->domainerror) {
826*4936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
827*4936397dSBarry Smith       PetscFunctionReturn(0);
82819717074SBarry Smith     }
829cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
83043e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8315ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
832ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
83306259719SBarry Smith       break;
8345e42470aSBarry Smith     }
8355e42470aSBarry Smith     count++;
8365e42470aSBarry Smith   }
837ad922ac9SBarry Smith   theend2:
8388f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8393c632250SBarry Smith   if (neP->postcheckstep) {
8403c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8413c632250SBarry Smith     if (changed_y) {
84279f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8433c632250SBarry Smith     }
8443c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8453c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
846*4936397dSBarry Smith       if (snes->domainerror) {
847*4936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
848*4936397dSBarry Smith         PetscFunctionReturn(0);
84919717074SBarry Smith       }
8508f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8513c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8528f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8533c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8543c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8558f99978dSLois Curfman McInnes     }
8568f99978dSLois Curfman McInnes   }
857d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8583a40ed3dSBarry Smith   PetscFunctionReturn(0);
8595e76c431SBarry Smith }
8602343ba6eSBarry Smith 
86104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8624a2ae208SSatish Balay #undef __FUNCT__
86317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
864c9e6a524SLois Curfman McInnes /*@C
86517bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8664b27c08aSLois Curfman McInnes    by the method SNESLS.
8675e76c431SBarry Smith 
8685e76c431SBarry Smith    Input Parameters:
869c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
870af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
871c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8725e76c431SBarry Smith 
873fee21e36SBarry Smith    Collective on SNES
874fee21e36SBarry Smith 
875c4a48953SLois Curfman McInnes    Available Routines:
87617bae607SBarry Smith +  SNESLineSearchCubic() - default line search
87717bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
87817bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
87917bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8805e76c431SBarry Smith 
881c4a48953SLois Curfman McInnes     Options Database Keys:
8824b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8834b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8844b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8854b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8863304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8873304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
888c4a48953SLois Curfman McInnes 
8895e76c431SBarry Smith    Calling sequence of func:
890bd208895SLois Curfman McInnes .vb
891af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
892a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
893bd208895SLois Curfman McInnes .ve
8945e76c431SBarry Smith 
8955e76c431SBarry Smith     Input parameters for func:
896c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
897af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8985e76c431SBarry Smith .   x - current iterate
8995e76c431SBarry Smith .   f - residual evaluated at x
9003c632250SBarry Smith .   y - search direction
9015e76c431SBarry Smith .   w - work vector
902c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9035e76c431SBarry Smith 
9045e76c431SBarry Smith     Output parameters for func:
905c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9063c632250SBarry Smith .   w - new iterate
9075e76c431SBarry Smith .   gnorm - 2-norm of g
9085e76c431SBarry Smith .   ynorm - 2-norm of search length
909a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
910f59ffdedSLois Curfman McInnes 
91136851e7fSLois Curfman McInnes     Level: advanced
91236851e7fSLois Curfman McInnes 
913f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
914f59ffdedSLois Curfman McInnes 
91517bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
91617bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9175e76c431SBarry Smith @*/
91817bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9195e76c431SBarry Smith {
920a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
92182bf6240SBarry Smith 
9223a40ed3dSBarry Smith   PetscFunctionBegin;
92317bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
92482bf6240SBarry Smith   if (f) {
925af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92682bf6240SBarry Smith   }
9273a40ed3dSBarry Smith   PetscFunctionReturn(0);
9285e76c431SBarry Smith }
9298e019c35SBarry Smith 
930a7cc72afSBarry 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*/
93104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
932fb2e594dSBarry Smith EXTERN_C_BEGIN
9334a2ae208SSatish Balay #undef __FUNCT__
93417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
93517bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
93682bf6240SBarry Smith {
93782bf6240SBarry Smith   PetscFunctionBegin;
9384b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9394b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
94082bf6240SBarry Smith   PetscFunctionReturn(0);
94182bf6240SBarry Smith }
942fb2e594dSBarry Smith EXTERN_C_END
94304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9444a2ae208SSatish Balay #undef __FUNCT__
9453c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
946c8dd0c0dSLois Curfman McInnes /*@C
9473c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9484b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
949c8dd0c0dSLois Curfman McInnes 
950c8dd0c0dSLois Curfman McInnes    Input Parameters:
951c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9523c632250SBarry Smith .  func - pointer to function
953c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
954c8dd0c0dSLois Curfman McInnes 
955c8dd0c0dSLois Curfman McInnes    Collective on SNES
956c8dd0c0dSLois Curfman McInnes 
957c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
958c8dd0c0dSLois Curfman McInnes .vb
9593c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
960c8dd0c0dSLois Curfman McInnes .ve
961b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
962b1ae27eaSLois Curfman McInnes    on failure.
963c8dd0c0dSLois Curfman McInnes 
964c8dd0c0dSLois Curfman McInnes    Input parameters for func:
965c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
966c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9673c632250SBarry Smith .  x - previous iterate
9683c632250SBarry Smith .  y - new search direction and length
9693c632250SBarry Smith -  w - current candidate iterate
970c8dd0c0dSLois Curfman McInnes 
971c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9723c632250SBarry Smith +  y - search direction (possibly changed)
9733c632250SBarry Smith .  w - current iterate (possibly modified)
9743c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9753c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
976c8dd0c0dSLois Curfman McInnes 
977c8dd0c0dSLois Curfman McInnes    Level: advanced
978c8dd0c0dSLois Curfman McInnes 
9799e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9809e247f21SBarry Smith 
9813c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9823c632250SBarry Smith 
9833c632250SBarry Smith    On input w = x + y
9843c632250SBarry Smith 
98517bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
986b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
987ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9888f99978dSLois Curfman McInnes 
98917bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
990ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
991ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
992ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9939e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
994b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9958f99978dSLois Curfman McInnes 
996c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
997c8dd0c0dSLois Curfman McInnes 
99817bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
999c8dd0c0dSLois Curfman McInnes @*/
10003c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
1001c8dd0c0dSLois Curfman McInnes {
10023c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
1003c8dd0c0dSLois Curfman McInnes 
1004c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10053c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1006c8dd0c0dSLois Curfman McInnes   if (f) {
1007c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1008c8dd0c0dSLois Curfman McInnes   }
1009c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1010c8dd0c0dSLois Curfman McInnes }
10119c3ca13aSBarry Smith #undef __FUNCT__
10129c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10139c3ca13aSBarry Smith /*@C
10149c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10157e4bb74cSBarry Smith          before the line search is called.
10169c3ca13aSBarry Smith 
10179c3ca13aSBarry Smith    Input Parameters:
10189c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10199c3ca13aSBarry Smith .  func - pointer to function
10209c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10219c3ca13aSBarry Smith 
10229c3ca13aSBarry Smith    Collective on SNES
10239c3ca13aSBarry Smith 
10249c3ca13aSBarry Smith    Calling sequence of func:
10259c3ca13aSBarry Smith .vb
10261e3347e8SBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
10279c3ca13aSBarry Smith .ve
10289c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10299c3ca13aSBarry Smith    on failure.
10309c3ca13aSBarry Smith 
10319c3ca13aSBarry Smith    Input parameters for func:
10329c3ca13aSBarry Smith +  snes - nonlinear context
10339c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10349c3ca13aSBarry Smith .  x - previous iterate
10359c3ca13aSBarry Smith -  y - new search direction and length
10369c3ca13aSBarry Smith 
10379c3ca13aSBarry Smith    Output parameters for func:
10389c3ca13aSBarry Smith +  y - search direction (possibly changed)
10399c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10409c3ca13aSBarry Smith 
10419c3ca13aSBarry Smith    Level: advanced
10429c3ca13aSBarry Smith 
10439c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10449c3ca13aSBarry Smith 
10459c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10469c3ca13aSBarry Smith 
10477e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10489c3ca13aSBarry Smith @*/
10499c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10509c3ca13aSBarry Smith {
10519c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10529c3ca13aSBarry Smith 
10539c3ca13aSBarry Smith   PetscFunctionBegin;
10549c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10559c3ca13aSBarry Smith   if (f) {
10569c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10579c3ca13aSBarry Smith   }
10589c3ca13aSBarry Smith   PetscFunctionReturn(0);
10599c3ca13aSBarry Smith }
10609c3ca13aSBarry Smith 
1061c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10629c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10639c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1064c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10654a2ae208SSatish Balay #undef __FUNCT__
10663c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10679c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1068c8dd0c0dSLois Curfman McInnes {
1069c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10703c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10713c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1072c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1073c8dd0c0dSLois Curfman McInnes }
1074c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10759c3ca13aSBarry Smith 
10769c3ca13aSBarry Smith EXTERN_C_BEGIN
10779c3ca13aSBarry Smith #undef __FUNCT__
10789c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10799c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10809c3ca13aSBarry Smith {
10819c3ca13aSBarry Smith   PetscFunctionBegin;
10829c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10839c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10849c3ca13aSBarry Smith   PetscFunctionReturn(0);
10859c3ca13aSBarry Smith }
10869c3ca13aSBarry Smith EXTERN_C_END
1087329e7e40SMatthew Knepley 
1088329e7e40SMatthew Knepley /*
10894b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
109004d965bbSLois Curfman McInnes 
109104d965bbSLois Curfman McInnes    Input Parameters:
109204d965bbSLois Curfman McInnes .  SNES - the SNES context
109304d965bbSLois Curfman McInnes .  viewer - visualization context
109404d965bbSLois Curfman McInnes 
109504d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
109604d965bbSLois Curfman McInnes */
10974a2ae208SSatish Balay #undef __FUNCT__
10984b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10996849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1100a935fc98SLois Curfman McInnes {
11014b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11022fc52814SBarry Smith   const char     *cstr;
1103dfbe8321SBarry Smith   PetscErrorCode ierr;
110432077d6dSBarry Smith   PetscTruth     iascii;
1105a935fc98SLois Curfman McInnes 
11063a40ed3dSBarry Smith   PetscFunctionBegin;
110732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
110832077d6dSBarry Smith   if (iascii) {
110917bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
111017bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
111117bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
111219bcc07fSBarry Smith     else                                                cstr = "unknown";
1113b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1114a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11155cd90555SBarry Smith   } else {
111679a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
111719bcc07fSBarry Smith   }
11183a40ed3dSBarry Smith   PetscFunctionReturn(0);
1119a935fc98SLois Curfman McInnes }
112004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
112104d965bbSLois Curfman McInnes /*
11224b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
112304d965bbSLois Curfman McInnes 
112404d965bbSLois Curfman McInnes    Input Parameter:
112504d965bbSLois Curfman McInnes .  snes - the SNES context
112604d965bbSLois Curfman McInnes 
112704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
112804d965bbSLois Curfman McInnes */
11294a2ae208SSatish Balay #undef __FUNCT__
11304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11316849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11325e42470aSBarry Smith {
11334b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1134e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1135dfbe8321SBarry Smith   PetscErrorCode ierr;
1136a7cc72afSBarry Smith   PetscInt       indx;
1137f1af5d2fSBarry Smith   PetscTruth     flg;
11385e42470aSBarry Smith 
11393a40ed3dSBarry Smith   PetscFunctionBegin;
1140b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11414b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11424b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11434b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1144186905e3SBarry Smith 
114517bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
114625cdf11fSBarry Smith     if (flg) {
114722e36657SBarry Smith       switch (indx) {
1148b49fd9e1SBarry Smith       case 0:
114917bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1150b49fd9e1SBarry Smith         break;
1151b49fd9e1SBarry Smith       case 1:
115217bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1153b49fd9e1SBarry Smith         break;
1154b49fd9e1SBarry Smith       case 2:
115517bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1156b49fd9e1SBarry Smith         break;
1157b49fd9e1SBarry Smith       case 3:
115817bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1159b49fd9e1SBarry Smith         break;
11605e42470aSBarry Smith       }
11615e42470aSBarry Smith     }
1162b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11633a40ed3dSBarry Smith   PetscFunctionReturn(0);
11645e42470aSBarry Smith }
116504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1166ebe3b25bSBarry Smith /*MC
1167ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
116804d965bbSLois Curfman McInnes 
1169ebe3b25bSBarry Smith    Options Database:
1170ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1171ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1172ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1173ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1174ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1175ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
117604d965bbSLois Curfman McInnes 
1177ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
117804d965bbSLois Curfman McInnes 
1179ee3001cbSBarry Smith    Level: beginner
1180ee3001cbSBarry Smith 
118117bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
118217bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
118317bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1184ebe3b25bSBarry Smith 
1185ebe3b25bSBarry Smith M*/
1186fb2e594dSBarry Smith EXTERN_C_BEGIN
11874a2ae208SSatish Balay #undef __FUNCT__
11884b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
118963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
11905e42470aSBarry Smith {
1191dfbe8321SBarry Smith   PetscErrorCode ierr;
11924b27c08aSLois Curfman McInnes   SNES_LS        *neP;
11935e42470aSBarry Smith 
11943a40ed3dSBarry Smith   PetscFunctionBegin;
1195e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1196e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1197e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1198e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1199e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
12005e42470aSBarry Smith 
120138f2d2fdSLisandro Dalcin   ierr                  = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
12025e42470aSBarry Smith   snes->data    	= (void*)neP;
12035e42470aSBarry Smith   neP->alpha		= 1.e-4;
12045e42470aSBarry Smith   neP->maxstep		= 1.e8;
12055e42470aSBarry Smith   neP->steptol		= 1.e-12;
120617bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1207c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12083c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12093c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12103c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12113c632250SBarry Smith   neP->precheck         = PETSC_NULL;
121282bf6240SBarry Smith 
1213557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1214557d3f75SLisandro Dalcin 					   "SNESLineSearchSet_LS",
1215557d3f75SLisandro Dalcin 					   SNESLineSearchSet_LS);CHKERRQ(ierr);
1216557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1217557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPostCheck_LS",
1218557d3f75SLisandro Dalcin 					   SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
1219557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1220557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPreCheck_LS",
1221557d3f75SLisandro Dalcin 					   SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
122282bf6240SBarry Smith 
12233a40ed3dSBarry Smith   PetscFunctionReturn(0);
12245e42470aSBarry Smith }
1225fb2e594dSBarry Smith EXTERN_C_END
122682bf6240SBarry Smith 
122782bf6240SBarry Smith 
122882bf6240SBarry Smith 
1229