xref: /petsc/src/snes/impls/ls/ls.c (revision 728696cf4dbfc82d2e79ed325f353b9daffcc864)
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 */
10d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckLocalMin_Private(SNES snes, Mat A, Vec F, PetscReal fnorm, PetscBool *ismin)
11d71ae5a4SJacob Faibussowitsch {
128a5d9ee4SBarry Smith   PetscReal a1;
13ace3abfcSBarry Smith   PetscBool hastranspose;
14f6c5ca78SBarry Smith   Vec       W;
15*728696cfSStefano Zampini   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
168a5d9ee4SBarry Smith 
178a5d9ee4SBarry Smith   PetscFunctionBegin;
188a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
19*728696cfSStefano Zampini   PetscCall(SNESGetObjective(snes, &objective, NULL));
20*728696cfSStefano Zampini   if (!objective) {
219566063dSJacob Faibussowitsch     PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
229566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W));
2336669109SBarry Smith     if (hastranspose) {
248a5d9ee4SBarry Smith       /* Compute || J^T F|| */
259566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(A, F, W));
269566063dSJacob Faibussowitsch       PetscCall(VecNorm(W, NORM_2, &a1));
279566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "|| J^T F|| %14.12e near zero implies found a local minimum\n", (double)(a1 / fnorm)));
288a5d9ee4SBarry Smith       if (a1 / fnorm < 1.e-4) *ismin = PETSC_TRUE;
2936669109SBarry Smith     } else {
3036669109SBarry Smith       Vec         work;
31ea709b57SSatish Balay       PetscScalar result;
3236669109SBarry Smith       PetscReal   wnorm;
3336669109SBarry Smith 
349566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(W, NULL));
359566063dSJacob Faibussowitsch       PetscCall(VecNorm(W, NORM_2, &wnorm));
369566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(W, &work));
379566063dSJacob Faibussowitsch       PetscCall(MatMult(A, W, work));
389566063dSJacob Faibussowitsch       PetscCall(VecDot(F, work, &result));
399566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&work));
4036669109SBarry Smith       a1 = PetscAbsScalar(result) / (fnorm * wnorm);
419566063dSJacob Faibussowitsch       PetscCall(PetscInfo(snes, "(F^T J random)/(|| F ||*||random|| %14.12e near zero implies found a local minimum\n", (double)a1));
4236669109SBarry Smith       if (a1 < 1.e-4) *ismin = PETSC_TRUE;
4336669109SBarry Smith     }
449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W));
45*728696cfSStefano Zampini   }
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
478a5d9ee4SBarry Smith }
488a5d9ee4SBarry Smith 
4974637425SBarry Smith /*
505ed2d596SBarry Smith      Checks if J^T(F - J*X) = 0
5174637425SBarry Smith */
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNEWTONLSCheckResidual_Private(SNES snes, Mat A, Vec F, Vec X)
53d71ae5a4SJacob Faibussowitsch {
5474637425SBarry Smith   PetscReal a1, a2;
55ace3abfcSBarry Smith   PetscBool hastranspose;
5674637425SBarry Smith 
5774637425SBarry Smith   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
5974637425SBarry Smith   if (hastranspose) {
60f6c5ca78SBarry Smith     Vec W1, W2;
61f6c5ca78SBarry Smith 
629566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W1));
639566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W2));
649566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, W1));
659566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1, -1.0, F));
6674637425SBarry Smith 
6774637425SBarry Smith     /* Compute || J^T W|| */
689566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, W1, W2));
699566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1, NORM_2, &a1));
709566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2, NORM_2, &a2));
7148a46eb9SPierre 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)));
729566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W1));
739566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W2));
7474637425SBarry Smith   }
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7674637425SBarry Smith }
7774637425SBarry Smith 
78f6dfbefdSBarry Smith /*
7904d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
8094b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
8104d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8204d965bbSLois Curfman McInnes      respectively.
8304d965bbSLois Curfman McInnes 
84fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8504d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8604d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
87fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
8804d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
8904d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
9004d7464bSBarry Smith      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
914b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
9204d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9304d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9404d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9504d965bbSLois Curfman McInnes      for all nonlinear solvers.
9604d965bbSLois Curfman McInnes 
9704d965bbSLois Curfman McInnes      Another key routine is:
9804d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
9904d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
10004d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
10104d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10204d965bbSLois Curfman McInnes 
10304d965bbSLois Curfman McInnes      Additional basic routines are:
10404d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10504d965bbSLois Curfman McInnes                                       have actually been used.
10604d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
107186905e3SBarry Smith      SNESView().
10804d965bbSLois Curfman McInnes 
10904d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
11004d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
11104d965bbSLois Curfman McInnes      above description applies to these categories also.
11204d965bbSLois Curfman McInnes 
113f6dfbefdSBarry Smith */
1145e42470aSBarry Smith /*
11504d7464bSBarry Smith    SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
11604d965bbSLois Curfman McInnes    method with a line search.
1175e76c431SBarry Smith 
11804d965bbSLois Curfman McInnes    Input Parameters:
11904d965bbSLois Curfman McInnes .  snes - the SNES context
1205e76c431SBarry Smith 
12104d965bbSLois Curfman McInnes    Output Parameter:
122c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12304d965bbSLois Curfman McInnes 
12404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1255e76c431SBarry Smith 
1265e76c431SBarry Smith */
127d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
128d71ae5a4SJacob Faibussowitsch {
129a7cc72afSBarry Smith   PetscInt             maxits, i, lits;
130422a814eSBarry Smith   SNESLineSearchReason lssucceed;
1314279555eSSatish Balay   PetscReal            fnorm, xnorm, ynorm;
132f6c5ca78SBarry Smith   Vec                  Y, X, F;
133f1c6b773SPeter Brune   SNESLineSearch       linesearch;
134c40d0f55SPeter Brune   SNESConvergedReason  reason;
1354279555eSSatish Balay #if defined(PETSC_USE_INFO)
1364279555eSSatish Balay   PetscReal gnorm;
1374279555eSSatish Balay #endif
1385e76c431SBarry Smith 
1393a40ed3dSBarry Smith   PetscFunctionBegin;
1400b121fc5SBarry 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);
141c579b300SPatrick Farrell 
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 */
149d12b9ff3SPeter Brune   Y      = snes->vec_sol_update; /* newton step */
1505e76c431SBarry Smith 
1519566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1524c49b128SBarry Smith   snes->iter = 0;
15375567043SBarry Smith   snes->norm = 0.0;
1549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1559566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
156ce5a860aSPeter Brune 
157ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
158efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1599566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1609566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
161b7281c8aSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
162b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1633ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
164b7281c8aSPeter Brune     }
165b7281c8aSPeter Brune 
1669566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(F, NORM_2, &fnorm));
1679566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(F, NORM_2, &fnorm));
168ce5a860aSPeter Brune   } else {
169e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1709566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1711aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
172c1c75074SPeter Brune   }
1731aa26658SKarl Rupp 
1749566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
175422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
1769566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1775e42470aSBarry Smith   snes->norm = fnorm;
1789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1799566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1809566063dSJacob Faibussowitsch   PetscCall(SNESMonitor(snes, 0, fnorm));
1813f149594SLisandro Dalcin 
18285385478SLisandro Dalcin   /* test convergence */
183dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
1843ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
185d034289bSBarry Smith 
1865e76c431SBarry Smith   for (i = 0; i < maxits; i++) {
187329e7e40SMatthew Knepley     /* Call general purpose update function */
188dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
189329e7e40SMatthew Knepley 
190ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
191efd4aadfSBarry Smith     if (snes->npc) {
192efd4aadfSBarry Smith       if (snes->npcside == PC_RIGHT) {
1939566063dSJacob Faibussowitsch         PetscCall(SNESSetInitialFunction(snes->npc, F));
1949566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1959566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1969566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1979566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
198c40d0f55SPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
199c40d0f55SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
2003ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
201c40d0f55SPeter Brune         }
2029566063dSJacob Faibussowitsch         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
203efd4aadfSBarry Smith       } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2049566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, F));
2059566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
206b7281c8aSPeter Brune         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
207b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
2083ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
209b7281c8aSPeter Brune         }
210ce5a860aSPeter Brune       }
211c40d0f55SPeter Brune     }
212c40d0f55SPeter Brune 
213ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2149566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
21507b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
2169566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2179566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, Y));
218422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2199566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
22063a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
22174637425SBarry Smith 
2221baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y));
22374637425SBarry Smith 
2244279555eSSatish Balay #if defined(PETSC_USE_INFO)
2254279555eSSatish Balay     gnorm = fnorm;
2264279555eSSatish Balay #endif
227ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
228d12b9ff3SPeter Brune          X <- X - lambda*Y
229d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
230ea4d3ed3SLois Curfman McInnes     */
2319566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y));
2329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed));
2339566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
2349566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed));
2357e49c775SBarry Smith     if (snes->reason) break;
236422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
237422a814eSBarry Smith     if (lssucceed) {
238c21ba15cSPeter Brune       if (snes->stol * xnorm > ynorm) {
239c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
2403ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
241c21ba15cSPeter Brune       }
24250ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
243ace3abfcSBarry Smith         PetscBool ismin;
244647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2459566063dSJacob Faibussowitsch         PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin));
2468a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2473add74b1SBarry Smith         if (snes->errorifnotconverged && snes->reason) {
2483add74b1SBarry Smith           PetscViewer monitor;
2499566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
25028b400f6SJacob 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]);
251f7d195e4SLawrence Mitchell           SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]);
2523add74b1SBarry Smith         }
2533505fcc1SBarry Smith         break;
2543505fcc1SBarry Smith       }
25550ffb88aSMatthew Knepley     }
25685385478SLisandro Dalcin     /* Monitor convergence */
2579566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
25885385478SLisandro Dalcin     snes->iter  = i + 1;
25985385478SLisandro Dalcin     snes->norm  = fnorm;
260c1e67a49SFande Kong     snes->ynorm = ynorm;
261c1e67a49SFande Kong     snes->xnorm = xnorm;
2629566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2639566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
2649566063dSJacob Faibussowitsch     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2652a997342SPeter Brune     /* Test for convergence */
266dbbe0bcdSBarry Smith     PetscUseTypeMethod(snes, converged, snes->iter, xnorm, ynorm, fnorm, &snes->reason, snes->cnvP);
2673f149594SLisandro Dalcin     if (snes->reason) break;
26816c95cacSBarry Smith   }
26952392280SLois Curfman McInnes   if (i == maxits) {
27063a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", maxits));
27185385478SLisandro Dalcin     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
27252392280SLois Curfman McInnes   }
2733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2745e76c431SBarry Smith }
275f6dfbefdSBarry Smith 
27604d965bbSLois Curfman McInnes /*
27704d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
27804d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
27904d965bbSLois Curfman McInnes 
28004d965bbSLois Curfman McInnes    Input Parameter:
28104d965bbSLois Curfman McInnes .  snes - the SNES context
28204d965bbSLois Curfman McInnes .  x - the solution vector
28304d965bbSLois Curfman McInnes 
28404d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
28504d965bbSLois Curfman McInnes 
28604d965bbSLois Curfman McInnes  */
287d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
288d71ae5a4SJacob Faibussowitsch {
2893a40ed3dSBarry Smith   PetscFunctionBegin;
2909566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
291efd4aadfSBarry Smith   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2935e76c431SBarry Smith }
2946b8b9a38SLisandro Dalcin 
295d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONLS(SNES snes)
296d71ae5a4SJacob Faibussowitsch {
2976b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2996b8b9a38SLisandro Dalcin }
3006b8b9a38SLisandro Dalcin 
30104d965bbSLois Curfman McInnes /*
30204d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
30304d7464bSBarry Smith    with SNESCreate_NEWTONLS().
30404d965bbSLois Curfman McInnes 
30504d965bbSLois Curfman McInnes    Input Parameter:
30604d965bbSLois Curfman McInnes .  snes - the SNES context
30704d965bbSLois Curfman McInnes 
308de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30904d965bbSLois Curfman McInnes  */
310d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
311d71ae5a4SJacob Faibussowitsch {
3123a40ed3dSBarry Smith   PetscFunctionBegin;
3139566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONLS(snes));
3149566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3165e76c431SBarry Smith }
31794298653SBarry Smith 
318329e7e40SMatthew Knepley /*
31904d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
32004d965bbSLois Curfman McInnes 
32104d965bbSLois Curfman McInnes    Input Parameters:
32204d965bbSLois Curfman McInnes .  SNES - the SNES context
32304d965bbSLois Curfman McInnes .  viewer - visualization context
32404d965bbSLois Curfman McInnes 
32504d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
32604d965bbSLois Curfman McInnes */
327d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer)
328d71ae5a4SJacob Faibussowitsch {
329ace3abfcSBarry Smith   PetscBool iascii;
330a935fc98SLois Curfman McInnes 
3313a40ed3dSBarry Smith   PetscFunctionBegin;
3329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3339371c9d4SSatish Balay   if (iascii) { }
3343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33550f6de3fSJed Brown }
33650f6de3fSJed Brown 
33704d965bbSLois Curfman McInnes /*
33804d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
33904d965bbSLois Curfman McInnes 
34004d965bbSLois Curfman McInnes    Input Parameter:
34104d965bbSLois Curfman McInnes .  snes - the SNES context
34204d965bbSLois Curfman McInnes 
34304d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
34404d965bbSLois Curfman McInnes */
345d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject)
346d71ae5a4SJacob Faibussowitsch {
3473a40ed3dSBarry Smith   PetscFunctionBegin;
3483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3495e42470aSBarry Smith }
3504827ddcaSBarry Smith 
351ebe3b25bSBarry Smith /*MC
35204d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
35304d965bbSLois Curfman McInnes 
3543c7db156SBarry Smith    Options Database Keys:
3555a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
3565a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
357f6dfbefdSBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`)
35882d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
35982d26c24SPeter 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)
3605a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
36182d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
36282d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
363acbee50cSBarry Smith 
364f6dfbefdSBarry Smith     Note:
365f6dfbefdSBarry Smith     This is the default nonlinear solver in `SNES`
36604d965bbSLois Curfman McInnes 
367ee3001cbSBarry Smith    Level: beginner
368ee3001cbSBarry Smith 
369db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()`
370f6dfbefdSBarry Smith           `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()`
371ebe3b25bSBarry Smith M*/
372d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
373d71ae5a4SJacob Faibussowitsch {
37404d7464bSBarry Smith   SNES_NEWTONLS *neP;
375d8d34be6SBarry Smith   SNESLineSearch linesearch;
3765e42470aSBarry Smith 
3773a40ed3dSBarry Smith   PetscFunctionBegin;
37804d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
37904d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
38004d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
38104d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
38204d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
38304d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
3845e42470aSBarry Smith 
385efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
38642f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
387efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
3884fc747eaSLawrence Mitchell 
3899566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
39048a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
391d8d34be6SBarry Smith 
3924fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
3934fc747eaSLawrence Mitchell 
3944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
3955e42470aSBarry Smith   snes->data = (void *)neP;
3963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3975e42470aSBarry Smith }
398