xref: /petsc/src/snes/impls/ls/ls.c (revision de8cb20032df8a6c332461b20fe56d14f0f1e0ae)
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 */
170*de8cb200SMatthew Knepley     if (snes->update != PETSC_NULL) {
171329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
172*de8cb200SMatthew 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;
2013505fcc1SBarry Smith       snes->nfailures++;
2023505fcc1SBarry Smith       snes->reason = SNES_DIVERGED_LS_FAILURE;
2038a5d9ee4SBarry Smith       ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2048a5d9ee4SBarry Smith       if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2053505fcc1SBarry Smith       break;
2063505fcc1SBarry Smith     }
2075e76c431SBarry Smith 
2080f5bd95cSBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
20984c569c9SBarry Smith     snes->iter = i+1;
2105e42470aSBarry Smith     snes->norm = fnorm;
2110f5bd95cSBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
21284c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
21394a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2145e76c431SBarry Smith 
2155e76c431SBarry Smith     /* Test for convergence */
21629e0b56fSBarry Smith     if (snes->converged) {
21729e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2183505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2193505fcc1SBarry Smith       if (snes->reason) {
22016c95cacSBarry Smith         break;
22116c95cacSBarry Smith       }
22216c95cacSBarry Smith     }
22329e0b56fSBarry Smith   }
22439e2f89bSBarry Smith   if (X != snes->vec_sol) {
225393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
226e15f7bb5SBarry Smith   }
227e15f7bb5SBarry Smith   if (F != snes->vec_func) {
228e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
229e15f7bb5SBarry Smith   }
23039e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
23139e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
23252392280SLois Curfman McInnes   if (i == maxits) {
233b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
23452392280SLois Curfman McInnes     i--;
2353505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
23652392280SLois Curfman McInnes   }
2370f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2380f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2395e42470aSBarry Smith   *outits = i+1;
2403a40ed3dSBarry Smith   PetscFunctionReturn(0);
2415e76c431SBarry Smith }
24204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
24304d965bbSLois Curfman McInnes /*
24404d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
2456831982aSBarry Smith    of the SNESEQLS nonlinear solver.
24604d965bbSLois Curfman McInnes 
24704d965bbSLois Curfman McInnes    Input Parameter:
24804d965bbSLois Curfman McInnes .  snes - the SNES context
24904d965bbSLois Curfman McInnes .  x - the solution vector
25004d965bbSLois Curfman McInnes 
25104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25204d965bbSLois Curfman McInnes 
25304d965bbSLois Curfman McInnes    Notes:
25404d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
25504d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
25604d965bbSLois Curfman McInnes    the call to SNESSolve().
25704d965bbSLois Curfman McInnes  */
2584a2ae208SSatish Balay #undef __FUNCT__
2594a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp_EQ_LS"
260f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
2615e76c431SBarry Smith {
2625e42470aSBarry Smith   int ierr;
2633a40ed3dSBarry Smith 
2643a40ed3dSBarry Smith   PetscFunctionBegin;
26581b6cf68SLois Curfman McInnes   snes->nwork = 4;
266d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
267b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
26881b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2693a40ed3dSBarry Smith   PetscFunctionReturn(0);
2705e76c431SBarry Smith }
27104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27204d965bbSLois Curfman McInnes /*
2736831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
27404d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
27504d965bbSLois Curfman McInnes 
27604d965bbSLois Curfman McInnes    Input Parameter:
27704d965bbSLois Curfman McInnes .  snes - the SNES context
27804d965bbSLois Curfman McInnes 
279de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28004d965bbSLois Curfman McInnes  */
2814a2ae208SSatish Balay #undef __FUNCT__
2824a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy_EQ_LS"
283e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
2845e76c431SBarry Smith {
285393d2d9aSLois Curfman McInnes   int  ierr;
2863a40ed3dSBarry Smith 
2873a40ed3dSBarry Smith   PetscFunctionBegin;
2885baf8537SBarry Smith   if (snes->nwork) {
2894b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29021c89e3eSBarry Smith   }
2915c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2923a40ed3dSBarry Smith   PetscFunctionReturn(0);
2935e76c431SBarry Smith }
29404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2954a2ae208SSatish Balay #undef __FUNCT__
2964a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
29782bf6240SBarry Smith 
2984b828684SBarry Smith /*@C
2995e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3005e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3015e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3025e76c431SBarry Smith 
303c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
304c7afd0dbSLois Curfman McInnes 
3055e76c431SBarry Smith    Input Parameters:
306c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
307af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3085e76c431SBarry Smith .  x - current iterate
3095e76c431SBarry Smith .  f - residual evaluated at x
3105e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3115e76c431SBarry Smith .  w - work vector
312c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3135e76c431SBarry Smith 
314c4a48953SLois Curfman McInnes    Output Parameters:
315c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3165e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3175e76c431SBarry Smith .  gnorm - 2-norm of g
3185e76c431SBarry Smith .  ynorm - 2-norm of search length
319c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
320fee21e36SBarry Smith 
321c4a48953SLois Curfman McInnes    Options Database Key:
322c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
323c4a48953SLois Curfman McInnes 
32436851e7fSLois Curfman McInnes    Level: advanced
32536851e7fSLois Curfman McInnes 
32628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
32728ae5a21SLois Curfman McInnes 
328f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
329af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3305e76c431SBarry Smith @*/
3315d1a10b1SBarry Smith int SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
3325e76c431SBarry Smith {
3335e42470aSBarry Smith   int           ierr;
334ea709b57SSatish Balay   PetscScalar   mone = -1.0;
3356831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
3368f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
3375334005bSBarry Smith 
3383a40ed3dSBarry Smith   PetscFunctionBegin;
339761aaf1bSLois Curfman McInnes   *flag = 0;
340b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
341a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
342ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3438f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3448f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3458f99978dSLois Curfman McInnes   }
346ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
347a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
348b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3493a40ed3dSBarry Smith   PetscFunctionReturn(0);
3505e76c431SBarry Smith }
35104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3524a2ae208SSatish Balay #undef __FUNCT__
3534a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
35482bf6240SBarry Smith 
35529e0b56fSBarry Smith /*@C
35629e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
35729e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
35829e0b56fSBarry Smith    even compute the norm of the function or search direction; this
35929e0b56fSBarry Smith    is intended only when you know the full step is fine and are
36029e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
361c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
362c7afd0dbSLois Curfman McInnes 
363c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
36429e0b56fSBarry Smith 
36529e0b56fSBarry Smith    Input Parameters:
366c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
367af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
36829e0b56fSBarry Smith .  x - current iterate
36929e0b56fSBarry Smith .  f - residual evaluated at x
37029e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
37129e0b56fSBarry Smith .  w - work vector
372c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
37329e0b56fSBarry Smith 
37429e0b56fSBarry Smith    Output Parameters:
375c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
37629e0b56fSBarry Smith .  gnorm - not changed
37729e0b56fSBarry Smith .  ynorm - not changed
378c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
379fee21e36SBarry Smith 
38029e0b56fSBarry Smith    Options Database Key:
381c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
38229e0b56fSBarry Smith 
3838cbba510SBarry Smith    Notes:
384ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
385ea56f5baSLois Curfman McInnes    either the options
386ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
387ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
388329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
389329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
390329f5518SBarry Smith 
391329f5518SBarry Smith    During the final iteration this will not evaluate the function at
392329f5518SBarry Smith    the solution point. This is to save a function evaluation while
393329f5518SBarry Smith    using pseudo-timestepping.
3948cbba510SBarry Smith 
395ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
396ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
397ea56f5baSLois Curfman McInnes    correct, since they are not computed.
398ea56f5baSLois Curfman McInnes 
399ea56f5baSLois Curfman McInnes    Level: advanced
4008cbba510SBarry Smith 
40129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
40229e0b56fSBarry Smith 
40329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
40429e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
40529e0b56fSBarry Smith @*/
4065d1a10b1SBarry Smith int SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
40729e0b56fSBarry Smith {
40829e0b56fSBarry Smith   int           ierr;
409ea709b57SSatish Balay   PetscScalar   mone = -1.0;
4106831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
4118f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
41229e0b56fSBarry Smith 
4133a40ed3dSBarry Smith   PetscFunctionBegin;
4148cbba510SBarry Smith   *flag = 0;
415b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
41629e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4178f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4188f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4198f99978dSLois Curfman McInnes   }
420329f5518SBarry Smith 
421329f5518SBarry Smith   /* don't evaluate function the last time through */
422329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
42329e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
424329f5518SBarry Smith   }
425b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
42729e0b56fSBarry Smith }
42804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4294a2ae208SSatish Balay #undef __FUNCT__
4304a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4314b828684SBarry Smith /*@C
432f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4335e76c431SBarry Smith 
434c7afd0dbSLois Curfman McInnes    Collective on SNES
435c7afd0dbSLois Curfman McInnes 
4365e76c431SBarry Smith    Input Parameters:
437c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
438af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4395e76c431SBarry Smith .  x - current iterate
4405e76c431SBarry Smith .  f - residual evaluated at x
4415e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4425e76c431SBarry Smith .  w - work vector
443c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4445e76c431SBarry Smith 
445393d2d9aSLois Curfman McInnes    Output Parameters:
446c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4475e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4485e76c431SBarry Smith .  gnorm - 2-norm of g
4495e76c431SBarry Smith .  ynorm - 2-norm of search length
450c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
451fee21e36SBarry Smith 
452c4a48953SLois Curfman McInnes    Options Database Key:
453c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
454c4a48953SLois Curfman McInnes 
4555e76c431SBarry Smith    Notes:
4565e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4575e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4585e76c431SBarry Smith 
45936851e7fSLois Curfman McInnes    Level: advanced
46036851e7fSLois Curfman McInnes 
46128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
46228ae5a21SLois Curfman McInnes 
463af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4645e76c431SBarry Smith @*/
4655d1a10b1SBarry Smith int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,int *flag)
4665e76c431SBarry Smith {
467406556e6SLois Curfman McInnes   /*
468406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
469406556e6SLois Curfman McInnes      minimization problem:
470406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
471406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
472406556e6SLois Curfman McInnes    */
473406556e6SLois Curfman McInnes 
474329f5518SBarry Smith   PetscReal     steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2;
475329f5518SBarry Smith   PetscReal     maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
476aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
477ea709b57SSatish Balay   PetscScalar   cinitslope,clambda;
4786b5873e3SBarry Smith #endif
4795e42470aSBarry Smith   int           ierr,count;
4806831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
481ea709b57SSatish Balay   PetscScalar   mone = -1.0,scale;
4828f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
4835e76c431SBarry Smith 
4843a40ed3dSBarry Smith   PetscFunctionBegin;
485b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
486761aaf1bSLois Curfman McInnes   *flag   = 0;
4875e76c431SBarry Smith   alpha   = neP->alpha;
4885e76c431SBarry Smith   maxstep = neP->maxstep;
4895e76c431SBarry Smith   steptol = neP->steptol;
4905e76c431SBarry Smith 
491cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
492a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
493b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
494a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
495a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
496a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
497ad922ac9SBarry Smith     goto theend1;
49894a9d846SBarry Smith   }
4995e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5005e42470aSBarry Smith     scale = maxstep/(*ynorm);
501aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
502b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
5036b5873e3SBarry Smith #else
504b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
5056b5873e3SBarry Smith #endif
506393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5075e76c431SBarry Smith     *ynorm = maxstep;
5085e76c431SBarry Smith   }
5095e76c431SBarry Smith   minlambda = steptol/(*ynorm);
510a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
511aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
512a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
513329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5145e42470aSBarry Smith #else
515a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5165e42470aSBarry Smith #endif
5175e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5185e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5195e76c431SBarry Smith 
520393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5215334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
52278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
523313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5245d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
525393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
526b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
527ad922ac9SBarry Smith     goto theend1;
5285e76c431SBarry Smith   }
5295e76c431SBarry Smith 
5305e76c431SBarry Smith   /* Fit points with quadratic */
531313b5538SBarry Smith   lambda = 1.0;
532a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5335e76c431SBarry Smith   lambdaprev = lambda;
5345e76c431SBarry Smith   gnormprev = *gnorm;
535329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
536ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
537ddd12767SBarry Smith   else                         lambda = lambdatemp;
538393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
539ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
540aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
541ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5425e42470aSBarry Smith #else
543ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5445e42470aSBarry Smith #endif
54578b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
546cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5475d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
548393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
549b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
550ad922ac9SBarry Smith     goto theend1;
5515e76c431SBarry Smith   }
5525e76c431SBarry Smith 
5535e76c431SBarry Smith   /* Fit points with cubic */
5545e76c431SBarry Smith   count = 1;
5555e76c431SBarry Smith   while (1) {
5565e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
557b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
558b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
559393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
560761aaf1bSLois Curfman McInnes       *flag = -1; break;
5615e76c431SBarry Smith     }
5625d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
5635d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
564ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5652b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5665e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5675e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5685e76c431SBarry Smith     if (a == 0.0) {
5695e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5705e76c431SBarry Smith     } else {
5715e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
5725e76c431SBarry Smith     }
5735e76c431SBarry Smith     lambdaprev = lambda;
5745e76c431SBarry Smith     gnormprev  = *gnorm;
575329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
576329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
5775e76c431SBarry Smith     else                         lambda     = lambdatemp;
578393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
579ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
580aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
581ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
582393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5835e42470aSBarry Smith #else
584ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5855e42470aSBarry Smith #endif
58678b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
587cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5885d1a10b1SBarry Smith     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
589393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
590b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
5912715565aSLois Curfman McInnes       goto theend1;
5922b022350SLois Curfman McInnes     } else {
593b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
5945e76c431SBarry Smith     }
5955e76c431SBarry Smith     count++;
5965e76c431SBarry Smith   }
597ad922ac9SBarry Smith   theend1:
5988f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
5998f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6008f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6018f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6028f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6038f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6048f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6058f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6068f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6078f99978dSLois Curfman McInnes     }
6088f99978dSLois Curfman McInnes   }
609b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6103a40ed3dSBarry Smith   PetscFunctionReturn(0);
6115e76c431SBarry Smith }
61204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6134a2ae208SSatish Balay #undef __FUNCT__
6144a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6154b828684SBarry Smith /*@C
616f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6175e76c431SBarry Smith 
618c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
619c7afd0dbSLois Curfman McInnes 
6205e42470aSBarry Smith    Input Parameters:
621c7afd0dbSLois Curfman McInnes +  snes - the SNES context
622af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6235e42470aSBarry Smith .  x - current iterate
6245e42470aSBarry Smith .  f - residual evaluated at x
6255e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6265e42470aSBarry Smith .  w - work vector
627c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6285e42470aSBarry Smith 
629c4a48953SLois Curfman McInnes    Output Parameters:
630c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6315e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6325e42470aSBarry Smith .  gnorm - 2-norm of g
6335e42470aSBarry Smith .  ynorm - 2-norm of search length
634c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
635fee21e36SBarry Smith 
636c4a48953SLois Curfman McInnes    Options Database Key:
637c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
638c4a48953SLois Curfman McInnes 
6395e42470aSBarry Smith    Notes:
6406831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
6415e42470aSBarry Smith 
64236851e7fSLois Curfman McInnes    Level: advanced
64336851e7fSLois Curfman McInnes 
644f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
645f59ffdedSLois Curfman McInnes 
646af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6475e42470aSBarry Smith @*/
6485d1a10b1SBarry 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)
6495e76c431SBarry Smith {
650406556e6SLois Curfman McInnes   /*
651406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
652406556e6SLois Curfman McInnes      minimization problem:
653406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
654406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
655406556e6SLois Curfman McInnes    */
656329f5518SBarry Smith   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
657aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
658ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
6596b5873e3SBarry Smith #endif
6605e42470aSBarry Smith   int        ierr,count;
6616831982aSBarry Smith   SNES_EQ_LS     *neP = (SNES_EQ_LS*)snes->data;
662ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
6638f99978dSLois Curfman McInnes   PetscTruth     change_y = PETSC_FALSE;
6645e76c431SBarry Smith 
6653a40ed3dSBarry Smith   PetscFunctionBegin;
666b0a32e0cSBarry Smith   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
667761aaf1bSLois Curfman McInnes   *flag   = 0;
6685e42470aSBarry Smith   alpha   = neP->alpha;
6695e42470aSBarry Smith   maxstep = neP->maxstep;
6705e42470aSBarry Smith   steptol = neP->steptol;
6715e76c431SBarry Smith 
6723505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
673b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
674b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
675b37302c6SLois Curfman McInnes     *gnorm = fnorm;
676b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
677b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
678ad922ac9SBarry Smith     goto theend2;
67994a9d846SBarry Smith   }
6805e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
6815e42470aSBarry Smith     scale = maxstep/(*ynorm);
682393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
6835e42470aSBarry Smith     *ynorm = maxstep;
6845e76c431SBarry Smith   }
6855e42470aSBarry Smith   minlambda = steptol/(*ynorm);
686a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
687aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
688a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
689329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
6905e42470aSBarry Smith #else
691a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
6925e42470aSBarry Smith #endif
6935e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
6945e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
6955e42470aSBarry Smith 
696393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
6975334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
69878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
699cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7005d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
701393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
702b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
703ad922ac9SBarry Smith     goto theend2;
7045e42470aSBarry Smith   }
7055e42470aSBarry Smith 
7065e42470aSBarry Smith   /* Fit points with quadratic */
707313b5538SBarry Smith   lambda = 1.0;
7085e42470aSBarry Smith   count = 1;
7095e42470aSBarry Smith   while (1) {
7105e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
711b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
712b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
713393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
714761aaf1bSLois Curfman McInnes       *flag = -1; break;
7155e42470aSBarry Smith     }
716a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
717329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
718329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
719329f5518SBarry Smith     else                         lambda     = lambdatemp;
720393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7213505fcc1SBarry Smith     lambdaneg = -lambda;
722aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7233505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7245e42470aSBarry Smith #else
7253505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7265e42470aSBarry Smith #endif
72778b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
728cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7295d1a10b1SBarry Smith     if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
730393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
731b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
73206259719SBarry Smith       break;
7335e42470aSBarry Smith     }
7345e42470aSBarry Smith     count++;
7355e42470aSBarry Smith   }
736ad922ac9SBarry Smith   theend2:
7378f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7388f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7398f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7408f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7418f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7428f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7438f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7448f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7458f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7468f99978dSLois Curfman McInnes     }
7478f99978dSLois Curfman McInnes   }
748b0a32e0cSBarry Smith   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7493a40ed3dSBarry Smith   PetscFunctionReturn(0);
7505e76c431SBarry Smith }
75104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7524a2ae208SSatish Balay #undef __FUNCT__
7534a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
754c9e6a524SLois Curfman McInnes /*@C
75577c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7566831982aSBarry Smith    by the method SNESEQLS.
7575e76c431SBarry Smith 
7585e76c431SBarry Smith    Input Parameters:
759c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
760af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
761c7afd0dbSLois Curfman McInnes -  func - pointer to int function
7625e76c431SBarry Smith 
763fee21e36SBarry Smith    Collective on SNES
764fee21e36SBarry Smith 
765c4a48953SLois Curfman McInnes    Available Routines:
766c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
767f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
768f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
769af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
7705e76c431SBarry Smith 
771c4a48953SLois Curfman McInnes     Options Database Keys:
772af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
773c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
774c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
775c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
776c4a48953SLois Curfman McInnes 
7775e76c431SBarry Smith    Calling sequence of func:
778bd208895SLois Curfman McInnes .vb
779af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
780329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
781bd208895SLois Curfman McInnes .ve
7825e76c431SBarry Smith 
7835e76c431SBarry Smith     Input parameters for func:
784c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
785af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
7865e76c431SBarry Smith .   x - current iterate
7875e76c431SBarry Smith .   f - residual evaluated at x
7885e76c431SBarry Smith .   y - search direction (contains new iterate on output)
7895e76c431SBarry Smith .   w - work vector
790c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
7915e76c431SBarry Smith 
7925e76c431SBarry Smith     Output parameters for func:
793c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
7945e76c431SBarry Smith .   y - new iterate (contains search direction on input)
7955e76c431SBarry Smith .   gnorm - 2-norm of g
7965e76c431SBarry Smith .   ynorm - 2-norm of search length
797c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
798761aaf1bSLois Curfman McInnes            on failure.
799f59ffdedSLois Curfman McInnes 
80036851e7fSLois Curfman McInnes     Level: advanced
80136851e7fSLois Curfman McInnes 
802f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
803f59ffdedSLois Curfman McInnes 
8045d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8055d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8065e76c431SBarry Smith @*/
807329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
8085e76c431SBarry Smith {
809329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
81082bf6240SBarry Smith 
8113a40ed3dSBarry Smith   PetscFunctionBegin;
812b9617806SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)())&f);CHKERRQ(ierr);
81382bf6240SBarry Smith   if (f) {
814af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
81582bf6240SBarry Smith   }
8163a40ed3dSBarry Smith   PetscFunctionReturn(0);
8175e76c431SBarry Smith }
81804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
819fb2e594dSBarry Smith EXTERN_C_BEGIN
8204a2ae208SSatish Balay #undef __FUNCT__
8214a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
822af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
82387828ca2SBarry Smith                          PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
82482bf6240SBarry Smith {
82582bf6240SBarry Smith   PetscFunctionBegin;
8266831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
8276831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
82882bf6240SBarry Smith   PetscFunctionReturn(0);
82982bf6240SBarry Smith }
830fb2e594dSBarry Smith EXTERN_C_END
83104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8324a2ae208SSatish Balay #undef __FUNCT__
8334a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
834c8dd0c0dSLois Curfman McInnes /*@C
835530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8366831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
837c8dd0c0dSLois Curfman McInnes 
838c8dd0c0dSLois Curfman McInnes    Input Parameters:
839c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
840c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
841c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
842c8dd0c0dSLois Curfman McInnes 
843c8dd0c0dSLois Curfman McInnes    Collective on SNES
844c8dd0c0dSLois Curfman McInnes 
845c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
846c8dd0c0dSLois Curfman McInnes .vb
847b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
848c8dd0c0dSLois Curfman McInnes .ve
849b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
850b1ae27eaSLois Curfman McInnes    on failure.
851c8dd0c0dSLois Curfman McInnes 
852c8dd0c0dSLois Curfman McInnes    Input parameters for func:
853c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
854c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8558f99978dSLois Curfman McInnes -  x - current candidate iterate
856c8dd0c0dSLois Curfman McInnes 
857c8dd0c0dSLois Curfman McInnes    Output parameters for func:
858c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
859c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
860c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
861c8dd0c0dSLois Curfman McInnes 
862c8dd0c0dSLois Curfman McInnes    Level: advanced
863c8dd0c0dSLois Curfman McInnes 
8648f99978dSLois Curfman McInnes    Notes:
865b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
866b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
867b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
868b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
869ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
8708f99978dSLois Curfman McInnes 
871b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
872b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
873b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
874ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
875ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
876ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
877ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
878b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
8798f99978dSLois Curfman McInnes 
880c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
881c8dd0c0dSLois Curfman McInnes 
882c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
883c8dd0c0dSLois Curfman McInnes @*/
8848f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
885c8dd0c0dSLois Curfman McInnes {
8868f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
887c8dd0c0dSLois Curfman McInnes 
888c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
889b9617806SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)())&f);CHKERRQ(ierr);
890c8dd0c0dSLois Curfman McInnes   if (f) {
891c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
892c8dd0c0dSLois Curfman McInnes   }
893c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
894c8dd0c0dSLois Curfman McInnes }
895c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
896c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
8974a2ae208SSatish Balay #undef __FUNCT__
8984a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
8998f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
900c8dd0c0dSLois Curfman McInnes {
901c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9026831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
9036831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
904c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
905c8dd0c0dSLois Curfman McInnes }
906c8dd0c0dSLois Curfman McInnes EXTERN_C_END
907c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
90804d965bbSLois Curfman McInnes /*
909329e7e40SMatthew Knepley    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
910329e7e40SMatthew Knepley 
911329e7e40SMatthew Knepley    Input Parameter:
912329e7e40SMatthew Knepley .  snes - the SNES context
913329e7e40SMatthew Knepley 
914329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
915329e7e40SMatthew Knepley */
916329e7e40SMatthew Knepley #undef __FUNCT__
917329e7e40SMatthew Knepley #define __FUNCT__ "SNESPrintHelp_EQ_LS"
918329e7e40SMatthew Knepley static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
919329e7e40SMatthew Knepley {
920329e7e40SMatthew Knepley   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
921329e7e40SMatthew Knepley 
922329e7e40SMatthew Knepley   PetscFunctionBegin;
923329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
924329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
925329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
926329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
927329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
928329e7e40SMatthew Knepley   PetscFunctionReturn(0);
929329e7e40SMatthew Knepley }
930329e7e40SMatthew Knepley 
931329e7e40SMatthew Knepley /*
9326831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
93304d965bbSLois Curfman McInnes 
93404d965bbSLois Curfman McInnes    Input Parameters:
93504d965bbSLois Curfman McInnes .  SNES - the SNES context
93604d965bbSLois Curfman McInnes .  viewer - visualization context
93704d965bbSLois Curfman McInnes 
93804d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
93904d965bbSLois Curfman McInnes */
9404a2ae208SSatish Balay #undef __FUNCT__
9414a2ae208SSatish Balay #define __FUNCT__ "SNESView_EQ_LS"
942b0a32e0cSBarry Smith static int SNESView_EQ_LS(SNES snes,PetscViewer viewer)
943a935fc98SLois Curfman McInnes {
9446831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
94519bcc07fSBarry Smith   char       *cstr;
94651695f54SLois Curfman McInnes   int        ierr;
9476831982aSBarry Smith   PetscTruth isascii;
948a935fc98SLois Curfman McInnes 
9493a40ed3dSBarry Smith   PetscFunctionBegin;
950b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
9510f5bd95cSBarry Smith   if (isascii) {
95219bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
95319bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
95419bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
95519bcc07fSBarry Smith     else                                                cstr = "unknown";
956b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
957b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9585cd90555SBarry Smith   } else {
95929bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
96019bcc07fSBarry Smith   }
9613a40ed3dSBarry Smith   PetscFunctionReturn(0);
962a935fc98SLois Curfman McInnes }
96304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
96404d965bbSLois Curfman McInnes /*
9656831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
96604d965bbSLois Curfman McInnes 
96704d965bbSLois Curfman McInnes    Input Parameter:
96804d965bbSLois Curfman McInnes .  snes - the SNES context
96904d965bbSLois Curfman McInnes 
97004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
97104d965bbSLois Curfman McInnes */
9724a2ae208SSatish Balay #undef __FUNCT__
9734a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions_EQ_LS"
974f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
9755e42470aSBarry Smith {
9766831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
977186905e3SBarry Smith   char       ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
978f1af5d2fSBarry Smith   int        ierr;
979f1af5d2fSBarry Smith   PetscTruth flg;
9805e42470aSBarry Smith 
9813a40ed3dSBarry Smith   PetscFunctionBegin;
982b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
98387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
98487828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
98587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
986186905e3SBarry Smith 
987b0a32e0cSBarry Smith     ierr = PetscOptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);CHKERRQ(ierr);
98825cdf11fSBarry Smith     if (flg) {
989c38d4ed2SBarry Smith       PetscTruth isbasic,isnonorms,isquad,iscubic;
9900f5bd95cSBarry Smith 
991186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[0],&isbasic);CHKERRQ(ierr);
992186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[1],&isnonorms);CHKERRQ(ierr);
993186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[2],&isquad);CHKERRQ(ierr);
994186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[3],&iscubic);CHKERRQ(ierr);
9950f5bd95cSBarry Smith 
9960f5bd95cSBarry Smith       if (isbasic) {
997af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
9980f5bd95cSBarry Smith       } else if (isnonorms) {
999af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
10000f5bd95cSBarry Smith       } else if (isquad) {
1001af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
10020f5bd95cSBarry Smith       } else if (iscubic) {
1003af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
10045e42470aSBarry Smith       }
100529bbc08cSBarry Smith       else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
10065e42470aSBarry Smith     }
1007b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10083a40ed3dSBarry Smith   PetscFunctionReturn(0);
10095e42470aSBarry Smith }
101004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
101104d965bbSLois Curfman McInnes /*
10126831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
10136831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
101404d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
101504d965bbSLois Curfman McInnes 
101604d965bbSLois Curfman McInnes 
101704d965bbSLois Curfman McInnes    Input Parameter:
101804d965bbSLois Curfman McInnes .  snes - the SNES context
101904d965bbSLois Curfman McInnes 
102004d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
102104d965bbSLois Curfman McInnes  */
1022fb2e594dSBarry Smith EXTERN_C_BEGIN
10234a2ae208SSatish Balay #undef __FUNCT__
10244a2ae208SSatish Balay #define __FUNCT__ "SNESCreate_EQ_LS"
1025f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
10265e42470aSBarry Smith {
102782bf6240SBarry Smith   int        ierr;
10286831982aSBarry Smith   SNES_EQ_LS *neP;
10295e42470aSBarry Smith 
10303a40ed3dSBarry Smith   PetscFunctionBegin;
1031a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
103229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1033a8c6a408SBarry Smith   }
103482bf6240SBarry Smith 
1035f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
1036f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
1037f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
1038f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
1039329e7e40SMatthew Knepley   snes->printhelp       = SNESPrintHelp_EQ_LS;
1040f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
1041f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
10425baf8537SBarry Smith   snes->nwork           = 0;
10435e42470aSBarry Smith 
1044b0a32e0cSBarry Smith   ierr                  = PetscNew(SNES_EQ_LS,&neP);CHKERRQ(ierr);
1045b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_EQ_LS));
10465e42470aSBarry Smith   snes->data    	= (void*)neP;
10475e42470aSBarry Smith   neP->alpha		= 1.e-4;
10485e42470aSBarry Smith   neP->maxstep		= 1.e8;
10495e42470aSBarry Smith   neP->steptol		= 1.e-12;
10505e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1051c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1052c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1053c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
105482bf6240SBarry Smith 
10555d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
10565d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
105782bf6240SBarry Smith 
10583a40ed3dSBarry Smith   PetscFunctionReturn(0);
10595e42470aSBarry Smith }
1060fb2e594dSBarry Smith EXTERN_C_END
106182bf6240SBarry Smith 
106282bf6240SBarry Smith 
106382bf6240SBarry Smith 
1064