xref: /petsc/src/snes/impls/ls/ls.c (revision 63dd3a1af947c5d07342e150a17d503f381a847e)
1*63dd3a1aSKris Buschelman #define PETSCSNES_DLL
25e76c431SBarry Smith 
370f55243SBarry Smith #include "src/snes/impls/ls/ls.h"
45e42470aSBarry Smith 
58a5d9ee4SBarry Smith /*
68a5d9ee4SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the function,
78a5d9ee4SBarry Smith     but not a zero. In the case when one cannot compute J^T F we use the fact that
836669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
936669109SBarry Smith     for this trick.
108a5d9ee4SBarry Smith */
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
13dfbe8321SBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(Mat A,Vec F,Vec W,PetscReal fnorm,PetscTruth *ismin)
148a5d9ee4SBarry Smith {
158a5d9ee4SBarry Smith   PetscReal      a1;
16dfbe8321SBarry Smith   PetscErrorCode ierr;
1736669109SBarry Smith   PetscTruth     hastranspose;
188a5d9ee4SBarry Smith 
198a5d9ee4SBarry Smith   PetscFunctionBegin;
208a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2136669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2236669109SBarry Smith   if (hastranspose) {
238a5d9ee4SBarry Smith     /* Compute || J^T F|| */
248a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
258a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
2663ba0a88SBarry Smith     ierr = PetscLogInfo((0,"SNESSolve_LS: || J^T F|| %g near zero implies found a local minimum\n",a1/fnorm));CHKERRQ(ierr);
278a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2836669109SBarry Smith   } else {
2936669109SBarry Smith     Vec         work;
30ea709b57SSatish Balay     PetscScalar result;
3136669109SBarry Smith     PetscReal   wnorm;
3236669109SBarry Smith 
3336669109SBarry Smith     ierr = VecSetRandom(PETSC_NULL,W);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3736669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
3836669109SBarry Smith     ierr = VecDestroy(work);CHKERRQ(ierr);
3936669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
4063ba0a88SBarry Smith     ierr = PetscLogInfo((0,"SNESSolve_LS: (F^T J random)/(|| F ||*||random|| %g near zero implies found a local minimum\n",a1));CHKERRQ(ierr);
4136669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4236669109SBarry Smith   }
438a5d9ee4SBarry Smith   PetscFunctionReturn(0);
448a5d9ee4SBarry Smith }
458a5d9ee4SBarry Smith 
4674637425SBarry Smith /*
475ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4874637425SBarry Smith */
494a2ae208SSatish Balay #undef __FUNCT__
504a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
51dfbe8321SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(Mat A,Vec F,Vec X,Vec W1,Vec W2)
5274637425SBarry Smith {
5374637425SBarry Smith   PetscReal     a1,a2;
54dfbe8321SBarry Smith   PetscErrorCode ierr;
5574637425SBarry Smith   PetscTruth    hastranspose;
56ea709b57SSatish Balay   PetscScalar   mone = -1.0;
5774637425SBarry Smith 
5874637425SBarry Smith   PetscFunctionBegin;
5974637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
6074637425SBarry Smith   if (hastranspose) {
6174637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6274637425SBarry Smith     ierr = VecAXPY(&mone,F,W1);CHKERRQ(ierr);
6374637425SBarry Smith 
6474637425SBarry Smith     /* Compute || J^T W|| */
6574637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6674637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6774637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
689994e62eSSatish Balay     if (a1 != 0) {
6963ba0a88SBarry Smith       ierr = PetscLogInfo((0,"SNESSolve_LS: ||J^T(F-Ax)||/||F-AX|| %g near zero implies inconsistent rhs\n",a2/a1));CHKERRQ(ierr);
7085471664SBarry Smith     }
7174637425SBarry Smith   }
7274637425SBarry Smith   PetscFunctionReturn(0);
7374637425SBarry Smith }
7474637425SBarry Smith 
7504d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7604d965bbSLois Curfman McInnes 
7704d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7894b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7904d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8004d965bbSLois Curfman McInnes      respectively.
8104d965bbSLois Curfman McInnes 
82fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8304d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8404d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
85fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8604d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8704d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
884b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
894b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
9004d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9104d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9204d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9304d965bbSLois Curfman McInnes      for all nonlinear solvers.
9404d965bbSLois Curfman McInnes 
9504d965bbSLois Curfman McInnes      Another key routine is:
9604d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9704d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9804d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9904d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10004d965bbSLois Curfman McInnes 
10104d965bbSLois Curfman McInnes      Additional basic routines are:
10204d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10304d965bbSLois Curfman McInnes                                       have actually been used.
10404d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
105186905e3SBarry Smith      SNESView().
10604d965bbSLois Curfman McInnes 
10704d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10804d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10904d965bbSLois Curfman McInnes      above description applies to these categories also.
11004d965bbSLois Curfman McInnes 
11104d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1125e42470aSBarry Smith /*
1134b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11404d965bbSLois Curfman McInnes    method with a line search.
1155e76c431SBarry Smith 
11604d965bbSLois Curfman McInnes    Input Parameters:
11704d965bbSLois Curfman McInnes .  snes - the SNES context
1185e76c431SBarry Smith 
11904d965bbSLois Curfman McInnes    Output Parameter:
120c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12104d965bbSLois Curfman McInnes 
12204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1235e76c431SBarry Smith 
1245e76c431SBarry Smith    Notes:
1255e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1265e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1275e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1285e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
129393d2d9aSLois Curfman McInnes    and Schnabel.
1305e76c431SBarry Smith */
1314a2ae208SSatish Balay #undef __FUNCT__
1324b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
133dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1345e76c431SBarry Smith {
1354b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
1366849ba73SBarry Smith   PetscErrorCode ierr;
137a7cc72afSBarry Smith   PetscInt       maxits,i,lits;
138a7cc72afSBarry Smith   PetscTruth     lssucceed;
139112a2221SBarry Smith   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
140329f5518SBarry Smith   PetscReal      fnorm,gnorm,xnorm,ynorm;
1415e42470aSBarry Smith   Vec            Y,X,F,G,W,TMP;
142c293cc10SBarry Smith   KSP            ksp;
1435e76c431SBarry Smith 
1443a40ed3dSBarry Smith   PetscFunctionBegin;
14594b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
146184914b5SBarry Smith   snes->reason  = SNES_CONVERGED_ITERATING;
147184914b5SBarry Smith 
1485e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1495e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
15039e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1515e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1525e42470aSBarry Smith   G		= snes->work[1];
1535e42470aSBarry Smith   W		= snes->work[2];
1545e76c431SBarry Smith 
1554c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1564c49b128SBarry Smith   snes->iter = 0;
1574c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1585334005bSBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);  /*  F(X)      */
159cddf8d76SBarry Smith   ierr = VecNorm(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
16043e71028SBarry Smith   if (fnorm != fnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1610f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1625e42470aSBarry Smith   snes->norm = fnorm;
1630f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
16484c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
16594a424c1SBarry Smith   SNESMonitor(snes,0,fnorm);
1665e76c431SBarry Smith 
16770441072SBarry Smith   if (fnorm < snes->abstol) {snes->reason = SNES_CONVERGED_FNORM_ABS; PetscFunctionReturn(0);}
16894a9d846SBarry Smith 
169d034289bSBarry Smith   /* set parameter for default relative tolerance convergence test */
170d034289bSBarry Smith   snes->ttol = fnorm*snes->rtol;
171d034289bSBarry Smith 
1725e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1735e76c431SBarry Smith 
174329e7e40SMatthew Knepley     /* Call general purpose update function */
175abc0a331SBarry Smith     if (snes->update) {
176329e7e40SMatthew Knepley       ierr = (*snes->update)(snes, snes->iter);CHKERRQ(ierr);
177de8cb200SMatthew Knepley     }
178329e7e40SMatthew Knepley 
179ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1805334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18194b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18223ce1328SBarry Smith     ierr = KSPSolve(snes->ksp,F,Y);CHKERRQ(ierr);
183c293cc10SBarry Smith     ierr = KSPGetIterationNumber(ksp,&lits);CHKERRQ(ierr);
18474637425SBarry Smith 
185b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
18674637425SBarry Smith       ierr = SNESLSCheckResidual_Private(snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
18785471664SBarry Smith     }
18874637425SBarry Smith 
18990cbd584SBarry Smith     /* should check what happened to the linear solve? */
1903505fcc1SBarry Smith     snes->linear_its += lits;
19163ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: iter=%D, linear solve iterations=%D\n",snes->iter,lits));CHKERRQ(ierr);
192ea4d3ed3SLois Curfman McInnes 
193ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
194ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
195ea4d3ed3SLois Curfman McInnes        and evaluate G(Y) = function(Y))
196ea4d3ed3SLois Curfman McInnes     */
19781b6cf68SLois Curfman McInnes     ierr = VecCopy(Y,snes->vec_sol_update_always);CHKERRQ(ierr);
198a7cc72afSBarry Smith     ierr = (*neP->LineSearch)(snes,neP->lsP,X,F,G,Y,W,fnorm,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
19963ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",fnorm,gnorm,ynorm,(int)lssucceed));CHKERRQ(ierr);
20043e71028SBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
201a3b891d8SBarry Smith 
202a3b891d8SBarry Smith     TMP = F; F = G; snes->vec_func_always = F; G = TMP;
203a3b891d8SBarry Smith     TMP = X; X = Y; snes->vec_sol_always = X;  Y = TMP;
204a3b891d8SBarry Smith     fnorm = gnorm;
205a3b891d8SBarry Smith 
2065ed2d596SBarry Smith     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
2075ed2d596SBarry Smith     snes->iter = i+1;
2085ed2d596SBarry Smith     snes->norm = fnorm;
2095ed2d596SBarry Smith     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
2105ed2d596SBarry Smith     SNESLogConvHistory(snes,fnorm,lits);
2115ed2d596SBarry Smith     SNESMonitor(snes,i+1,fnorm);
2125ed2d596SBarry Smith 
213a7cc72afSBarry Smith     if (!lssucceed) {
2148a5d9ee4SBarry Smith       PetscTruth ismin;
21550ffb88aSMatthew Knepley 
21650ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
2173505fcc1SBarry Smith         snes->reason = SNES_DIVERGED_LS_FAILURE;
2188a5d9ee4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes->jacobian,F,W,fnorm,&ismin);CHKERRQ(ierr);
2198a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2203505fcc1SBarry Smith         break;
2213505fcc1SBarry Smith       }
22250ffb88aSMatthew Knepley     }
2235e76c431SBarry Smith 
2245e76c431SBarry Smith     /* Test for convergence */
22529e0b56fSBarry Smith     if (snes->converged) {
22629e0b56fSBarry Smith       ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm = || X || */
2273505fcc1SBarry Smith       ierr = (*snes->converged)(snes,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2283505fcc1SBarry Smith       if (snes->reason) {
22916c95cacSBarry Smith         break;
23016c95cacSBarry Smith       }
23116c95cacSBarry Smith     }
23229e0b56fSBarry Smith   }
23339e2f89bSBarry Smith   if (X != snes->vec_sol) {
234393d2d9aSLois Curfman McInnes     ierr = VecCopy(X,snes->vec_sol);CHKERRQ(ierr);
235e15f7bb5SBarry Smith   }
236e15f7bb5SBarry Smith   if (F != snes->vec_func) {
237e15f7bb5SBarry Smith     ierr = VecCopy(F,snes->vec_func);CHKERRQ(ierr);
238e15f7bb5SBarry Smith   }
23939e2f89bSBarry Smith   snes->vec_sol_always  = snes->vec_sol;
24039e2f89bSBarry Smith   snes->vec_func_always = snes->vec_func;
24152392280SLois Curfman McInnes   if (i == maxits) {
24263ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESSolve_LS: Maximum number of iterations has been reached: %D\n",maxits));CHKERRQ(ierr);
2433505fcc1SBarry Smith     snes->reason = SNES_DIVERGED_MAX_IT;
24452392280SLois Curfman McInnes   }
2453a40ed3dSBarry Smith   PetscFunctionReturn(0);
2465e76c431SBarry Smith }
24704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
24804d965bbSLois Curfman McInnes /*
2494b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2504b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25104d965bbSLois Curfman McInnes 
25204d965bbSLois Curfman McInnes    Input Parameter:
25304d965bbSLois Curfman McInnes .  snes - the SNES context
25404d965bbSLois Curfman McInnes .  x - the solution vector
25504d965bbSLois Curfman McInnes 
25604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
25704d965bbSLois Curfman McInnes 
25804d965bbSLois Curfman McInnes    Notes:
2594b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26004d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26104d965bbSLois Curfman McInnes    the call to SNESSolve().
26204d965bbSLois Curfman McInnes  */
2634a2ae208SSatish Balay #undef __FUNCT__
2644b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
265dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2665e76c431SBarry Smith {
267dfbe8321SBarry Smith   PetscErrorCode ierr;
2683a40ed3dSBarry Smith 
2693a40ed3dSBarry Smith   PetscFunctionBegin;
27081b6cf68SLois Curfman McInnes   snes->nwork = 4;
271d7e8b826SBarry Smith   ierr = VecDuplicateVecs(snes->vec_sol,snes->nwork,&snes->work);CHKERRQ(ierr);
272efee365bSSatish Balay   ierr = PetscLogObjectParents(snes,snes->nwork,snes->work);CHKERRQ(ierr);
27381b6cf68SLois Curfman McInnes   snes->vec_sol_update_always = snes->work[3];
2743a40ed3dSBarry Smith   PetscFunctionReturn(0);
2755e76c431SBarry Smith }
27604d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
27704d965bbSLois Curfman McInnes /*
2784b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2794b27c08aSLois Curfman McInnes    with SNESCreate_LS().
28004d965bbSLois Curfman McInnes 
28104d965bbSLois Curfman McInnes    Input Parameter:
28204d965bbSLois Curfman McInnes .  snes - the SNES context
28304d965bbSLois Curfman McInnes 
284de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
28504d965bbSLois Curfman McInnes  */
2864a2ae208SSatish Balay #undef __FUNCT__
2874b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
288dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
2895e76c431SBarry Smith {
290dfbe8321SBarry Smith   PetscErrorCode ierr;
2913a40ed3dSBarry Smith 
2923a40ed3dSBarry Smith   PetscFunctionBegin;
2935baf8537SBarry Smith   if (snes->nwork) {
2944b0e389bSBarry Smith     ierr = VecDestroyVecs(snes->work,snes->nwork);CHKERRQ(ierr);
29521c89e3eSBarry Smith   }
2965c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
2973a40ed3dSBarry Smith   PetscFunctionReturn(0);
2985e76c431SBarry Smith }
29904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3004a2ae208SSatish Balay #undef __FUNCT__
3014a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearch"
30282bf6240SBarry Smith 
3034b828684SBarry Smith /*@C
3045e42470aSBarry Smith    SNESNoLineSearch - This routine is not a line search at all;
3055e76c431SBarry Smith    it simply uses the full Newton step.  Thus, this routine is intended
3065e76c431SBarry Smith    to serve as a template and is not recommended for general use.
3075e76c431SBarry Smith 
308c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
309c7afd0dbSLois Curfman McInnes 
3105e76c431SBarry Smith    Input Parameters:
311c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
312af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
3135e76c431SBarry Smith .  x - current iterate
3145e76c431SBarry Smith .  f - residual evaluated at x
3155e76c431SBarry Smith .  y - search direction (contains new iterate on output)
3165e76c431SBarry Smith .  w - work vector
317c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
3185e76c431SBarry Smith 
319c4a48953SLois Curfman McInnes    Output Parameters:
320c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
3215e76c431SBarry Smith .  y - new iterate (contains search direction on input)
3225e76c431SBarry Smith .  gnorm - 2-norm of g
3235e76c431SBarry Smith .  ynorm - 2-norm of search length
324a7cc72afSBarry Smith -  flag - PETSC_TRUE on success, PETSC_FALSE on failure
325fee21e36SBarry Smith 
326c4a48953SLois Curfman McInnes    Options Database Key:
3274b27c08aSLois Curfman McInnes .  -snes_ls basic - Activates SNESNoLineSearch()
328c4a48953SLois Curfman McInnes 
32936851e7fSLois Curfman McInnes    Level: advanced
33036851e7fSLois Curfman McInnes 
33128ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
33228ae5a21SLois Curfman McInnes 
333f59ffdedSLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
334af2d14d2SLois Curfman McInnes           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
3355e76c431SBarry Smith @*/
336*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESNoLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
3375e76c431SBarry Smith {
338dfbe8321SBarry Smith   PetscErrorCode ierr;
339ea709b57SSatish Balay   PetscScalar    mone = -1.0;
3404b27c08aSLois Curfman McInnes   SNES_LS        *neP = (SNES_LS*)snes->data;
3418f99978dSLois Curfman McInnes   PetscTruth     change_y = PETSC_FALSE;
3425334005bSBarry Smith 
3433a40ed3dSBarry Smith   PetscFunctionBegin;
344a7cc72afSBarry Smith   *flag = PETSC_TRUE;
345d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
346a3b38805SLois Curfman McInnes   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);  /* ynorm = || y || */
347ea4d3ed3SLois Curfman McInnes   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
3488f99978dSLois Curfman McInnes   if (neP->CheckStep) {
3498f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
3508f99978dSLois Curfman McInnes   }
351ea4d3ed3SLois Curfman McInnes   ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
352a3b38805SLois Curfman McInnes   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);  /* gnorm = || g || */
35343e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
354d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
3553a40ed3dSBarry Smith   PetscFunctionReturn(0);
3565e76c431SBarry Smith }
35704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
3584a2ae208SSatish Balay #undef __FUNCT__
3594a2ae208SSatish Balay #define __FUNCT__ "SNESNoLineSearchNoNorms"
36082bf6240SBarry Smith 
36129e0b56fSBarry Smith /*@C
36229e0b56fSBarry Smith    SNESNoLineSearchNoNorms - This routine is not a line search at
36329e0b56fSBarry Smith    all; it simply uses the full Newton step. This version does not
36429e0b56fSBarry Smith    even compute the norm of the function or search direction; this
36529e0b56fSBarry Smith    is intended only when you know the full step is fine and are
36629e0b56fSBarry Smith    not checking for convergence of the nonlinear iteration (for
367c7afd0dbSLois Curfman McInnes    example, you are running always for a fixed number of Newton steps).
368c7afd0dbSLois Curfman McInnes 
369c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
37029e0b56fSBarry Smith 
37129e0b56fSBarry Smith    Input Parameters:
372c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
373af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
37429e0b56fSBarry Smith .  x - current iterate
37529e0b56fSBarry Smith .  f - residual evaluated at x
37629e0b56fSBarry Smith .  y - search direction (contains new iterate on output)
37729e0b56fSBarry Smith .  w - work vector
378c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
37929e0b56fSBarry Smith 
38029e0b56fSBarry Smith    Output Parameters:
381c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
38229e0b56fSBarry Smith .  gnorm - not changed
38329e0b56fSBarry Smith .  ynorm - not changed
384a7cc72afSBarry Smith -  flag - set to PETSC_TRUE indicating a successful line search
385fee21e36SBarry Smith 
38629e0b56fSBarry Smith    Options Database Key:
3874b27c08aSLois Curfman McInnes .  -snes_ls basicnonorms - Activates SNESNoLineSearchNoNorms()
38829e0b56fSBarry Smith 
3898cbba510SBarry Smith    Notes:
390ea56f5baSLois Curfman McInnes    SNESNoLineSearchNoNorms() must be used in conjunction with
391ea56f5baSLois Curfman McInnes    either the options
392ea56f5baSLois Curfman McInnes $     -snes_no_convergence_test -snes_max_it <its>
393ea56f5baSLois Curfman McInnes    or alternatively a user-defined custom test set via
394329f5518SBarry Smith    SNESSetConvergenceTest(); or a -snes_max_it of 1,
395329f5518SBarry Smith    otherwise, the SNES solver will generate an error.
396329f5518SBarry Smith 
397329f5518SBarry Smith    During the final iteration this will not evaluate the function at
398329f5518SBarry Smith    the solution point. This is to save a function evaluation while
399329f5518SBarry Smith    using pseudo-timestepping.
4008cbba510SBarry Smith 
401ea56f5baSLois Curfman McInnes    The residual norms printed by monitoring routines such as
402ea56f5baSLois Curfman McInnes    SNESDefaultMonitor() (as activated via -snes_monitor) will not be
403ea56f5baSLois Curfman McInnes    correct, since they are not computed.
404ea56f5baSLois Curfman McInnes 
405ea56f5baSLois Curfman McInnes    Level: advanced
4068cbba510SBarry Smith 
40729e0b56fSBarry Smith .keywords: SNES, nonlinear, line search, cubic
40829e0b56fSBarry Smith 
40929e0b56fSBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(),
41029e0b56fSBarry Smith           SNESSetLineSearch(), SNESNoLineSearch()
41129e0b56fSBarry Smith @*/
412*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESNoLineSearchNoNorms(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
41329e0b56fSBarry Smith {
414dfbe8321SBarry Smith   PetscErrorCode ierr;
415ea709b57SSatish Balay   PetscScalar mone = -1.0;
4164b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
4178f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
41829e0b56fSBarry Smith 
4193a40ed3dSBarry Smith   PetscFunctionBegin;
420a7cc72afSBarry Smith   *flag = PETSC_TRUE;
421d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
42229e0b56fSBarry Smith   ierr = VecAYPX(&mone,x,y);CHKERRQ(ierr);            /* y <- y - x      */
4238f99978dSLois Curfman McInnes   if (neP->CheckStep) {
4248f99978dSLois Curfman McInnes    ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
4258f99978dSLois Curfman McInnes   }
426329f5518SBarry Smith 
427329f5518SBarry Smith   /* don't evaluate function the last time through */
428329f5518SBarry Smith   if (snes->iter < snes->max_its-1) {
42929e0b56fSBarry Smith     ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr); /* Compute F(y)    */
430329f5518SBarry Smith   }
431d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
4323a40ed3dSBarry Smith   PetscFunctionReturn(0);
43329e0b56fSBarry Smith }
43404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
4354a2ae208SSatish Balay #undef __FUNCT__
4364a2ae208SSatish Balay #define __FUNCT__ "SNESCubicLineSearch"
4374b828684SBarry Smith /*@C
438f525115eSLois Curfman McInnes    SNESCubicLineSearch - Performs a cubic line search (default line search method).
4395e76c431SBarry Smith 
440c7afd0dbSLois Curfman McInnes    Collective on SNES
441c7afd0dbSLois Curfman McInnes 
4425e76c431SBarry Smith    Input Parameters:
443c7afd0dbSLois Curfman McInnes +  snes - nonlinear context
444af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
4455e76c431SBarry Smith .  x - current iterate
4465e76c431SBarry Smith .  f - residual evaluated at x
4475e76c431SBarry Smith .  y - search direction (contains new iterate on output)
4485e76c431SBarry Smith .  w - work vector
449c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
4505e76c431SBarry Smith 
451393d2d9aSLois Curfman McInnes    Output Parameters:
452c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
4535e76c431SBarry Smith .  y - new iterate (contains search direction on input)
4545e76c431SBarry Smith .  gnorm - 2-norm of g
4555e76c431SBarry Smith .  ynorm - 2-norm of search length
456a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
457fee21e36SBarry Smith 
458c4a48953SLois Curfman McInnes    Options Database Key:
4594b27c08aSLois Curfman McInnes $  -snes_ls cubic - Activates SNESCubicLineSearch()
460c4a48953SLois Curfman McInnes 
4615e76c431SBarry Smith    Notes:
4625e76c431SBarry Smith    This line search is taken from "Numerical Methods for Unconstrained
4635e76c431SBarry Smith    Optimization and Nonlinear Equations" by Dennis and Schnabel, page 325.
4645e76c431SBarry Smith 
46536851e7fSLois Curfman McInnes    Level: advanced
46636851e7fSLois Curfman McInnes 
46728ae5a21SLois Curfman McInnes .keywords: SNES, nonlinear, line search, cubic
46828ae5a21SLois Curfman McInnes 
469af2d14d2SLois Curfman McInnes .seealso: SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
4705e76c431SBarry Smith @*/
471*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCubicLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
4725e76c431SBarry Smith {
473406556e6SLois Curfman McInnes   /*
474406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
475406556e6SLois Curfman McInnes      minimization problem:
476406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
477406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm, and fnorm = || f ||_2.
478406556e6SLois Curfman McInnes    */
479406556e6SLois Curfman McInnes 
4805ca10a37SBarry Smith   PetscReal   steptol,initslope,lambdaprev,gnormprev,a,b,d,t1,t2,rellength;
481329f5518SBarry Smith   PetscReal   maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg;
482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
483ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
4846b5873e3SBarry Smith #endif
485dfbe8321SBarry Smith   PetscErrorCode ierr;
486a7cc72afSBarry Smith   PetscInt count;
4874b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
488ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
4898f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
4905e76c431SBarry Smith 
4913a40ed3dSBarry Smith   PetscFunctionBegin;
492d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
493a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
4945e76c431SBarry Smith   alpha   = neP->alpha;
4955e76c431SBarry Smith   maxstep = neP->maxstep;
4965e76c431SBarry Smith   steptol = neP->steptol;
4975e76c431SBarry Smith 
498cddf8d76SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
49963b9fa5eSBarry Smith   if (*ynorm == 0.0) {
50063ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
501a1c2b6eeSLois Curfman McInnes     *gnorm = fnorm;
502a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
503a1c2b6eeSLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
504a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
505ad922ac9SBarry Smith     goto theend1;
50694a9d846SBarry Smith   }
5075e76c431SBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
5085e42470aSBarry Smith     scale = maxstep/(*ynorm);
509aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
51063ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",PetscRealPart(scale),*ynorm));CHKERRQ(ierr);
5116b5873e3SBarry Smith #else
51263ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Scaling step by %g old ynorm %g\n",scale,*ynorm));CHKERRQ(ierr);
5136b5873e3SBarry Smith #endif
514393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
5155e76c431SBarry Smith     *ynorm = maxstep;
5165e76c431SBarry Smith   }
517ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
5185ca10a37SBarry Smith   minlambda = steptol/rellength;
519a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
520aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
521a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
522329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
5235e42470aSBarry Smith #else
524a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
5255e42470aSBarry Smith #endif
5265e76c431SBarry Smith   if (initslope > 0.0) initslope = -initslope;
5275e76c431SBarry Smith   if (initslope == 0.0) initslope = -1.0;
5285e76c431SBarry Smith 
529393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
5305334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
53143e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
53263ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
53343e71028SBarry Smith     ierr = VecCopy(w,y);CHKERRQ(ierr);
53443e71028SBarry Smith     *flag = PETSC_FALSE;
53543e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
53643e71028SBarry Smith     goto theend1;
53743e71028SBarry Smith   }
53878b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
539313b5538SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
54043e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5415d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
542393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
54363ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Using full step\n"));CHKERRQ(ierr);
544ad922ac9SBarry Smith     goto theend1;
5455e76c431SBarry Smith   }
5465e76c431SBarry Smith 
5475e76c431SBarry Smith   /* Fit points with quadratic */
548313b5538SBarry Smith   lambda = 1.0;
549a703fe33SLois Curfman McInnes   lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
5505e76c431SBarry Smith   lambdaprev = lambda;
55143e71028SBarry Smith   printf("tmp %g ddsd %g %g %g %g\n",lambdatemp,((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope),*gnorm,fnorm,initslope);
5525e76c431SBarry Smith   gnormprev = *gnorm;
553329f5518SBarry Smith   if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
554ddd12767SBarry Smith   if (lambdatemp <= .1*lambda) lambda = .1*lambda;
555ddd12767SBarry Smith   else                         lambda = lambdatemp;
556393d2d9aSLois Curfman McInnes   ierr   = VecCopy(x,w);CHKERRQ(ierr);
557ea4d3ed3SLois Curfman McInnes   lambdaneg = -lambda;
558aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
559ea4d3ed3SLois Curfman McInnes   clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
5605e42470aSBarry Smith #else
561ea4d3ed3SLois Curfman McInnes   ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
5625e42470aSBarry Smith #endif
56343e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
56463ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while attempting quadratic backtracking! %D \n"));CHKERRQ(ierr);
56543e71028SBarry Smith     ierr = VecCopy(w,y);CHKERRQ(ierr);
56643e71028SBarry Smith     *flag = PETSC_FALSE;
56743e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
56843e71028SBarry Smith     goto theend1;
56943e71028SBarry Smith   }
57078b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
571cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
57243e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
5735ed2d596SBarry Smith   if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
574393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
57563ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Quadratically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
576ad922ac9SBarry Smith     goto theend1;
5775e76c431SBarry Smith   }
5785e76c431SBarry Smith 
5795e76c431SBarry Smith   /* Fit points with cubic */
5805e76c431SBarry Smith   count = 1;
5818229bfc2SMatthew Knepley   while (count < 10000) {
5825e76c431SBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
58363ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
58463ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
585f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
58643e71028SBarry Smith       *flag = PETSC_FALSE;
58743e71028SBarry Smith       break;
5885e76c431SBarry Smith     }
58943e71028SBarry Smith     printf("lamd %g %g\n",lambda,lambdaprev);
5905d1a10b1SBarry Smith     t1 = .5*((*gnorm)*(*gnorm) - fnorm*fnorm) - lambda*initslope;
5915d1a10b1SBarry Smith     t2 = .5*(gnormprev*gnormprev  - fnorm*fnorm) - lambdaprev*initslope;
592ddd12767SBarry Smith     a  = (t1/(lambda*lambda) - t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5932b022350SLois Curfman McInnes     b  = (-lambdaprev*t1/(lambda*lambda) + lambda*t2/(lambdaprev*lambdaprev))/(lambda-lambdaprev);
5945e76c431SBarry Smith     d  = b*b - 3*a*initslope;
5955e76c431SBarry Smith     if (d < 0.0) d = 0.0;
5965e76c431SBarry Smith     if (a == 0.0) {
5975e76c431SBarry Smith       lambdatemp = -initslope/(2.0*b);
5985e76c431SBarry Smith     } else {
5995e76c431SBarry Smith       lambdatemp = (-b + sqrt(d))/(3.0*a);
6005e76c431SBarry Smith     }
6015e76c431SBarry Smith     lambdaprev = lambda;
6025e76c431SBarry Smith     gnormprev  = *gnorm;
603329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
604329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
6055e76c431SBarry Smith     else                         lambda     = lambdatemp;
606393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
607ea4d3ed3SLois Curfman McInnes     lambdaneg = -lambda;
608aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
609ea4d3ed3SLois Curfman McInnes     clambda = lambdaneg;
610393d2d9aSLois Curfman McInnes     ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
6115e42470aSBarry Smith #else
612ea4d3ed3SLois Curfman McInnes     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
6135e42470aSBarry Smith #endif
61443e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
61563ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
61663ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
61743e71028SBarry Smith       ierr = VecCopy(w,y);CHKERRQ(ierr);
618ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
61943e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
620ed98166eSMatthew Knepley       break;
621ed98166eSMatthew Knepley     }
62243e71028SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
623cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
62443e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
6255ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* is reduction enough? */
626393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
62763ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubically determined step, lambda=%18.16e\n",lambda));CHKERRQ(ierr);
62843e71028SBarry Smith       break;
6292b022350SLois Curfman McInnes     } else {
63063ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESCubicLineSearch: Cubic step no good, shrinking lambda,  lambda=%18.16e\n",lambda));CHKERRQ(ierr);
6315e76c431SBarry Smith     }
6325e76c431SBarry Smith     count++;
6335e76c431SBarry Smith   }
6348229bfc2SMatthew Knepley   if (count >= 10000) {
635cd7625f8SBarry Smith     SETERRQ(PETSC_ERR_LIB, "Lambda was decreased more than 10,000 times, so something is probably wrong with the function evaluation");
6368229bfc2SMatthew Knepley   }
637ad922ac9SBarry Smith   theend1:
6388f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
63943e71028SBarry Smith   if (neP->CheckStep && *flag) {
6408f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
641abc0a331SBarry Smith     if (change_y) { /* recompute the function if the step has changed */
6428f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
6438f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
64443e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
64543e71028SBarry Smith       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
6468f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
6478f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
6488f99978dSLois Curfman McInnes     }
6498f99978dSLois Curfman McInnes   }
650d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
6513a40ed3dSBarry Smith   PetscFunctionReturn(0);
6525e76c431SBarry Smith }
65304d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
6544a2ae208SSatish Balay #undef __FUNCT__
6554a2ae208SSatish Balay #define __FUNCT__ "SNESQuadraticLineSearch"
6564b828684SBarry Smith /*@C
657f525115eSLois Curfman McInnes    SNESQuadraticLineSearch - Performs a quadratic line search.
6585e76c431SBarry Smith 
659c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
660c7afd0dbSLois Curfman McInnes 
6615e42470aSBarry Smith    Input Parameters:
662c7afd0dbSLois Curfman McInnes +  snes - the SNES context
663af2d14d2SLois Curfman McInnes .  lsctx - optional context for line search (not used here)
6645e42470aSBarry Smith .  x - current iterate
6655e42470aSBarry Smith .  f - residual evaluated at x
6665e42470aSBarry Smith .  y - search direction (contains new iterate on output)
6675e42470aSBarry Smith .  w - work vector
668c7afd0dbSLois Curfman McInnes -  fnorm - 2-norm of f
6695e42470aSBarry Smith 
670c4a48953SLois Curfman McInnes    Output Parameters:
671c7afd0dbSLois Curfman McInnes +  g - residual evaluated at new iterate y
6725e42470aSBarry Smith .  y - new iterate (contains search direction on input)
6735e42470aSBarry Smith .  gnorm - 2-norm of g
6745e42470aSBarry Smith .  ynorm - 2-norm of search length
675a7cc72afSBarry Smith -  flag - PETSC_TRUE if line search succeeds; PETSC_FALSE on failure.
676fee21e36SBarry Smith 
677c4a48953SLois Curfman McInnes    Options Database Key:
6784b27c08aSLois Curfman McInnes .  -snes_ls quadratic - Activates SNESQuadraticLineSearch()
679c4a48953SLois Curfman McInnes 
6805e42470aSBarry Smith    Notes:
6814b27c08aSLois Curfman McInnes    Use SNESSetLineSearch() to set this routine within the SNESLS method.
6825e42470aSBarry Smith 
68336851e7fSLois Curfman McInnes    Level: advanced
68436851e7fSLois Curfman McInnes 
685f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, quadratic, line search
686f59ffdedSLois Curfman McInnes 
687af2d14d2SLois Curfman McInnes .seealso: SNESCubicLineSearch(), SNESNoLineSearch(), SNESSetLineSearch(), SNESNoLineSearchNoNorms()
6885e42470aSBarry Smith @*/
689*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESQuadraticLineSearch(SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
6905e76c431SBarry Smith {
691406556e6SLois Curfman McInnes   /*
692406556e6SLois Curfman McInnes      Note that for line search purposes we work with with the related
693406556e6SLois Curfman McInnes      minimization problem:
694406556e6SLois Curfman McInnes         min  z(x):  R^n -> R,
695406556e6SLois Curfman McInnes      where z(x) = .5 * fnorm*fnorm,and fnorm = || f ||_2.
696406556e6SLois Curfman McInnes    */
6975ca10a37SBarry Smith   PetscReal   steptol,initslope,maxstep,minlambda,alpha,lambda,lambdatemp,lambdaneg,rellength;
698aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
699ea709b57SSatish Balay   PetscScalar cinitslope,clambda;
7006b5873e3SBarry Smith #endif
701dfbe8321SBarry Smith   PetscErrorCode ierr;
702a7cc72afSBarry Smith   PetscInt count;
7034b27c08aSLois Curfman McInnes   SNES_LS     *neP = (SNES_LS*)snes->data;
704ea709b57SSatish Balay   PetscScalar mone = -1.0,scale;
7058f99978dSLois Curfman McInnes   PetscTruth  change_y = PETSC_FALSE;
7065e76c431SBarry Smith 
7073a40ed3dSBarry Smith   PetscFunctionBegin;
708d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
709a7cc72afSBarry Smith   *flag   = PETSC_TRUE;
7105e42470aSBarry Smith   alpha   = neP->alpha;
7115e42470aSBarry Smith   maxstep = neP->maxstep;
7125e42470aSBarry Smith   steptol = neP->steptol;
7135e76c431SBarry Smith 
7143505fcc1SBarry Smith   ierr = VecNorm(y,NORM_2,ynorm);CHKERRQ(ierr);
71563b9fa5eSBarry Smith   if (*ynorm == 0.0) {
71663ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Search direction and size is 0\n"));CHKERRQ(ierr);
717b37302c6SLois Curfman McInnes     *gnorm = fnorm;
718b37302c6SLois Curfman McInnes     ierr   = VecCopy(x,y);CHKERRQ(ierr);
719b37302c6SLois Curfman McInnes     ierr   = VecCopy(f,g);CHKERRQ(ierr);
720a7cc72afSBarry Smith     *flag  = PETSC_FALSE;
721ad922ac9SBarry Smith     goto theend2;
72294a9d846SBarry Smith   }
7235e42470aSBarry Smith   if (*ynorm > maxstep) {	/* Step too big, so scale back */
7245e42470aSBarry Smith     scale = maxstep/(*ynorm);
725393d2d9aSLois Curfman McInnes     ierr = VecScale(&scale,y);CHKERRQ(ierr);
7265e42470aSBarry Smith     *ynorm = maxstep;
7275e76c431SBarry Smith   }
728ba6a83e5SMatthew Knepley   ierr      = VecMaxPointwiseDivide(y,x,&rellength);CHKERRQ(ierr);
7295ca10a37SBarry Smith   minlambda = steptol/rellength;
730a703fe33SLois Curfman McInnes   ierr = MatMult(snes->jacobian,y,w);CHKERRQ(ierr);
731aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
732a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&cinitslope);CHKERRQ(ierr);
733329f5518SBarry Smith   initslope = PetscRealPart(cinitslope);
7345e42470aSBarry Smith #else
735a703fe33SLois Curfman McInnes   ierr = VecDot(f,w,&initslope);CHKERRQ(ierr);
7365e42470aSBarry Smith #endif
7375e42470aSBarry Smith   if (initslope > 0.0) initslope = -initslope;
7385e42470aSBarry Smith   if (initslope == 0.0) initslope = -1.0;
7395e42470aSBarry Smith 
740393d2d9aSLois Curfman McInnes   ierr = VecCopy(y,w);CHKERRQ(ierr);
7415334005bSBarry Smith   ierr = VecAYPX(&mone,x,w);CHKERRQ(ierr);
74243e71028SBarry Smith   if (snes->nfuncs >= snes->max_funcs) {
74363ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while checking full step length!\n"));CHKERRQ(ierr);
74443e71028SBarry Smith     ierr = VecCopy(w,y);CHKERRQ(ierr);
74543e71028SBarry Smith     *flag = PETSC_FALSE;
74643e71028SBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
74743e71028SBarry Smith     goto theend2;
74843e71028SBarry Smith   }
74978b31e54SBarry Smith   ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
750cddf8d76SBarry Smith   ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
75143e71028SBarry Smith   if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7525d1a10b1SBarry Smith   if (.5*(*gnorm)*(*gnorm) <= .5*fnorm*fnorm + alpha*initslope) { /* Sufficient reduction */
753393d2d9aSLois Curfman McInnes     ierr = VecCopy(w,y);CHKERRQ(ierr);
75463ba0a88SBarry Smith     ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch: Using full step\n"));CHKERRQ(ierr);
755ad922ac9SBarry Smith     goto theend2;
7565e42470aSBarry Smith   }
7575e42470aSBarry Smith 
7585e42470aSBarry Smith   /* Fit points with quadratic */
759313b5538SBarry Smith   lambda = 1.0;
7605e42470aSBarry Smith   count = 1;
7615ca10a37SBarry Smith   while (PETSC_TRUE) {
7625e42470aSBarry Smith     if (lambda <= minlambda) { /* bad luck; use full step */
76363ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Unable to find good step length! %D \n",count));CHKERRQ(ierr);
76463ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%g, gnorm=%g, ynorm=%g, lambda=%g, initial slope=%g\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
765f168f371SBarry Smith       ierr = VecCopy(x,y);CHKERRQ(ierr);
76643e71028SBarry Smith       *flag = PETSC_FALSE;
76743e71028SBarry Smith       break;
7685e42470aSBarry Smith     }
769a703fe33SLois Curfman McInnes     lambdatemp = -initslope/((*gnorm)*(*gnorm) - fnorm*fnorm - 2.0*initslope);
770329f5518SBarry Smith     if (lambdatemp > .5*lambda)  lambdatemp = .5*lambda;
771329f5518SBarry Smith     if (lambdatemp <= .1*lambda) lambda     = .1*lambda;
772329f5518SBarry Smith     else                         lambda     = lambdatemp;
773393d2d9aSLois Curfman McInnes     ierr = VecCopy(x,w);CHKERRQ(ierr);
7743505fcc1SBarry Smith     lambdaneg = -lambda;
775aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX)
7763505fcc1SBarry Smith     clambda = lambdaneg; ierr = VecAXPY(&clambda,y,w);CHKERRQ(ierr);
7775e42470aSBarry Smith #else
7783505fcc1SBarry Smith     ierr = VecAXPY(&lambdaneg,y,w);CHKERRQ(ierr);
7795e42470aSBarry Smith #endif
78043e71028SBarry Smith     if (snes->nfuncs >= snes->max_funcs) {
78163ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Exceeded maximum function evaluations, while looking for good step length! %D \n",count));CHKERRQ(ierr);
78263ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lambda=%18.16e, initial slope=%18.16e\n",fnorm,*gnorm,*ynorm,lambda,initslope));CHKERRQ(ierr);
78343e71028SBarry Smith       ierr = VecCopy(w,y);CHKERRQ(ierr);
784ed98166eSMatthew Knepley       *flag = PETSC_FALSE;
78543e71028SBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
786ed98166eSMatthew Knepley       break;
787ed98166eSMatthew Knepley     }
78843e71028SBarry Smith     ierr = SNESComputeFunction(snes,w,g);CHKERRQ(ierr);
789cddf8d76SBarry Smith     ierr = VecNorm(g,NORM_2,gnorm);CHKERRQ(ierr);
79043e71028SBarry Smith     if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
7915ed2d596SBarry Smith     if (.5*(*gnorm)*(*gnorm) < .5*fnorm*fnorm + lambda*alpha*initslope) { /* sufficient reduction */
792393d2d9aSLois Curfman McInnes       ierr = VecCopy(w,y);CHKERRQ(ierr);
79363ba0a88SBarry Smith       ierr = PetscLogInfo((snes,"SNESQuadraticLineSearch:Quadratically determined step, lambda=%g\n",lambda));CHKERRQ(ierr);
79406259719SBarry Smith       break;
7955e42470aSBarry Smith     }
7965e42470aSBarry Smith     count++;
7975e42470aSBarry Smith   }
798ad922ac9SBarry Smith   theend2:
7998f99978dSLois Curfman McInnes   /* Optional user-defined check for line search step validity */
8008f99978dSLois Curfman McInnes   if (neP->CheckStep) {
8018f99978dSLois Curfman McInnes     ierr = (*neP->CheckStep)(snes,neP->checkP,y,&change_y);CHKERRQ(ierr);
802abc0a331SBarry Smith     if (change_y) { /* recompute the function if the step has changed */
8038f99978dSLois Curfman McInnes       ierr = SNESComputeFunction(snes,y,g);CHKERRQ(ierr);
8048f99978dSLois Curfman McInnes       ierr = VecNormBegin(g,NORM_2,gnorm);CHKERRQ(ierr);
80543e71028SBarry Smith       if (*gnorm != *gnorm) SETERRQ(PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
80643e71028SBarry Smith       ierr = VecNormBegin(y,NORM_2,ynorm);CHKERRQ(ierr);
8078f99978dSLois Curfman McInnes       ierr = VecNormEnd(y,NORM_2,ynorm);CHKERRQ(ierr);
8088f99978dSLois Curfman McInnes       ierr = VecNormEnd(g,NORM_2,gnorm);CHKERRQ(ierr);
8098f99978dSLois Curfman McInnes     }
8108f99978dSLois Curfman McInnes   }
811d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_LineSearch,snes,x,f,g);CHKERRQ(ierr);
8123a40ed3dSBarry Smith   PetscFunctionReturn(0);
8135e76c431SBarry Smith }
8142343ba6eSBarry Smith 
81504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8164a2ae208SSatish Balay #undef __FUNCT__
8174a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch"
818c9e6a524SLois Curfman McInnes /*@C
81977c4ece6SBarry Smith    SNESSetLineSearch - Sets the line search routine to be used
8204b27c08aSLois Curfman McInnes    by the method SNESLS.
8215e76c431SBarry Smith 
8225e76c431SBarry Smith    Input Parameters:
823c7afd0dbSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
824af2d14d2SLois Curfman McInnes .  lsctx - optional user-defined context for use by line search
825c7afd0dbSLois Curfman McInnes -  func - pointer to int function
8265e76c431SBarry Smith 
827fee21e36SBarry Smith    Collective on SNES
828fee21e36SBarry Smith 
829c4a48953SLois Curfman McInnes    Available Routines:
830c7afd0dbSLois Curfman McInnes +  SNESCubicLineSearch() - default line search
831f59ffdedSLois Curfman McInnes .  SNESQuadraticLineSearch() - quadratic line search
832f59ffdedSLois Curfman McInnes .  SNESNoLineSearch() - the full Newton step (actually not a line search)
833af2d14d2SLois Curfman McInnes -  SNESNoLineSearchNoNorms() - the full Newton step (calculating no norms; faster in parallel)
8345e76c431SBarry Smith 
835c4a48953SLois Curfman McInnes     Options Database Keys:
8364b27c08aSLois Curfman McInnes +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
8374b27c08aSLois Curfman McInnes .   -snes_ls_alpha <alpha> - Sets alpha
8384b27c08aSLois Curfman McInnes .   -snes_ls_maxstep <max> - Sets maxstep
8394b27c08aSLois Curfman McInnes -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
8403304466cSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
8413304466cSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
842c4a48953SLois Curfman McInnes 
8435e76c431SBarry Smith    Calling sequence of func:
844bd208895SLois Curfman McInnes .vb
845af2d14d2SLois Curfman McInnes    func (SNES snes,void *lsctx,Vec x,Vec f,Vec g,Vec y,Vec w,
846a7cc72afSBarry Smith          PetscReal fnorm,PetscReal *ynorm,PetscReal *gnorm,PetscTruth *flag)
847bd208895SLois Curfman McInnes .ve
8485e76c431SBarry Smith 
8495e76c431SBarry Smith     Input parameters for func:
850c7afd0dbSLois Curfman McInnes +   snes - nonlinear context
851af2d14d2SLois Curfman McInnes .   lsctx - optional user-defined context for line search
8525e76c431SBarry Smith .   x - current iterate
8535e76c431SBarry Smith .   f - residual evaluated at x
8545e76c431SBarry Smith .   y - search direction (contains new iterate on output)
8555e76c431SBarry Smith .   w - work vector
856c7afd0dbSLois Curfman McInnes -   fnorm - 2-norm of f
8575e76c431SBarry Smith 
8585e76c431SBarry Smith     Output parameters for func:
859c7afd0dbSLois Curfman McInnes +   g - residual evaluated at new iterate y
8605e76c431SBarry Smith .   y - new iterate (contains search direction on input)
8615e76c431SBarry Smith .   gnorm - 2-norm of g
8625e76c431SBarry Smith .   ynorm - 2-norm of search length
863a7cc72afSBarry Smith -   flag - set to PETSC_TRUE if the line search succeeds; PETSC_FALSE on failure.
864f59ffdedSLois Curfman McInnes 
86536851e7fSLois Curfman McInnes     Level: advanced
86636851e7fSLois Curfman McInnes 
867f59ffdedSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search, routine
868f59ffdedSLois Curfman McInnes 
8695d1a10b1SBarry Smith .seealso: SNESCubicLineSearch(), SNESQuadraticLineSearch(), SNESNoLineSearch(), SNESNoLineSearchNoNorms(),
8705d1a10b1SBarry Smith           SNESSetLineSearchCheck(), SNESSetLineSearchParams(), SNESGetLineSearchParams()
8715e76c431SBarry Smith @*/
872*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void *lsctx)
8735e76c431SBarry Smith {
874a7cc72afSBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*),void*);
87582bf6240SBarry Smith 
8763a40ed3dSBarry Smith   PetscFunctionBegin;
877c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearch_C",(void (**)(void))&f);CHKERRQ(ierr);
87882bf6240SBarry Smith   if (f) {
879af2d14d2SLois Curfman McInnes     ierr = (*f)(snes,func,lsctx);CHKERRQ(ierr);
88082bf6240SBarry Smith   }
8813a40ed3dSBarry Smith   PetscFunctionReturn(0);
8825e76c431SBarry Smith }
8838e019c35SBarry Smith 
884a7cc72afSBarry Smith typedef PetscErrorCode (*FCN2)(SNES,void*,Vec,Vec,Vec,Vec,Vec,PetscReal,PetscReal*,PetscReal*,PetscTruth*); /* force argument to next function to not be extern C*/
88504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
886fb2e594dSBarry Smith EXTERN_C_BEGIN
8874a2ae208SSatish Balay #undef __FUNCT__
8884a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearch_LS"
889*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearch_LS(SNES snes,FCN2 func,void *lsctx)
89082bf6240SBarry Smith {
89182bf6240SBarry Smith   PetscFunctionBegin;
8924b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->LineSearch = func;
8934b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->lsP        = lsctx;
89482bf6240SBarry Smith   PetscFunctionReturn(0);
89582bf6240SBarry Smith }
896fb2e594dSBarry Smith EXTERN_C_END
89704d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
8984a2ae208SSatish Balay #undef __FUNCT__
8994a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck"
900c8dd0c0dSLois Curfman McInnes /*@C
901530e4296SLois Curfman McInnes    SNESSetLineSearchCheck - Sets a routine to check the validity of new iterate computed
9024b27c08aSLois Curfman McInnes    by the line search routine in the Newton-based method SNESLS.
903c8dd0c0dSLois Curfman McInnes 
904c8dd0c0dSLois Curfman McInnes    Input Parameters:
905c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context obtained from SNESCreate()
906c8dd0c0dSLois Curfman McInnes .  func - pointer to int function
907c8dd0c0dSLois Curfman McInnes -  checkctx - optional user-defined context for use by step checking routine
908c8dd0c0dSLois Curfman McInnes 
909c8dd0c0dSLois Curfman McInnes    Collective on SNES
910c8dd0c0dSLois Curfman McInnes 
911c8dd0c0dSLois Curfman McInnes    Calling sequence of func:
912c8dd0c0dSLois Curfman McInnes .vb
913b1ae27eaSLois Curfman McInnes    int func (SNES snes, void *checkctx, Vec x, PetscTruth *flag)
914c8dd0c0dSLois Curfman McInnes .ve
915b1ae27eaSLois Curfman McInnes    where func returns an error code of 0 on success and a nonzero
916b1ae27eaSLois Curfman McInnes    on failure.
917c8dd0c0dSLois Curfman McInnes 
918c8dd0c0dSLois Curfman McInnes    Input parameters for func:
919c8dd0c0dSLois Curfman McInnes +  snes - nonlinear context
920c8dd0c0dSLois Curfman McInnes .  checkctx - optional user-defined context for use by step checking routine
9218f99978dSLois Curfman McInnes -  x - current candidate iterate
922c8dd0c0dSLois Curfman McInnes 
923c8dd0c0dSLois Curfman McInnes    Output parameters for func:
924c8dd0c0dSLois Curfman McInnes +  x - current iterate (possibly modified)
925c8dd0c0dSLois Curfman McInnes -  flag - flag indicating whether x has been modified (either
926c8dd0c0dSLois Curfman McInnes            PETSC_TRUE of PETSC_FALSE)
927c8dd0c0dSLois Curfman McInnes 
928c8dd0c0dSLois Curfman McInnes    Level: advanced
929c8dd0c0dSLois Curfman McInnes 
9308f99978dSLois Curfman McInnes    Notes:
931b1ae27eaSLois Curfman McInnes    SNESNoLineSearch() and SNESNoLineSearchNoNorms() accept the new
932b1ae27eaSLois Curfman McInnes    iterate computed by the line search checking routine.  In particular,
933b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1}, (2) pass control
934b1ae27eaSLois Curfman McInnes    to the checking routine, and then (3) compute the corresponding nonlinear
935ea56f5baSLois Curfman McInnes    function f(u_{i+1}) with the (possibly altered) iterate u_{i+1}.
9368f99978dSLois Curfman McInnes 
937b1ae27eaSLois Curfman McInnes    SNESQuadraticLineSearch() and SNESCubicLineSearch() also accept the
938b1ae27eaSLois Curfman McInnes    new iterate computed by the line search checking routine.  In particular,
939b1ae27eaSLois Curfman McInnes    these routines (1) compute a candidate iterate u_{i+1} as well as a
940ea56f5baSLois Curfman McInnes    candidate nonlinear function f(u_{i+1}), (2) pass control to the checking
941ea56f5baSLois Curfman McInnes    routine, and then (3) force a re-evaluation of f(u_{i+1}) if any changes
942ea56f5baSLois Curfman McInnes    were made to the candidate iterate in the checking routine (as indicated
943ea56f5baSLois Curfman McInnes    by flag=PETSC_TRUE).  The overhead of this function re-evaluation can be
944b1ae27eaSLois Curfman McInnes    very costly, so use this feature with caution!
9458f99978dSLois Curfman McInnes 
946c8dd0c0dSLois Curfman McInnes .keywords: SNES, nonlinear, set, line search check, step check, routine
947c8dd0c0dSLois Curfman McInnes 
948c8dd0c0dSLois Curfman McInnes .seealso: SNESSetLineSearch()
949c8dd0c0dSLois Curfman McInnes @*/
950*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearchCheck(SNES snes,PetscErrorCode (*func)(SNES,void*,Vec,PetscTruth*),void *checkctx)
951c8dd0c0dSLois Curfman McInnes {
9526849ba73SBarry Smith   PetscErrorCode ierr,(*f)(SNES,PetscErrorCode (*)(SNES,void*,Vec,PetscTruth*),void*);
953c8dd0c0dSLois Curfman McInnes 
954c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
955c134de8dSSatish Balay   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESSetLineSearchCheck_C",(void (**)(void))&f);CHKERRQ(ierr);
956c8dd0c0dSLois Curfman McInnes   if (f) {
957c8dd0c0dSLois Curfman McInnes     ierr = (*f)(snes,func,checkctx);CHKERRQ(ierr);
958c8dd0c0dSLois Curfman McInnes   }
959c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
960c8dd0c0dSLois Curfman McInnes }
961c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
9626849ba73SBarry Smith typedef PetscErrorCode (*FCN)(SNES,void*,Vec,PetscTruth*); /* force argument to next function to not be extern C*/
963c8dd0c0dSLois Curfman McInnes EXTERN_C_BEGIN
9644a2ae208SSatish Balay #undef __FUNCT__
9654a2ae208SSatish Balay #define __FUNCT__ "SNESSetLineSearchCheck_LS"
966*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLineSearchCheck_LS(SNES snes,FCN func,void *checkctx)
967c8dd0c0dSLois Curfman McInnes {
968c8dd0c0dSLois Curfman McInnes   PetscFunctionBegin;
9694b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->CheckStep = func;
9704b27c08aSLois Curfman McInnes   ((SNES_LS *)(snes->data))->checkP    = checkctx;
971c8dd0c0dSLois Curfman McInnes   PetscFunctionReturn(0);
972c8dd0c0dSLois Curfman McInnes }
973c8dd0c0dSLois Curfman McInnes EXTERN_C_END
974c8dd0c0dSLois Curfman McInnes /* -------------------------------------------------------------------------- */
97504d965bbSLois Curfman McInnes /*
9764b27c08aSLois Curfman McInnes    SNESPrintHelp_LS - Prints all options for the SNES_LS method.
977329e7e40SMatthew Knepley 
978329e7e40SMatthew Knepley    Input Parameter:
979329e7e40SMatthew Knepley .  snes - the SNES context
980329e7e40SMatthew Knepley 
981329e7e40SMatthew Knepley    Application Interface Routine: SNESPrintHelp()
982329e7e40SMatthew Knepley */
983329e7e40SMatthew Knepley #undef __FUNCT__
9844b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESPrintHelp_LS"
9856849ba73SBarry Smith static PetscErrorCode SNESPrintHelp_LS(SNES snes,char *p)
986329e7e40SMatthew Knepley {
9874b27c08aSLois Curfman McInnes   SNES_LS *ls = (SNES_LS *)snes->data;
988329e7e40SMatthew Knepley 
989329e7e40SMatthew Knepley   PetscFunctionBegin;
9904b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm," method SNES_LS (ls) for systems of nonlinear equations:\n");
9914b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls [cubic,quadratic,basic,basicnonorms,...]\n",p);
9924b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_alpha <alpha> (default %g)\n",p,ls->alpha);
9934b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_maxstep <max> (default %g)\n",p,ls->maxstep);
9944b27c08aSLois Curfman McInnes   (*PetscHelpPrintf)(snes->comm,"   %ssnes_ls_steptol <tol> (default %g)\n",p,ls->steptol);
995329e7e40SMatthew Knepley   PetscFunctionReturn(0);
996329e7e40SMatthew Knepley }
997329e7e40SMatthew Knepley 
998329e7e40SMatthew Knepley /*
9994b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
100004d965bbSLois Curfman McInnes 
100104d965bbSLois Curfman McInnes    Input Parameters:
100204d965bbSLois Curfman McInnes .  SNES - the SNES context
100304d965bbSLois Curfman McInnes .  viewer - visualization context
100404d965bbSLois Curfman McInnes 
100504d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
100604d965bbSLois Curfman McInnes */
10074a2ae208SSatish Balay #undef __FUNCT__
10084b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
10096849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
1010a935fc98SLois Curfman McInnes {
10114b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
10122fc52814SBarry Smith   const char *cstr;
1013dfbe8321SBarry Smith   PetscErrorCode ierr;
101432077d6dSBarry Smith   PetscTruth iascii;
1015a935fc98SLois Curfman McInnes 
10163a40ed3dSBarry Smith   PetscFunctionBegin;
101732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
101832077d6dSBarry Smith   if (iascii) {
101919bcc07fSBarry Smith     if (ls->LineSearch == SNESNoLineSearch)             cstr = "SNESNoLineSearch";
102019bcc07fSBarry Smith     else if (ls->LineSearch == SNESQuadraticLineSearch) cstr = "SNESQuadraticLineSearch";
102119bcc07fSBarry Smith     else if (ls->LineSearch == SNESCubicLineSearch)     cstr = "SNESCubicLineSearch";
102219bcc07fSBarry Smith     else                                                cstr = "unknown";
1023b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
1024b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%g, maxstep=%g, steptol=%g\n",ls->alpha,ls->maxstep,ls->steptol);CHKERRQ(ierr);
10255cd90555SBarry Smith   } else {
102679a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES EQ LS",((PetscObject)viewer)->type_name);
102719bcc07fSBarry Smith   }
10283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1029a935fc98SLois Curfman McInnes }
103004d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
103104d965bbSLois Curfman McInnes /*
10324b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
103304d965bbSLois Curfman McInnes 
103404d965bbSLois Curfman McInnes    Input Parameter:
103504d965bbSLois Curfman McInnes .  snes - the SNES context
103604d965bbSLois Curfman McInnes 
103704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
103804d965bbSLois Curfman McInnes */
10394a2ae208SSatish Balay #undef __FUNCT__
10404b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
10416849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
10425e42470aSBarry Smith {
10434b27c08aSLois Curfman McInnes   SNES_LS    *ls = (SNES_LS *)snes->data;
1044e5999256SBarry Smith   const char *lses[] = {"basic","basicnonorms","quadratic","cubic"};
1045dfbe8321SBarry Smith   PetscErrorCode ierr;
1046a7cc72afSBarry Smith   PetscInt indx;
1047f1af5d2fSBarry Smith   PetscTruth flg;
10485e42470aSBarry Smith 
10493a40ed3dSBarry Smith   PetscFunctionBegin;
1050b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
10514b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_alpha","Function norm must decrease by","None",ls->alpha,&ls->alpha,0);CHKERRQ(ierr);
10524b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_maxstep","Step must be less than","None",ls->maxstep,&ls->maxstep,0);CHKERRQ(ierr);
10534b27c08aSLois Curfman McInnes     ierr = PetscOptionsReal("-snes_ls_steptol","Step must be greater than","None",ls->steptol,&ls->steptol,0);CHKERRQ(ierr);
1054186905e3SBarry Smith 
105522e36657SBarry Smith     ierr = PetscOptionsEList("-snes_ls","Line search used","SNESSetLineSearch",lses,4,"cubic",&indx,&flg);CHKERRQ(ierr);
105625cdf11fSBarry Smith     if (flg) {
105722e36657SBarry Smith       switch (indx) {
1058b49fd9e1SBarry Smith       case 0:
1059af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearch,PETSC_NULL);CHKERRQ(ierr);
1060b49fd9e1SBarry Smith         break;
1061b49fd9e1SBarry Smith       case 1:
1062af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESNoLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
1063b49fd9e1SBarry Smith         break;
1064b49fd9e1SBarry Smith       case 2:
1065af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESQuadraticLineSearch,PETSC_NULL);CHKERRQ(ierr);
1066b49fd9e1SBarry Smith         break;
1067b49fd9e1SBarry Smith       case 3:
1068af2d14d2SLois Curfman McInnes         ierr = SNESSetLineSearch(snes,SNESCubicLineSearch,PETSC_NULL);CHKERRQ(ierr);
1069b49fd9e1SBarry Smith         break;
10705e42470aSBarry Smith       }
10715e42470aSBarry Smith     }
1072b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
10733a40ed3dSBarry Smith   PetscFunctionReturn(0);
10745e42470aSBarry Smith }
107504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
1076ebe3b25bSBarry Smith /*MC
1077ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
107804d965bbSLois Curfman McInnes 
1079ebe3b25bSBarry Smith    Options Database:
1080ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
1081ebe3b25bSBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha
1082ebe3b25bSBarry Smith .   -snes_ls_maxstep <max> - Sets maxstep
1083ebe3b25bSBarry Smith -   -snes_ls_steptol <steptol> - Sets steptol, this is the minimum step size that the line search code
1084ebe3b25bSBarry Smith                    will accept; min p[i]/x[i] < steptol. The -snes_stol <stol> is the minimum step length
1085ebe3b25bSBarry Smith                    the default convergence test will use and is based on 2-norm(p) < stol*2-norm(x)
108604d965bbSLois Curfman McInnes 
1087ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
108804d965bbSLois Curfman McInnes 
1089ee3001cbSBarry Smith    Level: beginner
1090ee3001cbSBarry Smith 
1091ebe3b25bSBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESSetLineSearch(),
1092ebe3b25bSBarry Smith            SNESSetLineSearchCheck(), SNESNoLineSearch(), SNESCubicLineSearch(), SNESQuadraticLineSearch(),
1093ebe3b25bSBarry Smith           SNESSetLineSearch(), SNESNoLineSearchNoNorms()
1094ebe3b25bSBarry Smith 
1095ebe3b25bSBarry Smith M*/
1096fb2e594dSBarry Smith EXTERN_C_BEGIN
10974a2ae208SSatish Balay #undef __FUNCT__
10984b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
1099*63dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate_LS(SNES snes)
11005e42470aSBarry Smith {
1101dfbe8321SBarry Smith   PetscErrorCode ierr;
11024b27c08aSLois Curfman McInnes   SNES_LS *neP;
11035e42470aSBarry Smith 
11043a40ed3dSBarry Smith   PetscFunctionBegin;
11054b27c08aSLois Curfman McInnes   snes->setup		= SNESSetUp_LS;
11064b27c08aSLois Curfman McInnes   snes->solve		= SNESSolve_LS;
11074b27c08aSLois Curfman McInnes   snes->destroy		= SNESDestroy_LS;
11084b27c08aSLois Curfman McInnes   snes->converged	= SNESConverged_LS;
11094b27c08aSLois Curfman McInnes   snes->printhelp       = SNESPrintHelp_LS;
11104b27c08aSLois Curfman McInnes   snes->setfromoptions  = SNESSetFromOptions_LS;
11114b27c08aSLois Curfman McInnes   snes->view            = SNESView_LS;
11125baf8537SBarry Smith   snes->nwork           = 0;
11135e42470aSBarry Smith 
11144b27c08aSLois Curfman McInnes   ierr                  = PetscNew(SNES_LS,&neP);CHKERRQ(ierr);
111552e6d16bSBarry Smith   ierr = PetscLogObjectMemory(snes,sizeof(SNES_LS));CHKERRQ(ierr);
11165e42470aSBarry Smith   snes->data    	= (void*)neP;
11175e42470aSBarry Smith   neP->alpha		= 1.e-4;
11185e42470aSBarry Smith   neP->maxstep		= 1.e8;
11195e42470aSBarry Smith   neP->steptol		= 1.e-12;
11205e42470aSBarry Smith   neP->LineSearch       = SNESCubicLineSearch;
1121c8dd0c0dSLois Curfman McInnes   neP->lsP              = PETSC_NULL;
1122c8dd0c0dSLois Curfman McInnes   neP->CheckStep        = PETSC_NULL;
1123c8dd0c0dSLois Curfman McInnes   neP->checkP           = PETSC_NULL;
112482bf6240SBarry Smith 
11255d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearch_C","SNESSetLineSearch_LS",SNESSetLineSearch_LS);CHKERRQ(ierr);
11265d1a10b1SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESSetLineSearchCheck_C","SNESSetLineSearchCheck_LS",SNESSetLineSearchCheck_LS);CHKERRQ(ierr);
112782bf6240SBarry Smith 
11283a40ed3dSBarry Smith   PetscFunctionReturn(0);
11295e42470aSBarry Smith }
1130fb2e594dSBarry Smith EXTERN_C_END
113182bf6240SBarry Smith 
113282bf6240SBarry Smith 
113382bf6240SBarry Smith 
1134