xref: /petsc/src/snes/impls/ls/ls.c (revision 4827ddca42da2d0ec0d1a0d9b26e4c5ddc271c44)
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);
251e2582c4SBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %G near zero implies found a local minimum\n",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);
3736669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
3836669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
391e2582c4SBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",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) {
671e2582c4SBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",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);
17294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
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);
219ae15b995SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,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);
24485385478SLisandro Dalcin     SNESMonitor(snes,snes->iter,snes->norm);
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;
27985385478SLisandro Dalcin   if (!snes->vec_sol_update) {
28085385478SLisandro Dalcin     ierr = VecDuplicate(snes->vec_sol,&snes->vec_sol_update);CHKERRQ(ierr);
28185385478SLisandro Dalcin     ierr = PetscLogObjectParent(snes,snes->vec_sol_update);CHKERRQ(ierr);
28285385478SLisandro Dalcin   }
283e74804ceSBarry Smith   if (!snes->work) {
28485385478SLisandro Dalcin     snes->nwork = 3;
28585385478SLisandro Dalcin     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
286efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
287e74804ceSBarry Smith   }
2883a40ed3dSBarry Smith   PetscFunctionReturn(0);
2895e76c431SBarry Smith }
29004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
29104d965bbSLois Curfman McInnes /*
2924b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2934b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29404d965bbSLois Curfman McInnes 
29504d965bbSLois Curfman McInnes    Input Parameter:
29604d965bbSLois Curfman McInnes .  snes - the SNES context
29704d965bbSLois Curfman McInnes 
298de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29904d965bbSLois Curfman McInnes  */
3004a2ae208SSatish Balay #undef __FUNCT__
3014b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
302dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3035e76c431SBarry Smith {
30494298653SBarry Smith   SNES_LS        *ls = (SNES_LS*) snes->data;
305dfbe8321SBarry Smith   PetscErrorCode ierr;
3063a40ed3dSBarry Smith 
3073a40ed3dSBarry Smith   PetscFunctionBegin;
30885385478SLisandro Dalcin   if (snes->vec_sol_update) {
30985385478SLisandro Dalcin     ierr = VecDestroy(snes->vec_sol_update);CHKERRQ(ierr);
31085385478SLisandro Dalcin     snes->vec_sol_update = PETSC_NULL;
31185385478SLisandro Dalcin   }
312574deadeSBarry Smith   if (snes->work) {
313574deadeSBarry Smith     ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);
31421c89e3eSBarry Smith   }
31594298653SBarry Smith   if (ls->monitor) {
31694298653SBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr);
31794298653SBarry Smith   }
3185c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
319557d3f75SLisandro Dalcin 
320557d3f75SLisandro Dalcin   /* clear composed functions */
321557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
322557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
323557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
324557d3f75SLisandro Dalcin 
3253a40ed3dSBarry Smith   PetscFunctionReturn(0);
3265e76c431SBarry Smith }
32704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3284a2ae208SSatish Balay #undef __FUNCT__
32917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
33082bf6240SBarry Smith 
3314b828684SBarry Smith /*@C
33217bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3335e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3345e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3355e76c431SBarry Smith 
3363f9fe445SBarry Smith    Logically Collective on SNES and Vec
337c7afd0dbSLois Curfman McInnes 
3385e76c431SBarry Smith    Input Parameters:
339c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
340af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3415e76c431SBarry Smith .  x - current iterate
3425e76c431SBarry Smith .  f - residual evaluated at x
3433c632250SBarry Smith .  y - search direction
344dc357ad2SBarry Smith .  fnorm - 2-norm of f
345dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
3465e76c431SBarry Smith 
347c4a48953SLois Curfman McInnes    Output Parameters:
348c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3493c632250SBarry Smith .  w - new iterate
3505e76c431SBarry Smith .  gnorm - 2-norm of g
3515e76c431SBarry Smith .  ynorm - 2-norm of search length
352a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
353fee21e36SBarry Smith 
354c4a48953SLois Curfman McInnes    Options Database Key:
35517bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
356c4a48953SLois Curfman McInnes 
35736851e7fSLois Curfman McInnes    Level: advanced
35836851e7fSLois Curfman McInnes 
35928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
36028ae5a21SLois Curfman McInnes 
36117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
36217bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3635e76c431SBarry Smith @*/
3647087cfbeSBarry 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)
3655e76c431SBarry Smith {
366dfbe8321SBarry Smith   PetscErrorCode ierr;
3674b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
368ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3695334005bSBarry Smith 
3703a40ed3dSBarry Smith   PetscFunctionBegin;
371a7cc72afSBarry Smith   *flag = PETSC_TRUE;
372d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
373a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
37479f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3753c632250SBarry Smith   if (neP->postcheckstep) {
3763c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3778f99978dSLois Curfman McInnes   }
3783c632250SBarry Smith   if (changed_y) {
37979f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3803c632250SBarry Smith   }
3814936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
3824936397dSBarry Smith   if (!snes->domainerror) {
383a3b38805SLois Curfman McInnes     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
38462cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
3854936397dSBarry Smith   }
386d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3873a40ed3dSBarry Smith   PetscFunctionReturn(0);
3885e76c431SBarry Smith }
38904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3904a2ae208SSatish Balay #undef __FUNCT__
39117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
39282bf6240SBarry Smith 
39329e0b56fSBarry Smith /*@C
39417bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
39529e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
39629e0b56fSBarry Smith    even compute the norm of the function or search direction; this
39729e0b56fSBarry Smith    is intended only when you know the full step is fine and are
39829e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
399c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
400c7afd0dbSLois Curfman McInnes 
4013f9fe445SBarry Smith    Logically Collective on SNES and Vec
40229e0b56fSBarry Smith 
40329e0b56fSBarry Smith    Input Parameters:
404c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
405af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
40629e0b56fSBarry Smith .  x - current iterate
40729e0b56fSBarry Smith .  f - residual evaluated at x
4083c632250SBarry Smith .  y - search direction
40929e0b56fSBarry Smith .  w - work vector
410dc357ad2SBarry Smith .  fnorm - 2-norm of f
411dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
41229e0b56fSBarry Smith 
41329e0b56fSBarry Smith    Output Parameters:
414c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4153c632250SBarry Smith .  w - new iterate
41629e0b56fSBarry Smith .  gnorm - not changed
41729e0b56fSBarry Smith .  ynorm - not changed
418a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
419fee21e36SBarry Smith 
42029e0b56fSBarry Smith    Options Database Key:
42117bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
42229e0b56fSBarry Smith 
4238cbba510SBarry Smith    Notes:
42417bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
425ea56f5baSLois Curfman McInnes    either the options
426ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
427ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
428329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
429329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
430329f5518SBarry Smith 
431329f5518SBarry Smith    During the final iteration this will not evaluate the function at
432329f5518SBarry Smith    the solution point. This is to save a function evaluation while
433329f5518SBarry Smith    using pseudo-timestepping.
4348cbba510SBarry Smith 
435ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
436a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
437ea56f5baSLois Curfman McInnes    correct, since they are not computed.
438ea56f5baSLois Curfman McInnes 
439ea56f5baSLois Curfman McInnes    Level: advanced
4408cbba510SBarry Smith 
44129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
44229e0b56fSBarry Smith 
44317bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
44417bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
44529e0b56fSBarry Smith @*/
4467087cfbeSBarry 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)
44729e0b56fSBarry Smith {
448dfbe8321SBarry Smith   PetscErrorCode ierr;
4494b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
450ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
45129e0b56fSBarry Smith 
4523a40ed3dSBarry Smith   PetscFunctionBegin;
453a7cc72afSBarry Smith   *flag = PETSC_TRUE;
454d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
45579f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4563c632250SBarry Smith   if (neP->postcheckstep) {
4573c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4583c632250SBarry Smith   }
4593c632250SBarry Smith   if (changed_y) {
46079f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4618f99978dSLois Curfman McInnes   }
462329f5518SBarry Smith 
463329f5518SBarry Smith   /* don't evaluate function the last time through */
464329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4654936397dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
466329f5518SBarry Smith   }
467d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4683a40ed3dSBarry Smith   PetscFunctionReturn(0);
46929e0b56fSBarry Smith }
47004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4714a2ae208SSatish Balay #undef __FUNCT__
47217bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4734b828684SBarry Smith /*@C
47417bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4755e76c431SBarry Smith 
476c7afd0dbSLois Curfman McInnes    Collective on SNES
477c7afd0dbSLois Curfman McInnes 
4785e76c431SBarry Smith    Input Parameters:
479c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
480af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4815e76c431SBarry Smith .  x - current iterate
4825e76c431SBarry Smith .  f - residual evaluated at x
4833c632250SBarry Smith .  y - search direction
4845e76c431SBarry Smith .  w - work vector
485dc357ad2SBarry Smith .  fnorm - 2-norm of f
486dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
4875e76c431SBarry Smith 
488393d2d9aSLois Curfman McInnes    Output Parameters:
489c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4903c632250SBarry Smith .  w - new iterate
4915e76c431SBarry Smith .  gnorm - 2-norm of g
4925e76c431SBarry Smith .  ynorm - 2-norm of search length
493a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
494fee21e36SBarry Smith 
495c4a48953SLois Curfman McInnes    Options Database Key:
496e106eecfSBarry Smith +  -snes_ls cubic - Activates SNESLineSearchCubic()
497e106eecfSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
498e106eecfSBarry 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)
499e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
500e106eecfSBarry Smith 
501c4a48953SLois Curfman McInnes 
5025e76c431SBarry Smith    Notes:
5035e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
5045e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
5055e76c431SBarry Smith 
50636851e7fSLois Curfman McInnes    Level: advanced
50736851e7fSLois Curfman McInnes 
50828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
50928ae5a21SLois Curfman McInnes 
51017bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5115e76c431SBarry Smith @*/
5127087cfbeSBarry 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)
5135e76c431SBarry Smith {
514406556e6SLois Curfman McInnes   /*
515406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
516406556e6SLois Curfman McInnes      minimization problem:
517406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
518406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
519406556e6SLois Curfman McInnes    */
520406556e6SLois Curfman McInnes 
521e106eecfSBarry Smith   PetscReal      initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
522e106eecfSBarry Smith   PetscReal      minlambda,lambda,lambdatemp;
523aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
524f5b6597dSBarry Smith   PetscScalar    cinitslope;
5256b5873e3SBarry Smith #endif
526dfbe8321SBarry Smith   PetscErrorCode ierr;
527a7cc72afSBarry Smith   PetscInt       count;
5284b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
529ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
53094298653SBarry Smith   MPI_Comm       comm;
5315e76c431SBarry Smith 
5323a40ed3dSBarry Smith   PetscFunctionBegin;
53394298653SBarry Smith   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
534d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
535a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5365e76c431SBarry Smith 
537cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
53875567043SBarry Smith   if (*ynorm == 0.0) {
53994298653SBarry Smith     if (neP->monitor) {
54094298653SBarry Smith       ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line search: Initial direction and size is 0\n");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) {
55094298653SBarry Smith       ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line search: Scaling step by %G old ynorm %G\n",neP->maxstep/(*ynorm),*ynorm);CHKERRQ(ierr);
55194298653SBarry Smith     }
552e106eecfSBarry Smith     ierr = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
553e106eecfSBarry Smith     *ynorm = neP->maxstep;
5545e76c431SBarry Smith   }
555ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
556e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
557a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
558aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
559a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
560329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5615e42470aSBarry Smith #else
562a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5635e42470aSBarry Smith #endif
5645e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5655e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5665e76c431SBarry Smith 
567e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
56843e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
569ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
57043e71028SBarry Smith     *flag = PETSC_FALSE;
57143e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
57243e71028SBarry Smith     goto theend1;
57343e71028SBarry Smith   }
5744936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
5754936397dSBarry Smith   if (snes->domainerror) {
5764936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
5774936397dSBarry Smith     PetscFunctionReturn(0);
57819717074SBarry Smith   }
579313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
58062cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
581f5b6597dSBarry Smith   ierr = PetscInfo2(snes,"Initial fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
582e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
58394298653SBarry Smith     if (neP->monitor) {
58494298653SBarry Smith       ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line search: Using full step: fnorm %G gnorm %G\n",fnorm,*gnorm);CHKERRQ(ierr);
58594298653SBarry Smith     }
586ad922ac9SBarry Smith     goto theend1;
5875e76c431SBarry Smith   }
5885e76c431SBarry Smith 
5895e76c431SBarry Smith   /* Fit points with quadratic */
590313b5538SBarry Smith   lambda     = 1.0;
591a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5925e76c431SBarry Smith   lambdaprev = lambda;
5935e76c431SBarry Smith   gnormprev  = *gnorm;
594329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
595ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
596ddd12767SBarry Smith   else                         lambda = lambdatemp;
597275d25c3SBarry Smith 
598e68848bdSBarry Smith   ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
59943e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
600ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
60143e71028SBarry Smith     *flag = PETSC_FALSE;
60243e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
60343e71028SBarry Smith     goto theend1;
60443e71028SBarry Smith   }
605f5b6597dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6064936397dSBarry Smith   if (snes->domainerror) {
6074936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6084936397dSBarry Smith     PetscFunctionReturn(0);
60919717074SBarry Smith   }
610cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
61162cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
61294298653SBarry Smith   if (neP->monitor) {
61394298653SBarry Smith     ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line search: gnorm after quadratic fit %G\n",*gnorm);CHKERRQ(ierr);
61494298653SBarry Smith   }
615e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
61694298653SBarry Smith     if (neP->monitor) {
61794298653SBarry Smith       ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line search: Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
61894298653SBarry Smith     }
619ad922ac9SBarry Smith     goto theend1;
6205e76c431SBarry Smith   }
6215e76c431SBarry Smith 
6225e76c431SBarry Smith   /* Fit points with cubic */
6235e76c431SBarry Smith   count = 1;
624bb9195b6SBarry Smith   while (PETSC_TRUE) {
625dc357ad2SBarry Smith     if (lambda <= minlambda) {
62694298653SBarry Smith       if (neP->monitor) {
62794298653SBarry Smith 	ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line search: unable to find good step length! After %D tries \n",count);CHKERRQ(ierr);
628d9822059SBarry Smith 	ierr = PetscViewerASCIIMonitorPrintf(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);
62994298653SBarry Smith       }
63043e71028SBarry Smith       *flag = PETSC_FALSE;
63143e71028SBarry Smith       break;
6325e76c431SBarry Smith     }
6335d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6345d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
635ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6362b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6375e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6385e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6395e76c431SBarry Smith     if (a == 0.0) {
6405e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6415e76c431SBarry Smith     } else {
6425e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6435e76c431SBarry Smith     }
6445e76c431SBarry Smith     lambdaprev = lambda;
6455e76c431SBarry Smith     gnormprev  = *gnorm;
646329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
647329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6485e76c431SBarry Smith     else                         lambda     = lambdatemp;
649e68848bdSBarry Smith     ierr  = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
65043e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
651ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
652ae15b995SBarry 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);
653ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
65443e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
655ed98166eSMatthew Knepley       break;
656ed98166eSMatthew Knepley     }
657f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6584936397dSBarry Smith     if (snes->domainerror) {
6594936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6604936397dSBarry Smith       PetscFunctionReturn(0);
66119717074SBarry Smith     }
662cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
66362cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
664e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* is reduction enough? */
66594298653SBarry Smith       if (neP->monitor) {
666d9822059SBarry Smith 	ierr = PetscPrintf(comm,"    Line search: Cubically determined step, current gnorm %G lambda=%18.16e\n",*gnorm,(double)lambda);CHKERRQ(ierr);
66794298653SBarry Smith       }
66843e71028SBarry Smith       break;
6692b022350SLois Curfman McInnes     } else {
67094298653SBarry Smith       if (neP->monitor) {
671d9822059SBarry Smith         ierr = PetscPrintf(comm,"    Line search: Cubic step no good, shrinking lambda, current gnorm %G lambda=%18.16e\n",*gnorm,(double)lambda);CHKERRQ(ierr);
67294298653SBarry Smith       }
6735e76c431SBarry Smith     }
6745e76c431SBarry Smith     count++;
6755e76c431SBarry Smith   }
676ad922ac9SBarry Smith   theend1:
6778f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6783c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6793c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6803c632250SBarry Smith     if (changed_y) {
68179f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6823c632250SBarry Smith     }
6833c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
684f5b6597dSBarry Smith       ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
6854936397dSBarry Smith       if (snes->domainerror) {
6864936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6874936397dSBarry Smith         PetscFunctionReturn(0);
68819717074SBarry Smith       }
6898f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
69062cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
691a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6928f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
693a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6948f99978dSLois Curfman McInnes     }
6958f99978dSLois Curfman McInnes   }
696d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6973a40ed3dSBarry Smith   PetscFunctionReturn(0);
6985e76c431SBarry Smith }
69904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7004a2ae208SSatish Balay #undef __FUNCT__
70117bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
7024b828684SBarry Smith /*@C
70317bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
7045e76c431SBarry Smith 
705c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
706c7afd0dbSLois Curfman McInnes 
7075e42470aSBarry Smith    Input Parameters:
708c7afd0dbSLois Curfman McInnes +  snes - the SNES context
709af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
7105e42470aSBarry Smith .  x - current iterate
7115e42470aSBarry Smith .  f - residual evaluated at x
7123c632250SBarry Smith .  y - search direction
7135e42470aSBarry Smith .  w - work vector
714dc357ad2SBarry Smith .  fnorm - 2-norm of f
715dc357ad2SBarry Smith -  xnorm - norm of x if known, otherwise 0
7165e42470aSBarry Smith 
717c4a48953SLois Curfman McInnes    Output Parameters:
7187f3332b4SBarry Smith +  g - residual evaluated at new iterate w
719e106eecfSBarry Smith .  w - new iterate (x + lambda*y)
7205e42470aSBarry Smith .  gnorm - 2-norm of g
7215e42470aSBarry Smith .  ynorm - 2-norm of search length
722a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
723fee21e36SBarry Smith 
724ce9499c7SBarry Smith    Options Database Keys:
725ce9499c7SBarry Smith +  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
726ce9499c7SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
727e106eecfSBarry 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)
728e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use minlambda/ max_i ( y[i]/x[i] )
729e106eecfSBarry Smith 
7305e42470aSBarry Smith    Notes:
73117bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7325e42470aSBarry Smith 
73336851e7fSLois Curfman McInnes    Level: advanced
73436851e7fSLois Curfman McInnes 
735f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
736f59ffdedSLois Curfman McInnes 
73717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7385e42470aSBarry Smith @*/
7397087cfbeSBarry 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)
7405e76c431SBarry Smith {
741406556e6SLois Curfman McInnes   /*
742406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
743406556e6SLois Curfman McInnes      minimization problem:
744406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
745406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
746406556e6SLois Curfman McInnes    */
747e106eecfSBarry Smith   PetscReal      initslope,minlambda,lambda,lambdatemp,rellength;
748aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
749f5b6597dSBarry Smith   PetscScalar    cinitslope;
7506b5873e3SBarry Smith #endif
751dfbe8321SBarry Smith   PetscErrorCode ierr;
752a7cc72afSBarry Smith   PetscInt       count;
7534b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
754ace3abfcSBarry Smith   PetscBool      changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7555e76c431SBarry Smith 
7563a40ed3dSBarry Smith   PetscFunctionBegin;
757d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
758a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7595e76c431SBarry Smith 
7603505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
76163b9fa5eSBarry Smith   if (*ynorm == 0.0) {
76294298653SBarry Smith     if (neP->monitor) {
76394298653SBarry Smith       ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Direction and size is 0\n");CHKERRQ(ierr);
76494298653SBarry Smith     }
765b37302c6SLois Curfman McInnes     *gnorm = fnorm;
766e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
767b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
768a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
769ad922ac9SBarry Smith     goto theend2;
77094a9d846SBarry Smith   }
771e106eecfSBarry Smith   if (*ynorm > neP->maxstep) {	/* Step too big, so scale back */
772e106eecfSBarry Smith     ierr   = VecScale(y,neP->maxstep/(*ynorm));CHKERRQ(ierr);
773e106eecfSBarry Smith     *ynorm = neP->maxstep;
7745e76c431SBarry Smith   }
775ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
776e106eecfSBarry Smith   minlambda = neP->minlambda/rellength;
777a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
778aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
779a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
780329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7815e42470aSBarry Smith #else
782a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7835e42470aSBarry Smith #endif
7845e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7855e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
786f8833479SBarry Smith   ierr = PetscInfo1(snes,"Initslope %G \n",initslope);CHKERRQ(ierr);
7875e42470aSBarry Smith 
788e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
78943e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
790ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
79143e71028SBarry Smith     *flag = PETSC_FALSE;
79243e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
79343e71028SBarry Smith     goto theend2;
79443e71028SBarry Smith   }
7954936397dSBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
7964936397dSBarry Smith   if (snes->domainerror) {
7974936397dSBarry Smith     ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7984936397dSBarry Smith     PetscFunctionReturn(0);
79919717074SBarry Smith   }
800cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
80162cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
802e106eecfSBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + neP->alpha*initslope) { /* Sufficient reduction */
80394298653SBarry Smith     if (neP->monitor) {
80494298653SBarry Smith       ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line Search: Using full step\n");CHKERRQ(ierr);
80594298653SBarry Smith     }
806ad922ac9SBarry Smith     goto theend2;
8075e42470aSBarry Smith   }
8085e42470aSBarry Smith 
8095e42470aSBarry Smith   /* Fit points with quadratic */
810313b5538SBarry Smith   lambda = 1.0;
8115e42470aSBarry Smith   count = 1;
8125ca10a37SBarry Smith   while (PETSC_TRUE) {
8135e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
81494298653SBarry Smith       if (neP->monitor) {
81594298653SBarry Smith         ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: Unable to find good step length! %D \n",count);CHKERRQ(ierr);
81694298653SBarry Smith         ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"Line search: fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
81794298653SBarry Smith       }
818e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
81943e71028SBarry Smith       *flag = PETSC_FALSE;
82043e71028SBarry Smith       break;
8215e42470aSBarry Smith     }
822a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
823329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
824329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
825329f5518SBarry Smith     else                         lambda     = lambdatemp;
826275d25c3SBarry Smith 
827e68848bdSBarry Smith     ierr = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
82843e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
829ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
830d9822059SBarry 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);
831ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
83243e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
833ed98166eSMatthew Knepley       break;
834ed98166eSMatthew Knepley     }
835f5b6597dSBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
8364936397dSBarry Smith     if (snes->domainerror) {
8374936397dSBarry Smith       ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8384936397dSBarry Smith       PetscFunctionReturn(0);
83919717074SBarry Smith     }
840cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
84162cbcd01SMatthew G Knepley     if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
842e106eecfSBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*neP->alpha*initslope) { /* sufficient reduction */
84394298653SBarry Smith       if (neP->monitor) {
84494298653SBarry Smith         ierr = PetscViewerASCIIMonitorPrintf(neP->monitor,"    Line Search: Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
84594298653SBarry Smith       }
84606259719SBarry Smith       break;
8475e42470aSBarry Smith     }
8485e42470aSBarry Smith     count++;
8495e42470aSBarry Smith   }
850ad922ac9SBarry Smith   theend2:
8518f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8523c632250SBarry Smith   if (neP->postcheckstep) {
8533c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8543c632250SBarry Smith     if (changed_y) {
85579f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8563c632250SBarry Smith     }
8573c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8583c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
8594936397dSBarry Smith       if (snes->domainerror) {
8604936397dSBarry Smith         ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8614936397dSBarry Smith         PetscFunctionReturn(0);
86219717074SBarry Smith       }
8638f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
864a9f58f88SMatthew Knepley       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
8658f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
866a9f58f88SMatthew Knepley       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
86762cbcd01SMatthew G Knepley       if (PetscIsInfOrNanReal(*gnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8688f99978dSLois Curfman McInnes     }
8698f99978dSLois Curfman McInnes   }
870d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8713a40ed3dSBarry Smith   PetscFunctionReturn(0);
8725e76c431SBarry Smith }
8732343ba6eSBarry Smith 
87404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8754a2ae208SSatish Balay #undef __FUNCT__
87617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
877c9e6a524SLois Curfman McInnes /*@C
87817bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8794b27c08aSLois Curfman McInnes    by the method SNESLS.
8805e76c431SBarry Smith 
8815e76c431SBarry Smith    Input Parameters:
882c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
883af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
884c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8855e76c431SBarry Smith 
8863f9fe445SBarry Smith    Logically Collective on SNES
887fee21e36SBarry Smith 
888c4a48953SLois Curfman McInnes    Available Routines:
88917bae607SBarry Smith +  SNESLineSearchCubic() - default line search
89017bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
89117bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
89217bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8935e76c431SBarry Smith 
894c4a48953SLois Curfman McInnes     Options Database Keys:
8954b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8964b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
897e106eecfSBarry 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)
898e106eecfSBarry Smith -   -snes_ls_minlambda <minlambda> - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
899c4a48953SLois Curfman McInnes 
9005e76c431SBarry Smith    Calling sequence of func:
901bd208895SLois Curfman McInnes .vb
902ace3abfcSBarry 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)
903bd208895SLois Curfman McInnes .ve
9045e76c431SBarry Smith 
9055e76c431SBarry Smith     Input parameters for func:
906c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
907af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
9085e76c431SBarry Smith .   x - current iterate
9095e76c431SBarry Smith .   f - residual evaluated at x
9103c632250SBarry Smith .   y - search direction
911c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
9125e76c431SBarry Smith 
9135e76c431SBarry Smith     Output parameters for func:
914c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9153c632250SBarry Smith .   w - new iterate
9165e76c431SBarry Smith .   gnorm - 2-norm of g
9175e76c431SBarry Smith .   ynorm - 2-norm of search length
918a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
919f59ffdedSLois Curfman McInnes 
92036851e7fSLois Curfman McInnes     Level: advanced
92136851e7fSLois Curfman McInnes 
922f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
923f59ffdedSLois Curfman McInnes 
92417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
92517bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9265e76c431SBarry Smith @*/
9277087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal*,PetscReal*,PetscBool *),void *lsctx)
9285e76c431SBarry Smith {
9294ac538c5SBarry Smith   PetscErrorCode ierr;
93082bf6240SBarry Smith 
9313a40ed3dSBarry Smith   PetscFunctionBegin;
9324ac538c5SBarry 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);
9333a40ed3dSBarry Smith   PetscFunctionReturn(0);
9345e76c431SBarry Smith }
9358e019c35SBarry Smith 
936ace3abfcSBarry 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*/
93704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
938fb2e594dSBarry Smith EXTERN_C_BEGIN
9394a2ae208SSatish Balay #undef __FUNCT__
94017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
9417087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
94282bf6240SBarry Smith {
94382bf6240SBarry Smith   PetscFunctionBegin;
9444b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9454b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
94682bf6240SBarry Smith   PetscFunctionReturn(0);
94782bf6240SBarry Smith }
948fb2e594dSBarry Smith EXTERN_C_END
94904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9504a2ae208SSatish Balay #undef __FUNCT__
9513c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
952c8dd0c0dSLois Curfman McInnes /*@C
9533c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9544b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
955c8dd0c0dSLois Curfman McInnes 
956c8dd0c0dSLois Curfman McInnes    Input Parameters:
957c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9583c632250SBarry Smith .  func - pointer to function
959c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
960c8dd0c0dSLois Curfman McInnes 
9613f9fe445SBarry Smith    Logically Collective on SNES
962c8dd0c0dSLois Curfman McInnes 
963c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
964c8dd0c0dSLois Curfman McInnes .vb
965ace3abfcSBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscBool  *changed_y,PetscBool  *changed_w)
966c8dd0c0dSLois Curfman McInnes .ve
967b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
968b1ae27eaSLois Curfman McInnes    on failure.
969c8dd0c0dSLois Curfman McInnes 
970c8dd0c0dSLois Curfman McInnes    Input parameters for func:
971c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
972c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9733c632250SBarry Smith .  x - previous iterate
9743c632250SBarry Smith .  y - new search direction and length
9753c632250SBarry Smith -  w - current candidate iterate
976c8dd0c0dSLois Curfman McInnes 
977c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9783c632250SBarry Smith +  y - search direction (possibly changed)
9793c632250SBarry Smith .  w - current iterate (possibly modified)
9803c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9813c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
982c8dd0c0dSLois Curfman McInnes 
983c8dd0c0dSLois Curfman McInnes    Level: advanced
984c8dd0c0dSLois Curfman McInnes 
9859e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9869e247f21SBarry Smith 
9873c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9883c632250SBarry Smith 
9893c632250SBarry Smith    On input w = x + y
9903c632250SBarry Smith 
99117bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
992b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
993ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9948f99978dSLois Curfman McInnes 
99517bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
996ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
997ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
998ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9999e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
1000b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
10018f99978dSLois Curfman McInnes 
1002c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
1003c8dd0c0dSLois Curfman McInnes 
100417bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
1005c8dd0c0dSLois Curfman McInnes @*/
10067087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void *checkctx)
1007c8dd0c0dSLois Curfman McInnes {
10084ac538c5SBarry Smith   PetscErrorCode ierr;
1009c8dd0c0dSLois Curfman McInnes 
1010c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10114ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetPostCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr);
1012c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1013c8dd0c0dSLois Curfman McInnes }
101494298653SBarry Smith 
10159c3ca13aSBarry Smith #undef __FUNCT__
10169c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10179c3ca13aSBarry Smith /*@C
10189c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10197e4bb74cSBarry Smith          before the line search is called.
10209c3ca13aSBarry Smith 
10219c3ca13aSBarry Smith    Input Parameters:
10229c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10239c3ca13aSBarry Smith .  func - pointer to function
10249c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10259c3ca13aSBarry Smith 
10263f9fe445SBarry Smith    Logically Collective on SNES
10279c3ca13aSBarry Smith 
10289c3ca13aSBarry Smith    Calling sequence of func:
10299c3ca13aSBarry Smith .vb
1030ace3abfcSBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscBool  *changed_y)
10319c3ca13aSBarry Smith .ve
10329c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10339c3ca13aSBarry Smith    on failure.
10349c3ca13aSBarry Smith 
10359c3ca13aSBarry Smith    Input parameters for func:
10369c3ca13aSBarry Smith +  snes - nonlinear context
10379c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10389c3ca13aSBarry Smith .  x - previous iterate
10399c3ca13aSBarry Smith -  y - new search direction and length
10409c3ca13aSBarry Smith 
10419c3ca13aSBarry Smith    Output parameters for func:
10429c3ca13aSBarry Smith +  y - search direction (possibly changed)
10439c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10449c3ca13aSBarry Smith 
10459c3ca13aSBarry Smith    Level: advanced
10469c3ca13aSBarry Smith 
10479c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10489c3ca13aSBarry Smith 
10499c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10509c3ca13aSBarry Smith 
10517e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10529c3ca13aSBarry Smith @*/
10537087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscBool *),void *checkctx)
10549c3ca13aSBarry Smith {
10554ac538c5SBarry Smith   PetscErrorCode ierr;
10569c3ca13aSBarry Smith 
10579c3ca13aSBarry Smith   PetscFunctionBegin;
10584ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetPreCheck_C",(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscBool *),void*),(snes,func,checkctx));CHKERRQ(ierr);
10599c3ca13aSBarry Smith   PetscFunctionReturn(0);
10609c3ca13aSBarry Smith }
10619c3ca13aSBarry Smith 
106294298653SBarry Smith #undef __FUNCT__
106394298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor"
106494298653SBarry Smith /*@C
106594298653SBarry Smith    SNESLineSearchSetMonitor - Prints information about the progress or lack of progress of the line search
106694298653SBarry Smith 
106794298653SBarry Smith    Input Parameters:
106894298653SBarry Smith +  snes - nonlinear context obtained from SNESCreate()
106994298653SBarry Smith -  flg - PETSC_TRUE to monitor the line search
107094298653SBarry Smith 
10713f9fe445SBarry Smith    Logically Collective on SNES
107294298653SBarry Smith 
107394298653SBarry Smith    Options Database:
107494298653SBarry Smith .   -snes_ls_monitor
107594298653SBarry Smith 
107694298653SBarry Smith    Level: intermediate
107794298653SBarry Smith 
107894298653SBarry Smith 
107994298653SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
108094298653SBarry Smith @*/
10817087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor(SNES snes,PetscBool  flg)
108294298653SBarry Smith {
10834ac538c5SBarry Smith   PetscErrorCode ierr;
108494298653SBarry Smith 
108594298653SBarry Smith   PetscFunctionBegin;
10864ac538c5SBarry Smith   ierr = PetscTryMethod(snes,"SNESLineSearchSetMonitor_C",(SNES,PetscBool),(snes,flg));CHKERRQ(ierr);
108794298653SBarry Smith   PetscFunctionReturn(0);
108894298653SBarry Smith }
108994298653SBarry Smith 
1090c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1091ace3abfcSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscBool *,PetscBool *); /* force argument to next function to not be extern C*/
1092ace3abfcSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscBool *);                 /* force argument to next function to not be extern C*/
1093c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10944a2ae208SSatish Balay #undef __FUNCT__
10953c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10967087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1097c8dd0c0dSLois Curfman McInnes {
1098c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10993c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
11003c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1101c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1102c8dd0c0dSLois Curfman McInnes }
1103c8dd0c0dSLois Curfman McInnes EXTERN_C_END
11049c3ca13aSBarry Smith 
11059c3ca13aSBarry Smith EXTERN_C_BEGIN
11069c3ca13aSBarry Smith #undef __FUNCT__
11079c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
11087087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
11099c3ca13aSBarry Smith {
11109c3ca13aSBarry Smith   PetscFunctionBegin;
11119c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
11129c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
11139c3ca13aSBarry Smith   PetscFunctionReturn(0);
11149c3ca13aSBarry Smith }
11159c3ca13aSBarry Smith EXTERN_C_END
1116329e7e40SMatthew Knepley 
111794298653SBarry Smith EXTERN_C_BEGIN
111894298653SBarry Smith #undef __FUNCT__
111994298653SBarry Smith #define __FUNCT__ "SNESLineSearchSetMonitor_LS"
11207087cfbeSBarry Smith PetscErrorCode  SNESLineSearchSetMonitor_LS(SNES snes,PetscBool  flg)
112194298653SBarry Smith {
112294298653SBarry Smith   SNES_LS        *ls = (SNES_LS*)snes->data;
112394298653SBarry Smith   PetscErrorCode ierr;
112494298653SBarry Smith 
112594298653SBarry Smith   PetscFunctionBegin;
112694298653SBarry Smith   if (flg && !ls->monitor) {
112794298653SBarry Smith     ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,"stdout",((PetscObject)snes)->tablevel,&ls->monitor);CHKERRQ(ierr);
112894298653SBarry Smith   } else if (!flg && ls->monitor) {
112994298653SBarry Smith     ierr = PetscViewerASCIIMonitorDestroy(ls->monitor);CHKERRQ(ierr);
113094298653SBarry Smith   }
113194298653SBarry Smith   PetscFunctionReturn(0);
113294298653SBarry Smith }
113394298653SBarry Smith EXTERN_C_END
113494298653SBarry Smith 
1135329e7e40SMatthew Knepley /*
11364b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
113704d965bbSLois Curfman McInnes 
113804d965bbSLois Curfman McInnes    Input Parameters:
113904d965bbSLois Curfman McInnes .  SNES - the SNES context
114004d965bbSLois Curfman McInnes .  viewer - visualization context
114104d965bbSLois Curfman McInnes 
114204d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
114304d965bbSLois Curfman McInnes */
11444a2ae208SSatish Balay #undef __FUNCT__
11454b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11466849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1147a935fc98SLois Curfman McInnes {
11484b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11492fc52814SBarry Smith   const char     *cstr;
1150dfbe8321SBarry Smith   PetscErrorCode ierr;
1151ace3abfcSBarry Smith   PetscBool      iascii;
1152a935fc98SLois Curfman McInnes 
11533a40ed3dSBarry Smith   PetscFunctionBegin;
11542692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
115532077d6dSBarry Smith   if (iascii) {
115617bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
115717bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
115817bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
115919bcc07fSBarry Smith     else                                                cstr = "unknown";
1160b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1161e106eecfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, minlambda=%G\n",ls->alpha,ls->maxstep,ls->minlambda);CHKERRQ(ierr);
11625cd90555SBarry Smith   } else {
1163e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
116419bcc07fSBarry Smith   }
11653a40ed3dSBarry Smith   PetscFunctionReturn(0);
1166a935fc98SLois Curfman McInnes }
116704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
116804d965bbSLois Curfman McInnes /*
11694b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
117004d965bbSLois Curfman McInnes 
117104d965bbSLois Curfman McInnes    Input Parameter:
117204d965bbSLois Curfman McInnes .  snes - the SNES context
117304d965bbSLois Curfman McInnes 
117404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
117504d965bbSLois Curfman McInnes */
11764a2ae208SSatish Balay #undef __FUNCT__
11774b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11786849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11795e42470aSBarry Smith {
11804b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1181e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1182dfbe8321SBarry Smith   PetscErrorCode ierr;
1183a7cc72afSBarry Smith   PetscInt       indx;
1184ace3abfcSBarry Smith   PetscBool      flg,set;
11855e42470aSBarry Smith 
11863a40ed3dSBarry Smith   PetscFunctionBegin;
1187b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11884b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11894b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
1190e106eecfSBarry Smith     ierr = PetscOptionsReal("-snes_ls_minlambda","Minimum lambda allowed","None",ls->minlambda,&ls->minlambda,0);CHKERRQ(ierr);
1191acfcf0e5SJed Brown     ierr = PetscOptionsBool("-snes_ls_monitor","Print progress of line searches","SNESLineSearchSetMonitor",ls->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
119294298653SBarry Smith     if (set) {ierr = SNESLineSearchSetMonitor(snes,flg);CHKERRQ(ierr);}
1193186905e3SBarry Smith 
119417bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
119525cdf11fSBarry Smith     if (flg) {
119622e36657SBarry Smith       switch (indx) {
1197b49fd9e1SBarry Smith       case 0:
119817bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1199b49fd9e1SBarry Smith         break;
1200b49fd9e1SBarry Smith       case 1:
120117bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1202b49fd9e1SBarry Smith         break;
1203b49fd9e1SBarry Smith       case 2:
120417bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1205b49fd9e1SBarry Smith         break;
1206b49fd9e1SBarry Smith       case 3:
120717bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1208b49fd9e1SBarry Smith         break;
12095e42470aSBarry Smith       }
12105e42470aSBarry Smith     }
1211b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
12123a40ed3dSBarry Smith   PetscFunctionReturn(0);
12135e42470aSBarry Smith }
1214*4827ddcaSBarry Smith 
1215*4827ddcaSBarry Smith EXTERN_C_BEGIN
1216*4827ddcaSBarry Smith extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
1217*4827ddcaSBarry Smith EXTERN_C_END
1218*4827ddcaSBarry Smith 
121904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1220ebe3b25bSBarry Smith /*MC
1221ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
122204d965bbSLois Curfman McInnes 
1223ebe3b25bSBarry Smith    Options Database:
1224ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1225ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1226e106eecfSBarry 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)
1227acbee50cSBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
1228acbee50cSBarry Smith -   -snes_ls_monitor - print information about progress of line searches
1229acbee50cSBarry Smith 
123004d965bbSLois Curfman McInnes 
1231ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
123204d965bbSLois Curfman McInnes 
1233ee3001cbSBarry Smith    Level: beginner
1234ee3001cbSBarry Smith 
123517bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
123617bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
1237b3dd4ab5SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
1238ebe3b25bSBarry Smith 
1239ebe3b25bSBarry Smith M*/
1240fb2e594dSBarry Smith EXTERN_C_BEGIN
12414a2ae208SSatish Balay #undef __FUNCT__
12424b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
12437087cfbeSBarry Smith PetscErrorCode  SNESCreate_LS(SNES snes)
12445e42470aSBarry Smith {
1245dfbe8321SBarry Smith   PetscErrorCode ierr;
12464b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12475e42470aSBarry Smith 
12483a40ed3dSBarry Smith   PetscFunctionBegin;
1249e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1250e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1251e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1252e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1253e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
12545e42470aSBarry Smith 
125538f2d2fdSLisandro Dalcin   ierr                  = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
12565e42470aSBarry Smith   snes->data    	= (void*)neP;
12575e42470aSBarry Smith   neP->alpha		= 1.e-4;
12585e42470aSBarry Smith   neP->maxstep		= 1.e8;
1259e106eecfSBarry Smith   neP->minlambda        = 1.e-12;
126017bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1261c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12623c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12633c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12643c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12653c632250SBarry Smith   neP->precheck         = PETSC_NULL;
126682bf6240SBarry Smith 
126794298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetMonitor_C","SNESLineSearchSetMonitor_LS",SNESLineSearchSetMonitor_LS);CHKERRQ(ierr);
1268*4827ddcaSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetParams_C","SNESLineSearchSetParams_LS",SNESLineSearchSetParams_LS);CHKERRQ(ierr);
126994298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
127094298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
127194298653SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
127282bf6240SBarry Smith 
12733a40ed3dSBarry Smith   PetscFunctionReturn(0);
12745e42470aSBarry Smith }
1275fb2e594dSBarry Smith EXTERN_C_END
127682bf6240SBarry Smith 
127782bf6240SBarry Smith 
127882bf6240SBarry Smith 
1279