xref: /petsc/src/snes/impls/ls/ls.c (revision 43e71028e8898ddee9dd80069fd00b0172b84ee3)
15e76c431SBarry Smith 
270f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
35e42470aSBarry Smith 
48a5d9ee4SBarry Smith /*
58a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
68a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
736669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
836669109SBarry Smith     for this trick.
98a5d9ee4SBarry Smith */
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal      a1;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
1636669109SBarry Smith   PetscTruth     hastranspose;
178a5d9ee4SBarry Smith 
188a5d9ee4SBarry Smith   PetscFunctionBegin;
198a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2036669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2136669109SBarry Smith   if (hastranspose) {
228a5d9ee4SBarry Smith     /* Compute || J^T F|| */
238a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
248a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
254b27c08aSLois Curfman McInnes     PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm);
268a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2736669109SBarry Smith   } else {
2836669109SBarry Smith     Vec         work;
29ea709b57SSatish Balay     PetscScalar result;
3036669109SBarry Smith     PetscReal   wnorm;
3136669109SBarry Smith 
3236669109SBarry Smith     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
3336669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
3736669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
3836669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
394b27c08aSLois Curfman McInnes     PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1);
4036669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4136669109SBarry Smith   }
428a5d9ee4SBarry Smith   PetscFunctionReturn(0);
438a5d9ee4SBarry Smith }
448a5d9ee4SBarry Smith 
4574637425SBarry Smith /*
465ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4774637425SBarry Smith */
484a2ae208SSatish Balay #undef __FUNCT__
494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
50dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5174637425SBarry Smith {
5274637425SBarry Smith   PetscReal     a1,a2;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
5474637425SBarry Smith   PetscTruth    hastranspose;
55ea709b57SSatish Balay   PetscScalar   mone = -1.0;
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);
6174637425SBarry Smith     ierr = VecAXPY(&mone,F,W1);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);
679994e62eSSatish Balay     if (a1 != 0) {
684b27c08aSLois Curfman McInnes       PetscLogInfo(0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1);
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;
141c293cc10SBarry Smith   KSP            ksp;
1425e76c431SBarry Smith 
1433a40ed3dSBarry Smith   PetscFunctionBegin;
14494b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
145184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
146184914b5SBarry Smith 
1475e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1485e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1505e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1515e42470aSBarry Smith   G		= snes->work[1];
1525e42470aSBarry Smith   W		= snes->work[2];
1535e76c431SBarry Smith 
1544c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1554c49b128SBarry Smith   snes->iter = 0;
1564c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1575334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
158cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
159*43e71028SBarry Smith   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1600f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1615e42470aSBarry Smith   snes->norm = fnorm;
1620f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16384c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16494a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1655e76c431SBarry Smith 
16670441072SBarry Smith   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16794a9d846SBarry Smith 
168d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
169d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
170d034289bSBarry Smith 
1715e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1725e76c431SBarry Smith 
173329e7e40SMatthew Knepley     /* Call general purpose update function */
174abc0a331SBarry Smith     if (snes->update) {
175329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
176de8cb200SMatthew Knepley     }
177329e7e40SMatthew Knepley 
178ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1795334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18094b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18123ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
182c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
18374637425SBarry Smith 
184b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
18574637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
18685471664SBarry Smith     }
18774637425SBarry Smith 
18890cbd584SBarry Smith     /* should check what happened to the linear solve? */
1893505fcc1SBarry Smith     snes->linear_its += lits;
19077431f27SBarry Smith     PetscLogInfo(snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits);
191ea4d3ed3SLois Curfman McInnes 
192ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
193ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
194ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
195ea4d3ed3SLois Curfman McInnes     */
19681b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
197a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
198a7cc72afSBarry Smith     PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed);
199*43e71028SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
200a3b891d8SBarry Smith 
201a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
202a3b891d8SBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
203a3b891d8SBarry Smith     fnorm = gnorm;
204a3b891d8SBarry Smith 
2055ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2065ed2d596SBarry Smith     snes->iter = i+1;
2075ed2d596SBarry Smith     snes->norm = fnorm;
2085ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2095ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2105ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2115ed2d596SBarry Smith 
212a7cc72afSBarry Smith     if (!lssucceed) {
2138a5d9ee4SBarry Smith       PetscTruth ismin;
21450ffb88aSMatthew Knepley 
21550ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2163505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2178a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2188a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2193505fcc1SBarry Smith         break;
2203505fcc1SBarry Smith       }
22150ffb88aSMatthew Knepley     }
2225e76c431SBarry Smith 
2235e76c431SBarry Smith     /* Test for convergence */
22429e0b56fSBarry Smith     if (snes->converged) {
22529e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2263505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2273505fcc1SBarry Smith       if (snes->reason) {
22816c95cacSBarry Smith         break;
22916c95cacSBarry Smith       }
23016c95cacSBarry Smith     }
23129e0b56fSBarry Smith   }
23239e2f89bSBarry Smith   if (X != snes->vec_sol) {
233393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
234e15f7bb5SBarry Smith   }
235e15f7bb5SBarry Smith   if (F != snes->vec_func) {
236e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
237e15f7bb5SBarry Smith   }
23839e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
23939e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
24052392280SLois Curfman McInnes   if (i == maxits) {
24177431f27SBarry Smith     PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits);
2423505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
24352392280SLois Curfman McInnes   }
2443a40ed3dSBarry Smith   PetscFunctionReturn(0);
2455e76c431SBarry Smith }
24604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
24704d965bbSLois Curfman McInnes /*
2484b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2494b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25004d965bbSLois Curfman McInnes 
25104d965bbSLois Curfman McInnes    Input Parameter:
25204d965bbSLois Curfman McInnes .  snes - the SNES context
25304d965bbSLois Curfman McInnes .  x - the solution vector
25404d965bbSLois Curfman McInnes 
25504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25604d965bbSLois Curfman McInnes 
25704d965bbSLois Curfman McInnes    Notes:
2584b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
25904d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26004d965bbSLois Curfman McInnes    the call to SNESSolve().
26104d965bbSLois Curfman McInnes  */
2624a2ae208SSatish Balay #undef __FUNCT__
2634b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
264dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2655e76c431SBarry Smith {
266dfbe8321SBarry Smith   PetscErrorCode ierr;
2673a40ed3dSBarry Smith 
2683a40ed3dSBarry Smith   PetscFunctionBegin;
26981b6cf68SLois Curfman McInnes   snes->nwork = 4;
270d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
271efee365bSSatish Balay   ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
27281b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2733a40ed3dSBarry Smith   PetscFunctionReturn(0);
2745e76c431SBarry Smith }
27504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27604d965bbSLois Curfman McInnes /*
2774b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2784b27c08aSLois Curfman McInnes    with SNESCreate_LS().
27904d965bbSLois Curfman McInnes 
28004d965bbSLois Curfman McInnes    Input Parameter:
28104d965bbSLois Curfman McInnes .  snes - the SNES context
28204d965bbSLois Curfman McInnes 
283de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28404d965bbSLois Curfman McInnes  */
2854a2ae208SSatish Balay #undef __FUNCT__
2864b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
287dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2885e76c431SBarry Smith {
289dfbe8321SBarry Smith   PetscErrorCode ierr;
2903a40ed3dSBarry Smith 
2913a40ed3dSBarry Smith   PetscFunctionBegin;
2925baf8537SBarry Smith   if (snes->nwork) {
2934b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29421c89e3eSBarry Smith   }
2955c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2963a40ed3dSBarry Smith   PetscFunctionReturn(0);
2975e76c431SBarry Smith }
29804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2994a2ae208SSatish Balay #undef __FUNCT__
3004a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
30182bf6240SBarry Smith 
3024b828684SBarry Smith /*@C
3035e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3045e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3055e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3065e76c431SBarry Smith 
307c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
308c7afd0dbSLois Curfman McInnes 
3095e76c431SBarry Smith    Input Parameters:
310c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
311af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3125e76c431SBarry Smith .  x - current iterate
3135e76c431SBarry Smith .  f - residual evaluated at x
3145e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3155e76c431SBarry Smith .  w - work vector
316c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3175e76c431SBarry Smith 
318c4a48953SLois Curfman McInnes    Output Parameters:
319c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3205e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3215e76c431SBarry Smith .  gnorm - 2-norm of g
3225e76c431SBarry Smith .  ynorm - 2-norm of search length
323a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
324fee21e36SBarry Smith 
325c4a48953SLois Curfman McInnes    Options Database Key:
3264b27c08aSLois Curfman McInnes .  -snes_ls basic - Activates SNESNoLineSearch()
327c4a48953SLois Curfman McInnes 
32836851e7fSLois Curfman McInnes    Level: advanced
32936851e7fSLois Curfman McInnes 
33028ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33128ae5a21SLois Curfman McInnes 
332f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
333af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3345e76c431SBarry Smith @*/
335a7cc72afSBarry Smith PetscErrorCode SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
3365e76c431SBarry Smith {
337dfbe8321SBarry Smith   PetscErrorCode ierr;
338ea709b57SSatish Balay   PetscScalar    mone = -1.0;
3394b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3408f99978dSLois Curfman McInnes   PetscTruth     change_y = PETSC_FALSE;
3415334005bSBarry Smith 
3423a40ed3dSBarry Smith   PetscFunctionBegin;
343a7cc72afSBarry Smith   *flag = PETSC_TRUE;
344d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
345a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
346ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3478f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3488f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3498f99978dSLois Curfman McInnes   }
350ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
351a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
352*43e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
353d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3543a40ed3dSBarry Smith   PetscFunctionReturn(0);
3555e76c431SBarry Smith }
35604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3574a2ae208SSatish Balay #undef __FUNCT__
3584a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
35982bf6240SBarry Smith 
36029e0b56fSBarry Smith /*@C
36129e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
36229e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
36329e0b56fSBarry Smith    even compute the norm of the function or search direction; this
36429e0b56fSBarry Smith    is intended only when you know the full step is fine and are
36529e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
366c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
367c7afd0dbSLois Curfman McInnes 
368c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
36929e0b56fSBarry Smith 
37029e0b56fSBarry Smith    Input Parameters:
371c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
372af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
37329e0b56fSBarry Smith .  x - current iterate
37429e0b56fSBarry Smith .  f - residual evaluated at x
37529e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
37629e0b56fSBarry Smith .  w - work vector
377c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
37829e0b56fSBarry Smith 
37929e0b56fSBarry Smith    Output Parameters:
380c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
38129e0b56fSBarry Smith .  gnorm - not changed
38229e0b56fSBarry Smith .  ynorm - not changed
383a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
384fee21e36SBarry Smith 
38529e0b56fSBarry Smith    Options Database Key:
3864b27c08aSLois Curfman McInnes .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
38729e0b56fSBarry Smith 
3888cbba510SBarry Smith    Notes:
389ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
390ea56f5baSLois Curfman McInnes    either the options
391ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
392ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
393329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
394329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
395329f5518SBarry Smith 
396329f5518SBarry Smith    During the final iteration this will not evaluate the function at
397329f5518SBarry Smith    the solution point. This is to save a function evaluation while
398329f5518SBarry Smith    using pseudo-timestepping.
3998cbba510SBarry Smith 
400ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
401ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
402ea56f5baSLois Curfman McInnes    correct, since they are not computed.
403ea56f5baSLois Curfman McInnes 
404ea56f5baSLois Curfman McInnes    Level: advanced
4058cbba510SBarry Smith 
40629e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
40729e0b56fSBarry Smith 
40829e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
40929e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
41029e0b56fSBarry Smith @*/
411a7cc72afSBarry Smith PetscErrorCode SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
41229e0b56fSBarry Smith {
413dfbe8321SBarry Smith   PetscErrorCode ierr;
414ea709b57SSatish Balay   PetscScalar mone = -1.0;
4154b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
4168f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
41729e0b56fSBarry Smith 
4183a40ed3dSBarry Smith   PetscFunctionBegin;
419a7cc72afSBarry Smith   *flag = PETSC_TRUE;
420d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
42129e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4228f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4238f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4248f99978dSLois Curfman McInnes   }
425329f5518SBarry Smith 
426329f5518SBarry Smith   /* don't evaluate function the last time through */
427329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
42829e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
429329f5518SBarry Smith   }
430d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4313a40ed3dSBarry Smith   PetscFunctionReturn(0);
43229e0b56fSBarry Smith }
43304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4344a2ae208SSatish Balay #undef __FUNCT__
4354a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4364b828684SBarry Smith /*@C
437f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4385e76c431SBarry Smith 
439c7afd0dbSLois Curfman McInnes    Collective on SNES
440c7afd0dbSLois Curfman McInnes 
4415e76c431SBarry Smith    Input Parameters:
442c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
443af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4445e76c431SBarry Smith .  x - current iterate
4455e76c431SBarry Smith .  f - residual evaluated at x
4465e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4475e76c431SBarry Smith .  w - work vector
448c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4495e76c431SBarry Smith 
450393d2d9aSLois Curfman McInnes    Output Parameters:
451c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4525e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4535e76c431SBarry Smith .  gnorm - 2-norm of g
4545e76c431SBarry Smith .  ynorm - 2-norm of search length
455a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
456fee21e36SBarry Smith 
457c4a48953SLois Curfman McInnes    Options Database Key:
4584b27c08aSLois Curfman McInnes $  -snes_ls cubic - Activates SNESCubicLineSearch()
459c4a48953SLois Curfman McInnes 
4605e76c431SBarry Smith    Notes:
4615e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4625e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4635e76c431SBarry Smith 
46436851e7fSLois Curfman McInnes    Level: advanced
46536851e7fSLois Curfman McInnes 
46628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
46728ae5a21SLois Curfman McInnes 
468af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4695e76c431SBarry Smith @*/
470a7cc72afSBarry Smith PetscErrorCode SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
4715e76c431SBarry Smith {
472406556e6SLois Curfman McInnes   /*
473406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
474406556e6SLois Curfman McInnes      minimization problem:
475406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
476406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
477406556e6SLois Curfman McInnes    */
478406556e6SLois Curfman McInnes 
4795ca10a37SBarry Smith   PetscReal   steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
480329f5518SBarry Smith   PetscReal   maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
481aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
482ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
4836b5873e3SBarry Smith #endif
484dfbe8321SBarry Smith   PetscErrorCode ierr;
485a7cc72afSBarry Smith   PetscInt count;
4864b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
487ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
4888f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
4895e76c431SBarry Smith 
4903a40ed3dSBarry Smith   PetscFunctionBegin;
491d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
492a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
4935e76c431SBarry Smith   alpha   = neP->alpha;
4945e76c431SBarry Smith   maxstep = neP->maxstep;
4955e76c431SBarry Smith   steptol = neP->steptol;
4965e76c431SBarry Smith 
497cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
49863b9fa5eSBarry Smith   if (*ynorm == 0.0) {
49963b9fa5eSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n");
500a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
501a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
502a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
503a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
504ad922ac9SBarry Smith     goto theend1;
50594a9d846SBarry Smith   }
5065e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5075e42470aSBarry Smith     scale = maxstep/(*ynorm);
508aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
509b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
5106b5873e3SBarry Smith #else
511b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
5126b5873e3SBarry Smith #endif
513393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5145e76c431SBarry Smith     *ynorm = maxstep;
5155e76c431SBarry Smith   }
516ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5175ca10a37SBarry Smith   minlambda = steptol/rellength;
518a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
519aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
520a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
521329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5225e42470aSBarry Smith #else
523a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5245e42470aSBarry Smith #endif
5255e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5265e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5275e76c431SBarry Smith 
528393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5295334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
530*43e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
531*43e71028SBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n");
532*43e71028SBarry Smith     ierr = VecCopy(w,y);CHKERRQ(ierr);
533*43e71028SBarry Smith     *flag = PETSC_FALSE;
534*43e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
535*43e71028SBarry Smith     goto theend1;
536*43e71028SBarry Smith   }
53778b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
538313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
539*43e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5405d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
541393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
542b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
543ad922ac9SBarry Smith     goto theend1;
5445e76c431SBarry Smith   }
5455e76c431SBarry Smith 
5465e76c431SBarry Smith   /* Fit points with quadratic */
547313b5538SBarry Smith   lambda = 1.0;
548a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5495e76c431SBarry Smith   lambdaprev = lambda;
550*43e71028SBarry Smith   printf("tmp %g ddsd %g %g %g %g\n",lambdatemp,((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope),*gnorm,fnorm,initslope);
5515e76c431SBarry Smith   gnormprev = *gnorm;
552329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
553ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
554ddd12767SBarry Smith   else                         lambda = lambdatemp;
555393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
556ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
557aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
558ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5595e42470aSBarry Smith #else
560ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5615e42470aSBarry Smith #endif
562*43e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
563*43e71028SBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n");
564*43e71028SBarry Smith     ierr = VecCopy(w,y);CHKERRQ(ierr);
565*43e71028SBarry Smith     *flag = PETSC_FALSE;
566*43e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
567*43e71028SBarry Smith     goto theend1;
568*43e71028SBarry Smith   }
56978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
570cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
571*43e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5725ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
573393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
5745ed2d596SBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda);
575ad922ac9SBarry Smith     goto theend1;
5765e76c431SBarry Smith   }
5775e76c431SBarry Smith 
5785e76c431SBarry Smith   /* Fit points with cubic */
5795e76c431SBarry Smith   count = 1;
5808229bfc2SMatthew Knepley   while (count < 10000) {
5815e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
58277431f27SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count);
5835ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
584f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
585*43e71028SBarry Smith       *flag = PETSC_FALSE;
586*43e71028SBarry Smith       break;
5875e76c431SBarry Smith     }
588*43e71028SBarry Smith     printf("lamd %g %g\n",lambda,lambdaprev);
5895d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
5905d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
591ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5922b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5935e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5945e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5955e76c431SBarry Smith     if (a == 0.0) {
5965e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5975e76c431SBarry Smith     } else {
5985e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5995e76c431SBarry Smith     }
6005e76c431SBarry Smith     lambdaprev = lambda;
6015e76c431SBarry Smith     gnormprev  = *gnorm;
602329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
603329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6045e76c431SBarry Smith     else                         lambda     = lambdatemp;
605393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
606ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
607aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
608ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
609393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6105e42470aSBarry Smith #else
611ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6125e42470aSBarry Smith #endif
613*43e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
614*43e71028SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
615ed98166eSMatthew Knepley       PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
616*43e71028SBarry Smith       ierr = VecCopy(w,y);CHKERRQ(ierr);
617ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
618*43e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
619ed98166eSMatthew Knepley       break;
620ed98166eSMatthew Knepley     }
621*43e71028SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
622cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
623*43e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6245ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
625393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
6265ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda);
627*43e71028SBarry Smith       break;
6282b022350SLois Curfman McInnes     } else {
6295ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);
6305e76c431SBarry Smith     }
6315e76c431SBarry Smith     count++;
6325e76c431SBarry Smith   }
6338229bfc2SMatthew Knepley   if (count >= 10000) {
634cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6358229bfc2SMatthew Knepley   }
636ad922ac9SBarry Smith   theend1:
6378f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
638*43e71028SBarry Smith   if (neP->CheckStep && *flag) {
6398f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
640abc0a331SBarry Smith     if (change_y) { /* recompute the function if the step has changed */
6418f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6428f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
643*43e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
644*43e71028SBarry Smith       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6458f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6468f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6478f99978dSLois Curfman McInnes     }
6488f99978dSLois Curfman McInnes   }
649d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6503a40ed3dSBarry Smith   PetscFunctionReturn(0);
6515e76c431SBarry Smith }
65204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6534a2ae208SSatish Balay #undef __FUNCT__
6544a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6554b828684SBarry Smith /*@C
656f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6575e76c431SBarry Smith 
658c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
659c7afd0dbSLois Curfman McInnes 
6605e42470aSBarry Smith    Input Parameters:
661c7afd0dbSLois Curfman McInnes +  snes - the SNES context
662af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6635e42470aSBarry Smith .  x - current iterate
6645e42470aSBarry Smith .  f - residual evaluated at x
6655e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6665e42470aSBarry Smith .  w - work vector
667c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6685e42470aSBarry Smith 
669c4a48953SLois Curfman McInnes    Output Parameters:
670c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6715e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6725e42470aSBarry Smith .  gnorm - 2-norm of g
6735e42470aSBarry Smith .  ynorm - 2-norm of search length
674a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
675fee21e36SBarry Smith 
676c4a48953SLois Curfman McInnes    Options Database Key:
6774b27c08aSLois Curfman McInnes .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
678c4a48953SLois Curfman McInnes 
6795e42470aSBarry Smith    Notes:
6804b27c08aSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNESLS method.
6815e42470aSBarry Smith 
68236851e7fSLois Curfman McInnes    Level: advanced
68336851e7fSLois Curfman McInnes 
684f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
685f59ffdedSLois Curfman McInnes 
686af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6875e42470aSBarry Smith @*/
688a7cc72afSBarry Smith PetscErrorCode SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
6895e76c431SBarry Smith {
690406556e6SLois Curfman McInnes   /*
691406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
692406556e6SLois Curfman McInnes      minimization problem:
693406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
694406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
695406556e6SLois Curfman McInnes    */
6965ca10a37SBarry Smith   PetscReal   steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
697aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
698ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
6996b5873e3SBarry Smith #endif
700dfbe8321SBarry Smith   PetscErrorCode ierr;
701a7cc72afSBarry Smith   PetscInt count;
7024b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
703ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
7048f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
7055e76c431SBarry Smith 
7063a40ed3dSBarry Smith   PetscFunctionBegin;
707d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
708a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7095e42470aSBarry Smith   alpha   = neP->alpha;
7105e42470aSBarry Smith   maxstep = neP->maxstep;
7115e42470aSBarry Smith   steptol = neP->steptol;
7125e76c431SBarry Smith 
7133505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
71463b9fa5eSBarry Smith   if (*ynorm == 0.0) {
715b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
716b37302c6SLois Curfman McInnes     *gnorm = fnorm;
717b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
718b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
719a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
720ad922ac9SBarry Smith     goto theend2;
72194a9d846SBarry Smith   }
7225e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7235e42470aSBarry Smith     scale = maxstep/(*ynorm);
724393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
7255e42470aSBarry Smith     *ynorm = maxstep;
7265e76c431SBarry Smith   }
727ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7285ca10a37SBarry Smith   minlambda = steptol/rellength;
729a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
730aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
731a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
732329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7335e42470aSBarry Smith #else
734a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7355e42470aSBarry Smith #endif
7365e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
7375e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7385e42470aSBarry Smith 
739393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7405334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
741*43e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
742*43e71028SBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n");
743*43e71028SBarry Smith     ierr = VecCopy(w,y);CHKERRQ(ierr);
744*43e71028SBarry Smith     *flag = PETSC_FALSE;
745*43e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
746*43e71028SBarry Smith     goto theend2;
747*43e71028SBarry Smith   }
74878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
749cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
750*43e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7515d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
752393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
753b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
754ad922ac9SBarry Smith     goto theend2;
7555e42470aSBarry Smith   }
7565e42470aSBarry Smith 
7575e42470aSBarry Smith   /* Fit points with quadratic */
758313b5538SBarry Smith   lambda = 1.0;
7595e42470aSBarry Smith   count = 1;
7605ca10a37SBarry Smith   while (PETSC_TRUE) {
7615e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
76277431f27SBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count);
763b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
764f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
765*43e71028SBarry Smith       *flag = PETSC_FALSE;
766*43e71028SBarry Smith       break;
7675e42470aSBarry Smith     }
768a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
769329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
770329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
771329f5518SBarry Smith     else                         lambda     = lambdatemp;
772393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7733505fcc1SBarry Smith     lambdaneg = -lambda;
774aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7753505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7765e42470aSBarry Smith #else
7773505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7785e42470aSBarry Smith #endif
779*43e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
780*43e71028SBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count);
781*43e71028SBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
782*43e71028SBarry Smith       ierr = VecCopy(w,y);CHKERRQ(ierr);
783ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
784*43e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
785ed98166eSMatthew Knepley       break;
786ed98166eSMatthew Knepley     }
787*43e71028SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
788cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
789*43e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7905ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
791393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
792b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
79306259719SBarry Smith       break;
7945e42470aSBarry Smith     }
7955e42470aSBarry Smith     count++;
7965e42470aSBarry Smith   }
797ad922ac9SBarry Smith   theend2:
7988f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7998f99978dSLois Curfman McInnes   if (neP->CheckStep) {
8008f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
801abc0a331SBarry Smith     if (change_y) { /* recompute the function if the step has changed */
8028f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
8038f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
804*43e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
805*43e71028SBarry Smith       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
8068f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
8078f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8088f99978dSLois Curfman McInnes     }
8098f99978dSLois Curfman McInnes   }
810d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8113a40ed3dSBarry Smith   PetscFunctionReturn(0);
8125e76c431SBarry Smith }
8132343ba6eSBarry Smith 
81404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8154a2ae208SSatish Balay #undef __FUNCT__
8164a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
817c9e6a524SLois Curfman McInnes /*@C
81877c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
8194b27c08aSLois Curfman McInnes    by the method SNESLS.
8205e76c431SBarry Smith 
8215e76c431SBarry Smith    Input Parameters:
822c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
823af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
824c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8255e76c431SBarry Smith 
826fee21e36SBarry Smith    Collective on SNES
827fee21e36SBarry Smith 
828c4a48953SLois Curfman McInnes    Available Routines:
829c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
830f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
831f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
832af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8335e76c431SBarry Smith 
834c4a48953SLois Curfman McInnes     Options Database Keys:
8354b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8364b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8374b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8384b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8393304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8403304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
841c4a48953SLois Curfman McInnes 
8425e76c431SBarry Smith    Calling sequence of func:
843bd208895SLois Curfman McInnes .vb
844af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
845a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
846bd208895SLois Curfman McInnes .ve
8475e76c431SBarry Smith 
8485e76c431SBarry Smith     Input parameters for func:
849c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
850af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8515e76c431SBarry Smith .   x - current iterate
8525e76c431SBarry Smith .   f - residual evaluated at x
8535e76c431SBarry Smith .   y - search direction (contains new iterate on output)
8545e76c431SBarry Smith .   w - work vector
855c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8565e76c431SBarry Smith 
8575e76c431SBarry Smith     Output parameters for func:
858c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8595e76c431SBarry Smith .   y - new iterate (contains search direction on input)
8605e76c431SBarry Smith .   gnorm - 2-norm of g
8615e76c431SBarry Smith .   ynorm - 2-norm of search length
862a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
863f59ffdedSLois Curfman McInnes 
86436851e7fSLois Curfman McInnes     Level: advanced
86536851e7fSLois Curfman McInnes 
866f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
867f59ffdedSLois Curfman McInnes 
8685d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8695d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8705e76c431SBarry Smith @*/
871a7cc72afSBarry Smith PetscErrorCode SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
8725e76c431SBarry Smith {
873a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
87482bf6240SBarry Smith 
8753a40ed3dSBarry Smith   PetscFunctionBegin;
876c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
87782bf6240SBarry Smith   if (f) {
878af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
87982bf6240SBarry Smith   }
8803a40ed3dSBarry Smith   PetscFunctionReturn(0);
8815e76c431SBarry Smith }
8828e019c35SBarry Smith 
883a7cc72afSBarry 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*/
88404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
885fb2e594dSBarry Smith EXTERN_C_BEGIN
8864a2ae208SSatish Balay #undef __FUNCT__
8874a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
888dfbe8321SBarry Smith PetscErrorCode SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
88982bf6240SBarry Smith {
89082bf6240SBarry Smith   PetscFunctionBegin;
8914b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
8924b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
89382bf6240SBarry Smith   PetscFunctionReturn(0);
89482bf6240SBarry Smith }
895fb2e594dSBarry Smith EXTERN_C_END
89604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8974a2ae208SSatish Balay #undef __FUNCT__
8984a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
899c8dd0c0dSLois Curfman McInnes /*@C
900530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
9014b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
902c8dd0c0dSLois Curfman McInnes 
903c8dd0c0dSLois Curfman McInnes    Input Parameters:
904c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
905c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
906c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
907c8dd0c0dSLois Curfman McInnes 
908c8dd0c0dSLois Curfman McInnes    Collective on SNES
909c8dd0c0dSLois Curfman McInnes 
910c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
911c8dd0c0dSLois Curfman McInnes .vb
912b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
913c8dd0c0dSLois Curfman McInnes .ve
914b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
915b1ae27eaSLois Curfman McInnes    on failure.
916c8dd0c0dSLois Curfman McInnes 
917c8dd0c0dSLois Curfman McInnes    Input parameters for func:
918c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
919c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9208f99978dSLois Curfman McInnes -  x - current candidate iterate
921c8dd0c0dSLois Curfman McInnes 
922c8dd0c0dSLois Curfman McInnes    Output parameters for func:
923c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
924c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
925c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
926c8dd0c0dSLois Curfman McInnes 
927c8dd0c0dSLois Curfman McInnes    Level: advanced
928c8dd0c0dSLois Curfman McInnes 
9298f99978dSLois Curfman McInnes    Notes:
930b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
931b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
932b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
933b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
934ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9358f99978dSLois Curfman McInnes 
936b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
937b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
938b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
939ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
940ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
941ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
942ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
943b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9448f99978dSLois Curfman McInnes 
945c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
946c8dd0c0dSLois Curfman McInnes 
947c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
948c8dd0c0dSLois Curfman McInnes @*/
9496849ba73SBarry Smith PetscErrorCode SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
950c8dd0c0dSLois Curfman McInnes {
9516849ba73SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*);
952c8dd0c0dSLois Curfman McInnes 
953c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
954c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
955c8dd0c0dSLois Curfman McInnes   if (f) {
956c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
957c8dd0c0dSLois Curfman McInnes   }
958c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
959c8dd0c0dSLois Curfman McInnes }
960c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9616849ba73SBarry Smith typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/
962c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
9634a2ae208SSatish Balay #undef __FUNCT__
9644a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
965dfbe8321SBarry Smith PetscErrorCode SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx)
966c8dd0c0dSLois Curfman McInnes {
967c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9684b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
9694b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP    = checkctx;
970c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
971c8dd0c0dSLois Curfman McInnes }
972c8dd0c0dSLois Curfman McInnes EXTERN_C_END
973c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
97404d965bbSLois Curfman McInnes /*
9754b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
976329e7e40SMatthew Knepley 
977329e7e40SMatthew Knepley    Input Parameter:
978329e7e40SMatthew Knepley .  snes - the SNES context
979329e7e40SMatthew Knepley 
980329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
981329e7e40SMatthew Knepley */
982329e7e40SMatthew Knepley #undef __FUNCT__
9834b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
9846849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
985329e7e40SMatthew Knepley {
9864b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
987329e7e40SMatthew Knepley 
988329e7e40SMatthew Knepley   PetscFunctionBegin;
9894b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
9904b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
9914b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
9924b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
9934b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
994329e7e40SMatthew Knepley   PetscFunctionReturn(0);
995329e7e40SMatthew Knepley }
996329e7e40SMatthew Knepley 
997329e7e40SMatthew Knepley /*
9984b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
99904d965bbSLois Curfman McInnes 
100004d965bbSLois Curfman McInnes    Input Parameters:
100104d965bbSLois Curfman McInnes .  SNES - the SNES context
100204d965bbSLois Curfman McInnes .  viewer - visualization context
100304d965bbSLois Curfman McInnes 
100404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
100504d965bbSLois Curfman McInnes */
10064a2ae208SSatish Balay #undef __FUNCT__
10074b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10086849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1009a935fc98SLois Curfman McInnes {
10104b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
10112fc52814SBarry Smith   const char *cstr;
1012dfbe8321SBarry Smith   PetscErrorCode ierr;
101332077d6dSBarry Smith   PetscTruth iascii;
1014a935fc98SLois Curfman McInnes 
10153a40ed3dSBarry Smith   PetscFunctionBegin;
101632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
101732077d6dSBarry Smith   if (iascii) {
101819bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
101919bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
102019bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
102119bcc07fSBarry Smith     else                                                cstr = "unknown";
1022b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1023b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
10245cd90555SBarry Smith   } else {
102579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
102619bcc07fSBarry Smith   }
10273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1028a935fc98SLois Curfman McInnes }
102904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
103004d965bbSLois Curfman McInnes /*
10314b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
103204d965bbSLois Curfman McInnes 
103304d965bbSLois Curfman McInnes    Input Parameter:
103404d965bbSLois Curfman McInnes .  snes - the SNES context
103504d965bbSLois Curfman McInnes 
103604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
103704d965bbSLois Curfman McInnes */
10384a2ae208SSatish Balay #undef __FUNCT__
10394b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
10406849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
10415e42470aSBarry Smith {
10424b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
1043e5999256SBarry Smith   const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1044dfbe8321SBarry Smith   PetscErrorCode ierr;
1045a7cc72afSBarry Smith   PetscInt indx;
1046f1af5d2fSBarry Smith   PetscTruth flg;
10475e42470aSBarry Smith 
10483a40ed3dSBarry Smith   PetscFunctionBegin;
1049b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
10504b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
10514b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
10524b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1053186905e3SBarry Smith 
105422e36657SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
105525cdf11fSBarry Smith     if (flg) {
105622e36657SBarry Smith       switch (indx) {
1057b49fd9e1SBarry Smith       case 0:
1058af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1059b49fd9e1SBarry Smith         break;
1060b49fd9e1SBarry Smith       case 1:
1061af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1062b49fd9e1SBarry Smith         break;
1063b49fd9e1SBarry Smith       case 2:
1064af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1065b49fd9e1SBarry Smith         break;
1066b49fd9e1SBarry Smith       case 3:
1067af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1068b49fd9e1SBarry Smith         break;
10695e42470aSBarry Smith       }
10705e42470aSBarry Smith     }
1071b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10723a40ed3dSBarry Smith   PetscFunctionReturn(0);
10735e42470aSBarry Smith }
107404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1075ebe3b25bSBarry Smith /*MC
1076ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
107704d965bbSLois Curfman McInnes 
1078ebe3b25bSBarry Smith    Options Database:
1079ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1080ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1081ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1082ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1083ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1084ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
108504d965bbSLois Curfman McInnes 
1086ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
108704d965bbSLois Curfman McInnes 
1088ee3001cbSBarry Smith    Level: beginner
1089ee3001cbSBarry Smith 
1090ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1091ebe3b25bSBarry Smith            SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1092ebe3b25bSBarry Smith           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1093ebe3b25bSBarry Smith 
1094ebe3b25bSBarry Smith M*/
1095fb2e594dSBarry Smith EXTERN_C_BEGIN
10964a2ae208SSatish Balay #undef __FUNCT__
10974b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
1098dfbe8321SBarry Smith PetscErrorCode SNESCreate_LS(SNES snes)
10995e42470aSBarry Smith {
1100dfbe8321SBarry Smith   PetscErrorCode ierr;
11014b27c08aSLois Curfman McInnes   SNES_LS *neP;
11025e42470aSBarry Smith 
11033a40ed3dSBarry Smith   PetscFunctionBegin;
11044b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
11054b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
11064b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
11074b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
11084b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
11094b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
11104b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
11115baf8537SBarry Smith   snes->nwork           = 0;
11125e42470aSBarry Smith 
11134b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
111452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
11155e42470aSBarry Smith   snes->data    	= (void*)neP;
11165e42470aSBarry Smith   neP->alpha		= 1.e-4;
11175e42470aSBarry Smith   neP->maxstep		= 1.e8;
11185e42470aSBarry Smith   neP->steptol		= 1.e-12;
11195e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1120c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1121c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1122c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
112382bf6240SBarry Smith 
11245d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
11255d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
112682bf6240SBarry Smith 
11273a40ed3dSBarry Smith   PetscFunctionReturn(0);
11285e42470aSBarry Smith }
1129fb2e594dSBarry Smith EXTERN_C_END
113082bf6240SBarry Smith 
113182bf6240SBarry Smith 
113282bf6240SBarry Smith 
1133