xref: /petsc/src/snes/impls/ls/ls.c (revision 55d4ea136f85ef5893a626c7e78d838f07b3f308)
15e76c431SBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h>
35e42470aSBarry Smith 
421167be9SJed Brown const char *const SNESLineSearchTypes[] = {"BASIC","BASICNONORMS","QUADRATIC","CUBIC","SNESLineSearchType","SNES_LS_",0};
521167be9SJed Brown 
621167be9SJed Brown const char *SNESLineSearchTypeName(SNESLineSearchType type)
721167be9SJed Brown {
821167be9SJed Brown   return (SNES_LS_BASIC <= type && type <= SNES_LS_CUBIC) ? SNESLineSearchTypes[type] : "unknown";
921167be9SJed Brown }
10490ca623SBarry Smith 
118a5d9ee4SBarry Smith /*
1220f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
1320f52e01SBarry 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
1436669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
1520f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
168a5d9ee4SBarry Smith */
174a2ae208SSatish Balay #undef __FUNCT__
184a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
19ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
208a5d9ee4SBarry Smith {
218a5d9ee4SBarry Smith   PetscReal      a1;
22dfbe8321SBarry Smith   PetscErrorCode ierr;
23ace3abfcSBarry Smith   PetscBool      hastranspose;
248a5d9ee4SBarry Smith 
258a5d9ee4SBarry Smith   PetscFunctionBegin;
268a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2736669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2836669109SBarry Smith   if (hastranspose) {
298a5d9ee4SBarry Smith     /* Compute || J^T F|| */
308a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
318a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
328f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
338a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
3436669109SBarry Smith   } else {
3536669109SBarry Smith     Vec         work;
36ea709b57SSatish Balay     PetscScalar result;
3736669109SBarry Smith     PetscReal   wnorm;
3836669109SBarry Smith 
39abb0e124SMatthew Knepley     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
4036669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
4136669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
4236669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
4336669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
446bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
4536669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
468f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
4736669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4836669109SBarry Smith   }
498a5d9ee4SBarry Smith   PetscFunctionReturn(0);
508a5d9ee4SBarry Smith }
518a5d9ee4SBarry Smith 
5274637425SBarry Smith /*
535ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
5474637425SBarry Smith */
554a2ae208SSatish Balay #undef __FUNCT__
564a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
571e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
5874637425SBarry Smith {
5974637425SBarry Smith   PetscReal      a1,a2;
60dfbe8321SBarry Smith   PetscErrorCode ierr;
61ace3abfcSBarry Smith   PetscBool      hastranspose;
6274637425SBarry Smith 
6374637425SBarry Smith   PetscFunctionBegin;
6474637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
6574637425SBarry Smith   if (hastranspose) {
6674637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6779f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6874637425SBarry Smith 
6974637425SBarry Smith     /* Compute || J^T W|| */
7074637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
7174637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
7274637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
7375567043SBarry Smith     if (a1 != 0.0) {
748f1a2a5eSBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
7585471664SBarry Smith     }
7674637425SBarry Smith   }
7774637425SBarry Smith   PetscFunctionReturn(0);
7874637425SBarry Smith }
7974637425SBarry Smith 
8004d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
8104d965bbSLois Curfman McInnes 
8204d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
8394b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
8404d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8504d965bbSLois Curfman McInnes      respectively.
8604d965bbSLois Curfman McInnes 
87fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8804d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8904d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
90fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
9104d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
9204d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
934b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
944b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
9504d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9604d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9704d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9804d965bbSLois Curfman McInnes      for all nonlinear solvers.
9904d965bbSLois Curfman McInnes 
10004d965bbSLois Curfman McInnes      Another key routine is:
10104d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
10204d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
10304d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
10404d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10504d965bbSLois Curfman McInnes 
10604d965bbSLois Curfman McInnes      Additional basic routines are:
10704d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10804d965bbSLois Curfman McInnes                                       have actually been used.
10904d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
110186905e3SBarry Smith      SNESView().
11104d965bbSLois Curfman McInnes 
11204d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
11304d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
11404d965bbSLois Curfman McInnes      above description applies to these categories also.
11504d965bbSLois Curfman McInnes 
11604d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1175e42470aSBarry Smith /*
1184b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11904d965bbSLois Curfman McInnes    method with a line search.
1205e76c431SBarry Smith 
12104d965bbSLois Curfman McInnes    Input Parameters:
12204d965bbSLois Curfman McInnes .  snes - the SNES context
1235e76c431SBarry Smith 
12404d965bbSLois Curfman McInnes    Output Parameter:
125c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12604d965bbSLois Curfman McInnes 
12704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1285e76c431SBarry Smith 
1295e76c431SBarry Smith    Notes:
1305e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1315e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1325e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1335e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
134393d2d9aSLois Curfman McInnes    and Schnabel.
1355e76c431SBarry Smith */
1364a2ae208SSatish Balay #undef __FUNCT__
1374b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
138dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1395e76c431SBarry Smith {
1404b27c08aSLois Curfman McInnes   SNES_LS            *neP = (SNES_LS*)snes->data;
1416849ba73SBarry Smith   PetscErrorCode     ierr;
142a7cc72afSBarry Smith   PetscInt           maxits,i,lits;
143ace3abfcSBarry Smith   PetscBool          lssucceed;
144112a2221SBarry Smith   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
14585385478SLisandro Dalcin   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
14685385478SLisandro Dalcin   Vec                Y,X,F,G,W;
1473d4c4710SBarry Smith   KSPConvergedReason kspreason;
1485e76c431SBarry Smith 
1493a40ed3dSBarry Smith   PetscFunctionBegin;
1503d4c4710SBarry Smith   snes->numFailures            = 0;
1513d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
152184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
153184914b5SBarry Smith 
1545e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1555e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
15639e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1575e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1585e42470aSBarry Smith   G		= snes->work[1];
1595e42470aSBarry Smith   W		= snes->work[2];
1605e76c431SBarry Smith 
1614c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1624c49b128SBarry Smith   snes->iter = 0;
16375567043SBarry Smith   snes->norm = 0.0;
1644c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16519717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1664936397dSBarry Smith   if (snes->domainerror) {
1674936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1684936397dSBarry Smith     PetscFunctionReturn(0);
1694936397dSBarry Smith   }
1702613ca53SBarry Smith   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1712613ca53SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
1722613ca53SBarry Smith   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
1732613ca53SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
17462cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1750f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1765e42470aSBarry Smith   snes->norm = fnorm;
1770f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17884c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
1797a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1803f149594SLisandro Dalcin 
1813f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
1823f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
18385385478SLisandro Dalcin   /* test convergence */
18485385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
18506ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
186d034289bSBarry Smith 
1875e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1885e76c431SBarry Smith 
189329e7e40SMatthew Knepley     /* Call general purpose update function */
190e7788613SBarry Smith     if (snes->ops->update) {
191e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
192de8cb200SMatthew Knepley     }
193329e7e40SMatthew Knepley 
194ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1955334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
19694b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
19771f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
1983d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1993d4c4710SBarry Smith     if (kspreason < 0) {
2003d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
201ef998cc9SBarry 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);
2023d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
203ab81109fSSatish Balay         break;
2043d4c4710SBarry Smith       }
2053d4c4710SBarry Smith     }
2063d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2073f149594SLisandro Dalcin     snes->linear_its += lits;
2083f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
20974637425SBarry Smith 
2109c3ca13aSBarry Smith     if (neP->precheckstep) {
211ace3abfcSBarry Smith       PetscBool  changed_y = PETSC_FALSE;
2129c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
2139c3ca13aSBarry Smith     }
2149c3ca13aSBarry Smith 
215b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
2161e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
21785471664SBarry Smith     }
21874637425SBarry Smith 
219ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
220ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
221e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
222ea4d3ed3SLois Curfman McInnes     */
22385385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
2243f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
225dc357ad2SBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
2268f1a2a5eSBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
2274a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2284936397dSBarry Smith     if (snes->domainerror) {
2294936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2304936397dSBarry Smith       PetscFunctionReturn(0);
2314936397dSBarry Smith     }
232a7cc72afSBarry Smith     if (!lssucceed) {
23350ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
234ace3abfcSBarry Smith 	PetscBool  ismin;
235647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2361e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
2378a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2383505fcc1SBarry Smith         break;
2393505fcc1SBarry Smith       }
24050ffb88aSMatthew Knepley     }
24185385478SLisandro Dalcin     /* Update function and solution vectors */
24285385478SLisandro Dalcin     fnorm = gnorm;
24385385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
24485385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
24585385478SLisandro Dalcin     /* Monitor convergence */
24685385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
24785385478SLisandro Dalcin     snes->iter = i+1;
24885385478SLisandro Dalcin     snes->norm = fnorm;
24985385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
25085385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
2517a03ce2fSLisandro Dalcin     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
25285385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
25385385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
254e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2553f149594SLisandro Dalcin     if (snes->reason) break;
25616c95cacSBarry Smith   }
25752392280SLois Curfman McInnes   if (i == maxits) {
258ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
25985385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
26052392280SLois Curfman McInnes   }
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
2625e76c431SBarry Smith }
26304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
26404d965bbSLois Curfman McInnes /*
2654b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2664b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
26704d965bbSLois Curfman McInnes 
26804d965bbSLois Curfman McInnes    Input Parameter:
26904d965bbSLois Curfman McInnes .  snes - the SNES context
27004d965bbSLois Curfman McInnes .  x - the solution vector
27104d965bbSLois Curfman McInnes 
27204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
27304d965bbSLois Curfman McInnes 
27404d965bbSLois Curfman McInnes    Notes:
2754b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
27604d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
27704d965bbSLois Curfman McInnes    the call to SNESSolve().
27804d965bbSLois Curfman McInnes  */
2794a2ae208SSatish Balay #undef __FUNCT__
2804b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
281dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2825e76c431SBarry Smith {
283dfbe8321SBarry Smith   PetscErrorCode ierr;
2843a40ed3dSBarry Smith 
2853a40ed3dSBarry Smith   PetscFunctionBegin;
28658c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
2873a40ed3dSBarry Smith   PetscFunctionReturn(0);
2885e76c431SBarry Smith }
28904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2906b8b9a38SLisandro Dalcin 
2916b8b9a38SLisandro Dalcin #undef __FUNCT__
2926b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS"
2936b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes)
2946b8b9a38SLisandro Dalcin {
2956b8b9a38SLisandro Dalcin   PetscErrorCode ierr;
2966b8b9a38SLisandro Dalcin 
2976b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2986b8b9a38SLisandro Dalcin   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2996b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
3006b8b9a38SLisandro Dalcin }
3016b8b9a38SLisandro Dalcin 
30204d965bbSLois Curfman McInnes /*
3034b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
3044b27c08aSLois Curfman McInnes    with SNESCreate_LS().
30504d965bbSLois Curfman McInnes 
30604d965bbSLois Curfman McInnes    Input Parameter:
30704d965bbSLois Curfman McInnes .  snes - the SNES context
30804d965bbSLois Curfman McInnes 
309de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
31004d965bbSLois Curfman McInnes  */
3114a2ae208SSatish Balay #undef __FUNCT__
3124b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
313dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3145e76c431SBarry Smith {
31594298653SBarry Smith   SNES_LS        *ls = (SNES_LS*) snes->data;
316dfbe8321SBarry Smith   PetscErrorCode ierr;
3173a40ed3dSBarry Smith 
3183a40ed3dSBarry Smith   PetscFunctionBegin;
3196b8b9a38SLisandro Dalcin   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
320649052a6SBarry Smith   ierr = PetscViewerDestroy(&ls->monitor);CHKERRQ(ierr);
3215c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
322557d3f75SLisandro Dalcin 
323557d3f75SLisandro Dalcin   /* clear composed functions */
32458c9b817SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
32558c9b817SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","",PETSC_NULL);CHKERRQ(ierr);
326557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
327557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
328557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
329557d3f75SLisandro Dalcin 
3303a40ed3dSBarry Smith   PetscFunctionReturn(0);
3315e76c431SBarry Smith }
33204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3334a2ae208SSatish Balay #undef __FUNCT__
33417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
33582bf6240SBarry Smith 
3364b828684SBarry Smith /*@C
33717bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3385e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3395e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3405e76c431SBarry Smith 
3413f9fe445SBarry Smith    Logically Collective on SNES and Vec
342c7afd0dbSLois Curfman McInnes 
3435e76c431SBarry Smith    Input Parameters:
344c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
345af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3465e76c431SBarry Smith .  x - current iterate
3475e76c431SBarry Smith .  f - residual evaluated at x
3483c632250SBarry Smith .  y - search direction
349dc357ad2SBarry Smith .  fnorm - 2-norm of f
350dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
3515e76c431SBarry Smith 
352c4a48953SLois Curfman McInnes    Output Parameters:
353c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3543c632250SBarry Smith .  w - new iterate
3555e76c431SBarry Smith .  gnorm - 2-norm of g
3565e76c431SBarry Smith .  ynorm - 2-norm of search length
357a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
358fee21e36SBarry Smith 
359c4a48953SLois Curfman McInnes    Options Database Key:
36017bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
361c4a48953SLois Curfman McInnes 
36236851e7fSLois Curfman McInnes    Level: advanced
36336851e7fSLois Curfman McInnes 
36428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
36528ae5a21SLois Curfman McInnes 
36617bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
36717bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3685e76c431SBarry Smith @*/
3697087cfbeSBarry Smith PetscErrorCode  SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
3705e76c431SBarry Smith {
371dfbe8321SBarry Smith   PetscErrorCode ierr;
3724b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
373ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3745334005bSBarry Smith 
3753a40ed3dSBarry Smith   PetscFunctionBegin;
376a7cc72afSBarry Smith   *flag = PETSC_TRUE;
377d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
378a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
379*55d4ea13SBarry Smith   ierr = VecWAXPY(w,-neP->damping,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3803c632250SBarry Smith   if (neP->postcheckstep) {
3813c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3828f99978dSLois Curfman McInnes   }
3833c632250SBarry Smith   if (changed_y) {
384*55d4ea13SBarry Smith     ierr = VecWAXPY(w,-neP->damping,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3853c632250SBarry Smith   }
3864936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
3874936397dSBarry Smith   if (!snes->domainerror) {
388a3b38805SLois Curfman McInnes     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
38962cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3904936397dSBarry Smith   }
391d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3923a40ed3dSBarry Smith   PetscFunctionReturn(0);
3935e76c431SBarry Smith }
39404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3954a2ae208SSatish Balay #undef __FUNCT__
39617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
39782bf6240SBarry Smith 
39829e0b56fSBarry Smith /*@C
39917bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
40029e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
40129e0b56fSBarry Smith    even compute the norm of the function or search direction; this
40229e0b56fSBarry Smith    is intended only when you know the full step is fine and are
40329e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
404c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
405c7afd0dbSLois Curfman McInnes 
4063f9fe445SBarry Smith    Logically Collective on SNES and Vec
40729e0b56fSBarry Smith 
40829e0b56fSBarry Smith    Input Parameters:
409c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
410af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
41129e0b56fSBarry Smith .  x - current iterate
41229e0b56fSBarry Smith .  f - residual evaluated at x
4133c632250SBarry Smith .  y - search direction
41429e0b56fSBarry Smith .  w - work vector
415dc357ad2SBarry Smith .  fnorm - 2-norm of f
416dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
41729e0b56fSBarry Smith 
41829e0b56fSBarry Smith    Output Parameters:
419c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4203c632250SBarry Smith .  w - new iterate
42129e0b56fSBarry Smith .  gnorm - not changed
42229e0b56fSBarry Smith .  ynorm - not changed
423a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
424fee21e36SBarry Smith 
42529e0b56fSBarry Smith    Options Database Key:
42617bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
42729e0b56fSBarry Smith 
4288cbba510SBarry Smith    Notes:
42917bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
430ea56f5baSLois Curfman McInnes    either the options
431ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
432ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
433329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
434329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
435329f5518SBarry Smith 
436329f5518SBarry Smith    During the final iteration this will not evaluate the function at
437329f5518SBarry Smith    the solution point. This is to save a function evaluation while
438329f5518SBarry Smith    using pseudo-timestepping.
4398cbba510SBarry Smith 
440ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
441a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
442ea56f5baSLois Curfman McInnes    correct, since they are not computed.
443ea56f5baSLois Curfman McInnes 
444ea56f5baSLois Curfman McInnes    Level: advanced
4458cbba510SBarry Smith 
44629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
44729e0b56fSBarry Smith 
44817bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
44917bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
45029e0b56fSBarry Smith @*/
4517087cfbeSBarry Smith PetscErrorCode  SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
45229e0b56fSBarry Smith {
453dfbe8321SBarry Smith   PetscErrorCode ierr;
4544b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
455ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
45629e0b56fSBarry Smith 
4573a40ed3dSBarry Smith   PetscFunctionBegin;
458a7cc72afSBarry Smith   *flag = PETSC_TRUE;
459d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
460*55d4ea13SBarry Smith   ierr = VecWAXPY(w,-neP->damping,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4613c632250SBarry Smith   if (neP->postcheckstep) {
4623c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4633c632250SBarry Smith   }
4643c632250SBarry Smith   if (changed_y) {
465*55d4ea13SBarry Smith     ierr = VecWAXPY(w,-neP->damping,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4668f99978dSLois Curfman McInnes   }
467329f5518SBarry Smith 
468329f5518SBarry Smith   /* don't evaluate function the last time through */
469329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4704936397dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
471329f5518SBarry Smith   }
472d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4733a40ed3dSBarry Smith   PetscFunctionReturn(0);
47429e0b56fSBarry Smith }
47504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4764a2ae208SSatish Balay #undef __FUNCT__
47717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4784b828684SBarry Smith /*@C
47917bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4805e76c431SBarry Smith 
481c7afd0dbSLois Curfman McInnes    Collective on SNES
482c7afd0dbSLois Curfman McInnes 
4835e76c431SBarry Smith    Input Parameters:
484c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
485af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4865e76c431SBarry Smith .  x - current iterate
4875e76c431SBarry Smith .  f - residual evaluated at x
4883c632250SBarry Smith .  y - search direction
4895e76c431SBarry Smith .  w - work vector
490dc357ad2SBarry Smith .  fnorm - 2-norm of f
491dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
4925e76c431SBarry Smith 
493393d2d9aSLois Curfman McInnes    Output Parameters:
494c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4953c632250SBarry Smith .  w - new iterate
4965e76c431SBarry Smith .  gnorm - 2-norm of g
4975e76c431SBarry Smith .  ynorm - 2-norm of search length
498a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
499fee21e36SBarry Smith 
500c4a48953SLois Curfman McInnes    Options Database Key:
501e106eecfSBarry Smith +  -snes_ls cubic - Activates SNESLineSearchCubic()
502e106eecfSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
503e106eecfSBarry 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)
504e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
505e106eecfSBarry Smith 
506c4a48953SLois Curfman McInnes 
5075e76c431SBarry Smith    Notes:
5085e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
5095e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
5105e76c431SBarry Smith 
51136851e7fSLois Curfman McInnes    Level: advanced
51236851e7fSLois Curfman McInnes 
51328ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
51428ae5a21SLois Curfman McInnes 
51517bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5165e76c431SBarry Smith @*/
5177087cfbeSBarry Smith PetscErrorCode  SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
5185e76c431SBarry Smith {
519406556e6SLois Curfman McInnes   /*
520406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
521406556e6SLois Curfman McInnes      minimization problem:
522406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
523406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
524406556e6SLois Curfman McInnes    */
525406556e6SLois Curfman McInnes 
526e106eecfSBarry Smith   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
527e106eecfSBarry Smith   PetscReal      minlambda,lambda,lambdatemp;
528aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
529f5b6597dSBarry Smith   PetscScalar    cinitslope;
5306b5873e3SBarry Smith #endif
531dfbe8321SBarry Smith   PetscErrorCode ierr;
532a7cc72afSBarry Smith   PetscInt       count;
5334b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
534ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
53594298653SBarry Smith   MPI_Comm       comm;
5365e76c431SBarry Smith 
5373a40ed3dSBarry Smith   PetscFunctionBegin;
53894298653SBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
539d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
540a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5415e76c431SBarry Smith 
542cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
54375567043SBarry Smith   if (*ynorm == 0.0) {
54494298653SBarry Smith     if (neP->monitor) {
545649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
546649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
547649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
54894298653SBarry Smith     }
549a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5503c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
551a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
552a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
553ad922ac9SBarry Smith     goto theend1;
55494a9d846SBarry Smith   }
555e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
55694298653SBarry Smith     if (neP->monitor) {
557649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5588f1a2a5eSBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Scaling step by %14.12e old ynorm %14.12e\n",(double)(neP->maxstep/(*ynorm)),(double)(*ynorm));CHKERRQ(ierr);
559649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
56094298653SBarry Smith     }
561e106eecfSBarry Smith     ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
562e106eecfSBarry Smith     *ynorm = neP->maxstep;
5635e76c431SBarry Smith   }
564ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
565e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
566a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
567aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
568a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
569329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5705e42470aSBarry Smith #else
571a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5725e42470aSBarry Smith #endif
5735e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5745e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5755e76c431SBarry Smith 
576e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
57743e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
578ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
57943e71028SBarry Smith     *flag = PETSC_FALSE;
58043e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
58143e71028SBarry Smith     goto theend1;
58243e71028SBarry Smith   }
5834936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5844936397dSBarry Smith   if (snes->domainerror) {
5854936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5864936397dSBarry Smith     PetscFunctionReturn(0);
58719717074SBarry Smith   }
588313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
58962cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5908f1a2a5eSBarry Smith   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
591e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
59294298653SBarry Smith     if (neP->monitor) {
593649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
5948f1a2a5eSBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
595649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
59694298653SBarry Smith     }
597ad922ac9SBarry Smith     goto theend1;
5985e76c431SBarry Smith   }
5995e76c431SBarry Smith 
6005e76c431SBarry Smith   /* Fit points with quadratic */
601313b5538SBarry Smith   lambda     = 1.0;
602a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
6035e76c431SBarry Smith   lambdaprev = lambda;
6045e76c431SBarry Smith   gnormprev  = *gnorm;
605329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
606ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
607ddd12767SBarry Smith   else                         lambda = lambdatemp;
608275d25c3SBarry Smith 
609e68848bdSBarry Smith   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
61043e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
611ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
61243e71028SBarry Smith     *flag = PETSC_FALSE;
61343e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
61443e71028SBarry Smith     goto theend1;
61543e71028SBarry Smith   }
616f5b6597dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6174936397dSBarry Smith   if (snes->domainerror) {
6184936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6194936397dSBarry Smith     PetscFunctionReturn(0);
62019717074SBarry Smith   }
621cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
62262cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
62394298653SBarry Smith   if (neP->monitor) {
624649052a6SBarry Smith     ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
6258f1a2a5eSBarry Smith     ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr);
626649052a6SBarry Smith     ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
62794298653SBarry Smith   }
628e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
62994298653SBarry Smith     if (neP->monitor) {
630649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
631649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
632649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
63394298653SBarry Smith     }
634ad922ac9SBarry Smith     goto theend1;
6355e76c431SBarry Smith   }
6365e76c431SBarry Smith 
6375e76c431SBarry Smith   /* Fit points with cubic */
6385e76c431SBarry Smith   count = 1;
639bb9195b6SBarry Smith   while (PETSC_TRUE) {
640dc357ad2SBarry Smith     if (lambda <= minlambda) {
64194298653SBarry Smith       if (neP->monitor) {
642649052a6SBarry Smith         ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
643649052a6SBarry Smith 	ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
644649052a6SBarry Smith 	ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, minlambda=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)minlambda,(double)lambda,(double)initslope);CHKERRQ(ierr);
645649052a6SBarry Smith         ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
64694298653SBarry Smith       }
64743e71028SBarry Smith       *flag = PETSC_FALSE;
64843e71028SBarry Smith       break;
6495e76c431SBarry Smith     }
6505d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6515d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
652ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6532b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6545e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6555e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6565e76c431SBarry Smith     if (a == 0.0) {
6575e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6585e76c431SBarry Smith     } else {
6595e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6605e76c431SBarry Smith     }
6615e76c431SBarry Smith     lambdaprev = lambda;
6625e76c431SBarry Smith     gnormprev  = *gnorm;
663329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
664329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6655e76c431SBarry Smith     else                         lambda     = lambdatemp;
666e68848bdSBarry Smith     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
66743e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
668ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
669ae15b995SBarry 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);
670ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
67143e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
672ed98166eSMatthew Knepley       break;
673ed98166eSMatthew Knepley     }
674f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6754936397dSBarry Smith     if (snes->domainerror) {
6764936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6774936397dSBarry Smith       PetscFunctionReturn(0);
67819717074SBarry Smith     }
679cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
68062cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
681e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */
68294298653SBarry Smith       if (neP->monitor) {
6838f1a2a5eSBarry Smith 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
68494298653SBarry Smith       }
68543e71028SBarry Smith       break;
6862b022350SLois Curfman McInnes     } else {
68794298653SBarry Smith       if (neP->monitor) {
6888f1a2a5eSBarry Smith         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %12.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
68994298653SBarry Smith       }
6905e76c431SBarry Smith     }
6915e76c431SBarry Smith     count++;
6925e76c431SBarry Smith   }
693ad922ac9SBarry Smith   theend1:
6948f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6953c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6963c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6973c632250SBarry Smith     if (changed_y) {
69879f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6993c632250SBarry Smith     }
7003c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
701f5b6597dSBarry Smith       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
7024936397dSBarry Smith       if (snes->domainerror) {
7034936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7044936397dSBarry Smith         PetscFunctionReturn(0);
70519717074SBarry Smith       }
7068f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
70762cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
708a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7098f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
710a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7118f99978dSLois Curfman McInnes     }
7128f99978dSLois Curfman McInnes   }
713d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7143a40ed3dSBarry Smith   PetscFunctionReturn(0);
7155e76c431SBarry Smith }
71604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7174a2ae208SSatish Balay #undef __FUNCT__
71817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
7194b828684SBarry Smith /*@C
72017bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
7215e76c431SBarry Smith 
722c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
723c7afd0dbSLois Curfman McInnes 
7245e42470aSBarry Smith    Input Parameters:
725c7afd0dbSLois Curfman McInnes +  snes - the SNES context
726af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
7275e42470aSBarry Smith .  x - current iterate
7285e42470aSBarry Smith .  f - residual evaluated at x
7293c632250SBarry Smith .  y - search direction
7305e42470aSBarry Smith .  w - work vector
731dc357ad2SBarry Smith .  fnorm - 2-norm of f
732dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
7335e42470aSBarry Smith 
734c4a48953SLois Curfman McInnes    Output Parameters:
7357f3332b4SBarry Smith +  g - residual evaluated at new iterate w
736e106eecfSBarry Smith .  w - new iterate (x + lambda*y)
7375e42470aSBarry Smith .  gnorm - 2-norm of g
7385e42470aSBarry Smith .  ynorm - 2-norm of search length
739a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
740fee21e36SBarry Smith 
741ce9499c7SBarry Smith    Options Database Keys:
742ce9499c7SBarry Smith +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
743ce9499c7SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
744e106eecfSBarry 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)
745e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
746e106eecfSBarry Smith 
7475e42470aSBarry Smith    Notes:
74817bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7495e42470aSBarry Smith 
75036851e7fSLois Curfman McInnes    Level: advanced
75136851e7fSLois Curfman McInnes 
752f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
753f59ffdedSLois Curfman McInnes 
75417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7555e42470aSBarry Smith @*/
7567087cfbeSBarry Smith PetscErrorCode  SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
7575e76c431SBarry Smith {
758406556e6SLois Curfman McInnes   /*
759406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
760406556e6SLois Curfman McInnes      minimization problem:
761406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
762406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
763406556e6SLois Curfman McInnes    */
764e106eecfSBarry Smith   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
765aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
766f5b6597dSBarry Smith   PetscScalar    cinitslope;
7676b5873e3SBarry Smith #endif
768dfbe8321SBarry Smith   PetscErrorCode ierr;
769a7cc72afSBarry Smith   PetscInt       count;
7704b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
771ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7725e76c431SBarry Smith 
7733a40ed3dSBarry Smith   PetscFunctionBegin;
774d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
775a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7765e76c431SBarry Smith 
7773505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
77863b9fa5eSBarry Smith   if (*ynorm == 0.0) {
77994298653SBarry Smith     if (neP->monitor) {
780649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
781649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
782649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
78394298653SBarry Smith     }
784b37302c6SLois Curfman McInnes     *gnorm = fnorm;
785e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
786b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
787a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
788ad922ac9SBarry Smith     goto theend2;
78994a9d846SBarry Smith   }
790e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
791e106eecfSBarry Smith     ierr   = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
792e106eecfSBarry Smith     *ynorm = neP->maxstep;
7935e76c431SBarry Smith   }
794ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
795e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
796a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
797aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
798a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
799329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
8005e42470aSBarry Smith #else
801a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
8025e42470aSBarry Smith #endif
8035e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
8045e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
8058f1a2a5eSBarry Smith   ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr);
8065e42470aSBarry Smith 
807e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
80843e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
809ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
81043e71028SBarry Smith     *flag = PETSC_FALSE;
81143e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
81243e71028SBarry Smith     goto theend2;
81343e71028SBarry Smith   }
8144936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8154936397dSBarry Smith   if (snes->domainerror) {
8164936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8174936397dSBarry Smith     PetscFunctionReturn(0);
81819717074SBarry Smith   }
819cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82062cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
821e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
82294298653SBarry Smith     if (neP->monitor) {
823649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
824649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
825649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
82694298653SBarry Smith     }
827ad922ac9SBarry Smith     goto theend2;
8285e42470aSBarry Smith   }
8295e42470aSBarry Smith 
8305e42470aSBarry Smith   /* Fit points with quadratic */
831313b5538SBarry Smith   lambda = 1.0;
8325e42470aSBarry Smith   count = 1;
8335ca10a37SBarry Smith   while (PETSC_TRUE) {
8345e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
83594298653SBarry Smith       if (neP->monitor) {
836649052a6SBarry Smith         ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
837649052a6SBarry Smith         ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
8388f1a2a5eSBarry Smith         ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: fnorm=%14.12e, gnorm=%14.12e, ynorm=%14.12e, lambda=%14.12e, initial slope=%14.12e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
839649052a6SBarry Smith         ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
84094298653SBarry Smith       }
841e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
84243e71028SBarry Smith       *flag = PETSC_FALSE;
84343e71028SBarry Smith       break;
8445e42470aSBarry Smith     }
845a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
846329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
847329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
848329f5518SBarry Smith     else                         lambda     = lambdatemp;
849275d25c3SBarry Smith 
850e68848bdSBarry Smith     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
85143e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
852ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
853d9822059SBarry Smith       ierr  = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",(double)fnorm,(double)*gnorm,(double)*ynorm,(double)lambda,(double)initslope);CHKERRQ(ierr);
854ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
85543e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
856ed98166eSMatthew Knepley       break;
857ed98166eSMatthew Knepley     }
858f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8594936397dSBarry Smith     if (snes->domainerror) {
8604936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8614936397dSBarry Smith       PetscFunctionReturn(0);
86219717074SBarry Smith     }
863cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
86462cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
865e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
86694298653SBarry Smith       if (neP->monitor) {
867649052a6SBarry Smith         ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
8688f1a2a5eSBarry Smith         ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr);
869649052a6SBarry Smith         ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
87094298653SBarry Smith       }
87106259719SBarry Smith       break;
8725e42470aSBarry Smith     }
8735e42470aSBarry Smith     count++;
8745e42470aSBarry Smith   }
875ad922ac9SBarry Smith   theend2:
8768f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8773c632250SBarry Smith   if (neP->postcheckstep) {
8783c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8793c632250SBarry Smith     if (changed_y) {
88079f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8813c632250SBarry Smith     }
8823c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8833c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
8844936397dSBarry Smith       if (snes->domainerror) {
8854936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8864936397dSBarry Smith         PetscFunctionReturn(0);
88719717074SBarry Smith       }
8888f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
889a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
8908f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
891a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
89262cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8938f99978dSLois Curfman McInnes     }
8948f99978dSLois Curfman McInnes   }
895d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8963a40ed3dSBarry Smith   PetscFunctionReturn(0);
8975e76c431SBarry Smith }
8982343ba6eSBarry Smith 
89904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9004a2ae208SSatish Balay #undef __FUNCT__
90117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
902c9e6a524SLois Curfman McInnes /*@C
90317bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
9044b27c08aSLois Curfman McInnes    by the method SNESLS.
9055e76c431SBarry Smith 
9065e76c431SBarry Smith    Input Parameters:
907c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
908af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
909c7afd0dbSLois Curfman McInnes -  func - pointer to int function
9105e76c431SBarry Smith 
9113f9fe445SBarry Smith    Logically Collective on SNES
912fee21e36SBarry Smith 
913c4a48953SLois Curfman McInnes    Available Routines:
91417bae607SBarry Smith +  SNESLineSearchCubic() - default line search
91517bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
91617bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
91717bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
9185e76c431SBarry Smith 
919c4a48953SLois Curfman McInnes     Options Database Keys:
9204b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
921*55d4ea13SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
922e106eecfSBarry 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)
923e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
924c4a48953SLois Curfman McInnes 
9255e76c431SBarry Smith    Calling sequence of func:
926bd208895SLois Curfman McInnes .vb
927ace3abfcSBarry Smith    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal xnorm,PetscReal *ynorm,PetscReal *gnorm,PetscBool  *flag)
928bd208895SLois Curfman McInnes .ve
9295e76c431SBarry Smith 
9305e76c431SBarry Smith     Input parameters for func:
931c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
932af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
9335e76c431SBarry Smith .   x - current iterate
9345e76c431SBarry Smith .   f - residual evaluated at x
9353c632250SBarry Smith .   y - search direction
936c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9375e76c431SBarry Smith 
9385e76c431SBarry Smith     Output parameters for func:
939c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9403c632250SBarry Smith .   w - new iterate
9415e76c431SBarry Smith .   gnorm - 2-norm of g
9425e76c431SBarry Smith .   ynorm - 2-norm of search length
943a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
944f59ffdedSLois Curfman McInnes 
94536851e7fSLois Curfman McInnes     Level: advanced
94636851e7fSLois Curfman McInnes 
947f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
948f59ffdedSLois Curfman McInnes 
94917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
95017bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9515e76c431SBarry Smith @*/
9527087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx)
9535e76c431SBarry Smith {
9544ac538c5SBarry Smith   PetscErrorCode ierr;
95582bf6240SBarry Smith 
9563a40ed3dSBarry Smith   PetscFunctionBegin;
9574ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSet_C",(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void*),(snes,func,lsctx));CHKERRQ(ierr);
9583a40ed3dSBarry Smith   PetscFunctionReturn(0);
9595e76c431SBarry Smith }
9608e019c35SBarry Smith 
961ace3abfcSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *); /* force argument to next function to not be extern C*/
96204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
963fb2e594dSBarry Smith EXTERN_C_BEGIN
9644a2ae208SSatish Balay #undef __FUNCT__
96517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
9667087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
96782bf6240SBarry Smith {
96882bf6240SBarry Smith   PetscFunctionBegin;
9694b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9704b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
97182bf6240SBarry Smith   PetscFunctionReturn(0);
97282bf6240SBarry Smith }
973fb2e594dSBarry Smith EXTERN_C_END
97404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9754a2ae208SSatish Balay #undef __FUNCT__
9763c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
977c8dd0c0dSLois Curfman McInnes /*@C
9783c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9794b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
980c8dd0c0dSLois Curfman McInnes 
981c8dd0c0dSLois Curfman McInnes    Input Parameters:
982c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9833c632250SBarry Smith .  func - pointer to function
984c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
985c8dd0c0dSLois Curfman McInnes 
9863f9fe445SBarry Smith    Logically Collective on SNES
987c8dd0c0dSLois Curfman McInnes 
988c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
989c8dd0c0dSLois Curfman McInnes .vb
990ace3abfcSBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool  *changed_y,PetscBool  *changed_w)
991c8dd0c0dSLois Curfman McInnes .ve
992b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
993b1ae27eaSLois Curfman McInnes    on failure.
994c8dd0c0dSLois Curfman McInnes 
995c8dd0c0dSLois Curfman McInnes    Input parameters for func:
996c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
997c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9983c632250SBarry Smith .  x - previous iterate
9993c632250SBarry Smith .  y - new search direction and length
10003c632250SBarry Smith -  w - current candidate iterate
1001c8dd0c0dSLois Curfman McInnes 
1002c8dd0c0dSLois Curfman McInnes    Output parameters for func:
10033c632250SBarry Smith +  y - search direction (possibly changed)
10043c632250SBarry Smith .  w - current iterate (possibly modified)
10053c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
10063c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
1007c8dd0c0dSLois Curfman McInnes 
1008c8dd0c0dSLois Curfman McInnes    Level: advanced
1009c8dd0c0dSLois Curfman McInnes 
10109e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10119e247f21SBarry Smith 
10123c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
10133c632250SBarry Smith 
1014481b4421SBarry Smith    On input w = x - y
10153c632250SBarry Smith 
101617bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
1017b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
1018ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
10198f99978dSLois Curfman McInnes 
102017bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
1021ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
1022ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
1023ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
10249e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
1025b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
10268f99978dSLois Curfman McInnes 
1027c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
1028c8dd0c0dSLois Curfman McInnes 
102917bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1030c8dd0c0dSLois Curfman McInnes @*/
10317087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx)
1032c8dd0c0dSLois Curfman McInnes {
10334ac538c5SBarry Smith   PetscErrorCode ierr;
1034c8dd0c0dSLois Curfman McInnes 
1035c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10364ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr);
1037c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1038c8dd0c0dSLois Curfman McInnes }
103994298653SBarry Smith 
10409c3ca13aSBarry Smith #undef __FUNCT__
10419c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10429c3ca13aSBarry Smith /*@C
10439c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10447e4bb74cSBarry Smith          before the line search is called.
10459c3ca13aSBarry Smith 
10469c3ca13aSBarry Smith    Input Parameters:
10479c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10489c3ca13aSBarry Smith .  func - pointer to function
10499c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10509c3ca13aSBarry Smith 
10513f9fe445SBarry Smith    Logically Collective on SNES
10529c3ca13aSBarry Smith 
10539c3ca13aSBarry Smith    Calling sequence of func:
10549c3ca13aSBarry Smith .vb
1055ace3abfcSBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool  *changed_y)
10569c3ca13aSBarry Smith .ve
10579c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10589c3ca13aSBarry Smith    on failure.
10599c3ca13aSBarry Smith 
10609c3ca13aSBarry Smith    Input parameters for func:
10619c3ca13aSBarry Smith +  snes - nonlinear context
10629c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10639c3ca13aSBarry Smith .  x - previous iterate
10649c3ca13aSBarry Smith -  y - new search direction and length
10659c3ca13aSBarry Smith 
10669c3ca13aSBarry Smith    Output parameters for func:
10679c3ca13aSBarry Smith +  y - search direction (possibly changed)
10689c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10699c3ca13aSBarry Smith 
10709c3ca13aSBarry Smith    Level: advanced
10719c3ca13aSBarry Smith 
10729c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10739c3ca13aSBarry Smith 
10749c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10759c3ca13aSBarry Smith 
10767e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10779c3ca13aSBarry Smith @*/
10787087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx)
10799c3ca13aSBarry Smith {
10804ac538c5SBarry Smith   PetscErrorCode ierr;
10819c3ca13aSBarry Smith 
10829c3ca13aSBarry Smith   PetscFunctionBegin;
10834ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr);
10849c3ca13aSBarry Smith   PetscFunctionReturn(0);
10859c3ca13aSBarry Smith }
10869c3ca13aSBarry Smith 
108794298653SBarry Smith #undef __FUNCT__
108894298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor"
108994298653SBarry Smith /*@C
109094298653SBarry Smith    SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search
109194298653SBarry Smith 
109294298653SBarry Smith    Input Parameters:
109394298653SBarry Smith +  snes - nonlinear context obtained from SNESCreate()
109494298653SBarry Smith -  flg - PETSC_TRUE to monitor the line search
109594298653SBarry Smith 
10963f9fe445SBarry Smith    Logically Collective on SNES
109794298653SBarry Smith 
109894298653SBarry Smith    Options Database:
109994298653SBarry Smith .   -snes_ls_monitor
110094298653SBarry Smith 
110194298653SBarry Smith    Level: intermediate
110294298653SBarry Smith 
110394298653SBarry Smith 
110494298653SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
110594298653SBarry Smith @*/
11067087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor(SNES snes,PetscBool  flg)
110794298653SBarry Smith {
11084ac538c5SBarry Smith   PetscErrorCode ierr;
110994298653SBarry Smith 
111094298653SBarry Smith   PetscFunctionBegin;
11114ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr);
111294298653SBarry Smith   PetscFunctionReturn(0);
111394298653SBarry Smith }
111494298653SBarry Smith 
1115c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1116ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/
1117ace3abfcSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *);                 /* force argument to next function to not be extern C*/
1118c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
11194a2ae208SSatish Balay #undef __FUNCT__
11203c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
11217087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1122c8dd0c0dSLois Curfman McInnes {
1123c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
11243c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
11253c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1126c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1127c8dd0c0dSLois Curfman McInnes }
1128c8dd0c0dSLois Curfman McInnes EXTERN_C_END
11299c3ca13aSBarry Smith 
11309c3ca13aSBarry Smith EXTERN_C_BEGIN
11319c3ca13aSBarry Smith #undef __FUNCT__
11329c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
11337087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
11349c3ca13aSBarry Smith {
11359c3ca13aSBarry Smith   PetscFunctionBegin;
11369c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
11379c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
11389c3ca13aSBarry Smith   PetscFunctionReturn(0);
11399c3ca13aSBarry Smith }
11409c3ca13aSBarry Smith EXTERN_C_END
1141329e7e40SMatthew Knepley 
114294298653SBarry Smith EXTERN_C_BEGIN
114394298653SBarry Smith #undef __FUNCT__
114494298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor_LS"
11457087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor_LS(SNES snes,PetscBool  flg)
114694298653SBarry Smith {
114794298653SBarry Smith   SNES_LS        *ls = (SNES_LS*)snes->data;
114894298653SBarry Smith   PetscErrorCode ierr;
114994298653SBarry Smith 
115094298653SBarry Smith   PetscFunctionBegin;
115194298653SBarry Smith   if (flg && !ls->monitor) {
1152649052a6SBarry Smith     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&ls->monitor);CHKERRQ(ierr);
115394298653SBarry Smith   } else if (!flg && ls->monitor) {
1154649052a6SBarry Smith     ierr = PetscViewerDestroy(&ls->monitor);CHKERRQ(ierr);
115594298653SBarry Smith   }
115694298653SBarry Smith   PetscFunctionReturn(0);
115794298653SBarry Smith }
115894298653SBarry Smith EXTERN_C_END
115994298653SBarry Smith 
1160329e7e40SMatthew Knepley /*
11614b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
116204d965bbSLois Curfman McInnes 
116304d965bbSLois Curfman McInnes    Input Parameters:
116404d965bbSLois Curfman McInnes .  SNES - the SNES context
116504d965bbSLois Curfman McInnes .  viewer - visualization context
116604d965bbSLois Curfman McInnes 
116704d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
116804d965bbSLois Curfman McInnes */
11694a2ae208SSatish Balay #undef __FUNCT__
11704b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11716849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1172a935fc98SLois Curfman McInnes {
11734b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11742fc52814SBarry Smith   const char     *cstr;
1175dfbe8321SBarry Smith   PetscErrorCode ierr;
1176ace3abfcSBarry Smith   PetscBool      iascii;
1177a935fc98SLois Curfman McInnes 
11783a40ed3dSBarry Smith   PetscFunctionBegin;
11792692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
118032077d6dSBarry Smith   if (iascii) {
118117bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
118221167be9SJed Brown     if (ls->LineSearch == SNESLineSearchNoNorms)        cstr = "SNESLineSearchNoNorms";
118317bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
118417bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
118519bcc07fSBarry Smith     else                                                cstr = "unknown";
1186b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
11878f1a2a5eSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)ls->alpha,(double)ls->maxstep,(double)ls->minlambda);CHKERRQ(ierr);
1188*55d4ea13SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  damping factor=%14.12e\n",(double)ls->damping);CHKERRQ(ierr);
118919bcc07fSBarry Smith   }
11903a40ed3dSBarry Smith   PetscFunctionReturn(0);
1191a935fc98SLois Curfman McInnes }
119204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
119304d965bbSLois Curfman McInnes /*
11944b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
119504d965bbSLois Curfman McInnes 
119604d965bbSLois Curfman McInnes    Input Parameter:
119704d965bbSLois Curfman McInnes .  snes - the SNES context
119804d965bbSLois Curfman McInnes 
119904d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
120004d965bbSLois Curfman McInnes */
12014a2ae208SSatish Balay #undef __FUNCT__
12024b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
12036849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
12045e42470aSBarry Smith {
12054b27c08aSLois Curfman McInnes   SNES_LS            *ls = (SNES_LS *)snes->data;
1206dfbe8321SBarry Smith   PetscErrorCode     ierr;
1207490ca623SBarry Smith   SNESLineSearchType indx;
1208ace3abfcSBarry Smith   PetscBool          flg,set;
12095e42470aSBarry Smith 
12103a40ed3dSBarry Smith   PetscFunctionBegin;
1211b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
12124b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
12134b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1214e106eecfSBarry Smith     ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr);
1215*55d4ea13SBarry Smith     ierr = PetscOptionsReal("-snes_ls_damping","Damping parameter","SNES",ls->damping,&ls->damping,&flg);CHKERRQ(ierr);
1216acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
121794298653SBarry Smith     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
1218186905e3SBarry Smith 
1219490ca623SBarry Smith     ierr = PetscOptionsEnum("-snes_ls","Line search used","SNESLineSearchSet",SNESLineSearchTypes,(PetscEnum)SNES_LS_CUBIC,(PetscEnum*)&indx,&flg);CHKERRQ(ierr);
122025cdf11fSBarry Smith     if (flg) {
122122e36657SBarry Smith       switch (indx) {
1222490ca623SBarry Smith       case SNES_LS_BASIC:
122317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1224b49fd9e1SBarry Smith         break;
1225490ca623SBarry Smith       case SNES_LS_BASIC_NONORMS:
122617bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1227b49fd9e1SBarry Smith         break;
1228490ca623SBarry Smith       case SNES_LS_QUADRATIC:
122917bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1230b49fd9e1SBarry Smith         break;
1231490ca623SBarry Smith       case SNES_LS_CUBIC:
123217bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1233b49fd9e1SBarry Smith         break;
12345e42470aSBarry Smith       }
12355e42470aSBarry Smith     }
1236b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
12373a40ed3dSBarry Smith   PetscFunctionReturn(0);
12385e42470aSBarry Smith }
12394827ddcaSBarry Smith 
12404827ddcaSBarry Smith EXTERN_C_BEGIN
12414827ddcaSBarry Smith extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
12424827ddcaSBarry Smith EXTERN_C_END
12434827ddcaSBarry Smith 
124404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1245ebe3b25bSBarry Smith /*MC
1246ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
124704d965bbSLois Curfman McInnes 
1248ebe3b25bSBarry Smith    Options Database:
1249ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1250*55d4ea13SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
1251e106eecfSBarry 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)
1252acbee50cSBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1253*55d4ea13SBarry Smith .   -snes_ls_monitor - print information about progress of line searches
1254*55d4ea13SBarry Smith -   -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms
1255acbee50cSBarry Smith 
125604d965bbSLois Curfman McInnes 
1257ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
125804d965bbSLois Curfman McInnes 
1259ee3001cbSBarry Smith    Level: beginner
1260ee3001cbSBarry Smith 
126117bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
126217bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1263b3dd4ab5SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
1264ebe3b25bSBarry Smith 
1265ebe3b25bSBarry Smith M*/
1266fb2e594dSBarry Smith EXTERN_C_BEGIN
12674a2ae208SSatish Balay #undef __FUNCT__
12684b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
12697087cfbeSBarry Smith PetscErrorCode  SNESCreate_LS(SNES snes)
12705e42470aSBarry Smith {
1271dfbe8321SBarry Smith   PetscErrorCode ierr;
12724b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12735e42470aSBarry Smith 
12743a40ed3dSBarry Smith   PetscFunctionBegin;
1275e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1276e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1277e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1278e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1279e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
12806b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_LS;
12815e42470aSBarry Smith 
128238f2d2fdSLisandro Dalcin   ierr                  = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
12835e42470aSBarry Smith   snes->data    	= (void*)neP;
12845e42470aSBarry Smith   neP->alpha		= 1.e-4;
12855e42470aSBarry Smith   neP->maxstep		= 1.e8;
1286e106eecfSBarry Smith   neP->minlambda        = 1.e-12;
1287*55d4ea13SBarry Smith   neP->damping          = 1.0;
128817bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1289c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12903c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12913c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12923c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12933c632250SBarry Smith   neP->precheck         = PETSC_NULL;
129482bf6240SBarry Smith 
129594298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr);
12964827ddcaSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);CHKERRQ(ierr);
129794298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
129894298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
129994298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
130082bf6240SBarry Smith 
13013a40ed3dSBarry Smith   PetscFunctionReturn(0);
13025e42470aSBarry Smith }
1303fb2e594dSBarry Smith EXTERN_C_END
130482bf6240SBarry Smith 
130582bf6240SBarry Smith 
130682bf6240SBarry Smith 
1307