xref: /petsc/src/snes/impls/ls/ls.c (revision 32077d6df18512e94702c86e7d562481ed0cebd0)
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"
128a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal  a1;
158a5d9ee4SBarry Smith   int        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"
5074637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5174637425SBarry Smith {
5274637425SBarry Smith   PetscReal     a1,a2;
5374637425SBarry Smith   int           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"
132c9780b6fSBarry Smith int SNESSolve_LS(SNES snes)
1335e76c431SBarry Smith {
1344b27c08aSLois Curfman McInnes   SNES_LS      *neP = (SNES_LS*)snes->data;
13584c569c9SBarry Smith   int          maxits,i,ierr,lits,lsfail;
136112a2221SBarry Smith   MatStructure flg = DIFFERENT_NONZERO_PATTERN;
137329f5518SBarry Smith   PetscReal    fnorm,gnorm,xnorm,ynorm;
1385e42470aSBarry Smith   Vec          Y,X,F,G,W,TMP;
139c293cc10SBarry Smith   KSP          ksp;
1405e76c431SBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
14294b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
143184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
144184914b5SBarry Smith 
1455e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1465e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14739e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1485e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1495e42470aSBarry Smith   G		= snes->work[1];
1505e42470aSBarry Smith   W		= snes->work[2];
1515e76c431SBarry Smith 
1524c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1534c49b128SBarry Smith   snes->iter = 0;
1544c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1555334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
156cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1570f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1585e42470aSBarry Smith   snes->norm = fnorm;
1590f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16084c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16194a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1625e76c431SBarry Smith 
163c9780b6fSBarry Smith   if (fnorm < snes->atol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16494a9d846SBarry Smith 
165d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
166d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
167d034289bSBarry Smith 
1685e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1695e76c431SBarry Smith 
170329e7e40SMatthew Knepley     /* Call general purpose update function */
171de8cb200SMatthew Knepley     if (snes->update != PETSC_NULL) {
172329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
173de8cb200SMatthew Knepley     }
174329e7e40SMatthew Knepley 
175ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1765334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
17794b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
17823ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
179c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
18074637425SBarry Smith 
181b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
18274637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
18385471664SBarry Smith     }
18474637425SBarry Smith 
18590cbd584SBarry Smith     /* should check what happened to the linear solve? */
1863505fcc1SBarry Smith     snes->linear_its += lits;
1874b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
188ea4d3ed3SLois Curfman McInnes 
189ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
190ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
191ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
192ea4d3ed3SLois Curfman McInnes     */
19381b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
194af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
1954b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
196a3b891d8SBarry Smith 
197a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
198a3b891d8SBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
199a3b891d8SBarry Smith     fnorm = gnorm;
200a3b891d8SBarry Smith 
2015ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2025ed2d596SBarry Smith     snes->iter = i+1;
2035ed2d596SBarry Smith     snes->norm = fnorm;
2045ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2055ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2065ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2075ed2d596SBarry Smith 
2083505fcc1SBarry Smith     if (lsfail) {
2098a5d9ee4SBarry Smith       PetscTruth ismin;
21050ffb88aSMatthew Knepley 
21150ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2123505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2138a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2148a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2153505fcc1SBarry Smith         break;
2163505fcc1SBarry Smith       }
21750ffb88aSMatthew Knepley     }
2185e76c431SBarry Smith 
2195e76c431SBarry Smith     /* Test for convergence */
22029e0b56fSBarry Smith     if (snes->converged) {
22129e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2223505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2233505fcc1SBarry Smith       if (snes->reason) {
22416c95cacSBarry Smith         break;
22516c95cacSBarry Smith       }
22616c95cacSBarry Smith     }
22729e0b56fSBarry Smith   }
22839e2f89bSBarry Smith   if (X != snes->vec_sol) {
229393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
230e15f7bb5SBarry Smith   }
231e15f7bb5SBarry Smith   if (F != snes->vec_func) {
232e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
233e15f7bb5SBarry Smith   }
23439e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
23539e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
23652392280SLois Curfman McInnes   if (i == maxits) {
2374b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %d\n",maxits);
2383505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
23952392280SLois Curfman McInnes   }
2403a40ed3dSBarry Smith   PetscFunctionReturn(0);
2415e76c431SBarry Smith }
24204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
24304d965bbSLois Curfman McInnes /*
2444b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2454b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
24604d965bbSLois Curfman McInnes 
24704d965bbSLois Curfman McInnes    Input Parameter:
24804d965bbSLois Curfman McInnes .  snes - the SNES context
24904d965bbSLois Curfman McInnes .  x - the solution vector
25004d965bbSLois Curfman McInnes 
25104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25204d965bbSLois Curfman McInnes 
25304d965bbSLois Curfman McInnes    Notes:
2544b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
25504d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
25604d965bbSLois Curfman McInnes    the call to SNESSolve().
25704d965bbSLois Curfman McInnes  */
2584a2ae208SSatish Balay #undef __FUNCT__
2594b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
2604b27c08aSLois Curfman McInnes int SNESSetUp_LS(SNES snes)
2615e76c431SBarry Smith {
2625e42470aSBarry Smith   int ierr;
2633a40ed3dSBarry Smith 
2643a40ed3dSBarry Smith   PetscFunctionBegin;
26581b6cf68SLois Curfman McInnes   snes->nwork = 4;
266d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
267b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
26881b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2693a40ed3dSBarry Smith   PetscFunctionReturn(0);
2705e76c431SBarry Smith }
27104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27204d965bbSLois Curfman McInnes /*
2734b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2744b27c08aSLois Curfman McInnes    with SNESCreate_LS().
27504d965bbSLois Curfman McInnes 
27604d965bbSLois Curfman McInnes    Input Parameter:
27704d965bbSLois Curfman McInnes .  snes - the SNES context
27804d965bbSLois Curfman McInnes 
279de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28004d965bbSLois Curfman McInnes  */
2814a2ae208SSatish Balay #undef __FUNCT__
2824b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
2834b27c08aSLois Curfman McInnes int SNESDestroy_LS(SNES snes)
2845e76c431SBarry Smith {
285393d2d9aSLois Curfman McInnes   int  ierr;
2863a40ed3dSBarry Smith 
2873a40ed3dSBarry Smith   PetscFunctionBegin;
2885baf8537SBarry Smith   if (snes->nwork) {
2894b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29021c89e3eSBarry Smith   }
2915c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2923a40ed3dSBarry Smith   PetscFunctionReturn(0);
2935e76c431SBarry Smith }
29404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
29782bf6240SBarry Smith 
2984b828684SBarry Smith /*@C
2995e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3005e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3015e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3025e76c431SBarry Smith 
303c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
304c7afd0dbSLois Curfman McInnes 
3055e76c431SBarry Smith    Input Parameters:
306c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
307af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3085e76c431SBarry Smith .  x - current iterate
3095e76c431SBarry Smith .  f - residual evaluated at x
3105e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3115e76c431SBarry Smith .  w - work vector
312c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3135e76c431SBarry Smith 
314c4a48953SLois Curfman McInnes    Output Parameters:
315c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3165e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3175e76c431SBarry Smith .  gnorm - 2-norm of g
3185e76c431SBarry Smith .  ynorm - 2-norm of search length
319c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
320fee21e36SBarry Smith 
321c4a48953SLois Curfman McInnes    Options Database Key:
3224b27c08aSLois Curfman McInnes .  -snes_ls basic - Activates SNESNoLineSearch()
323c4a48953SLois Curfman McInnes 
32436851e7fSLois Curfman McInnes    Level: advanced
32536851e7fSLois Curfman McInnes 
32628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
32728ae5a21SLois Curfman McInnes 
328f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
329af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3305e76c431SBarry Smith @*/
3315d1a10b1SBarry Smith int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
3325e76c431SBarry Smith {
3335e42470aSBarry Smith   int         ierr;
334ea709b57SSatish Balay   PetscScalar mone = -1.0;
3354b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
3368f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
3375334005bSBarry Smith 
3383a40ed3dSBarry Smith   PetscFunctionBegin;
339761aaf1bSLois Curfman McInnes   *flag = 0;
340d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
341a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
342ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3438f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3448f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3458f99978dSLois Curfman McInnes   }
346ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
347a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
348d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3493a40ed3dSBarry Smith   PetscFunctionReturn(0);
3505e76c431SBarry Smith }
35104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3524a2ae208SSatish Balay #undef __FUNCT__
3534a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
35482bf6240SBarry Smith 
35529e0b56fSBarry Smith /*@C
35629e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
35729e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
35829e0b56fSBarry Smith    even compute the norm of the function or search direction; this
35929e0b56fSBarry Smith    is intended only when you know the full step is fine and are
36029e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
361c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
362c7afd0dbSLois Curfman McInnes 
363c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
36429e0b56fSBarry Smith 
36529e0b56fSBarry Smith    Input Parameters:
366c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
367af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
36829e0b56fSBarry Smith .  x - current iterate
36929e0b56fSBarry Smith .  f - residual evaluated at x
37029e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
37129e0b56fSBarry Smith .  w - work vector
372c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
37329e0b56fSBarry Smith 
37429e0b56fSBarry Smith    Output Parameters:
375c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
37629e0b56fSBarry Smith .  gnorm - not changed
37729e0b56fSBarry Smith .  ynorm - not changed
378c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
379fee21e36SBarry Smith 
38029e0b56fSBarry Smith    Options Database Key:
3814b27c08aSLois Curfman McInnes .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
38229e0b56fSBarry Smith 
3838cbba510SBarry Smith    Notes:
384ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
385ea56f5baSLois Curfman McInnes    either the options
386ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
387ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
388329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
389329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
390329f5518SBarry Smith 
391329f5518SBarry Smith    During the final iteration this will not evaluate the function at
392329f5518SBarry Smith    the solution point. This is to save a function evaluation while
393329f5518SBarry Smith    using pseudo-timestepping.
3948cbba510SBarry Smith 
395ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
396ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
397ea56f5baSLois Curfman McInnes    correct, since they are not computed.
398ea56f5baSLois Curfman McInnes 
399ea56f5baSLois Curfman McInnes    Level: advanced
4008cbba510SBarry Smith 
40129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
40229e0b56fSBarry Smith 
40329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
40429e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
40529e0b56fSBarry Smith @*/
4065d1a10b1SBarry Smith int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
40729e0b56fSBarry Smith {
40829e0b56fSBarry Smith   int         ierr;
409ea709b57SSatish Balay   PetscScalar mone = -1.0;
4104b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
4118f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
41229e0b56fSBarry Smith 
4133a40ed3dSBarry Smith   PetscFunctionBegin;
4148cbba510SBarry Smith   *flag = 0;
415d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
41629e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4178f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4188f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4198f99978dSLois Curfman McInnes   }
420329f5518SBarry Smith 
421329f5518SBarry Smith   /* don't evaluate function the last time through */
422329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
42329e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
424329f5518SBarry Smith   }
425d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
42729e0b56fSBarry Smith }
42804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4294a2ae208SSatish Balay #undef __FUNCT__
4304a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4314b828684SBarry Smith /*@C
432f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4335e76c431SBarry Smith 
434c7afd0dbSLois Curfman McInnes    Collective on SNES
435c7afd0dbSLois Curfman McInnes 
4365e76c431SBarry Smith    Input Parameters:
437c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
438af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4395e76c431SBarry Smith .  x - current iterate
4405e76c431SBarry Smith .  f - residual evaluated at x
4415e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4425e76c431SBarry Smith .  w - work vector
443c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4445e76c431SBarry Smith 
445393d2d9aSLois Curfman McInnes    Output Parameters:
446c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4475e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4485e76c431SBarry Smith .  gnorm - 2-norm of g
4495e76c431SBarry Smith .  ynorm - 2-norm of search length
450c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
451fee21e36SBarry Smith 
452c4a48953SLois Curfman McInnes    Options Database Key:
4534b27c08aSLois Curfman McInnes $  -snes_ls cubic - Activates SNESCubicLineSearch()
454c4a48953SLois Curfman McInnes 
4555e76c431SBarry Smith    Notes:
4565e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4575e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4585e76c431SBarry Smith 
45936851e7fSLois Curfman McInnes    Level: advanced
46036851e7fSLois Curfman McInnes 
46128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
46228ae5a21SLois Curfman McInnes 
463af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4645e76c431SBarry Smith @*/
4655d1a10b1SBarry Smith int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
4665e76c431SBarry Smith {
467406556e6SLois Curfman McInnes   /*
468406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
469406556e6SLois Curfman McInnes      minimization problem:
470406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
471406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
472406556e6SLois Curfman McInnes    */
473406556e6SLois Curfman McInnes 
4745ca10a37SBarry Smith   PetscReal   steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
475329f5518SBarry Smith   PetscReal   maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
476aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
477ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
4786b5873e3SBarry Smith #endif
4795e42470aSBarry Smith   int         ierr,count;
4804b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
481ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
4828f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
4835e76c431SBarry Smith 
4843a40ed3dSBarry Smith   PetscFunctionBegin;
485d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
486761aaf1bSLois Curfman McInnes   *flag   = 0;
4875e76c431SBarry Smith   alpha   = neP->alpha;
4885e76c431SBarry Smith   maxstep = neP->maxstep;
4895e76c431SBarry Smith   steptol = neP->steptol;
4905e76c431SBarry Smith 
491cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
49263b9fa5eSBarry Smith   if (*ynorm == 0.0) {
49363b9fa5eSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n");
494a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
495a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
496a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
49763b9fa5eSBarry Smith     *flag  = -1;
498ad922ac9SBarry Smith     goto theend1;
49994a9d846SBarry Smith   }
5005e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5015e42470aSBarry Smith     scale = maxstep/(*ynorm);
502aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
503b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
5046b5873e3SBarry Smith #else
505b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
5066b5873e3SBarry Smith #endif
507393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5085e76c431SBarry Smith     *ynorm = maxstep;
5095e76c431SBarry Smith   }
510ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5115ca10a37SBarry Smith   minlambda = steptol/rellength;
512a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
513aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
514a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
515329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5165e42470aSBarry Smith #else
517a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5185e42470aSBarry Smith #endif
5195e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5205e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5215e76c431SBarry Smith 
522393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5235334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
52478b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
525313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5265d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
527393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
528b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
529ad922ac9SBarry Smith     goto theend1;
5305e76c431SBarry Smith   }
5315e76c431SBarry Smith 
5325e76c431SBarry Smith   /* Fit points with quadratic */
533313b5538SBarry Smith   lambda = 1.0;
534a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5355e76c431SBarry Smith   lambdaprev = lambda;
5365e76c431SBarry Smith   gnormprev = *gnorm;
537329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
538ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
539ddd12767SBarry Smith   else                         lambda = lambdatemp;
540393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
541ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
542aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
543ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5445e42470aSBarry Smith #else
545ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5465e42470aSBarry Smith #endif
54778b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
548cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5495ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
550393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
5515ed2d596SBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda);
552ad922ac9SBarry Smith     goto theend1;
5535e76c431SBarry Smith   }
5545e76c431SBarry Smith 
5555e76c431SBarry Smith   /* Fit points with cubic */
5565e76c431SBarry Smith   count = 1;
5578229bfc2SMatthew Knepley   while (count < 10000) {
5585e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
559b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
5605ed2d596SBarry 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);
561f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
562761aaf1bSLois Curfman McInnes       *flag = -1; break;
5635e76c431SBarry Smith     }
5645d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
5655d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
566ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5672b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5685e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5695e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5705e76c431SBarry Smith     if (a == 0.0) {
5715e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5725e76c431SBarry Smith     } else {
5735e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5745e76c431SBarry Smith     }
5755e76c431SBarry Smith     lambdaprev = lambda;
5765e76c431SBarry Smith     gnormprev  = *gnorm;
577329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
578329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
5795e76c431SBarry Smith     else                         lambda     = lambdatemp;
580393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
581ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
582aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
583ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
584393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5855e42470aSBarry Smith #else
586ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5875e42470aSBarry Smith #endif
58878b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
589cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5905ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
591393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
5925ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda);
5932715565aSLois Curfman McInnes       goto theend1;
5942b022350SLois Curfman McInnes     } else {
5955ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);
5965e76c431SBarry Smith     }
5975e76c431SBarry Smith     count++;
5985e76c431SBarry Smith   }
5998229bfc2SMatthew Knepley   if (count >= 10000) {
600cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6018229bfc2SMatthew Knepley   }
602ad922ac9SBarry Smith   theend1:
6038f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6048f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6058f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6068f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6078f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6088f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6098f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6108f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6118f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6128f99978dSLois Curfman McInnes     }
6138f99978dSLois Curfman McInnes   }
614d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6153a40ed3dSBarry Smith   PetscFunctionReturn(0);
6165e76c431SBarry Smith }
61704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6184a2ae208SSatish Balay #undef __FUNCT__
6194a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6204b828684SBarry Smith /*@C
621f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6225e76c431SBarry Smith 
623c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
624c7afd0dbSLois Curfman McInnes 
6255e42470aSBarry Smith    Input Parameters:
626c7afd0dbSLois Curfman McInnes +  snes - the SNES context
627af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6285e42470aSBarry Smith .  x - current iterate
6295e42470aSBarry Smith .  f - residual evaluated at x
6305e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6315e42470aSBarry Smith .  w - work vector
632c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6335e42470aSBarry Smith 
634c4a48953SLois Curfman McInnes    Output Parameters:
635c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6365e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6375e42470aSBarry Smith .  gnorm - 2-norm of g
6385e42470aSBarry Smith .  ynorm - 2-norm of search length
639c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
640fee21e36SBarry Smith 
641c4a48953SLois Curfman McInnes    Options Database Key:
6424b27c08aSLois Curfman McInnes .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
643c4a48953SLois Curfman McInnes 
6445e42470aSBarry Smith    Notes:
6454b27c08aSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNESLS method.
6465e42470aSBarry Smith 
64736851e7fSLois Curfman McInnes    Level: advanced
64836851e7fSLois Curfman McInnes 
649f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
650f59ffdedSLois Curfman McInnes 
651af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6525e42470aSBarry Smith @*/
6535d1a10b1SBarry Smith int SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
6545e76c431SBarry Smith {
655406556e6SLois Curfman McInnes   /*
656406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
657406556e6SLois Curfman McInnes      minimization problem:
658406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
659406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
660406556e6SLois Curfman McInnes    */
6615ca10a37SBarry Smith   PetscReal   steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
662aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
663ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
6646b5873e3SBarry Smith #endif
6655e42470aSBarry Smith   int         ierr,count;
6664b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
667ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
6688f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
6695e76c431SBarry Smith 
6703a40ed3dSBarry Smith   PetscFunctionBegin;
671d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
672761aaf1bSLois Curfman McInnes   *flag   = 0;
6735e42470aSBarry Smith   alpha   = neP->alpha;
6745e42470aSBarry Smith   maxstep = neP->maxstep;
6755e42470aSBarry Smith   steptol = neP->steptol;
6765e76c431SBarry Smith 
6773505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
67863b9fa5eSBarry Smith   if (*ynorm == 0.0) {
679b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
680b37302c6SLois Curfman McInnes     *gnorm = fnorm;
681b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
682b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
68363b9fa5eSBarry Smith     *flag  = -1;
684ad922ac9SBarry Smith     goto theend2;
68594a9d846SBarry Smith   }
6865e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6875e42470aSBarry Smith     scale = maxstep/(*ynorm);
688393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6895e42470aSBarry Smith     *ynorm = maxstep;
6905e76c431SBarry Smith   }
691ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
6925ca10a37SBarry Smith   minlambda = steptol/rellength;
693a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
694aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
695a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
696329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
6975e42470aSBarry Smith #else
698a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6995e42470aSBarry Smith #endif
7005e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
7015e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7025e42470aSBarry Smith 
703393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7045334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
70578b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
706cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7075d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
708393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
709b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
710ad922ac9SBarry Smith     goto theend2;
7115e42470aSBarry Smith   }
7125e42470aSBarry Smith 
7135e42470aSBarry Smith   /* Fit points with quadratic */
714313b5538SBarry Smith   lambda = 1.0;
7155e42470aSBarry Smith   count = 1;
7165ca10a37SBarry Smith   while (PETSC_TRUE) {
7175e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
718b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
719b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
720f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
721761aaf1bSLois Curfman McInnes       *flag = -1; break;
7225e42470aSBarry Smith     }
723a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
724329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
725329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
726329f5518SBarry Smith     else                         lambda     = lambdatemp;
727393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7283505fcc1SBarry Smith     lambdaneg = -lambda;
729aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7303505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7315e42470aSBarry Smith #else
7323505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7335e42470aSBarry Smith #endif
73478b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
735cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7365ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
737393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
738b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
73906259719SBarry Smith       break;
7405e42470aSBarry Smith     }
7415e42470aSBarry Smith     count++;
7425e42470aSBarry Smith   }
743ad922ac9SBarry Smith   theend2:
7448f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7458f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7468f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7478f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7488f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7498f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7508f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7518f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7528f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7538f99978dSLois Curfman McInnes     }
7548f99978dSLois Curfman McInnes   }
755d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7563a40ed3dSBarry Smith   PetscFunctionReturn(0);
7575e76c431SBarry Smith }
7582343ba6eSBarry Smith 
75904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7604a2ae208SSatish Balay #undef __FUNCT__
7614a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
762c9e6a524SLois Curfman McInnes /*@C
76377c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7644b27c08aSLois Curfman McInnes    by the method SNESLS.
7655e76c431SBarry Smith 
7665e76c431SBarry Smith    Input Parameters:
767c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
768af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
769c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7705e76c431SBarry Smith 
771fee21e36SBarry Smith    Collective on SNES
772fee21e36SBarry Smith 
773c4a48953SLois Curfman McInnes    Available Routines:
774c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
775f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
776f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
777af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7785e76c431SBarry Smith 
779c4a48953SLois Curfman McInnes     Options Database Keys:
7804b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
7814b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
7824b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
7834b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
7843304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
7853304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
786c4a48953SLois Curfman McInnes 
7875e76c431SBarry Smith    Calling sequence of func:
788bd208895SLois Curfman McInnes .vb
789af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
790329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
791bd208895SLois Curfman McInnes .ve
7925e76c431SBarry Smith 
7935e76c431SBarry Smith     Input parameters for func:
794c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
795af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
7965e76c431SBarry Smith .   x - current iterate
7975e76c431SBarry Smith .   f - residual evaluated at x
7985e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7995e76c431SBarry Smith .   w - work vector
800c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8015e76c431SBarry Smith 
8025e76c431SBarry Smith     Output parameters for func:
803c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8045e76c431SBarry Smith .   y - new iterate (contains search direction on input)
8055e76c431SBarry Smith .   gnorm - 2-norm of g
8065e76c431SBarry Smith .   ynorm - 2-norm of search length
807c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
808761aaf1bSLois Curfman McInnes            on failure.
809f59ffdedSLois Curfman McInnes 
81036851e7fSLois Curfman McInnes     Level: advanced
81136851e7fSLois Curfman McInnes 
812f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
813f59ffdedSLois Curfman McInnes 
8145d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8155d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8165e76c431SBarry Smith @*/
817329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
8185e76c431SBarry Smith {
819329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
82082bf6240SBarry Smith 
8213a40ed3dSBarry Smith   PetscFunctionBegin;
822c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
82382bf6240SBarry Smith   if (f) {
824af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
82582bf6240SBarry Smith   }
8263a40ed3dSBarry Smith   PetscFunctionReturn(0);
8275e76c431SBarry Smith }
8288e019c35SBarry Smith 
8298e019c35SBarry Smith typedef int (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*); /* force argument to next function to not be extern C*/
83004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
831fb2e594dSBarry Smith EXTERN_C_BEGIN
8324a2ae208SSatish Balay #undef __FUNCT__
8334a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
8348e019c35SBarry Smith int SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
83582bf6240SBarry Smith {
83682bf6240SBarry Smith   PetscFunctionBegin;
8374b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
8384b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
83982bf6240SBarry Smith   PetscFunctionReturn(0);
84082bf6240SBarry Smith }
841fb2e594dSBarry Smith EXTERN_C_END
84204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8434a2ae208SSatish Balay #undef __FUNCT__
8444a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
845c8dd0c0dSLois Curfman McInnes /*@C
846530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8474b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
848c8dd0c0dSLois Curfman McInnes 
849c8dd0c0dSLois Curfman McInnes    Input Parameters:
850c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
851c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
852c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
853c8dd0c0dSLois Curfman McInnes 
854c8dd0c0dSLois Curfman McInnes    Collective on SNES
855c8dd0c0dSLois Curfman McInnes 
856c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
857c8dd0c0dSLois Curfman McInnes .vb
858b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
859c8dd0c0dSLois Curfman McInnes .ve
860b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
861b1ae27eaSLois Curfman McInnes    on failure.
862c8dd0c0dSLois Curfman McInnes 
863c8dd0c0dSLois Curfman McInnes    Input parameters for func:
864c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
865c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8668f99978dSLois Curfman McInnes -  x - current candidate iterate
867c8dd0c0dSLois Curfman McInnes 
868c8dd0c0dSLois Curfman McInnes    Output parameters for func:
869c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
870c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
871c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
872c8dd0c0dSLois Curfman McInnes 
873c8dd0c0dSLois Curfman McInnes    Level: advanced
874c8dd0c0dSLois Curfman McInnes 
8758f99978dSLois Curfman McInnes    Notes:
876b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
877b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
878b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
879b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
880ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8818f99978dSLois Curfman McInnes 
882b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
883b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
884b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
885ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
886ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
887ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
888ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
889b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8908f99978dSLois Curfman McInnes 
891c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
892c8dd0c0dSLois Curfman McInnes 
893c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
894c8dd0c0dSLois Curfman McInnes @*/
8958f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
896c8dd0c0dSLois Curfman McInnes {
8978f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
898c8dd0c0dSLois Curfman McInnes 
899c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
900c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
901c8dd0c0dSLois Curfman McInnes   if (f) {
902c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
903c8dd0c0dSLois Curfman McInnes   }
904c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
905c8dd0c0dSLois Curfman McInnes }
906c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9072343ba6eSBarry Smith typedef int (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/
908c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
9094a2ae208SSatish Balay #undef __FUNCT__
9104a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
9112343ba6eSBarry Smith int SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx)
912c8dd0c0dSLois Curfman McInnes {
913c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9144b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
9154b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP    = checkctx;
916c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
917c8dd0c0dSLois Curfman McInnes }
918c8dd0c0dSLois Curfman McInnes EXTERN_C_END
919c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
92004d965bbSLois Curfman McInnes /*
9214b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
922329e7e40SMatthew Knepley 
923329e7e40SMatthew Knepley    Input Parameter:
924329e7e40SMatthew Knepley .  snes - the SNES context
925329e7e40SMatthew Knepley 
926329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
927329e7e40SMatthew Knepley */
928329e7e40SMatthew Knepley #undef __FUNCT__
9294b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
9304b27c08aSLois Curfman McInnes static int SNESPrintHelp_LS(SNES snes,char *p)
931329e7e40SMatthew Knepley {
9324b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
933329e7e40SMatthew Knepley 
934329e7e40SMatthew Knepley   PetscFunctionBegin;
9354b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
9364b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
9374b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
9384b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
9394b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
940329e7e40SMatthew Knepley   PetscFunctionReturn(0);
941329e7e40SMatthew Knepley }
942329e7e40SMatthew Knepley 
943329e7e40SMatthew Knepley /*
9444b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
94504d965bbSLois Curfman McInnes 
94604d965bbSLois Curfman McInnes    Input Parameters:
94704d965bbSLois Curfman McInnes .  SNES - the SNES context
94804d965bbSLois Curfman McInnes .  viewer - visualization context
94904d965bbSLois Curfman McInnes 
95004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
95104d965bbSLois Curfman McInnes */
9524a2ae208SSatish Balay #undef __FUNCT__
9534b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
9544b27c08aSLois Curfman McInnes static int SNESView_LS(SNES snes,PetscViewer viewer)
955a935fc98SLois Curfman McInnes {
9564b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
9572fc52814SBarry Smith   const char *cstr;
95851695f54SLois Curfman McInnes   int        ierr;
959*32077d6dSBarry Smith   PetscTruth iascii;
960a935fc98SLois Curfman McInnes 
9613a40ed3dSBarry Smith   PetscFunctionBegin;
962*32077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
963*32077d6dSBarry Smith   if (iascii) {
96419bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
96519bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
96619bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
96719bcc07fSBarry Smith     else                                                cstr = "unknown";
968b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
969b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9705cd90555SBarry Smith   } else {
97129bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
97219bcc07fSBarry Smith   }
9733a40ed3dSBarry Smith   PetscFunctionReturn(0);
974a935fc98SLois Curfman McInnes }
97504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
97604d965bbSLois Curfman McInnes /*
9774b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
97804d965bbSLois Curfman McInnes 
97904d965bbSLois Curfman McInnes    Input Parameter:
98004d965bbSLois Curfman McInnes .  snes - the SNES context
98104d965bbSLois Curfman McInnes 
98204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
98304d965bbSLois Curfman McInnes */
9844a2ae208SSatish Balay #undef __FUNCT__
9854b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
9864b27c08aSLois Curfman McInnes static int SNESSetFromOptions_LS(SNES snes)
9875e42470aSBarry Smith {
9884b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
989e5999256SBarry Smith   const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
99022e36657SBarry Smith   int        ierr,indx;
991f1af5d2fSBarry Smith   PetscTruth flg;
9925e42470aSBarry Smith 
9933a40ed3dSBarry Smith   PetscFunctionBegin;
994b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
9954b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
9964b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
9974b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
998186905e3SBarry Smith 
99922e36657SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
100025cdf11fSBarry Smith     if (flg) {
100122e36657SBarry Smith       switch (indx) {
1002b49fd9e1SBarry Smith       case 0:
1003af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1004b49fd9e1SBarry Smith         break;
1005b49fd9e1SBarry Smith       case 1:
1006af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1007b49fd9e1SBarry Smith         break;
1008b49fd9e1SBarry Smith       case 2:
1009af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1010b49fd9e1SBarry Smith         break;
1011b49fd9e1SBarry Smith       case 3:
1012af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1013b49fd9e1SBarry Smith         break;
10145e42470aSBarry Smith       }
10155e42470aSBarry Smith     }
1016b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10173a40ed3dSBarry Smith   PetscFunctionReturn(0);
10185e42470aSBarry Smith }
101904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1020ebe3b25bSBarry Smith /*MC
1021ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
102204d965bbSLois Curfman McInnes 
1023ebe3b25bSBarry Smith    Options Database:
1024ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1025ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1026ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1027ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1028ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1029ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
103004d965bbSLois Curfman McInnes 
1031ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
103204d965bbSLois Curfman McInnes 
1033ee3001cbSBarry Smith    Level: beginner
1034ee3001cbSBarry Smith 
1035ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1036ebe3b25bSBarry Smith            SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1037ebe3b25bSBarry Smith           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1038ebe3b25bSBarry Smith 
1039ebe3b25bSBarry Smith M*/
1040fb2e594dSBarry Smith EXTERN_C_BEGIN
10414a2ae208SSatish Balay #undef __FUNCT__
10424b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
10434b27c08aSLois Curfman McInnes int SNESCreate_LS(SNES snes)
10445e42470aSBarry Smith {
104582bf6240SBarry Smith   int     ierr;
10464b27c08aSLois Curfman McInnes   SNES_LS *neP;
10475e42470aSBarry Smith 
10483a40ed3dSBarry Smith   PetscFunctionBegin;
10494b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
10504b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
10514b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
10524b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
10534b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
10544b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
10554b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
10565baf8537SBarry Smith   snes->nwork           = 0;
10575e42470aSBarry Smith 
10584b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
10594b27c08aSLois Curfman McInnes   PetscLogObjectMemory(snes,sizeof(SNES_LS));
10605e42470aSBarry Smith   snes->data    	= (void*)neP;
10615e42470aSBarry Smith   neP->alpha		= 1.e-4;
10625e42470aSBarry Smith   neP->maxstep		= 1.e8;
10635e42470aSBarry Smith   neP->steptol		= 1.e-12;
10645e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1065c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1066c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1067c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
106882bf6240SBarry Smith 
10695d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
10705d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
107182bf6240SBarry Smith 
10723a40ed3dSBarry Smith   PetscFunctionReturn(0);
10735e42470aSBarry Smith }
1074fb2e594dSBarry Smith EXTERN_C_END
107582bf6240SBarry Smith 
107682bf6240SBarry Smith 
107782bf6240SBarry Smith 
1078