xref: /petsc/src/snes/impls/ls/ls.c (revision 50ffb88a32740ec1e28141e8b16795f505e93f38)
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);
26b0a32e0cSBarry Smith     PetscLogInfo(0,"SNESSolve_EQ_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);
40b0a32e0cSBarry Smith     PetscLogInfo(0,"SNESSolve_EQ_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 /*
4774637425SBarry Smith      Checks if J^T(F - AX) = 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) {
69b0a32e0cSBarry Smith       PetscLogInfo(0,"SNESSolve_EQ_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,
7804d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, 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
8804d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
8904d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) 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 /*
11304d965bbSLois Curfman McInnes    SNESSolve_EQ_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__
1324a2ae208SSatish Balay #define __FUNCT__ "SNESSolve_EQ_LS"
133f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
1345e76c431SBarry Smith {
1356831982aSBarry Smith   SNES_EQ_LS          *neP = (SNES_EQ_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;
1405e76c431SBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
142184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
143184914b5SBarry Smith 
1445e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1455e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14639e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1475e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1485e42470aSBarry Smith   G		= snes->work[1];
1495e42470aSBarry Smith   W		= snes->work[2];
1505e76c431SBarry Smith 
1514c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1524c49b128SBarry Smith   snes->iter = 0;
1534c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1545334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
155cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1560f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1575e42470aSBarry Smith   snes->norm = fnorm;
1580f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15984c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16094a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1615e76c431SBarry Smith 
162184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16394a9d846SBarry Smith 
164d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
165d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
166d034289bSBarry Smith 
1675e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1685e76c431SBarry Smith 
169329e7e40SMatthew Knepley     /* Call general purpose update function */
170de8cb200SMatthew Knepley     if (snes->update != PETSC_NULL) {
171329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
172de8cb200SMatthew Knepley     }
173329e7e40SMatthew Knepley 
174ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1755334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
1765334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
17778b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
17874637425SBarry Smith 
179b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
18074637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
18185471664SBarry Smith     }
18274637425SBarry Smith 
18390cbd584SBarry Smith     /* should check what happened to the linear solve? */
1843505fcc1SBarry Smith     snes->linear_its += lits;
185b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
186ea4d3ed3SLois Curfman McInnes 
187ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
188ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
189ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
190ea4d3ed3SLois Curfman McInnes     */
19181b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
192af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
193b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
194a3b891d8SBarry Smith 
195a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
196a3b891d8SBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
197a3b891d8SBarry Smith     fnorm = gnorm;
198a3b891d8SBarry Smith 
1993505fcc1SBarry Smith     if (lsfail) {
2008a5d9ee4SBarry Smith       PetscTruth ismin;
201*50ffb88aSMatthew Knepley 
202*50ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2033505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2048a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2058a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2063505fcc1SBarry Smith         break;
2073505fcc1SBarry Smith       }
208*50ffb88aSMatthew Knepley     }
2095e76c431SBarry Smith 
2100f5bd95cSBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
21184c569c9SBarry Smith     snes->iter = i+1;
2125e42470aSBarry Smith     snes->norm = fnorm;
2130f5bd95cSBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
21484c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
21594a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2165e76c431SBarry Smith 
2175e76c431SBarry Smith     /* Test for convergence */
21829e0b56fSBarry Smith     if (snes->converged) {
21929e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2203505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2213505fcc1SBarry Smith       if (snes->reason) {
22216c95cacSBarry Smith         break;
22316c95cacSBarry Smith       }
22416c95cacSBarry Smith     }
22529e0b56fSBarry Smith   }
22639e2f89bSBarry Smith   if (X != snes->vec_sol) {
227393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
228e15f7bb5SBarry Smith   }
229e15f7bb5SBarry Smith   if (F != snes->vec_func) {
230e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
231e15f7bb5SBarry Smith   }
23239e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
23339e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
23452392280SLois Curfman McInnes   if (i == maxits) {
235b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
23652392280SLois Curfman McInnes     i--;
2373505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
23852392280SLois Curfman McInnes   }
2390f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2400f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2415e42470aSBarry Smith   *outits = i+1;
2423a40ed3dSBarry Smith   PetscFunctionReturn(0);
2435e76c431SBarry Smith }
24404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
24504d965bbSLois Curfman McInnes /*
24604d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
2476831982aSBarry Smith    of the SNESEQLS nonlinear solver.
24804d965bbSLois Curfman McInnes 
24904d965bbSLois Curfman McInnes    Input Parameter:
25004d965bbSLois Curfman McInnes .  snes - the SNES context
25104d965bbSLois Curfman McInnes .  x - the solution vector
25204d965bbSLois Curfman McInnes 
25304d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25404d965bbSLois Curfman McInnes 
25504d965bbSLois Curfman McInnes    Notes:
25604d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
25704d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
25804d965bbSLois Curfman McInnes    the call to SNESSolve().
25904d965bbSLois Curfman McInnes  */
2604a2ae208SSatish Balay #undef __FUNCT__
2614a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp_EQ_LS"
262f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
2635e76c431SBarry Smith {
2645e42470aSBarry Smith   int ierr;
2653a40ed3dSBarry Smith 
2663a40ed3dSBarry Smith   PetscFunctionBegin;
26781b6cf68SLois Curfman McInnes   snes->nwork = 4;
268d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
269b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
27081b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2713a40ed3dSBarry Smith   PetscFunctionReturn(0);
2725e76c431SBarry Smith }
27304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27404d965bbSLois Curfman McInnes /*
2756831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
27604d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
27704d965bbSLois Curfman McInnes 
27804d965bbSLois Curfman McInnes    Input Parameter:
27904d965bbSLois Curfman McInnes .  snes - the SNES context
28004d965bbSLois Curfman McInnes 
281de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28204d965bbSLois Curfman McInnes  */
2834a2ae208SSatish Balay #undef __FUNCT__
2844a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy_EQ_LS"
285e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
2865e76c431SBarry Smith {
287393d2d9aSLois Curfman McInnes   int  ierr;
2883a40ed3dSBarry Smith 
2893a40ed3dSBarry Smith   PetscFunctionBegin;
2905baf8537SBarry Smith   if (snes->nwork) {
2914b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29221c89e3eSBarry Smith   }
2935c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2943a40ed3dSBarry Smith   PetscFunctionReturn(0);
2955e76c431SBarry Smith }
29604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2974a2ae208SSatish Balay #undef __FUNCT__
2984a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
29982bf6240SBarry Smith 
3004b828684SBarry Smith /*@C
3015e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3025e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3035e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3045e76c431SBarry Smith 
305c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
306c7afd0dbSLois Curfman McInnes 
3075e76c431SBarry Smith    Input Parameters:
308c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
309af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3105e76c431SBarry Smith .  x - current iterate
3115e76c431SBarry Smith .  f - residual evaluated at x
3125e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3135e76c431SBarry Smith .  w - work vector
314c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3155e76c431SBarry Smith 
316c4a48953SLois Curfman McInnes    Output Parameters:
317c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3185e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3195e76c431SBarry Smith .  gnorm - 2-norm of g
3205e76c431SBarry Smith .  ynorm - 2-norm of search length
321c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
322fee21e36SBarry Smith 
323c4a48953SLois Curfman McInnes    Options Database Key:
324c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
325c4a48953SLois Curfman McInnes 
32636851e7fSLois Curfman McInnes    Level: advanced
32736851e7fSLois Curfman McInnes 
32828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
32928ae5a21SLois Curfman McInnes 
330f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
331af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3325e76c431SBarry Smith @*/
3335d1a10b1SBarry 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)
3345e76c431SBarry Smith {
3355e42470aSBarry Smith   int           ierr;
336ea709b57SSatish Balay   PetscScalar   mone = -1.0;
3376831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
3388f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
3395334005bSBarry Smith 
3403a40ed3dSBarry Smith   PetscFunctionBegin;
341761aaf1bSLois Curfman McInnes   *flag = 0;
342d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
343a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
344ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3458f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3468f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3478f99978dSLois Curfman McInnes   }
348ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
349a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
350d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3513a40ed3dSBarry Smith   PetscFunctionReturn(0);
3525e76c431SBarry Smith }
35304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3544a2ae208SSatish Balay #undef __FUNCT__
3554a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
35682bf6240SBarry Smith 
35729e0b56fSBarry Smith /*@C
35829e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
35929e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
36029e0b56fSBarry Smith    even compute the norm of the function or search direction; this
36129e0b56fSBarry Smith    is intended only when you know the full step is fine and are
36229e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
363c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
364c7afd0dbSLois Curfman McInnes 
365c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
36629e0b56fSBarry Smith 
36729e0b56fSBarry Smith    Input Parameters:
368c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
369af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
37029e0b56fSBarry Smith .  x - current iterate
37129e0b56fSBarry Smith .  f - residual evaluated at x
37229e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
37329e0b56fSBarry Smith .  w - work vector
374c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
37529e0b56fSBarry Smith 
37629e0b56fSBarry Smith    Output Parameters:
377c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
37829e0b56fSBarry Smith .  gnorm - not changed
37929e0b56fSBarry Smith .  ynorm - not changed
380c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
381fee21e36SBarry Smith 
38229e0b56fSBarry Smith    Options Database Key:
383c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
38429e0b56fSBarry Smith 
3858cbba510SBarry Smith    Notes:
386ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
387ea56f5baSLois Curfman McInnes    either the options
388ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
389ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
390329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
391329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
392329f5518SBarry Smith 
393329f5518SBarry Smith    During the final iteration this will not evaluate the function at
394329f5518SBarry Smith    the solution point. This is to save a function evaluation while
395329f5518SBarry Smith    using pseudo-timestepping.
3968cbba510SBarry Smith 
397ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
398ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
399ea56f5baSLois Curfman McInnes    correct, since they are not computed.
400ea56f5baSLois Curfman McInnes 
401ea56f5baSLois Curfman McInnes    Level: advanced
4028cbba510SBarry Smith 
40329e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
40429e0b56fSBarry Smith 
40529e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
40629e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
40729e0b56fSBarry Smith @*/
4085d1a10b1SBarry 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)
40929e0b56fSBarry Smith {
41029e0b56fSBarry Smith   int           ierr;
411ea709b57SSatish Balay   PetscScalar   mone = -1.0;
4126831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
4138f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
41429e0b56fSBarry Smith 
4153a40ed3dSBarry Smith   PetscFunctionBegin;
4168cbba510SBarry Smith   *flag = 0;
417d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
41829e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4198f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4208f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4218f99978dSLois Curfman McInnes   }
422329f5518SBarry Smith 
423329f5518SBarry Smith   /* don't evaluate function the last time through */
424329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
42529e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
426329f5518SBarry Smith   }
427d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4283a40ed3dSBarry Smith   PetscFunctionReturn(0);
42929e0b56fSBarry Smith }
43004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4314a2ae208SSatish Balay #undef __FUNCT__
4324a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4334b828684SBarry Smith /*@C
434f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4355e76c431SBarry Smith 
436c7afd0dbSLois Curfman McInnes    Collective on SNES
437c7afd0dbSLois Curfman McInnes 
4385e76c431SBarry Smith    Input Parameters:
439c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
440af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4415e76c431SBarry Smith .  x - current iterate
4425e76c431SBarry Smith .  f - residual evaluated at x
4435e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4445e76c431SBarry Smith .  w - work vector
445c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4465e76c431SBarry Smith 
447393d2d9aSLois Curfman McInnes    Output Parameters:
448c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4495e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4505e76c431SBarry Smith .  gnorm - 2-norm of g
4515e76c431SBarry Smith .  ynorm - 2-norm of search length
452c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
453fee21e36SBarry Smith 
454c4a48953SLois Curfman McInnes    Options Database Key:
455c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
456c4a48953SLois Curfman McInnes 
4575e76c431SBarry Smith    Notes:
4585e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4595e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4605e76c431SBarry Smith 
46136851e7fSLois Curfman McInnes    Level: advanced
46236851e7fSLois Curfman McInnes 
46328ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
46428ae5a21SLois Curfman McInnes 
465af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4665e76c431SBarry Smith @*/
4675d1a10b1SBarry 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)
4685e76c431SBarry Smith {
469406556e6SLois Curfman McInnes   /*
470406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
471406556e6SLois Curfman McInnes      minimization problem:
472406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
473406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
474406556e6SLois Curfman McInnes    */
475406556e6SLois Curfman McInnes 
476329f5518SBarry Smith   PetscReal     steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
477329f5518SBarry Smith   PetscReal     maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
478aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
479ea709b57SSatish Balay   PetscScalar   cinitslope,clambda;
4806b5873e3SBarry Smith #endif
4815e42470aSBarry Smith   int           ierr,count;
4826831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
483ea709b57SSatish Balay   PetscScalar   mone = -1.0,scale;
4848f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
4855e76c431SBarry Smith 
4863a40ed3dSBarry Smith   PetscFunctionBegin;
487d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
488761aaf1bSLois Curfman McInnes   *flag   = 0;
4895e76c431SBarry Smith   alpha   = neP->alpha;
4905e76c431SBarry Smith   maxstep = neP->maxstep;
4915e76c431SBarry Smith   steptol = neP->steptol;
4925e76c431SBarry Smith 
493cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
494a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
495b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
496a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
497a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
498a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
499ad922ac9SBarry Smith     goto theend1;
50094a9d846SBarry Smith   }
5015e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5025e42470aSBarry Smith     scale = maxstep/(*ynorm);
503aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
504b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
5056b5873e3SBarry Smith #else
506b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
5076b5873e3SBarry Smith #endif
508393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5095e76c431SBarry Smith     *ynorm = maxstep;
5105e76c431SBarry Smith   }
5115e76c431SBarry Smith   minlambda = steptol/(*ynorm);
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);
5495d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
550393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
551b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
552ad922ac9SBarry Smith     goto theend1;
5535e76c431SBarry Smith   }
5545e76c431SBarry Smith 
5555e76c431SBarry Smith   /* Fit points with cubic */
5565e76c431SBarry Smith   count = 1;
5575e76c431SBarry Smith   while (1) {
5585e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
559b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
560b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
561393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,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);
5905d1a10b1SBarry Smith     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
591393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
592b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
5932715565aSLois Curfman McInnes       goto theend1;
5942b022350SLois Curfman McInnes     } else {
595b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
5965e76c431SBarry Smith     }
5975e76c431SBarry Smith     count++;
5985e76c431SBarry Smith   }
599ad922ac9SBarry Smith   theend1:
6008f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6018f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6028f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6038f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6048f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6058f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6068f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6078f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6088f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6098f99978dSLois Curfman McInnes     }
6108f99978dSLois Curfman McInnes   }
611d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6123a40ed3dSBarry Smith   PetscFunctionReturn(0);
6135e76c431SBarry Smith }
61404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6154a2ae208SSatish Balay #undef __FUNCT__
6164a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6174b828684SBarry Smith /*@C
618f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6195e76c431SBarry Smith 
620c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
621c7afd0dbSLois Curfman McInnes 
6225e42470aSBarry Smith    Input Parameters:
623c7afd0dbSLois Curfman McInnes +  snes - the SNES context
624af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6255e42470aSBarry Smith .  x - current iterate
6265e42470aSBarry Smith .  f - residual evaluated at x
6275e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6285e42470aSBarry Smith .  w - work vector
629c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6305e42470aSBarry Smith 
631c4a48953SLois Curfman McInnes    Output Parameters:
632c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6335e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6345e42470aSBarry Smith .  gnorm - 2-norm of g
6355e42470aSBarry Smith .  ynorm - 2-norm of search length
636c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
637fee21e36SBarry Smith 
638c4a48953SLois Curfman McInnes    Options Database Key:
639c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
640c4a48953SLois Curfman McInnes 
6415e42470aSBarry Smith    Notes:
6426831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
6435e42470aSBarry Smith 
64436851e7fSLois Curfman McInnes    Level: advanced
64536851e7fSLois Curfman McInnes 
646f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
647f59ffdedSLois Curfman McInnes 
648af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6495e42470aSBarry Smith @*/
6505d1a10b1SBarry 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)
6515e76c431SBarry Smith {
652406556e6SLois Curfman McInnes   /*
653406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
654406556e6SLois Curfman McInnes      minimization problem:
655406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
656406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
657406556e6SLois Curfman McInnes    */
658329f5518SBarry Smith   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
659aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
660ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
6616b5873e3SBarry Smith #endif
6625e42470aSBarry Smith   int        ierr,count;
6636831982aSBarry Smith   SNES_EQ_LS     *neP = (SNES_EQ_LS*)snes->data;
664ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
6658f99978dSLois Curfman McInnes   PetscTruth     change_y = PETSC_FALSE;
6665e76c431SBarry Smith 
6673a40ed3dSBarry Smith   PetscFunctionBegin;
668d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
669761aaf1bSLois Curfman McInnes   *flag   = 0;
6705e42470aSBarry Smith   alpha   = neP->alpha;
6715e42470aSBarry Smith   maxstep = neP->maxstep;
6725e42470aSBarry Smith   steptol = neP->steptol;
6735e76c431SBarry Smith 
6743505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
675b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
676b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
677b37302c6SLois Curfman McInnes     *gnorm = fnorm;
678b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
679b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
680ad922ac9SBarry Smith     goto theend2;
68194a9d846SBarry Smith   }
6825e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6835e42470aSBarry Smith     scale = maxstep/(*ynorm);
684393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6855e42470aSBarry Smith     *ynorm = maxstep;
6865e76c431SBarry Smith   }
6875e42470aSBarry Smith   minlambda = steptol/(*ynorm);
688a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
689aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
690a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
691329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
6925e42470aSBarry Smith #else
693a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6945e42470aSBarry Smith #endif
6955e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
6965e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6975e42470aSBarry Smith 
698393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6995334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
70078b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
701cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7025d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
703393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
704b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
705ad922ac9SBarry Smith     goto theend2;
7065e42470aSBarry Smith   }
7075e42470aSBarry Smith 
7085e42470aSBarry Smith   /* Fit points with quadratic */
709313b5538SBarry Smith   lambda = 1.0;
7105e42470aSBarry Smith   count = 1;
7115e42470aSBarry Smith   while (1) {
7125e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
713b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
714b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
715393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
716761aaf1bSLois Curfman McInnes       *flag = -1; break;
7175e42470aSBarry Smith     }
718a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
719329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
720329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
721329f5518SBarry Smith     else                         lambda     = lambdatemp;
722393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7233505fcc1SBarry Smith     lambdaneg = -lambda;
724aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7253505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7265e42470aSBarry Smith #else
7273505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7285e42470aSBarry Smith #endif
72978b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
730cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7315d1a10b1SBarry Smith     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
732393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
733b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
73406259719SBarry Smith       break;
7355e42470aSBarry Smith     }
7365e42470aSBarry Smith     count++;
7375e42470aSBarry Smith   }
738ad922ac9SBarry Smith   theend2:
7398f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7408f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7418f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7428f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7438f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7448f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7458f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7468f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7478f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7488f99978dSLois Curfman McInnes     }
7498f99978dSLois Curfman McInnes   }
750d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7513a40ed3dSBarry Smith   PetscFunctionReturn(0);
7525e76c431SBarry Smith }
75304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7544a2ae208SSatish Balay #undef __FUNCT__
7554a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
756c9e6a524SLois Curfman McInnes /*@C
75777c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7586831982aSBarry Smith    by the method SNESEQLS.
7595e76c431SBarry Smith 
7605e76c431SBarry Smith    Input Parameters:
761c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
762af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
763c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7645e76c431SBarry Smith 
765fee21e36SBarry Smith    Collective on SNES
766fee21e36SBarry Smith 
767c4a48953SLois Curfman McInnes    Available Routines:
768c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
769f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
770f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
771af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7725e76c431SBarry Smith 
773c4a48953SLois Curfman McInnes     Options Database Keys:
774af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
775c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
776c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
777c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
778c4a48953SLois Curfman McInnes 
7795e76c431SBarry Smith    Calling sequence of func:
780bd208895SLois Curfman McInnes .vb
781af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
782329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
783bd208895SLois Curfman McInnes .ve
7845e76c431SBarry Smith 
7855e76c431SBarry Smith     Input parameters for func:
786c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
787af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
7885e76c431SBarry Smith .   x - current iterate
7895e76c431SBarry Smith .   f - residual evaluated at x
7905e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7915e76c431SBarry Smith .   w - work vector
792c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
7935e76c431SBarry Smith 
7945e76c431SBarry Smith     Output parameters for func:
795c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7965e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7975e76c431SBarry Smith .   gnorm - 2-norm of g
7985e76c431SBarry Smith .   ynorm - 2-norm of search length
799c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
800761aaf1bSLois Curfman McInnes            on failure.
801f59ffdedSLois Curfman McInnes 
80236851e7fSLois Curfman McInnes     Level: advanced
80336851e7fSLois Curfman McInnes 
804f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
805f59ffdedSLois Curfman McInnes 
8065d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8075d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8085e76c431SBarry Smith @*/
809329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
8105e76c431SBarry Smith {
811329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
81282bf6240SBarry Smith 
8133a40ed3dSBarry Smith   PetscFunctionBegin;
814c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
81582bf6240SBarry Smith   if (f) {
816af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
81782bf6240SBarry Smith   }
8183a40ed3dSBarry Smith   PetscFunctionReturn(0);
8195e76c431SBarry Smith }
82004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
821fb2e594dSBarry Smith EXTERN_C_BEGIN
8224a2ae208SSatish Balay #undef __FUNCT__
8234a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
824af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
82587828ca2SBarry Smith                          PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
82682bf6240SBarry Smith {
82782bf6240SBarry Smith   PetscFunctionBegin;
8286831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
8296831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
83082bf6240SBarry Smith   PetscFunctionReturn(0);
83182bf6240SBarry Smith }
832fb2e594dSBarry Smith EXTERN_C_END
83304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8344a2ae208SSatish Balay #undef __FUNCT__
8354a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
836c8dd0c0dSLois Curfman McInnes /*@C
837530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8386831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
839c8dd0c0dSLois Curfman McInnes 
840c8dd0c0dSLois Curfman McInnes    Input Parameters:
841c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
842c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
843c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
844c8dd0c0dSLois Curfman McInnes 
845c8dd0c0dSLois Curfman McInnes    Collective on SNES
846c8dd0c0dSLois Curfman McInnes 
847c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
848c8dd0c0dSLois Curfman McInnes .vb
849b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
850c8dd0c0dSLois Curfman McInnes .ve
851b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
852b1ae27eaSLois Curfman McInnes    on failure.
853c8dd0c0dSLois Curfman McInnes 
854c8dd0c0dSLois Curfman McInnes    Input parameters for func:
855c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
856c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8578f99978dSLois Curfman McInnes -  x - current candidate iterate
858c8dd0c0dSLois Curfman McInnes 
859c8dd0c0dSLois Curfman McInnes    Output parameters for func:
860c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
861c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
862c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
863c8dd0c0dSLois Curfman McInnes 
864c8dd0c0dSLois Curfman McInnes    Level: advanced
865c8dd0c0dSLois Curfman McInnes 
8668f99978dSLois Curfman McInnes    Notes:
867b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
868b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
869b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
870b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
871ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8728f99978dSLois Curfman McInnes 
873b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
874b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
875b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
876ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
877ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
878ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
879ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
880b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8818f99978dSLois Curfman McInnes 
882c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
883c8dd0c0dSLois Curfman McInnes 
884c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
885c8dd0c0dSLois Curfman McInnes @*/
8868f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
887c8dd0c0dSLois Curfman McInnes {
8888f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
889c8dd0c0dSLois Curfman McInnes 
890c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
891c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
892c8dd0c0dSLois Curfman McInnes   if (f) {
893c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
894c8dd0c0dSLois Curfman McInnes   }
895c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
896c8dd0c0dSLois Curfman McInnes }
897c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
898c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
8994a2ae208SSatish Balay #undef __FUNCT__
9004a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
9018f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
902c8dd0c0dSLois Curfman McInnes {
903c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9046831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
9056831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
906c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
907c8dd0c0dSLois Curfman McInnes }
908c8dd0c0dSLois Curfman McInnes EXTERN_C_END
909c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
91004d965bbSLois Curfman McInnes /*
911329e7e40SMatthew Knepley    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
912329e7e40SMatthew Knepley 
913329e7e40SMatthew Knepley    Input Parameter:
914329e7e40SMatthew Knepley .  snes - the SNES context
915329e7e40SMatthew Knepley 
916329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
917329e7e40SMatthew Knepley */
918329e7e40SMatthew Knepley #undef __FUNCT__
919329e7e40SMatthew Knepley #define __FUNCT__ "SNESPrintHelp_EQ_LS"
920329e7e40SMatthew Knepley static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
921329e7e40SMatthew Knepley {
922329e7e40SMatthew Knepley   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
923329e7e40SMatthew Knepley 
924329e7e40SMatthew Knepley   PetscFunctionBegin;
925329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
926329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
927329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
928329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
929329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
930329e7e40SMatthew Knepley   PetscFunctionReturn(0);
931329e7e40SMatthew Knepley }
932329e7e40SMatthew Knepley 
933329e7e40SMatthew Knepley /*
9346831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
93504d965bbSLois Curfman McInnes 
93604d965bbSLois Curfman McInnes    Input Parameters:
93704d965bbSLois Curfman McInnes .  SNES - the SNES context
93804d965bbSLois Curfman McInnes .  viewer - visualization context
93904d965bbSLois Curfman McInnes 
94004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
94104d965bbSLois Curfman McInnes */
9424a2ae208SSatish Balay #undef __FUNCT__
9434a2ae208SSatish Balay #define __FUNCT__ "SNESView_EQ_LS"
944b0a32e0cSBarry Smith static int SNESView_EQ_LS(SNES snes,PetscViewer viewer)
945a935fc98SLois Curfman McInnes {
9466831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
94719bcc07fSBarry Smith   char       *cstr;
94851695f54SLois Curfman McInnes   int        ierr;
9496831982aSBarry Smith   PetscTruth isascii;
950a935fc98SLois Curfman McInnes 
9513a40ed3dSBarry Smith   PetscFunctionBegin;
952b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
9530f5bd95cSBarry Smith   if (isascii) {
95419bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
95519bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
95619bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
95719bcc07fSBarry Smith     else                                                cstr = "unknown";
958b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
959b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9605cd90555SBarry Smith   } else {
96129bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
96219bcc07fSBarry Smith   }
9633a40ed3dSBarry Smith   PetscFunctionReturn(0);
964a935fc98SLois Curfman McInnes }
96504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
96604d965bbSLois Curfman McInnes /*
9676831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
96804d965bbSLois Curfman McInnes 
96904d965bbSLois Curfman McInnes    Input Parameter:
97004d965bbSLois Curfman McInnes .  snes - the SNES context
97104d965bbSLois Curfman McInnes 
97204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
97304d965bbSLois Curfman McInnes */
9744a2ae208SSatish Balay #undef __FUNCT__
9754a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions_EQ_LS"
976f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
9775e42470aSBarry Smith {
9786831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
979186905e3SBarry Smith   char       ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
980f1af5d2fSBarry Smith   int        ierr;
981f1af5d2fSBarry Smith   PetscTruth flg;
9825e42470aSBarry Smith 
9833a40ed3dSBarry Smith   PetscFunctionBegin;
984b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
98587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
98687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
98787828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
988186905e3SBarry Smith 
989b0a32e0cSBarry Smith     ierr = PetscOptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);CHKERRQ(ierr);
99025cdf11fSBarry Smith     if (flg) {
991c38d4ed2SBarry Smith       PetscTruth isbasic,isnonorms,isquad,iscubic;
9920f5bd95cSBarry Smith 
993186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[0],&isbasic);CHKERRQ(ierr);
994186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[1],&isnonorms);CHKERRQ(ierr);
995186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[2],&isquad);CHKERRQ(ierr);
996186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[3],&iscubic);CHKERRQ(ierr);
9970f5bd95cSBarry Smith 
9980f5bd95cSBarry Smith       if (isbasic) {
999af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
10000f5bd95cSBarry Smith       } else if (isnonorms) {
1001af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
10020f5bd95cSBarry Smith       } else if (isquad) {
1003af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
10040f5bd95cSBarry Smith       } else if (iscubic) {
1005af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
10065e42470aSBarry Smith       }
100729bbc08cSBarry Smith       else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
10085e42470aSBarry Smith     }
1009b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10103a40ed3dSBarry Smith   PetscFunctionReturn(0);
10115e42470aSBarry Smith }
101204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
101304d965bbSLois Curfman McInnes /*
10146831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
10156831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
101604d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
101704d965bbSLois Curfman McInnes 
101804d965bbSLois Curfman McInnes 
101904d965bbSLois Curfman McInnes    Input Parameter:
102004d965bbSLois Curfman McInnes .  snes - the SNES context
102104d965bbSLois Curfman McInnes 
102204d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
102304d965bbSLois Curfman McInnes  */
1024fb2e594dSBarry Smith EXTERN_C_BEGIN
10254a2ae208SSatish Balay #undef __FUNCT__
10264a2ae208SSatish Balay #define __FUNCT__ "SNESCreate_EQ_LS"
1027f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
10285e42470aSBarry Smith {
102982bf6240SBarry Smith   int        ierr;
10306831982aSBarry Smith   SNES_EQ_LS *neP;
10315e42470aSBarry Smith 
10323a40ed3dSBarry Smith   PetscFunctionBegin;
1033a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
103429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1035a8c6a408SBarry Smith   }
103682bf6240SBarry Smith 
1037f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
1038f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
1039f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
1040f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
1041329e7e40SMatthew Knepley   snes->printhelp       = SNESPrintHelp_EQ_LS;
1042f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
1043f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
10445baf8537SBarry Smith   snes->nwork           = 0;
10455e42470aSBarry Smith 
1046b0a32e0cSBarry Smith   ierr                  = PetscNew(SNES_EQ_LS,&neP);CHKERRQ(ierr);
1047b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_EQ_LS));
10485e42470aSBarry Smith   snes->data    	= (void*)neP;
10495e42470aSBarry Smith   neP->alpha		= 1.e-4;
10505e42470aSBarry Smith   neP->maxstep		= 1.e8;
10515e42470aSBarry Smith   neP->steptol		= 1.e-12;
10525e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1053c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1054c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1055c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
105682bf6240SBarry Smith 
10575d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
10585d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
105982bf6240SBarry Smith 
10603a40ed3dSBarry Smith   PetscFunctionReturn(0);
10615e42470aSBarry Smith }
1062fb2e594dSBarry Smith EXTERN_C_END
106382bf6240SBarry Smith 
106482bf6240SBarry Smith 
106582bf6240SBarry Smith 
1066