xref: /petsc/src/snes/impls/ls/ls.c (revision c8dd0c0d1465a0be5d606d0075fe0d07d9b5cf5f)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*c8dd0c0dSLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.125 1999/03/14 17:17:00 curfman Exp curfman $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
570f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
65e42470aSBarry Smith 
704d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
804d965bbSLois Curfman McInnes 
904d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
1004d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
1104d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
1204d965bbSLois Curfman McInnes      respectively.
1304d965bbSLois Curfman McInnes 
14fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
1504d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
1604d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
17fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
1804d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
1904d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
2004d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
2104d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
2204d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
2304d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
2404d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
2504d965bbSLois Curfman McInnes      for all nonlinear solvers.
2604d965bbSLois Curfman McInnes 
2704d965bbSLois Curfman McInnes      Another key routine is:
2804d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
2904d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
3004d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
3104d965bbSLois Curfman McInnes      SNESSolve() if necessary.
3204d965bbSLois Curfman McInnes 
3304d965bbSLois Curfman McInnes      Additional basic routines are:
3404d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
3504d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
3604d965bbSLois Curfman McInnes                                       have actually been used.
3704d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
3804d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
3904d965bbSLois Curfman McInnes 
4004d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
4104d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
4204d965bbSLois Curfman McInnes      above description applies to these categories also.
4304d965bbSLois Curfman McInnes 
4404d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
455e42470aSBarry Smith /*
4604d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
4704d965bbSLois Curfman McInnes    method with a line search.
485e76c431SBarry Smith 
4904d965bbSLois Curfman McInnes    Input Parameters:
5004d965bbSLois Curfman McInnes .  snes - the SNES context
515e76c431SBarry Smith 
5204d965bbSLois Curfman McInnes    Output Parameter:
53c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
5404d965bbSLois Curfman McInnes 
5504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
565e76c431SBarry Smith 
575e76c431SBarry Smith    Notes:
585e76c431SBarry Smith    This implements essentially a truncated Newton method with a
595e76c431SBarry Smith    line search.  By default a cubic backtracking line search
605e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
615e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
62393d2d9aSLois Curfman McInnes    and Schnabel.
635e76c431SBarry Smith */
645615d1e5SSatish Balay #undef __FUNC__
655615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS"
66f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
675e76c431SBarry Smith {
685e42470aSBarry Smith   SNES_LS       *neP = (SNES_LS *) snes->data;
6984c569c9SBarry Smith   int           maxits, i, ierr, lits, lsfail;
70112a2221SBarry Smith   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
7184c569c9SBarry Smith   double        fnorm, gnorm, xnorm, ynorm;
725e42470aSBarry Smith   Vec           Y, X, F, G, W, TMP;
735e76c431SBarry Smith 
743a40ed3dSBarry Smith   PetscFunctionBegin;
755e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
765e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
7739e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
785e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
795e42470aSBarry Smith   G		= snes->work[1];
805e42470aSBarry Smith   W		= snes->work[2];
815e76c431SBarry Smith 
825334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
83cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
8406d0569aSBarry Smith   PetscAMSTakeAccess(snes);
8584c569c9SBarry Smith   snes->iter = 0;
865e42470aSBarry Smith   snes->norm = fnorm;
8706d0569aSBarry Smith   PetscAMSGrantAccess(snes);
8884c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
8994a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
905e76c431SBarry Smith 
913a40ed3dSBarry Smith   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
9294a9d846SBarry Smith 
93d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
94d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
95d034289bSBarry Smith 
965e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
975e76c431SBarry Smith 
98ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
995334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
1005334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
10178b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
1027a00f4a9SLois Curfman McInnes     snes->linear_its += PetscAbsInt(lits);
103af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
104ea4d3ed3SLois Curfman McInnes 
105ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
106ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
107ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
108ea4d3ed3SLois Curfman McInnes     */
10981b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
110af2d14d2SLois Curfman McInnes     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
111af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
112761aaf1bSLois Curfman McInnes     if (lsfail) snes->nfailures++;
1135e76c431SBarry Smith 
11439e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
11539e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1165e76c431SBarry Smith     fnorm = gnorm;
1175e76c431SBarry Smith 
11806d0569aSBarry Smith     PetscAMSTakeAccess(snes);
11984c569c9SBarry Smith     snes->iter = i+1;
1205e42470aSBarry Smith     snes->norm = fnorm;
12106d0569aSBarry Smith     PetscAMSGrantAccess(snes);
12284c569c9SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
12394a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1245e76c431SBarry Smith 
1255e76c431SBarry Smith     /* Test for convergence */
12629e0b56fSBarry Smith     if (snes->converged) {
12729e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
128bbb6d6a8SBarry Smith       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
12916c95cacSBarry Smith         break;
13016c95cacSBarry Smith       }
13116c95cacSBarry Smith     }
13229e0b56fSBarry Smith   }
13339e2f89bSBarry Smith   if (X != snes->vec_sol) {
134393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
13539e2f89bSBarry Smith     snes->vec_sol_always  = snes->vec_sol;
13639e2f89bSBarry Smith     snes->vec_func_always = snes->vec_func;
13739e2f89bSBarry Smith   }
13852392280SLois Curfman McInnes   if (i == maxits) {
139981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
14052392280SLois Curfman McInnes     i--;
14152392280SLois Curfman McInnes   }
1425e42470aSBarry Smith   *outits = i+1;
1433a40ed3dSBarry Smith   PetscFunctionReturn(0);
1445e76c431SBarry Smith }
14504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
14604d965bbSLois Curfman McInnes /*
14704d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
14804d965bbSLois Curfman McInnes    of the SNES_EQ_LS nonlinear solver.
14904d965bbSLois Curfman McInnes 
15004d965bbSLois Curfman McInnes    Input Parameter:
15104d965bbSLois Curfman McInnes .  snes - the SNES context
15204d965bbSLois Curfman McInnes .  x - the solution vector
15304d965bbSLois Curfman McInnes 
15404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
15504d965bbSLois Curfman McInnes 
15604d965bbSLois Curfman McInnes    Notes:
15704d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
15804d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
15904d965bbSLois Curfman McInnes    the call to SNESSolve().
16004d965bbSLois Curfman McInnes  */
1615615d1e5SSatish Balay #undef __FUNC__
1625615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
163f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
1645e76c431SBarry Smith {
1655e42470aSBarry Smith   int ierr;
1663a40ed3dSBarry Smith 
1673a40ed3dSBarry Smith   PetscFunctionBegin;
16881b6cf68SLois Curfman McInnes   snes->nwork = 4;
169d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1705e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
17181b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1723a40ed3dSBarry Smith   PetscFunctionReturn(0);
1735e76c431SBarry Smith }
17404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
17504d965bbSLois Curfman McInnes /*
17604d965bbSLois Curfman McInnes    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
17704d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
17804d965bbSLois Curfman McInnes 
17904d965bbSLois Curfman McInnes    Input Parameter:
18004d965bbSLois Curfman McInnes .  snes - the SNES context
18104d965bbSLois Curfman McInnes 
182de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
18304d965bbSLois Curfman McInnes  */
1845615d1e5SSatish Balay #undef __FUNC__
185d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
186e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
1875e76c431SBarry Smith {
188393d2d9aSLois Curfman McInnes   int  ierr;
1893a40ed3dSBarry Smith 
1903a40ed3dSBarry Smith   PetscFunctionBegin;
1915baf8537SBarry Smith   if (snes->nwork) {
1924b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
19321c89e3eSBarry Smith   }
1940452661fSBarry Smith   PetscFree(snes->data);
1953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1965e76c431SBarry Smith }
19704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1985615d1e5SSatish Balay #undef __FUNC__
1995615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
20082bf6240SBarry Smith 
2014b828684SBarry Smith /*@C
2025e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2035e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2045e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2055e76c431SBarry Smith 
206c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
207c7afd0dbSLois Curfman McInnes 
2085e76c431SBarry Smith    Input Parameters:
209c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
210af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
2115e76c431SBarry Smith .  x - current iterate
2125e76c431SBarry Smith .  f - residual evaluated at x
2135e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2145e76c431SBarry Smith .  w - work vector
215c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
2165e76c431SBarry Smith 
217c4a48953SLois Curfman McInnes    Output Parameters:
218c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
2195e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2205e76c431SBarry Smith .  gnorm - 2-norm of g
2215e76c431SBarry Smith .  ynorm - 2-norm of search length
222c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
223fee21e36SBarry Smith 
224c4a48953SLois Curfman McInnes    Options Database Key:
225c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
226c4a48953SLois Curfman McInnes 
22736851e7fSLois Curfman McInnes    Level: advanced
22836851e7fSLois Curfman McInnes 
22928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
23028ae5a21SLois Curfman McInnes 
231f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
232af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
2335e76c431SBarry Smith @*/
234af2d14d2SLois Curfman McInnes int SNESNoLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
235761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
2365e76c431SBarry Smith {
2375e42470aSBarry Smith   int    ierr;
2385334005bSBarry Smith   Scalar mone = -1.0;
2395334005bSBarry Smith 
2403a40ed3dSBarry Smith   PetscFunctionBegin;
241761aaf1bSLois Curfman McInnes   *flag = 0;
2427857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
243cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
244ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
245ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
246cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
2477857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
2495e76c431SBarry Smith }
25004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2515615d1e5SSatish Balay #undef __FUNC__
25229e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
25382bf6240SBarry Smith 
25429e0b56fSBarry Smith /*@C
25529e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
25629e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
25729e0b56fSBarry Smith    even compute the norm of the function or search direction; this
25829e0b56fSBarry Smith    is intended only when you know the full step is fine and are
25929e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
260c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
261c7afd0dbSLois Curfman McInnes 
262c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
26329e0b56fSBarry Smith 
26429e0b56fSBarry Smith    Input Parameters:
265c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
266af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
26729e0b56fSBarry Smith .  x - current iterate
26829e0b56fSBarry Smith .  f - residual evaluated at x
26929e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
27029e0b56fSBarry Smith .  w - work vector
271c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
27229e0b56fSBarry Smith 
27329e0b56fSBarry Smith    Output Parameters:
274c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
27529e0b56fSBarry Smith .  gnorm - not changed
27629e0b56fSBarry Smith .  ynorm - not changed
277c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
278fee21e36SBarry Smith 
27929e0b56fSBarry Smith    Options Database Key:
280c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
28129e0b56fSBarry Smith 
28236851e7fSLois Curfman McInnes    Level: advanced
28336851e7fSLois Curfman McInnes 
28429e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
28529e0b56fSBarry Smith 
28629e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
28729e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
28829e0b56fSBarry Smith @*/
289af2d14d2SLois Curfman McInnes int SNESNoLineSearchNoNorms(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
29029e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
29129e0b56fSBarry Smith {
29229e0b56fSBarry Smith   int    ierr;
29329e0b56fSBarry Smith   Scalar mone = -1.0;
29429e0b56fSBarry Smith 
2953a40ed3dSBarry Smith   PetscFunctionBegin;
29629e0b56fSBarry Smith   *flag = 0;
29729e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
29829e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
29929e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
30029e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3013a40ed3dSBarry Smith   PetscFunctionReturn(0);
30229e0b56fSBarry Smith }
30304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30429e0b56fSBarry Smith #undef __FUNC__
3055615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3064b828684SBarry Smith /*@C
307f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3085e76c431SBarry Smith 
309c7afd0dbSLois Curfman McInnes    Collective on SNES
310c7afd0dbSLois Curfman McInnes 
3115e76c431SBarry Smith    Input Parameters:
312c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
313af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3145e76c431SBarry Smith .  x - current iterate
3155e76c431SBarry Smith .  f - residual evaluated at x
3165e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3175e76c431SBarry Smith .  w - work vector
318c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3195e76c431SBarry Smith 
320393d2d9aSLois Curfman McInnes    Output Parameters:
321c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3225e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3235e76c431SBarry Smith .  gnorm - 2-norm of g
3245e76c431SBarry Smith .  ynorm - 2-norm of search length
325c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
326fee21e36SBarry Smith 
327c4a48953SLois Curfman McInnes    Options Database Key:
328c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
329c4a48953SLois Curfman McInnes 
3305e76c431SBarry Smith    Notes:
3315e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3325e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3335e76c431SBarry Smith 
33436851e7fSLois Curfman McInnes    Level: advanced
33536851e7fSLois Curfman McInnes 
33628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33728ae5a21SLois Curfman McInnes 
338af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3395e76c431SBarry Smith @*/
340af2d14d2SLois Curfman McInnes int SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
341761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3425e76c431SBarry Smith {
343406556e6SLois Curfman McInnes   /*
344406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
345406556e6SLois Curfman McInnes      minimization problem:
346406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
347406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
348406556e6SLois Curfman McInnes    */
349406556e6SLois Curfman McInnes 
350ddd12767SBarry Smith   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
351ea4d3ed3SLois Curfman McInnes   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
3523a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
3535e42470aSBarry Smith   Scalar  cinitslope, clambda;
3546b5873e3SBarry Smith #endif
3555e42470aSBarry Smith   int     ierr, count;
3565e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
3575334005bSBarry Smith   Scalar  mone = -1.0,scale;
3585e76c431SBarry Smith 
3593a40ed3dSBarry Smith   PetscFunctionBegin;
3607857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
361761aaf1bSLois Curfman McInnes   *flag   = 0;
3625e76c431SBarry Smith   alpha   = neP->alpha;
3635e76c431SBarry Smith   maxstep = neP->maxstep;
3645e76c431SBarry Smith   steptol = neP->steptol;
3655e76c431SBarry Smith 
366cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
367a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
368a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
369a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
370a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
371a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
372ad922ac9SBarry Smith     goto theend1;
37394a9d846SBarry Smith   }
3745e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3755e42470aSBarry Smith     scale = maxstep/(*ynorm);
3763a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
37792318cfeSSatish Balay     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
3786b5873e3SBarry Smith #else
37994a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
3806b5873e3SBarry Smith #endif
381393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
3825e76c431SBarry Smith     *ynorm = maxstep;
3835e76c431SBarry Smith   }
3845e76c431SBarry Smith   minlambda = steptol/(*ynorm);
385a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
3863a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
387a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
38892318cfeSSatish Balay   initslope = PetscReal(cinitslope);
3895e42470aSBarry Smith #else
390a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
3915e42470aSBarry Smith #endif
3925e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3935e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3945e76c431SBarry Smith 
395393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
3965334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
39778b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
398cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);
399406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
400393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
40194a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
402ad922ac9SBarry Smith     goto theend1;
4035e76c431SBarry Smith   }
4045e76c431SBarry Smith 
4055e76c431SBarry Smith   /* Fit points with quadratic */
4065e76c431SBarry Smith   lambda = 1.0; count = 0;
407a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4085e76c431SBarry Smith   lambdaprev = lambda;
4095e76c431SBarry Smith   gnormprev = *gnorm;
410ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
411ddd12767SBarry Smith   else lambda = lambdatemp;
412393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w); CHKERRQ(ierr);
413ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
4143a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
415ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4165e42470aSBarry Smith #else
417ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4185e42470aSBarry Smith #endif
41978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
420cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
421406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
422393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
423ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
424ad922ac9SBarry Smith     goto theend1;
4255e76c431SBarry Smith   }
4265e76c431SBarry Smith 
4275e76c431SBarry Smith   /* Fit points with cubic */
4285e76c431SBarry Smith   count = 1;
4295e76c431SBarry Smith   while (1) {
4305e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4312b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4322b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
433ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
434393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
435761aaf1bSLois Curfman McInnes       *flag = -1; break;
4365e76c431SBarry Smith     }
437406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
438406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
439ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4402b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4415e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4425e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4435e76c431SBarry Smith     if (a == 0.0) {
4445e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4455e76c431SBarry Smith     } else {
4465e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4475e76c431SBarry Smith     }
4485e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4495e76c431SBarry Smith       lambdatemp = .5*lambda;
4505e76c431SBarry Smith     }
4515e76c431SBarry Smith     lambdaprev = lambda;
4525e76c431SBarry Smith     gnormprev = *gnorm;
4535e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4545e76c431SBarry Smith       lambda = .1*lambda;
4555e76c431SBarry Smith     }
4565e76c431SBarry Smith     else lambda = lambdatemp;
457393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w ); CHKERRQ(ierr);
458ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
4593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
460ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
461393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4625e42470aSBarry Smith #else
463ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4645e42470aSBarry Smith #endif
46578b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
466cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
467406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
468393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
469ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
4702715565aSLois Curfman McInnes       goto theend1;
4712b022350SLois Curfman McInnes     } else {
4722b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
4735e76c431SBarry Smith     }
4745e76c431SBarry Smith     count++;
4755e76c431SBarry Smith   }
476ad922ac9SBarry Smith   theend1:
4777857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4783a40ed3dSBarry Smith   PetscFunctionReturn(0);
4795e76c431SBarry Smith }
48004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4815615d1e5SSatish Balay #undef __FUNC__
4825615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
4834b828684SBarry Smith /*@C
484f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
4855e76c431SBarry Smith 
486c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
487c7afd0dbSLois Curfman McInnes 
4885e42470aSBarry Smith    Input Parameters:
489c7afd0dbSLois Curfman McInnes +  snes - the SNES context
490af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4915e42470aSBarry Smith .  x - current iterate
4925e42470aSBarry Smith .  f - residual evaluated at x
4935e42470aSBarry Smith .  y - search direction (contains new iterate on output)
4945e42470aSBarry Smith .  w - work vector
495c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4965e42470aSBarry Smith 
497c4a48953SLois Curfman McInnes    Output Parameters:
498c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4995e42470aSBarry Smith .  y - new iterate (contains search direction on input)
5005e42470aSBarry Smith .  gnorm - 2-norm of g
5015e42470aSBarry Smith .  ynorm - 2-norm of search length
502c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
503fee21e36SBarry Smith 
504c4a48953SLois Curfman McInnes    Options Database Key:
505c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
506c4a48953SLois Curfman McInnes 
5075e42470aSBarry Smith    Notes:
508c7afd0dbSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
5095e42470aSBarry Smith 
51036851e7fSLois Curfman McInnes    Level: advanced
51136851e7fSLois Curfman McInnes 
512f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
513f59ffdedSLois Curfman McInnes 
514af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
5155e42470aSBarry Smith @*/
516af2d14d2SLois Curfman McInnes int SNESQuadraticLineSearch(SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
517761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5185e76c431SBarry Smith {
519406556e6SLois Curfman McInnes   /*
520406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
521406556e6SLois Curfman McInnes      minimization problem:
522406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
523406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
524406556e6SLois Curfman McInnes    */
52540011551SBarry Smith   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
5263a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
5275e42470aSBarry Smith   Scalar  cinitslope,clambda;
5286b5873e3SBarry Smith #endif
5295e42470aSBarry Smith   int     ierr,count;
5305e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
5315334005bSBarry Smith   Scalar  mone = -1.0,scale;
5325e76c431SBarry Smith 
5333a40ed3dSBarry Smith   PetscFunctionBegin;
5347857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
535761aaf1bSLois Curfman McInnes   *flag   = 0;
5365e42470aSBarry Smith   alpha   = neP->alpha;
5375e42470aSBarry Smith   maxstep = neP->maxstep;
5385e42470aSBarry Smith   steptol = neP->steptol;
5395e76c431SBarry Smith 
540cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
541b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
54294a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
543b37302c6SLois Curfman McInnes     *gnorm = fnorm;
544b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
545b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
546ad922ac9SBarry Smith     goto theend2;
54794a9d846SBarry Smith   }
5485e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5495e42470aSBarry Smith     scale = maxstep/(*ynorm);
550393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
5515e42470aSBarry Smith     *ynorm = maxstep;
5525e76c431SBarry Smith   }
5535e42470aSBarry Smith   minlambda = steptol/(*ynorm);
554a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
5553a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
556a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
55792318cfeSSatish Balay   initslope = PetscReal(cinitslope);
5585e42470aSBarry Smith #else
559a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
5605e42470aSBarry Smith #endif
5615e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
5625e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
5635e42470aSBarry Smith 
564393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
5655334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
56678b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
567cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
568406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
569393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
57094a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
571ad922ac9SBarry Smith     goto theend2;
5725e42470aSBarry Smith   }
5735e42470aSBarry Smith 
5745e42470aSBarry Smith   /* Fit points with quadratic */
5755e42470aSBarry Smith   lambda = 1.0; count = 0;
5765e42470aSBarry Smith   count = 1;
5775e42470aSBarry Smith   while (1) {
5785e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
579981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
580981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
581ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
582393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
583761aaf1bSLois Curfman McInnes       *flag = -1; break;
5845e42470aSBarry Smith     }
585a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5865e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
5875e42470aSBarry Smith       lambda = .1*lambda;
5885e42470aSBarry Smith     } else lambda = lambdatemp;
589393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w); CHKERRQ(ierr);
5905334005bSBarry Smith     lambda = -lambda;
5913a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
592393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
5935e42470aSBarry Smith #else
594393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
5955e42470aSBarry Smith #endif
59678b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
597cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
598406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
599393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
600981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
60106259719SBarry Smith       break;
6025e42470aSBarry Smith     }
6035e42470aSBarry Smith     count++;
6045e42470aSBarry Smith   }
605ad922ac9SBarry Smith   theend2:
6067857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
6073a40ed3dSBarry Smith   PetscFunctionReturn(0);
6085e76c431SBarry Smith }
60904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6105615d1e5SSatish Balay #undef __FUNC__
611d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
612c9e6a524SLois Curfman McInnes /*@C
61377c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
614f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
6155e76c431SBarry Smith 
6165e76c431SBarry Smith    Input Parameters:
617c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
618af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
619c7afd0dbSLois Curfman McInnes -  func - pointer to int function
6205e76c431SBarry Smith 
621fee21e36SBarry Smith    Collective on SNES
622fee21e36SBarry Smith 
623c4a48953SLois Curfman McInnes    Available Routines:
624c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
625f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
626f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
627af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
6285e76c431SBarry Smith 
629c4a48953SLois Curfman McInnes     Options Database Keys:
630af2d14d2SLois Curfman McInnes +   -snes_eq_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
631c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
632c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
633c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
634c4a48953SLois Curfman McInnes 
6355e76c431SBarry Smith    Calling sequence of func:
636bd208895SLois Curfman McInnes .vb
637af2d14d2SLois Curfman McInnes    func (SNES snes, void *lsctx, Vec x, Vec f, Vec g, Vec y, Vec w,
638bd208895SLois Curfman McInnes          double fnorm, double *ynorm, double *gnorm, *flag)
639bd208895SLois Curfman McInnes .ve
6405e76c431SBarry Smith 
6415e76c431SBarry Smith     Input parameters for func:
642c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
643af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
6445e76c431SBarry Smith .   x - current iterate
6455e76c431SBarry Smith .   f - residual evaluated at x
6465e76c431SBarry Smith .   y - search direction (contains new iterate on output)
6475e76c431SBarry Smith .   w - work vector
648c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
6495e76c431SBarry Smith 
6505e76c431SBarry Smith     Output parameters for func:
651c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
6525e76c431SBarry Smith .   y - new iterate (contains search direction on input)
6535e76c431SBarry Smith .   gnorm - 2-norm of g
6545e76c431SBarry Smith .   ynorm - 2-norm of search length
655c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
656761aaf1bSLois Curfman McInnes            on failure.
657f59ffdedSLois Curfman McInnes 
65836851e7fSLois Curfman McInnes     Level: advanced
65936851e7fSLois Curfman McInnes 
660f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
661f59ffdedSLois Curfman McInnes 
662af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms()
6635e76c431SBarry Smith @*/
664af2d14d2SLois Curfman McInnes int SNESSetLineSearch(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void *lsctx)
6655e76c431SBarry Smith {
666af2d14d2SLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*),void*);
66782bf6240SBarry Smith 
6683a40ed3dSBarry Smith   PetscFunctionBegin;
66994d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
67082bf6240SBarry Smith   if (f) {
671af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
67282bf6240SBarry Smith   }
6733a40ed3dSBarry Smith   PetscFunctionReturn(0);
6745e76c431SBarry Smith }
67504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
676fb2e594dSBarry Smith EXTERN_C_BEGIN
67782bf6240SBarry Smith #undef __FUNC__
67882bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
679af2d14d2SLois Curfman McInnes int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,
680af2d14d2SLois Curfman McInnes                          double,double*,double*,int*),void *lsctx)
68182bf6240SBarry Smith {
68282bf6240SBarry Smith   PetscFunctionBegin;
68382bf6240SBarry Smith   ((SNES_LS *)(snes->data))->LineSearch = func;
684af2d14d2SLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP = lsctx;
68582bf6240SBarry Smith   PetscFunctionReturn(0);
68682bf6240SBarry Smith }
687fb2e594dSBarry Smith EXTERN_C_END
68804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
689*c8dd0c0dSLois Curfman McInnes #undef __FUNC__
690*c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck"
691*c8dd0c0dSLois Curfman McInnes /*@C
692*c8dd0c0dSLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the step computed by the
693*c8dd0c0dSLois Curfman McInnes    line search routine in the Newton-based method SNES_EQ_LS.
694*c8dd0c0dSLois Curfman McInnes 
695*c8dd0c0dSLois Curfman McInnes    Input Parameters:
696*c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
697*c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
698*c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
699*c8dd0c0dSLois Curfman McInnes 
700*c8dd0c0dSLois Curfman McInnes    Collective on SNES
701*c8dd0c0dSLois Curfman McInnes 
702*c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
703*c8dd0c0dSLois Curfman McInnes .vb
704*c8dd0c0dSLois Curfman McInnes    func (SNES snes, void *lsctx, Vec x, int *flag)
705*c8dd0c0dSLois Curfman McInnes .ve
706*c8dd0c0dSLois Curfman McInnes 
707*c8dd0c0dSLois Curfman McInnes     Input parameters for func:
708*c8dd0c0dSLois Curfman McInnes +   snes - nonlinear context
709*c8dd0c0dSLois Curfman McInnes .   checkctx - optional user-defined context for use by step checking routine
710*c8dd0c0dSLois Curfman McInnes -   x - current iterate
711*c8dd0c0dSLois Curfman McInnes 
712*c8dd0c0dSLois Curfman McInnes     Output parameters for func:
713*c8dd0c0dSLois Curfman McInnes +   x - current iterate (possibly modified)
714*c8dd0c0dSLois Curfman McInnes -   flag - flag indicating whether x has been modified (either
715*c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
716*c8dd0c0dSLois Curfman McInnes 
717*c8dd0c0dSLois Curfman McInnes     Level: advanced
718*c8dd0c0dSLois Curfman McInnes 
719*c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
720*c8dd0c0dSLois Curfman McInnes 
721*c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
722*c8dd0c0dSLois Curfman McInnes @*/
723*c8dd0c0dSLois Curfman McInnes int SNESSetLineSearchCheck(SNES snes,int (*func)(SNES,void*,Vec,int*),void *checkctx)
724*c8dd0c0dSLois Curfman McInnes {
725*c8dd0c0dSLois Curfman McInnes   int ierr, (*f)(SNES,int (*)(SNES,void*,Vec,int*),void*);
726*c8dd0c0dSLois Curfman McInnes 
727*c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
728*c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void **)&f);CHKERRQ(ierr);
729*c8dd0c0dSLois Curfman McInnes   if (f) {
730*c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
731*c8dd0c0dSLois Curfman McInnes   }
732*c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
733*c8dd0c0dSLois Curfman McInnes }
734*c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
735*c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
736*c8dd0c0dSLois Curfman McInnes #undef __FUNC__
737*c8dd0c0dSLois Curfman McInnes #define __FUNC__ "SNESSetLineSearchCheck_LS"
738*c8dd0c0dSLois Curfman McInnes int SNESSetLineSearchCheck_LS(SNES snes,int (*func)(SNES,void*,Vec,int*),void *checkctx)
739*c8dd0c0dSLois Curfman McInnes {
740*c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
741*c8dd0c0dSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
742*c8dd0c0dSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP = checkctx;
743*c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
744*c8dd0c0dSLois Curfman McInnes }
745*c8dd0c0dSLois Curfman McInnes EXTERN_C_END
746*c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
74704d965bbSLois Curfman McInnes /*
74804d965bbSLois Curfman McInnes    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
74982bf6240SBarry Smith 
75004d965bbSLois Curfman McInnes    Input Parameter:
75104d965bbSLois Curfman McInnes .  snes - the SNES context
75204d965bbSLois Curfman McInnes 
75304d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
75404d965bbSLois Curfman McInnes */
7555615d1e5SSatish Balay #undef __FUNC__
756d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
757f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
7585e42470aSBarry Smith {
7595e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
7606daaf66cSBarry Smith 
7613a40ed3dSBarry Smith   PetscFunctionBegin;
76276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
763af2d14d2SLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
76476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
76576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
76676be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
7673a40ed3dSBarry Smith   PetscFunctionReturn(0);
7685e42470aSBarry Smith }
76904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
77004d965bbSLois Curfman McInnes /*
771de49b36dSLois Curfman McInnes    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
77204d965bbSLois Curfman McInnes 
77304d965bbSLois Curfman McInnes    Input Parameters:
77404d965bbSLois Curfman McInnes .  SNES - the SNES context
77504d965bbSLois Curfman McInnes .  viewer - visualization context
77604d965bbSLois Curfman McInnes 
77704d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
77804d965bbSLois Curfman McInnes */
7795615d1e5SSatish Balay #undef __FUNC__
780d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
781e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
782a935fc98SLois Curfman McInnes {
783a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
78419bcc07fSBarry Smith   char       *cstr;
78551695f54SLois Curfman McInnes   int        ierr;
78619bcc07fSBarry Smith   ViewerType vtype;
787a935fc98SLois Curfman McInnes 
7883a40ed3dSBarry Smith   PetscFunctionBegin;
78919bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
7903f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
79119bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
79219bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
79319bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
79419bcc07fSBarry Smith     else                                                cstr = "unknown";
7950ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);
7960ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);
7975cd90555SBarry Smith   } else {
7985cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
79919bcc07fSBarry Smith   }
8003a40ed3dSBarry Smith   PetscFunctionReturn(0);
801a935fc98SLois Curfman McInnes }
80204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
80304d965bbSLois Curfman McInnes /*
80404d965bbSLois Curfman McInnes    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
80504d965bbSLois Curfman McInnes 
80604d965bbSLois Curfman McInnes    Input Parameter:
80704d965bbSLois Curfman McInnes .  snes - the SNES context
80804d965bbSLois Curfman McInnes 
80904d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
81004d965bbSLois Curfman McInnes */
8115615d1e5SSatish Balay #undef __FUNC__
8125615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
813f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
8145e42470aSBarry Smith {
8155e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
8165e42470aSBarry Smith   char    ver[16];
8175e42470aSBarry Smith   double  tmp;
81825cdf11fSBarry Smith   int     ierr,flg;
8195e42470aSBarry Smith 
8203a40ed3dSBarry Smith   PetscFunctionBegin;
82109d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
82225cdf11fSBarry Smith   if (flg) {
8235e42470aSBarry Smith     ls->alpha = tmp;
8245e42470aSBarry Smith   }
82509d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
82625cdf11fSBarry Smith   if (flg) {
8275e42470aSBarry Smith     ls->maxstep = tmp;
8285e42470aSBarry Smith   }
82909d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
83025cdf11fSBarry Smith   if (flg) {
8315e42470aSBarry Smith     ls->steptol = tmp;
8325e42470aSBarry Smith   }
83309d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
83425cdf11fSBarry Smith   if (flg) {
83548d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
836af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
837b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"basicnonorms")) {
838af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
839b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"quadratic")) {
840af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
841b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"cubic")) {
842af2d14d2SLois Curfman McInnes       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
8435e42470aSBarry Smith     }
844a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
8455e42470aSBarry Smith   }
8463a40ed3dSBarry Smith   PetscFunctionReturn(0);
8475e42470aSBarry Smith }
84804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
84904d965bbSLois Curfman McInnes /*
85004d965bbSLois Curfman McInnes    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
85104d965bbSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
85204d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
85304d965bbSLois Curfman McInnes 
85404d965bbSLois Curfman McInnes 
85504d965bbSLois Curfman McInnes    Input Parameter:
85604d965bbSLois Curfman McInnes .  snes - the SNES context
85704d965bbSLois Curfman McInnes 
85804d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
85904d965bbSLois Curfman McInnes  */
860fb2e594dSBarry Smith EXTERN_C_BEGIN
8615615d1e5SSatish Balay #undef __FUNC__
8625615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
863f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
8645e42470aSBarry Smith {
86582bf6240SBarry Smith   int     ierr;
8665e42470aSBarry Smith   SNES_LS *neP;
8675e42470aSBarry Smith 
8683a40ed3dSBarry Smith   PetscFunctionBegin;
869a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
870a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
871a8c6a408SBarry Smith   }
87282bf6240SBarry Smith 
873f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
874f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
875f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
876f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
877f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
878f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
879f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
8805baf8537SBarry Smith   snes->nwork           = 0;
8815e42470aSBarry Smith 
8820452661fSBarry Smith   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
883ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
8845e42470aSBarry Smith   snes->data    	= (void *) neP;
8855e42470aSBarry Smith   neP->alpha		= 1.e-4;
8865e42470aSBarry Smith   neP->maxstep		= 1.e8;
8875e42470aSBarry Smith   neP->steptol		= 1.e-12;
8885e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
889*c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
890*c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
891*c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
89282bf6240SBarry Smith 
89394d884c6SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
894e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
895*c8dd0c0dSLois Curfman McInnes   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",
896*c8dd0c0dSLois Curfman McInnes                     (void*)SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
89782bf6240SBarry Smith 
8983a40ed3dSBarry Smith   PetscFunctionReturn(0);
8995e42470aSBarry Smith }
900fb2e594dSBarry Smith EXTERN_C_END
90182bf6240SBarry Smith 
90282bf6240SBarry Smith 
90382bf6240SBarry Smith 
904