xref: /petsc/src/snes/impls/ls/ls.c (revision 8f1a2a5e609bf3d0ec1c7e7bf09507165b68d274)
15e76c431SBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h>
35e42470aSBarry Smith 
48a5d9ee4SBarry Smith /*
520f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
620f52e01SBarry 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
736669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
820f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
98a5d9ee4SBarry Smith */
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal      a1;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
16ace3abfcSBarry Smith   PetscBool      hastranspose;
178a5d9ee4SBarry Smith 
188a5d9ee4SBarry Smith   PetscFunctionBegin;
198a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2036669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2136669109SBarry Smith   if (hastranspose) {
228a5d9ee4SBarry Smith     /* Compute || J^T F|| */
238a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
248a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
25*8f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
268a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2736669109SBarry Smith   } else {
2836669109SBarry Smith     Vec         work;
29ea709b57SSatish Balay     PetscScalar result;
3036669109SBarry Smith     PetscReal   wnorm;
3136669109SBarry Smith 
32abb0e124SMatthew Knepley     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
3336669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
376bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
3836669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
39*8f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
4036669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4136669109SBarry Smith   }
428a5d9ee4SBarry Smith   PetscFunctionReturn(0);
438a5d9ee4SBarry Smith }
448a5d9ee4SBarry Smith 
4574637425SBarry Smith /*
465ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4774637425SBarry Smith */
484a2ae208SSatish Balay #undef __FUNCT__
494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
501e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
5174637425SBarry Smith {
5274637425SBarry Smith   PetscReal      a1,a2;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
54ace3abfcSBarry Smith   PetscBool      hastranspose;
5574637425SBarry Smith 
5674637425SBarry Smith   PetscFunctionBegin;
5774637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5874637425SBarry Smith   if (hastranspose) {
5974637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6079f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6174637425SBarry Smith 
6274637425SBarry Smith     /* Compute || J^T W|| */
6374637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6474637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6574637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
6675567043SBarry Smith     if (a1 != 0.0) {
67*8f1a2a5eSBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
6885471664SBarry Smith     }
6974637425SBarry Smith   }
7074637425SBarry Smith   PetscFunctionReturn(0);
7174637425SBarry Smith }
7274637425SBarry Smith 
7304d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7404d965bbSLois Curfman McInnes 
7504d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7694b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7704d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
7804d965bbSLois Curfman McInnes      respectively.
7904d965bbSLois Curfman McInnes 
80fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8104d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8204d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
83fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8404d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8504d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
864b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
874b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
8804d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
8904d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9004d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9104d965bbSLois Curfman McInnes      for all nonlinear solvers.
9204d965bbSLois Curfman McInnes 
9304d965bbSLois Curfman McInnes      Another key routine is:
9404d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9504d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9604d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9704d965bbSLois Curfman McInnes      SNESSolve() if necessary.
9804d965bbSLois Curfman McInnes 
9904d965bbSLois Curfman McInnes      Additional basic routines are:
10004d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10104d965bbSLois Curfman McInnes                                       have actually been used.
10204d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
103186905e3SBarry Smith      SNESView().
10404d965bbSLois Curfman McInnes 
10504d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10604d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10704d965bbSLois Curfman McInnes      above description applies to these categories also.
10804d965bbSLois Curfman McInnes 
10904d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1105e42470aSBarry Smith /*
1114b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11204d965bbSLois Curfman McInnes    method with a line search.
1135e76c431SBarry Smith 
11404d965bbSLois Curfman McInnes    Input Parameters:
11504d965bbSLois Curfman McInnes .  snes - the SNES context
1165e76c431SBarry Smith 
11704d965bbSLois Curfman McInnes    Output Parameter:
118c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
11904d965bbSLois Curfman McInnes 
12004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1215e76c431SBarry Smith 
1225e76c431SBarry Smith    Notes:
1235e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1245e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1255e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1265e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
127393d2d9aSLois Curfman McInnes    and Schnabel.
1285e76c431SBarry Smith */
1294a2ae208SSatish Balay #undef __FUNCT__
1304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
131dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1325e76c431SBarry Smith {
1334b27c08aSLois Curfman McInnes   SNES_LS            *neP = (SNES_LS*)snes->data;
1346849ba73SBarry Smith   PetscErrorCode     ierr;
135a7cc72afSBarry Smith   PetscInt           maxits,i,lits;
136ace3abfcSBarry Smith   PetscBool          lssucceed;
137112a2221SBarry Smith   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
13885385478SLisandro Dalcin   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
13985385478SLisandro Dalcin   Vec                Y,X,F,G,W;
1403d4c4710SBarry Smith   KSPConvergedReason kspreason;
1415e76c431SBarry Smith 
1423a40ed3dSBarry Smith   PetscFunctionBegin;
1433d4c4710SBarry Smith   snes->numFailures            = 0;
1443d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
145184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
146184914b5SBarry Smith 
1475e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1485e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1505e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1515e42470aSBarry Smith   G		= snes->work[1];
1525e42470aSBarry Smith   W		= snes->work[2];
1535e76c431SBarry Smith 
1544c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1554c49b128SBarry Smith   snes->iter = 0;
15675567043SBarry Smith   snes->norm = 0.0;
1574c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15819717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1594936397dSBarry Smith   if (snes->domainerror) {
1604936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1614936397dSBarry Smith     PetscFunctionReturn(0);
1624936397dSBarry Smith   }
1632613ca53SBarry Smith   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1642613ca53SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
1652613ca53SBarry Smith   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
1662613ca53SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
16762cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1680f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1695e42470aSBarry Smith   snes->norm = fnorm;
1700f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17184c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
1727a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1733f149594SLisandro Dalcin 
1743f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
1753f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
17685385478SLisandro Dalcin   /* test convergence */
17785385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
17806ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
179d034289bSBarry Smith 
1805e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1815e76c431SBarry Smith 
182329e7e40SMatthew Knepley     /* Call general purpose update function */
183e7788613SBarry Smith     if (snes->ops->update) {
184e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
185de8cb200SMatthew Knepley     }
186329e7e40SMatthew Knepley 
187ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1885334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18994b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
19071f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
1913d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1923d4c4710SBarry Smith     if (kspreason < 0) {
1933d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
194ef998cc9SBarry 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);
1953d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
196ab81109fSSatish Balay         break;
1973d4c4710SBarry Smith       }
1983d4c4710SBarry Smith     }
1993d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
2003f149594SLisandro Dalcin     snes->linear_its += lits;
2013f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
20274637425SBarry Smith 
2039c3ca13aSBarry Smith     if (neP->precheckstep) {
204ace3abfcSBarry Smith       PetscBool  changed_y = PETSC_FALSE;
2059c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
2069c3ca13aSBarry Smith     }
2079c3ca13aSBarry Smith 
208b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
2091e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
21085471664SBarry Smith     }
21174637425SBarry Smith 
212ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
213ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
214e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
215ea4d3ed3SLois Curfman McInnes     */
21685385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
2173f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
218dc357ad2SBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,xnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
219*8f1a2a5eSBarry 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);
2204a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2214936397dSBarry Smith     if (snes->domainerror) {
2224936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2234936397dSBarry Smith       PetscFunctionReturn(0);
2244936397dSBarry Smith     }
225a7cc72afSBarry Smith     if (!lssucceed) {
22650ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
227ace3abfcSBarry Smith 	PetscBool  ismin;
228647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2291e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
2308a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2313505fcc1SBarry Smith         break;
2323505fcc1SBarry Smith       }
23350ffb88aSMatthew Knepley     }
23485385478SLisandro Dalcin     /* Update function and solution vectors */
23585385478SLisandro Dalcin     fnorm = gnorm;
23685385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
23785385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
23885385478SLisandro Dalcin     /* Monitor convergence */
23985385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
24085385478SLisandro Dalcin     snes->iter = i+1;
24185385478SLisandro Dalcin     snes->norm = fnorm;
24285385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
24385385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
2447a03ce2fSLisandro Dalcin     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
24585385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
24685385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
247e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2483f149594SLisandro Dalcin     if (snes->reason) break;
24916c95cacSBarry Smith   }
25052392280SLois Curfman McInnes   if (i == maxits) {
251ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
25285385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
25352392280SLois Curfman McInnes   }
2543a40ed3dSBarry Smith   PetscFunctionReturn(0);
2555e76c431SBarry Smith }
25604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25704d965bbSLois Curfman McInnes /*
2584b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2594b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
26004d965bbSLois Curfman McInnes 
26104d965bbSLois Curfman McInnes    Input Parameter:
26204d965bbSLois Curfman McInnes .  snes - the SNES context
26304d965bbSLois Curfman McInnes .  x - the solution vector
26404d965bbSLois Curfman McInnes 
26504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26604d965bbSLois Curfman McInnes 
26704d965bbSLois Curfman McInnes    Notes:
2684b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26904d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
27004d965bbSLois Curfman McInnes    the call to SNESSolve().
27104d965bbSLois Curfman McInnes  */
2724a2ae208SSatish Balay #undef __FUNCT__
2734b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
274dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2755e76c431SBarry Smith {
276dfbe8321SBarry Smith   PetscErrorCode ierr;
2773a40ed3dSBarry Smith 
2783a40ed3dSBarry Smith   PetscFunctionBegin;
27958c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
2803a40ed3dSBarry Smith   PetscFunctionReturn(0);
2815e76c431SBarry Smith }
28204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2836b8b9a38SLisandro Dalcin 
2846b8b9a38SLisandro Dalcin #undef __FUNCT__
2856b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS"
2866b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes)
2876b8b9a38SLisandro Dalcin {
2886b8b9a38SLisandro Dalcin   PetscErrorCode ierr;
2896b8b9a38SLisandro Dalcin 
2906b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2916b8b9a38SLisandro Dalcin   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2926b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
2936b8b9a38SLisandro Dalcin }
2946b8b9a38SLisandro Dalcin 
29504d965bbSLois Curfman McInnes /*
2964b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2974b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29804d965bbSLois Curfman McInnes 
29904d965bbSLois Curfman McInnes    Input Parameter:
30004d965bbSLois Curfman McInnes .  snes - the SNES context
30104d965bbSLois Curfman McInnes 
302de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30304d965bbSLois Curfman McInnes  */
3044a2ae208SSatish Balay #undef __FUNCT__
3054b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
306dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3075e76c431SBarry Smith {
30894298653SBarry Smith   SNES_LS        *ls = (SNES_LS*) snes->data;
309dfbe8321SBarry Smith   PetscErrorCode ierr;
3103a40ed3dSBarry Smith 
3113a40ed3dSBarry Smith   PetscFunctionBegin;
3126b8b9a38SLisandro Dalcin   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
313649052a6SBarry Smith   ierr = PetscViewerDestroy(&ls->monitor);CHKERRQ(ierr);
3145c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
315557d3f75SLisandro Dalcin 
316557d3f75SLisandro Dalcin   /* clear composed functions */
31758c9b817SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","",PETSC_NULL);CHKERRQ(ierr);
31858c9b817SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","",PETSC_NULL);CHKERRQ(ierr);
319557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
320557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
321557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
322557d3f75SLisandro Dalcin 
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
3245e76c431SBarry Smith }
32504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3264a2ae208SSatish Balay #undef __FUNCT__
32717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
32882bf6240SBarry Smith 
3294b828684SBarry Smith /*@C
33017bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3315e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3325e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3335e76c431SBarry Smith 
3343f9fe445SBarry Smith    Logically Collective on SNES and Vec
335c7afd0dbSLois Curfman McInnes 
3365e76c431SBarry Smith    Input Parameters:
337c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
338af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3395e76c431SBarry Smith .  x - current iterate
3405e76c431SBarry Smith .  f - residual evaluated at x
3413c632250SBarry Smith .  y - search direction
342dc357ad2SBarry Smith .  fnorm - 2-norm of f
343dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
3445e76c431SBarry Smith 
345c4a48953SLois Curfman McInnes    Output Parameters:
346c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3473c632250SBarry Smith .  w - new iterate
3485e76c431SBarry Smith .  gnorm - 2-norm of g
3495e76c431SBarry Smith .  ynorm - 2-norm of search length
350a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
351fee21e36SBarry Smith 
352c4a48953SLois Curfman McInnes    Options Database Key:
35317bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
354c4a48953SLois Curfman McInnes 
35536851e7fSLois Curfman McInnes    Level: advanced
35636851e7fSLois Curfman McInnes 
35728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
35828ae5a21SLois Curfman McInnes 
35917bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
36017bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3615e76c431SBarry Smith @*/
3627087cfbeSBarry 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)
3635e76c431SBarry Smith {
364dfbe8321SBarry Smith   PetscErrorCode ierr;
3654b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
366ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3675334005bSBarry Smith 
3683a40ed3dSBarry Smith   PetscFunctionBegin;
369a7cc72afSBarry Smith   *flag = PETSC_TRUE;
370d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
371a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
37279f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3733c632250SBarry Smith   if (neP->postcheckstep) {
3743c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3758f99978dSLois Curfman McInnes   }
3763c632250SBarry Smith   if (changed_y) {
37779f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3783c632250SBarry Smith   }
3794936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
3804936397dSBarry Smith   if (!snes->domainerror) {
381a3b38805SLois Curfman McInnes     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
38262cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3834936397dSBarry Smith   }
384d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3853a40ed3dSBarry Smith   PetscFunctionReturn(0);
3865e76c431SBarry Smith }
38704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3884a2ae208SSatish Balay #undef __FUNCT__
38917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
39082bf6240SBarry Smith 
39129e0b56fSBarry Smith /*@C
39217bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
39329e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
39429e0b56fSBarry Smith    even compute the norm of the function or search direction; this
39529e0b56fSBarry Smith    is intended only when you know the full step is fine and are
39629e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
397c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
398c7afd0dbSLois Curfman McInnes 
3993f9fe445SBarry Smith    Logically Collective on SNES and Vec
40029e0b56fSBarry Smith 
40129e0b56fSBarry Smith    Input Parameters:
402c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
403af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
40429e0b56fSBarry Smith .  x - current iterate
40529e0b56fSBarry Smith .  f - residual evaluated at x
4063c632250SBarry Smith .  y - search direction
40729e0b56fSBarry Smith .  w - work vector
408dc357ad2SBarry Smith .  fnorm - 2-norm of f
409dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
41029e0b56fSBarry Smith 
41129e0b56fSBarry Smith    Output Parameters:
412c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4133c632250SBarry Smith .  w - new iterate
41429e0b56fSBarry Smith .  gnorm - not changed
41529e0b56fSBarry Smith .  ynorm - not changed
416a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
417fee21e36SBarry Smith 
41829e0b56fSBarry Smith    Options Database Key:
41917bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
42029e0b56fSBarry Smith 
4218cbba510SBarry Smith    Notes:
42217bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
423ea56f5baSLois Curfman McInnes    either the options
424ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
425ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
426329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
427329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
428329f5518SBarry Smith 
429329f5518SBarry Smith    During the final iteration this will not evaluate the function at
430329f5518SBarry Smith    the solution point. This is to save a function evaluation while
431329f5518SBarry Smith    using pseudo-timestepping.
4328cbba510SBarry Smith 
433ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
434a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
435ea56f5baSLois Curfman McInnes    correct, since they are not computed.
436ea56f5baSLois Curfman McInnes 
437ea56f5baSLois Curfman McInnes    Level: advanced
4388cbba510SBarry Smith 
43929e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
44029e0b56fSBarry Smith 
44117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
44217bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
44329e0b56fSBarry Smith @*/
4447087cfbeSBarry 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)
44529e0b56fSBarry Smith {
446dfbe8321SBarry Smith   PetscErrorCode ierr;
4474b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
448ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
44929e0b56fSBarry Smith 
4503a40ed3dSBarry Smith   PetscFunctionBegin;
451a7cc72afSBarry Smith   *flag = PETSC_TRUE;
452d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
45379f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4543c632250SBarry Smith   if (neP->postcheckstep) {
4553c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4563c632250SBarry Smith   }
4573c632250SBarry Smith   if (changed_y) {
45879f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4598f99978dSLois Curfman McInnes   }
460329f5518SBarry Smith 
461329f5518SBarry Smith   /* don't evaluate function the last time through */
462329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4634936397dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
464329f5518SBarry Smith   }
465d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4663a40ed3dSBarry Smith   PetscFunctionReturn(0);
46729e0b56fSBarry Smith }
46804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4694a2ae208SSatish Balay #undef __FUNCT__
47017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4714b828684SBarry Smith /*@C
47217bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4735e76c431SBarry Smith 
474c7afd0dbSLois Curfman McInnes    Collective on SNES
475c7afd0dbSLois Curfman McInnes 
4765e76c431SBarry Smith    Input Parameters:
477c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
478af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4795e76c431SBarry Smith .  x - current iterate
4805e76c431SBarry Smith .  f - residual evaluated at x
4813c632250SBarry Smith .  y - search direction
4825e76c431SBarry Smith .  w - work vector
483dc357ad2SBarry Smith .  fnorm - 2-norm of f
484dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
4855e76c431SBarry Smith 
486393d2d9aSLois Curfman McInnes    Output Parameters:
487c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4883c632250SBarry Smith .  w - new iterate
4895e76c431SBarry Smith .  gnorm - 2-norm of g
4905e76c431SBarry Smith .  ynorm - 2-norm of search length
491a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
492fee21e36SBarry Smith 
493c4a48953SLois Curfman McInnes    Options Database Key:
494e106eecfSBarry Smith +  -snes_ls cubic - Activates SNESLineSearchCubic()
495e106eecfSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
496e106eecfSBarry 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)
497e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
498e106eecfSBarry Smith 
499c4a48953SLois Curfman McInnes 
5005e76c431SBarry Smith    Notes:
5015e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
5025e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
5035e76c431SBarry Smith 
50436851e7fSLois Curfman McInnes    Level: advanced
50536851e7fSLois Curfman McInnes 
50628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
50728ae5a21SLois Curfman McInnes 
50817bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5095e76c431SBarry Smith @*/
5107087cfbeSBarry 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)
5115e76c431SBarry Smith {
512406556e6SLois Curfman McInnes   /*
513406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
514406556e6SLois Curfman McInnes      minimization problem:
515406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
516406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
517406556e6SLois Curfman McInnes    */
518406556e6SLois Curfman McInnes 
519e106eecfSBarry Smith   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
520e106eecfSBarry Smith   PetscReal      minlambda,lambda,lambdatemp;
521aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
522f5b6597dSBarry Smith   PetscScalar    cinitslope;
5236b5873e3SBarry Smith #endif
524dfbe8321SBarry Smith   PetscErrorCode ierr;
525a7cc72afSBarry Smith   PetscInt       count;
5264b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
527ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
52894298653SBarry Smith   MPI_Comm       comm;
5295e76c431SBarry Smith 
5303a40ed3dSBarry Smith   PetscFunctionBegin;
53194298653SBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
532d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
533a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5345e76c431SBarry Smith 
535cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
53675567043SBarry Smith   if (*ynorm == 0.0) {
53794298653SBarry Smith     if (neP->monitor) {
538649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
539649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
540649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
54194298653SBarry Smith     }
542a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5433c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
544a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
545a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
546ad922ac9SBarry Smith     goto theend1;
54794a9d846SBarry Smith   }
548e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
54994298653SBarry Smith     if (neP->monitor) {
550649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
551*8f1a2a5eSBarry 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);
552649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
55394298653SBarry Smith     }
554e106eecfSBarry Smith     ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
555e106eecfSBarry Smith     *ynorm = neP->maxstep;
5565e76c431SBarry Smith   }
557ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
558e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
559a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
560aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
561a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
562329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5635e42470aSBarry Smith #else
564a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5655e42470aSBarry Smith #endif
5665e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5675e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5685e76c431SBarry Smith 
569e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
57043e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
571ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
57243e71028SBarry Smith     *flag = PETSC_FALSE;
57343e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
57443e71028SBarry Smith     goto theend1;
57543e71028SBarry Smith   }
5764936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5774936397dSBarry Smith   if (snes->domainerror) {
5784936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5794936397dSBarry Smith     PetscFunctionReturn(0);
58019717074SBarry Smith   }
581313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
58262cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
583*8f1a2a5eSBarry Smith   ierr = PetscInfo2(snes,"Initial fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
584e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
58594298653SBarry Smith     if (neP->monitor) {
586649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
587*8f1a2a5eSBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Using full step: fnorm %14.12e gnorm %14.12e\n",(double)fnorm,(double)*gnorm);CHKERRQ(ierr);
588649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
58994298653SBarry Smith     }
590ad922ac9SBarry Smith     goto theend1;
5915e76c431SBarry Smith   }
5925e76c431SBarry Smith 
5935e76c431SBarry Smith   /* Fit points with quadratic */
594313b5538SBarry Smith   lambda     = 1.0;
595a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5965e76c431SBarry Smith   lambdaprev = lambda;
5975e76c431SBarry Smith   gnormprev  = *gnorm;
598329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
599ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
600ddd12767SBarry Smith   else                         lambda = lambdatemp;
601275d25c3SBarry Smith 
602e68848bdSBarry Smith   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
60343e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
604ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
60543e71028SBarry Smith     *flag = PETSC_FALSE;
60643e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
60743e71028SBarry Smith     goto theend1;
60843e71028SBarry Smith   }
609f5b6597dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6104936397dSBarry Smith   if (snes->domainerror) {
6114936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6124936397dSBarry Smith     PetscFunctionReturn(0);
61319717074SBarry Smith   }
614cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
61562cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
61694298653SBarry Smith   if (neP->monitor) {
617649052a6SBarry Smith     ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
618*8f1a2a5eSBarry Smith     ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: gnorm after quadratic fit %14.12e\n",(double)*gnorm);CHKERRQ(ierr);
619649052a6SBarry Smith     ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
62094298653SBarry Smith   }
621e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
62294298653SBarry Smith     if (neP->monitor) {
623649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
624649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
625649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
62694298653SBarry Smith     }
627ad922ac9SBarry Smith     goto theend1;
6285e76c431SBarry Smith   }
6295e76c431SBarry Smith 
6305e76c431SBarry Smith   /* Fit points with cubic */
6315e76c431SBarry Smith   count = 1;
632bb9195b6SBarry Smith   while (PETSC_TRUE) {
633dc357ad2SBarry Smith     if (lambda <= minlambda) {
63494298653SBarry Smith       if (neP->monitor) {
635649052a6SBarry Smith         ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
636649052a6SBarry Smith 	ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
637649052a6SBarry 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);
638649052a6SBarry Smith         ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
63994298653SBarry Smith       }
64043e71028SBarry Smith       *flag = PETSC_FALSE;
64143e71028SBarry Smith       break;
6425e76c431SBarry Smith     }
6435d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6445d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
645ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6462b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6475e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6485e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6495e76c431SBarry Smith     if (a == 0.0) {
6505e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6515e76c431SBarry Smith     } else {
6525e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6535e76c431SBarry Smith     }
6545e76c431SBarry Smith     lambdaprev = lambda;
6555e76c431SBarry Smith     gnormprev  = *gnorm;
656329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
657329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6585e76c431SBarry Smith     else                         lambda     = lambdatemp;
659e68848bdSBarry Smith     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
66043e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
661ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
662ae15b995SBarry 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);
663ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
66443e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
665ed98166eSMatthew Knepley       break;
666ed98166eSMatthew Knepley     }
667f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6684936397dSBarry Smith     if (snes->domainerror) {
6694936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6704936397dSBarry Smith       PetscFunctionReturn(0);
67119717074SBarry Smith     }
672cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
67362cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
674e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */
67594298653SBarry Smith       if (neP->monitor) {
676*8f1a2a5eSBarry Smith 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %14.12e lambda=%18.16e\n",(double)*gnorm,(double)lambda);CHKERRQ(ierr);
67794298653SBarry Smith       }
67843e71028SBarry Smith       break;
6792b022350SLois Curfman McInnes     } else {
68094298653SBarry Smith       if (neP->monitor) {
681*8f1a2a5eSBarry 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);
68294298653SBarry Smith       }
6835e76c431SBarry Smith     }
6845e76c431SBarry Smith     count++;
6855e76c431SBarry Smith   }
686ad922ac9SBarry Smith   theend1:
6878f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6883c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6893c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6903c632250SBarry Smith     if (changed_y) {
69179f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6923c632250SBarry Smith     }
6933c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
694f5b6597dSBarry Smith       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6954936397dSBarry Smith       if (snes->domainerror) {
6964936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6974936397dSBarry Smith         PetscFunctionReturn(0);
69819717074SBarry Smith       }
6998f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
70062cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
701a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7028f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
703a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7048f99978dSLois Curfman McInnes     }
7058f99978dSLois Curfman McInnes   }
706d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7073a40ed3dSBarry Smith   PetscFunctionReturn(0);
7085e76c431SBarry Smith }
70904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7104a2ae208SSatish Balay #undef __FUNCT__
71117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
7124b828684SBarry Smith /*@C
71317bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
7145e76c431SBarry Smith 
715c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
716c7afd0dbSLois Curfman McInnes 
7175e42470aSBarry Smith    Input Parameters:
718c7afd0dbSLois Curfman McInnes +  snes - the SNES context
719af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
7205e42470aSBarry Smith .  x - current iterate
7215e42470aSBarry Smith .  f - residual evaluated at x
7223c632250SBarry Smith .  y - search direction
7235e42470aSBarry Smith .  w - work vector
724dc357ad2SBarry Smith .  fnorm - 2-norm of f
725dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
7265e42470aSBarry Smith 
727c4a48953SLois Curfman McInnes    Output Parameters:
7287f3332b4SBarry Smith +  g - residual evaluated at new iterate w
729e106eecfSBarry Smith .  w - new iterate (x + lambda*y)
7305e42470aSBarry Smith .  gnorm - 2-norm of g
7315e42470aSBarry Smith .  ynorm - 2-norm of search length
732a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
733fee21e36SBarry Smith 
734ce9499c7SBarry Smith    Options Database Keys:
735ce9499c7SBarry Smith +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
736ce9499c7SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
737e106eecfSBarry 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)
738e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
739e106eecfSBarry Smith 
7405e42470aSBarry Smith    Notes:
74117bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7425e42470aSBarry Smith 
74336851e7fSLois Curfman McInnes    Level: advanced
74436851e7fSLois Curfman McInnes 
745f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
746f59ffdedSLois Curfman McInnes 
74717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7485e42470aSBarry Smith @*/
7497087cfbeSBarry 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)
7505e76c431SBarry Smith {
751406556e6SLois Curfman McInnes   /*
752406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
753406556e6SLois Curfman McInnes      minimization problem:
754406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
755406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
756406556e6SLois Curfman McInnes    */
757e106eecfSBarry Smith   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
758aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
759f5b6597dSBarry Smith   PetscScalar    cinitslope;
7606b5873e3SBarry Smith #endif
761dfbe8321SBarry Smith   PetscErrorCode ierr;
762a7cc72afSBarry Smith   PetscInt       count;
7634b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
764ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7655e76c431SBarry Smith 
7663a40ed3dSBarry Smith   PetscFunctionBegin;
767d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
768a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7695e76c431SBarry Smith 
7703505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
77163b9fa5eSBarry Smith   if (*ynorm == 0.0) {
77294298653SBarry Smith     if (neP->monitor) {
773649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
774649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
775649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
77694298653SBarry Smith     }
777b37302c6SLois Curfman McInnes     *gnorm = fnorm;
778e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
779b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
780a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
781ad922ac9SBarry Smith     goto theend2;
78294a9d846SBarry Smith   }
783e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
784e106eecfSBarry Smith     ierr   = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
785e106eecfSBarry Smith     *ynorm = neP->maxstep;
7865e76c431SBarry Smith   }
787ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
788e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
789a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
790aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
791a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
792329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7935e42470aSBarry Smith #else
794a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7955e42470aSBarry Smith #endif
7965e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7975e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
798*8f1a2a5eSBarry Smith   ierr = PetscInfo1(snes,"Initslope %14.12e \n",(double)initslope);CHKERRQ(ierr);
7995e42470aSBarry Smith 
800e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
80143e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
802ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
80343e71028SBarry Smith     *flag = PETSC_FALSE;
80443e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
80543e71028SBarry Smith     goto theend2;
80643e71028SBarry Smith   }
8074936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8084936397dSBarry Smith   if (snes->domainerror) {
8094936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8104936397dSBarry Smith     PetscFunctionReturn(0);
81119717074SBarry Smith   }
812cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
81362cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
814e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
81594298653SBarry Smith     if (neP->monitor) {
816649052a6SBarry Smith       ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
817649052a6SBarry Smith       ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
818649052a6SBarry Smith       ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
81994298653SBarry Smith     }
820ad922ac9SBarry Smith     goto theend2;
8215e42470aSBarry Smith   }
8225e42470aSBarry Smith 
8235e42470aSBarry Smith   /* Fit points with quadratic */
824313b5538SBarry Smith   lambda = 1.0;
8255e42470aSBarry Smith   count = 1;
8265ca10a37SBarry Smith   while (PETSC_TRUE) {
8275e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
82894298653SBarry Smith       if (neP->monitor) {
829649052a6SBarry Smith         ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
830649052a6SBarry Smith         ierr = PetscViewerASCIIPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
831*8f1a2a5eSBarry 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);
832649052a6SBarry Smith         ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
83394298653SBarry Smith       }
834e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
83543e71028SBarry Smith       *flag = PETSC_FALSE;
83643e71028SBarry Smith       break;
8375e42470aSBarry Smith     }
838a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
839329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
840329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
841329f5518SBarry Smith     else                         lambda     = lambdatemp;
842275d25c3SBarry Smith 
843e68848bdSBarry Smith     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
84443e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
845ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
846d9822059SBarry 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);
847ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
84843e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
849ed98166eSMatthew Knepley       break;
850ed98166eSMatthew Knepley     }
851f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8524936397dSBarry Smith     if (snes->domainerror) {
8534936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8544936397dSBarry Smith       PetscFunctionReturn(0);
85519717074SBarry Smith     }
856cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
85762cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
858e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
85994298653SBarry Smith       if (neP->monitor) {
860649052a6SBarry Smith         ierr = PetscViewerASCIIAddTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
861*8f1a2a5eSBarry Smith         ierr = PetscViewerASCIIPrintf(neP->monitor,"    Line Search: Quadratically determined step, lambda=%14.12e\n",(double)lambda);CHKERRQ(ierr);
862649052a6SBarry Smith         ierr = PetscViewerASCIISubtractTab(neP->monitor,((PetscObject)snes)->tablevel);CHKERRQ(ierr);
86394298653SBarry Smith       }
86406259719SBarry Smith       break;
8655e42470aSBarry Smith     }
8665e42470aSBarry Smith     count++;
8675e42470aSBarry Smith   }
868ad922ac9SBarry Smith   theend2:
8698f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8703c632250SBarry Smith   if (neP->postcheckstep) {
8713c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8723c632250SBarry Smith     if (changed_y) {
87379f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8743c632250SBarry Smith     }
8753c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8763c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
8774936397dSBarry Smith       if (snes->domainerror) {
8784936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8794936397dSBarry Smith         PetscFunctionReturn(0);
88019717074SBarry Smith       }
8818f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
882a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
8838f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
884a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
88562cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8868f99978dSLois Curfman McInnes     }
8878f99978dSLois Curfman McInnes   }
888d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8893a40ed3dSBarry Smith   PetscFunctionReturn(0);
8905e76c431SBarry Smith }
8912343ba6eSBarry Smith 
89204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8934a2ae208SSatish Balay #undef __FUNCT__
89417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
895c9e6a524SLois Curfman McInnes /*@C
89617bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8974b27c08aSLois Curfman McInnes    by the method SNESLS.
8985e76c431SBarry Smith 
8995e76c431SBarry Smith    Input Parameters:
900c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
901af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
902c7afd0dbSLois Curfman McInnes -  func - pointer to int function
9035e76c431SBarry Smith 
9043f9fe445SBarry Smith    Logically Collective on SNES
905fee21e36SBarry Smith 
906c4a48953SLois Curfman McInnes    Available Routines:
90717bae607SBarry Smith +  SNESLineSearchCubic() - default line search
90817bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
90917bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
91017bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
9115e76c431SBarry Smith 
912c4a48953SLois Curfman McInnes     Options Database Keys:
9134b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
9144b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
915e106eecfSBarry 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)
916e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
917c4a48953SLois Curfman McInnes 
9185e76c431SBarry Smith    Calling sequence of func:
919bd208895SLois Curfman McInnes .vb
920ace3abfcSBarry 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)
921bd208895SLois Curfman McInnes .ve
9225e76c431SBarry Smith 
9235e76c431SBarry Smith     Input parameters for func:
924c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
925af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
9265e76c431SBarry Smith .   x - current iterate
9275e76c431SBarry Smith .   f - residual evaluated at x
9283c632250SBarry Smith .   y - search direction
929c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9305e76c431SBarry Smith 
9315e76c431SBarry Smith     Output parameters for func:
932c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9333c632250SBarry Smith .   w - new iterate
9345e76c431SBarry Smith .   gnorm - 2-norm of g
9355e76c431SBarry Smith .   ynorm - 2-norm of search length
936a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
937f59ffdedSLois Curfman McInnes 
93836851e7fSLois Curfman McInnes     Level: advanced
93936851e7fSLois Curfman McInnes 
940f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
941f59ffdedSLois Curfman McInnes 
94217bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
94317bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9445e76c431SBarry Smith @*/
9457087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx)
9465e76c431SBarry Smith {
9474ac538c5SBarry Smith   PetscErrorCode ierr;
94882bf6240SBarry Smith 
9493a40ed3dSBarry Smith   PetscFunctionBegin;
9504ac538c5SBarry 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);
9513a40ed3dSBarry Smith   PetscFunctionReturn(0);
9525e76c431SBarry Smith }
9538e019c35SBarry Smith 
954ace3abfcSBarry 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*/
95504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
956fb2e594dSBarry Smith EXTERN_C_BEGIN
9574a2ae208SSatish Balay #undef __FUNCT__
95817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
9597087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
96082bf6240SBarry Smith {
96182bf6240SBarry Smith   PetscFunctionBegin;
9624b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9634b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
96482bf6240SBarry Smith   PetscFunctionReturn(0);
96582bf6240SBarry Smith }
966fb2e594dSBarry Smith EXTERN_C_END
96704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9684a2ae208SSatish Balay #undef __FUNCT__
9693c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
970c8dd0c0dSLois Curfman McInnes /*@C
9713c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9724b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
973c8dd0c0dSLois Curfman McInnes 
974c8dd0c0dSLois Curfman McInnes    Input Parameters:
975c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9763c632250SBarry Smith .  func - pointer to function
977c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
978c8dd0c0dSLois Curfman McInnes 
9793f9fe445SBarry Smith    Logically Collective on SNES
980c8dd0c0dSLois Curfman McInnes 
981c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
982c8dd0c0dSLois Curfman McInnes .vb
983ace3abfcSBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool  *changed_y,PetscBool  *changed_w)
984c8dd0c0dSLois Curfman McInnes .ve
985b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
986b1ae27eaSLois Curfman McInnes    on failure.
987c8dd0c0dSLois Curfman McInnes 
988c8dd0c0dSLois Curfman McInnes    Input parameters for func:
989c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
990c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9913c632250SBarry Smith .  x - previous iterate
9923c632250SBarry Smith .  y - new search direction and length
9933c632250SBarry Smith -  w - current candidate iterate
994c8dd0c0dSLois Curfman McInnes 
995c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9963c632250SBarry Smith +  y - search direction (possibly changed)
9973c632250SBarry Smith .  w - current iterate (possibly modified)
9983c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9993c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
1000c8dd0c0dSLois Curfman McInnes 
1001c8dd0c0dSLois Curfman McInnes    Level: advanced
1002c8dd0c0dSLois Curfman McInnes 
10039e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10049e247f21SBarry Smith 
10053c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
10063c632250SBarry Smith 
1007481b4421SBarry Smith    On input w = x - y
10083c632250SBarry Smith 
100917bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
1010b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
1011ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
10128f99978dSLois Curfman McInnes 
101317bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
1014ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
1015ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
1016ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
10179e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
1018b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
10198f99978dSLois Curfman McInnes 
1020c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
1021c8dd0c0dSLois Curfman McInnes 
102217bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1023c8dd0c0dSLois Curfman McInnes @*/
10247087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx)
1025c8dd0c0dSLois Curfman McInnes {
10264ac538c5SBarry Smith   PetscErrorCode ierr;
1027c8dd0c0dSLois Curfman McInnes 
1028c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10294ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr);
1030c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1031c8dd0c0dSLois Curfman McInnes }
103294298653SBarry Smith 
10339c3ca13aSBarry Smith #undef __FUNCT__
10349c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10359c3ca13aSBarry Smith /*@C
10369c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10377e4bb74cSBarry Smith          before the line search is called.
10389c3ca13aSBarry Smith 
10399c3ca13aSBarry Smith    Input Parameters:
10409c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10419c3ca13aSBarry Smith .  func - pointer to function
10429c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10439c3ca13aSBarry Smith 
10443f9fe445SBarry Smith    Logically Collective on SNES
10459c3ca13aSBarry Smith 
10469c3ca13aSBarry Smith    Calling sequence of func:
10479c3ca13aSBarry Smith .vb
1048ace3abfcSBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool  *changed_y)
10499c3ca13aSBarry Smith .ve
10509c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10519c3ca13aSBarry Smith    on failure.
10529c3ca13aSBarry Smith 
10539c3ca13aSBarry Smith    Input parameters for func:
10549c3ca13aSBarry Smith +  snes - nonlinear context
10559c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10569c3ca13aSBarry Smith .  x - previous iterate
10579c3ca13aSBarry Smith -  y - new search direction and length
10589c3ca13aSBarry Smith 
10599c3ca13aSBarry Smith    Output parameters for func:
10609c3ca13aSBarry Smith +  y - search direction (possibly changed)
10619c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10629c3ca13aSBarry Smith 
10639c3ca13aSBarry Smith    Level: advanced
10649c3ca13aSBarry Smith 
10659c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10669c3ca13aSBarry Smith 
10679c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10689c3ca13aSBarry Smith 
10697e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10709c3ca13aSBarry Smith @*/
10717087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx)
10729c3ca13aSBarry Smith {
10734ac538c5SBarry Smith   PetscErrorCode ierr;
10749c3ca13aSBarry Smith 
10759c3ca13aSBarry Smith   PetscFunctionBegin;
10764ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr);
10779c3ca13aSBarry Smith   PetscFunctionReturn(0);
10789c3ca13aSBarry Smith }
10799c3ca13aSBarry Smith 
108094298653SBarry Smith #undef __FUNCT__
108194298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor"
108294298653SBarry Smith /*@C
108394298653SBarry Smith    SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search
108494298653SBarry Smith 
108594298653SBarry Smith    Input Parameters:
108694298653SBarry Smith +  snes - nonlinear context obtained from SNESCreate()
108794298653SBarry Smith -  flg - PETSC_TRUE to monitor the line search
108894298653SBarry Smith 
10893f9fe445SBarry Smith    Logically Collective on SNES
109094298653SBarry Smith 
109194298653SBarry Smith    Options Database:
109294298653SBarry Smith .   -snes_ls_monitor
109394298653SBarry Smith 
109494298653SBarry Smith    Level: intermediate
109594298653SBarry Smith 
109694298653SBarry Smith 
109794298653SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
109894298653SBarry Smith @*/
10997087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor(SNES snes,PetscBool  flg)
110094298653SBarry Smith {
11014ac538c5SBarry Smith   PetscErrorCode ierr;
110294298653SBarry Smith 
110394298653SBarry Smith   PetscFunctionBegin;
11044ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr);
110594298653SBarry Smith   PetscFunctionReturn(0);
110694298653SBarry Smith }
110794298653SBarry Smith 
1108c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1109ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/
1110ace3abfcSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *);                 /* force argument to next function to not be extern C*/
1111c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
11124a2ae208SSatish Balay #undef __FUNCT__
11133c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
11147087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1115c8dd0c0dSLois Curfman McInnes {
1116c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
11173c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
11183c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1119c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1120c8dd0c0dSLois Curfman McInnes }
1121c8dd0c0dSLois Curfman McInnes EXTERN_C_END
11229c3ca13aSBarry Smith 
11239c3ca13aSBarry Smith EXTERN_C_BEGIN
11249c3ca13aSBarry Smith #undef __FUNCT__
11259c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
11267087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
11279c3ca13aSBarry Smith {
11289c3ca13aSBarry Smith   PetscFunctionBegin;
11299c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
11309c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
11319c3ca13aSBarry Smith   PetscFunctionReturn(0);
11329c3ca13aSBarry Smith }
11339c3ca13aSBarry Smith EXTERN_C_END
1134329e7e40SMatthew Knepley 
113594298653SBarry Smith EXTERN_C_BEGIN
113694298653SBarry Smith #undef __FUNCT__
113794298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor_LS"
11387087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor_LS(SNES snes,PetscBool  flg)
113994298653SBarry Smith {
114094298653SBarry Smith   SNES_LS        *ls = (SNES_LS*)snes->data;
114194298653SBarry Smith   PetscErrorCode ierr;
114294298653SBarry Smith 
114394298653SBarry Smith   PetscFunctionBegin;
114494298653SBarry Smith   if (flg && !ls->monitor) {
1145649052a6SBarry Smith     ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,"stdout",&ls->monitor);CHKERRQ(ierr);
114694298653SBarry Smith   } else if (!flg && ls->monitor) {
1147649052a6SBarry Smith     ierr = PetscViewerDestroy(&ls->monitor);CHKERRQ(ierr);
114894298653SBarry Smith   }
114994298653SBarry Smith   PetscFunctionReturn(0);
115094298653SBarry Smith }
115194298653SBarry Smith EXTERN_C_END
115294298653SBarry Smith 
1153329e7e40SMatthew Knepley /*
11544b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
115504d965bbSLois Curfman McInnes 
115604d965bbSLois Curfman McInnes    Input Parameters:
115704d965bbSLois Curfman McInnes .  SNES - the SNES context
115804d965bbSLois Curfman McInnes .  viewer - visualization context
115904d965bbSLois Curfman McInnes 
116004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
116104d965bbSLois Curfman McInnes */
11624a2ae208SSatish Balay #undef __FUNCT__
11634b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11646849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1165a935fc98SLois Curfman McInnes {
11664b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11672fc52814SBarry Smith   const char     *cstr;
1168dfbe8321SBarry Smith   PetscErrorCode ierr;
1169ace3abfcSBarry Smith   PetscBool      iascii;
1170a935fc98SLois Curfman McInnes 
11713a40ed3dSBarry Smith   PetscFunctionBegin;
11722692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
117332077d6dSBarry Smith   if (iascii) {
117417bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
117517bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
117617bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
117719bcc07fSBarry Smith     else                                                cstr = "unknown";
1178b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1179*8f1a2a5eSBarry 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);
118019bcc07fSBarry Smith   }
11813a40ed3dSBarry Smith   PetscFunctionReturn(0);
1182a935fc98SLois Curfman McInnes }
118304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
118404d965bbSLois Curfman McInnes /*
11854b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
118604d965bbSLois Curfman McInnes 
118704d965bbSLois Curfman McInnes    Input Parameter:
118804d965bbSLois Curfman McInnes .  snes - the SNES context
118904d965bbSLois Curfman McInnes 
119004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
119104d965bbSLois Curfman McInnes */
11924a2ae208SSatish Balay #undef __FUNCT__
11934b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11946849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11955e42470aSBarry Smith {
11964b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1197e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1198dfbe8321SBarry Smith   PetscErrorCode ierr;
1199a7cc72afSBarry Smith   PetscInt       indx;
1200ace3abfcSBarry Smith   PetscBool      flg,set;
12015e42470aSBarry Smith 
12023a40ed3dSBarry Smith   PetscFunctionBegin;
1203b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
12044b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
12054b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1206e106eecfSBarry Smith     ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr);
1207acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
120894298653SBarry Smith     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
1209186905e3SBarry Smith 
121017bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
121125cdf11fSBarry Smith     if (flg) {
121222e36657SBarry Smith       switch (indx) {
1213b49fd9e1SBarry Smith       case 0:
121417bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1215b49fd9e1SBarry Smith         break;
1216b49fd9e1SBarry Smith       case 1:
121717bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1218b49fd9e1SBarry Smith         break;
1219b49fd9e1SBarry Smith       case 2:
122017bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1221b49fd9e1SBarry Smith         break;
1222b49fd9e1SBarry Smith       case 3:
122317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1224b49fd9e1SBarry Smith         break;
12255e42470aSBarry Smith       }
12265e42470aSBarry Smith     }
1227b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
12283a40ed3dSBarry Smith   PetscFunctionReturn(0);
12295e42470aSBarry Smith }
12304827ddcaSBarry Smith 
12314827ddcaSBarry Smith EXTERN_C_BEGIN
12324827ddcaSBarry Smith extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
12334827ddcaSBarry Smith EXTERN_C_END
12344827ddcaSBarry Smith 
123504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1236ebe3b25bSBarry Smith /*MC
1237ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
123804d965bbSLois Curfman McInnes 
1239ebe3b25bSBarry Smith    Options Database:
1240ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1241ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1242e106eecfSBarry 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)
1243acbee50cSBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1244acbee50cSBarry Smith -   -snes_ls_monitor - print information about progress of line searches
1245acbee50cSBarry Smith 
124604d965bbSLois Curfman McInnes 
1247ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
124804d965bbSLois Curfman McInnes 
1249ee3001cbSBarry Smith    Level: beginner
1250ee3001cbSBarry Smith 
125117bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
125217bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1253b3dd4ab5SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
1254ebe3b25bSBarry Smith 
1255ebe3b25bSBarry Smith M*/
1256fb2e594dSBarry Smith EXTERN_C_BEGIN
12574a2ae208SSatish Balay #undef __FUNCT__
12584b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
12597087cfbeSBarry Smith PetscErrorCode  SNESCreate_LS(SNES snes)
12605e42470aSBarry Smith {
1261dfbe8321SBarry Smith   PetscErrorCode ierr;
12624b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12635e42470aSBarry Smith 
12643a40ed3dSBarry Smith   PetscFunctionBegin;
1265e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1266e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1267e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1268e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1269e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
12706b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_LS;
12715e42470aSBarry Smith 
127238f2d2fdSLisandro Dalcin   ierr                  = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
12735e42470aSBarry Smith   snes->data    	= (void*)neP;
12745e42470aSBarry Smith   neP->alpha		= 1.e-4;
12755e42470aSBarry Smith   neP->maxstep		= 1.e8;
1276e106eecfSBarry Smith   neP->minlambda        = 1.e-12;
127717bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1278c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12793c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12803c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12813c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12823c632250SBarry Smith   neP->precheck         = PETSC_NULL;
128382bf6240SBarry Smith 
128494298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr);
12854827ddcaSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);CHKERRQ(ierr);
128694298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
128794298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
128894298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
128982bf6240SBarry Smith 
12903a40ed3dSBarry Smith   PetscFunctionReturn(0);
12915e42470aSBarry Smith }
1292fb2e594dSBarry Smith EXTERN_C_END
129382bf6240SBarry Smith 
129482bf6240SBarry Smith 
129582bf6240SBarry Smith 
1296