xref: /petsc/src/snes/impls/ls/ls.c (revision 15f5eeead57292b8b561bfc988cadbf25c9e1d6e)
15e76c431SBarry Smith 
2c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h>
35e42470aSBarry Smith 
48a5d9ee4SBarry Smith /*
520f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
620f52e01SBarry Smith     || F(u) ||_2 but not a zero, F(u) = 0. In the case when one cannot compute J^T F we use the fact that
736669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
820f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
98a5d9ee4SBarry Smith */
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckLocalMin_Private"
12ace3abfcSBarry Smith PetscErrorCode SNESLSCheckLocalMin_Private(SNES snes,Mat A,Vec F,Vec W,PetscReal fnorm,PetscBool  *ismin)
138a5d9ee4SBarry Smith {
148a5d9ee4SBarry Smith   PetscReal      a1;
15dfbe8321SBarry Smith   PetscErrorCode ierr;
16ace3abfcSBarry Smith   PetscBool      hastranspose;
178a5d9ee4SBarry Smith 
188a5d9ee4SBarry Smith   PetscFunctionBegin;
198a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
2036669109SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
2136669109SBarry Smith   if (hastranspose) {
228a5d9ee4SBarry Smith     /* Compute || J^T F|| */
238a5d9ee4SBarry Smith     ierr = MatMultTranspose(A,F,W);CHKERRQ(ierr);
248a5d9ee4SBarry Smith     ierr = VecNorm(W,NORM_2,&a1);CHKERRQ(ierr);
258f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"|| J^T F|| %14.12e near zero implies found a local minimum\n",(double)(a1/fnorm));CHKERRQ(ierr);
268a5d9ee4SBarry Smith     if (a1/fnorm < 1.e-4) *ismin = PETSC_TRUE;
2736669109SBarry Smith   } else {
2836669109SBarry Smith     Vec         work;
29ea709b57SSatish Balay     PetscScalar result;
3036669109SBarry Smith     PetscReal   wnorm;
3136669109SBarry Smith 
32abb0e124SMatthew Knepley     ierr = VecSetRandom(W,PETSC_NULL);CHKERRQ(ierr);
3336669109SBarry Smith     ierr = VecNorm(W,NORM_2,&wnorm);CHKERRQ(ierr);
3436669109SBarry Smith     ierr = VecDuplicate(W,&work);CHKERRQ(ierr);
3536669109SBarry Smith     ierr = MatMult(A,W,work);CHKERRQ(ierr);
3636669109SBarry Smith     ierr = VecDot(F,work,&result);CHKERRQ(ierr);
376bf464f9SBarry Smith     ierr = VecDestroy(&work);CHKERRQ(ierr);
3836669109SBarry Smith     a1   = PetscAbsScalar(result)/(fnorm*wnorm);
398f1a2a5eSBarry Smith     ierr = PetscInfo1(snes,"(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n",(double)a1);CHKERRQ(ierr);
4036669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4136669109SBarry Smith   }
428a5d9ee4SBarry Smith   PetscFunctionReturn(0);
438a5d9ee4SBarry Smith }
448a5d9ee4SBarry Smith 
4574637425SBarry Smith /*
465ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4774637425SBarry Smith */
484a2ae208SSatish Balay #undef __FUNCT__
494a2ae208SSatish Balay #define __FUNCT__ "SNESLSCheckResidual_Private"
501e2582c4SBarry Smith PetscErrorCode SNESLSCheckResidual_Private(SNES snes,Mat A,Vec F,Vec X,Vec W1,Vec W2)
5174637425SBarry Smith {
5274637425SBarry Smith   PetscReal      a1,a2;
53dfbe8321SBarry Smith   PetscErrorCode ierr;
54ace3abfcSBarry Smith   PetscBool      hastranspose;
5574637425SBarry Smith 
5674637425SBarry Smith   PetscFunctionBegin;
5774637425SBarry Smith   ierr = MatHasOperation(A,MATOP_MULT_TRANSPOSE,&hastranspose);CHKERRQ(ierr);
5874637425SBarry Smith   if (hastranspose) {
5974637425SBarry Smith     ierr = MatMult(A,X,W1);CHKERRQ(ierr);
6079f36460SBarry Smith     ierr = VecAXPY(W1,-1.0,F);CHKERRQ(ierr);
6174637425SBarry Smith 
6274637425SBarry Smith     /* Compute || J^T W|| */
6374637425SBarry Smith     ierr = MatMultTranspose(A,W1,W2);CHKERRQ(ierr);
6474637425SBarry Smith     ierr = VecNorm(W1,NORM_2,&a1);CHKERRQ(ierr);
6574637425SBarry Smith     ierr = VecNorm(W2,NORM_2,&a2);CHKERRQ(ierr);
6675567043SBarry Smith     if (a1 != 0.0) {
678f1a2a5eSBarry Smith       ierr = PetscInfo1(snes,"||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n",(double)(a2/a1));CHKERRQ(ierr);
6885471664SBarry Smith     }
6974637425SBarry Smith   }
7074637425SBarry Smith   PetscFunctionReturn(0);
7174637425SBarry Smith }
7274637425SBarry Smith 
7304d965bbSLois Curfman McInnes /*  --------------------------------------------------------------------
7404d965bbSLois Curfman McInnes 
7504d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7694b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7704d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
7804d965bbSLois Curfman McInnes      respectively.
7904d965bbSLois Curfman McInnes 
80fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8104d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8204d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
83fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8404d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8504d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
864b27c08aSLois Curfman McInnes      we use _LS (e.g., SNESCreate_LS, SNESSolve_LS) for solving
874b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
8804d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
8904d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9004d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9104d965bbSLois Curfman McInnes      for all nonlinear solvers.
9204d965bbSLois Curfman McInnes 
9304d965bbSLois Curfman McInnes      Another key routine is:
9404d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9504d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9604d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9704d965bbSLois Curfman McInnes      SNESSolve() if necessary.
9804d965bbSLois Curfman McInnes 
9904d965bbSLois Curfman McInnes      Additional basic routines are:
10004d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10104d965bbSLois Curfman McInnes                                       have actually been used.
10204d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
103186905e3SBarry Smith      SNESView().
10404d965bbSLois Curfman McInnes 
10504d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10604d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10704d965bbSLois Curfman McInnes      above description applies to these categories also.
10804d965bbSLois Curfman McInnes 
10904d965bbSLois Curfman McInnes     -------------------------------------------------------------------- */
1105e42470aSBarry Smith /*
1114b27c08aSLois Curfman McInnes    SNESSolve_LS - Solves a nonlinear system with a truncated Newton
11204d965bbSLois Curfman McInnes    method with a line search.
1135e76c431SBarry Smith 
11404d965bbSLois Curfman McInnes    Input Parameters:
11504d965bbSLois Curfman McInnes .  snes - the SNES context
1165e76c431SBarry Smith 
11704d965bbSLois Curfman McInnes    Output Parameter:
118c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
11904d965bbSLois Curfman McInnes 
12004d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1215e76c431SBarry Smith 
1225e76c431SBarry Smith    Notes:
1235e76c431SBarry Smith    This implements essentially a truncated Newton method with a
1245e76c431SBarry Smith    line search.  By default a cubic backtracking line search
1255e76c431SBarry Smith    is employed, as described in the text "Numerical Methods for
1265e76c431SBarry Smith    Unconstrained Optimization and Nonlinear Equations" by Dennis
127393d2d9aSLois Curfman McInnes    and Schnabel.
1285e76c431SBarry Smith */
1294a2ae208SSatish Balay #undef __FUNCT__
1304b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSolve_LS"
131dfbe8321SBarry Smith PetscErrorCode SNESSolve_LS(SNES snes)
1325e76c431SBarry Smith {
1336849ba73SBarry Smith   PetscErrorCode     ierr;
134a7cc72afSBarry Smith   PetscInt           maxits,i,lits;
135ace3abfcSBarry Smith   PetscBool          lssucceed;
136112a2221SBarry Smith   MatStructure       flg = DIFFERENT_NONZERO_PATTERN;
13785385478SLisandro Dalcin   PetscReal          fnorm,gnorm,xnorm=0,ynorm;
13885385478SLisandro Dalcin   Vec                Y,X,F,G,W;
1393d4c4710SBarry Smith   KSPConvergedReason kspreason;
1405e76c431SBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
1423d4c4710SBarry Smith   snes->numFailures            = 0;
1433d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
144184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
145184914b5SBarry Smith 
1465e42470aSBarry Smith   maxits	= snes->max_its;	/* maximum number of iterations */
1475e42470aSBarry Smith   X		= snes->vec_sol;	/* solution vector */
14839e2f89bSBarry Smith   F		= snes->vec_func;	/* residual vector */
1495e42470aSBarry Smith   Y		= snes->work[0];	/* work vectors */
1505e42470aSBarry Smith   G		= snes->work[1];
1515e42470aSBarry Smith   W		= snes->work[2];
1525e76c431SBarry Smith 
1534c49b128SBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1544c49b128SBarry Smith   snes->iter = 0;
15575567043SBarry Smith   snes->norm = 0.0;
1564c49b128SBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
15719717074SBarry Smith   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1584936397dSBarry Smith   if (snes->domainerror) {
1594936397dSBarry Smith     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1604936397dSBarry Smith     PetscFunctionReturn(0);
1614936397dSBarry Smith   }
1622613ca53SBarry Smith   ierr = VecNormBegin(F,NORM_2,&fnorm);CHKERRQ(ierr);	/* fnorm <- ||F||  */
1632613ca53SBarry Smith   ierr = VecNormBegin(X,NORM_2,&xnorm);CHKERRQ(ierr);	/* xnorm <- ||x||  */
1642613ca53SBarry Smith   ierr = VecNormEnd(F,NORM_2,&fnorm);CHKERRQ(ierr);
1652613ca53SBarry Smith   ierr = VecNormEnd(X,NORM_2,&xnorm);CHKERRQ(ierr);
16662cbcd01SMatthew G Knepley   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"User provided compute function generated a Not-a-Number");
1670f5bd95cSBarry Smith   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1685e42470aSBarry Smith   snes->norm = fnorm;
1690f5bd95cSBarry Smith   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
17084c569c9SBarry Smith   SNESLogConvHistory(snes,fnorm,0);
1717a03ce2fSLisandro Dalcin   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1723f149594SLisandro Dalcin 
1733f149594SLisandro Dalcin   /* set parameter for default relative tolerance convergence test */
1743f149594SLisandro Dalcin   snes->ttol = fnorm*snes->rtol;
17585385478SLisandro Dalcin   /* test convergence */
17685385478SLisandro Dalcin   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
17706ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
178d034289bSBarry Smith 
1795e76c431SBarry Smith   for (i=0; i<maxits; i++) {
1805e76c431SBarry Smith 
181329e7e40SMatthew Knepley     /* Call general purpose update function */
182e7788613SBarry Smith     if (snes->ops->update) {
183e7788613SBarry Smith       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
184de8cb200SMatthew Knepley     }
185329e7e40SMatthew Knepley 
186ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
1875334005bSBarry Smith     ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
18894b7f48cSBarry Smith     ierr = KSPSetOperators(snes->ksp,snes->jacobian,snes->jacobian_pre,flg);CHKERRQ(ierr);
18971f87433Sdalcinl     ierr = SNES_KSPSolve(snes,snes->ksp,F,Y);CHKERRQ(ierr);
1903d4c4710SBarry Smith     ierr = KSPGetConvergedReason(snes->ksp,&kspreason);CHKERRQ(ierr);
1913d4c4710SBarry Smith     if (kspreason < 0) {
1923d4c4710SBarry Smith       if (++snes->numLinearSolveFailures >= snes->maxLinearSolveFailures) {
193ef998cc9SBarry Smith         ierr = PetscInfo2(snes,"iter=%D, number linear solve failures %D greater than current SNES allowed, stopping solve\n",snes->iter,snes->numLinearSolveFailures);CHKERRQ(ierr);
1943d4c4710SBarry Smith         snes->reason = SNES_DIVERGED_LINEAR_SOLVE;
195ab81109fSSatish Balay         break;
1963d4c4710SBarry Smith       }
1973d4c4710SBarry Smith     }
1983d4c4710SBarry Smith     ierr = KSPGetIterationNumber(snes->ksp,&lits);CHKERRQ(ierr);
1993f149594SLisandro Dalcin     snes->linear_its += lits;
2003f149594SLisandro Dalcin     ierr = PetscInfo2(snes,"iter=%D, linear solve iterations=%D\n",snes->iter,lits);CHKERRQ(ierr);
20174637425SBarry Smith 
202ea630c6eSPeter Brune     if (snes->ops->precheckstep) {
203ace3abfcSBarry Smith       PetscBool  changed_y = PETSC_FALSE;
204ea630c6eSPeter Brune       ierr = (*snes->ops->precheckstep)(snes,X,Y,snes->precheck,&changed_y);CHKERRQ(ierr);
2059c3ca13aSBarry Smith     }
2069c3ca13aSBarry Smith 
207b0a32e0cSBarry Smith     if (PetscLogPrintInfo){
2081e2582c4SBarry Smith       ierr = SNESLSCheckResidual_Private(snes,snes->jacobian,F,Y,G,W);CHKERRQ(ierr);
20985471664SBarry Smith     }
21074637425SBarry Smith 
211ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
212ea4d3ed3SLois Curfman McInnes          Y <- X - lambda*Y
213e68848bdSBarry Smith        and evaluate G = function(Y) (depends on the line search).
214ea4d3ed3SLois Curfman McInnes     */
21585385478SLisandro Dalcin     ierr = VecCopy(Y,snes->vec_sol_update);CHKERRQ(ierr);
2163f149594SLisandro Dalcin     ynorm = 1; gnorm = fnorm;
217ea630c6eSPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
2188f1a2a5eSBarry Smith     ierr = PetscInfo4(snes,"fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n",(double)fnorm,(double)gnorm,(double)ynorm,(int)lssucceed);CHKERRQ(ierr);
2194a93e4ddSBarry Smith     if (snes->reason == SNES_DIVERGED_FUNCTION_COUNT) break;
2204936397dSBarry Smith     if (snes->domainerror) {
2214936397dSBarry Smith       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
2224936397dSBarry Smith       PetscFunctionReturn(0);
2234936397dSBarry Smith     }
224a7cc72afSBarry Smith     if (!lssucceed) {
22550ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
226ace3abfcSBarry Smith 	PetscBool  ismin;
227647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2281e2582c4SBarry Smith         ierr = SNESLSCheckLocalMin_Private(snes,snes->jacobian,G,W,gnorm,&ismin);CHKERRQ(ierr);
2298a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2303505fcc1SBarry Smith         break;
2313505fcc1SBarry Smith       }
23250ffb88aSMatthew Knepley     }
23385385478SLisandro Dalcin     /* Update function and solution vectors */
23485385478SLisandro Dalcin     fnorm = gnorm;
23585385478SLisandro Dalcin     ierr = VecCopy(G,F);CHKERRQ(ierr);
23685385478SLisandro Dalcin     ierr = VecCopy(W,X);CHKERRQ(ierr);
23785385478SLisandro Dalcin     /* Monitor convergence */
23885385478SLisandro Dalcin     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
23985385478SLisandro Dalcin     snes->iter = i+1;
24085385478SLisandro Dalcin     snes->norm = fnorm;
24185385478SLisandro Dalcin     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
24285385478SLisandro Dalcin     SNESLogConvHistory(snes,snes->norm,lits);
2437a03ce2fSLisandro Dalcin     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
24485385478SLisandro Dalcin     /* Test for convergence, xnorm = || X || */
24585385478SLisandro Dalcin     if (snes->ops->converged != SNESSkipConverged) { ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); }
246e7788613SBarry Smith     ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,ynorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
2473f149594SLisandro Dalcin     if (snes->reason) break;
24816c95cacSBarry Smith   }
24952392280SLois Curfman McInnes   if (i == maxits) {
250ae15b995SBarry Smith     ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",maxits);CHKERRQ(ierr);
25185385478SLisandro Dalcin     if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
25252392280SLois Curfman McInnes   }
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
2545e76c431SBarry Smith }
25504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
25604d965bbSLois Curfman McInnes /*
2574b27c08aSLois Curfman McInnes    SNESSetUp_LS - Sets up the internal data structures for the later use
2584b27c08aSLois Curfman McInnes    of the SNESLS nonlinear solver.
25904d965bbSLois Curfman McInnes 
26004d965bbSLois Curfman McInnes    Input Parameter:
26104d965bbSLois Curfman McInnes .  snes - the SNES context
26204d965bbSLois Curfman McInnes .  x - the solution vector
26304d965bbSLois Curfman McInnes 
26404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
26504d965bbSLois Curfman McInnes 
26604d965bbSLois Curfman McInnes    Notes:
2674b27c08aSLois Curfman McInnes    For basic use of the SNES solvers, the user need not explicitly call
26804d965bbSLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
26904d965bbSLois Curfman McInnes    the call to SNESSolve().
27004d965bbSLois Curfman McInnes  */
2714a2ae208SSatish Balay #undef __FUNCT__
2724b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetUp_LS"
273dfbe8321SBarry Smith PetscErrorCode SNESSetUp_LS(SNES snes)
2745e76c431SBarry Smith {
275dfbe8321SBarry Smith   PetscErrorCode ierr;
2763a40ed3dSBarry Smith 
2773a40ed3dSBarry Smith   PetscFunctionBegin;
27858c9b817SLisandro Dalcin   ierr = SNESDefaultGetWork(snes,3);CHKERRQ(ierr);
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
2805e76c431SBarry Smith }
28104d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
2826b8b9a38SLisandro Dalcin 
2836b8b9a38SLisandro Dalcin #undef __FUNCT__
2846b8b9a38SLisandro Dalcin #define __FUNCT__ "SNESReset_LS"
2856b8b9a38SLisandro Dalcin PetscErrorCode SNESReset_LS(SNES snes)
2866b8b9a38SLisandro Dalcin {
2876b8b9a38SLisandro Dalcin   PetscErrorCode ierr;
2886b8b9a38SLisandro Dalcin 
2896b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2906b8b9a38SLisandro Dalcin   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
2916b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
2926b8b9a38SLisandro Dalcin }
2936b8b9a38SLisandro Dalcin 
29404d965bbSLois Curfman McInnes /*
2954b27c08aSLois Curfman McInnes    SNESDestroy_LS - Destroys the private SNES_LS context that was created
2964b27c08aSLois Curfman McInnes    with SNESCreate_LS().
29704d965bbSLois Curfman McInnes 
29804d965bbSLois Curfman McInnes    Input Parameter:
29904d965bbSLois Curfman McInnes .  snes - the SNES context
30004d965bbSLois Curfman McInnes 
301de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30204d965bbSLois Curfman McInnes  */
3034a2ae208SSatish Balay #undef __FUNCT__
3044b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESDestroy_LS"
305dfbe8321SBarry Smith PetscErrorCode SNESDestroy_LS(SNES snes)
3065e76c431SBarry Smith {
307dfbe8321SBarry Smith   PetscErrorCode ierr;
3083a40ed3dSBarry Smith 
3093a40ed3dSBarry Smith   PetscFunctionBegin;
3106b8b9a38SLisandro Dalcin   ierr = SNESReset_LS(snes);CHKERRQ(ierr);
3115c704376SSatish Balay   ierr = PetscFree(snes->data);CHKERRQ(ierr);
312557d3f75SLisandro Dalcin 
3133a40ed3dSBarry Smith   PetscFunctionReturn(0);
3145e76c431SBarry Smith }
31504d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
31694298653SBarry Smith 
317329e7e40SMatthew Knepley /*
3184b27c08aSLois Curfman McInnes    SNESView_LS - Prints info from the SNESLS data structure.
31904d965bbSLois Curfman McInnes 
32004d965bbSLois Curfman McInnes    Input Parameters:
32104d965bbSLois Curfman McInnes .  SNES - the SNES context
32204d965bbSLois Curfman McInnes .  viewer - visualization context
32304d965bbSLois Curfman McInnes 
32404d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
32504d965bbSLois Curfman McInnes */
3264a2ae208SSatish Balay #undef __FUNCT__
3274b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESView_LS"
3286849ba73SBarry Smith static PetscErrorCode SNESView_LS(SNES snes,PetscViewer viewer)
329a935fc98SLois Curfman McInnes {
3302fc52814SBarry Smith   const char     *cstr;
331dfbe8321SBarry Smith   PetscErrorCode ierr;
332ace3abfcSBarry Smith   PetscBool      iascii;
333a935fc98SLois Curfman McInnes 
3343a40ed3dSBarry Smith   PetscFunctionBegin;
3352692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
33632077d6dSBarry Smith   if (iascii) {
337*15f5eeeaSPeter Brune     cstr = SNESLineSearchTypeName(snes->ls_type);
338b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  line search variant: %s\n",cstr);CHKERRQ(ierr);
339ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  alpha=%14.12e, maxstep=%14.12e, minlambda=%14.12e\n",(double)snes->ls_alpha,(double)snes->maxstep,(double)snes->steptol);CHKERRQ(ierr);
340ea630c6eSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"  damping factor=%14.12e\n",(double)snes->damping);CHKERRQ(ierr);
34119bcc07fSBarry Smith   }
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
343a935fc98SLois Curfman McInnes }
34404d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
34504d965bbSLois Curfman McInnes /*
3464b27c08aSLois Curfman McInnes    SNESSetFromOptions_LS - Sets various parameters for the SNESLS method.
34704d965bbSLois Curfman McInnes 
34804d965bbSLois Curfman McInnes    Input Parameter:
34904d965bbSLois Curfman McInnes .  snes - the SNES context
35004d965bbSLois Curfman McInnes 
35104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
35204d965bbSLois Curfman McInnes */
3534a2ae208SSatish Balay #undef __FUNCT__
3544b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESSetFromOptions_LS"
3556849ba73SBarry Smith static PetscErrorCode SNESSetFromOptions_LS(SNES snes)
3565e42470aSBarry Smith {
357dfbe8321SBarry Smith   PetscErrorCode     ierr;
3585e42470aSBarry Smith 
3593a40ed3dSBarry Smith   PetscFunctionBegin;
360b0a32e0cSBarry Smith   ierr = PetscOptionsHead("SNES Line search options");CHKERRQ(ierr);
361b0a32e0cSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3623a40ed3dSBarry Smith   PetscFunctionReturn(0);
3635e42470aSBarry Smith }
3644827ddcaSBarry Smith 
3654827ddcaSBarry Smith EXTERN_C_BEGIN
3664827ddcaSBarry Smith extern PetscErrorCode  SNESLineSearchSetParams_LS(SNES,PetscReal,PetscReal,PetscReal);
3674827ddcaSBarry Smith EXTERN_C_END
3684827ddcaSBarry Smith 
36904d965bbSLois Curfman McInnes /* -------------------------------------------------------------------------- */
370ebe3b25bSBarry Smith /*MC
371ebe3b25bSBarry Smith       SNESLS - Newton based nonlinear solver that uses a line search
37204d965bbSLois Curfman McInnes 
373ebe3b25bSBarry Smith    Options Database:
374ebe3b25bSBarry Smith +   -snes_ls [cubic,quadratic,basic,basicnonorms] - Selects line search
37555d4ea13SBarry Smith .   -snes_ls_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
376e106eecfSBarry Smith .   -snes_ls_maxstep <maxstep> - Sets the maximum stepsize the line search will use (if the 2-norm(y) > maxstep then scale y to be y = (maxstep/2-norm(y)) *y)
377acbee50cSBarry Smith .   -snes_ls_minlambda <minlambda>  - Sets the minimum lambda the line search will use  minlambda / max_i ( y[i]/x[i] )
37855d4ea13SBarry Smith .   -snes_ls_monitor - print information about progress of line searches
37955d4ea13SBarry Smith -   -snes_ls_damping - damping factor used if -snes_ls is basic or basicnonorms
380acbee50cSBarry Smith 
38104d965bbSLois Curfman McInnes 
382ebe3b25bSBarry Smith     Notes: This is the default nonlinear solver in SNES
38304d965bbSLois Curfman McInnes 
384ee3001cbSBarry Smith    Level: beginner
385ee3001cbSBarry Smith 
38617bae607SBarry Smith .seealso:  SNESCreate(), SNES, SNESSetType(), SNESTR, SNESLineSearchSet(),
38717bae607SBarry Smith            SNESLineSearchSetPostCheck(), SNESLineSearchNo(), SNESLineSearchCubic(), SNESLineSearchQuadratic(),
388b3dd4ab5SBarry Smith           SNESLineSearchSet(), SNESLineSearchNoNorms(), SNESLineSearchSetPreCheck(), SNESLineSearchSetParams(), SNESLineSearchGetParams()
389ebe3b25bSBarry Smith 
390ebe3b25bSBarry Smith M*/
391fb2e594dSBarry Smith EXTERN_C_BEGIN
3924a2ae208SSatish Balay #undef __FUNCT__
3934b27c08aSLois Curfman McInnes #define __FUNCT__ "SNESCreate_LS"
3947087cfbeSBarry Smith PetscErrorCode  SNESCreate_LS(SNES snes)
3955e42470aSBarry Smith {
396dfbe8321SBarry Smith   PetscErrorCode ierr;
3974b27c08aSLois Curfman McInnes   SNES_LS        *neP;
3985e42470aSBarry Smith 
3993a40ed3dSBarry Smith   PetscFunctionBegin;
400e7788613SBarry Smith   snes->ops->setup           = SNESSetUp_LS;
401e7788613SBarry Smith   snes->ops->solve           = SNESSolve_LS;
402e7788613SBarry Smith   snes->ops->destroy         = SNESDestroy_LS;
403e7788613SBarry Smith   snes->ops->setfromoptions  = SNESSetFromOptions_LS;
404e7788613SBarry Smith   snes->ops->view            = SNESView_LS;
4056b8b9a38SLisandro Dalcin   snes->ops->reset           = SNESReset_LS;
4065e42470aSBarry Smith 
40742f4f86dSBarry Smith   snes->usesksp                      = PETSC_TRUE;
40842f4f86dSBarry Smith   snes->usespc                       = PETSC_FALSE;
40938f2d2fdSLisandro Dalcin   ierr                               = PetscNewLog(snes,SNES_LS,&neP);CHKERRQ(ierr);
4105e42470aSBarry Smith   snes->data                         = (void*)neP;
411*15f5eeeaSPeter Brune   snes->ops->linesearchno            = SNESLineSearchNo;
412*15f5eeeaSPeter Brune   snes->ops->linesearchnonorms       = SNESLineSearchNoNorms;
413*15f5eeeaSPeter Brune   snes->ops->linesearchquadratic     = SNESLineSearchQuadratic;
414*15f5eeeaSPeter Brune   snes->ops->linesearchcubic         = SNESLineSearchCubic;
415ea630c6eSPeter Brune   snes->ops->linesearchexact         = PETSC_NULL;
416ea630c6eSPeter Brune   snes->ops->linesearchtest          = PETSC_NULL;
417ea630c6eSPeter Brune   snes->lsP                          = PETSC_NULL;
418ea630c6eSPeter Brune   snes->ops->postcheckstep           = PETSC_NULL;
419ea630c6eSPeter Brune   snes->postcheck                    = PETSC_NULL;
420ea630c6eSPeter Brune   snes->ops->precheckstep            = PETSC_NULL;
421ea630c6eSPeter Brune   snes->precheck                     = PETSC_NULL;
422*15f5eeeaSPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_CUBIC);CHKERRQ(ierr);
4233a40ed3dSBarry Smith   PetscFunctionReturn(0);
4245e42470aSBarry Smith }
425fb2e594dSBarry Smith EXTERN_C_END
426