xref: /petsc/src/snes/impls/ls/ls.c (revision 2fc52814d27bf1f4e71021c1c3ebb532b583ed60)
173f4d377SMatthew Knepley /*$Id: ls.c,v 1.172 2001/08/07 03:04:11 balay Exp $*/
25e76c431SBarry Smith 
370f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
45e42470aSBarry Smith 
58a5d9ee4SBarry Smith /*
68a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
78a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
836669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
936669109SBarry Smith     for this trick.
108a5d9ee4SBarry Smith */
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
138a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
148a5d9ee4SBarry Smith {
158a5d9ee4SBarry Smith   PetscReal  a1;
168a5d9ee4SBarry Smith   int        ierr;
1736669109SBarry Smith   PetscTruth hastranspose;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2136669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2236669109SBarry Smith   if (hastranspose) {
238a5d9ee4SBarry Smith     /* Compute || J^T F|| */
248a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
258a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
264b27c08aSLois Curfman McInnes     PetscLogInfo(0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm);
278a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2836669109SBarry Smith   } else {
2936669109SBarry Smith     Vec       work;
30ea709b57SSatish Balay     PetscScalar    result;
3136669109SBarry Smith     PetscReal wnorm;
3236669109SBarry Smith 
3336669109SBarry Smith     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3736669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
3836669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
3936669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
404b27c08aSLois Curfman McInnes     PetscLogInfo(0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1);
4136669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4236669109SBarry Smith   }
438a5d9ee4SBarry Smith   PetscFunctionReturn(0);
448a5d9ee4SBarry Smith }
458a5d9ee4SBarry Smith 
4674637425SBarry Smith /*
475ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4874637425SBarry Smith */
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
5174637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5274637425SBarry Smith {
5374637425SBarry Smith   PetscReal     a1,a2;
5474637425SBarry Smith   int           ierr;
5574637425SBarry Smith   PetscTruth    hastranspose;
56ea709b57SSatish Balay   PetscScalar   mone = -1.0;
5774637425SBarry Smith 
5874637425SBarry Smith   PetscFunctionBegin;
5974637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
6074637425SBarry Smith   if (hastranspose) {
6174637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6274637425SBarry Smith     ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr);
6374637425SBarry Smith 
6474637425SBarry Smith     /* Compute || J^T W|| */
6574637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6674637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6774637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
689994e62eSSatish Balay     if (a1 != 0) {
694b27c08aSLois Curfman McInnes       PetscLogInfo(0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1);
7085471664SBarry Smith     }
7174637425SBarry Smith   }
7274637425SBarry Smith   PetscFunctionReturn(0);
7374637425SBarry Smith }
7474637425SBarry Smith 
7504d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7604d965bbSLois Curfman McInnes 
7704d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7894b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7904d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8004d965bbSLois Curfman McInnes      respectively.
8104d965bbSLois Curfman McInnes 
82fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8304d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8404d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
85fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8604d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8704d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
884b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
894b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
9004d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9104d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9204d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9304d965bbSLois Curfman McInnes      for all nonlinear solvers.
9404d965bbSLois Curfman McInnes 
9504d965bbSLois Curfman McInnes      Another key routine is:
9604d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9704d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9804d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9904d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10004d965bbSLois Curfman McInnes 
10104d965bbSLois Curfman McInnes      Additional basic routines are:
10204d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10304d965bbSLois Curfman McInnes                                       have actually been used.
10404d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
105186905e3SBarry Smith      SNESView().
10604d965bbSLois Curfman McInnes 
10704d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10804d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10904d965bbSLois Curfman McInnes      above description applies to these categories also.
11004d965bbSLois Curfman McInnes 
11104d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1125e42470aSBarry Smith /*
1134b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11404d965bbSLois Curfman McInnes    method with a line search.
1155e76c431SBarry Smith 
11604d965bbSLois Curfman McInnes    Input Parameters:
11704d965bbSLois Curfman McInnes .  snes - the SNES context
1185e76c431SBarry Smith 
11904d965bbSLois Curfman McInnes    Output Parameter:
120c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12104d965bbSLois Curfman McInnes 
12204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1235e76c431SBarry Smith 
1245e76c431SBarry Smith    Notes:
1255e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1265e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1275e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1285e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
129393d2d9aSLois Curfman McInnes    and Schnabel.
1305e76c431SBarry Smith */
1314a2ae208SSatish Balay #undef __FUNCT__
1324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
1334b27c08aSLois Curfman McInnes int SNESSolve_LS(SNES snes,int *outits)
1345e76c431SBarry Smith {
1354b27c08aSLois Curfman McInnes   SNES_LS      *neP = (SNES_LS*)snes->data;
13684c569c9SBarry Smith   int          maxits,i,ierr,lits,lsfail;
137112a2221SBarry Smith   MatStructure flg = DIFFERENT_NONZERO_PATTERN;
138329f5518SBarry Smith   PetscReal    fnorm,gnorm,xnorm,ynorm;
1395e42470aSBarry Smith   Vec          Y,X,F,G,W,TMP;
140c293cc10SBarry Smith   KSP          ksp;
1415e76c431SBarry Smith 
1423a40ed3dSBarry Smith   PetscFunctionBegin;
14394b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
144184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
145184914b5SBarry Smith 
1465e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1475e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14839e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1495e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1505e42470aSBarry Smith   G		= snes->work[1];
1515e42470aSBarry Smith   W		= snes->work[2];
1525e76c431SBarry Smith 
1534c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1544c49b128SBarry Smith   snes->iter = 0;
1554c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1565334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
157cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1580f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1595e42470aSBarry Smith   snes->norm = fnorm;
1600f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16184c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1635e76c431SBarry Smith 
164184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16594a9d846SBarry Smith 
166d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
167d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
168d034289bSBarry Smith 
1695e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1705e76c431SBarry Smith 
171329e7e40SMatthew Knepley     /* Call general purpose update function */
172de8cb200SMatthew Knepley     if (snes->update != PETSC_NULL) {
173329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
174de8cb200SMatthew Knepley     }
175329e7e40SMatthew Knepley 
176ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1775334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
17894b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
17994b7f48cSBarry Smith     ierr = KSPSetRhs(snes->ksp,F);CHKERRQ(ierr);
18094b7f48cSBarry Smith     ierr = KSPSetSolution(snes->ksp,Y);CHKERRQ(ierr);
18194b7f48cSBarry Smith     ierr = KSPSolve(snes->ksp);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;
1904b27c08aSLois Curfman McInnes     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);
197af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
1984b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
199a3b891d8SBarry Smith 
200a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
201a3b891d8SBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
202a3b891d8SBarry Smith     fnorm = gnorm;
203a3b891d8SBarry Smith 
2045ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2055ed2d596SBarry Smith     snes->iter = i+1;
2065ed2d596SBarry Smith     snes->norm = fnorm;
2075ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2085ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2095ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2105ed2d596SBarry Smith 
2113505fcc1SBarry Smith     if (lsfail) {
2128a5d9ee4SBarry Smith       PetscTruth ismin;
21350ffb88aSMatthew Knepley 
21450ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2153505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2168a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2178a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2183505fcc1SBarry Smith         break;
2193505fcc1SBarry Smith       }
22050ffb88aSMatthew Knepley     }
2215e76c431SBarry Smith 
2225e76c431SBarry Smith     /* Test for convergence */
22329e0b56fSBarry Smith     if (snes->converged) {
22429e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2253505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2263505fcc1SBarry Smith       if (snes->reason) {
22716c95cacSBarry Smith         break;
22816c95cacSBarry Smith       }
22916c95cacSBarry Smith     }
23029e0b56fSBarry Smith   }
23139e2f89bSBarry Smith   if (X != snes->vec_sol) {
232393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
233e15f7bb5SBarry Smith   }
234e15f7bb5SBarry Smith   if (F != snes->vec_func) {
235e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
236e15f7bb5SBarry Smith   }
23739e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
23839e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
23952392280SLois Curfman McInnes   if (i == maxits) {
2404b27c08aSLois Curfman McInnes     PetscLogInfo(snes,"SNESSolve_LS: Maximum number of iterations has been reached: %d\n",maxits);
24152392280SLois Curfman McInnes     i--;
2423505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
24352392280SLois Curfman McInnes   }
2440f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2450f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2465e42470aSBarry Smith   *outits = i+1;
2473a40ed3dSBarry Smith   PetscFunctionReturn(0);
2485e76c431SBarry Smith }
24904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25004d965bbSLois Curfman McInnes /*
2514b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2524b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25304d965bbSLois Curfman McInnes 
25404d965bbSLois Curfman McInnes    Input Parameter:
25504d965bbSLois Curfman McInnes .  snes - the SNES context
25604d965bbSLois Curfman McInnes .  x - the solution vector
25704d965bbSLois Curfman McInnes 
25804d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25904d965bbSLois Curfman McInnes 
26004d965bbSLois Curfman McInnes    Notes:
2614b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26204d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26304d965bbSLois Curfman McInnes    the call to SNESSolve().
26404d965bbSLois Curfman McInnes  */
2654a2ae208SSatish Balay #undef __FUNCT__
2664b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
2674b27c08aSLois Curfman McInnes int SNESSetUp_LS(SNES snes)
2685e76c431SBarry Smith {
2695e42470aSBarry Smith   int ierr;
2703a40ed3dSBarry Smith 
2713a40ed3dSBarry Smith   PetscFunctionBegin;
27281b6cf68SLois Curfman McInnes   snes->nwork = 4;
273d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
274b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
27581b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2763a40ed3dSBarry Smith   PetscFunctionReturn(0);
2775e76c431SBarry Smith }
27804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27904d965bbSLois Curfman McInnes /*
2804b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2814b27c08aSLois Curfman McInnes    with SNESCreate_LS().
28204d965bbSLois Curfman McInnes 
28304d965bbSLois Curfman McInnes    Input Parameter:
28404d965bbSLois Curfman McInnes .  snes - the SNES context
28504d965bbSLois Curfman McInnes 
286de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28704d965bbSLois Curfman McInnes  */
2884a2ae208SSatish Balay #undef __FUNCT__
2894b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
2904b27c08aSLois Curfman McInnes int SNESDestroy_LS(SNES snes)
2915e76c431SBarry Smith {
292393d2d9aSLois Curfman McInnes   int  ierr;
2933a40ed3dSBarry Smith 
2943a40ed3dSBarry Smith   PetscFunctionBegin;
2955baf8537SBarry Smith   if (snes->nwork) {
2964b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29721c89e3eSBarry Smith   }
2985c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
3005e76c431SBarry Smith }
30104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3024a2ae208SSatish Balay #undef __FUNCT__
3034a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
30482bf6240SBarry Smith 
3054b828684SBarry Smith /*@C
3065e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3075e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3085e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3095e76c431SBarry Smith 
310c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
311c7afd0dbSLois Curfman McInnes 
3125e76c431SBarry Smith    Input Parameters:
313c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
314af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3155e76c431SBarry Smith .  x - current iterate
3165e76c431SBarry Smith .  f - residual evaluated at x
3175e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3185e76c431SBarry Smith .  w - work vector
319c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3205e76c431SBarry Smith 
321c4a48953SLois Curfman McInnes    Output Parameters:
322c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3235e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3245e76c431SBarry Smith .  gnorm - 2-norm of g
3255e76c431SBarry Smith .  ynorm - 2-norm of search length
326c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
327fee21e36SBarry Smith 
328c4a48953SLois Curfman McInnes    Options Database Key:
3294b27c08aSLois Curfman McInnes .  -snes_ls basic - Activates SNESNoLineSearch()
330c4a48953SLois Curfman McInnes 
33136851e7fSLois Curfman McInnes    Level: advanced
33236851e7fSLois Curfman McInnes 
33328ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33428ae5a21SLois Curfman McInnes 
335f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
336af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3375e76c431SBarry Smith @*/
3385d1a10b1SBarry 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)
3395e76c431SBarry Smith {
3405e42470aSBarry Smith   int         ierr;
341ea709b57SSatish Balay   PetscScalar mone = -1.0;
3424b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
3438f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
3445334005bSBarry Smith 
3453a40ed3dSBarry Smith   PetscFunctionBegin;
346761aaf1bSLois Curfman McInnes   *flag = 0;
347d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
348a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
349ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3508f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3518f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3528f99978dSLois Curfman McInnes   }
353ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
354a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
355d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3563a40ed3dSBarry Smith   PetscFunctionReturn(0);
3575e76c431SBarry Smith }
35804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3594a2ae208SSatish Balay #undef __FUNCT__
3604a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
36182bf6240SBarry Smith 
36229e0b56fSBarry Smith /*@C
36329e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
36429e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
36529e0b56fSBarry Smith    even compute the norm of the function or search direction; this
36629e0b56fSBarry Smith    is intended only when you know the full step is fine and are
36729e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
368c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
369c7afd0dbSLois Curfman McInnes 
370c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
37129e0b56fSBarry Smith 
37229e0b56fSBarry Smith    Input Parameters:
373c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
374af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
37529e0b56fSBarry Smith .  x - current iterate
37629e0b56fSBarry Smith .  f - residual evaluated at x
37729e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
37829e0b56fSBarry Smith .  w - work vector
379c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
38029e0b56fSBarry Smith 
38129e0b56fSBarry Smith    Output Parameters:
382c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
38329e0b56fSBarry Smith .  gnorm - not changed
38429e0b56fSBarry Smith .  ynorm - not changed
385c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
386fee21e36SBarry Smith 
38729e0b56fSBarry Smith    Options Database Key:
3884b27c08aSLois Curfman McInnes .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
38929e0b56fSBarry Smith 
3908cbba510SBarry Smith    Notes:
391ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
392ea56f5baSLois Curfman McInnes    either the options
393ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
394ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
395329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
396329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
397329f5518SBarry Smith 
398329f5518SBarry Smith    During the final iteration this will not evaluate the function at
399329f5518SBarry Smith    the solution point. This is to save a function evaluation while
400329f5518SBarry Smith    using pseudo-timestepping.
4018cbba510SBarry Smith 
402ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
403ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
404ea56f5baSLois Curfman McInnes    correct, since they are not computed.
405ea56f5baSLois Curfman McInnes 
406ea56f5baSLois Curfman McInnes    Level: advanced
4078cbba510SBarry Smith 
40829e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
40929e0b56fSBarry Smith 
41029e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
41129e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
41229e0b56fSBarry Smith @*/
4135d1a10b1SBarry 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)
41429e0b56fSBarry Smith {
41529e0b56fSBarry Smith   int         ierr;
416ea709b57SSatish Balay   PetscScalar mone = -1.0;
4174b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
4188f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
41929e0b56fSBarry Smith 
4203a40ed3dSBarry Smith   PetscFunctionBegin;
4218cbba510SBarry Smith   *flag = 0;
422d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
42329e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4248f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4258f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4268f99978dSLois Curfman McInnes   }
427329f5518SBarry Smith 
428329f5518SBarry Smith   /* don't evaluate function the last time through */
429329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
43029e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
431329f5518SBarry Smith   }
432d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4333a40ed3dSBarry Smith   PetscFunctionReturn(0);
43429e0b56fSBarry Smith }
43504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4364a2ae208SSatish Balay #undef __FUNCT__
4374a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4384b828684SBarry Smith /*@C
439f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4405e76c431SBarry Smith 
441c7afd0dbSLois Curfman McInnes    Collective on SNES
442c7afd0dbSLois Curfman McInnes 
4435e76c431SBarry Smith    Input Parameters:
444c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
445af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4465e76c431SBarry Smith .  x - current iterate
4475e76c431SBarry Smith .  f - residual evaluated at x
4485e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4495e76c431SBarry Smith .  w - work vector
450c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4515e76c431SBarry Smith 
452393d2d9aSLois Curfman McInnes    Output Parameters:
453c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4545e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4555e76c431SBarry Smith .  gnorm - 2-norm of g
4565e76c431SBarry Smith .  ynorm - 2-norm of search length
457c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
458fee21e36SBarry Smith 
459c4a48953SLois Curfman McInnes    Options Database Key:
4604b27c08aSLois Curfman McInnes $  -snes_ls cubic - Activates SNESCubicLineSearch()
461c4a48953SLois Curfman McInnes 
4625e76c431SBarry Smith    Notes:
4635e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4645e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4655e76c431SBarry Smith 
46636851e7fSLois Curfman McInnes    Level: advanced
46736851e7fSLois Curfman McInnes 
46828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
46928ae5a21SLois Curfman McInnes 
470af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4715e76c431SBarry Smith @*/
4725d1a10b1SBarry 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)
4735e76c431SBarry Smith {
474406556e6SLois Curfman McInnes   /*
475406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
476406556e6SLois Curfman McInnes      minimization problem:
477406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
478406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
479406556e6SLois Curfman McInnes    */
480406556e6SLois Curfman McInnes 
4815ca10a37SBarry Smith   PetscReal   steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
482329f5518SBarry Smith   PetscReal   maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
483aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
484ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
4856b5873e3SBarry Smith #endif
4865e42470aSBarry Smith   int         ierr,count;
4874b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
488ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
4898f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
4905e76c431SBarry Smith 
4913a40ed3dSBarry Smith   PetscFunctionBegin;
492d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
493761aaf1bSLois Curfman McInnes   *flag   = 0;
4945e76c431SBarry Smith   alpha   = neP->alpha;
4955e76c431SBarry Smith   maxstep = neP->maxstep;
4965e76c431SBarry Smith   steptol = neP->steptol;
4975e76c431SBarry Smith 
498cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
49963b9fa5eSBarry Smith   if (*ynorm == 0.0) {
50063b9fa5eSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size is 0\n");
501a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
502a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
503a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
50463b9fa5eSBarry Smith     *flag  = -1;
505ad922ac9SBarry Smith     goto theend1;
50694a9d846SBarry Smith   }
5075e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5085e42470aSBarry Smith     scale = maxstep/(*ynorm);
509aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
510b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
5116b5873e3SBarry Smith #else
512b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
5136b5873e3SBarry Smith #endif
514393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5155e76c431SBarry Smith     *ynorm = maxstep;
5165e76c431SBarry Smith   }
517ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5185ca10a37SBarry Smith   minlambda = steptol/rellength;
519a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
520aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
521a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
522329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5235e42470aSBarry Smith #else
524a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5255e42470aSBarry Smith #endif
5265e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5275e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5285e76c431SBarry Smith 
529393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5305334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
53178b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
532313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5335d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
534393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
535b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
536ad922ac9SBarry Smith     goto theend1;
5375e76c431SBarry Smith   }
5385e76c431SBarry Smith 
5395e76c431SBarry Smith   /* Fit points with quadratic */
540313b5538SBarry Smith   lambda = 1.0;
541a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5425e76c431SBarry Smith   lambdaprev = lambda;
5435e76c431SBarry Smith   gnormprev = *gnorm;
544329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
545ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
546ddd12767SBarry Smith   else                         lambda = lambdatemp;
547393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
548ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
549aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
550ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5515e42470aSBarry Smith #else
552ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5535e42470aSBarry Smith #endif
55478b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
555cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5565ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
557393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
5585ed2d596SBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda);
559ad922ac9SBarry Smith     goto theend1;
5605e76c431SBarry Smith   }
5615e76c431SBarry Smith 
5625e76c431SBarry Smith   /* Fit points with cubic */
5635e76c431SBarry Smith   count = 1;
5648229bfc2SMatthew Knepley   while (count < 10000) {
5655e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
566b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
5675ed2d596SBarry 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);
568f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
569761aaf1bSLois Curfman McInnes       *flag = -1; break;
5705e76c431SBarry Smith     }
5715d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
5725d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
573ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5742b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5755e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5765e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5775e76c431SBarry Smith     if (a == 0.0) {
5785e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5795e76c431SBarry Smith     } else {
5805e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5815e76c431SBarry Smith     }
5825e76c431SBarry Smith     lambdaprev = lambda;
5835e76c431SBarry Smith     gnormprev  = *gnorm;
584329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
585329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
5865e76c431SBarry Smith     else                         lambda     = lambdatemp;
587393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
588ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
589aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
590ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
591393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5925e42470aSBarry Smith #else
593ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5945e42470aSBarry Smith #endif
59578b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
596cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5975ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
598393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
5995ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda);
6002715565aSLois Curfman McInnes       goto theend1;
6012b022350SLois Curfman McInnes     } else {
6025ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);
6035e76c431SBarry Smith     }
6045e76c431SBarry Smith     count++;
6055e76c431SBarry Smith   }
6068229bfc2SMatthew Knepley   if (count >= 10000) {
607cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6088229bfc2SMatthew Knepley   }
609ad922ac9SBarry Smith   theend1:
6108f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6118f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6128f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6138f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6148f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6158f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6168f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6178f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6188f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6198f99978dSLois Curfman McInnes     }
6208f99978dSLois Curfman McInnes   }
621d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6223a40ed3dSBarry Smith   PetscFunctionReturn(0);
6235e76c431SBarry Smith }
62404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6254a2ae208SSatish Balay #undef __FUNCT__
6264a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6274b828684SBarry Smith /*@C
628f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6295e76c431SBarry Smith 
630c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
631c7afd0dbSLois Curfman McInnes 
6325e42470aSBarry Smith    Input Parameters:
633c7afd0dbSLois Curfman McInnes +  snes - the SNES context
634af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6355e42470aSBarry Smith .  x - current iterate
6365e42470aSBarry Smith .  f - residual evaluated at x
6375e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6385e42470aSBarry Smith .  w - work vector
639c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6405e42470aSBarry Smith 
641c4a48953SLois Curfman McInnes    Output Parameters:
642c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6435e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6445e42470aSBarry Smith .  gnorm - 2-norm of g
6455e42470aSBarry Smith .  ynorm - 2-norm of search length
646c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
647fee21e36SBarry Smith 
648c4a48953SLois Curfman McInnes    Options Database Key:
6494b27c08aSLois Curfman McInnes .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
650c4a48953SLois Curfman McInnes 
6515e42470aSBarry Smith    Notes:
6524b27c08aSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNESLS method.
6535e42470aSBarry Smith 
65436851e7fSLois Curfman McInnes    Level: advanced
65536851e7fSLois Curfman McInnes 
656f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
657f59ffdedSLois Curfman McInnes 
658af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6595e42470aSBarry Smith @*/
6605d1a10b1SBarry 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)
6615e76c431SBarry Smith {
662406556e6SLois Curfman McInnes   /*
663406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
664406556e6SLois Curfman McInnes      minimization problem:
665406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
666406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
667406556e6SLois Curfman McInnes    */
6685ca10a37SBarry Smith   PetscReal   steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
669aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
670ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
6716b5873e3SBarry Smith #endif
6725e42470aSBarry Smith   int         ierr,count;
6734b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
674ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
6758f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
6765e76c431SBarry Smith 
6773a40ed3dSBarry Smith   PetscFunctionBegin;
678d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
679761aaf1bSLois Curfman McInnes   *flag   = 0;
6805e42470aSBarry Smith   alpha   = neP->alpha;
6815e42470aSBarry Smith   maxstep = neP->maxstep;
6825e42470aSBarry Smith   steptol = neP->steptol;
6835e76c431SBarry Smith 
6843505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
68563b9fa5eSBarry Smith   if (*ynorm == 0.0) {
686b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
687b37302c6SLois Curfman McInnes     *gnorm = fnorm;
688b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
689b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
69063b9fa5eSBarry Smith     *flag  = -1;
691ad922ac9SBarry Smith     goto theend2;
69294a9d846SBarry Smith   }
6935e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6945e42470aSBarry Smith     scale = maxstep/(*ynorm);
695393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6965e42470aSBarry Smith     *ynorm = maxstep;
6975e76c431SBarry Smith   }
698ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
6995ca10a37SBarry Smith   minlambda = steptol/rellength;
700a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
701aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
702a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
703329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7045e42470aSBarry Smith #else
705a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7065e42470aSBarry Smith #endif
7075e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
7085e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7095e42470aSBarry Smith 
710393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7115334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
71278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
713cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7145d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
715393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
716b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
717ad922ac9SBarry Smith     goto theend2;
7185e42470aSBarry Smith   }
7195e42470aSBarry Smith 
7205e42470aSBarry Smith   /* Fit points with quadratic */
721313b5538SBarry Smith   lambda = 1.0;
7225e42470aSBarry Smith   count = 1;
7235ca10a37SBarry Smith   while (PETSC_TRUE) {
7245e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
725b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
726b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
727f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
728761aaf1bSLois Curfman McInnes       *flag = -1; break;
7295e42470aSBarry Smith     }
730a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
731329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
732329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
733329f5518SBarry Smith     else                         lambda     = lambdatemp;
734393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7353505fcc1SBarry Smith     lambdaneg = -lambda;
736aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7373505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7385e42470aSBarry Smith #else
7393505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7405e42470aSBarry Smith #endif
74178b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
742cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7435ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
744393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
745b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
74606259719SBarry Smith       break;
7475e42470aSBarry Smith     }
7485e42470aSBarry Smith     count++;
7495e42470aSBarry Smith   }
750ad922ac9SBarry Smith   theend2:
7518f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7528f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7538f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7548f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7558f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7568f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7578f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7588f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7598f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7608f99978dSLois Curfman McInnes     }
7618f99978dSLois Curfman McInnes   }
762d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7633a40ed3dSBarry Smith   PetscFunctionReturn(0);
7645e76c431SBarry Smith }
76504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7664a2ae208SSatish Balay #undef __FUNCT__
7674a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
768c9e6a524SLois Curfman McInnes /*@C
76977c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7704b27c08aSLois Curfman McInnes    by the method SNESLS.
7715e76c431SBarry Smith 
7725e76c431SBarry Smith    Input Parameters:
773c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
774af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
775c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7765e76c431SBarry Smith 
777fee21e36SBarry Smith    Collective on SNES
778fee21e36SBarry Smith 
779c4a48953SLois Curfman McInnes    Available Routines:
780c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
781f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
782f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
783af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7845e76c431SBarry Smith 
785c4a48953SLois Curfman McInnes     Options Database Keys:
7864b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
7874b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
7884b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
7894b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
7903304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
7913304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
792c4a48953SLois Curfman McInnes 
7935e76c431SBarry Smith    Calling sequence of func:
794bd208895SLois Curfman McInnes .vb
795af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
796329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
797bd208895SLois Curfman McInnes .ve
7985e76c431SBarry Smith 
7995e76c431SBarry Smith     Input parameters for func:
800c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
801af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8025e76c431SBarry Smith .   x - current iterate
8035e76c431SBarry Smith .   f - residual evaluated at x
8045e76c431SBarry Smith .   y - search direction (contains new iterate on output)
8055e76c431SBarry Smith .   w - work vector
806c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8075e76c431SBarry Smith 
8085e76c431SBarry Smith     Output parameters for func:
809c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8105e76c431SBarry Smith .   y - new iterate (contains search direction on input)
8115e76c431SBarry Smith .   gnorm - 2-norm of g
8125e76c431SBarry Smith .   ynorm - 2-norm of search length
813c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
814761aaf1bSLois Curfman McInnes            on failure.
815f59ffdedSLois Curfman McInnes 
81636851e7fSLois Curfman McInnes     Level: advanced
81736851e7fSLois Curfman McInnes 
818f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
819f59ffdedSLois Curfman McInnes 
8205d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8215d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8225e76c431SBarry Smith @*/
823329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
8245e76c431SBarry Smith {
825329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
82682bf6240SBarry Smith 
8273a40ed3dSBarry Smith   PetscFunctionBegin;
828c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
82982bf6240SBarry Smith   if (f) {
830af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
83182bf6240SBarry Smith   }
8323a40ed3dSBarry Smith   PetscFunctionReturn(0);
8335e76c431SBarry Smith }
83404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
835fb2e594dSBarry Smith EXTERN_C_BEGIN
8364a2ae208SSatish Balay #undef __FUNCT__
8374a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
838af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
83987828ca2SBarry Smith                          PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
84082bf6240SBarry Smith {
84182bf6240SBarry Smith   PetscFunctionBegin;
8424b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
8434b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
84482bf6240SBarry Smith   PetscFunctionReturn(0);
84582bf6240SBarry Smith }
846fb2e594dSBarry Smith EXTERN_C_END
84704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8484a2ae208SSatish Balay #undef __FUNCT__
8494a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
850c8dd0c0dSLois Curfman McInnes /*@C
851530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8524b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
853c8dd0c0dSLois Curfman McInnes 
854c8dd0c0dSLois Curfman McInnes    Input Parameters:
855c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
856c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
857c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
858c8dd0c0dSLois Curfman McInnes 
859c8dd0c0dSLois Curfman McInnes    Collective on SNES
860c8dd0c0dSLois Curfman McInnes 
861c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
862c8dd0c0dSLois Curfman McInnes .vb
863b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
864c8dd0c0dSLois Curfman McInnes .ve
865b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
866b1ae27eaSLois Curfman McInnes    on failure.
867c8dd0c0dSLois Curfman McInnes 
868c8dd0c0dSLois Curfman McInnes    Input parameters for func:
869c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
870c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8718f99978dSLois Curfman McInnes -  x - current candidate iterate
872c8dd0c0dSLois Curfman McInnes 
873c8dd0c0dSLois Curfman McInnes    Output parameters for func:
874c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
875c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
876c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
877c8dd0c0dSLois Curfman McInnes 
878c8dd0c0dSLois Curfman McInnes    Level: advanced
879c8dd0c0dSLois Curfman McInnes 
8808f99978dSLois Curfman McInnes    Notes:
881b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
882b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
883b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
884b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
885ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8868f99978dSLois Curfman McInnes 
887b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
888b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
889b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
890ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
891ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
892ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
893ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
894b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8958f99978dSLois Curfman McInnes 
896c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
897c8dd0c0dSLois Curfman McInnes 
898c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
899c8dd0c0dSLois Curfman McInnes @*/
9008f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
901c8dd0c0dSLois Curfman McInnes {
9028f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
903c8dd0c0dSLois Curfman McInnes 
904c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
905c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
906c8dd0c0dSLois Curfman McInnes   if (f) {
907c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
908c8dd0c0dSLois Curfman McInnes   }
909c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
910c8dd0c0dSLois Curfman McInnes }
911c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
912c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
9134a2ae208SSatish Balay #undef __FUNCT__
9144a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
9158f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
916c8dd0c0dSLois Curfman McInnes {
917c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9184b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
9194b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP    = checkctx;
920c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
921c8dd0c0dSLois Curfman McInnes }
922c8dd0c0dSLois Curfman McInnes EXTERN_C_END
923c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
92404d965bbSLois Curfman McInnes /*
9254b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
926329e7e40SMatthew Knepley 
927329e7e40SMatthew Knepley    Input Parameter:
928329e7e40SMatthew Knepley .  snes - the SNES context
929329e7e40SMatthew Knepley 
930329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
931329e7e40SMatthew Knepley */
932329e7e40SMatthew Knepley #undef __FUNCT__
9334b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
9344b27c08aSLois Curfman McInnes static int SNESPrintHelp_LS(SNES snes,char *p)
935329e7e40SMatthew Knepley {
9364b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
937329e7e40SMatthew Knepley 
938329e7e40SMatthew Knepley   PetscFunctionBegin;
9394b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
9404b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
9414b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
9424b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
9434b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
944329e7e40SMatthew Knepley   PetscFunctionReturn(0);
945329e7e40SMatthew Knepley }
946329e7e40SMatthew Knepley 
947329e7e40SMatthew Knepley /*
9484b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
94904d965bbSLois Curfman McInnes 
95004d965bbSLois Curfman McInnes    Input Parameters:
95104d965bbSLois Curfman McInnes .  SNES - the SNES context
95204d965bbSLois Curfman McInnes .  viewer - visualization context
95304d965bbSLois Curfman McInnes 
95404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
95504d965bbSLois Curfman McInnes */
9564a2ae208SSatish Balay #undef __FUNCT__
9574b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
9584b27c08aSLois Curfman McInnes static int SNESView_LS(SNES snes,PetscViewer viewer)
959a935fc98SLois Curfman McInnes {
9604b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
961*2fc52814SBarry Smith   const char *cstr;
96251695f54SLois Curfman McInnes   int        ierr;
9636831982aSBarry Smith   PetscTruth isascii;
964a935fc98SLois Curfman McInnes 
9653a40ed3dSBarry Smith   PetscFunctionBegin;
966b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
9670f5bd95cSBarry Smith   if (isascii) {
96819bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
96919bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
97019bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
97119bcc07fSBarry Smith     else                                                cstr = "unknown";
972b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
973b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9745cd90555SBarry Smith   } else {
97529bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
97619bcc07fSBarry Smith   }
9773a40ed3dSBarry Smith   PetscFunctionReturn(0);
978a935fc98SLois Curfman McInnes }
97904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
98004d965bbSLois Curfman McInnes /*
9814b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
98204d965bbSLois Curfman McInnes 
98304d965bbSLois Curfman McInnes    Input Parameter:
98404d965bbSLois Curfman McInnes .  snes - the SNES context
98504d965bbSLois Curfman McInnes 
98604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
98704d965bbSLois Curfman McInnes */
9884a2ae208SSatish Balay #undef __FUNCT__
9894b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
9904b27c08aSLois Curfman McInnes static int SNESSetFromOptions_LS(SNES snes)
9915e42470aSBarry Smith {
9924b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
993e5999256SBarry Smith   const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
99422e36657SBarry Smith   int        ierr,indx;
995f1af5d2fSBarry Smith   PetscTruth flg;
9965e42470aSBarry Smith 
9973a40ed3dSBarry Smith   PetscFunctionBegin;
998b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
9994b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
10004b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
10014b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1002186905e3SBarry Smith 
100322e36657SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
100425cdf11fSBarry Smith     if (flg) {
100522e36657SBarry Smith       switch (indx) {
1006b49fd9e1SBarry Smith       case 0:
1007af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1008b49fd9e1SBarry Smith         break;
1009b49fd9e1SBarry Smith       case 1:
1010af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1011b49fd9e1SBarry Smith         break;
1012b49fd9e1SBarry Smith       case 2:
1013af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1014b49fd9e1SBarry Smith         break;
1015b49fd9e1SBarry Smith       case 3:
1016af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1017b49fd9e1SBarry Smith         break;
10185e42470aSBarry Smith       }
10195e42470aSBarry Smith     }
1020b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10213a40ed3dSBarry Smith   PetscFunctionReturn(0);
10225e42470aSBarry Smith }
102304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
102404d965bbSLois Curfman McInnes /*
10254b27c08aSLois Curfman McInnes    SNESCreate_LS - Creates a nonlinear solver context for the SNESLS method,
10264b27c08aSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
102704d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
102804d965bbSLois Curfman McInnes 
102904d965bbSLois Curfman McInnes 
103004d965bbSLois Curfman McInnes    Input Parameter:
103104d965bbSLois Curfman McInnes .  snes - the SNES context
103204d965bbSLois Curfman McInnes 
103304d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
103404d965bbSLois Curfman McInnes  */
1035fb2e594dSBarry Smith EXTERN_C_BEGIN
10364a2ae208SSatish Balay #undef __FUNCT__
10374b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
10384b27c08aSLois Curfman McInnes int SNESCreate_LS(SNES snes)
10395e42470aSBarry Smith {
104082bf6240SBarry Smith   int     ierr;
10414b27c08aSLois Curfman McInnes   SNES_LS *neP;
10425e42470aSBarry Smith 
10433a40ed3dSBarry Smith   PetscFunctionBegin;
10444b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
10454b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
10464b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
10474b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
10484b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
10494b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
10504b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
10515baf8537SBarry Smith   snes->nwork           = 0;
10525e42470aSBarry Smith 
10534b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
10544b27c08aSLois Curfman McInnes   PetscLogObjectMemory(snes,sizeof(SNES_LS));
10555e42470aSBarry Smith   snes->data    	= (void*)neP;
10565e42470aSBarry Smith   neP->alpha		= 1.e-4;
10575e42470aSBarry Smith   neP->maxstep		= 1.e8;
10585e42470aSBarry Smith   neP->steptol		= 1.e-12;
10595e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1060c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1061c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1062c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
106382bf6240SBarry Smith 
10645d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
10655d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
106682bf6240SBarry Smith 
10673a40ed3dSBarry Smith   PetscFunctionReturn(0);
10685e42470aSBarry Smith }
1069fb2e594dSBarry Smith EXTERN_C_END
107082bf6240SBarry Smith 
107182bf6240SBarry Smith 
107282bf6240SBarry Smith 
1073