xref: /petsc/src/snes/impls/ls/ls.c (revision 3d4c47100ce0d1dc158aed3628ce92315b1883ec)
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;
141*3d4c4710SBarry Smith   KSPConvergedReason kspreason;
1425e76c431SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
144*3d4c4710SBarry Smith   snes->numFailures            = 0;
145*3d4c4710SBarry 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);
1665e76c431SBarry Smith 
16770441072SBarry Smith   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16894a9d846SBarry Smith 
169d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
170d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
171d034289bSBarry Smith 
1725e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1735e76c431SBarry Smith 
174329e7e40SMatthew Knepley     /* Call general purpose update function */
175abc0a331SBarry Smith     if (snes->update) {
176329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
177de8cb200SMatthew Knepley     }
178329e7e40SMatthew Knepley 
179ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1805334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18194b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18223ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
183*3d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
184*3d4c4710SBarry Smith     if (kspreason < 0) {
185*3d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
186*3d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
187*3d4c4710SBarry Smith         PetscFunctionReturn(0);
188*3d4c4710SBarry Smith       }
189*3d4c4710SBarry Smith     }
190*3d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
19174637425SBarry Smith 
1929c3ca13aSBarry Smith     if (neP->precheckstep) {
1939c3ca13aSBarry Smith       PetscTruth changed_y = PETSC_FALSE;
1949c3ca13aSBarry Smith       ierr = (*neP->precheckstep)(snes,X,Y,neP->precheck,&changed_y);CHKERRQ(ierr);
1959c3ca13aSBarry Smith     }
1969c3ca13aSBarry Smith 
197b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
19874637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
19985471664SBarry Smith     }
20074637425SBarry Smith 
20190cbd584SBarry Smith     /* should check what happened to the linear solve? */
2023505fcc1SBarry Smith     snes->linear_its += lits;
203ae15b995SBarry Smith     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
204ea4d3ed3SLois Curfman McInnes 
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);
210a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
211ae15b995SBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);CHKERRQ(ierr);
212a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
2133c632250SBarry Smith     TMP = X; X = W; snes->vec_sol_always = X;  W = TMP;
214a3b891d8SBarry Smith     fnorm = gnorm;
2154a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
216a3b891d8SBarry Smith 
2175ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2185ed2d596SBarry Smith     snes->iter = i+1;
2195ed2d596SBarry Smith     snes->norm = fnorm;
2205ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2215ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2225ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2235ed2d596SBarry Smith 
224a7cc72afSBarry Smith     if (!lssucceed) {
2258a5d9ee4SBarry Smith       PetscTruth ismin;
22650ffb88aSMatthew Knepley 
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 */
23629e0b56fSBarry Smith     if (snes->converged) {
23729e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2383505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2393505fcc1SBarry Smith       if (snes->reason) {
24016c95cacSBarry Smith         break;
24116c95cacSBarry Smith       }
24216c95cacSBarry Smith     }
24329e0b56fSBarry Smith   }
24439e2f89bSBarry Smith   if (X != snes->vec_sol) {
245393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
246e15f7bb5SBarry Smith   }
247e15f7bb5SBarry Smith   if (F != snes->vec_func) {
248e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
249e15f7bb5SBarry Smith   }
25039e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
25139e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
25252392280SLois Curfman McInnes   if (i == maxits) {
253ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
2543505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
25552392280SLois Curfman McInnes   }
2563a40ed3dSBarry Smith   PetscFunctionReturn(0);
2575e76c431SBarry Smith }
25804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25904d965bbSLois Curfman McInnes /*
2604b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2614b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
26204d965bbSLois Curfman McInnes 
26304d965bbSLois Curfman McInnes    Input Parameter:
26404d965bbSLois Curfman McInnes .  snes - the SNES context
26504d965bbSLois Curfman McInnes .  x - the solution vector
26604d965bbSLois Curfman McInnes 
26704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26804d965bbSLois Curfman McInnes 
26904d965bbSLois Curfman McInnes    Notes:
2704b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
27104d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
27204d965bbSLois Curfman McInnes    the call to SNESSolve().
27304d965bbSLois Curfman McInnes  */
2744a2ae208SSatish Balay #undef __FUNCT__
2754b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
276dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2775e76c431SBarry Smith {
278dfbe8321SBarry Smith   PetscErrorCode ierr;
2793a40ed3dSBarry Smith 
2803a40ed3dSBarry Smith   PetscFunctionBegin;
281e74804ceSBarry Smith   if (!snes->work) {
28281b6cf68SLois Curfman McInnes     snes->nwork = 4;
283d7e8b826SBarry Smith     ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
284efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
28581b6cf68SLois Curfman McInnes     snes->vec_sol_update_always = snes->work[3];
286e74804ceSBarry Smith   }
2873a40ed3dSBarry Smith   PetscFunctionReturn(0);
2885e76c431SBarry Smith }
28904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
29004d965bbSLois Curfman McInnes /*
2914b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2924b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29304d965bbSLois Curfman McInnes 
29404d965bbSLois Curfman McInnes    Input Parameter:
29504d965bbSLois Curfman McInnes .  snes - the SNES context
29604d965bbSLois Curfman McInnes 
297de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29804d965bbSLois Curfman McInnes  */
2994a2ae208SSatish Balay #undef __FUNCT__
3004b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
301dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3025e76c431SBarry Smith {
303dfbe8321SBarry Smith   PetscErrorCode ierr;
3043a40ed3dSBarry Smith 
3053a40ed3dSBarry Smith   PetscFunctionBegin;
3065baf8537SBarry Smith   if (snes->nwork) {
3074b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
30821c89e3eSBarry Smith   }
3095c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3103a40ed3dSBarry Smith   PetscFunctionReturn(0);
3115e76c431SBarry Smith }
31204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3134a2ae208SSatish Balay #undef __FUNCT__
31417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNo"
31582bf6240SBarry Smith 
3164b828684SBarry Smith /*@C
31717bae607SBarry Smith    SNESLineSearchNo - This routine is not a line search at all;
3185e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3195e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3205e76c431SBarry Smith 
321c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
322c7afd0dbSLois Curfman McInnes 
3235e76c431SBarry Smith    Input Parameters:
324c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
325af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3265e76c431SBarry Smith .  x - current iterate
3275e76c431SBarry Smith .  f - residual evaluated at x
3283c632250SBarry Smith .  y - search direction
3295e76c431SBarry Smith .  w - work vector
330c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3315e76c431SBarry Smith 
332c4a48953SLois Curfman McInnes    Output Parameters:
333c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3343c632250SBarry Smith .  w - new iterate
3355e76c431SBarry Smith .  gnorm - 2-norm of g
3365e76c431SBarry Smith .  ynorm - 2-norm of search length
337a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
338fee21e36SBarry Smith 
339c4a48953SLois Curfman McInnes    Options Database Key:
34017bae607SBarry Smith .  -snes_ls basic - Activates SNESLineSearchNo()
341c4a48953SLois Curfman McInnes 
34236851e7fSLois Curfman McInnes    Level: advanced
34336851e7fSLois Curfman McInnes 
34428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
34528ae5a21SLois Curfman McInnes 
34617bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
34717bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms()
3485e76c431SBarry Smith @*/
34917bae607SBarry 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)
3505e76c431SBarry Smith {
351dfbe8321SBarry Smith   PetscErrorCode ierr;
3524b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3533c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
3545334005bSBarry Smith 
3553a40ed3dSBarry Smith   PetscFunctionBegin;
356a7cc72afSBarry Smith   *flag = PETSC_TRUE;
357d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
358a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);         /* ynorm = || y || */
35979f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3603c632250SBarry Smith   if (neP->postcheckstep) {
3613c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
3628f99978dSLois Curfman McInnes   }
3633c632250SBarry Smith   if (changed_y) {
36479f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
3653c632250SBarry Smith   }
3663c632250SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
367d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
36819717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
36919717074SBarry Smith   }
370d5e45103SBarry Smith   CHKERRQ(ierr);
371d5e45103SBarry Smith 
372a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
37343e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
374d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3753a40ed3dSBarry Smith   PetscFunctionReturn(0);
3765e76c431SBarry Smith }
37704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3784a2ae208SSatish Balay #undef __FUNCT__
37917bae607SBarry Smith #define __FUNCT__ "SNESLineSearchNoNorms"
38082bf6240SBarry Smith 
38129e0b56fSBarry Smith /*@C
38217bae607SBarry Smith    SNESLineSearchNoNorms - This routine is not a line search at
38329e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
38429e0b56fSBarry Smith    even compute the norm of the function or search direction; this
38529e0b56fSBarry Smith    is intended only when you know the full step is fine and are
38629e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
387c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
388c7afd0dbSLois Curfman McInnes 
389c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
39029e0b56fSBarry Smith 
39129e0b56fSBarry Smith    Input Parameters:
392c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
393af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
39429e0b56fSBarry Smith .  x - current iterate
39529e0b56fSBarry Smith .  f - residual evaluated at x
3963c632250SBarry Smith .  y - search direction
39729e0b56fSBarry Smith .  w - work vector
398c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
39929e0b56fSBarry Smith 
40029e0b56fSBarry Smith    Output Parameters:
401c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4023c632250SBarry Smith .  w - new iterate
40329e0b56fSBarry Smith .  gnorm - not changed
40429e0b56fSBarry Smith .  ynorm - not changed
405a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
406fee21e36SBarry Smith 
40729e0b56fSBarry Smith    Options Database Key:
40817bae607SBarry Smith .  -snes_ls basicnonorms - Activates SNESLineSearchNoNorms()
40929e0b56fSBarry Smith 
4108cbba510SBarry Smith    Notes:
41117bae607SBarry Smith    SNESLineSearchNoNorms() must be used in conjunction with
412ea56f5baSLois Curfman McInnes    either the options
413ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
414ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
415329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
416329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
417329f5518SBarry Smith 
418329f5518SBarry Smith    During the final iteration this will not evaluate the function at
419329f5518SBarry Smith    the solution point. This is to save a function evaluation while
420329f5518SBarry Smith    using pseudo-timestepping.
4218cbba510SBarry Smith 
422ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
423ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
424ea56f5baSLois Curfman McInnes    correct, since they are not computed.
425ea56f5baSLois Curfman McInnes 
426ea56f5baSLois Curfman McInnes    Level: advanced
4278cbba510SBarry Smith 
42829e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
42929e0b56fSBarry Smith 
43017bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(),
43117bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNo()
43229e0b56fSBarry Smith @*/
43317bae607SBarry 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)
43429e0b56fSBarry Smith {
435dfbe8321SBarry Smith   PetscErrorCode ierr;
4364b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
4373c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
43829e0b56fSBarry Smith 
4393a40ed3dSBarry Smith   PetscFunctionBegin;
440a7cc72afSBarry Smith   *flag = PETSC_TRUE;
441d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
44279f36460SBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y      */
4433c632250SBarry Smith   if (neP->postcheckstep) {
4443c632250SBarry Smith    ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
4453c632250SBarry Smith   }
4463c632250SBarry Smith   if (changed_y) {
44779f36460SBarry Smith     ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);            /* w <- x - y   */
4488f99978dSLois Curfman McInnes   }
449329f5518SBarry Smith 
450329f5518SBarry Smith   /* don't evaluate function the last time through */
451329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
4523c632250SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
453d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
45419717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
45519717074SBarry Smith     }
456d5e45103SBarry Smith     CHKERRQ(ierr);
457329f5518SBarry Smith   }
458d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4593a40ed3dSBarry Smith   PetscFunctionReturn(0);
46029e0b56fSBarry Smith }
46104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4624a2ae208SSatish Balay #undef __FUNCT__
46317bae607SBarry Smith #define __FUNCT__ "SNESLineSearchCubic"
4644b828684SBarry Smith /*@C
46517bae607SBarry Smith    SNESLineSearchCubic - Performs a cubic line search (default line search method).
4665e76c431SBarry Smith 
467c7afd0dbSLois Curfman McInnes    Collective on SNES
468c7afd0dbSLois Curfman McInnes 
4695e76c431SBarry Smith    Input Parameters:
470c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
471af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4725e76c431SBarry Smith .  x - current iterate
4735e76c431SBarry Smith .  f - residual evaluated at x
4743c632250SBarry Smith .  y - search direction
4755e76c431SBarry Smith .  w - work vector
476c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4775e76c431SBarry Smith 
478393d2d9aSLois Curfman McInnes    Output Parameters:
479c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4803c632250SBarry Smith .  w - new iterate
4815e76c431SBarry Smith .  gnorm - 2-norm of g
4825e76c431SBarry Smith .  ynorm - 2-norm of search length
483a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
484fee21e36SBarry Smith 
485c4a48953SLois Curfman McInnes    Options Database Key:
48617bae607SBarry Smith $  -snes_ls cubic - Activates SNESLineSearchCubic()
487c4a48953SLois Curfman McInnes 
4885e76c431SBarry Smith    Notes:
4895e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4905e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4915e76c431SBarry Smith 
49236851e7fSLois Curfman McInnes    Level: advanced
49336851e7fSLois Curfman McInnes 
49428ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49528ae5a21SLois Curfman McInnes 
49617bae607SBarry Smith .seealso: SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
4975e76c431SBarry Smith @*/
49817bae607SBarry 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)
4995e76c431SBarry Smith {
500406556e6SLois Curfman McInnes   /*
501406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
502406556e6SLois Curfman McInnes      minimization problem:
503406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
504406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
505406556e6SLois Curfman McInnes    */
506406556e6SLois Curfman McInnes 
5075ca10a37SBarry Smith   PetscReal      steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
508275d25c3SBarry Smith   PetscReal      maxstep,minlambda,alpha,lambda,lambdatemp;
509aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
510ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
5116b5873e3SBarry Smith #endif
512dfbe8321SBarry Smith   PetscErrorCode ierr;
513a7cc72afSBarry Smith   PetscInt       count;
5144b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
51579f36460SBarry Smith   PetscScalar    scale;
5163c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
5175e76c431SBarry Smith 
5183a40ed3dSBarry Smith   PetscFunctionBegin;
519d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
520a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
5215e76c431SBarry Smith   alpha   = neP->alpha;
5225e76c431SBarry Smith   maxstep = neP->maxstep;
5235e76c431SBarry Smith   steptol = neP->steptol;
5245e76c431SBarry Smith 
525cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
5269e247f21SBarry Smith   if (!*ynorm) {
527ae15b995SBarry Smith     ierr = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
528a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
5293c632250SBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
530a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
531a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
532ad922ac9SBarry Smith     goto theend1;
53394a9d846SBarry Smith   }
5345e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5355e42470aSBarry Smith     scale = maxstep/(*ynorm);
536ae15b995SBarry Smith     ierr = PetscInfo2(snes,"Scaling step by %G old ynorm %G\n",PetscRealPart(scale),*ynorm);CHKERRQ(ierr);
5372dcb1b2aSMatthew Knepley     ierr = VecScale(y,scale);CHKERRQ(ierr);
5385e76c431SBarry Smith     *ynorm = maxstep;
5395e76c431SBarry Smith   }
540ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5415ca10a37SBarry Smith   minlambda = steptol/rellength;
542a703fe33SLois Curfman McInnes   ierr      = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
543aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
544a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
545329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5465e42470aSBarry Smith #else
547a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&initslope);CHKERRQ(ierr);
5485e42470aSBarry Smith #endif
5495e76c431SBarry Smith   if (initslope > 0.0)  initslope = -initslope;
5505e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5515e76c431SBarry Smith 
552e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
55343e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
554ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
55543e71028SBarry Smith     *flag = PETSC_FALSE;
55643e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
55743e71028SBarry Smith     goto theend1;
55843e71028SBarry Smith   }
55919717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
560d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
56119717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
56219717074SBarry Smith   }
563d5e45103SBarry Smith   CHKERRQ(ierr);
564313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
56543e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5665d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
567ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
568ad922ac9SBarry Smith     goto theend1;
5695e76c431SBarry Smith   }
5705e76c431SBarry Smith 
5715e76c431SBarry Smith   /* Fit points with quadratic */
572313b5538SBarry Smith   lambda     = 1.0;
573a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5745e76c431SBarry Smith   lambdaprev = lambda;
5755e76c431SBarry Smith   gnormprev  = *gnorm;
576329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
577ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
578ddd12767SBarry Smith   else                         lambda = lambdatemp;
579275d25c3SBarry Smith 
580aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
581a23f76dfSSatish Balay   clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
5825e42470aSBarry Smith #else
583e68848bdSBarry Smith   ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
5845e42470aSBarry Smith #endif
58543e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
586ae15b995SBarry Smith     ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n",snes->nfuncs);CHKERRQ(ierr);
58743e71028SBarry Smith     *flag = PETSC_FALSE;
58843e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
58943e71028SBarry Smith     goto theend1;
59043e71028SBarry Smith   }
59119717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
592d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
59319717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
59419717074SBarry Smith   }
595d5e45103SBarry Smith   CHKERRQ(ierr);
596cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
59743e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5985ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
599ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
600ad922ac9SBarry Smith     goto theend1;
6015e76c431SBarry Smith   }
6025e76c431SBarry Smith 
6035e76c431SBarry Smith   /* Fit points with cubic */
6045e76c431SBarry Smith   count = 1;
6058229bfc2SMatthew Knepley   while (count < 10000) {
6065e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
607ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
608ae15b995SBarry 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);
60943e71028SBarry Smith       *flag = PETSC_FALSE;
61043e71028SBarry Smith       break;
6115e76c431SBarry Smith     }
6125d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6135d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
614ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6152b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6165e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6175e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6185e76c431SBarry Smith     if (a == 0.0) {
6195e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6205e76c431SBarry Smith     } else {
6215e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6225e76c431SBarry Smith     }
6235e76c431SBarry Smith     lambdaprev = lambda;
6245e76c431SBarry Smith     gnormprev  = *gnorm;
625329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
626329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6275e76c431SBarry Smith     else                         lambda     = lambdatemp;
628aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
629275d25c3SBarry Smith     clambda   = lambda;
630e68848bdSBarry Smith     ierr      = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
6315e42470aSBarry Smith #else
632e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
6335e42470aSBarry Smith #endif
63443e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
635ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
636ae15b995SBarry 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);
637ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
63843e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
639ed98166eSMatthew Knepley       break;
640ed98166eSMatthew Knepley     }
64119717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
642d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
64319717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
64419717074SBarry Smith     }
645d5e45103SBarry Smith     CHKERRQ(ierr);
646cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
64743e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6485ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
649ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubically determined step, lambda=%18.16e\n",lambda);CHKERRQ(ierr);
65043e71028SBarry Smith       break;
6512b022350SLois Curfman McInnes     } else {
652ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);CHKERRQ(ierr);
6535e76c431SBarry Smith     }
6545e76c431SBarry Smith     count++;
6555e76c431SBarry Smith   }
6568229bfc2SMatthew Knepley   if (count >= 10000) {
657cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6588229bfc2SMatthew Knepley   }
659ad922ac9SBarry Smith   theend1:
6608f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6613c632250SBarry Smith   if (neP->postcheckstep && *flag) {
6623c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
6633c632250SBarry Smith     if (changed_y) {
66479f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
6653c632250SBarry Smith     }
6663c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
6673c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
668d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
66919717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
67019717074SBarry Smith       }
671d5e45103SBarry Smith       CHKERRQ(ierr);
6728f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
67343e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6743c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
6758f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6763c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
6778f99978dSLois Curfman McInnes     }
6788f99978dSLois Curfman McInnes   }
679d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6803a40ed3dSBarry Smith   PetscFunctionReturn(0);
6815e76c431SBarry Smith }
68204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6834a2ae208SSatish Balay #undef __FUNCT__
68417bae607SBarry Smith #define __FUNCT__ "SNESLineSearchQuadratic"
6854b828684SBarry Smith /*@C
68617bae607SBarry Smith    SNESLineSearchQuadratic - Performs a quadratic line search.
6875e76c431SBarry Smith 
688c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
689c7afd0dbSLois Curfman McInnes 
6905e42470aSBarry Smith    Input Parameters:
691c7afd0dbSLois Curfman McInnes +  snes - the SNES context
692af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6935e42470aSBarry Smith .  x - current iterate
6945e42470aSBarry Smith .  f - residual evaluated at x
6953c632250SBarry Smith .  y - search direction
6965e42470aSBarry Smith .  w - work vector
697c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6985e42470aSBarry Smith 
699c4a48953SLois Curfman McInnes    Output Parameters:
7007f3332b4SBarry Smith +  g - residual evaluated at new iterate w
7017f3332b4SBarry Smith .  w - new iterate (x + alpha*y)
7025e42470aSBarry Smith .  gnorm - 2-norm of g
7035e42470aSBarry Smith .  ynorm - 2-norm of search length
704a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
705fee21e36SBarry Smith 
706c4a48953SLois Curfman McInnes    Options Database Key:
70717bae607SBarry Smith .  -snes_ls quadratic - Activates SNESLineSearchQuadratic()
708c4a48953SLois Curfman McInnes 
7095e42470aSBarry Smith    Notes:
71017bae607SBarry Smith    Use SNESLineSearchSet() to set this routine within the SNESLS method.
7115e42470aSBarry Smith 
71236851e7fSLois Curfman McInnes    Level: advanced
71336851e7fSLois Curfman McInnes 
714f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
715f59ffdedSLois Curfman McInnes 
71617bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchNo(), SNESLineSearchSet(), SNESLineSearchNoNorms()
7175e42470aSBarry Smith @*/
71817bae607SBarry 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)
7195e76c431SBarry Smith {
720406556e6SLois Curfman McInnes   /*
721406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
722406556e6SLois Curfman McInnes      minimization problem:
723406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
724406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
725406556e6SLois Curfman McInnes    */
726275d25c3SBarry Smith   PetscReal      steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,rellength;
727aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
728ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
7296b5873e3SBarry Smith #endif
730dfbe8321SBarry Smith   PetscErrorCode ierr;
731a7cc72afSBarry Smith   PetscInt       count;
7324b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
73379f36460SBarry Smith   PetscScalar    scale;
7343c632250SBarry Smith   PetscTruth     changed_w = PETSC_FALSE,changed_y = PETSC_FALSE;
7355e76c431SBarry Smith 
7363a40ed3dSBarry Smith   PetscFunctionBegin;
737d5ba7fb7SMatthew Knepley   ierr    = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
738a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7395e42470aSBarry Smith   alpha   = neP->alpha;
7405e42470aSBarry Smith   maxstep = neP->maxstep;
7415e42470aSBarry Smith   steptol = neP->steptol;
7425e76c431SBarry Smith 
7433505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
74463b9fa5eSBarry Smith   if (*ynorm == 0.0) {
745ae15b995SBarry Smith     ierr   = PetscInfo(snes,"Search direction and size is 0\n");CHKERRQ(ierr);
746b37302c6SLois Curfman McInnes     *gnorm = fnorm;
747e68848bdSBarry Smith     ierr   = VecCopy(x,w);CHKERRQ(ierr);
748b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
749a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
750ad922ac9SBarry Smith     goto theend2;
75194a9d846SBarry Smith   }
7525e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7535e42470aSBarry Smith     scale  = maxstep/(*ynorm);
7542dcb1b2aSMatthew Knepley     ierr   = VecScale(y,scale);CHKERRQ(ierr);
7555e42470aSBarry Smith     *ynorm = maxstep;
7565e76c431SBarry Smith   }
757ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7585ca10a37SBarry Smith   minlambda = steptol/rellength;
759a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
760aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
761a703fe33SLois Curfman McInnes   ierr      = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
762329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7635e42470aSBarry Smith #else
764a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7655e42470aSBarry Smith #endif
7665e42470aSBarry Smith   if (initslope > 0.0)  initslope = -initslope;
7675e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7685e42470aSBarry Smith 
769e68848bdSBarry Smith   ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
77043e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
771ae15b995SBarry Smith     ierr  = PetscInfo(snes,"Exceeded maximum function evaluations, while checking full step length!\n");CHKERRQ(ierr);
77243e71028SBarry Smith     *flag = PETSC_FALSE;
77343e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
77443e71028SBarry Smith     goto theend2;
77543e71028SBarry Smith   }
77619717074SBarry Smith   ierr = SNESComputeFunction(snes,w,g);
777d5e45103SBarry Smith   if (PetscExceptionValue(ierr)) {
77819717074SBarry Smith     PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
77919717074SBarry Smith   }
780d5e45103SBarry Smith   CHKERRQ(ierr);
781cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
78243e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7835d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
784ae15b995SBarry Smith     ierr = PetscInfo(snes,"Using full step\n");CHKERRQ(ierr);
785ad922ac9SBarry Smith     goto theend2;
7865e42470aSBarry Smith   }
7875e42470aSBarry Smith 
7885e42470aSBarry Smith   /* Fit points with quadratic */
789313b5538SBarry Smith   lambda = 1.0;
7905e42470aSBarry Smith   count = 1;
7915ca10a37SBarry Smith   while (PETSC_TRUE) {
7925e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
793ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Unable to find good step length! %D \n",count);CHKERRQ(ierr);
794ae15b995SBarry Smith       ierr = PetscInfo5(snes,"fnorm=%G, gnorm=%G, ynorm=%G, lambda=%G, initial slope=%G\n",fnorm,*gnorm,*ynorm,lambda,initslope);CHKERRQ(ierr);
795e68848bdSBarry Smith       ierr = VecCopy(x,w);CHKERRQ(ierr);
79643e71028SBarry Smith       *flag = PETSC_FALSE;
79743e71028SBarry Smith       break;
7985e42470aSBarry Smith     }
799a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
800329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
801329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
802329f5518SBarry Smith     else                         lambda     = lambdatemp;
803275d25c3SBarry Smith 
804aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
805e68848bdSBarry Smith     clambda   = lambda; ierr = VecWAXPY(w,-clambda,y,x);CHKERRQ(ierr);
8065e42470aSBarry Smith #else
807e68848bdSBarry Smith     ierr      = VecWAXPY(w,-lambda,y,x);CHKERRQ(ierr);
8085e42470aSBarry Smith #endif
80943e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
810ae15b995SBarry Smith       ierr  = PetscInfo1(snes,"Exceeded maximum function evaluations, while looking for good step length! %D \n",count);CHKERRQ(ierr);
811ae15b995SBarry 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);
812ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
81343e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
814ed98166eSMatthew Knepley       break;
815ed98166eSMatthew Knepley     }
81619717074SBarry Smith     ierr = SNESComputeFunction(snes,w,g);
817d5e45103SBarry Smith     if (PetscExceptionValue(ierr)) {
81819717074SBarry Smith       PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
81919717074SBarry Smith     }
820d5e45103SBarry Smith     CHKERRQ(ierr);
821cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
82243e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8235ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
824ae15b995SBarry Smith       ierr = PetscInfo1(snes,"Quadratically determined step, lambda=%G\n",lambda);CHKERRQ(ierr);
82506259719SBarry Smith       break;
8265e42470aSBarry Smith     }
8275e42470aSBarry Smith     count++;
8285e42470aSBarry Smith   }
829ad922ac9SBarry Smith   theend2:
8308f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8313c632250SBarry Smith   if (neP->postcheckstep) {
8323c632250SBarry Smith     ierr = (*neP->postcheckstep)(snes,x,y,w,neP->postcheck,&changed_y,&changed_w);CHKERRQ(ierr);
8333c632250SBarry Smith     if (changed_y) {
83479f36460SBarry Smith       ierr = VecWAXPY(w,-1.0,y,x);CHKERRQ(ierr);
8353c632250SBarry Smith     }
8363c632250SBarry Smith     if (changed_y || changed_w) { /* recompute the function if the step has changed */
8373c632250SBarry Smith       ierr = SNESComputeFunction(snes,w,g);
838d5e45103SBarry Smith       if (PetscExceptionValue(ierr)) {
83919717074SBarry Smith         PetscErrorCode pierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(pierr);
84019717074SBarry Smith       }
841d5e45103SBarry Smith       CHKERRQ(ierr);
8428f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
8433c632250SBarry Smith       ierr = VecNormBegin(w,NORM_2,ynorm);CHKERRQ(ierr);
8448f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8453c632250SBarry Smith       ierr = VecNormEnd(w,NORM_2,ynorm);CHKERRQ(ierr);
8463c632250SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
8478f99978dSLois Curfman McInnes     }
8488f99978dSLois Curfman McInnes   }
849d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8503a40ed3dSBarry Smith   PetscFunctionReturn(0);
8515e76c431SBarry Smith }
8522343ba6eSBarry Smith 
85304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8544a2ae208SSatish Balay #undef __FUNCT__
85517bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet"
856c9e6a524SLois Curfman McInnes /*@C
85717bae607SBarry Smith    SNESLineSearchSet - Sets the line search routine to be used
8584b27c08aSLois Curfman McInnes    by the method SNESLS.
8595e76c431SBarry Smith 
8605e76c431SBarry Smith    Input Parameters:
861c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
862af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
863c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8645e76c431SBarry Smith 
865fee21e36SBarry Smith    Collective on SNES
866fee21e36SBarry Smith 
867c4a48953SLois Curfman McInnes    Available Routines:
86817bae607SBarry Smith +  SNESLineSearchCubic() - default line search
86917bae607SBarry Smith .  SNESLineSearchQuadratic() - quadratic line search
87017bae607SBarry Smith .  SNESLineSearchNo() - the full Newton step (actually not a line search)
87117bae607SBarry Smith -  SNESLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8725e76c431SBarry Smith 
873c4a48953SLois Curfman McInnes     Options Database Keys:
8744b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8754b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8764b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8774b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8783304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8793304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
880c4a48953SLois Curfman McInnes 
8815e76c431SBarry Smith    Calling sequence of func:
882bd208895SLois Curfman McInnes .vb
883af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
884a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
885bd208895SLois Curfman McInnes .ve
8865e76c431SBarry Smith 
8875e76c431SBarry Smith     Input parameters for func:
888c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
889af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8905e76c431SBarry Smith .   x - current iterate
8915e76c431SBarry Smith .   f - residual evaluated at x
8923c632250SBarry Smith .   y - search direction
8935e76c431SBarry Smith .   w - work vector
894c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8955e76c431SBarry Smith 
8965e76c431SBarry Smith     Output parameters for func:
897c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8983c632250SBarry Smith .   w - new iterate
8995e76c431SBarry Smith .   gnorm - 2-norm of g
9005e76c431SBarry Smith .   ynorm - 2-norm of search length
901a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
902f59ffdedSLois Curfman McInnes 
90336851e7fSLois Curfman McInnes     Level: advanced
90436851e7fSLois Curfman McInnes 
905f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
906f59ffdedSLois Curfman McInnes 
90717bae607SBarry Smith .seealso: SNESLineSearchCubic(), SNESLineSearchQuadratic(), SNESLineSearchNo(), SNESLineSearchNoNorms(),
90817bae607SBarry Smith           SNESLineSearchSetPostCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams(), SNESLineSearchSetPreCheck()
9095e76c431SBarry Smith @*/
91017bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
9115e76c431SBarry Smith {
912a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
91382bf6240SBarry Smith 
9143a40ed3dSBarry Smith   PetscFunctionBegin;
91517bae607SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSet_C",(void (**)(void))&f);CHKERRQ(ierr);
91682bf6240SBarry Smith   if (f) {
917af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
91882bf6240SBarry Smith   }
9193a40ed3dSBarry Smith   PetscFunctionReturn(0);
9205e76c431SBarry Smith }
9218e019c35SBarry Smith 
922a7cc72afSBarry 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*/
92304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
924fb2e594dSBarry Smith EXTERN_C_BEGIN
9254a2ae208SSatish Balay #undef __FUNCT__
92617bae607SBarry Smith #define __FUNCT__ "SNESLineSearchSet_LS"
92717bae607SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSet_LS(SNES snes,FCN2 func,void *lsctx)
92882bf6240SBarry Smith {
92982bf6240SBarry Smith   PetscFunctionBegin;
9304b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
9314b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
93282bf6240SBarry Smith   PetscFunctionReturn(0);
93382bf6240SBarry Smith }
934fb2e594dSBarry Smith EXTERN_C_END
93504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9364a2ae208SSatish Balay #undef __FUNCT__
9373c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck"
938c8dd0c0dSLois Curfman McInnes /*@C
9393c632250SBarry Smith    SNESLineSearchSetPostCheck - Sets a routine to check the validity of new iterate computed
9404b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
941c8dd0c0dSLois Curfman McInnes 
942c8dd0c0dSLois Curfman McInnes    Input Parameters:
943c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
9443c632250SBarry Smith .  func - pointer to function
945c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
946c8dd0c0dSLois Curfman McInnes 
947c8dd0c0dSLois Curfman McInnes    Collective on SNES
948c8dd0c0dSLois Curfman McInnes 
949c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
950c8dd0c0dSLois Curfman McInnes .vb
9513c632250SBarry Smith    int func (SNES snes, Vec x,Vec y,Vec w,void *checkctx, PetscTruth *changed_y,PetscTruth *changed_w)
952c8dd0c0dSLois Curfman McInnes .ve
953b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
954b1ae27eaSLois Curfman McInnes    on failure.
955c8dd0c0dSLois Curfman McInnes 
956c8dd0c0dSLois Curfman McInnes    Input parameters for func:
957c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
958c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9593c632250SBarry Smith .  x - previous iterate
9603c632250SBarry Smith .  y - new search direction and length
9613c632250SBarry Smith -  w - current candidate iterate
962c8dd0c0dSLois Curfman McInnes 
963c8dd0c0dSLois Curfman McInnes    Output parameters for func:
9643c632250SBarry Smith +  y - search direction (possibly changed)
9653c632250SBarry Smith .  w - current iterate (possibly modified)
9663c632250SBarry Smith .  changed_y - indicates search direction was changed by this routine
9673c632250SBarry Smith -  changed_w - indicates current iterate was changed by this routine
968c8dd0c0dSLois Curfman McInnes 
969c8dd0c0dSLois Curfman McInnes    Level: advanced
970c8dd0c0dSLois Curfman McInnes 
9719e247f21SBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
9729e247f21SBarry Smith 
9733c632250SBarry Smith    Only one of changed_y and changed_w can  be PETSC_TRUE
9743c632250SBarry Smith 
9753c632250SBarry Smith    On input w = x + y
9763c632250SBarry Smith 
97717bae607SBarry Smith    SNESLineSearchNo() and SNESLineSearchNoNorms() (1) compute a candidate iterate u_{i+1}, (2) pass control
978b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
979ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9808f99978dSLois Curfman McInnes 
98117bae607SBarry Smith    SNESLineSearchQuadratic() and SNESLineSearchCubic() (1) compute a candidate iterate u_{i+1} as well as a
982ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
983ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
984ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
9859e247f21SBarry Smith    by flag=PETSC_TRUE).  The overhead of this extra function re-evaluation can be
986b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9878f99978dSLois Curfman McInnes 
988c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
989c8dd0c0dSLois Curfman McInnes 
99017bae607SBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPreCheck()
991c8dd0c0dSLois Curfman McInnes @*/
9923c632250SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void *checkctx)
993c8dd0c0dSLois Curfman McInnes {
9943c632250SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*),void*);
995c8dd0c0dSLois Curfman McInnes 
996c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9973c632250SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPostCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
998c8dd0c0dSLois Curfman McInnes   if (f) {
999c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
1000c8dd0c0dSLois Curfman McInnes   }
1001c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1002c8dd0c0dSLois Curfman McInnes }
10039c3ca13aSBarry Smith #undef __FUNCT__
10049c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck"
10059c3ca13aSBarry Smith /*@C
10069c3ca13aSBarry Smith    SNESLineSearchSetPreCheck - Sets a routine to check the validity of a new direction given by the linear solve
10077e4bb74cSBarry Smith          before the line search is called.
10089c3ca13aSBarry Smith 
10099c3ca13aSBarry Smith    Input Parameters:
10109c3ca13aSBarry Smith +  snes - nonlinear context obtained from SNESCreate()
10119c3ca13aSBarry Smith .  func - pointer to function
10129c3ca13aSBarry Smith -  checkctx - optional user-defined context for use by step checking routine
10139c3ca13aSBarry Smith 
10149c3ca13aSBarry Smith    Collective on SNES
10159c3ca13aSBarry Smith 
10169c3ca13aSBarry Smith    Calling sequence of func:
10179c3ca13aSBarry Smith .vb
10181e3347e8SBarry Smith    int func (SNES snes, Vec x,Vec y,void *checkctx, PetscTruth *changed_y)
10199c3ca13aSBarry Smith .ve
10209c3ca13aSBarry Smith    where func returns an error code of 0 on success and a nonzero
10219c3ca13aSBarry Smith    on failure.
10229c3ca13aSBarry Smith 
10239c3ca13aSBarry Smith    Input parameters for func:
10249c3ca13aSBarry Smith +  snes - nonlinear context
10259c3ca13aSBarry Smith .  checkctx - optional user-defined context for use by step checking routine
10269c3ca13aSBarry Smith .  x - previous iterate
10279c3ca13aSBarry Smith -  y - new search direction and length
10289c3ca13aSBarry Smith 
10299c3ca13aSBarry Smith    Output parameters for func:
10309c3ca13aSBarry Smith +  y - search direction (possibly changed)
10319c3ca13aSBarry Smith -  changed_y - indicates search direction was changed by this routine
10329c3ca13aSBarry Smith 
10339c3ca13aSBarry Smith    Level: advanced
10349c3ca13aSBarry Smith 
10359c3ca13aSBarry Smith    Notes: All line searches accept the new iterate computed by the line search checking routine.
10369c3ca13aSBarry Smith 
10379c3ca13aSBarry Smith .keywords: SNES, nonlinear, set, line search check, step check, routine
10389c3ca13aSBarry Smith 
10397e4bb74cSBarry Smith .seealso: SNESLineSearchSet(), SNESLineSearchSetPostCheck(), SNESSetUpdate()
10409c3ca13aSBarry Smith @*/
10419c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck(SNES snes,PetscErrorCode (*func)(SNES,Vec,Vec,void*,PetscTruth*),void *checkctx)
10429c3ca13aSBarry Smith {
10439c3ca13aSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,Vec,Vec,void*,PetscTruth*),void*);
10449c3ca13aSBarry Smith 
10459c3ca13aSBarry Smith   PetscFunctionBegin;
10469c3ca13aSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESLineSearchSetPreCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
10479c3ca13aSBarry Smith   if (f) {
10489c3ca13aSBarry Smith     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
10499c3ca13aSBarry Smith   }
10509c3ca13aSBarry Smith   PetscFunctionReturn(0);
10519c3ca13aSBarry Smith }
10529c3ca13aSBarry Smith 
1053c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
10549c3ca13aSBarry Smith typedef PetscErrorCode (*FCN1)(SNES,Vec,Vec,Vec,void*,PetscTruth*,PetscTruth*); /* force argument to next function to not be extern C*/
10559c3ca13aSBarry Smith typedef PetscErrorCode (*FCN3)(SNES,Vec,Vec,void*,PetscTruth*);                 /* force argument to next function to not be extern C*/
1056c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
10574a2ae208SSatish Balay #undef __FUNCT__
10583c632250SBarry Smith #define __FUNCT__ "SNESLineSearchSetPostCheck_LS"
10599c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPostCheck_LS(SNES snes,FCN1 func,void *checkctx)
1060c8dd0c0dSLois Curfman McInnes {
1061c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
10623c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheckstep = func;
10633c632250SBarry Smith   ((SNES_LS *)(snes->data))->postcheck     = checkctx;
1064c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
1065c8dd0c0dSLois Curfman McInnes }
1066c8dd0c0dSLois Curfman McInnes EXTERN_C_END
10679c3ca13aSBarry Smith 
10689c3ca13aSBarry Smith EXTERN_C_BEGIN
10699c3ca13aSBarry Smith #undef __FUNCT__
10709c3ca13aSBarry Smith #define __FUNCT__ "SNESLineSearchSetPreCheck_LS"
10719c3ca13aSBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESLineSearchSetPreCheck_LS(SNES snes,FCN3 func,void *checkctx)
10729c3ca13aSBarry Smith {
10739c3ca13aSBarry Smith   PetscFunctionBegin;
10749c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheckstep = func;
10759c3ca13aSBarry Smith   ((SNES_LS *)(snes->data))->precheck     = checkctx;
10769c3ca13aSBarry Smith   PetscFunctionReturn(0);
10779c3ca13aSBarry Smith }
10789c3ca13aSBarry Smith EXTERN_C_END
1079c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
108004d965bbSLois Curfman McInnes /*
10814b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
1082329e7e40SMatthew Knepley 
1083329e7e40SMatthew Knepley    Input Parameter:
1084329e7e40SMatthew Knepley .  snes - the SNES context
1085329e7e40SMatthew Knepley 
1086329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
1087329e7e40SMatthew Knepley */
1088329e7e40SMatthew Knepley #undef __FUNCT__
10894b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
10906849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
1091329e7e40SMatthew Knepley {
10924b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
1093329e7e40SMatthew Knepley 
1094329e7e40SMatthew Knepley   PetscFunctionBegin;
10954b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
10964b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
1097a83599f4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %G)\n",p,ls->alpha);
1098a83599f4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %G)\n",p,ls->maxstep);
1099a83599f4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %G)\n",p,ls->steptol);
1100329e7e40SMatthew Knepley   PetscFunctionReturn(0);
1101329e7e40SMatthew Knepley }
1102329e7e40SMatthew Knepley 
1103329e7e40SMatthew Knepley /*
11044b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
110504d965bbSLois Curfman McInnes 
110604d965bbSLois Curfman McInnes    Input Parameters:
110704d965bbSLois Curfman McInnes .  SNES - the SNES context
110804d965bbSLois Curfman McInnes .  viewer - visualization context
110904d965bbSLois Curfman McInnes 
111004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
111104d965bbSLois Curfman McInnes */
11124a2ae208SSatish Balay #undef __FUNCT__
11134b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
11146849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1115a935fc98SLois Curfman McInnes {
11164b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
11172fc52814SBarry Smith   const char     *cstr;
1118dfbe8321SBarry Smith   PetscErrorCode ierr;
111932077d6dSBarry Smith   PetscTruth     iascii;
1120a935fc98SLois Curfman McInnes 
11213a40ed3dSBarry Smith   PetscFunctionBegin;
112232077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
112332077d6dSBarry Smith   if (iascii) {
112417bae607SBarry Smith     if (ls->LineSearch == SNESLineSearchNo)             cstr = "SNESLineSearchNo";
112517bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchQuadratic) cstr = "SNESLineSearchQuadratic";
112617bae607SBarry Smith     else if (ls->LineSearch == SNESLineSearchCubic)     cstr = "SNESLineSearchCubic";
112719bcc07fSBarry Smith     else                                                cstr = "unknown";
1128b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1129a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%G, maxstep=%G, steptol=%G\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
11305cd90555SBarry Smith   } else {
113179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
113219bcc07fSBarry Smith   }
11333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1134a935fc98SLois Curfman McInnes }
113504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
113604d965bbSLois Curfman McInnes /*
11374b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
113804d965bbSLois Curfman McInnes 
113904d965bbSLois Curfman McInnes    Input Parameter:
114004d965bbSLois Curfman McInnes .  snes - the SNES context
114104d965bbSLois Curfman McInnes 
114204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
114304d965bbSLois Curfman McInnes */
11444a2ae208SSatish Balay #undef __FUNCT__
11454b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
11466849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
11475e42470aSBarry Smith {
11484b27c08aSLois Curfman McInnes   SNES_LS        *ls = (SNES_LS *)snes->data;
1149e5999256SBarry Smith   const char     *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1150dfbe8321SBarry Smith   PetscErrorCode ierr;
1151a7cc72afSBarry Smith   PetscInt       indx;
1152f1af5d2fSBarry Smith   PetscTruth     flg;
11535e42470aSBarry Smith 
11543a40ed3dSBarry Smith   PetscFunctionBegin;
1155b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
11564b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
11574b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
11584b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1159186905e3SBarry Smith 
116017bae607SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESLineSearchSet",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
116125cdf11fSBarry Smith     if (flg) {
116222e36657SBarry Smith       switch (indx) {
1163b49fd9e1SBarry Smith       case 0:
116417bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
1165b49fd9e1SBarry Smith         break;
1166b49fd9e1SBarry Smith       case 1:
116717bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1168b49fd9e1SBarry Smith         break;
1169b49fd9e1SBarry Smith       case 2:
117017bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchQuadratic,PETSC_NULL);CHKERRQ(ierr);
1171b49fd9e1SBarry Smith         break;
1172b49fd9e1SBarry Smith       case 3:
117317bae607SBarry Smith         ierr = SNESLineSearchSet(snes,SNESLineSearchCubic,PETSC_NULL);CHKERRQ(ierr);
1174b49fd9e1SBarry Smith         break;
11755e42470aSBarry Smith       }
11765e42470aSBarry Smith     }
1177b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
11783a40ed3dSBarry Smith   PetscFunctionReturn(0);
11795e42470aSBarry Smith }
118004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1181ebe3b25bSBarry Smith /*MC
1182ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
118304d965bbSLois Curfman McInnes 
1184ebe3b25bSBarry Smith    Options Database:
1185ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1186ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1187ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1188ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1189ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1190ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
119104d965bbSLois Curfman McInnes 
1192ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
119304d965bbSLois Curfman McInnes 
1194ee3001cbSBarry Smith    Level: beginner
1195ee3001cbSBarry Smith 
119617bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
119717bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
119817bae607SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck()
1199ebe3b25bSBarry Smith 
1200ebe3b25bSBarry Smith M*/
1201fb2e594dSBarry Smith EXTERN_C_BEGIN
12024a2ae208SSatish Balay #undef __FUNCT__
12034b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
120463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
12055e42470aSBarry Smith {
1206dfbe8321SBarry Smith   PetscErrorCode ierr;
12074b27c08aSLois Curfman McInnes   SNES_LS        *neP;
12085e42470aSBarry Smith 
12093a40ed3dSBarry Smith   PetscFunctionBegin;
12104b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
12114b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
12124b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
12134b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
12144b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
12154b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
12164b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
12175baf8537SBarry Smith   snes->nwork           = 0;
12185e42470aSBarry Smith 
12194b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
122052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
12215e42470aSBarry Smith   snes->data    	= (void*)neP;
12225e42470aSBarry Smith   neP->alpha		= 1.e-4;
12235e42470aSBarry Smith   neP->maxstep		= 1.e8;
12245e42470aSBarry Smith   neP->steptol		= 1.e-12;
122517bae607SBarry Smith   neP->LineSearch       = SNESLineSearchCubic;
1226c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
12273c632250SBarry Smith   neP->postcheckstep    = PETSC_NULL;
12283c632250SBarry Smith   neP->postcheck        = PETSC_NULL;
12293c632250SBarry Smith   neP->precheckstep     = PETSC_NULL;
12303c632250SBarry Smith   neP->precheck         = PETSC_NULL;
123182bf6240SBarry Smith 
123217bae607SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSet_C","SNESLineSearchSet_LS",SNESLineSearchSet_LS);CHKERRQ(ierr);
12333c632250SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPostCheck_C","SNESLineSearchSetPostCheck_LS",SNESLineSearchSetPostCheck_LS);CHKERRQ(ierr);
12349c3ca13aSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetPreCheck_C","SNESLineSearchSetPreCheck_LS",SNESLineSearchSetPreCheck_LS);CHKERRQ(ierr);
123582bf6240SBarry Smith 
12363a40ed3dSBarry Smith   PetscFunctionReturn(0);
12375e42470aSBarry Smith }
1238fb2e594dSBarry Smith EXTERN_C_END
123982bf6240SBarry Smith 
124082bf6240SBarry Smith 
124182bf6240SBarry Smith 
1242