xref: /petsc/src/snes/impls/ls/ls.c (revision 420bcc1b905230dede3c88f397d1a4e60493adde)
1c6db04a5SJed Brown #include <../src/snes/impls/ls/lsimpl.h>
25e42470aSBarry Smith 
38a5d9ee4SBarry Smith /*
4*420bcc1bSBarry Smith      This file implements a truncated Newton method with a line search,
5*420bcc1bSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
6*420bcc1bSBarry Smith      and Mat interfaces for linear solvers, vectors, and matrices,
7*420bcc1bSBarry Smith      respectively.
8*420bcc1bSBarry Smith 
9*420bcc1bSBarry Smith      The following basic routines are required for each nonlinear solver:
10*420bcc1bSBarry Smith           SNESCreate_XXX()          - Creates a nonlinear solver context
11*420bcc1bSBarry Smith           SNESSetFromOptions_XXX()  - Sets runtime options
12*420bcc1bSBarry Smith           SNESSolve_XXX()           - Solves the nonlinear system
13*420bcc1bSBarry Smith           SNESDestroy_XXX()         - Destroys the nonlinear solver context
14*420bcc1bSBarry Smith      The suffix "_XXX" denotes a particular implementation, in this case
15*420bcc1bSBarry Smith      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
16*420bcc1bSBarry Smith      systems of nonlinear equations with a line search (LS) method.
17*420bcc1bSBarry Smith      These routines are actually called via the common user interface
18*420bcc1bSBarry Smith      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
19*420bcc1bSBarry Smith      SNESDestroy(), so the application code interface remains identical
20*420bcc1bSBarry Smith      for all nonlinear solvers.
21*420bcc1bSBarry Smith 
22*420bcc1bSBarry Smith      Another key routine is:
23*420bcc1bSBarry Smith           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
24*420bcc1bSBarry Smith      by setting data structures and options.   The interface routine SNESSetUp()
25*420bcc1bSBarry Smith      is not usually called directly by the user, but instead is called by
26*420bcc1bSBarry Smith      SNESSolve() if necessary.
27*420bcc1bSBarry Smith 
28*420bcc1bSBarry Smith      Additional basic routines are:
29*420bcc1bSBarry Smith           SNESView_XXX()            - Prints details of runtime options that
30*420bcc1bSBarry Smith                                       have actually been used.
31*420bcc1bSBarry Smith      These are called by application codes via the interface routines
32*420bcc1bSBarry Smith      SNESView().
33*420bcc1bSBarry Smith 
34*420bcc1bSBarry Smith      The various types of solvers (preconditioners, Krylov subspace methods,
35*420bcc1bSBarry Smith      nonlinear solvers, timesteppers) are all organized similarly, so the
36*420bcc1bSBarry Smith      above description applies to these categories also.
37*420bcc1bSBarry Smith 
38*420bcc1bSBarry Smith */
39*420bcc1bSBarry Smith 
40*420bcc1bSBarry Smith /*
4120f52e01SBarry Smith      Checks if J^T F = 0 which implies we've found a local minimum of the norm of the function,
4220f52e01SBarry 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
4336669109SBarry Smith     0 = (J^T F)^T W = F^T J W iff W not in the null space of J. Thanks for Jorge More
4420f52e01SBarry Smith     for this trick. One assumes that the probability that W is in the null space of J is very, very small.
458a5d9ee4SBarry Smith */
46d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin)
47d71ae5a4SJacob Faibussowitsch {
488a5d9ee4SBarry Smith   PetscReal a1;
49ace3abfcSBarry Smith   PetscBool hastranspose;
50f6c5ca78SBarry Smith   Vec       W;
51728696cfSStefano Zampini   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
528a5d9ee4SBarry Smith 
538a5d9ee4SBarry Smith   PetscFunctionBegin;
548a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
55728696cfSStefano Zampini   PetscCall(SNESGetObjective(snes, &objective, NULL));
56728696cfSStefano Zampini   if (!objective) {
579566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
589566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W));
5936669109SBarry Smith     if (hastranspose) {
608a5d9ee4SBarry Smith       /* Compute || J^T F|| */
619566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A, F, W));
629566063dSJacob Faibussowitsch       PetscCall(VecNorm(W, NORM_2, &a1));
639566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm)));
648a5d9ee4SBarry Smith       if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
6536669109SBarry Smith     } else {
6636669109SBarry Smith       Vec         work;
67ea709b57SSatish Balay       PetscScalar result;
6836669109SBarry Smith       PetscReal   wnorm;
6936669109SBarry Smith 
709566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(W, NULL));
719566063dSJacob Faibussowitsch       PetscCall(VecNorm(W, NORM_2, &wnorm));
729566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(W, &work));
739566063dSJacob Faibussowitsch       PetscCall(MatMult(A, W, work));
749566063dSJacob Faibussowitsch       PetscCall(VecDot(F, work, &result));
759566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work));
7636669109SBarry Smith       a1 = PetscAbsScalar(result) / (fnorm * wnorm);
779566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1));
7836669109SBarry Smith       if (a1 < 1.e-4) *ismin = PETSC_TRUE;
7936669109SBarry Smith     }
809566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W));
81728696cfSStefano Zampini   }
823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
838a5d9ee4SBarry Smith }
848a5d9ee4SBarry Smith 
8574637425SBarry Smith /*
865ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
8774637425SBarry Smith */
88d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X)
89d71ae5a4SJacob Faibussowitsch {
9074637425SBarry Smith   PetscReal a1, a2;
91ace3abfcSBarry Smith   PetscBool hastranspose;
92b13f2d03SStefano Zampini   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
9374637425SBarry Smith 
9474637425SBarry Smith   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
96b13f2d03SStefano Zampini   PetscCall(SNESGetObjective(snes, &objective, NULL));
97b13f2d03SStefano Zampini   if (hastranspose && !objective) {
98f6c5ca78SBarry Smith     Vec W1, W2;
99f6c5ca78SBarry Smith 
1009566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W1));
1019566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W2));
1029566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, W1));
1039566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1, -1.0, F));
10474637425SBarry Smith 
10574637425SBarry Smith     /* Compute || J^T W|| */
1069566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, W1, W2));
1079566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1, NORM_2, &a1));
1089566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2, NORM_2, &a2));
10948a46eb9SPierre 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)));
1109566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W1));
1119566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W2));
11274637425SBarry Smith   }
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11474637425SBarry Smith }
11574637425SBarry Smith 
116ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars
1175e42470aSBarry Smith /*
11804d7464bSBarry Smith   SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
11904d965bbSLois Curfman McInnes   method with a line search.
1205e76c431SBarry Smith 
1212fe279fdSBarry Smith   Input Parameter:
12204d965bbSLois Curfman McInnes . snes - the SNES context
1235e76c431SBarry Smith 
1245e76c431SBarry Smith */
12566976f2fSJacob Faibussowitsch static PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
126d71ae5a4SJacob Faibussowitsch {
127a7cc72afSBarry Smith   PetscInt             maxits, i, lits;
128422a814eSBarry Smith   SNESLineSearchReason lssucceed;
1294279555eSSatish Balay   PetscReal            fnorm, xnorm, ynorm;
130f6c5ca78SBarry Smith   Vec                  Y, X, F;
131f1c6b773SPeter Brune   SNESLineSearch       linesearch;
132c40d0f55SPeter Brune   SNESConvergedReason  reason;
1334279555eSSatish Balay #if defined(PETSC_USE_INFO)
1344279555eSSatish Balay   PetscReal gnorm;
1354279555eSSatish Balay #endif
1365e76c431SBarry Smith 
1373a40ed3dSBarry Smith   PetscFunctionBegin;
1380b121fc5SBarry 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);
139c579b300SPatrick Farrell 
1403d4c4710SBarry Smith   snes->numFailures            = 0;
1413d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
142184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
143184914b5SBarry Smith 
1445e42470aSBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
1455e42470aSBarry Smith   X      = snes->vec_sol;        /* solution vector */
14639e2f89bSBarry Smith   F      = snes->vec_func;       /* residual vector */
147d12b9ff3SPeter Brune   Y      = snes->vec_sol_update; /* newton step */
1485e76c431SBarry Smith 
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1504c49b128SBarry Smith   snes->iter = 0;
15175567043SBarry Smith   snes->norm = 0.0;
1529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1539566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
154ce5a860aSPeter Brune 
155ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
156efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1579566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1589566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
159eee3b80eSStefano Zampini     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
1605f3c5e7aSBarry Smith       PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
1613ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
162b7281c8aSPeter Brune     }
163b7281c8aSPeter Brune 
1649566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(F, NORM_2, &fnorm));
1659566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(F, NORM_2, &fnorm));
166ce5a860aSPeter Brune   } else {
167e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1689566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1691aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
170c1c75074SPeter Brune   }
1711aa26658SKarl Rupp 
1729566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
173422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
1749566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1755e42470aSBarry Smith   snes->norm = fnorm;
1769566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1779566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1783f149594SLisandro Dalcin 
17985385478SLisandro Dalcin   /* test convergence */
180d76a863bSStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
1812d157150SStefano Zampini   PetscCall(SNESMonitor(snes, 0, fnorm));
1823ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
183d034289bSBarry Smith 
1845e76c431SBarry Smith   for (i = 0; i < maxits; i++) {
185329e7e40SMatthew Knepley     /* Call general purpose update function */
186dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
187329e7e40SMatthew Knepley 
188ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
189efd4aadfSBarry Smith     if (snes->npc) {
190efd4aadfSBarry Smith       if (snes->npcside == PC_RIGHT) {
1919566063dSJacob Faibussowitsch         PetscCall(SNESSetInitialFunction(snes->npc, F));
1929566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1939566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1949566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1959566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
196eee3b80eSStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
1975f3c5e7aSBarry Smith           PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
1983ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
199c40d0f55SPeter Brune         }
2009566063dSJacob Faibussowitsch         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
201efd4aadfSBarry Smith       } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2029566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, F));
2039566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
204eee3b80eSStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
2055f3c5e7aSBarry Smith           PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
2063ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
207b7281c8aSPeter Brune         }
208ce5a860aSPeter Brune       }
209c40d0f55SPeter Brune     }
210c40d0f55SPeter Brune 
211ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2129566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
21307b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
2149566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2159566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, Y));
216422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2179566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
21863a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
21974637425SBarry Smith 
2201baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y));
22174637425SBarry Smith 
2224279555eSSatish Balay #if defined(PETSC_USE_INFO)
2234279555eSSatish Balay     gnorm = fnorm;
2244279555eSSatish Balay #endif
225ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
226d12b9ff3SPeter Brune          X <- X - lambda*Y
227d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
228ea4d3ed3SLois Curfman McInnes     */
2299566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y));
2309566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed));
2319566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
2329566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed));
2337e49c775SBarry Smith     if (snes->reason) break;
234422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
235422a814eSBarry Smith     if (lssucceed) {
236c21ba15cSPeter Brune       if (snes->stol * xnorm > ynorm) {
237c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
2383ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
239c21ba15cSPeter Brune       }
24050ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
241ace3abfcSBarry Smith         PetscBool ismin;
242647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2439566063dSJacob Faibussowitsch         PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin));
2448a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2453add74b1SBarry Smith         if (snes->errorifnotconverged && snes->reason) {
2463add74b1SBarry Smith           PetscViewer monitor;
2479566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
24828b400f6SJacob 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]);
249f7d195e4SLawrence Mitchell           SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]);
2503add74b1SBarry Smith         }
2513505fcc1SBarry Smith         break;
2523505fcc1SBarry Smith       }
25350ffb88aSMatthew Knepley     }
25485385478SLisandro Dalcin     /* Monitor convergence */
2559566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
25685385478SLisandro Dalcin     snes->iter  = i + 1;
25785385478SLisandro Dalcin     snes->norm  = fnorm;
258c1e67a49SFande Kong     snes->ynorm = ynorm;
259c1e67a49SFande Kong     snes->xnorm = xnorm;
2609566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2619566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
2622a997342SPeter Brune     /* Test for convergence */
2632d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
2642d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2653f149594SLisandro Dalcin     if (snes->reason) break;
26616c95cacSBarry Smith   }
2673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2685e76c431SBarry Smith }
269f6dfbefdSBarry Smith 
27004d965bbSLois Curfman McInnes /*
27104d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
27204d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
27304d965bbSLois Curfman McInnes 
27404d965bbSLois Curfman McInnes    Input Parameter:
27504d965bbSLois Curfman McInnes .  snes - the SNES context
27604d965bbSLois Curfman McInnes .  x - the solution vector
27704d965bbSLois Curfman McInnes 
27804d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
27904d965bbSLois Curfman McInnes 
28004d965bbSLois Curfman McInnes  */
28166976f2fSJacob Faibussowitsch static PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
282d71ae5a4SJacob Faibussowitsch {
2833a40ed3dSBarry Smith   PetscFunctionBegin;
2849566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
285efd4aadfSBarry Smith   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2875e76c431SBarry Smith }
2886b8b9a38SLisandro Dalcin 
28966976f2fSJacob Faibussowitsch static PetscErrorCode SNESReset_NEWTONLS(SNES snes)
290d71ae5a4SJacob Faibussowitsch {
2916b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2936b8b9a38SLisandro Dalcin }
2946b8b9a38SLisandro Dalcin 
29504d965bbSLois Curfman McInnes /*
29604d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
29704d7464bSBarry Smith    with SNESCreate_NEWTONLS().
29804d965bbSLois Curfman McInnes 
29904d965bbSLois Curfman McInnes    Input Parameter:
30004d965bbSLois Curfman McInnes .  snes - the SNES context
30104d965bbSLois Curfman McInnes 
302de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30304d965bbSLois Curfman McInnes  */
30466976f2fSJacob Faibussowitsch static PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
305d71ae5a4SJacob Faibussowitsch {
3063a40ed3dSBarry Smith   PetscFunctionBegin;
3079566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONLS(snes));
3089566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3105e76c431SBarry Smith }
31194298653SBarry Smith 
312329e7e40SMatthew Knepley /*
31304d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
31404d965bbSLois Curfman McInnes 
31504d965bbSLois Curfman McInnes    Input Parameters:
31604d965bbSLois Curfman McInnes .  SNES - the SNES context
31704d965bbSLois Curfman McInnes .  viewer - visualization context
31804d965bbSLois Curfman McInnes 
31904d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
32004d965bbSLois Curfman McInnes */
321d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer)
322d71ae5a4SJacob Faibussowitsch {
323ace3abfcSBarry Smith   PetscBool iascii;
324a935fc98SLois Curfman McInnes 
3253a40ed3dSBarry Smith   PetscFunctionBegin;
3269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3279371c9d4SSatish Balay   if (iascii) { }
3283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32950f6de3fSJed Brown }
33050f6de3fSJed Brown 
33104d965bbSLois Curfman McInnes /*
33204d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
33304d965bbSLois Curfman McInnes 
33404d965bbSLois Curfman McInnes    Input Parameter:
33504d965bbSLois Curfman McInnes .  snes - the SNES context
33604d965bbSLois Curfman McInnes 
33704d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
33804d965bbSLois Curfman McInnes */
339d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject)
340d71ae5a4SJacob Faibussowitsch {
3413a40ed3dSBarry Smith   PetscFunctionBegin;
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3435e42470aSBarry Smith }
3444827ddcaSBarry Smith 
345ebe3b25bSBarry Smith /*MC
34604d7464bSBarry Smith    SNESNEWTONLS - Newton based nonlinear solver that uses a line search
34704d965bbSLois Curfman McInnes 
3483c7db156SBarry Smith    Options Database Keys:
3495a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
3505a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
351f6dfbefdSBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`)
35282d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
35382d26c24SPeter 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)
3545a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
35582d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
35682d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
357acbee50cSBarry Smith 
358*420bcc1bSBarry Smith    Level: beginner
359*420bcc1bSBarry Smith 
360f6dfbefdSBarry Smith    Note:
361f6dfbefdSBarry Smith    This is the default nonlinear solver in `SNES`
36204d965bbSLois Curfman McInnes 
363*420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()`
364f6dfbefdSBarry Smith           `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()`
365ebe3b25bSBarry Smith M*/
366d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
367d71ae5a4SJacob Faibussowitsch {
36804d7464bSBarry Smith   SNES_NEWTONLS *neP;
369d8d34be6SBarry Smith   SNESLineSearch linesearch;
3705e42470aSBarry Smith 
3713a40ed3dSBarry Smith   PetscFunctionBegin;
37204d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
37304d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
37404d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
37504d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
37604d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
37704d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
3785e42470aSBarry Smith 
379efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
38042f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
381efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
3824fc747eaSLawrence Mitchell 
3839566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
38448a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
385d8d34be6SBarry Smith 
3864fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
3874fc747eaSLawrence Mitchell 
3884dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
3895e42470aSBarry Smith   snes->data = (void *)neP;
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3915e42470aSBarry Smith }
392