xref: /petsc/src/snes/impls/ls/ls.c (revision c2a78f1a4e367c9a8834b085c861515732dcb22d)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*c2a78f1aSLois Curfman McInnes static char vcid[] = "$Id: ls.c,v 1.106 1998/04/21 23:44:46 curfman 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 
904d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
1004d965bbSLois Curfman McInnes 
1104d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
1204d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
1304d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
1404d965bbSLois Curfman McInnes      respectively.
1504d965bbSLois Curfman McInnes 
1604d965bbSLois Curfman McInnes      The following basic routines are required for each nonlinear solver
1704d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
1804d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
1904d965bbSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear solver
2004d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
2104d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
2204d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
2304d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
2404d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
2504d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
2604d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
2704d965bbSLois Curfman McInnes      for all nonlinear solvers.
2804d965bbSLois Curfman McInnes 
2904d965bbSLois Curfman McInnes      Another key routine is:
3004d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
3104d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
3204d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
3304d965bbSLois Curfman McInnes      SNESSolve() if necessary.
3404d965bbSLois Curfman McInnes 
3504d965bbSLois Curfman McInnes      Additional basic routines are:
3604d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
3704d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
3804d965bbSLois Curfman McInnes                                     have actually been used.
3904d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
4004d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
4104d965bbSLois Curfman McInnes 
4204d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
4304d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
4404d965bbSLois Curfman McInnes      above description applies to these categories also.
4504d965bbSLois Curfman McInnes 
4604d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
475e42470aSBarry Smith /*
4804d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
4904d965bbSLois Curfman McInnes    method with a line search.
505e76c431SBarry Smith 
5104d965bbSLois Curfman McInnes    Input Parameters:
5204d965bbSLois Curfman McInnes .  snes - the SNES context
535e76c431SBarry Smith 
5404d965bbSLois Curfman McInnes    Output Parameter:
55*c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
5604d965bbSLois Curfman McInnes 
5704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
585e76c431SBarry Smith 
595e76c431SBarry Smith    Notes:
605e76c431SBarry Smith    This implements essentially a truncated Newton method with a
615e76c431SBarry Smith    line search.  By default a cubic backtracking line search
625e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
635e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
64393d2d9aSLois Curfman McInnes    and Schnabel.
655e76c431SBarry Smith */
665615d1e5SSatish Balay #undef __FUNC__
675615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS"
68f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
695e76c431SBarry Smith {
705e42470aSBarry Smith   SNES_LS       *neP = (SNES_LS *) snes->data;
71761aaf1bSLois Curfman McInnes   int           maxits, i, history_len, ierr, lits, lsfail;
72112a2221SBarry Smith   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
736b5873e3SBarry Smith   double        fnorm, gnorm, xnorm, ynorm, *history;
745e42470aSBarry Smith   Vec           Y, X, F, G, W, TMP;
755e76c431SBarry Smith 
763a40ed3dSBarry Smith   PetscFunctionBegin;
775e42470aSBarry Smith   history	= snes->conv_hist;	/* convergence history */
7851979daaSLois Curfman McInnes   history_len	= snes->conv_hist_size;	/* convergence history length */
795e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
805e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
8139e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
825e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
835e42470aSBarry Smith   G		= snes->work[1];
845e42470aSBarry Smith   W		= snes->work[2];
855e76c431SBarry Smith 
8625ed9cd1SBarry Smith   snes->iter = 0;
875334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
88cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
895e42470aSBarry Smith   snes->norm = fnorm;
9051979daaSLois Curfman McInnes   if (history) history[0] = fnorm;
9194a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
925e76c431SBarry Smith 
933a40ed3dSBarry Smith   if (fnorm < snes->atol) {*outits = 0; PetscFunctionReturn(0);}
9494a9d846SBarry Smith 
95d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
96d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
97d034289bSBarry Smith 
985e76c431SBarry Smith   for ( i=0; i<maxits; i++ ) {
995e42470aSBarry Smith     snes->iter = i+1;
1005e76c431SBarry Smith 
101ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1025334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg); CHKERRQ(ierr);
1035334005bSBarry Smith     ierr = SLESSetOperators(snes->sles,snes->jacobian,snes->jacobian_pre,flg); CHKERRQ(ierr);
10478b31e54SBarry Smith     ierr = SLESSolve(snes->sles,F,Y,&lits); CHKERRQ(ierr);
1057a00f4a9SLois Curfman McInnes     snes->linear_its += PetscAbsInt(lits);
106af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: iter=%d, linear solve iterations=%d\n",snes->iter,lits);
107ea4d3ed3SLois Curfman McInnes 
108ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
109ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
110ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
111ea4d3ed3SLois Curfman McInnes     */
11281b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always); CHKERRQ(ierr);
113ddd12767SBarry Smith     ierr = (*neP->LineSearch)(snes,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lsfail); CHKERRQ(ierr);
114af1ccdebSLois Curfman McInnes     PLogInfo(snes,"SNESSolve_EQ_LS: fnorm=%g, gnorm=%g, ynorm=%g, lsfail=%d\n",fnorm,gnorm,ynorm,lsfail);
115761aaf1bSLois Curfman McInnes     if (lsfail) snes->nfailures++;
1165e76c431SBarry Smith 
11739e2f89bSBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
11839e2f89bSBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
1195e76c431SBarry Smith     fnorm = gnorm;
1205e76c431SBarry Smith 
1215e42470aSBarry Smith     snes->norm = fnorm;
1225e76c431SBarry Smith     if (history && history_len > i+1) history[i+1] = fnorm;
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   }
14251979daaSLois Curfman McInnes   if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1;
1435e42470aSBarry Smith   *outits = i+1;
1443a40ed3dSBarry Smith   PetscFunctionReturn(0);
1455e76c431SBarry Smith }
14604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
14704d965bbSLois Curfman McInnes /*
14804d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
14904d965bbSLois Curfman McInnes    of the SNES_EQ_LS nonlinear solver.
15004d965bbSLois Curfman McInnes 
15104d965bbSLois Curfman McInnes    Input Parameter:
15204d965bbSLois Curfman McInnes .  snes - the SNES context
15304d965bbSLois Curfman McInnes .  x - the solution vector
15404d965bbSLois Curfman McInnes 
15504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
15604d965bbSLois Curfman McInnes 
15704d965bbSLois Curfman McInnes    Notes:
15804d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
15904d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
16004d965bbSLois Curfman McInnes    the call to SNESSolve().
16104d965bbSLois Curfman McInnes  */
1625615d1e5SSatish Balay #undef __FUNC__
1635615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
164f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
1655e76c431SBarry Smith {
1665e42470aSBarry Smith   int ierr;
1673a40ed3dSBarry Smith 
1683a40ed3dSBarry Smith   PetscFunctionBegin;
16981b6cf68SLois Curfman McInnes   snes->nwork = 4;
170d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1715e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
17281b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1733a40ed3dSBarry Smith   PetscFunctionReturn(0);
1745e76c431SBarry Smith }
17504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
17604d965bbSLois Curfman McInnes /*
17704d965bbSLois Curfman McInnes    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
17804d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
17904d965bbSLois Curfman McInnes 
18004d965bbSLois Curfman McInnes    Input Parameter:
18104d965bbSLois Curfman McInnes .  snes - the SNES context
18204d965bbSLois Curfman McInnes 
18304d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetDestroy()
18404d965bbSLois Curfman McInnes  */
1855615d1e5SSatish Balay #undef __FUNC__
186d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
187e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
1885e76c431SBarry Smith {
189393d2d9aSLois Curfman McInnes   int  ierr;
1903a40ed3dSBarry Smith 
1913a40ed3dSBarry Smith   PetscFunctionBegin;
1925baf8537SBarry Smith   if (snes->nwork) {
1934b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
19421c89e3eSBarry Smith   }
1950452661fSBarry Smith   PetscFree(snes->data);
1963a40ed3dSBarry Smith   PetscFunctionReturn(0);
1975e76c431SBarry Smith }
19804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1995615d1e5SSatish Balay #undef __FUNC__
2005615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
20182bf6240SBarry Smith 
2024b828684SBarry Smith /*@C
2035e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2045e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2055e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2065e76c431SBarry Smith 
2075e76c431SBarry Smith    Input Parameters:
2085e42470aSBarry Smith .  snes - nonlinear context
2095e76c431SBarry Smith .  x - current iterate
2105e76c431SBarry Smith .  f - residual evaluated at x
2115e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2125e76c431SBarry Smith .  w - work vector
2135e76c431SBarry Smith .  fnorm - 2-norm of f
2145e76c431SBarry Smith 
215c4a48953SLois Curfman McInnes    Output Parameters:
2165e76c431SBarry Smith .  g - residual evaluated at new iterate y
2175e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2185e76c431SBarry Smith .  gnorm - 2-norm of g
2195e76c431SBarry Smith .  ynorm - 2-norm of search length
220761aaf1bSLois Curfman McInnes .  flag - set to 0, indicating a successful line search
2215e76c431SBarry Smith 
222fee21e36SBarry Smith    Collective on SNES and Vec
223fee21e36SBarry Smith 
224c4a48953SLois Curfman McInnes    Options Database Key:
22509d61ba7SLois Curfman McInnes $  -snes_eq_ls basic
226c4a48953SLois Curfman McInnes 
22728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
22828ae5a21SLois Curfman McInnes 
229f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
23077c4ece6SBarry Smith           SNESSetLineSearch()
2315e76c431SBarry Smith @*/
2325e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
233761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
2345e76c431SBarry Smith {
2355e42470aSBarry Smith   int    ierr;
2365334005bSBarry Smith   Scalar mone = -1.0;
2375334005bSBarry Smith 
2383a40ed3dSBarry Smith   PetscFunctionBegin;
239761aaf1bSLois Curfman McInnes   *flag = 0;
2407857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
241cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
242ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
243ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
244cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
2457857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2463a40ed3dSBarry Smith   PetscFunctionReturn(0);
2475e76c431SBarry Smith }
24804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2495615d1e5SSatish Balay #undef __FUNC__
25029e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
25182bf6240SBarry Smith 
25229e0b56fSBarry Smith /*@C
25329e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
25429e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
25529e0b56fSBarry Smith    even compute the norm of the function or search direction; this
25629e0b56fSBarry Smith    is intended only when you know the full step is fine and are
25729e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
25829e0b56fSBarry Smith    example, you are running always for a fixed number of Newton
25929e0b56fSBarry Smith    steps).
26029e0b56fSBarry Smith 
26129e0b56fSBarry Smith    Input Parameters:
26229e0b56fSBarry Smith .  snes - nonlinear context
26329e0b56fSBarry Smith .  x - current iterate
26429e0b56fSBarry Smith .  f - residual evaluated at x
26529e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
26629e0b56fSBarry Smith .  w - work vector
26729e0b56fSBarry Smith .  fnorm - 2-norm of f
26829e0b56fSBarry Smith 
26929e0b56fSBarry Smith    Output Parameters:
27029e0b56fSBarry Smith .  g - residual evaluated at new iterate y
27129e0b56fSBarry Smith .  gnorm - not changed
27229e0b56fSBarry Smith .  ynorm - not changed
27329e0b56fSBarry Smith .  flag - set to 0, indicating a successful line search
27429e0b56fSBarry Smith 
275fee21e36SBarry Smith    Collective on SNES and Vec
276fee21e36SBarry Smith 
27729e0b56fSBarry Smith    Options Database Key:
27829e0b56fSBarry Smith $  -snes_eq_ls basicnonorms
27929e0b56fSBarry Smith 
28029e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
28129e0b56fSBarry Smith 
28229e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
28329e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
28429e0b56fSBarry Smith @*/
28529e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
28629e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
28729e0b56fSBarry Smith {
28829e0b56fSBarry Smith   int    ierr;
28929e0b56fSBarry Smith   Scalar mone = -1.0;
29029e0b56fSBarry Smith 
2913a40ed3dSBarry Smith   PetscFunctionBegin;
29229e0b56fSBarry Smith   *flag = 0;
29329e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
29429e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
29529e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
29629e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
29829e0b56fSBarry Smith }
29904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30029e0b56fSBarry Smith #undef __FUNC__
3015615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3024b828684SBarry Smith /*@C
303f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3045e76c431SBarry Smith 
3055e76c431SBarry Smith    Input Parameters:
3065e42470aSBarry Smith .  snes - nonlinear context
3075e76c431SBarry Smith .  x - current iterate
3085e76c431SBarry Smith .  f - residual evaluated at x
3095e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3105e76c431SBarry Smith .  w - work vector
3115e76c431SBarry Smith .  fnorm - 2-norm of f
3125e76c431SBarry Smith 
313393d2d9aSLois Curfman McInnes    Output Parameters:
3145e76c431SBarry Smith .  g - residual evaluated at new iterate y
3155e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3165e76c431SBarry Smith .  gnorm - 2-norm of g
3175e76c431SBarry Smith .  ynorm - 2-norm of search length
318761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
3195e76c431SBarry Smith 
320fee21e36SBarry Smith    Collective on SNES
321fee21e36SBarry Smith 
322c4a48953SLois Curfman McInnes    Options Database Key:
32309d61ba7SLois Curfman McInnes $  -snes_eq_ls cubic
324c4a48953SLois Curfman McInnes 
3255e76c431SBarry Smith    Notes:
3265e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3275e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3285e76c431SBarry Smith 
32928ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33028ae5a21SLois Curfman McInnes 
33177c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
3325e76c431SBarry Smith @*/
3335e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
334761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3355e76c431SBarry Smith {
336406556e6SLois Curfman McInnes   /*
337406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
338406556e6SLois Curfman McInnes      minimization problem:
339406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
340406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
341406556e6SLois Curfman McInnes    */
342406556e6SLois Curfman McInnes 
343ddd12767SBarry Smith   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
344ea4d3ed3SLois Curfman McInnes   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
3453a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
3465e42470aSBarry Smith   Scalar  cinitslope, clambda;
3476b5873e3SBarry Smith #endif
3485e42470aSBarry Smith   int     ierr, count;
3495e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
3505334005bSBarry Smith   Scalar  mone = -1.0,scale;
3515e76c431SBarry Smith 
3523a40ed3dSBarry Smith   PetscFunctionBegin;
3537857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
354761aaf1bSLois Curfman McInnes   *flag   = 0;
3555e76c431SBarry Smith   alpha   = neP->alpha;
3565e76c431SBarry Smith   maxstep = neP->maxstep;
3575e76c431SBarry Smith   steptol = neP->steptol;
3585e76c431SBarry Smith 
359cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
360a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
361a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
362a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
363a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
364a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
365ad922ac9SBarry Smith     goto theend1;
36694a9d846SBarry Smith   }
3675e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3685e42470aSBarry Smith     scale = maxstep/(*ynorm);
3693a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
37094a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",real(scale));
3716b5873e3SBarry Smith #else
37294a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
3736b5873e3SBarry Smith #endif
374393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
3755e76c431SBarry Smith     *ynorm = maxstep;
3765e76c431SBarry Smith   }
3775e76c431SBarry Smith   minlambda = steptol/(*ynorm);
378a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
3793a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
380a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
3815e42470aSBarry Smith   initslope = real(cinitslope);
3825e42470aSBarry Smith #else
383a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
3845e42470aSBarry Smith #endif
3855e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3865e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3875e76c431SBarry Smith 
388393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
3895334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
39078b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
391cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);
392406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
393393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
39494a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
395ad922ac9SBarry Smith     goto theend1;
3965e76c431SBarry Smith   }
3975e76c431SBarry Smith 
3985e76c431SBarry Smith   /* Fit points with quadratic */
3995e76c431SBarry Smith   lambda = 1.0; count = 0;
400a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4015e76c431SBarry Smith   lambdaprev = lambda;
4025e76c431SBarry Smith   gnormprev = *gnorm;
403ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
404ddd12767SBarry Smith   else lambda = lambdatemp;
405393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w); CHKERRQ(ierr);
406ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
4073a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
408ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4095e42470aSBarry Smith #else
410ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4115e42470aSBarry Smith #endif
41278b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
413cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
414406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
415393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
416ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
417ad922ac9SBarry Smith     goto theend1;
4185e76c431SBarry Smith   }
4195e76c431SBarry Smith 
4205e76c431SBarry Smith   /* Fit points with cubic */
4215e76c431SBarry Smith   count = 1;
4225e76c431SBarry Smith   while (1) {
4235e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4242b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4252b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
426ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
427393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
428761aaf1bSLois Curfman McInnes       *flag = -1; break;
4295e76c431SBarry Smith     }
430406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
431406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
432ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4332b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4345e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4355e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4365e76c431SBarry Smith     if (a == 0.0) {
4375e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4385e76c431SBarry Smith     } else {
4395e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4405e76c431SBarry Smith     }
4415e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4425e76c431SBarry Smith       lambdatemp = .5*lambda;
4435e76c431SBarry Smith     }
4445e76c431SBarry Smith     lambdaprev = lambda;
4455e76c431SBarry Smith     gnormprev = *gnorm;
4465e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4475e76c431SBarry Smith       lambda = .1*lambda;
4485e76c431SBarry Smith     }
4495e76c431SBarry Smith     else lambda = lambdatemp;
450393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w ); CHKERRQ(ierr);
451ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
4523a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
453ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
454393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4555e42470aSBarry Smith #else
456ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4575e42470aSBarry Smith #endif
45878b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
459cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
460406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
461393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
462ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
4632715565aSLois Curfman McInnes       goto theend1;
4642b022350SLois Curfman McInnes     } else {
4652b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
4665e76c431SBarry Smith     }
4675e76c431SBarry Smith     count++;
4685e76c431SBarry Smith   }
469ad922ac9SBarry Smith   theend1:
4707857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4713a40ed3dSBarry Smith   PetscFunctionReturn(0);
4725e76c431SBarry Smith }
47304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4745615d1e5SSatish Balay #undef __FUNC__
4755615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
4764b828684SBarry Smith /*@C
477f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
4785e76c431SBarry Smith 
4795e42470aSBarry Smith    Input Parameters:
480f59ffdedSLois Curfman McInnes .  snes - the SNES context
4815e42470aSBarry Smith .  x - current iterate
4825e42470aSBarry Smith .  f - residual evaluated at x
4835e42470aSBarry Smith .  y - search direction (contains new iterate on output)
4845e42470aSBarry Smith .  w - work vector
4855e42470aSBarry Smith .  fnorm - 2-norm of f
4865e42470aSBarry Smith 
487c4a48953SLois Curfman McInnes    Output Parameters:
4885e42470aSBarry Smith .  g - residual evaluated at new iterate y
4895e42470aSBarry Smith .  y - new iterate (contains search direction on input)
4905e42470aSBarry Smith .  gnorm - 2-norm of g
4915e42470aSBarry Smith .  ynorm - 2-norm of search length
492761aaf1bSLois Curfman McInnes .  flag - 0 if line search succeeds; -1 on failure.
4935e42470aSBarry Smith 
494fee21e36SBarry Smith    Collective on SNES and Vec
495fee21e36SBarry Smith 
496c4a48953SLois Curfman McInnes    Options Database Key:
49709d61ba7SLois Curfman McInnes $  -snes_eq_ls quadratic
498c4a48953SLois Curfman McInnes 
4995e42470aSBarry Smith    Notes:
50077c4ece6SBarry Smith    Use SNESSetLineSearch()
501f63b844aSLois Curfman McInnes    to set this routine within the SNES_EQ_LS method.
5025e42470aSBarry Smith 
503f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
504f59ffdedSLois Curfman McInnes 
50577c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
5065e42470aSBarry Smith @*/
5075e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
508761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5095e76c431SBarry Smith {
510406556e6SLois Curfman McInnes   /*
511406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
512406556e6SLois Curfman McInnes      minimization problem:
513406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
514406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
515406556e6SLois Curfman McInnes    */
51640011551SBarry Smith   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
5173a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
5185e42470aSBarry Smith   Scalar  cinitslope,clambda;
5196b5873e3SBarry Smith #endif
5205e42470aSBarry Smith   int     ierr,count;
5215e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
5225334005bSBarry Smith   Scalar  mone = -1.0,scale;
5235e76c431SBarry Smith 
5243a40ed3dSBarry Smith   PetscFunctionBegin;
5257857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
526761aaf1bSLois Curfman McInnes   *flag   = 0;
5275e42470aSBarry Smith   alpha   = neP->alpha;
5285e42470aSBarry Smith   maxstep = neP->maxstep;
5295e42470aSBarry Smith   steptol = neP->steptol;
5305e76c431SBarry Smith 
531cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
532b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
53394a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
534b37302c6SLois Curfman McInnes     *gnorm = fnorm;
535b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
536b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
537ad922ac9SBarry Smith     goto theend2;
53894a9d846SBarry Smith   }
5395e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5405e42470aSBarry Smith     scale = maxstep/(*ynorm);
541393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
5425e42470aSBarry Smith     *ynorm = maxstep;
5435e76c431SBarry Smith   }
5445e42470aSBarry Smith   minlambda = steptol/(*ynorm);
545a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
5463a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
547a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
5485e42470aSBarry Smith   initslope = real(cinitslope);
5495e42470aSBarry Smith #else
550a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
5515e42470aSBarry Smith #endif
5525e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
5535e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
5545e42470aSBarry Smith 
555393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
5565334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
55778b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
558cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
559406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
560393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
56194a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
562ad922ac9SBarry Smith     goto theend2;
5635e42470aSBarry Smith   }
5645e42470aSBarry Smith 
5655e42470aSBarry Smith   /* Fit points with quadratic */
5665e42470aSBarry Smith   lambda = 1.0; count = 0;
5675e42470aSBarry Smith   count = 1;
5685e42470aSBarry Smith   while (1) {
5695e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
570981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
571981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
572ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
573393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
574761aaf1bSLois Curfman McInnes       *flag = -1; break;
5755e42470aSBarry Smith     }
576a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5775e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
5785e42470aSBarry Smith       lambda = .1*lambda;
5795e42470aSBarry Smith     } else lambda = lambdatemp;
580393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w); CHKERRQ(ierr);
5815334005bSBarry Smith     lambda = -lambda;
5823a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
583393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
5845e42470aSBarry Smith #else
585393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
5865e42470aSBarry Smith #endif
58778b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
588cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
589406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
590393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
591981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
59206259719SBarry Smith       break;
5935e42470aSBarry Smith     }
5945e42470aSBarry Smith     count++;
5955e42470aSBarry Smith   }
596ad922ac9SBarry Smith   theend2:
5977857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
5983a40ed3dSBarry Smith   PetscFunctionReturn(0);
5995e76c431SBarry Smith }
60004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6015615d1e5SSatish Balay #undef __FUNC__
602d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
603c9e6a524SLois Curfman McInnes /*@C
60477c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
605f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
6065e76c431SBarry Smith 
6075e76c431SBarry Smith    Input Parameters:
608eafb4bcbSBarry Smith .  snes - nonlinear context obtained from SNESCreate()
6095e76c431SBarry Smith .  func - pointer to int function
6105e76c431SBarry Smith 
611fee21e36SBarry Smith    Collective on SNES
612fee21e36SBarry Smith 
613c4a48953SLois Curfman McInnes    Available Routines:
614f59ffdedSLois Curfman McInnes .  SNESCubicLineSearch() - default line search
615f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
616f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
61782bf6240SBarry Smith .  SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel)
6185e76c431SBarry Smith 
619c4a48953SLois Curfman McInnes     Options Database Keys:
62009d61ba7SLois Curfman McInnes $   -snes_eq_ls [basic,quadratic,cubic]
621912ebf9aSLois Curfman McInnes $   -snes_eq_ls_alpha <alpha>
622912ebf9aSLois Curfman McInnes $   -snes_eq_ls_maxstep <max>
623912ebf9aSLois Curfman McInnes $   -snes_eq_ls_steptol <steptol>
624c4a48953SLois Curfman McInnes 
6255e76c431SBarry Smith    Calling sequence of func:
626f59ffdedSLois Curfman McInnes    func (SNES snes, Vec x, Vec f, Vec g, Vec y,
627761aaf1bSLois Curfman McInnes          Vec w, double fnorm, double *ynorm,
628761aaf1bSLois Curfman McInnes          double *gnorm, *flag)
6295e76c431SBarry Smith 
6305e76c431SBarry Smith     Input parameters for func:
6315e42470aSBarry Smith .   snes - nonlinear context
6325e76c431SBarry Smith .   x - current iterate
6335e76c431SBarry Smith .   f - residual evaluated at x
6345e76c431SBarry Smith .   y - search direction (contains new iterate on output)
6355e76c431SBarry Smith .   w - work vector
6365e76c431SBarry Smith .   fnorm - 2-norm of f
6375e76c431SBarry Smith 
6385e76c431SBarry Smith     Output parameters for func:
6395e76c431SBarry Smith .   g - residual evaluated at new iterate y
6405e76c431SBarry Smith .   y - new iterate (contains search direction on input)
6415e76c431SBarry Smith .   gnorm - 2-norm of g
6425e76c431SBarry Smith .   ynorm - 2-norm of search length
643761aaf1bSLois Curfman McInnes .   flag - set to 0 if the line search succeeds; a nonzero integer
644761aaf1bSLois Curfman McInnes            on failure.
645f59ffdedSLois Curfman McInnes 
646f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
647f59ffdedSLois Curfman McInnes 
648f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
6495e76c431SBarry Smith @*/
65077c4ece6SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
651761aaf1bSLois Curfman McInnes                              double,double*,double*,int*))
6525e76c431SBarry Smith {
65382bf6240SBarry Smith   int ierr, (*f)(SNES,int (*f)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*));
65482bf6240SBarry Smith 
6553a40ed3dSBarry Smith   PetscFunctionBegin;
656e1311b90SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch",(void **)&f);CHKERRQ(ierr);
65782bf6240SBarry Smith   if (f) {
65882bf6240SBarry Smith     ierr = (*f)(snes,func);CHKERRQ(ierr);
65982bf6240SBarry Smith   }
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
6615e76c431SBarry Smith }
66204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
66382bf6240SBarry Smith #undef __FUNC__
66482bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
66582bf6240SBarry Smith int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
66682bf6240SBarry Smith                          double,double*,double*,int*))
66782bf6240SBarry Smith {
66882bf6240SBarry Smith   PetscFunctionBegin;
66982bf6240SBarry Smith   ((SNES_LS *)(snes->data))->LineSearch = func;
67082bf6240SBarry Smith   PetscFunctionReturn(0);
67182bf6240SBarry Smith }
67204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
67304d965bbSLois Curfman McInnes /*
67404d965bbSLois Curfman McInnes    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
67582bf6240SBarry Smith 
67604d965bbSLois Curfman McInnes    Input Parameter:
67704d965bbSLois Curfman McInnes .  snes - the SNES context
67804d965bbSLois Curfman McInnes 
67904d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
68004d965bbSLois Curfman McInnes */
68104d965bbSLois Curfman McInnes 
6825615d1e5SSatish Balay #undef __FUNC__
683d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
684f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
6855e42470aSBarry Smith {
6865e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
6876daaf66cSBarry Smith 
6883a40ed3dSBarry Smith   PetscFunctionBegin;
68976be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
69076be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
69176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
69276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
69376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
6943a40ed3dSBarry Smith   PetscFunctionReturn(0);
6955e42470aSBarry Smith }
69604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
69704d965bbSLois Curfman McInnes /*
69804d965bbSLois Curfman McInnes    SNESView_EQ_LS - Prints the SNES_EQ_LS data structure.
69904d965bbSLois Curfman McInnes 
70004d965bbSLois Curfman McInnes    Input Parameters:
70104d965bbSLois Curfman McInnes .  SNES - the SNES context
70204d965bbSLois Curfman McInnes .  viewer - visualization context
70304d965bbSLois Curfman McInnes 
70404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
70504d965bbSLois Curfman McInnes */
7065615d1e5SSatish Balay #undef __FUNC__
707d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
708e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
709a935fc98SLois Curfman McInnes {
710a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
711d636dbe3SBarry Smith   FILE       *fd;
71219bcc07fSBarry Smith   char       *cstr;
71351695f54SLois Curfman McInnes   int        ierr;
71419bcc07fSBarry Smith   ViewerType vtype;
715a935fc98SLois Curfman McInnes 
7163a40ed3dSBarry Smith   PetscFunctionBegin;
71719bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
71819bcc07fSBarry Smith   if (vtype  == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
71990ace30eSBarry Smith     ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
72019bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch) cstr = "SNESNoLineSearch";
72119bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
72219bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch) cstr = "SNESCubicLineSearch";
72319bcc07fSBarry Smith     else cstr = "unknown";
72477c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"    line search variant: %s\n",cstr);
72577c4ece6SBarry Smith     PetscFPrintf(snes->comm,fd,"    alpha=%g, maxstep=%g, steptol=%g\n",
726a935fc98SLois Curfman McInnes                  ls->alpha,ls->maxstep,ls->steptol);
7275cd90555SBarry Smith   } else {
7285cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
72919bcc07fSBarry Smith   }
7303a40ed3dSBarry Smith   PetscFunctionReturn(0);
731a935fc98SLois Curfman McInnes }
73204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
73304d965bbSLois Curfman McInnes /*
73404d965bbSLois Curfman McInnes    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
73504d965bbSLois Curfman McInnes 
73604d965bbSLois Curfman McInnes    Input Parameter:
73704d965bbSLois Curfman McInnes .  snes - the SNES context
73804d965bbSLois Curfman McInnes 
73904d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
74004d965bbSLois Curfman McInnes */
7415615d1e5SSatish Balay #undef __FUNC__
7425615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
743f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
7445e42470aSBarry Smith {
7455e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
7465e42470aSBarry Smith   char    ver[16];
7475e42470aSBarry Smith   double  tmp;
74825cdf11fSBarry Smith   int     ierr,flg;
7495e42470aSBarry Smith 
7503a40ed3dSBarry Smith   PetscFunctionBegin;
75109d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
75225cdf11fSBarry Smith   if (flg) {
7535e42470aSBarry Smith     ls->alpha = tmp;
7545e42470aSBarry Smith   }
75509d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
75625cdf11fSBarry Smith   if (flg) {
7575e42470aSBarry Smith     ls->maxstep = tmp;
7585e42470aSBarry Smith   }
75909d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
76025cdf11fSBarry Smith   if (flg) {
7615e42470aSBarry Smith     ls->steptol = tmp;
7625e42470aSBarry Smith   }
76309d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
76425cdf11fSBarry Smith   if (flg) {
76548d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
76677c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESNoLineSearch);
7675e42470aSBarry Smith     }
76829e0b56fSBarry Smith     else if (!PetscStrcmp(ver,"basicnonorms")) {
76929e0b56fSBarry Smith       SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);
77029e0b56fSBarry Smith     }
77148d91487SBarry Smith     else if (!PetscStrcmp(ver,"quadratic")) {
77277c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESQuadraticLineSearch);
7735e42470aSBarry Smith     }
77448d91487SBarry Smith     else if (!PetscStrcmp(ver,"cubic")) {
77577c4ece6SBarry Smith       SNESSetLineSearch(snes,SNESCubicLineSearch);
7765e42470aSBarry Smith     }
777a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
7785e42470aSBarry Smith   }
7793a40ed3dSBarry Smith   PetscFunctionReturn(0);
7805e42470aSBarry Smith }
78104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
78204d965bbSLois Curfman McInnes /*
78304d965bbSLois Curfman McInnes    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
78404d965bbSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
78504d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
78604d965bbSLois Curfman McInnes 
78704d965bbSLois Curfman McInnes 
78804d965bbSLois Curfman McInnes    Input Parameter:
78904d965bbSLois Curfman McInnes .  snes - the SNES context
79004d965bbSLois Curfman McInnes 
79104d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
79204d965bbSLois Curfman McInnes  */
7935615d1e5SSatish Balay #undef __FUNC__
7945615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
795f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
7965e42470aSBarry Smith {
79782bf6240SBarry Smith   int     ierr;
7985e42470aSBarry Smith   SNES_LS *neP;
7995e42470aSBarry Smith 
8003a40ed3dSBarry Smith   PetscFunctionBegin;
801a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
802a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
803a8c6a408SBarry Smith   }
80482bf6240SBarry Smith 
805f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
806f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
807f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
808f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
809f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
810f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
811f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
8125baf8537SBarry Smith   snes->nwork           = 0;
8135e42470aSBarry Smith 
8140452661fSBarry Smith   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
815ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
8165e42470aSBarry Smith   snes->data    	= (void *) neP;
8175e42470aSBarry Smith   neP->alpha		= 1.e-4;
8185e42470aSBarry Smith   neP->maxstep		= 1.e8;
8195e42470aSBarry Smith   neP->steptol		= 1.e-12;
8205e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
82182bf6240SBarry Smith 
822e1311b90SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch","SNESSetLineSearch_LS",
823e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
82482bf6240SBarry Smith 
8253a40ed3dSBarry Smith   PetscFunctionReturn(0);
8265e42470aSBarry Smith }
8275e42470aSBarry Smith 
82882bf6240SBarry Smith 
82982bf6240SBarry Smith 
83082bf6240SBarry Smith 
831