xref: /petsc/src/snes/impls/ls/ls.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 */
109371c9d4SSatish Balay static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin) {
118a5d9ee4SBarry Smith   PetscReal a1;
12ace3abfcSBarry Smith   PetscBool hastranspose;
13f6c5ca78SBarry Smith   Vec       W;
148a5d9ee4SBarry Smith 
158a5d9ee4SBarry Smith   PetscFunctionBegin;
168a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
179566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(F, &W));
1936669109SBarry Smith   if (hastranspose) {
208a5d9ee4SBarry Smith     /* Compute || J^T F|| */
219566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, F, W));
229566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &a1));
239566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm)));
248a5d9ee4SBarry Smith     if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
2536669109SBarry Smith   } else {
2636669109SBarry Smith     Vec         work;
27ea709b57SSatish Balay     PetscScalar result;
2836669109SBarry Smith     PetscReal   wnorm;
2936669109SBarry Smith 
309566063dSJacob Faibussowitsch     PetscCall(VecSetRandom(W, NULL));
319566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &wnorm));
329566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(W, &work));
339566063dSJacob Faibussowitsch     PetscCall(MatMult(A, W, work));
349566063dSJacob Faibussowitsch     PetscCall(VecDot(F, work, &result));
359566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&work));
3636669109SBarry Smith     a1 = PetscAbsScalar(result) / (fnorm * wnorm);
379566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1));
3836669109SBarry Smith     if (a1 < 1.e-4) *ismin = PETSC_TRUE;
3936669109SBarry Smith   }
409566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&W));
418a5d9ee4SBarry Smith   PetscFunctionReturn(0);
428a5d9ee4SBarry Smith }
438a5d9ee4SBarry Smith 
4474637425SBarry Smith /*
455ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
4674637425SBarry Smith */
479371c9d4SSatish Balay static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X) {
4874637425SBarry Smith   PetscReal a1, a2;
49ace3abfcSBarry Smith   PetscBool hastranspose;
5074637425SBarry Smith 
5174637425SBarry Smith   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
5374637425SBarry Smith   if (hastranspose) {
54f6c5ca78SBarry Smith     Vec W1, W2;
55f6c5ca78SBarry Smith 
569566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W1));
579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W2));
589566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, W1));
599566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1, -1.0, F));
6074637425SBarry Smith 
6174637425SBarry Smith     /* Compute || J^T W|| */
629566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, W1, W2));
639566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1, NORM_2, &a1));
649566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2, NORM_2, &a2));
6548a46eb9SPierre Jolivet     if (a1 != 0.0) PetscCall(PetscInfo(snes, "||J^T(F-Ax)||/||F-AX|| %14.12e near zero implies inconsistent rhs\n", (double)(a2 / a1)));
669566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W1));
679566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W2));
6874637425SBarry Smith   }
6974637425SBarry Smith   PetscFunctionReturn(0);
7074637425SBarry Smith }
7174637425SBarry Smith 
72f6dfbefdSBarry Smith /*
7304d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
7494b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
7504d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
7604d965bbSLois Curfman McInnes      respectively.
7704d965bbSLois Curfman McInnes 
78fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
7904d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8004d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
81fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8204d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8304d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
8404d7464bSBarry Smith      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
854b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
8604d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
8704d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
8804d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
8904d965bbSLois Curfman McInnes      for all nonlinear solvers.
9004d965bbSLois Curfman McInnes 
9104d965bbSLois Curfman McInnes      Another key routine is:
9204d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9304d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
9404d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
9504d965bbSLois Curfman McInnes      SNESSolve() if necessary.
9604d965bbSLois Curfman McInnes 
9704d965bbSLois Curfman McInnes      Additional basic routines are:
9804d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
9904d965bbSLois Curfman McInnes                                       have actually been used.
10004d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
101186905e3SBarry Smith      SNESView().
10204d965bbSLois Curfman McInnes 
10304d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
10404d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
10504d965bbSLois Curfman McInnes      above description applies to these categories also.
10604d965bbSLois Curfman McInnes 
107f6dfbefdSBarry Smith */
1085e42470aSBarry Smith /*
10904d7464bSBarry Smith    SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
11004d965bbSLois Curfman McInnes    method with a line search.
1115e76c431SBarry Smith 
11204d965bbSLois Curfman McInnes    Input Parameters:
11304d965bbSLois Curfman McInnes .  snes - the SNES context
1145e76c431SBarry Smith 
11504d965bbSLois Curfman McInnes    Output Parameter:
116c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
11704d965bbSLois Curfman McInnes 
11804d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1195e76c431SBarry Smith 
1205e76c431SBarry Smith */
1219371c9d4SSatish Balay PetscErrorCode SNESSolve_NEWTONLS(SNES snes) {
122a7cc72afSBarry Smith   PetscInt             maxits, i, lits;
123422a814eSBarry Smith   SNESLineSearchReason lssucceed;
1241936c542SPeter Brune   PetscReal            fnorm, gnorm, xnorm, ynorm;
125f6c5ca78SBarry Smith   Vec                  Y, X, F;
126f1c6b773SPeter Brune   SNESLineSearch       linesearch;
127c40d0f55SPeter Brune   SNESConvergedReason  reason;
1285e76c431SBarry Smith 
1293a40ed3dSBarry Smith   PetscFunctionBegin;
1300b121fc5SBarry Smith   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
131c579b300SPatrick Farrell 
1323d4c4710SBarry Smith   snes->numFailures            = 0;
1333d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
134184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
135184914b5SBarry Smith 
1365e42470aSBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
1375e42470aSBarry Smith   X      = snes->vec_sol;        /* solution vector */
13839e2f89bSBarry Smith   F      = snes->vec_func;       /* residual vector */
139d12b9ff3SPeter Brune   Y      = snes->vec_sol_update; /* newton step */
1405e76c431SBarry Smith 
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1424c49b128SBarry Smith   snes->iter = 0;
14375567043SBarry Smith   snes->norm = 0.0;
1449566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1459566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
146ce5a860aSPeter Brune 
147ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
148efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1499566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1509566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
151b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
152b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
153b7281c8aSPeter Brune       PetscFunctionReturn(0);
154b7281c8aSPeter Brune     }
155b7281c8aSPeter Brune 
1569566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(F, NORM_2, &fnorm));
1579566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(F, NORM_2, &fnorm));
158ce5a860aSPeter Brune   } else {
159e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1609566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1611aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
162c1c75074SPeter Brune   }
1631aa26658SKarl Rupp 
1649566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
165422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
1669566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1675e42470aSBarry Smith   snes->norm = fnorm;
1689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1699566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1709566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
1713f149594SLisandro Dalcin 
17285385478SLisandro Dalcin   /* test convergence */
173dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
17406ee9f85SBarry Smith   if (snes->reason) PetscFunctionReturn(0);
175d034289bSBarry Smith 
1765e76c431SBarry Smith   for (i = 0; i < maxits; i++) {
177329e7e40SMatthew Knepley     /* Call general purpose update function */
178dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
179329e7e40SMatthew Knepley 
180ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
181efd4aadfSBarry Smith     if (snes->npc) {
182efd4aadfSBarry Smith       if (snes->npcside == PC_RIGHT) {
1839566063dSJacob Faibussowitsch         PetscCall(SNESSetInitialFunction(snes->npc, F));
1849566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1859566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1869566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1879566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
188c40d0f55SPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
189c40d0f55SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
190c40d0f55SPeter Brune           PetscFunctionReturn(0);
191c40d0f55SPeter Brune         }
1929566063dSJacob Faibussowitsch         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
193efd4aadfSBarry Smith       } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
1949566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, F));
1959566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
196b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
197b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
198b7281c8aSPeter Brune           PetscFunctionReturn(0);
199b7281c8aSPeter Brune         }
200ce5a860aSPeter Brune       }
201c40d0f55SPeter Brune     }
202c40d0f55SPeter Brune 
203ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2049566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
20507b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
2069566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2079566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, Y));
208422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2099566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
21063a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
21174637425SBarry Smith 
2121baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y));
21374637425SBarry Smith 
214ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
215d12b9ff3SPeter Brune          X <- X - lambda*Y
216d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
217ea4d3ed3SLois Curfman McInnes     */
2181936c542SPeter Brune     gnorm = fnorm;
2199566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y));
2209566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed));
2219566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
2229566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed));
2237e49c775SBarry Smith     if (snes->reason) break;
224422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
225422a814eSBarry Smith     if (lssucceed) {
226c21ba15cSPeter Brune       if (snes->stol * xnorm > ynorm) {
227c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
228c21ba15cSPeter Brune         PetscFunctionReturn(0);
229c21ba15cSPeter Brune       }
23050ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
231ace3abfcSBarry Smith         PetscBool ismin;
232647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2339566063dSJacob Faibussowitsch         PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin));
2348a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2353add74b1SBarry Smith         if (snes->errorifnotconverged && snes->reason) {
2363add74b1SBarry Smith           PetscViewer monitor;
2379566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
23828b400f6SJacob Faibussowitsch           PetscCheck(monitor, PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due to %s. Suggest running with -snes_linesearch_monitor", SNESConvergedReasons[snes->reason]);
239f7d195e4SLawrence Mitchell           SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]);
2403add74b1SBarry Smith         }
2413505fcc1SBarry Smith         break;
2423505fcc1SBarry Smith       }
24350ffb88aSMatthew Knepley     }
24485385478SLisandro Dalcin     /* Monitor convergence */
2459566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
24685385478SLisandro Dalcin     snes->iter  = i + 1;
24785385478SLisandro Dalcin     snes->norm  = fnorm;
248c1e67a49SFande Kong     snes->ynorm = ynorm;
249c1e67a49SFande Kong     snes->xnorm = xnorm;
2509566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2519566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
2529566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2532a997342SPeter Brune     /* Test for convergence */
254dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
2553f149594SLisandro Dalcin     if (snes->reason) break;
25616c95cacSBarry Smith   }
25752392280SLois Curfman McInnes   if (i == maxits) {
25863a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
25985385478SLisandro Dalcin     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
26052392280SLois Curfman McInnes   }
2613a40ed3dSBarry Smith   PetscFunctionReturn(0);
2625e76c431SBarry Smith }
263f6dfbefdSBarry Smith 
26404d965bbSLois Curfman McInnes /*
26504d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
26604d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
26704d965bbSLois Curfman McInnes 
26804d965bbSLois Curfman McInnes    Input Parameter:
26904d965bbSLois Curfman McInnes .  snes - the SNES context
27004d965bbSLois Curfman McInnes .  x - the solution vector
27104d965bbSLois Curfman McInnes 
27204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
27304d965bbSLois Curfman McInnes 
27404d965bbSLois Curfman McInnes  */
2759371c9d4SSatish Balay PetscErrorCode SNESSetUp_NEWTONLS(SNES snes) {
2763a40ed3dSBarry Smith   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
278efd4aadfSBarry Smith   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
2793a40ed3dSBarry Smith   PetscFunctionReturn(0);
2805e76c431SBarry Smith }
2816b8b9a38SLisandro Dalcin 
2829371c9d4SSatish Balay PetscErrorCode SNESReset_NEWTONLS(SNES snes) {
2836b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2846b8b9a38SLisandro Dalcin   PetscFunctionReturn(0);
2856b8b9a38SLisandro Dalcin }
2866b8b9a38SLisandro Dalcin 
28704d965bbSLois Curfman McInnes /*
28804d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
28904d7464bSBarry Smith    with SNESCreate_NEWTONLS().
29004d965bbSLois Curfman McInnes 
29104d965bbSLois Curfman McInnes    Input Parameter:
29204d965bbSLois Curfman McInnes .  snes - the SNES context
29304d965bbSLois Curfman McInnes 
294de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
29504d965bbSLois Curfman McInnes  */
2969371c9d4SSatish Balay PetscErrorCode SNESDestroy_NEWTONLS(SNES snes) {
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONLS(snes));
2999566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3003a40ed3dSBarry Smith   PetscFunctionReturn(0);
3015e76c431SBarry Smith }
30294298653SBarry Smith 
303329e7e40SMatthew Knepley /*
30404d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
30504d965bbSLois Curfman McInnes 
30604d965bbSLois Curfman McInnes    Input Parameters:
30704d965bbSLois Curfman McInnes .  SNES - the SNES context
30804d965bbSLois Curfman McInnes .  viewer - visualization context
30904d965bbSLois Curfman McInnes 
31004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
31104d965bbSLois Curfman McInnes */
3129371c9d4SSatish Balay static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer) {
313ace3abfcSBarry Smith   PetscBool iascii;
314a935fc98SLois Curfman McInnes 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
3169566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3179371c9d4SSatish Balay   if (iascii) { }
31850f6de3fSJed Brown   PetscFunctionReturn(0);
31950f6de3fSJed Brown }
32050f6de3fSJed Brown 
32104d965bbSLois Curfman McInnes /*
32204d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
32304d965bbSLois Curfman McInnes 
32404d965bbSLois Curfman McInnes    Input Parameter:
32504d965bbSLois Curfman McInnes .  snes - the SNES context
32604d965bbSLois Curfman McInnes 
32704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
32804d965bbSLois Curfman McInnes */
3299371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject) {
3303a40ed3dSBarry Smith   PetscFunctionBegin;
3313a40ed3dSBarry Smith   PetscFunctionReturn(0);
3325e42470aSBarry Smith }
3334827ddcaSBarry Smith 
334ebe3b25bSBarry Smith /*MC
33504d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
33604d965bbSLois Curfman McInnes 
337ebe3b25bSBarry Smith    Options Database:
3385a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
3395a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
340f6dfbefdSBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`)
34182d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
34282d26c24SPeter Brune .   -snes_linesearch_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)
3435a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
34482d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
34582d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
346acbee50cSBarry Smith 
347f6dfbefdSBarry Smith     Note:
348f6dfbefdSBarry Smith     This is the default nonlinear solver in `SNES`
34904d965bbSLois Curfman McInnes 
350ee3001cbSBarry Smith    Level: beginner
351ee3001cbSBarry Smith 
352db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()`
353f6dfbefdSBarry Smith           `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()`
354ebe3b25bSBarry Smith 
355ebe3b25bSBarry Smith M*/
3569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes) {
35704d7464bSBarry Smith   SNES_NEWTONLS *neP;
358d8d34be6SBarry Smith   SNESLineSearch linesearch;
3595e42470aSBarry Smith 
3603a40ed3dSBarry Smith   PetscFunctionBegin;
36104d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
36204d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
36304d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
36404d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
36504d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
36604d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
3675e42470aSBarry Smith 
368efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
36942f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
370efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
3714fc747eaSLawrence Mitchell 
3729566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
37348a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
374d8d34be6SBarry Smith 
3754fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
3764fc747eaSLawrence Mitchell 
377*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
3785e42470aSBarry Smith   snes->data = (void *)neP;
3793a40ed3dSBarry Smith   PetscFunctionReturn(0);
3805e42470aSBarry Smith }
381