xref: /petsc/src/snes/impls/ls/ls.c (revision 5ed2d596fd9914de76abbaed56c65f4792ae57ae)
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 
55ca10a37SBarry Smith 
65ca10a37SBarry Smith #undef __FUNCT__
75ca10a37SBarry Smith #define __FUNCT__ "VecMaxScale_SNES"
85ca10a37SBarry Smith /*
95ca10a37SBarry Smith             max { p[i]/x[i] }
105ca10a37SBarry Smith */
115ca10a37SBarry Smith int VecMaxScale_SNES(Vec p,Vec x,PetscReal *m)
125ca10a37SBarry Smith {
135ca10a37SBarry Smith   int         ierr,i,n;
145ca10a37SBarry Smith   PetscScalar *pa,*xa;
155ca10a37SBarry Smith   PetscReal   t;
165ca10a37SBarry Smith   MPI_Comm    comm;
175ca10a37SBarry Smith 
185ca10a37SBarry Smith   PetscFunctionBegin;
195ca10a37SBarry Smith   ierr = VecGetLocalSize(p,&n);CHKERRQ(ierr);
205ca10a37SBarry Smith   ierr = PetscObjectGetComm((PetscObject)p,&comm);CHKERRQ(ierr);
215ca10a37SBarry Smith 
225ca10a37SBarry Smith   ierr = VecGetArray(p,&pa);CHKERRQ(ierr);
235ca10a37SBarry Smith   ierr = VecGetArray(x,&xa);CHKERRQ(ierr);
245ca10a37SBarry Smith   t = 0.0;
255ca10a37SBarry Smith   for ( i=0; i<n; i++) {
265ca10a37SBarry Smith     if (xa[i] != 0.0) {
275ca10a37SBarry Smith       t = PetscMax(PetscAbsScalar(pa[i]/xa[i]),t);
285ca10a37SBarry Smith     } else {
295ca10a37SBarry Smith       t = PetscMax(PetscAbsScalar(pa[i]),t);
305ca10a37SBarry Smith     }
315ca10a37SBarry Smith   }
325ca10a37SBarry Smith   ierr = MPI_Allreduce(&t,m,1,MPI_DOUBLE,MPI_MAX,comm);CHKERRQ(ierr);
335ca10a37SBarry Smith   ierr = VecRestoreArray(p,&pa);CHKERRQ(ierr);
345ca10a37SBarry Smith   ierr = VecRestoreArray(x,&xa);CHKERRQ(ierr);
355ca10a37SBarry Smith   PetscFunctionReturn(0);
365ca10a37SBarry Smith }
375ca10a37SBarry Smith 
388a5d9ee4SBarry Smith /*
398a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
408a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
4136669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
4236669109SBarry Smith     for this trick.
438a5d9ee4SBarry Smith */
444a2ae208SSatish Balay #undef __FUNCT__
454a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
468a5d9ee4SBarry Smith int SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
478a5d9ee4SBarry Smith {
488a5d9ee4SBarry Smith   PetscReal  a1;
498a5d9ee4SBarry Smith   int        ierr;
5036669109SBarry Smith   PetscTruth hastranspose;
518a5d9ee4SBarry Smith 
528a5d9ee4SBarry Smith   PetscFunctionBegin;
538a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
5436669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5536669109SBarry Smith   if (hastranspose) {
568a5d9ee4SBarry Smith     /* Compute || J^T F|| */
578a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
588a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
59b0a32e0cSBarry Smith     PetscLogInfo(0,"SNESSolve_EQ_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm);
608a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
6136669109SBarry Smith   } else {
6236669109SBarry Smith     Vec       work;
63ea709b57SSatish Balay     PetscScalar    result;
6436669109SBarry Smith     PetscReal wnorm;
6536669109SBarry Smith 
6636669109SBarry Smith     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
6736669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
6836669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
6936669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
7036669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
7136669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
7236669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
73b0a32e0cSBarry Smith     PetscLogInfo(0,"SNESSolve_EQ_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1);
7436669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
7536669109SBarry Smith   }
768a5d9ee4SBarry Smith   PetscFunctionReturn(0);
778a5d9ee4SBarry Smith }
788a5d9ee4SBarry Smith 
7974637425SBarry Smith /*
80*5ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
8174637425SBarry Smith */
824a2ae208SSatish Balay #undef __FUNCT__
834a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
8474637425SBarry Smith int SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
8574637425SBarry Smith {
8674637425SBarry Smith   PetscReal     a1,a2;
8774637425SBarry Smith   int           ierr;
8874637425SBarry Smith   PetscTruth    hastranspose;
89ea709b57SSatish Balay   PetscScalar   mone = -1.0;
9074637425SBarry Smith 
9174637425SBarry Smith   PetscFunctionBegin;
9274637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
9374637425SBarry Smith   if (hastranspose) {
9474637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
9574637425SBarry Smith     ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr);
9674637425SBarry Smith 
9774637425SBarry Smith     /* Compute || J^T W|| */
9874637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
9974637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
10074637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
1019994e62eSSatish Balay     if (a1 != 0) {
102b0a32e0cSBarry Smith       PetscLogInfo(0,"SNESSolve_EQ_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1);
10385471664SBarry Smith     }
10474637425SBarry Smith   }
10574637425SBarry Smith   PetscFunctionReturn(0);
10674637425SBarry Smith }
10774637425SBarry Smith 
10804d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
10904d965bbSLois Curfman McInnes 
11004d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
11104d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
11204d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
11304d965bbSLois Curfman McInnes      respectively.
11404d965bbSLois Curfman McInnes 
115fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
11604d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
11704d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
118fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
11904d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
12004d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
12104d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
12204d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
12304d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
12404d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
12504d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
12604d965bbSLois Curfman McInnes      for all nonlinear solvers.
12704d965bbSLois Curfman McInnes 
12804d965bbSLois Curfman McInnes      Another key routine is:
12904d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
13004d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
13104d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
13204d965bbSLois Curfman McInnes      SNESSolve() if necessary.
13304d965bbSLois Curfman McInnes 
13404d965bbSLois Curfman McInnes      Additional basic routines are:
13504d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
13604d965bbSLois Curfman McInnes                                       have actually been used.
13704d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
138186905e3SBarry Smith      SNESView().
13904d965bbSLois Curfman McInnes 
14004d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
14104d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
14204d965bbSLois Curfman McInnes      above description applies to these categories also.
14304d965bbSLois Curfman McInnes 
14404d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1455e42470aSBarry Smith /*
14604d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
14704d965bbSLois Curfman McInnes    method with a line search.
1485e76c431SBarry Smith 
14904d965bbSLois Curfman McInnes    Input Parameters:
15004d965bbSLois Curfman McInnes .  snes - the SNES context
1515e76c431SBarry Smith 
15204d965bbSLois Curfman McInnes    Output Parameter:
153c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
15404d965bbSLois Curfman McInnes 
15504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1565e76c431SBarry Smith 
1575e76c431SBarry Smith    Notes:
1585e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1595e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1605e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1615e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
162393d2d9aSLois Curfman McInnes    and Schnabel.
1635e76c431SBarry Smith */
1644a2ae208SSatish Balay #undef __FUNCT__
1654a2ae208SSatish Balay #define __FUNCT__ "SNESSolve_EQ_LS"
166f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
1675e76c431SBarry Smith {
1686831982aSBarry Smith   SNES_EQ_LS          *neP = (SNES_EQ_LS*)snes->data;
16984c569c9SBarry Smith   int                 maxits,i,ierr,lits,lsfail;
170112a2221SBarry Smith   MatStructure        flg = DIFFERENT_NONZERO_PATTERN;
171329f5518SBarry Smith   PetscReal           fnorm,gnorm,xnorm,ynorm;
1725e42470aSBarry Smith   Vec                 Y,X,F,G,W,TMP;
1735e76c431SBarry Smith 
1743a40ed3dSBarry Smith   PetscFunctionBegin;
175184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
176184914b5SBarry Smith 
1775e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1785e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
17939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1805e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1815e42470aSBarry Smith   G		= snes->work[1];
1825e42470aSBarry Smith   W		= snes->work[2];
1835e76c431SBarry Smith 
1844c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1854c49b128SBarry Smith   snes->iter = 0;
1864c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1875334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
188cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1890f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1905e42470aSBarry Smith   snes->norm = fnorm;
1910f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
19284c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
19394a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1945e76c431SBarry Smith 
195184914b5SBarry Smith   if (fnorm < snes->atol) {*outits = 0; snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
19694a9d846SBarry Smith 
197d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
198d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
199d034289bSBarry Smith 
2005e76c431SBarry Smith   for (i=0; i<maxits; i++) {
2015e76c431SBarry Smith 
202329e7e40SMatthew Knepley     /* Call general purpose update function */
203de8cb200SMatthew Knepley     if (snes->update != PETSC_NULL) {
204329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
205de8cb200SMatthew Knepley     }
206329e7e40SMatthew Knepley 
207ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2085334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
2095334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
21078b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits);CHKERRQ(ierr);
21174637425SBarry Smith 
212b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
21374637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
21485471664SBarry Smith     }
21574637425SBarry Smith 
21690cbd584SBarry Smith     /* should check what happened to the linear solve? */
2173505fcc1SBarry Smith     snes->linear_its += lits;
218b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
219ea4d3ed3SLois Curfman McInnes 
220ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
221ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
222ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
223ea4d3ed3SLois Curfman McInnes     */
22481b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
225af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail);CHKERRQ(ierr);
226*5ed2d596SBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
227a3b891d8SBarry Smith 
228a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
229a3b891d8SBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
230a3b891d8SBarry Smith     fnorm = gnorm;
231a3b891d8SBarry Smith 
232*5ed2d596SBarry Smith 
233*5ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
234*5ed2d596SBarry Smith     snes->iter = i+1;
235*5ed2d596SBarry Smith     snes->norm = fnorm;
236*5ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
237*5ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
238*5ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
239*5ed2d596SBarry Smith 
2403505fcc1SBarry Smith     if (lsfail) {
2418a5d9ee4SBarry Smith       PetscTruth ismin;
24250ffb88aSMatthew Knepley 
24350ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2443505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2458a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2468a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2473505fcc1SBarry Smith         break;
2483505fcc1SBarry Smith       }
24950ffb88aSMatthew Knepley     }
2505e76c431SBarry Smith 
2515e76c431SBarry Smith     /* Test for convergence */
25229e0b56fSBarry Smith     if (snes->converged) {
25329e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2543505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2553505fcc1SBarry Smith       if (snes->reason) {
25616c95cacSBarry Smith         break;
25716c95cacSBarry Smith       }
25816c95cacSBarry Smith     }
25929e0b56fSBarry Smith   }
26039e2f89bSBarry Smith   if (X != snes->vec_sol) {
261393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
262e15f7bb5SBarry Smith   }
263e15f7bb5SBarry Smith   if (F != snes->vec_func) {
264e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
265e15f7bb5SBarry Smith   }
26639e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
26739e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
26852392280SLois Curfman McInnes   if (i == maxits) {
269b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
27052392280SLois Curfman McInnes     i--;
2713505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
27252392280SLois Curfman McInnes   }
2730f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2740f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2755e42470aSBarry Smith   *outits = i+1;
2763a40ed3dSBarry Smith   PetscFunctionReturn(0);
2775e76c431SBarry Smith }
27804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27904d965bbSLois Curfman McInnes /*
28004d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
2816831982aSBarry Smith    of the SNESEQLS nonlinear solver.
28204d965bbSLois Curfman McInnes 
28304d965bbSLois Curfman McInnes    Input Parameter:
28404d965bbSLois Curfman McInnes .  snes - the SNES context
28504d965bbSLois Curfman McInnes .  x - the solution vector
28604d965bbSLois Curfman McInnes 
28704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
28804d965bbSLois Curfman McInnes 
28904d965bbSLois Curfman McInnes    Notes:
29004d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
29104d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
29204d965bbSLois Curfman McInnes    the call to SNESSolve().
29304d965bbSLois Curfman McInnes  */
2944a2ae208SSatish Balay #undef __FUNCT__
2954a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp_EQ_LS"
296f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
2975e76c431SBarry Smith {
2985e42470aSBarry Smith   int ierr;
2993a40ed3dSBarry Smith 
3003a40ed3dSBarry Smith   PetscFunctionBegin;
30181b6cf68SLois Curfman McInnes   snes->nwork = 4;
302d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
303b0a32e0cSBarry Smith   PetscLogObjectParents(snes,snes->nwork,snes->work);
30481b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
3053a40ed3dSBarry Smith   PetscFunctionReturn(0);
3065e76c431SBarry Smith }
30704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30804d965bbSLois Curfman McInnes /*
3096831982aSBarry Smith    SNESDestroy_EQ_LS - Destroys the private SNESEQLS context that was created
31004d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
31104d965bbSLois Curfman McInnes 
31204d965bbSLois Curfman McInnes    Input Parameter:
31304d965bbSLois Curfman McInnes .  snes - the SNES context
31404d965bbSLois Curfman McInnes 
315de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
31604d965bbSLois Curfman McInnes  */
3174a2ae208SSatish Balay #undef __FUNCT__
3184a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy_EQ_LS"
319e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
3205e76c431SBarry Smith {
321393d2d9aSLois Curfman McInnes   int  ierr;
3223a40ed3dSBarry Smith 
3233a40ed3dSBarry Smith   PetscFunctionBegin;
3245baf8537SBarry Smith   if (snes->nwork) {
3254b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
32621c89e3eSBarry Smith   }
3275c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
3295e76c431SBarry Smith }
33004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3314a2ae208SSatish Balay #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
33382bf6240SBarry Smith 
3344b828684SBarry Smith /*@C
3355e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3365e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3375e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3385e76c431SBarry Smith 
339c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
340c7afd0dbSLois Curfman McInnes 
3415e76c431SBarry Smith    Input Parameters:
342c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
343af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3445e76c431SBarry Smith .  x - current iterate
3455e76c431SBarry Smith .  f - residual evaluated at x
3465e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3475e76c431SBarry Smith .  w - work vector
348c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3495e76c431SBarry Smith 
350c4a48953SLois Curfman McInnes    Output Parameters:
351c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3525e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3535e76c431SBarry Smith .  gnorm - 2-norm of g
3545e76c431SBarry Smith .  ynorm - 2-norm of search length
355c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
356fee21e36SBarry Smith 
357c4a48953SLois Curfman McInnes    Options Database Key:
358c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
359c4a48953SLois Curfman McInnes 
36036851e7fSLois Curfman McInnes    Level: advanced
36136851e7fSLois Curfman McInnes 
36228ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
36328ae5a21SLois Curfman McInnes 
364f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
365af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3665e76c431SBarry Smith @*/
3675d1a10b1SBarry 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)
3685e76c431SBarry Smith {
3695e42470aSBarry Smith   int           ierr;
370ea709b57SSatish Balay   PetscScalar   mone = -1.0;
3716831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
3728f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
3735334005bSBarry Smith 
3743a40ed3dSBarry Smith   PetscFunctionBegin;
375761aaf1bSLois Curfman McInnes   *flag = 0;
376d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
377a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
378ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3798f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3808f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3818f99978dSLois Curfman McInnes   }
382ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
383a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
384d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3853a40ed3dSBarry Smith   PetscFunctionReturn(0);
3865e76c431SBarry Smith }
38704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3884a2ae208SSatish Balay #undef __FUNCT__
3894a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
39082bf6240SBarry Smith 
39129e0b56fSBarry Smith /*@C
39229e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
39329e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
39429e0b56fSBarry Smith    even compute the norm of the function or search direction; this
39529e0b56fSBarry Smith    is intended only when you know the full step is fine and are
39629e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
397c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
398c7afd0dbSLois Curfman McInnes 
399c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
40029e0b56fSBarry Smith 
40129e0b56fSBarry Smith    Input Parameters:
402c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
403af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
40429e0b56fSBarry Smith .  x - current iterate
40529e0b56fSBarry Smith .  f - residual evaluated at x
40629e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
40729e0b56fSBarry Smith .  w - work vector
408c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
40929e0b56fSBarry Smith 
41029e0b56fSBarry Smith    Output Parameters:
411c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
41229e0b56fSBarry Smith .  gnorm - not changed
41329e0b56fSBarry Smith .  ynorm - not changed
414c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
415fee21e36SBarry Smith 
41629e0b56fSBarry Smith    Options Database Key:
417c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
41829e0b56fSBarry Smith 
4198cbba510SBarry Smith    Notes:
420ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
421ea56f5baSLois Curfman McInnes    either the options
422ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
423ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
424329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
425329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
426329f5518SBarry Smith 
427329f5518SBarry Smith    During the final iteration this will not evaluate the function at
428329f5518SBarry Smith    the solution point. This is to save a function evaluation while
429329f5518SBarry Smith    using pseudo-timestepping.
4308cbba510SBarry Smith 
431ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
432ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
433ea56f5baSLois Curfman McInnes    correct, since they are not computed.
434ea56f5baSLois Curfman McInnes 
435ea56f5baSLois Curfman McInnes    Level: advanced
4368cbba510SBarry Smith 
43729e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
43829e0b56fSBarry Smith 
43929e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
44029e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
44129e0b56fSBarry Smith @*/
4425d1a10b1SBarry 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)
44329e0b56fSBarry Smith {
44429e0b56fSBarry Smith   int           ierr;
445ea709b57SSatish Balay   PetscScalar   mone = -1.0;
4466831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
4478f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
44829e0b56fSBarry Smith 
4493a40ed3dSBarry Smith   PetscFunctionBegin;
4508cbba510SBarry Smith   *flag = 0;
451d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
45229e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4538f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4548f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4558f99978dSLois Curfman McInnes   }
456329f5518SBarry Smith 
457329f5518SBarry Smith   /* don't evaluate function the last time through */
458329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
45929e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
460329f5518SBarry Smith   }
461d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4623a40ed3dSBarry Smith   PetscFunctionReturn(0);
46329e0b56fSBarry Smith }
46404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4654a2ae208SSatish Balay #undef __FUNCT__
4664a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4674b828684SBarry Smith /*@C
468f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4695e76c431SBarry Smith 
470c7afd0dbSLois Curfman McInnes    Collective on SNES
471c7afd0dbSLois Curfman McInnes 
4725e76c431SBarry Smith    Input Parameters:
473c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
474af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4755e76c431SBarry Smith .  x - current iterate
4765e76c431SBarry Smith .  f - residual evaluated at x
4775e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4785e76c431SBarry Smith .  w - work vector
479c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4805e76c431SBarry Smith 
481393d2d9aSLois Curfman McInnes    Output Parameters:
482c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4835e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4845e76c431SBarry Smith .  gnorm - 2-norm of g
4855e76c431SBarry Smith .  ynorm - 2-norm of search length
486c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
487fee21e36SBarry Smith 
488c4a48953SLois Curfman McInnes    Options Database Key:
489c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
490c4a48953SLois Curfman McInnes 
4915e76c431SBarry Smith    Notes:
4925e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4935e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4945e76c431SBarry Smith 
49536851e7fSLois Curfman McInnes    Level: advanced
49636851e7fSLois Curfman McInnes 
49728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
49828ae5a21SLois Curfman McInnes 
499af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
5005e76c431SBarry Smith @*/
5015d1a10b1SBarry 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)
5025e76c431SBarry Smith {
503406556e6SLois Curfman McInnes   /*
504406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
505406556e6SLois Curfman McInnes      minimization problem:
506406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
507406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
508406556e6SLois Curfman McInnes    */
509406556e6SLois Curfman McInnes 
5105ca10a37SBarry Smith   PetscReal     steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
511329f5518SBarry Smith   PetscReal     maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
512aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
513ea709b57SSatish Balay   PetscScalar   cinitslope,clambda;
5146b5873e3SBarry Smith #endif
5155e42470aSBarry Smith   int           ierr,count;
5166831982aSBarry Smith   SNES_EQ_LS    *neP = (SNES_EQ_LS*)snes->data;
517ea709b57SSatish Balay   PetscScalar   mone = -1.0,scale;
5188f99978dSLois Curfman McInnes   PetscTruth    change_y = PETSC_FALSE;
5195e76c431SBarry Smith 
5203a40ed3dSBarry Smith   PetscFunctionBegin;
521d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
522761aaf1bSLois Curfman McInnes   *flag   = 0;
5235e76c431SBarry Smith   alpha   = neP->alpha;
5245e76c431SBarry Smith   maxstep = neP->maxstep;
5255e76c431SBarry Smith   steptol = neP->steptol;
5265e76c431SBarry Smith 
527cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
528a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
529b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
530a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
531a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
532a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
533ad922ac9SBarry Smith     goto theend1;
53494a9d846SBarry Smith   }
5355e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5365e42470aSBarry Smith     scale = maxstep/(*ynorm);
537aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
538b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm);
5396b5873e3SBarry Smith #else
540b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm);
5416b5873e3SBarry Smith #endif
542393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5435e76c431SBarry Smith     *ynorm = maxstep;
5445e76c431SBarry Smith   }
5455ca10a37SBarry Smith   ierr      = VecMaxScale_SNES(y,x,&rellength);CHKERRQ(ierr);
5465ca10a37SBarry Smith   minlambda = steptol/rellength;
547a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
548aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
549a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
550329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5515e42470aSBarry Smith #else
552a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5535e42470aSBarry Smith #endif
5545e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5555e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5565e76c431SBarry Smith 
557393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5585334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
55978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
560313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
5615d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
562393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
563b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
564ad922ac9SBarry Smith     goto theend1;
5655e76c431SBarry Smith   }
5665e76c431SBarry Smith 
5675e76c431SBarry Smith   /* Fit points with quadratic */
568313b5538SBarry Smith   lambda = 1.0;
569a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5705e76c431SBarry Smith   lambdaprev = lambda;
5715e76c431SBarry Smith   gnormprev = *gnorm;
572329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
573ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
574ddd12767SBarry Smith   else                         lambda = lambdatemp;
575393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
576ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
577aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
578ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5795e42470aSBarry Smith #else
580ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5815e42470aSBarry Smith #endif
58278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
583cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
584*5ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
585393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
586*5ed2d596SBarry Smith     PetscLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda);
587ad922ac9SBarry Smith     goto theend1;
5885e76c431SBarry Smith   }
5895e76c431SBarry Smith 
5905e76c431SBarry Smith   /* Fit points with cubic */
5915e76c431SBarry Smith   count = 1;
5925ca10a37SBarry Smith   while (PETSC_TRUE) {
5935e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
594b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
595*5ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope);
596f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
597761aaf1bSLois Curfman McInnes       *flag = -1; break;
5985e76c431SBarry Smith     }
5995d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
6005d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
601ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6022b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
6035e76c431SBarry Smith     d  = b*b - 3*a*initslope;
6045e76c431SBarry Smith     if (d < 0.0) d = 0.0;
6055e76c431SBarry Smith     if (a == 0.0) {
6065e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
6075e76c431SBarry Smith     } else {
6085e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6095e76c431SBarry Smith     }
6105e76c431SBarry Smith     lambdaprev = lambda;
6115e76c431SBarry Smith     gnormprev  = *gnorm;
612329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
613329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6145e76c431SBarry Smith     else                         lambda     = lambdatemp;
615393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
616ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
617aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
618ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
619393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6205e42470aSBarry Smith #else
621ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6225e42470aSBarry Smith #endif
62378b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
624cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
625*5ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
626393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
627*5ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda);
6282715565aSLois Curfman McInnes       goto theend1;
6292b022350SLois Curfman McInnes     } else {
630*5ed2d596SBarry Smith       PetscLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda);
6315e76c431SBarry Smith     }
6325e76c431SBarry Smith     count++;
6335e76c431SBarry Smith   }
634ad922ac9SBarry Smith   theend1:
6358f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
6368f99978dSLois Curfman McInnes   if (neP->CheckStep) {
6378f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
6388f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
6398f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6408f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6418f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
6428f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6438f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6448f99978dSLois Curfman McInnes     }
6458f99978dSLois Curfman McInnes   }
646d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6473a40ed3dSBarry Smith   PetscFunctionReturn(0);
6485e76c431SBarry Smith }
64904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6504a2ae208SSatish Balay #undef __FUNCT__
6514a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6524b828684SBarry Smith /*@C
653f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6545e76c431SBarry Smith 
655c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
656c7afd0dbSLois Curfman McInnes 
6575e42470aSBarry Smith    Input Parameters:
658c7afd0dbSLois Curfman McInnes +  snes - the SNES context
659af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6605e42470aSBarry Smith .  x - current iterate
6615e42470aSBarry Smith .  f - residual evaluated at x
6625e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6635e42470aSBarry Smith .  w - work vector
664c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6655e42470aSBarry Smith 
666c4a48953SLois Curfman McInnes    Output Parameters:
667c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6685e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6695e42470aSBarry Smith .  gnorm - 2-norm of g
6705e42470aSBarry Smith .  ynorm - 2-norm of search length
671c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
672fee21e36SBarry Smith 
673c4a48953SLois Curfman McInnes    Options Database Key:
674c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
675c4a48953SLois Curfman McInnes 
6765e42470aSBarry Smith    Notes:
6776831982aSBarry Smith    Use SNESSetLineSearch() to set this routine within the SNESEQLS method.
6785e42470aSBarry Smith 
67936851e7fSLois Curfman McInnes    Level: advanced
68036851e7fSLois Curfman McInnes 
681f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
682f59ffdedSLois Curfman McInnes 
683af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6845e42470aSBarry Smith @*/
6855d1a10b1SBarry 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)
6865e76c431SBarry Smith {
687406556e6SLois Curfman McInnes   /*
688406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
689406556e6SLois Curfman McInnes      minimization problem:
690406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
691406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
692406556e6SLois Curfman McInnes    */
6935ca10a37SBarry Smith   PetscReal  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
694aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
695ea709b57SSatish Balay   PetscScalar    cinitslope,clambda;
6966b5873e3SBarry Smith #endif
6975e42470aSBarry Smith   int        ierr,count;
6986831982aSBarry Smith   SNES_EQ_LS     *neP = (SNES_EQ_LS*)snes->data;
699ea709b57SSatish Balay   PetscScalar    mone = -1.0,scale;
7008f99978dSLois Curfman McInnes   PetscTruth     change_y = PETSC_FALSE;
7015e76c431SBarry Smith 
7023a40ed3dSBarry Smith   PetscFunctionBegin;
703d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
704761aaf1bSLois Curfman McInnes   *flag   = 0;
7055e42470aSBarry Smith   alpha   = neP->alpha;
7065e42470aSBarry Smith   maxstep = neP->maxstep;
7075e42470aSBarry Smith   steptol = neP->steptol;
7085e76c431SBarry Smith 
7093505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
710b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
711b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
712b37302c6SLois Curfman McInnes     *gnorm = fnorm;
713b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y);CHKERRQ(ierr);
714b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g);CHKERRQ(ierr);
715ad922ac9SBarry Smith     goto theend2;
71694a9d846SBarry Smith   }
7175e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7185e42470aSBarry Smith     scale = maxstep/(*ynorm);
719393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
7205e42470aSBarry Smith     *ynorm = maxstep;
7215e76c431SBarry Smith   }
7225ca10a37SBarry Smith   ierr      = VecMaxScale_SNES(y,x,&rellength);CHKERRQ(ierr);
7235ca10a37SBarry Smith   minlambda = steptol/rellength;
724a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
725aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
726a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
727329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7285e42470aSBarry Smith #else
729a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7305e42470aSBarry Smith #endif
7315e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
7325e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7335e42470aSBarry Smith 
734393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7355334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
73678b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
737cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
7385d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
739393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
740b0a32e0cSBarry Smith     PetscLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
741ad922ac9SBarry Smith     goto theend2;
7425e42470aSBarry Smith   }
7435e42470aSBarry Smith 
7445e42470aSBarry Smith   /* Fit points with quadratic */
745313b5538SBarry Smith   lambda = 1.0;
7465e42470aSBarry Smith   count = 1;
7475ca10a37SBarry Smith   while (PETSC_TRUE) {
7485e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
749b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
750b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope);
751f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
752761aaf1bSLois Curfman McInnes       *flag = -1; break;
7535e42470aSBarry Smith     }
754a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
755329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
756329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
757329f5518SBarry Smith     else                         lambda     = lambdatemp;
758393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7593505fcc1SBarry Smith     lambdaneg = -lambda;
760aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7613505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7625e42470aSBarry Smith #else
7633505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7645e42470aSBarry Smith #endif
76578b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
766cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
767*5ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
768393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
769b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
77006259719SBarry Smith       break;
7715e42470aSBarry Smith     }
7725e42470aSBarry Smith     count++;
7735e42470aSBarry Smith   }
774ad922ac9SBarry Smith   theend2:
7758f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
7768f99978dSLois Curfman McInnes   if (neP->CheckStep) {
7778f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
7788f99978dSLois Curfman McInnes     if (change_y == PETSC_TRUE) { /* recompute the function if the step has changed */
7798f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
7808f99978dSLois Curfman McInnes       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
7818f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
7828f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
7838f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
7848f99978dSLois Curfman McInnes     }
7858f99978dSLois Curfman McInnes   }
786d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
7873a40ed3dSBarry Smith   PetscFunctionReturn(0);
7885e76c431SBarry Smith }
78904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
7904a2ae208SSatish Balay #undef __FUNCT__
7914a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
792c9e6a524SLois Curfman McInnes /*@C
79377c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
7946831982aSBarry Smith    by the method SNESEQLS.
7955e76c431SBarry Smith 
7965e76c431SBarry Smith    Input Parameters:
797c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
798af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
799c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8005e76c431SBarry Smith 
801fee21e36SBarry Smith    Collective on SNES
802fee21e36SBarry Smith 
803c4a48953SLois Curfman McInnes    Available Routines:
804c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
805f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
806f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
807af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8085e76c431SBarry Smith 
809c4a48953SLois Curfman McInnes     Options Database Keys:
810af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
811c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
812c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
8133304466cSBarry Smith -   -snes_eq_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8143304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8153304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
816c4a48953SLois Curfman McInnes 
8175e76c431SBarry Smith    Calling sequence of func:
818bd208895SLois Curfman McInnes .vb
819af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
820329f5518SBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,*flag)
821bd208895SLois Curfman McInnes .ve
8225e76c431SBarry Smith 
8235e76c431SBarry Smith     Input parameters for func:
824c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
825af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8265e76c431SBarry Smith .   x - current iterate
8275e76c431SBarry Smith .   f - residual evaluated at x
8285e76c431SBarry Smith .   y - search direction (contains new iterate on output)
8295e76c431SBarry Smith .   w - work vector
830c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8315e76c431SBarry Smith 
8325e76c431SBarry Smith     Output parameters for func:
833c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8345e76c431SBarry Smith .   y - new iterate (contains search direction on input)
8355e76c431SBarry Smith .   gnorm - 2-norm of g
8365e76c431SBarry Smith .   ynorm - 2-norm of search length
837c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
838761aaf1bSLois Curfman McInnes            on failure.
839f59ffdedSLois Curfman McInnes 
84036851e7fSLois Curfman McInnes     Level: advanced
84136851e7fSLois Curfman McInnes 
842f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
843f59ffdedSLois Curfman McInnes 
8445d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8455d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8465e76c431SBarry Smith @*/
847329f5518SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
8485e76c431SBarry Smith {
849329f5518SBarry Smith   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,int*),void*);
85082bf6240SBarry Smith 
8513a40ed3dSBarry Smith   PetscFunctionBegin;
852c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
85382bf6240SBarry Smith   if (f) {
854af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
85582bf6240SBarry Smith   }
8563a40ed3dSBarry Smith   PetscFunctionReturn(0);
8575e76c431SBarry Smith }
85804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
859fb2e594dSBarry Smith EXTERN_C_BEGIN
8604a2ae208SSatish Balay #undef __FUNCT__
8614a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
862af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
86387828ca2SBarry Smith                          PetscReal,PetscReal*,PetscReal*,int*),void *lsctx)
86482bf6240SBarry Smith {
86582bf6240SBarry Smith   PetscFunctionBegin;
8666831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->LineSearch = func;
8676831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->lsP        = lsctx;
86882bf6240SBarry Smith   PetscFunctionReturn(0);
86982bf6240SBarry Smith }
870fb2e594dSBarry Smith EXTERN_C_END
87104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8724a2ae208SSatish Balay #undef __FUNCT__
8734a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
874c8dd0c0dSLois Curfman McInnes /*@C
875530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
8766831982aSBarry Smith    by the line search routine in the Newton-based method SNESEQLS.
877c8dd0c0dSLois Curfman McInnes 
878c8dd0c0dSLois Curfman McInnes    Input Parameters:
879c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
880c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
881c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
882c8dd0c0dSLois Curfman McInnes 
883c8dd0c0dSLois Curfman McInnes    Collective on SNES
884c8dd0c0dSLois Curfman McInnes 
885c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
886c8dd0c0dSLois Curfman McInnes .vb
887b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
888c8dd0c0dSLois Curfman McInnes .ve
889b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
890b1ae27eaSLois Curfman McInnes    on failure.
891c8dd0c0dSLois Curfman McInnes 
892c8dd0c0dSLois Curfman McInnes    Input parameters for func:
893c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
894c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
8958f99978dSLois Curfman McInnes -  x - current candidate iterate
896c8dd0c0dSLois Curfman McInnes 
897c8dd0c0dSLois Curfman McInnes    Output parameters for func:
898c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
899c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
900c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
901c8dd0c0dSLois Curfman McInnes 
902c8dd0c0dSLois Curfman McInnes    Level: advanced
903c8dd0c0dSLois Curfman McInnes 
9048f99978dSLois Curfman McInnes    Notes:
905b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
906b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
907b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
908b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
909ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9108f99978dSLois Curfman McInnes 
911b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
912b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
913b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
914ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
915ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
916ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
917ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
918b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9198f99978dSLois Curfman McInnes 
920c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
921c8dd0c0dSLois Curfman McInnes 
922c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
923c8dd0c0dSLois Curfman McInnes @*/
9248f99978dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
925c8dd0c0dSLois Curfman McInnes {
9268f99978dSLois Curfman McInnes   int ierr,(*f)(SNES,int (*)(SNES,void*,Vec,PetscTruth*),void*);
927c8dd0c0dSLois Curfman McInnes 
928c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
929c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
930c8dd0c0dSLois Curfman McInnes   if (f) {
931c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
932c8dd0c0dSLois Curfman McInnes   }
933c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
934c8dd0c0dSLois Curfman McInnes }
935c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
936c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
9374a2ae208SSatish Balay #undef __FUNCT__
9384a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
9398f99978dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
940c8dd0c0dSLois Curfman McInnes {
941c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9426831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->CheckStep = func;
9436831982aSBarry Smith   ((SNES_EQ_LS *)(snes->data))->checkP    = checkctx;
944c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
945c8dd0c0dSLois Curfman McInnes }
946c8dd0c0dSLois Curfman McInnes EXTERN_C_END
947c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
94804d965bbSLois Curfman McInnes /*
949329e7e40SMatthew Knepley    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
950329e7e40SMatthew Knepley 
951329e7e40SMatthew Knepley    Input Parameter:
952329e7e40SMatthew Knepley .  snes - the SNES context
953329e7e40SMatthew Knepley 
954329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
955329e7e40SMatthew Knepley */
956329e7e40SMatthew Knepley #undef __FUNCT__
957329e7e40SMatthew Knepley #define __FUNCT__ "SNESPrintHelp_EQ_LS"
958329e7e40SMatthew Knepley static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
959329e7e40SMatthew Knepley {
960329e7e40SMatthew Knepley   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
961329e7e40SMatthew Knepley 
962329e7e40SMatthew Knepley   PetscFunctionBegin;
963329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
964329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
965329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
966329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
967329e7e40SMatthew Knepley   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
968329e7e40SMatthew Knepley   PetscFunctionReturn(0);
969329e7e40SMatthew Knepley }
970329e7e40SMatthew Knepley 
971329e7e40SMatthew Knepley /*
9726831982aSBarry Smith    SNESView_EQ_LS - Prints info from the SNESEQLS data structure.
97304d965bbSLois Curfman McInnes 
97404d965bbSLois Curfman McInnes    Input Parameters:
97504d965bbSLois Curfman McInnes .  SNES - the SNES context
97604d965bbSLois Curfman McInnes .  viewer - visualization context
97704d965bbSLois Curfman McInnes 
97804d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
97904d965bbSLois Curfman McInnes */
9804a2ae208SSatish Balay #undef __FUNCT__
9814a2ae208SSatish Balay #define __FUNCT__ "SNESView_EQ_LS"
982b0a32e0cSBarry Smith static int SNESView_EQ_LS(SNES snes,PetscViewer viewer)
983a935fc98SLois Curfman McInnes {
9846831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
98519bcc07fSBarry Smith   char       *cstr;
98651695f54SLois Curfman McInnes   int        ierr;
9876831982aSBarry Smith   PetscTruth isascii;
988a935fc98SLois Curfman McInnes 
9893a40ed3dSBarry Smith   PetscFunctionBegin;
990b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
9910f5bd95cSBarry Smith   if (isascii) {
99219bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
99319bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
99419bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
99519bcc07fSBarry Smith     else                                                cstr = "unknown";
996b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
997b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
9985cd90555SBarry Smith   } else {
99929bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
100019bcc07fSBarry Smith   }
10013a40ed3dSBarry Smith   PetscFunctionReturn(0);
1002a935fc98SLois Curfman McInnes }
100304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
100404d965bbSLois Curfman McInnes /*
10056831982aSBarry Smith    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNESEQLS method.
100604d965bbSLois Curfman McInnes 
100704d965bbSLois Curfman McInnes    Input Parameter:
100804d965bbSLois Curfman McInnes .  snes - the SNES context
100904d965bbSLois Curfman McInnes 
101004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
101104d965bbSLois Curfman McInnes */
10124a2ae208SSatish Balay #undef __FUNCT__
10134a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions_EQ_LS"
1014f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
10155e42470aSBarry Smith {
10166831982aSBarry Smith   SNES_EQ_LS *ls = (SNES_EQ_LS *)snes->data;
1017186905e3SBarry Smith   char       ver[16],*lses[] = {"basic","basicnonorms","quadratic","cubic"};
1018f1af5d2fSBarry Smith   int        ierr;
1019f1af5d2fSBarry Smith   PetscTruth flg;
10205e42470aSBarry Smith 
10213a40ed3dSBarry Smith   PetscFunctionBegin;
1022b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
102387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
102487828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
102587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_eq_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1026186905e3SBarry Smith 
1027b0a32e0cSBarry Smith     ierr = PetscOptionsEList("-snes_eq_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",ver,16,&flg);CHKERRQ(ierr);
102825cdf11fSBarry Smith     if (flg) {
1029c38d4ed2SBarry Smith       PetscTruth isbasic,isnonorms,isquad,iscubic;
10300f5bd95cSBarry Smith 
1031186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[0],&isbasic);CHKERRQ(ierr);
1032186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[1],&isnonorms);CHKERRQ(ierr);
1033186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[2],&isquad);CHKERRQ(ierr);
1034186905e3SBarry Smith       ierr = PetscStrcmp(ver,lses[3],&iscubic);CHKERRQ(ierr);
10350f5bd95cSBarry Smith 
10360f5bd95cSBarry Smith       if (isbasic) {
1037af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
10380f5bd95cSBarry Smith       } else if (isnonorms) {
1039af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
10400f5bd95cSBarry Smith       } else if (isquad) {
1041af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
10420f5bd95cSBarry Smith       } else if (iscubic) {
1043af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
10445e42470aSBarry Smith       }
104529bbc08cSBarry Smith       else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown line search");}
10465e42470aSBarry Smith     }
1047b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10483a40ed3dSBarry Smith   PetscFunctionReturn(0);
10495e42470aSBarry Smith }
105004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
105104d965bbSLois Curfman McInnes /*
10526831982aSBarry Smith    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNESEQLS method,
10536831982aSBarry Smith    SNES_EQ_LS, and sets this as the private data within the generic nonlinear solver
105404d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
105504d965bbSLois Curfman McInnes 
105604d965bbSLois Curfman McInnes 
105704d965bbSLois Curfman McInnes    Input Parameter:
105804d965bbSLois Curfman McInnes .  snes - the SNES context
105904d965bbSLois Curfman McInnes 
106004d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
106104d965bbSLois Curfman McInnes  */
1062fb2e594dSBarry Smith EXTERN_C_BEGIN
10634a2ae208SSatish Balay #undef __FUNCT__
10644a2ae208SSatish Balay #define __FUNCT__ "SNESCreate_EQ_LS"
1065f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
10665e42470aSBarry Smith {
106782bf6240SBarry Smith   int        ierr;
10686831982aSBarry Smith   SNES_EQ_LS *neP;
10695e42470aSBarry Smith 
10703a40ed3dSBarry Smith   PetscFunctionBegin;
1071a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
107229bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1073a8c6a408SBarry Smith   }
107482bf6240SBarry Smith 
1075f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
1076f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
1077f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
1078f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
1079329e7e40SMatthew Knepley   snes->printhelp       = SNESPrintHelp_EQ_LS;
1080f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
1081f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
10825baf8537SBarry Smith   snes->nwork           = 0;
10835e42470aSBarry Smith 
1084b0a32e0cSBarry Smith   ierr                  = PetscNew(SNES_EQ_LS,&neP);CHKERRQ(ierr);
1085b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_EQ_LS));
10865e42470aSBarry Smith   snes->data    	= (void*)neP;
10875e42470aSBarry Smith   neP->alpha		= 1.e-4;
10885e42470aSBarry Smith   neP->maxstep		= 1.e8;
10895e42470aSBarry Smith   neP->steptol		= 1.e-12;
10905e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1091c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1092c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1093c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
109482bf6240SBarry Smith 
10955d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
10965d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
109782bf6240SBarry Smith 
10983a40ed3dSBarry Smith   PetscFunctionReturn(0);
10995e42470aSBarry Smith }
1100fb2e594dSBarry Smith EXTERN_C_END
110182bf6240SBarry Smith 
110282bf6240SBarry Smith 
110382bf6240SBarry Smith 
1104