xref: /petsc/src/snes/impls/ls/ls.c (revision 06d0569a1a49123e2b3e747d9c02bd1779fb4f7e)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*06d0569aSBarry Smith static char vcid[] = "$Id: ls.c,v 1.121 1999/02/01 03:26:42 curfman Exp bsmith $";
35e76c431SBarry Smith #endif
45e76c431SBarry Smith 
570f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
65e42470aSBarry Smith 
704d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
804d965bbSLois Curfman McInnes 
904d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
1004d965bbSLois Curfman McInnes      for solving a system of nonlinear equations, using the SLES, Vec,
1104d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
1204d965bbSLois Curfman McInnes      respectively.
1304d965bbSLois Curfman McInnes 
14fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
1504d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
1604d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
17fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
1804d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
1904d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
2004d965bbSLois Curfman McInnes      we use _EQ_LS (e.g., SNESCreate_EQ_LS, SNESSolve_EQ_LS) for solving
2104d965bbSLois Curfman McInnes      systems of nonlinear equations (EQ) with a line search (LS) method.
2204d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
2304d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
2404d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
2504d965bbSLois Curfman McInnes      for all nonlinear solvers.
2604d965bbSLois Curfman McInnes 
2704d965bbSLois Curfman McInnes      Another key routine is:
2804d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
2904d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
3004d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
3104d965bbSLois Curfman McInnes      SNESSolve() if necessary.
3204d965bbSLois Curfman McInnes 
3304d965bbSLois Curfman McInnes      Additional basic routines are:
3404d965bbSLois Curfman McInnes           SNESPrintHelp_XXX()       - Prints nonlinear solver runtime options
3504d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
3604d965bbSLois Curfman McInnes                                       have actually been used.
3704d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
3804d965bbSLois Curfman McInnes      SNESPrintHelp() and SNESView().
3904d965bbSLois Curfman McInnes 
4004d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
4104d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
4204d965bbSLois Curfman McInnes      above description applies to these categories also.
4304d965bbSLois Curfman McInnes 
4404d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
455e42470aSBarry Smith /*
4604d965bbSLois Curfman McInnes    SNESSolve_EQ_LS - Solves a nonlinear system with a truncated Newton
4704d965bbSLois Curfman McInnes    method with a line search.
485e76c431SBarry Smith 
4904d965bbSLois Curfman McInnes    Input Parameters:
5004d965bbSLois Curfman McInnes .  snes - the SNES context
515e76c431SBarry Smith 
5204d965bbSLois Curfman McInnes    Output Parameter:
53c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
5404d965bbSLois Curfman McInnes 
5504d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
565e76c431SBarry Smith 
575e76c431SBarry Smith    Notes:
585e76c431SBarry Smith    This implements essentially a truncated Newton method with a
595e76c431SBarry Smith    line search.  By default a cubic backtracking line search
605e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
615e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
62393d2d9aSLois Curfman McInnes    and Schnabel.
635e76c431SBarry Smith */
645615d1e5SSatish Balay #undef __FUNC__
655615d1e5SSatish Balay #define __FUNC__ "SNESSolve_EQ_LS"
66f63b844aSLois Curfman McInnes int SNESSolve_EQ_LS(SNES snes,int *outits)
675e76c431SBarry Smith {
685e42470aSBarry Smith   SNES_LS       *neP = (SNES_LS *) snes->data;
69761aaf1bSLois Curfman McInnes   int           maxits, i, history_len, ierr, lits, lsfail;
70112a2221SBarry Smith   MatStructure  flg = DIFFERENT_NONZERO_PATTERN;
716b5873e3SBarry Smith   double        fnorm, gnorm, xnorm, ynorm, *history;
725e42470aSBarry Smith   Vec           Y, X, F, G, W, TMP;
735e76c431SBarry Smith 
743a40ed3dSBarry Smith   PetscFunctionBegin;
755e42470aSBarry Smith   history	= snes->conv_hist;	/* convergence history */
7651979daaSLois Curfman McInnes   history_len	= snes->conv_hist_size;	/* convergence history length */
775e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
785e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
7939e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
805e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
815e42470aSBarry Smith   G		= snes->work[1];
825e42470aSBarry Smith   W		= snes->work[2];
835e76c431SBarry Smith 
8425ed9cd1SBarry Smith   snes->iter = 0;
855334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F); CHKERRQ(ierr);  /*  F(X)      */
86cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm); CHKERRQ(ierr);	/* fnorm <- ||F||  */
87*06d0569aSBarry Smith   PetscAMSTakeAccess(snes);
885e42470aSBarry Smith   snes->norm = fnorm;
89*06d0569aSBarry Smith   PetscAMSGrantAccess(snes);
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 
121*06d0569aSBarry Smith     PetscAMSTakeAccess(snes);
1225e42470aSBarry Smith     snes->norm = fnorm;
123*06d0569aSBarry Smith     PetscAMSGrantAccess(snes);
1245e76c431SBarry Smith     if (history && history_len > i+1) history[i+1] = fnorm;
12594a424c1SBarry Smith     SNESMonitor(snes,i+1,fnorm);
1265e76c431SBarry Smith 
1275e76c431SBarry Smith     /* Test for convergence */
12829e0b56fSBarry Smith     if (snes->converged) {
12929e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm); CHKERRQ(ierr);	/* xnorm = || X || */
130bbb6d6a8SBarry Smith       if ((*snes->converged)(snes,xnorm,ynorm,fnorm,snes->cnvP)) {
13116c95cacSBarry Smith         break;
13216c95cacSBarry Smith       }
13316c95cacSBarry Smith     }
13429e0b56fSBarry Smith   }
13539e2f89bSBarry Smith   if (X != snes->vec_sol) {
136393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol); CHKERRQ(ierr);
13739e2f89bSBarry Smith     snes->vec_sol_always  = snes->vec_sol;
13839e2f89bSBarry Smith     snes->vec_func_always = snes->vec_func;
13939e2f89bSBarry Smith   }
14052392280SLois Curfman McInnes   if (i == maxits) {
141981c4779SBarry Smith     PLogInfo(snes,"SNESSolve_EQ_LS: Maximum number of iterations has been reached: %d\n",maxits);
14252392280SLois Curfman McInnes     i--;
14352392280SLois Curfman McInnes   }
14451979daaSLois Curfman McInnes   if (history) snes->conv_act_size = (history_len < i+1) ? history_len : i+1;
1455e42470aSBarry Smith   *outits = i+1;
1463a40ed3dSBarry Smith   PetscFunctionReturn(0);
1475e76c431SBarry Smith }
14804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
14904d965bbSLois Curfman McInnes /*
15004d965bbSLois Curfman McInnes    SNESSetUp_EQ_LS - Sets up the internal data structures for the later use
15104d965bbSLois Curfman McInnes    of the SNES_EQ_LS nonlinear solver.
15204d965bbSLois Curfman McInnes 
15304d965bbSLois Curfman McInnes    Input Parameter:
15404d965bbSLois Curfman McInnes .  snes - the SNES context
15504d965bbSLois Curfman McInnes .  x - the solution vector
15604d965bbSLois Curfman McInnes 
15704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
15804d965bbSLois Curfman McInnes 
15904d965bbSLois Curfman McInnes    Notes:
16004d965bbSLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
16104d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
16204d965bbSLois Curfman McInnes    the call to SNESSolve().
16304d965bbSLois Curfman McInnes  */
1645615d1e5SSatish Balay #undef __FUNC__
1655615d1e5SSatish Balay #define __FUNC__ "SNESSetUp_EQ_LS"
166f63b844aSLois Curfman McInnes int SNESSetUp_EQ_LS(SNES snes)
1675e76c431SBarry Smith {
1685e42470aSBarry Smith   int ierr;
1693a40ed3dSBarry Smith 
1703a40ed3dSBarry Smith   PetscFunctionBegin;
17181b6cf68SLois Curfman McInnes   snes->nwork = 4;
172d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
1735e42470aSBarry Smith   PLogObjectParents(snes,snes->nwork,snes->work);
17481b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
1753a40ed3dSBarry Smith   PetscFunctionReturn(0);
1765e76c431SBarry Smith }
17704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
17804d965bbSLois Curfman McInnes /*
17904d965bbSLois Curfman McInnes    SNESDestroy_EQ_LS - Destroys the private SNES_LS context that was created
18004d965bbSLois Curfman McInnes    with SNESCreate_EQ_LS().
18104d965bbSLois Curfman McInnes 
18204d965bbSLois Curfman McInnes    Input Parameter:
18304d965bbSLois Curfman McInnes .  snes - the SNES context
18404d965bbSLois Curfman McInnes 
185de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
18604d965bbSLois Curfman McInnes  */
1875615d1e5SSatish Balay #undef __FUNC__
188d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy_EQ_LS"
189e1311b90SBarry Smith int SNESDestroy_EQ_LS(SNES snes)
1905e76c431SBarry Smith {
191393d2d9aSLois Curfman McInnes   int  ierr;
1923a40ed3dSBarry Smith 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
1945baf8537SBarry Smith   if (snes->nwork) {
1954b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork); CHKERRQ(ierr);
19621c89e3eSBarry Smith   }
1970452661fSBarry Smith   PetscFree(snes->data);
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
1995e76c431SBarry Smith }
20004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2015615d1e5SSatish Balay #undef __FUNC__
2025615d1e5SSatish Balay #define __FUNC__ "SNESNoLineSearch"
20382bf6240SBarry Smith 
2044b828684SBarry Smith /*@C
2055e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
2065e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
2075e76c431SBarry Smith    to serve as a template and is not recommended for general use.
2085e76c431SBarry Smith 
209c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
210c7afd0dbSLois Curfman McInnes 
2115e76c431SBarry Smith    Input Parameters:
212c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
2135e76c431SBarry Smith .  x - current iterate
2145e76c431SBarry Smith .  f - residual evaluated at x
2155e76c431SBarry Smith .  y - search direction (contains new iterate on output)
2165e76c431SBarry Smith .  w - work vector
217c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
2185e76c431SBarry Smith 
219c4a48953SLois Curfman McInnes    Output Parameters:
220c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
2215e76c431SBarry Smith .  y - new iterate (contains search direction on input)
2225e76c431SBarry Smith .  gnorm - 2-norm of g
2235e76c431SBarry Smith .  ynorm - 2-norm of search length
224c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
225fee21e36SBarry Smith 
226c4a48953SLois Curfman McInnes    Options Database Key:
227c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basic - Activates SNESNoLineSearch()
228c4a48953SLois Curfman McInnes 
22936851e7fSLois Curfman McInnes    Level: advanced
23036851e7fSLois Curfman McInnes 
23128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
23228ae5a21SLois Curfman McInnes 
233f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
23477c4ece6SBarry Smith           SNESSetLineSearch()
2355e76c431SBarry Smith @*/
2365e42470aSBarry Smith int SNESNoLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
237761aaf1bSLois Curfman McInnes                      double fnorm, double *ynorm, double *gnorm,int *flag )
2385e76c431SBarry Smith {
2395e42470aSBarry Smith   int    ierr;
2405334005bSBarry Smith   Scalar mone = -1.0;
2415334005bSBarry Smith 
2423a40ed3dSBarry Smith   PetscFunctionBegin;
243761aaf1bSLois Curfman McInnes   *flag = 0;
2447857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
245cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);       /* ynorm = || y || */
246ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
247ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
248cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);       /* gnorm = || g || */
2497857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
2503a40ed3dSBarry Smith   PetscFunctionReturn(0);
2515e76c431SBarry Smith }
25204d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2535615d1e5SSatish Balay #undef __FUNC__
25429e0b56fSBarry Smith #define __FUNC__ "SNESNoLineSearchNoNorms"
25582bf6240SBarry Smith 
25629e0b56fSBarry Smith /*@C
25729e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
25829e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
25929e0b56fSBarry Smith    even compute the norm of the function or search direction; this
26029e0b56fSBarry Smith    is intended only when you know the full step is fine and are
26129e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
262c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
263c7afd0dbSLois Curfman McInnes 
264c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
26529e0b56fSBarry Smith 
26629e0b56fSBarry Smith    Input Parameters:
267c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
26829e0b56fSBarry Smith .  x - current iterate
26929e0b56fSBarry Smith .  f - residual evaluated at x
27029e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
27129e0b56fSBarry Smith .  w - work vector
272c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
27329e0b56fSBarry Smith 
27429e0b56fSBarry Smith    Output Parameters:
275c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
27629e0b56fSBarry Smith .  gnorm - not changed
27729e0b56fSBarry Smith .  ynorm - not changed
278c7afd0dbSLois Curfman McInnes -  flag - set to 0, indicating a successful line search
279fee21e36SBarry Smith 
28029e0b56fSBarry Smith    Options Database Key:
281c7afd0dbSLois Curfman McInnes .  -snes_eq_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
28229e0b56fSBarry Smith 
28336851e7fSLois Curfman McInnes    Level: advanced
28436851e7fSLois Curfman McInnes 
28529e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
28629e0b56fSBarry Smith 
28729e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
28829e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
28929e0b56fSBarry Smith @*/
29029e0b56fSBarry Smith int SNESNoLineSearchNoNorms(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
29129e0b56fSBarry Smith                      double fnorm, double *ynorm, double *gnorm,int *flag )
29229e0b56fSBarry Smith {
29329e0b56fSBarry Smith   int    ierr;
29429e0b56fSBarry Smith   Scalar mone = -1.0;
29529e0b56fSBarry Smith 
2963a40ed3dSBarry Smith   PetscFunctionBegin;
29729e0b56fSBarry Smith   *flag = 0;
29829e0b56fSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
29929e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y); CHKERRQ(ierr);            /* y <- y - x      */
30029e0b56fSBarry Smith   ierr = SNESComputeFunction(snes,y,g); CHKERRQ(ierr); /* Compute F(y)    */
30129e0b56fSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
3023a40ed3dSBarry Smith   PetscFunctionReturn(0);
30329e0b56fSBarry Smith }
30404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
30529e0b56fSBarry Smith #undef __FUNC__
3065615d1e5SSatish Balay #define __FUNC__ "SNESCubicLineSearch"
3074b828684SBarry Smith /*@C
308f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
3095e76c431SBarry Smith 
310c7afd0dbSLois Curfman McInnes    Collective on SNES
311c7afd0dbSLois Curfman McInnes 
3125e76c431SBarry Smith    Input Parameters:
313c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
3145e76c431SBarry Smith .  x - current iterate
3155e76c431SBarry Smith .  f - residual evaluated at x
3165e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3175e76c431SBarry Smith .  w - work vector
318c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3195e76c431SBarry Smith 
320393d2d9aSLois Curfman McInnes    Output Parameters:
321c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3225e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3235e76c431SBarry Smith .  gnorm - 2-norm of g
3245e76c431SBarry Smith .  ynorm - 2-norm of search length
325c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
326fee21e36SBarry Smith 
327c4a48953SLois Curfman McInnes    Options Database Key:
328c7afd0dbSLois Curfman McInnes $  -snes_eq_ls cubic - Activates SNESCubicLineSearch()
329c4a48953SLois Curfman McInnes 
3305e76c431SBarry Smith    Notes:
3315e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
3325e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
3335e76c431SBarry Smith 
33436851e7fSLois Curfman McInnes    Level: advanced
33536851e7fSLois Curfman McInnes 
33628ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33728ae5a21SLois Curfman McInnes 
33877c4ece6SBarry Smith .seealso: SNESNoLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
3395e76c431SBarry Smith @*/
3405e42470aSBarry Smith int SNESCubicLineSearch(SNES snes,Vec x,Vec f,Vec g,Vec y,Vec w,
341761aaf1bSLois Curfman McInnes                         double fnorm,double *ynorm,double *gnorm,int *flag)
3425e76c431SBarry Smith {
343406556e6SLois Curfman McInnes   /*
344406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
345406556e6SLois Curfman McInnes      minimization problem:
346406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
347406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
348406556e6SLois Curfman McInnes    */
349406556e6SLois Curfman McInnes 
350ddd12767SBarry Smith   double  steptol, initslope, lambdaprev, gnormprev, a, b, d, t1, t2;
351ea4d3ed3SLois Curfman McInnes   double  maxstep, minlambda, alpha, lambda, lambdatemp, lambdaneg;
3523a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
3535e42470aSBarry Smith   Scalar  cinitslope, clambda;
3546b5873e3SBarry Smith #endif
3555e42470aSBarry Smith   int     ierr, count;
3565e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
3575334005bSBarry Smith   Scalar  mone = -1.0,scale;
3585e76c431SBarry Smith 
3593a40ed3dSBarry Smith   PetscFunctionBegin;
3607857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
361761aaf1bSLois Curfman McInnes   *flag   = 0;
3625e76c431SBarry Smith   alpha   = neP->alpha;
3635e76c431SBarry Smith   maxstep = neP->maxstep;
3645e76c431SBarry Smith   steptol = neP->steptol;
3655e76c431SBarry Smith 
366cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm); CHKERRQ(ierr);
367a1c2b6eeSLois Curfman McInnes   if (*ynorm < snes->atol) {
368a1c2b6eeSLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Search direction and size are nearly 0\n");
369a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
370a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
371a1c2b6eeSLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
372ad922ac9SBarry Smith     goto theend1;
37394a9d846SBarry Smith   }
3745e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
3755e42470aSBarry Smith     scale = maxstep/(*ynorm);
3763a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
37792318cfeSSatish Balay     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",PetscReal(scale));
3786b5873e3SBarry Smith #else
37994a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Scaling step by %g\n",scale);
3806b5873e3SBarry Smith #endif
381393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
3825e76c431SBarry Smith     *ynorm = maxstep;
3835e76c431SBarry Smith   }
3845e76c431SBarry Smith   minlambda = steptol/(*ynorm);
385a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
3863a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
387a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
38892318cfeSSatish Balay   initslope = PetscReal(cinitslope);
3895e42470aSBarry Smith #else
390a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
3915e42470aSBarry Smith #endif
3925e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
3935e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
3945e76c431SBarry Smith 
395393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
3965334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
39778b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
398cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);
399406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
400393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
40194a424c1SBarry Smith     PLogInfo(snes,"SNESCubicLineSearch: Using full step\n");
402ad922ac9SBarry Smith     goto theend1;
4035e76c431SBarry Smith   }
4045e76c431SBarry Smith 
4055e76c431SBarry Smith   /* Fit points with quadratic */
4065e76c431SBarry Smith   lambda = 1.0; count = 0;
407a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
4085e76c431SBarry Smith   lambdaprev = lambda;
4095e76c431SBarry Smith   gnormprev = *gnorm;
410ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
411ddd12767SBarry Smith   else lambda = lambdatemp;
412393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w); CHKERRQ(ierr);
413ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
4143a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
415ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4165e42470aSBarry Smith #else
417ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4185e42470aSBarry Smith #endif
41978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
420cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
421406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
422393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
423ea4d3ed3SLois Curfman McInnes     PLogInfo(snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%g\n",lambda);
424ad922ac9SBarry Smith     goto theend1;
4255e76c431SBarry Smith   }
4265e76c431SBarry Smith 
4275e76c431SBarry Smith   /* Fit points with cubic */
4285e76c431SBarry Smith   count = 1;
4295e76c431SBarry Smith   while (1) {
4305e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
4312b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:Unable to find good step length! %d \n",count);
4322b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
433ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
434393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
435761aaf1bSLois Curfman McInnes       *flag = -1; break;
4365e76c431SBarry Smith     }
437406556e6SLois Curfman McInnes     t1 = ((*gnorm)*(*gnorm) - fnorm*fnorm)*0.5 - lambda*initslope;
438406556e6SLois Curfman McInnes     t2 = (gnormprev*gnormprev  - fnorm*fnorm)*0.5 - lambdaprev*initslope;
439ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4402b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
4415e76c431SBarry Smith     d  = b*b - 3*a*initslope;
4425e76c431SBarry Smith     if (d < 0.0) d = 0.0;
4435e76c431SBarry Smith     if (a == 0.0) {
4445e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
4455e76c431SBarry Smith     } else {
4465e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
4475e76c431SBarry Smith     }
4485e76c431SBarry Smith     if (lambdatemp > .5*lambda) {
4495e76c431SBarry Smith       lambdatemp = .5*lambda;
4505e76c431SBarry Smith     }
4515e76c431SBarry Smith     lambdaprev = lambda;
4525e76c431SBarry Smith     gnormprev = *gnorm;
4535e76c431SBarry Smith     if (lambdatemp <= .1*lambda) {
4545e76c431SBarry Smith       lambda = .1*lambda;
4555e76c431SBarry Smith     }
4565e76c431SBarry Smith     else lambda = lambdatemp;
457393d2d9aSLois Curfman McInnes     ierr = VecCopy( x, w ); CHKERRQ(ierr);
458ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
4593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
460ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
461393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
4625e42470aSBarry Smith #else
463ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w); CHKERRQ(ierr);
4645e42470aSBarry Smith #endif
46578b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
466cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
467406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* is reduction enough? */
468393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
469ea4d3ed3SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubically determined step, lambda=%g\n",lambda);
4702715565aSLois Curfman McInnes       goto theend1;
4712b022350SLois Curfman McInnes     } else {
4722b022350SLois Curfman McInnes       PLogInfo(snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%g\n",lambda);
4735e76c431SBarry Smith     }
4745e76c431SBarry Smith     count++;
4755e76c431SBarry Smith   }
476ad922ac9SBarry Smith   theend1:
4777857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
4783a40ed3dSBarry Smith   PetscFunctionReturn(0);
4795e76c431SBarry Smith }
48004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4815615d1e5SSatish Balay #undef __FUNC__
4825615d1e5SSatish Balay #define __FUNC__ "SNESQuadraticLineSearch"
4834b828684SBarry Smith /*@C
484f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
4855e76c431SBarry Smith 
486c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
487c7afd0dbSLois Curfman McInnes 
4885e42470aSBarry Smith    Input Parameters:
489c7afd0dbSLois Curfman McInnes +  snes - the SNES context
4905e42470aSBarry Smith .  x - current iterate
4915e42470aSBarry Smith .  f - residual evaluated at x
4925e42470aSBarry Smith .  y - search direction (contains new iterate on output)
4935e42470aSBarry Smith .  w - work vector
494c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4955e42470aSBarry Smith 
496c4a48953SLois Curfman McInnes    Output Parameters:
497c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4985e42470aSBarry Smith .  y - new iterate (contains search direction on input)
4995e42470aSBarry Smith .  gnorm - 2-norm of g
5005e42470aSBarry Smith .  ynorm - 2-norm of search length
501c7afd0dbSLois Curfman McInnes -  flag - 0 if line search succeeds; -1 on failure.
502fee21e36SBarry Smith 
503c4a48953SLois Curfman McInnes    Options Database Key:
504c7afd0dbSLois Curfman McInnes .  -snes_eq_ls quadratic - Activates SNESQuadraticLineSearch()
505c4a48953SLois Curfman McInnes 
5065e42470aSBarry Smith    Notes:
507c7afd0dbSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNES_EQ_LS method.
5085e42470aSBarry Smith 
50936851e7fSLois Curfman McInnes    Level: advanced
51036851e7fSLois Curfman McInnes 
511f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
512f59ffdedSLois Curfman McInnes 
51377c4ece6SBarry Smith .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch()
5145e42470aSBarry Smith @*/
5155e42470aSBarry Smith int SNESQuadraticLineSearch(SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
516761aaf1bSLois Curfman McInnes                            double fnorm, double *ynorm, double *gnorm,int *flag)
5175e76c431SBarry Smith {
518406556e6SLois Curfman McInnes   /*
519406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
520406556e6SLois Curfman McInnes      minimization problem:
521406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
522406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
523406556e6SLois Curfman McInnes    */
52440011551SBarry Smith   double  steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp;
5253a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
5265e42470aSBarry Smith   Scalar  cinitslope,clambda;
5276b5873e3SBarry Smith #endif
5285e42470aSBarry Smith   int     ierr,count;
5295e42470aSBarry Smith   SNES_LS *neP = (SNES_LS *) snes->data;
5305334005bSBarry Smith   Scalar  mone = -1.0,scale;
5315e76c431SBarry Smith 
5323a40ed3dSBarry Smith   PetscFunctionBegin;
5337857610eSBarry Smith   PLogEventBegin(SNES_LineSearch,snes,x,f,g);
534761aaf1bSLois Curfman McInnes   *flag   = 0;
5355e42470aSBarry Smith   alpha   = neP->alpha;
5365e42470aSBarry Smith   maxstep = neP->maxstep;
5375e42470aSBarry Smith   steptol = neP->steptol;
5385e76c431SBarry Smith 
539cddf8d76SBarry Smith   VecNorm(y, NORM_2,ynorm );
540b37302c6SLois Curfman McInnes   if (*ynorm < snes->atol) {
54194a9d846SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Search direction and size is 0\n");
542b37302c6SLois Curfman McInnes     *gnorm = fnorm;
543b37302c6SLois Curfman McInnes     ierr = VecCopy(x,y); CHKERRQ(ierr);
544b37302c6SLois Curfman McInnes     ierr = VecCopy(f,g); CHKERRQ(ierr);
545ad922ac9SBarry Smith     goto theend2;
54694a9d846SBarry Smith   }
5475e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5485e42470aSBarry Smith     scale = maxstep/(*ynorm);
549393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y); CHKERRQ(ierr);
5505e42470aSBarry Smith     *ynorm = maxstep;
5515e76c431SBarry Smith   }
5525e42470aSBarry Smith   minlambda = steptol/(*ynorm);
553a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w); CHKERRQ(ierr);
5543a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
555a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope); CHKERRQ(ierr);
55692318cfeSSatish Balay   initslope = PetscReal(cinitslope);
5575e42470aSBarry Smith #else
558a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope); CHKERRQ(ierr);
5595e42470aSBarry Smith #endif
5605e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
5615e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
5625e42470aSBarry Smith 
563393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w); CHKERRQ(ierr);
5645334005bSBarry Smith   ierr = VecAYPX(&mone,x,w); CHKERRQ(ierr);
56578b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
566cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
567406556e6SLois Curfman McInnes   if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* Sufficient reduction */
568393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y); CHKERRQ(ierr);
56994a424c1SBarry Smith     PLogInfo(snes,"SNESQuadraticLineSearch: Using full step\n");
570ad922ac9SBarry Smith     goto theend2;
5715e42470aSBarry Smith   }
5725e42470aSBarry Smith 
5735e42470aSBarry Smith   /* Fit points with quadratic */
5745e42470aSBarry Smith   lambda = 1.0; count = 0;
5755e42470aSBarry Smith   count = 1;
5765e42470aSBarry Smith   while (1) {
5775e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
578981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Unable to find good step length! %d \n",count);
579981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",
580ea4d3ed3SLois Curfman McInnes                fnorm,*gnorm,*ynorm,lambda,initslope);
581393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
582761aaf1bSLois Curfman McInnes       *flag = -1; break;
5835e42470aSBarry Smith     }
584a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5855e42470aSBarry Smith     if (lambdatemp <= .1*lambda) {
5865e42470aSBarry Smith       lambda = .1*lambda;
5875e42470aSBarry Smith     } else lambda = lambdatemp;
588393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w); CHKERRQ(ierr);
5895334005bSBarry Smith     lambda = -lambda;
5903a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
591393d2d9aSLois Curfman McInnes     clambda = lambda; ierr = VecAXPY(&clambda,y,w); CHKERRQ(ierr);
5925e42470aSBarry Smith #else
593393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&lambda,y,w); CHKERRQ(ierr);
5945e42470aSBarry Smith #endif
59578b31e54SBarry Smith     ierr = SNESComputeFunction(snes,w,g); CHKERRQ(ierr);
596cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm); CHKERRQ(ierr);
597406556e6SLois Curfman McInnes     if ((*gnorm)*(*gnorm)*0.5 <= fnorm*fnorm*0.5 + alpha*initslope) { /* sufficient reduction */
598393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y); CHKERRQ(ierr);
599981c4779SBarry Smith       PLogInfo(snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda);
60006259719SBarry Smith       break;
6015e42470aSBarry Smith     }
6025e42470aSBarry Smith     count++;
6035e42470aSBarry Smith   }
604ad922ac9SBarry Smith   theend2:
6057857610eSBarry Smith   PLogEventEnd(SNES_LineSearch,snes,x,f,g);
6063a40ed3dSBarry Smith   PetscFunctionReturn(0);
6075e76c431SBarry Smith }
60804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6095615d1e5SSatish Balay #undef __FUNC__
610d4bb536fSBarry Smith #define __FUNC__ "SNESSetLineSearch"
611c9e6a524SLois Curfman McInnes /*@C
61277c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
613f63b844aSLois Curfman McInnes    by the method SNES_EQ_LS.
6145e76c431SBarry Smith 
6155e76c431SBarry Smith    Input Parameters:
616c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
617c7afd0dbSLois Curfman McInnes -  func - pointer to int function
6185e76c431SBarry Smith 
619fee21e36SBarry Smith    Collective on SNES
620fee21e36SBarry Smith 
621c4a48953SLois Curfman McInnes    Available Routines:
622c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
623f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
624f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
625c7afd0dbSLois Curfman McInnes -  SNESNoLineSearchNoNorms() - use the full Newton step (calculate no norms; faster in parallel)
6265e76c431SBarry Smith 
627c4a48953SLois Curfman McInnes     Options Database Keys:
628c7afd0dbSLois Curfman McInnes +   -snes_eq_ls [basic,quadratic,cubic] - Selects line search
629c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_alpha <alpha> - Sets alpha
630c7afd0dbSLois Curfman McInnes .   -snes_eq_ls_maxstep <max> - Sets maxstep
631c7afd0dbSLois Curfman McInnes -   -snes_eq_ls_steptol <steptol> - Sets steptol
632c4a48953SLois Curfman McInnes 
6335e76c431SBarry Smith    Calling sequence of func:
634bd208895SLois Curfman McInnes .vb
635bd208895SLois Curfman McInnes    func (SNES snes, Vec x, Vec f, Vec g, Vec y, Vec w,
636bd208895SLois Curfman McInnes          double fnorm, double *ynorm, double *gnorm, *flag)
637bd208895SLois Curfman McInnes .ve
6385e76c431SBarry Smith 
6395e76c431SBarry Smith     Input parameters for func:
640c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
6415e76c431SBarry Smith .   x - current iterate
6425e76c431SBarry Smith .   f - residual evaluated at x
6435e76c431SBarry Smith .   y - search direction (contains new iterate on output)
6445e76c431SBarry Smith .   w - work vector
645c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
6465e76c431SBarry Smith 
6475e76c431SBarry Smith     Output parameters for func:
648c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
6495e76c431SBarry Smith .   y - new iterate (contains search direction on input)
6505e76c431SBarry Smith .   gnorm - 2-norm of g
6515e76c431SBarry Smith .   ynorm - 2-norm of search length
652c7afd0dbSLois Curfman McInnes -   flag - set to 0 if the line search succeeds; a nonzero integer
653761aaf1bSLois Curfman McInnes            on failure.
654f59ffdedSLois Curfman McInnes 
65536851e7fSLois Curfman McInnes     Level: advanced
65636851e7fSLois Curfman McInnes 
657f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
658f59ffdedSLois Curfman McInnes 
659f59ffdedSLois Curfman McInnes .seealso: SNESNoLineSearch(), SNESQuadraticLineSearch(), SNESCubicLineSearch()
6605e76c431SBarry Smith @*/
66177c4ece6SBarry Smith int SNESSetLineSearch(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
662761aaf1bSLois Curfman McInnes                              double,double*,double*,int*))
6635e76c431SBarry Smith {
664415aef20SBarry Smith   int ierr, (*f)(SNES,int (*)(SNES,Vec,Vec,Vec,Vec,Vec,double,double*,double*,int*));
66582bf6240SBarry Smith 
6663a40ed3dSBarry Smith   PetscFunctionBegin;
66794d884c6SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void **)&f);CHKERRQ(ierr);
66882bf6240SBarry Smith   if (f) {
66982bf6240SBarry Smith     ierr = (*f)(snes,func);CHKERRQ(ierr);
67082bf6240SBarry Smith   }
6713a40ed3dSBarry Smith   PetscFunctionReturn(0);
6725e76c431SBarry Smith }
67304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
674fb2e594dSBarry Smith EXTERN_C_BEGIN
67582bf6240SBarry Smith #undef __FUNC__
67682bf6240SBarry Smith #define __FUNC__ "SNESSetLineSearch_LS"
67782bf6240SBarry Smith int SNESSetLineSearch_LS(SNES snes,int (*func)(SNES,Vec,Vec,Vec,Vec,Vec,
67882bf6240SBarry Smith                          double,double*,double*,int*))
67982bf6240SBarry Smith {
68082bf6240SBarry Smith   PetscFunctionBegin;
68182bf6240SBarry Smith   ((SNES_LS *)(snes->data))->LineSearch = func;
68282bf6240SBarry Smith   PetscFunctionReturn(0);
68382bf6240SBarry Smith }
684fb2e594dSBarry Smith EXTERN_C_END
68504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
68604d965bbSLois Curfman McInnes /*
68704d965bbSLois Curfman McInnes    SNESPrintHelp_EQ_LS - Prints all options for the SNES_EQ_LS method.
68882bf6240SBarry Smith 
68904d965bbSLois Curfman McInnes    Input Parameter:
69004d965bbSLois Curfman McInnes .  snes - the SNES context
69104d965bbSLois Curfman McInnes 
69204d965bbSLois Curfman McInnes    Application Interface Routine: SNESPrintHelp()
69304d965bbSLois Curfman McInnes */
6945615d1e5SSatish Balay #undef __FUNC__
695d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp_EQ_LS"
696f63b844aSLois Curfman McInnes static int SNESPrintHelp_EQ_LS(SNES snes,char *p)
6975e42470aSBarry Smith {
6985e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
6996daaf66cSBarry Smith 
7003a40ed3dSBarry Smith   PetscFunctionBegin;
70176be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm," method SNES_EQ_LS (ls) for systems of nonlinear equations:\n");
70276be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls [basic,quadratic,cubic]\n",p);
70376be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
70476be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
70576be9ce4SBarry Smith   (*PetscHelpPrintf)(snes->comm,"   %ssnes_eq_ls_steptol <tol> (default %g)\n",p,ls->steptol);
7063a40ed3dSBarry Smith   PetscFunctionReturn(0);
7075e42470aSBarry Smith }
70804d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
70904d965bbSLois Curfman McInnes /*
710de49b36dSLois Curfman McInnes    SNESView_EQ_LS - Prints info from the SNES_EQ_LS data structure.
71104d965bbSLois Curfman McInnes 
71204d965bbSLois Curfman McInnes    Input Parameters:
71304d965bbSLois Curfman McInnes .  SNES - the SNES context
71404d965bbSLois Curfman McInnes .  viewer - visualization context
71504d965bbSLois Curfman McInnes 
71604d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
71704d965bbSLois Curfman McInnes */
7185615d1e5SSatish Balay #undef __FUNC__
719d4bb536fSBarry Smith #define __FUNC__ "SNESView_EQ_LS"
720e1311b90SBarry Smith static int SNESView_EQ_LS(SNES snes,Viewer viewer)
721a935fc98SLois Curfman McInnes {
722a935fc98SLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
72319bcc07fSBarry Smith   char       *cstr;
72451695f54SLois Curfman McInnes   int        ierr;
72519bcc07fSBarry Smith   ViewerType vtype;
726a935fc98SLois Curfman McInnes 
7273a40ed3dSBarry Smith   PetscFunctionBegin;
72819bcc07fSBarry Smith   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
7293f1db9ecSBarry Smith   if (PetscTypeCompare(vtype,ASCII_VIEWER)) {
73019bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
73119bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
73219bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
73319bcc07fSBarry Smith     else                                                cstr = "unknown";
7340ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);
7350ef38995SBarry Smith     ViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);
7365cd90555SBarry Smith   } else {
7375cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
73819bcc07fSBarry Smith   }
7393a40ed3dSBarry Smith   PetscFunctionReturn(0);
740a935fc98SLois Curfman McInnes }
74104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
74204d965bbSLois Curfman McInnes /*
74304d965bbSLois Curfman McInnes    SNESSetFromOptions_EQ_LS - Sets various parameters for the SNES_EQ_LS method.
74404d965bbSLois Curfman McInnes 
74504d965bbSLois Curfman McInnes    Input Parameter:
74604d965bbSLois Curfman McInnes .  snes - the SNES context
74704d965bbSLois Curfman McInnes 
74804d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
74904d965bbSLois Curfman McInnes */
7505615d1e5SSatish Balay #undef __FUNC__
7515615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions_EQ_LS"
752f63b844aSLois Curfman McInnes static int SNESSetFromOptions_EQ_LS(SNES snes)
7535e42470aSBarry Smith {
7545e42470aSBarry Smith   SNES_LS *ls = (SNES_LS *)snes->data;
7555e42470aSBarry Smith   char    ver[16];
7565e42470aSBarry Smith   double  tmp;
75725cdf11fSBarry Smith   int     ierr,flg;
7585e42470aSBarry Smith 
7593a40ed3dSBarry Smith   PetscFunctionBegin;
76009d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_alpha",&tmp, &flg);CHKERRQ(ierr);
76125cdf11fSBarry Smith   if (flg) {
7625e42470aSBarry Smith     ls->alpha = tmp;
7635e42470aSBarry Smith   }
76409d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_maxstep",&tmp, &flg);CHKERRQ(ierr);
76525cdf11fSBarry Smith   if (flg) {
7665e42470aSBarry Smith     ls->maxstep = tmp;
7675e42470aSBarry Smith   }
76809d61ba7SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_eq_ls_steptol",&tmp, &flg);CHKERRQ(ierr);
76925cdf11fSBarry Smith   if (flg) {
7705e42470aSBarry Smith     ls->steptol = tmp;
7715e42470aSBarry Smith   }
77209d61ba7SLois Curfman McInnes   ierr = OptionsGetString(snes->prefix,"-snes_eq_ls",ver,16, &flg); CHKERRQ(ierr);
77325cdf11fSBarry Smith   if (flg) {
77448d91487SBarry Smith     if (!PetscStrcmp(ver,"basic")) {
775b2d88d52SBarry Smith       ierr = SNESSetLineSearch(snes,SNESNoLineSearch);CHKERRQ(ierr);
776b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"basicnonorms")) {
777b2d88d52SBarry Smith       ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms);CHKERRQ(ierr);
778b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"quadratic")) {
779b2d88d52SBarry Smith       ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch);CHKERRQ(ierr);
780b2d88d52SBarry Smith     } else if (!PetscStrcmp(ver,"cubic")) {
781b2d88d52SBarry Smith       ierr = SNESSetLineSearch(snes,SNESCubicLineSearch);CHKERRQ(ierr);
7825e42470aSBarry Smith     }
783a8c6a408SBarry Smith     else {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown line search");}
7845e42470aSBarry Smith   }
7853a40ed3dSBarry Smith   PetscFunctionReturn(0);
7865e42470aSBarry Smith }
78704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
78804d965bbSLois Curfman McInnes /*
78904d965bbSLois Curfman McInnes    SNESCreate_EQ_LS - Creates a nonlinear solver context for the SNES_EQ_LS method,
79004d965bbSLois Curfman McInnes    SNES_LS, and sets this as the private data within the generic nonlinear solver
79104d965bbSLois Curfman McInnes    context, SNES, that was created within SNESCreate().
79204d965bbSLois Curfman McInnes 
79304d965bbSLois Curfman McInnes 
79404d965bbSLois Curfman McInnes    Input Parameter:
79504d965bbSLois Curfman McInnes .  snes - the SNES context
79604d965bbSLois Curfman McInnes 
79704d965bbSLois Curfman McInnes    Application Interface Routine: SNESCreate()
79804d965bbSLois Curfman McInnes  */
799fb2e594dSBarry Smith EXTERN_C_BEGIN
8005615d1e5SSatish Balay #undef __FUNC__
8015615d1e5SSatish Balay #define __FUNC__ "SNESCreate_EQ_LS"
802f63b844aSLois Curfman McInnes int SNESCreate_EQ_LS(SNES snes)
8035e42470aSBarry Smith {
80482bf6240SBarry Smith   int     ierr;
8055e42470aSBarry Smith   SNES_LS *neP;
8065e42470aSBarry Smith 
8073a40ed3dSBarry Smith   PetscFunctionBegin;
808a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
809a8c6a408SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only");
810a8c6a408SBarry Smith   }
81182bf6240SBarry Smith 
812f63b844aSLois Curfman McInnes   snes->setup		= SNESSetUp_EQ_LS;
813f63b844aSLois Curfman McInnes   snes->solve		= SNESSolve_EQ_LS;
814f63b844aSLois Curfman McInnes   snes->destroy		= SNESDestroy_EQ_LS;
815f63b844aSLois Curfman McInnes   snes->converged	= SNESConverged_EQ_LS;
816f63b844aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_EQ_LS;
817f63b844aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_EQ_LS;
818f63b844aSLois Curfman McInnes   snes->view            = SNESView_EQ_LS;
8195baf8537SBarry Smith   snes->nwork           = 0;
8205e42470aSBarry Smith 
8210452661fSBarry Smith   neP			= PetscNew(SNES_LS);   CHKPTRQ(neP);
822ff782a9fSLois Curfman McInnes   PLogObjectMemory(snes,sizeof(SNES_LS));
8235e42470aSBarry Smith   snes->data    	= (void *) neP;
8245e42470aSBarry Smith   neP->alpha		= 1.e-4;
8255e42470aSBarry Smith   neP->maxstep		= 1.e8;
8265e42470aSBarry Smith   neP->steptol		= 1.e-12;
8275e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
82882bf6240SBarry Smith 
82994d884c6SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",
830e1311b90SBarry Smith                     (void*)SNESSetLineSearch_LS);CHKERRQ(ierr);
83182bf6240SBarry Smith 
8323a40ed3dSBarry Smith   PetscFunctionReturn(0);
8335e42470aSBarry Smith }
834fb2e594dSBarry Smith EXTERN_C_END
83582bf6240SBarry Smith 
83682bf6240SBarry Smith 
83782bf6240SBarry Smith 
838