xref: /petsc/src/snes/impls/ls/ls.c (revision 04d965bb56167bd2a0a34aea0cbc0c7064553e9a)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*04d965bbSLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.105 1998/04/13 17:56:04 bsmith Exp curfman $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
55e76c431SBarry Smith #include <math.h>
670f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
7b16a3bb1SBarry Smith #include "pinclude/pviewer.h"
85e42470aSBarry Smith 
9*04d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
10*04d965bbSLois Curfman McInnes 
11*04d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
12*04d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
13*04d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
14*04d965bbSLois Curfman McInnes      respectively.
15*04d965bbSLois Curfman McInnes 
16*04d965bbSLois Curfman McInnes      The following basic routines are required for each nonlinear solver
17*04d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
18*04d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
19*04d965bbSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear solver
20*04d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
21*04d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
22*04d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
23*04d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
24*04d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
25*04d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
26*04d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
27*04d965bbSLois Curfman McInnes      for all nonlinear solvers.
28*04d965bbSLois Curfman McInnes 
29*04d965bbSLois Curfman McInnes      Another key routine is:
30*04d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
31*04d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
32*04d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
33*04d965bbSLois Curfman McInnes      SNESSolve() if necessary.
34*04d965bbSLois Curfman McInnes 
35*04d965bbSLois Curfman McInnes      Additional basic routines are:
36*04d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
37*04d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
38*04d965bbSLois Curfman McInnes                                     have actually been used.
39*04d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
40*04d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
41*04d965bbSLois Curfman McInnes 
42*04d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
43*04d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
44*04d965bbSLois Curfman McInnes      above description applies to these categories also.
45*04d965bbSLois Curfman McInnes 
46*04d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
475e42470aSBarry Smith /*
48*04d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
49*04d965bbSLois Curfman McInnes    method with a line search.
505e76c431SBarry Smith 
51*04d965bbSLois Curfman McInnes    Input Parameters:
52*04d965bbSLois Curfman McInnes .  snes - the SNES context
53*04d965bbSLois Curfman McInnes .  x - the solution vector
545e76c431SBarry Smith 
55*04d965bbSLois Curfman McInnes    Output Parameter:
56*04d965bbSLois Curfman McInnes .  its - number of iterations until termination
57*04d965bbSLois Curfman McInnes 
58*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
595e76c431SBarry Smith 
605e76c431SBarry Smith    Notes:
615e76c431SBarry Smith    This implements essentially a truncated Newton method with a
625e76c431SBarry Smith    line search.  By default a cubic backtracking line search
635e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
645e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
65393d2d9aSLois Curfman McInnes    and Schnabel.
665e76c431SBarry Smith */
675615d1e5SSatish Balay #undef __FUNC__
685615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS"
69f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
705e76c431SBarry Smith {
715e42470aSBarry Smith   SNES_LS       *neP = (SNES_LS *) snes->data;
72761aaf1bSLois Curfman McInnes   int           maxits, i, history_len, ierr, lits, lsfail;
73112a2221SBarry Smith   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
746b5873e3SBarry Smith   double        fnorm, gnorm, xnorm, ynorm, *history;
755e42470aSBarry Smith   Vec           Y, X, F, G, W, TMP;
765e76c431SBarry Smith 
773a40ed3dSBarry Smith   PetscFunctionBegin;
785e42470aSBarry Smith   history	= snes->conv_hist;	/* convergence history */
7951979daaSLois Curfman McInnes   history_len	= snes->conv_hist_size;	/* convergence history length */
805e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
815e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
8239e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
835e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
845e42470aSBarry Smith   G		= snes->work[1];
855e42470aSBarry Smith   W		= snes->work[2];
865e76c431SBarry Smith 
8725ed9cd1SBarry Smith   snes->iter = 0;
885334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
89cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
905e42470aSBarry Smith   snes->norm = fnorm;
9151979daaSLois Curfman McInnes   if (history) history[0] = fnorm;
9294a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
935e76c431SBarry Smith 
943a40ed3dSBarry Smith   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
9594a9d846SBarry Smith 
96d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
97d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
98d034289bSBarry Smith 
995e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
1005e42470aSBarry Smith     snes->iter = i+1;
1015e76c431SBarry Smith 
102ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1035334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
1045334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
10578b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
1067a00f4a9SLois Curfman McInnes     snes->linear_its += PetscAbsInt(lits);
107af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
108ea4d3ed3SLois Curfman McInnes 
109ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
110ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
111ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
112ea4d3ed3SLois Curfman McInnes     */
11381b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
114ddd12767SBarry Smith     ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
115af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
116761aaf1bSLois Curfman McInnes     if (lsfail) snes->nfailures++;
1175e76c431SBarry Smith 
11839e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
11939e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1205e76c431SBarry Smith     fnorm = gnorm;
1215e76c431SBarry Smith 
1225e42470aSBarry Smith     snes->norm = fnorm;
1235e76c431SBarry Smith     if (history && history_len > i+1) history[i+1] = fnorm;
12494a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1255e76c431SBarry Smith 
1265e76c431SBarry Smith     /* Test for convergence */
12729e0b56fSBarry Smith     if (snes->converged) {
12829e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
129bbb6d6a8SBarry Smith       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
13016c95cacSBarry Smith         break;
13116c95cacSBarry Smith       }
13216c95cacSBarry Smith     }
13329e0b56fSBarry Smith   }
13439e2f89bSBarry Smith   if (X != snes->vec_sol) {
135393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
13639e2f89bSBarry Smith     snes->vec_sol_always  = snes->vec_sol;
13739e2f89bSBarry Smith     snes->vec_func_always = snes->vec_func;
13839e2f89bSBarry Smith   }
13952392280SLois Curfman McInnes   if (i == maxits) {
140981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
14152392280SLois Curfman McInnes     i--;
14252392280SLois Curfman McInnes   }
14351979daaSLois Curfman McInnes   if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1;
1445e42470aSBarry Smith   *outits = i+1;
1453a40ed3dSBarry Smith   PetscFunctionReturn(0);
1465e76c431SBarry Smith }
147*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
148*04d965bbSLois Curfman McInnes /*
149*04d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
150*04d965bbSLois Curfman McInnes    of the SNES_EQ_LS nonlinear solver.
151*04d965bbSLois Curfman McInnes 
152*04d965bbSLois Curfman McInnes    Input Parameter:
153*04d965bbSLois Curfman McInnes .  snes - the SNES context
154*04d965bbSLois Curfman McInnes .  x - the solution vector
155*04d965bbSLois Curfman McInnes 
156*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
157*04d965bbSLois Curfman McInnes 
158*04d965bbSLois Curfman McInnes    Notes:
159*04d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
160*04d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
161*04d965bbSLois Curfman McInnes    the call to SNESSolve().
162*04d965bbSLois Curfman McInnes  */
1635615d1e5SSatish Balay #undef __FUNC__
1645615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
165f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
1665e76c431SBarry Smith {
1675e42470aSBarry Smith   int ierr;
1683a40ed3dSBarry Smith 
1693a40ed3dSBarry Smith   PetscFunctionBegin;
17081b6cf68SLois Curfman McInnes   snes->nwork = 4;
171d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1725e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
17381b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1743a40ed3dSBarry Smith   PetscFunctionReturn(0);
1755e76c431SBarry Smith }
176*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
177*04d965bbSLois Curfman McInnes /*
178*04d965bbSLois Curfman McInnes    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
179*04d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
180*04d965bbSLois Curfman McInnes 
181*04d965bbSLois Curfman McInnes    Input Parameter:
182*04d965bbSLois Curfman McInnes .  snes - the SNES context
183*04d965bbSLois Curfman McInnes 
184*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetDestroy()
185*04d965bbSLois Curfman McInnes  */
1865615d1e5SSatish Balay #undef __FUNC__
187d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
188e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
1895e76c431SBarry Smith {
190393d2d9aSLois Curfman McInnes   int  ierr;
1913a40ed3dSBarry Smith 
1923a40ed3dSBarry Smith   PetscFunctionBegin;
1935baf8537SBarry Smith   if (snes->nwork) {
1944b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
19521c89e3eSBarry Smith   }
1960452661fSBarry Smith   PetscFree(snes->data);
1973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1985e76c431SBarry Smith }
199*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2005615d1e5SSatish Balay #undef __FUNC__
2015615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
20282bf6240SBarry Smith 
2034b828684SBarry Smith /*@C
2045e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2055e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2065e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2075e76c431SBarry Smith 
2085e76c431SBarry Smith    Input Parameters:
2095e42470aSBarry Smith .  snes - nonlinear context
2105e76c431SBarry Smith .  x - current iterate
2115e76c431SBarry Smith .  f - residual evaluated at x
2125e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2135e76c431SBarry Smith .  w - work vector
2145e76c431SBarry Smith .  fnorm - 2-norm of f
2155e76c431SBarry Smith 
216c4a48953SLois Curfman McInnes    Output Parameters:
2175e76c431SBarry Smith .  g - residual evaluated at new iterate y
2185e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2195e76c431SBarry Smith .  gnorm - 2-norm of g
2205e76c431SBarry Smith .  ynorm - 2-norm of search length
221761aaf1bSLois Curfman McInnes .  flag - set to 0, indicating a successful line search
2225e76c431SBarry Smith 
223fee21e36SBarry Smith    Collective on SNES and Vec
224fee21e36SBarry Smith 
225c4a48953SLois Curfman McInnes    Options Database Key:
22609d61ba7SLois Curfman McInnes $  -snes_eq_ls basic
227c4a48953SLois Curfman McInnes 
22828ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
22928ae5a21SLois Curfman McInnes 
230f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
23177c4ece6SBarry Smith           SNESSetLineSearch()
2325e76c431SBarry Smith @*/
2335e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
234761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
2355e76c431SBarry Smith {
2365e42470aSBarry Smith   int    ierr;
2375334005bSBarry Smith   Scalar mone = -1.0;
2385334005bSBarry Smith 
2393a40ed3dSBarry Smith   PetscFunctionBegin;
240761aaf1bSLois Curfman McInnes   *flag = 0;
2417857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
242cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
243ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
244ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
245cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
2467857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2473a40ed3dSBarry Smith   PetscFunctionReturn(0);
2485e76c431SBarry Smith }
249*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2505615d1e5SSatish Balay #undef __FUNC__
25129e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
25282bf6240SBarry Smith 
25329e0b56fSBarry Smith /*@C
25429e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
25529e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
25629e0b56fSBarry Smith    even compute the norm of the function or search direction; this
25729e0b56fSBarry Smith    is intended only when you know the full step is fine and are
25829e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
25929e0b56fSBarry Smith    example, you are running always for a fixed number of Newton
26029e0b56fSBarry Smith    steps).
26129e0b56fSBarry Smith 
26229e0b56fSBarry Smith    Input Parameters:
26329e0b56fSBarry Smith .  snes - nonlinear context
26429e0b56fSBarry Smith .  x - current iterate
26529e0b56fSBarry Smith .  f - residual evaluated at x
26629e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
26729e0b56fSBarry Smith .  w - work vector
26829e0b56fSBarry Smith .  fnorm - 2-norm of f
26929e0b56fSBarry Smith 
27029e0b56fSBarry Smith    Output Parameters:
27129e0b56fSBarry Smith .  g - residual evaluated at new iterate y
27229e0b56fSBarry Smith .  gnorm - not changed
27329e0b56fSBarry Smith .  ynorm - not changed
27429e0b56fSBarry Smith .  flag - set to 0, indicating a successful line search
27529e0b56fSBarry Smith 
276fee21e36SBarry Smith    Collective on SNES and Vec
277fee21e36SBarry Smith 
27829e0b56fSBarry Smith    Options Database Key:
27929e0b56fSBarry Smith $  -snes_eq_ls basicnonorms
28029e0b56fSBarry Smith 
28129e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
28229e0b56fSBarry Smith 
28329e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
28429e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
28529e0b56fSBarry Smith @*/
28629e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
28729e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
28829e0b56fSBarry Smith {
28929e0b56fSBarry Smith   int    ierr;
29029e0b56fSBarry Smith   Scalar mone = -1.0;
29129e0b56fSBarry Smith 
2923a40ed3dSBarry Smith   PetscFunctionBegin;
29329e0b56fSBarry Smith   *flag = 0;
29429e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
29529e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
29629e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
29729e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2983a40ed3dSBarry Smith   PetscFunctionReturn(0);
29929e0b56fSBarry Smith }
300*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30129e0b56fSBarry Smith #undef __FUNC__
3025615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3034b828684SBarry Smith /*@C
304f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3055e76c431SBarry Smith 
3065e76c431SBarry Smith    Input Parameters:
3075e42470aSBarry Smith .  snes - nonlinear context
3085e76c431SBarry Smith .  x - current iterate
3095e76c431SBarry Smith .  f - residual evaluated at x
3105e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3115e76c431SBarry Smith .  w - work vector
3125e76c431SBarry Smith .  fnorm - 2-norm of f
3135e76c431SBarry Smith 
314393d2d9aSLois Curfman McInnes    Output Parameters:
3155e76c431SBarry Smith .  g - residual evaluated at new iterate y
3165e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3175e76c431SBarry Smith .  gnorm - 2-norm of g
3185e76c431SBarry Smith .  ynorm - 2-norm of search length
319761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
3205e76c431SBarry Smith 
321fee21e36SBarry Smith    Collective on SNES
322fee21e36SBarry Smith 
323c4a48953SLois Curfman McInnes    Options Database Key:
32409d61ba7SLois Curfman McInnes $  -snes_eq_ls cubic
325c4a48953SLois Curfman McInnes 
3265e76c431SBarry Smith    Notes:
3275e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3285e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3295e76c431SBarry Smith 
33028ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33128ae5a21SLois Curfman McInnes 
33277c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
3335e76c431SBarry Smith @*/
3345e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
335761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3365e76c431SBarry Smith {
337406556e6SLois Curfman McInnes   /*
338406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
339406556e6SLois Curfman McInnes      minimization problem:
340406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
341406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
342406556e6SLois Curfman McInnes    */
343406556e6SLois Curfman McInnes 
344ddd12767SBarry Smith   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
345ea4d3ed3SLois Curfman McInnes   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
3463a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
3475e42470aSBarry Smith   Scalar  cinitslope, clambda;
3486b5873e3SBarry Smith #endif
3495e42470aSBarry Smith   int     ierr, count;
3505e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
3515334005bSBarry Smith   Scalar  mone = -1.0,scale;
3525e76c431SBarry Smith 
3533a40ed3dSBarry Smith   PetscFunctionBegin;
3547857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
355761aaf1bSLois Curfman McInnes   *flag   = 0;
3565e76c431SBarry Smith   alpha   = neP->alpha;
3575e76c431SBarry Smith   maxstep = neP->maxstep;
3585e76c431SBarry Smith   steptol = neP->steptol;
3595e76c431SBarry Smith 
360cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
361a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
362a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
363a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
364a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
365a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
366ad922ac9SBarry Smith     goto theend1;
36794a9d846SBarry Smith   }
3685e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3695e42470aSBarry Smith     scale = maxstep/(*ynorm);
3703a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
37194a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",real(scale));
3726b5873e3SBarry Smith #else
37394a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
3746b5873e3SBarry Smith #endif
375393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
3765e76c431SBarry Smith     *ynorm = maxstep;
3775e76c431SBarry Smith   }
3785e76c431SBarry Smith   minlambda = steptol/(*ynorm);
379a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
3803a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
381a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
3825e42470aSBarry Smith   initslope = real(cinitslope);
3835e42470aSBarry Smith #else
384a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
3855e42470aSBarry Smith #endif
3865e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3875e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3885e76c431SBarry Smith 
389393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
3905334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
39178b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
392cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);
393406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
394393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
39594a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
396ad922ac9SBarry Smith     goto theend1;
3975e76c431SBarry Smith   }
3985e76c431SBarry Smith 
3995e76c431SBarry Smith   /* Fit points with quadratic */
4005e76c431SBarry Smith   lambda = 1.0; count = 0;
401a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4025e76c431SBarry Smith   lambdaprev = lambda;
4035e76c431SBarry Smith   gnormprev = *gnorm;
404ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
405ddd12767SBarry Smith   else lambda = lambdatemp;
406393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w); CHKERRQ(ierr);
407ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
4083a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
409ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4105e42470aSBarry Smith #else
411ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4125e42470aSBarry Smith #endif
41378b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
414cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
415406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
416393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
417ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
418ad922ac9SBarry Smith     goto theend1;
4195e76c431SBarry Smith   }
4205e76c431SBarry Smith 
4215e76c431SBarry Smith   /* Fit points with cubic */
4225e76c431SBarry Smith   count = 1;
4235e76c431SBarry Smith   while (1) {
4245e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4252b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4262b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
427ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
428393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
429761aaf1bSLois Curfman McInnes       *flag = -1; break;
4305e76c431SBarry Smith     }
431406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
432406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
433ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4342b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4355e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4365e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4375e76c431SBarry Smith     if (a == 0.0) {
4385e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4395e76c431SBarry Smith     } else {
4405e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4415e76c431SBarry Smith     }
4425e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4435e76c431SBarry Smith       lambdatemp = .5*lambda;
4445e76c431SBarry Smith     }
4455e76c431SBarry Smith     lambdaprev = lambda;
4465e76c431SBarry Smith     gnormprev = *gnorm;
4475e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4485e76c431SBarry Smith       lambda = .1*lambda;
4495e76c431SBarry Smith     }
4505e76c431SBarry Smith     else lambda = lambdatemp;
451393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w ); CHKERRQ(ierr);
452ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
4533a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
454ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
455393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4565e42470aSBarry Smith #else
457ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4585e42470aSBarry Smith #endif
45978b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
460cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
461406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
462393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
463ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
4642715565aSLois Curfman McInnes       goto theend1;
4652b022350SLois Curfman McInnes     } else {
4662b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
4675e76c431SBarry Smith     }
4685e76c431SBarry Smith     count++;
4695e76c431SBarry Smith   }
470ad922ac9SBarry Smith   theend1:
4717857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4723a40ed3dSBarry Smith   PetscFunctionReturn(0);
4735e76c431SBarry Smith }
474*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4755615d1e5SSatish Balay #undef __FUNC__
4765615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
4774b828684SBarry Smith /*@C
478f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
4795e76c431SBarry Smith 
4805e42470aSBarry Smith    Input Parameters:
481f59ffdedSLois Curfman McInnes .  snes - the SNES context
4825e42470aSBarry Smith .  x - current iterate
4835e42470aSBarry Smith .  f - residual evaluated at x
4845e42470aSBarry Smith .  y - search direction (contains new iterate on output)
4855e42470aSBarry Smith .  w - work vector
4865e42470aSBarry Smith .  fnorm - 2-norm of f
4875e42470aSBarry Smith 
488c4a48953SLois Curfman McInnes    Output Parameters:
4895e42470aSBarry Smith .  g - residual evaluated at new iterate y
4905e42470aSBarry Smith .  y - new iterate (contains search direction on input)
4915e42470aSBarry Smith .  gnorm - 2-norm of g
4925e42470aSBarry Smith .  ynorm - 2-norm of search length
493761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
4945e42470aSBarry Smith 
495fee21e36SBarry Smith    Collective on SNES and Vec
496fee21e36SBarry Smith 
497c4a48953SLois Curfman McInnes    Options Database Key:
49809d61ba7SLois Curfman McInnes $  -snes_eq_ls quadratic
499c4a48953SLois Curfman McInnes 
5005e42470aSBarry Smith    Notes:
50177c4ece6SBarry Smith    Use SNESSetLineSearch()
502f63b844aSLois Curfman McInnes    to set this routine within the SNES_EQ_LS method.
5035e42470aSBarry Smith 
504f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
505f59ffdedSLois Curfman McInnes 
50677c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
5075e42470aSBarry Smith @*/
5085e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
509761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5105e76c431SBarry Smith {
511406556e6SLois Curfman McInnes   /*
512406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
513406556e6SLois Curfman McInnes      minimization problem:
514406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
515406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
516406556e6SLois Curfman McInnes    */
51740011551SBarry Smith   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
5183a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
5195e42470aSBarry Smith   Scalar  cinitslope,clambda;
5206b5873e3SBarry Smith #endif
5215e42470aSBarry Smith   int     ierr,count;
5225e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
5235334005bSBarry Smith   Scalar  mone = -1.0,scale;
5245e76c431SBarry Smith 
5253a40ed3dSBarry Smith   PetscFunctionBegin;
5267857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
527761aaf1bSLois Curfman McInnes   *flag   = 0;
5285e42470aSBarry Smith   alpha   = neP->alpha;
5295e42470aSBarry Smith   maxstep = neP->maxstep;
5305e42470aSBarry Smith   steptol = neP->steptol;
5315e76c431SBarry Smith 
532cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
533b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
53494a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
535b37302c6SLois Curfman McInnes     *gnorm = fnorm;
536b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
537b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
538ad922ac9SBarry Smith     goto theend2;
53994a9d846SBarry Smith   }
5405e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5415e42470aSBarry Smith     scale = maxstep/(*ynorm);
542393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
5435e42470aSBarry Smith     *ynorm = maxstep;
5445e76c431SBarry Smith   }
5455e42470aSBarry Smith   minlambda = steptol/(*ynorm);
546a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
5473a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
548a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
5495e42470aSBarry Smith   initslope = real(cinitslope);
5505e42470aSBarry Smith #else
551a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
5525e42470aSBarry Smith #endif
5535e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
5545e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
5555e42470aSBarry Smith 
556393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
5575334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
55878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
559cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
560406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
561393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
56294a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
563ad922ac9SBarry Smith     goto theend2;
5645e42470aSBarry Smith   }
5655e42470aSBarry Smith 
5665e42470aSBarry Smith   /* Fit points with quadratic */
5675e42470aSBarry Smith   lambda = 1.0; count = 0;
5685e42470aSBarry Smith   count = 1;
5695e42470aSBarry Smith   while (1) {
5705e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
571981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
572981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
573ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
574393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
575761aaf1bSLois Curfman McInnes       *flag = -1; break;
5765e42470aSBarry Smith     }
577a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5785e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
5795e42470aSBarry Smith       lambda = .1*lambda;
5805e42470aSBarry Smith     } else lambda = lambdatemp;
581393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w); CHKERRQ(ierr);
5825334005bSBarry Smith     lambda = -lambda;
5833a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
584393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
5855e42470aSBarry Smith #else
586393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
5875e42470aSBarry Smith #endif
58878b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
589cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
590406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
591393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
592981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
59306259719SBarry Smith       break;
5945e42470aSBarry Smith     }
5955e42470aSBarry Smith     count++;
5965e42470aSBarry Smith   }
597ad922ac9SBarry Smith   theend2:
5987857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5993a40ed3dSBarry Smith   PetscFunctionReturn(0);
6005e76c431SBarry Smith }
601*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6025615d1e5SSatish Balay #undef __FUNC__
603d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
604c9e6a524SLois Curfman McInnes /*@C
60577c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
606f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
6075e76c431SBarry Smith 
6085e76c431SBarry Smith    Input Parameters:
609eafb4bcbSBarry Smith .  snes - nonlinear context obtained from SNESCreate()
6105e76c431SBarry Smith .  func - pointer to int function
6115e76c431SBarry Smith 
612fee21e36SBarry Smith    Collective on SNES
613fee21e36SBarry Smith 
614c4a48953SLois Curfman McInnes    Available Routines:
615f59ffdedSLois Curfman McInnes .  SNESCubicLineSearch() - default line search
616f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
617f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
61882bf6240SBarry Smith .  SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel)
6195e76c431SBarry Smith 
620c4a48953SLois Curfman McInnes     Options Database Keys:
62109d61ba7SLois Curfman McInnes $   -snes_eq_ls [basic,quadratic,cubic]
622912ebf9aSLois Curfman McInnes $   -snes_eq_ls_alpha <alpha>
623912ebf9aSLois Curfman McInnes $   -snes_eq_ls_maxstep <max>
624912ebf9aSLois Curfman McInnes $   -snes_eq_ls_steptol <steptol>
625c4a48953SLois Curfman McInnes 
6265e76c431SBarry Smith    Calling sequence of func:
627f59ffdedSLois Curfman McInnes    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
628761aaf1bSLois Curfman McInnes          Vec w, double fnorm, double *ynorm,
629761aaf1bSLois Curfman McInnes          double *gnorm, *flag)
6305e76c431SBarry Smith 
6315e76c431SBarry Smith     Input parameters for func:
6325e42470aSBarry Smith .   snes - nonlinear context
6335e76c431SBarry Smith .   x - current iterate
6345e76c431SBarry Smith .   f - residual evaluated at x
6355e76c431SBarry Smith .   y - search direction (contains new iterate on output)
6365e76c431SBarry Smith .   w - work vector
6375e76c431SBarry Smith .   fnorm - 2-norm of f
6385e76c431SBarry Smith 
6395e76c431SBarry Smith     Output parameters for func:
6405e76c431SBarry Smith .   g - residual evaluated at new iterate y
6415e76c431SBarry Smith .   y - new iterate (contains search direction on input)
6425e76c431SBarry Smith .   gnorm - 2-norm of g
6435e76c431SBarry Smith .   ynorm - 2-norm of search length
644761aaf1bSLois Curfman McInnes .   flag - set to 0 if the line search succeeds; a nonzero integer
645761aaf1bSLois Curfman McInnes            on failure.
646f59ffdedSLois Curfman McInnes 
647f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
648f59ffdedSLois Curfman McInnes 
649f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
6505e76c431SBarry Smith @*/
65177c4ece6SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
652761aaf1bSLois Curfman McInnes                              double,double*,double*,int*))
6535e76c431SBarry Smith {
65482bf6240SBarry Smith   int ierr, (*f)(SNES,int (*f)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*));
65582bf6240SBarry Smith 
6563a40ed3dSBarry Smith   PetscFunctionBegin;
657e1311b90SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch",(void **)&f);CHKERRQ(ierr);
65882bf6240SBarry Smith   if (f) {
65982bf6240SBarry Smith     ierr = (*f)(snes,func);CHKERRQ(ierr);
66082bf6240SBarry Smith   }
6613a40ed3dSBarry Smith   PetscFunctionReturn(0);
6625e76c431SBarry Smith }
663*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
66482bf6240SBarry Smith #undef __FUNC__
66582bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
66682bf6240SBarry Smith int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
66782bf6240SBarry Smith                          double,double*,double*,int*))
66882bf6240SBarry Smith {
66982bf6240SBarry Smith   PetscFunctionBegin;
67082bf6240SBarry Smith   ((SNES_LS *)(snes->data))->LineSearch = func;
67182bf6240SBarry Smith   PetscFunctionReturn(0);
67282bf6240SBarry Smith }
673*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
674*04d965bbSLois Curfman McInnes /*
675*04d965bbSLois Curfman McInnes    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
67682bf6240SBarry Smith 
677*04d965bbSLois Curfman McInnes    Input Parameter:
678*04d965bbSLois Curfman McInnes .  snes - the SNES context
679*04d965bbSLois Curfman McInnes 
680*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
681*04d965bbSLois Curfman McInnes */
682*04d965bbSLois Curfman McInnes 
6835615d1e5SSatish Balay #undef __FUNC__
684d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
685f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
6865e42470aSBarry Smith {
6875e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
6886daaf66cSBarry Smith 
6893a40ed3dSBarry Smith   PetscFunctionBegin;
69076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
69176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
69276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
69376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
69476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
6953a40ed3dSBarry Smith   PetscFunctionReturn(0);
6965e42470aSBarry Smith }
697*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
698*04d965bbSLois Curfman McInnes /*
699*04d965bbSLois Curfman McInnes    SNESView_EQ_LS - Prints the SNES_EQ_LS data structure.
700*04d965bbSLois Curfman McInnes 
701*04d965bbSLois Curfman McInnes    Input Parameters:
702*04d965bbSLois Curfman McInnes .  SNES - the SNES context
703*04d965bbSLois Curfman McInnes .  viewer - visualization context
704*04d965bbSLois Curfman McInnes 
705*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
706*04d965bbSLois Curfman McInnes */
7075615d1e5SSatish Balay #undef __FUNC__
708d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
709e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
710a935fc98SLois Curfman McInnes {
711a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
712d636dbe3SBarry Smith   FILE       *fd;
71319bcc07fSBarry Smith   char       *cstr;
71451695f54SLois Curfman McInnes   int        ierr;
71519bcc07fSBarry Smith   ViewerType vtype;
716a935fc98SLois Curfman McInnes 
7173a40ed3dSBarry Smith   PetscFunctionBegin;
71819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
71919bcc07fSBarry Smith   if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
72090ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
72119bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
72219bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
72319bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
72419bcc07fSBarry Smith     else cstr = "unknown";
72577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"    line search variant: %s\n",cstr);
72677c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"    alpha=%g, maxstep=%g, steptol=%g\n",
727a935fc98SLois Curfman McInnes                  ls->alpha,ls->maxstep,ls->steptol);
7285cd90555SBarry Smith   } else {
7295cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
73019bcc07fSBarry Smith   }
7313a40ed3dSBarry Smith   PetscFunctionReturn(0);
732a935fc98SLois Curfman McInnes }
733*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
734*04d965bbSLois Curfman McInnes /*
735*04d965bbSLois Curfman McInnes    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
736*04d965bbSLois Curfman McInnes 
737*04d965bbSLois Curfman McInnes    Input Parameter:
738*04d965bbSLois Curfman McInnes .  snes - the SNES context
739*04d965bbSLois Curfman McInnes 
740*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
741*04d965bbSLois Curfman McInnes */
7425615d1e5SSatish Balay #undef __FUNC__
7435615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
744f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
7455e42470aSBarry Smith {
7465e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
7475e42470aSBarry Smith   char    ver[16];
7485e42470aSBarry Smith   double  tmp;
74925cdf11fSBarry Smith   int     ierr,flg;
7505e42470aSBarry Smith 
7513a40ed3dSBarry Smith   PetscFunctionBegin;
75209d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
75325cdf11fSBarry Smith   if (flg) {
7545e42470aSBarry Smith     ls->alpha = tmp;
7555e42470aSBarry Smith   }
75609d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
75725cdf11fSBarry Smith   if (flg) {
7585e42470aSBarry Smith     ls->maxstep = tmp;
7595e42470aSBarry Smith   }
76009d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
76125cdf11fSBarry Smith   if (flg) {
7625e42470aSBarry Smith     ls->steptol = tmp;
7635e42470aSBarry Smith   }
76409d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
76525cdf11fSBarry Smith   if (flg) {
76648d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
76777c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESNoLineSearch);
7685e42470aSBarry Smith     }
76929e0b56fSBarry Smith     else if (!PetscStrcmp(ver,"basicnonorms")) {
77029e0b56fSBarry Smith       SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);
77129e0b56fSBarry Smith     }
77248d91487SBarry Smith     else if (!PetscStrcmp(ver,"quadratic")) {
77377c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESQuadraticLineSearch);
7745e42470aSBarry Smith     }
77548d91487SBarry Smith     else if (!PetscStrcmp(ver,"cubic")) {
77677c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESCubicLineSearch);
7775e42470aSBarry Smith     }
778a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
7795e42470aSBarry Smith   }
7803a40ed3dSBarry Smith   PetscFunctionReturn(0);
7815e42470aSBarry Smith }
782*04d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
783*04d965bbSLois Curfman McInnes /*
784*04d965bbSLois Curfman McInnes    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
785*04d965bbSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
786*04d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
787*04d965bbSLois Curfman McInnes 
788*04d965bbSLois Curfman McInnes 
789*04d965bbSLois Curfman McInnes    Input Parameter:
790*04d965bbSLois Curfman McInnes .  snes - the SNES context
791*04d965bbSLois Curfman McInnes 
792*04d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
793*04d965bbSLois Curfman McInnes  */
7945615d1e5SSatish Balay #undef __FUNC__
7955615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
796f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
7975e42470aSBarry Smith {
79882bf6240SBarry Smith   int     ierr;
7995e42470aSBarry Smith   SNES_LS *neP;
8005e42470aSBarry Smith 
8013a40ed3dSBarry Smith   PetscFunctionBegin;
802a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
803a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
804a8c6a408SBarry Smith   }
80582bf6240SBarry Smith 
806f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
807f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
808f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
809f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
810f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
811f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
812f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
8135baf8537SBarry Smith   snes->nwork           = 0;
8145e42470aSBarry Smith 
8150452661fSBarry Smith   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
816ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
8175e42470aSBarry Smith   snes->data    	= (void *) neP;
8185e42470aSBarry Smith   neP->alpha		= 1.e-4;
8195e42470aSBarry Smith   neP->maxstep		= 1.e8;
8205e42470aSBarry Smith   neP->steptol		= 1.e-12;
8215e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
82282bf6240SBarry Smith 
823e1311b90SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch","SNESSetLineSearch_LS",
824e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
82582bf6240SBarry Smith 
8263a40ed3dSBarry Smith   PetscFunctionReturn(0);
8275e42470aSBarry Smith }
8285e42470aSBarry Smith 
82982bf6240SBarry Smith 
83082bf6240SBarry Smith 
83182bf6240SBarry Smith 
832