xref: /petsc/src/snes/impls/ls/ls.c (revision 3f1495945e0c17d25f364df1e0c912f77d59eb1e)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
25e76c431SBarry Smith 
370f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
45e42470aSBarry Smith 
58a5d9ee4SBarry Smith /*
68a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
78a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
836669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
936669109SBarry Smith     for this trick.
108a5d9ee4SBarry Smith */
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
13dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
148a5d9ee4SBarry Smith {
158a5d9ee4SBarry Smith   PetscReal      a1;
16dfbe8321SBarry Smith   PetscErrorCode ierr;
1736669109SBarry Smith   PetscTruth     hastranspose;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2136669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2236669109SBarry Smith   if (hastranspose) {
238a5d9ee4SBarry Smith     /* Compute || J^T F|| */
248a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
258a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
26ae15b995SBarry Smith     ierr = PetscInfo1(0,"|| J^T F|| %G near zero implies found a local minimum\n",a1/fnorm);CHKERRQ(ierr);
278a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2836669109SBarry Smith   } else {
2936669109SBarry Smith     Vec         work;
30ea709b57SSatish Balay     PetscScalar result;
3136669109SBarry Smith     PetscReal   wnorm;
3236669109SBarry Smith 
33abb0e124SMatthew Knepley     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3736669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
3836669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
3936669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
40ae15b995SBarry Smith     ierr = PetscInfo1(0,"(F^T J random)/(|| F ||*||random|| %G near zero implies found a local minimum\n",a1);CHKERRQ(ierr);
4136669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4236669109SBarry Smith   }
438a5d9ee4SBarry Smith   PetscFunctionReturn(0);
448a5d9ee4SBarry Smith }
458a5d9ee4SBarry Smith 
4674637425SBarry Smith /*
475ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4874637425SBarry Smith */
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
51dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5274637425SBarry Smith {
5374637425SBarry Smith   PetscReal      a1,a2;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
5574637425SBarry Smith   PetscTruth     hastranspose;
5674637425SBarry Smith 
5774637425SBarry Smith   PetscFunctionBegin;
5874637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5974637425SBarry Smith   if (hastranspose) {
6074637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6179f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6274637425SBarry Smith 
6374637425SBarry Smith     /* Compute || J^T W|| */
6474637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6574637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6674637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
679e247f21SBarry Smith     if (a1) {
68ae15b995SBarry Smith       ierr = PetscInfo1(0,"||J^T(F-Ax)||/||F-AX|| %G near zero implies inconsistent rhs\n",a2/a1);CHKERRQ(ierr);
6985471664SBarry Smith     }
7074637425SBarry Smith   }
7174637425SBarry Smith   PetscFunctionReturn(0);
7274637425SBarry Smith }
7374637425SBarry Smith 
7404d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7504d965bbSLois Curfman McInnes 
7604d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7794b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7804d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
7904d965bbSLois Curfman McInnes      respectively.
8004d965bbSLois Curfman McInnes 
81fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8204d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8304d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
84fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8504d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8604d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
874b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
884b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
8904d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9004d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9104d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9204d965bbSLois Curfman McInnes      for all nonlinear solvers.
9304d965bbSLois Curfman McInnes 
9404d965bbSLois Curfman McInnes      Another key routine is:
9504d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9604d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9704d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9804d965bbSLois Curfman McInnes      SNESSolve() if necessary.
9904d965bbSLois Curfman McInnes 
10004d965bbSLois Curfman McInnes      Additional basic routines are:
10104d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10204d965bbSLois Curfman McInnes                                       have actually been used.
10304d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
104186905e3SBarry Smith      SNESView().
10504d965bbSLois Curfman McInnes 
10604d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10704d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10804d965bbSLois Curfman McInnes      above description applies to these categories also.
10904d965bbSLois Curfman McInnes 
11004d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1115e42470aSBarry Smith /*
1124b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11304d965bbSLois Curfman McInnes    method with a line search.
1145e76c431SBarry Smith 
11504d965bbSLois Curfman McInnes    Input Parameters:
11604d965bbSLois Curfman McInnes .  snes - the SNES context
1175e76c431SBarry Smith 
11804d965bbSLois Curfman McInnes    Output Parameter:
119c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12004d965bbSLois Curfman McInnes 
12104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1225e76c431SBarry Smith 
1235e76c431SBarry Smith    Notes:
1245e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1255e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1265e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1275e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
128393d2d9aSLois Curfman McInnes    and Schnabel.
1295e76c431SBarry Smith */
1304a2ae208SSatish Balay #undef __FUNCT__
1314b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
132dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1335e76c431SBarry Smith {
1344b27c08aSLois Curfman McInnes   SNES_LS            *neP = (SNES_LS*)snes->data;
1356849ba73SBarry Smith   PetscErrorCode     ierr;
136a7cc72afSBarry Smith   PetscInt           maxits,i,lits;
137a7cc72afSBarry Smith   PetscTruth         lssucceed;
138112a2221SBarry Smith   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
139329f5518SBarry Smith   PetscReal          fnorm,gnorm,xnorm,ynorm;
1405e42470aSBarry Smith   Vec                Y,X,F,G,W,TMP;
1413d4c4710SBarry Smith   KSPConvergedReason kspreason;
1425e76c431SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
1443d4c4710SBarry Smith   snes->numFailures            = 0;
1453d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
146184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
147184914b5SBarry Smith 
1485e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1495e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
15039e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1515e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1525e42470aSBarry Smith   G		= snes->work[1];
1535e42470aSBarry Smith   W		= snes->work[2];
1545e76c431SBarry Smith 
1554c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1564c49b128SBarry Smith   snes->iter = 0;
1574c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15819717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
159cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
16043e71028SBarry Smith   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1610f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1625e42470aSBarry Smith   snes->norm = fnorm;
1630f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16484c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16594a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
166*3f149594SLisandro Dalcin 
167*3f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
168*3f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
169*3f149594SLisandro Dalcin   if (snes->ops->converged) {
170e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
171*3f149594SLisandro Dalcin   }
17206ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
173d034289bSBarry Smith 
1745e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1755e76c431SBarry Smith 
176329e7e40SMatthew Knepley     /* Call general purpose update function */
177e7788613SBarry Smith     if (snes->ops->update) {
178e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
179de8cb200SMatthew Knepley     }
180329e7e40SMatthew Knepley 
181ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1825334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18394b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18471f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
1853d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1863d4c4710SBarry Smith     if (kspreason < 0) {
1873d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
1883d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
1893d4c4710SBarry Smith         PetscFunctionReturn(0);
1903d4c4710SBarry Smith       }
1913d4c4710SBarry Smith     }
1923d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
193*3f149594SLisandro Dalcin     snes->linear_its += lits;
194*3f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
19574637425SBarry Smith 
1969c3ca13aSBarry Smith     if (neP->precheckstep) {
1979c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
1989c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
1999c3ca13aSBarry Smith     }
2009c3ca13aSBarry Smith 
201b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
20274637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
20385471664SBarry Smith     }
20474637425SBarry Smith 
205ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
206ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
207e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
208ea4d3ed3SLois Curfman McInnes     */
20981b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
210*3f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
211a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
212ae15b995SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
213a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
2143c632250SBarry Smith     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
215a3b891d8SBarry Smith     fnorm = gnorm;
2164a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
217a3b891d8SBarry Smith 
2185ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2195ed2d596SBarry Smith     snes->iter = i+1;
2205ed2d596SBarry Smith     snes->norm = fnorm;
2215ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2225ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2235ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2245ed2d596SBarry Smith 
225a7cc72afSBarry Smith     if (!lssucceed) {
2268a5d9ee4SBarry Smith       PetscTruth ismin;
22750ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2283505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2298a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2308a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2313505fcc1SBarry Smith         break;
2323505fcc1SBarry Smith       }
23350ffb88aSMatthew Knepley     }
2345e76c431SBarry Smith 
2355e76c431SBarry Smith     /* Test for convergence */
236e7788613SBarry Smith     if (snes->ops->converged) {
23729e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
238e7788613SBarry Smith       ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
239*3f149594SLisandro Dalcin       if (snes->reason) break;
24016c95cacSBarry Smith     }
24129e0b56fSBarry Smith   }
24239e2f89bSBarry Smith   if (X != snes->vec_sol) {
243393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
244e15f7bb5SBarry Smith   }
245e15f7bb5SBarry Smith   if (F != snes->vec_func) {
246e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
247e15f7bb5SBarry Smith   }
24839e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
24939e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
25052392280SLois Curfman McInnes   if (i == maxits) {
251ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
2523505fcc1SBarry Smith     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;
279e74804ceSBarry Smith   if (!snes->work) {
28081b6cf68SLois Curfman McInnes     snes->nwork = 4;
281d7e8b826SBarry Smith     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
282efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
28381b6cf68SLois Curfman McInnes     snes->vec_sol_update_always = snes->work[3];
284e74804ceSBarry Smith   }
2853a40ed3dSBarry Smith   PetscFunctionReturn(0);
2865e76c431SBarry Smith }
28704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
28804d965bbSLois Curfman McInnes /*
2894b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2904b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29104d965bbSLois Curfman McInnes 
29204d965bbSLois Curfman McInnes    Input Parameter:
29304d965bbSLois Curfman McInnes .  snes - the SNES context
29404d965bbSLois Curfman McInnes 
295de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29604d965bbSLois Curfman McInnes  */
2974a2ae208SSatish Balay #undef __FUNCT__
2984b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
299dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3005e76c431SBarry Smith {
301dfbe8321SBarry Smith   PetscErrorCode ierr;
3023a40ed3dSBarry Smith 
3033a40ed3dSBarry Smith   PetscFunctionBegin;
3045baf8537SBarry Smith   if (snes->nwork) {
3054b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
30621c89e3eSBarry Smith   }
3075c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
308557d3f75SLisandro Dalcin 
309557d3f75SLisandro Dalcin   /* clear composed functions */
310557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","",PETSC_NULL);CHKERRQ(ierr);
311557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","",PETSC_NULL);CHKERRQ(ierr);
312557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","",PETSC_NULL);CHKERRQ(ierr);
313557d3f75SLisandro Dalcin 
3143a40ed3dSBarry Smith   PetscFunctionReturn(0);
3155e76c431SBarry Smith }
31604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3174a2ae208SSatish Balay #undef __FUNCT__
31817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
31982bf6240SBarry Smith 
3204b828684SBarry Smith /*@C
32117bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3225e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3235e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3245e76c431SBarry Smith 
325c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
326c7afd0dbSLois Curfman McInnes 
3275e76c431SBarry Smith    Input Parameters:
328c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
329af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3305e76c431SBarry Smith .  x - current iterate
3315e76c431SBarry Smith .  f - residual evaluated at x
3323c632250SBarry Smith .  y - search direction
3335e76c431SBarry Smith .  w - work vector
334c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3355e76c431SBarry Smith 
336c4a48953SLois Curfman McInnes    Output Parameters:
337c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3383c632250SBarry Smith .  w - new iterate
3395e76c431SBarry Smith .  gnorm - 2-norm of g
3405e76c431SBarry Smith .  ynorm - 2-norm of search length
341a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
342fee21e36SBarry Smith 
343c4a48953SLois Curfman McInnes    Options Database Key:
34417bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
345c4a48953SLois Curfman McInnes 
34636851e7fSLois Curfman McInnes    Level: advanced
34736851e7fSLois Curfman McInnes 
34828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
34928ae5a21SLois Curfman McInnes 
35017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
35117bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3525e76c431SBarry Smith @*/
35317bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNo(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
3545e76c431SBarry Smith {
355dfbe8321SBarry Smith   PetscErrorCode ierr;
3564b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3573c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3585334005bSBarry Smith 
3593a40ed3dSBarry Smith   PetscFunctionBegin;
360a7cc72afSBarry Smith   *flag = PETSC_TRUE;
361d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
362a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
36379f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3643c632250SBarry Smith   if (neP->postcheckstep) {
3653c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3668f99978dSLois Curfman McInnes   }
3673c632250SBarry Smith   if (changed_y) {
36879f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3693c632250SBarry Smith   }
3703c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
371d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
37219717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
37319717074SBarry Smith   }
374d5e45103SBarry Smith   CHKERRQ(ierr);
375d5e45103SBarry Smith 
376a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
37743e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
378d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3793a40ed3dSBarry Smith   PetscFunctionReturn(0);
3805e76c431SBarry Smith }
38104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3824a2ae208SSatish Balay #undef __FUNCT__
38317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
38482bf6240SBarry Smith 
38529e0b56fSBarry Smith /*@C
38617bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
38729e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
38829e0b56fSBarry Smith    even compute the norm of the function or search direction; this
38929e0b56fSBarry Smith    is intended only when you know the full step is fine and are
39029e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
391c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
392c7afd0dbSLois Curfman McInnes 
393c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
39429e0b56fSBarry Smith 
39529e0b56fSBarry Smith    Input Parameters:
396c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
397af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
39829e0b56fSBarry Smith .  x - current iterate
39929e0b56fSBarry Smith .  f - residual evaluated at x
4003c632250SBarry Smith .  y - search direction
40129e0b56fSBarry Smith .  w - work vector
402c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
40329e0b56fSBarry Smith 
40429e0b56fSBarry Smith    Output Parameters:
405c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4063c632250SBarry Smith .  w - new iterate
40729e0b56fSBarry Smith .  gnorm - not changed
40829e0b56fSBarry Smith .  ynorm - not changed
409a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
410fee21e36SBarry Smith 
41129e0b56fSBarry Smith    Options Database Key:
41217bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
41329e0b56fSBarry Smith 
4148cbba510SBarry Smith    Notes:
41517bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
416ea56f5baSLois Curfman McInnes    either the options
417ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
418ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
419329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
420329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
421329f5518SBarry Smith 
422329f5518SBarry Smith    During the final iteration this will not evaluate the function at
423329f5518SBarry Smith    the solution point. This is to save a function evaluation while
424329f5518SBarry Smith    using pseudo-timestepping.
4258cbba510SBarry Smith 
426ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
427a6570f20SBarry Smith    SNESMonitorDefault() (as activated via -snes_monitor) will not be
428ea56f5baSLois Curfman McInnes    correct, since they are not computed.
429ea56f5baSLois Curfman McInnes 
430ea56f5baSLois Curfman McInnes    Level: advanced
4318cbba510SBarry Smith 
43229e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
43329e0b56fSBarry Smith 
43417bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
43517bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
43629e0b56fSBarry Smith @*/
43717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
43829e0b56fSBarry Smith {
439dfbe8321SBarry Smith   PetscErrorCode ierr;
4404b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4413c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
44229e0b56fSBarry Smith 
4433a40ed3dSBarry Smith   PetscFunctionBegin;
444a7cc72afSBarry Smith   *flag = PETSC_TRUE;
445d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
44679f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4473c632250SBarry Smith   if (neP->postcheckstep) {
4483c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4493c632250SBarry Smith   }
4503c632250SBarry Smith   if (changed_y) {
45179f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4528f99978dSLois Curfman McInnes   }
453329f5518SBarry Smith 
454329f5518SBarry Smith   /* don't evaluate function the last time through */
455329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4563c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
457d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
45819717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
45919717074SBarry Smith     }
460d5e45103SBarry Smith     CHKERRQ(ierr);
461329f5518SBarry Smith   }
462d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4633a40ed3dSBarry Smith   PetscFunctionReturn(0);
46429e0b56fSBarry Smith }
46504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4664a2ae208SSatish Balay #undef __FUNCT__
46717bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4684b828684SBarry Smith /*@C
46917bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4705e76c431SBarry Smith 
471c7afd0dbSLois Curfman McInnes    Collective on SNES
472c7afd0dbSLois Curfman McInnes 
4735e76c431SBarry Smith    Input Parameters:
474c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
475af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4765e76c431SBarry Smith .  x - current iterate
4775e76c431SBarry Smith .  f - residual evaluated at x
4783c632250SBarry Smith .  y - search direction
4795e76c431SBarry Smith .  w - work vector
480c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4815e76c431SBarry Smith 
482393d2d9aSLois Curfman McInnes    Output Parameters:
483c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4843c632250SBarry Smith .  w - new iterate
4855e76c431SBarry Smith .  gnorm - 2-norm of g
4865e76c431SBarry Smith .  ynorm - 2-norm of search length
487a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
488fee21e36SBarry Smith 
489c4a48953SLois Curfman McInnes    Options Database Key:
49017bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
491c4a48953SLois Curfman McInnes 
4925e76c431SBarry Smith    Notes:
4935e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4945e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4955e76c431SBarry Smith 
49636851e7fSLois Curfman McInnes    Level: advanced
49736851e7fSLois Curfman McInnes 
49828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49928ae5a21SLois Curfman McInnes 
50017bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
5015e76c431SBarry Smith @*/
50217bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchCubic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
5035e76c431SBarry Smith {
504406556e6SLois Curfman McInnes   /*
505406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
506406556e6SLois Curfman McInnes      minimization problem:
507406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
508406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
509406556e6SLois Curfman McInnes    */
510406556e6SLois Curfman McInnes 
5115ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
512275d25c3SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
513aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
514ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5156b5873e3SBarry Smith #endif
516dfbe8321SBarry Smith   PetscErrorCode ierr;
517a7cc72afSBarry Smith   PetscInt       count;
5184b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
51979f36460SBarry Smith   PetscScalar    scale;
5203c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5215e76c431SBarry Smith 
5223a40ed3dSBarry Smith   PetscFunctionBegin;
523d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
524a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5255e76c431SBarry Smith   alpha   = neP->alpha;
5265e76c431SBarry Smith   maxstep = neP->maxstep;
5275e76c431SBarry Smith   steptol = neP->steptol;
5285e76c431SBarry Smith 
529cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5309e247f21SBarry Smith   if (!*ynorm) {
531ae15b995SBarry Smith     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
532a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5333c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
534a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
535a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
536ad922ac9SBarry Smith     goto theend1;
53794a9d846SBarry Smith   }
5385e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5395e42470aSBarry Smith     scale = maxstep/(*ynorm);
540ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr);
5412dcb1b2aSMatthew Knepley     ierr = VecScale(y,scale);CHKERRQ(ierr);
5425e76c431SBarry Smith     *ynorm = maxstep;
5435e76c431SBarry Smith   }
544ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5455ca10a37SBarry Smith   minlambda = steptol/rellength;
546a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
547aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
548a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
549329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5505e42470aSBarry Smith #else
551a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5525e42470aSBarry Smith #endif
5535e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5545e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5555e76c431SBarry Smith 
556e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
55743e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
558ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
55943e71028SBarry Smith     *flag = PETSC_FALSE;
56043e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
56143e71028SBarry Smith     goto theend1;
56243e71028SBarry Smith   }
56319717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
564d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
56519717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
56619717074SBarry Smith   }
567d5e45103SBarry Smith   CHKERRQ(ierr);
568313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56943e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5705d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
571ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
572ad922ac9SBarry Smith     goto theend1;
5735e76c431SBarry Smith   }
5745e76c431SBarry Smith 
5755e76c431SBarry Smith   /* Fit points with quadratic */
576313b5538SBarry Smith   lambda     = 1.0;
577a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5785e76c431SBarry Smith   lambdaprev = lambda;
5795e76c431SBarry Smith   gnormprev  = *gnorm;
580329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
581ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
582ddd12767SBarry Smith   else                         lambda = lambdatemp;
583275d25c3SBarry Smith 
584aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
585a23f76dfSSatish Balay   clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
5865e42470aSBarry Smith #else
587e68848bdSBarry Smith   ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
5885e42470aSBarry Smith #endif
58943e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
590ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
59143e71028SBarry Smith     *flag = PETSC_FALSE;
59243e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
59343e71028SBarry Smith     goto theend1;
59443e71028SBarry Smith   }
59519717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
596d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59719717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59819717074SBarry Smith   }
599d5e45103SBarry Smith   CHKERRQ(ierr);
600cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
60143e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6025ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
603ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
604ad922ac9SBarry Smith     goto theend1;
6055e76c431SBarry Smith   }
6065e76c431SBarry Smith 
6075e76c431SBarry Smith   /* Fit points with cubic */
6085e76c431SBarry Smith   count = 1;
6098229bfc2SMatthew Knepley   while (count < 10000) {
6105e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
611ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
612ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
61343e71028SBarry Smith       *flag = PETSC_FALSE;
61443e71028SBarry Smith       break;
6155e76c431SBarry Smith     }
6165d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6175d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
618ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6192b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6205e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6215e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6225e76c431SBarry Smith     if (a == 0.0) {
6235e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6245e76c431SBarry Smith     } else {
6255e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6265e76c431SBarry Smith     }
6275e76c431SBarry Smith     lambdaprev = lambda;
6285e76c431SBarry Smith     gnormprev  = *gnorm;
629329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
630329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6315e76c431SBarry Smith     else                         lambda     = lambdatemp;
632aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
633275d25c3SBarry Smith     clambda   = lambda;
634e68848bdSBarry Smith     ierr      = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
6355e42470aSBarry Smith #else
636e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
6375e42470aSBarry Smith #endif
63843e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
639ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
640ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
641ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
64243e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
643ed98166eSMatthew Knepley       break;
644ed98166eSMatthew Knepley     }
64519717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
646d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64719717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64819717074SBarry Smith     }
649d5e45103SBarry Smith     CHKERRQ(ierr);
650cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
65143e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6525ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
653ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
65443e71028SBarry Smith       break;
6552b022350SLois Curfman McInnes     } else {
656ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);CHKERRQ(ierr);
6575e76c431SBarry Smith     }
6585e76c431SBarry Smith     count++;
6595e76c431SBarry Smith   }
6608229bfc2SMatthew Knepley   if (count >= 10000) {
661cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6628229bfc2SMatthew Knepley   }
663ad922ac9SBarry Smith   theend1:
6648f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6653c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6663c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6673c632250SBarry Smith     if (changed_y) {
66879f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6693c632250SBarry Smith     }
6703c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6713c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
672d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
67319717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
67419717074SBarry Smith       }
675d5e45103SBarry Smith       CHKERRQ(ierr);
6768f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67743e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6783c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6798f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6803c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6818f99978dSLois Curfman McInnes     }
6828f99978dSLois Curfman McInnes   }
683d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6843a40ed3dSBarry Smith   PetscFunctionReturn(0);
6855e76c431SBarry Smith }
68604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6874a2ae208SSatish Balay #undef __FUNCT__
68817bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6894b828684SBarry Smith /*@C
69017bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6915e76c431SBarry Smith 
692c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
693c7afd0dbSLois Curfman McInnes 
6945e42470aSBarry Smith    Input Parameters:
695c7afd0dbSLois Curfman McInnes +  snes - the SNES context
696af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6975e42470aSBarry Smith .  x - current iterate
6985e42470aSBarry Smith .  f - residual evaluated at x
6993c632250SBarry Smith .  y - search direction
7005e42470aSBarry Smith .  w - work vector
701c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
7025e42470aSBarry Smith 
703c4a48953SLois Curfman McInnes    Output Parameters:
7047f3332b4SBarry Smith +  g - residual evaluated at new iterate w
7057f3332b4SBarry Smith .  w - new iterate (x + alpha*y)
7065e42470aSBarry Smith .  gnorm - 2-norm of g
7075e42470aSBarry Smith .  ynorm - 2-norm of search length
708a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
709fee21e36SBarry Smith 
710c4a48953SLois Curfman McInnes    Options Database Key:
71117bae607SBarry Smith .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
712c4a48953SLois Curfman McInnes 
7135e42470aSBarry Smith    Notes:
71417bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7155e42470aSBarry Smith 
71636851e7fSLois Curfman McInnes    Level: advanced
71736851e7fSLois Curfman McInnes 
718f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
719f59ffdedSLois Curfman McInnes 
72017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7215e42470aSBarry Smith @*/
72217bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchQuadratic(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
7235e76c431SBarry Smith {
724406556e6SLois Curfman McInnes   /*
725406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
726406556e6SLois Curfman McInnes      minimization problem:
727406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
728406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
729406556e6SLois Curfman McInnes    */
730275d25c3SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
731aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
732ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7336b5873e3SBarry Smith #endif
734dfbe8321SBarry Smith   PetscErrorCode ierr;
735a7cc72afSBarry Smith   PetscInt       count;
7364b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
73779f36460SBarry Smith   PetscScalar    scale;
7383c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7395e76c431SBarry Smith 
7403a40ed3dSBarry Smith   PetscFunctionBegin;
741d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
742a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7435e42470aSBarry Smith   alpha   = neP->alpha;
7445e42470aSBarry Smith   maxstep = neP->maxstep;
7455e42470aSBarry Smith   steptol = neP->steptol;
7465e76c431SBarry Smith 
7473505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74863b9fa5eSBarry Smith   if (*ynorm == 0.0) {
749ae15b995SBarry Smith     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
750b37302c6SLois Curfman McInnes     *gnorm = fnorm;
751e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
752b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
753a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
754ad922ac9SBarry Smith     goto theend2;
75594a9d846SBarry Smith   }
7565e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7575e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7582dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7595e42470aSBarry Smith     *ynorm = maxstep;
7605e76c431SBarry Smith   }
761ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7625ca10a37SBarry Smith   minlambda = steptol/rellength;
763a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
764aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
765a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
766329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7675e42470aSBarry Smith #else
768a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7695e42470aSBarry Smith #endif
7705e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7715e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7725e42470aSBarry Smith 
773e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
77443e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
775ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
77643e71028SBarry Smith     *flag = PETSC_FALSE;
77743e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77843e71028SBarry Smith     goto theend2;
77943e71028SBarry Smith   }
78019717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
781d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
78219717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
78319717074SBarry Smith   }
784d5e45103SBarry Smith   CHKERRQ(ierr);
785cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
78643e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7875d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
788ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
789ad922ac9SBarry Smith     goto theend2;
7905e42470aSBarry Smith   }
7915e42470aSBarry Smith 
7925e42470aSBarry Smith   /* Fit points with quadratic */
793313b5538SBarry Smith   lambda = 1.0;
7945e42470aSBarry Smith   count = 1;
7955ca10a37SBarry Smith   while (PETSC_TRUE) {
7965e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
797ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
798ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
799e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
80043e71028SBarry Smith       *flag = PETSC_FALSE;
80143e71028SBarry Smith       break;
8025e42470aSBarry Smith     }
803a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
804329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
805329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
806329f5518SBarry Smith     else                         lambda     = lambdatemp;
807275d25c3SBarry Smith 
808aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
809e68848bdSBarry Smith     clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
8105e42470aSBarry Smith #else
811e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8125e42470aSBarry Smith #endif
81343e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
814ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
815ae15b995SBarry 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);
816ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
81743e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
818ed98166eSMatthew Knepley       break;
819ed98166eSMatthew Knepley     }
82019717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
821d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
82219717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
82319717074SBarry Smith     }
824d5e45103SBarry Smith     CHKERRQ(ierr);
825cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82643e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8275ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
828ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
82906259719SBarry Smith       break;
8305e42470aSBarry Smith     }
8315e42470aSBarry Smith     count++;
8325e42470aSBarry Smith   }
833ad922ac9SBarry Smith   theend2:
8348f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8353c632250SBarry Smith   if (neP->postcheckstep) {
8363c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8373c632250SBarry Smith     if (changed_y) {
83879f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8393c632250SBarry Smith     }
8403c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8413c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
842d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
84319717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
84419717074SBarry Smith       }
845d5e45103SBarry Smith       CHKERRQ(ierr);
8468f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8473c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8488f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8493c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8503c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8518f99978dSLois Curfman McInnes     }
8528f99978dSLois Curfman McInnes   }
853d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8543a40ed3dSBarry Smith   PetscFunctionReturn(0);
8555e76c431SBarry Smith }
8562343ba6eSBarry Smith 
85704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8584a2ae208SSatish Balay #undef __FUNCT__
85917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
860c9e6a524SLois Curfman McInnes /*@C
86117bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8624b27c08aSLois Curfman McInnes    by the method SNESLS.
8635e76c431SBarry Smith 
8645e76c431SBarry Smith    Input Parameters:
865c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
866af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
867c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8685e76c431SBarry Smith 
869fee21e36SBarry Smith    Collective on SNES
870fee21e36SBarry Smith 
871c4a48953SLois Curfman McInnes    Available Routines:
87217bae607SBarry Smith +  SNESLineSearchCubic() - default line search
87317bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
87417bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
87517bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8765e76c431SBarry Smith 
877c4a48953SLois Curfman McInnes     Options Database Keys:
8784b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8794b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8804b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8814b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8823304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8833304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
884c4a48953SLois Curfman McInnes 
8855e76c431SBarry Smith    Calling sequence of func:
886bd208895SLois Curfman McInnes .vb
887af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
888a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
889bd208895SLois Curfman McInnes .ve
8905e76c431SBarry Smith 
8915e76c431SBarry Smith     Input parameters for func:
892c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
893af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8945e76c431SBarry Smith .   x - current iterate
8955e76c431SBarry Smith .   f - residual evaluated at x
8963c632250SBarry Smith .   y - search direction
8975e76c431SBarry Smith .   w - work vector
898c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8995e76c431SBarry Smith 
9005e76c431SBarry Smith     Output parameters for func:
901c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
9023c632250SBarry Smith .   w - new iterate
9035e76c431SBarry Smith .   gnorm - 2-norm of g
9045e76c431SBarry Smith .   ynorm - 2-norm of search length
905a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
906f59ffdedSLois Curfman McInnes 
90736851e7fSLois Curfman McInnes     Level: advanced
90836851e7fSLois Curfman McInnes 
909f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
910f59ffdedSLois Curfman McInnes 
91117bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
91217bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9135e76c431SBarry Smith @*/
91417bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9155e76c431SBarry Smith {
916a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
91782bf6240SBarry Smith 
9183a40ed3dSBarry Smith   PetscFunctionBegin;
91917bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
92082bf6240SBarry Smith   if (f) {
921af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
92282bf6240SBarry Smith   }
9233a40ed3dSBarry Smith   PetscFunctionReturn(0);
9245e76c431SBarry Smith }
9258e019c35SBarry Smith 
926a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
92704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
928fb2e594dSBarry Smith EXTERN_C_BEGIN
9294a2ae208SSatish Balay #undef __FUNCT__
93017bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
93117bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
93282bf6240SBarry Smith {
93382bf6240SBarry Smith   PetscFunctionBegin;
9344b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9354b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
93682bf6240SBarry Smith   PetscFunctionReturn(0);
93782bf6240SBarry Smith }
938fb2e594dSBarry Smith EXTERN_C_END
93904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9404a2ae208SSatish Balay #undef __FUNCT__
9413c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
942c8dd0c0dSLois Curfman McInnes /*@C
9433c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9444b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
945c8dd0c0dSLois Curfman McInnes 
946c8dd0c0dSLois Curfman McInnes    Input Parameters:
947c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9483c632250SBarry Smith .  func - pointer to function
949c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
950c8dd0c0dSLois Curfman McInnes 
951c8dd0c0dSLois Curfman McInnes    Collective on SNES
952c8dd0c0dSLois Curfman McInnes 
953c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
954c8dd0c0dSLois Curfman McInnes .vb
9553c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
956c8dd0c0dSLois Curfman McInnes .ve
957b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
958b1ae27eaSLois Curfman McInnes    on failure.
959c8dd0c0dSLois Curfman McInnes 
960c8dd0c0dSLois Curfman McInnes    Input parameters for func:
961c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
962c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9633c632250SBarry Smith .  x - previous iterate
9643c632250SBarry Smith .  y - new search direction and length
9653c632250SBarry Smith -  w - current candidate iterate
966c8dd0c0dSLois Curfman McInnes 
967c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9683c632250SBarry Smith +  y - search direction (possibly changed)
9693c632250SBarry Smith .  w - current iterate (possibly modified)
9703c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9713c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
972c8dd0c0dSLois Curfman McInnes 
973c8dd0c0dSLois Curfman McInnes    Level: advanced
974c8dd0c0dSLois Curfman McInnes 
9759e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9769e247f21SBarry Smith 
9773c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9783c632250SBarry Smith 
9793c632250SBarry Smith    On input w = x + y
9803c632250SBarry Smith 
98117bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
982b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
983ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9848f99978dSLois Curfman McInnes 
98517bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
986ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
987ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
988ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9899e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
990b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9918f99978dSLois Curfman McInnes 
992c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
993c8dd0c0dSLois Curfman McInnes 
99417bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
995c8dd0c0dSLois Curfman McInnes @*/
9963c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
997c8dd0c0dSLois Curfman McInnes {
9983c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
999c8dd0c0dSLois Curfman McInnes 
1000c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10013c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
1002c8dd0c0dSLois Curfman McInnes   if (f) {
1003c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1004c8dd0c0dSLois Curfman McInnes   }
1005c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1006c8dd0c0dSLois Curfman McInnes }
10079c3ca13aSBarry Smith #undef __FUNCT__
10089c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10099c3ca13aSBarry Smith /*@C
10109c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10117e4bb74cSBarry Smith          before the line search is called.
10129c3ca13aSBarry Smith 
10139c3ca13aSBarry Smith    Input Parameters:
10149c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10159c3ca13aSBarry Smith .  func - pointer to function
10169c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10179c3ca13aSBarry Smith 
10189c3ca13aSBarry Smith    Collective on SNES
10199c3ca13aSBarry Smith 
10209c3ca13aSBarry Smith    Calling sequence of func:
10219c3ca13aSBarry Smith .vb
10221e3347e8SBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
10239c3ca13aSBarry Smith .ve
10249c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10259c3ca13aSBarry Smith    on failure.
10269c3ca13aSBarry Smith 
10279c3ca13aSBarry Smith    Input parameters for func:
10289c3ca13aSBarry Smith +  snes - nonlinear context
10299c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10309c3ca13aSBarry Smith .  x - previous iterate
10319c3ca13aSBarry Smith -  y - new search direction and length
10329c3ca13aSBarry Smith 
10339c3ca13aSBarry Smith    Output parameters for func:
10349c3ca13aSBarry Smith +  y - search direction (possibly changed)
10359c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10369c3ca13aSBarry Smith 
10379c3ca13aSBarry Smith    Level: advanced
10389c3ca13aSBarry Smith 
10399c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10409c3ca13aSBarry Smith 
10419c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10429c3ca13aSBarry Smith 
10437e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10449c3ca13aSBarry Smith @*/
10459c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10469c3ca13aSBarry Smith {
10479c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10489c3ca13aSBarry Smith 
10499c3ca13aSBarry Smith   PetscFunctionBegin;
10509c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10519c3ca13aSBarry Smith   if (f) {
10529c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10539c3ca13aSBarry Smith   }
10549c3ca13aSBarry Smith   PetscFunctionReturn(0);
10559c3ca13aSBarry Smith }
10569c3ca13aSBarry Smith 
1057c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10589c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10599c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1060c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10614a2ae208SSatish Balay #undef __FUNCT__
10623c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10639c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1064c8dd0c0dSLois Curfman McInnes {
1065c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10663c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10673c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1068c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1069c8dd0c0dSLois Curfman McInnes }
1070c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10719c3ca13aSBarry Smith 
10729c3ca13aSBarry Smith EXTERN_C_BEGIN
10739c3ca13aSBarry Smith #undef __FUNCT__
10749c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10759c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10769c3ca13aSBarry Smith {
10779c3ca13aSBarry Smith   PetscFunctionBegin;
10789c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10799c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10809c3ca13aSBarry Smith   PetscFunctionReturn(0);
10819c3ca13aSBarry Smith }
10829c3ca13aSBarry Smith EXTERN_C_END
1083329e7e40SMatthew Knepley 
1084329e7e40SMatthew Knepley /*
10854b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
108604d965bbSLois Curfman McInnes 
108704d965bbSLois Curfman McInnes    Input Parameters:
108804d965bbSLois Curfman McInnes .  SNES - the SNES context
108904d965bbSLois Curfman McInnes .  viewer - visualization context
109004d965bbSLois Curfman McInnes 
109104d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
109204d965bbSLois Curfman McInnes */
10934a2ae208SSatish Balay #undef __FUNCT__
10944b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10956849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1096a935fc98SLois Curfman McInnes {
10974b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
10982fc52814SBarry Smith   const char     *cstr;
1099dfbe8321SBarry Smith   PetscErrorCode ierr;
110032077d6dSBarry Smith   PetscTruth     iascii;
1101a935fc98SLois Curfman McInnes 
11023a40ed3dSBarry Smith   PetscFunctionBegin;
110332077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
110432077d6dSBarry Smith   if (iascii) {
110517bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
110617bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
110717bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
110819bcc07fSBarry Smith     else                                                cstr = "unknown";
1109b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1110a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11115cd90555SBarry Smith   } else {
111279a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
111319bcc07fSBarry Smith   }
11143a40ed3dSBarry Smith   PetscFunctionReturn(0);
1115a935fc98SLois Curfman McInnes }
111604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
111704d965bbSLois Curfman McInnes /*
11184b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
111904d965bbSLois Curfman McInnes 
112004d965bbSLois Curfman McInnes    Input Parameter:
112104d965bbSLois Curfman McInnes .  snes - the SNES context
112204d965bbSLois Curfman McInnes 
112304d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
112404d965bbSLois Curfman McInnes */
11254a2ae208SSatish Balay #undef __FUNCT__
11264b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11276849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11285e42470aSBarry Smith {
11294b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1130e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1131dfbe8321SBarry Smith   PetscErrorCode ierr;
1132a7cc72afSBarry Smith   PetscInt       indx;
1133f1af5d2fSBarry Smith   PetscTruth     flg;
11345e42470aSBarry Smith 
11353a40ed3dSBarry Smith   PetscFunctionBegin;
1136b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11374b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11384b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11394b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1140186905e3SBarry Smith 
114117bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
114225cdf11fSBarry Smith     if (flg) {
114322e36657SBarry Smith       switch (indx) {
1144b49fd9e1SBarry Smith       case 0:
114517bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1146b49fd9e1SBarry Smith         break;
1147b49fd9e1SBarry Smith       case 1:
114817bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1149b49fd9e1SBarry Smith         break;
1150b49fd9e1SBarry Smith       case 2:
115117bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1152b49fd9e1SBarry Smith         break;
1153b49fd9e1SBarry Smith       case 3:
115417bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1155b49fd9e1SBarry Smith         break;
11565e42470aSBarry Smith       }
11575e42470aSBarry Smith     }
1158b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11593a40ed3dSBarry Smith   PetscFunctionReturn(0);
11605e42470aSBarry Smith }
116104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1162*3f149594SLisandro Dalcin #undef __FUNCT__
1163*3f149594SLisandro Dalcin #define __FUNCT__ "SNESConverged_LS"
1164*3f149594SLisandro Dalcin /*@C
1165*3f149594SLisandro Dalcin    SNESConverged_LS - Monitors the convergence of the line search
1166*3f149594SLisandro Dalcin    method SNESLS for solving systems of nonlinear equations (default).
1167*3f149594SLisandro Dalcin 
1168*3f149594SLisandro Dalcin    Collective on SNES
1169*3f149594SLisandro Dalcin 
1170*3f149594SLisandro Dalcin    Input Parameters:
1171*3f149594SLisandro Dalcin +  snes - the SNES context
1172*3f149594SLisandro Dalcin .  it - the iteration (0 indicates before any Newton steps)
1173*3f149594SLisandro Dalcin .  xnorm - 2-norm of current iterate
1174*3f149594SLisandro Dalcin .  pnorm - 2-norm of current step
1175*3f149594SLisandro Dalcin .  fnorm - 2-norm of function at current iterate
1176*3f149594SLisandro Dalcin -  dummy - unused context
1177*3f149594SLisandro Dalcin 
1178*3f149594SLisandro Dalcin    Output Parameter:
1179*3f149594SLisandro Dalcin .   reason  - one of
1180*3f149594SLisandro Dalcin $  SNES_CONVERGED_FNORM_ABS       - (fnorm < abstol),
1181*3f149594SLisandro Dalcin $  SNES_CONVERGED_PNORM_RELATIVE  - (pnorm < xtol*xnorm),
1182*3f149594SLisandro Dalcin $  SNES_CONVERGED_FNORM_RELATIVE  - (fnorm < rtol*fnorm0),
1183*3f149594SLisandro Dalcin $  SNES_DIVERGED_FUNCTION_COUNT   - (nfct > maxf),
1184*3f149594SLisandro Dalcin $  SNES_DIVERGED_FNORM_NAN        - (fnorm == NaN),
1185*3f149594SLisandro Dalcin $  SNES_CONVERGED_ITERATING       - (otherwise),
1186*3f149594SLisandro Dalcin 
1187*3f149594SLisandro Dalcin    where
1188*3f149594SLisandro Dalcin +    maxf - maximum number of function evaluations,
1189*3f149594SLisandro Dalcin             set with SNESSetTolerances()
1190*3f149594SLisandro Dalcin .    nfct - number of function evaluations,
1191*3f149594SLisandro Dalcin .    abstol - absolute function norm tolerance,
1192*3f149594SLisandro Dalcin             set with SNESSetTolerances()
1193*3f149594SLisandro Dalcin -    rtol - relative function norm tolerance, set with SNESSetTolerances()
1194*3f149594SLisandro Dalcin 
1195*3f149594SLisandro Dalcin    Level: intermediate
1196*3f149594SLisandro Dalcin 
1197*3f149594SLisandro Dalcin .keywords: SNES, nonlinear, default, converged, convergence
1198*3f149594SLisandro Dalcin 
1199*3f149594SLisandro Dalcin .seealso: SNESSetConvergenceTest()
1200*3f149594SLisandro Dalcin @*/
1201*3f149594SLisandro Dalcin PetscErrorCode PETSCSNES_DLLEXPORT SNESConverged_LS(SNES snes,PetscInt it,PetscReal xnorm,PetscReal pnorm,PetscReal fnorm,SNESConvergedReason *reason,void *dummy)
1202*3f149594SLisandro Dalcin {
1203*3f149594SLisandro Dalcin   PetscErrorCode ierr;
1204*3f149594SLisandro Dalcin   PetscFunctionBegin;
1205*3f149594SLisandro Dalcin   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1206*3f149594SLisandro Dalcin   PetscValidType(snes,1);
1207*3f149594SLisandro Dalcin   PetscValidPointer(reason,6);
1208*3f149594SLisandro Dalcin   ierr = SNESDefaultConverged(snes,it,xnorm,pnorm,fnorm,reason,dummy);CHKERRQ(ierr);
1209*3f149594SLisandro Dalcin   PetscFunctionReturn(0);
1210*3f149594SLisandro Dalcin }
1211*3f149594SLisandro Dalcin /* -------------------------------------------------------------------------- */
1212ebe3b25bSBarry Smith /*MC
1213ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
121404d965bbSLois Curfman McInnes 
1215ebe3b25bSBarry Smith    Options Database:
1216ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1217ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1218ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1219ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1220ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1221ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
122204d965bbSLois Curfman McInnes 
1223ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
122404d965bbSLois Curfman McInnes 
1225ee3001cbSBarry Smith    Level: beginner
1226ee3001cbSBarry Smith 
122717bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
122817bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
122917bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1230ebe3b25bSBarry Smith 
1231ebe3b25bSBarry Smith M*/
1232fb2e594dSBarry Smith EXTERN_C_BEGIN
12334a2ae208SSatish Balay #undef __FUNCT__
12344b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
123563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
12365e42470aSBarry Smith {
1237dfbe8321SBarry Smith   PetscErrorCode ierr;
12384b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12395e42470aSBarry Smith 
12403a40ed3dSBarry Smith   PetscFunctionBegin;
1241e7788613SBarry Smith   snes->ops->setup	     = SNESSetUp_LS;
1242e7788613SBarry Smith   snes->ops->solve	     = SNESSolve_LS;
1243e7788613SBarry Smith   snes->ops->destroy	     = SNESDestroy_LS;
1244e7788613SBarry Smith   snes->ops->converged	     = SNESConverged_LS;
1245e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
1246e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
12475baf8537SBarry Smith   snes->nwork                = 0;
12485e42470aSBarry Smith 
12494b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
125052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
12515e42470aSBarry Smith   snes->data    	= (void*)neP;
12525e42470aSBarry Smith   neP->alpha		= 1.e-4;
12535e42470aSBarry Smith   neP->maxstep		= 1.e8;
12545e42470aSBarry Smith   neP->steptol		= 1.e-12;
125517bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1256c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12573c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12583c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12593c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12603c632250SBarry Smith   neP->precheck         = PETSC_NULL;
126182bf6240SBarry Smith 
1262557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C",
1263557d3f75SLisandro Dalcin 					   "SNESLineSearchSet_LS",
1264557d3f75SLisandro Dalcin 					   SNESLineSearchSet_LS);CHKERRQ(ierr);
1265557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C",
1266557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPostCheck_LS",
1267557d3f75SLisandro Dalcin 					   SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
1268557d3f75SLisandro Dalcin   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C",
1269557d3f75SLisandro Dalcin 					   "SNESLineSearchSetPreCheck_LS",
1270557d3f75SLisandro Dalcin 					   SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
127182bf6240SBarry Smith 
12723a40ed3dSBarry Smith   PetscFunctionReturn(0);
12735e42470aSBarry Smith }
1274fb2e594dSBarry Smith EXTERN_C_END
127582bf6240SBarry Smith 
127682bf6240SBarry Smith 
127782bf6240SBarry Smith 
1278