xref: /petsc/src/snes/impls/ls/ls.c (revision ceaaa4989964adb3f5eb130cb04b8f49c83e49c9)
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;
15728696cfSStefano Zampini   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
168a5d9ee4SBarry Smith 
178a5d9ee4SBarry Smith   PetscFunctionBegin;
188a5d9ee4SBarry Smith   *ismin = PETSC_FALSE;
19728696cfSStefano Zampini   PetscCall(SNESGetObjective(snes, &objective, NULL));
20728696cfSStefano 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));
45728696cfSStefano 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;
56b13f2d03SStefano Zampini   PetscErrorCode (*objective)(SNES, Vec, PetscReal *, void *);
5774637425SBarry Smith 
5874637425SBarry Smith   PetscFunctionBegin;
599566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_MULT_TRANSPOSE, &hastranspose));
60b13f2d03SStefano Zampini   PetscCall(SNESGetObjective(snes, &objective, NULL));
61b13f2d03SStefano Zampini   if (hastranspose && !objective) {
62f6c5ca78SBarry Smith     Vec W1, W2;
63f6c5ca78SBarry Smith 
649566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W1));
659566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(F, &W2));
669566063dSJacob Faibussowitsch     PetscCall(MatMult(A, X, W1));
679566063dSJacob Faibussowitsch     PetscCall(VecAXPY(W1, -1.0, F));
6874637425SBarry Smith 
6974637425SBarry Smith     /* Compute || J^T W|| */
709566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(A, W1, W2));
719566063dSJacob Faibussowitsch     PetscCall(VecNorm(W1, NORM_2, &a1));
729566063dSJacob Faibussowitsch     PetscCall(VecNorm(W2, NORM_2, &a2));
7348a46eb9SPierre 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)));
749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W1));
759566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&W2));
7674637425SBarry Smith   }
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7874637425SBarry Smith }
7974637425SBarry Smith 
80f6dfbefdSBarry Smith /*
8104d965bbSLois Curfman McInnes      This file implements a truncated Newton method with a line search,
8294b7f48cSBarry Smith      for solving a system of nonlinear equations, using the KSP, Vec,
8304d965bbSLois Curfman McInnes      and Mat interfaces for linear solvers, vectors, and matrices,
8404d965bbSLois Curfman McInnes      respectively.
8504d965bbSLois Curfman McInnes 
86fe6c6bacSLois Curfman McInnes      The following basic routines are required for each nonlinear solver:
8704d965bbSLois Curfman McInnes           SNESCreate_XXX()          - Creates a nonlinear solver context
8804d965bbSLois Curfman McInnes           SNESSetFromOptions_XXX()  - Sets runtime options
89fe6c6bacSLois Curfman McInnes           SNESSolve_XXX()           - Solves the nonlinear system
9004d965bbSLois Curfman McInnes           SNESDestroy_XXX()         - Destroys the nonlinear solver context
9104d965bbSLois Curfman McInnes      The suffix "_XXX" denotes a particular implementation, in this case
9204d7464bSBarry Smith      we use _NEWTONLS (e.g., SNESCreate_NEWTONLS, SNESSolve_NEWTONLS) for solving
934b27c08aSLois Curfman McInnes      systems of nonlinear equations with a line search (LS) method.
9404d965bbSLois Curfman McInnes      These routines are actually called via the common user interface
9504d965bbSLois Curfman McInnes      routines SNESCreate(), SNESSetFromOptions(), SNESSolve(), and
9604d965bbSLois Curfman McInnes      SNESDestroy(), so the application code interface remains identical
9704d965bbSLois Curfman McInnes      for all nonlinear solvers.
9804d965bbSLois Curfman McInnes 
9904d965bbSLois Curfman McInnes      Another key routine is:
10004d965bbSLois Curfman McInnes           SNESSetUp_XXX()           - Prepares for the use of a nonlinear solver
10104d965bbSLois Curfman McInnes      by setting data structures and options.   The interface routine SNESSetUp()
10204d965bbSLois Curfman McInnes      is not usually called directly by the user, but instead is called by
10304d965bbSLois Curfman McInnes      SNESSolve() if necessary.
10404d965bbSLois Curfman McInnes 
10504d965bbSLois Curfman McInnes      Additional basic routines are:
10604d965bbSLois Curfman McInnes           SNESView_XXX()            - Prints details of runtime options that
10704d965bbSLois Curfman McInnes                                       have actually been used.
10804d965bbSLois Curfman McInnes      These are called by application codes via the interface routines
109186905e3SBarry Smith      SNESView().
11004d965bbSLois Curfman McInnes 
11104d965bbSLois Curfman McInnes      The various types of solvers (preconditioners, Krylov subspace methods,
11204d965bbSLois Curfman McInnes      nonlinear solvers, timesteppers) are all organized similarly, so the
11304d965bbSLois Curfman McInnes      above description applies to these categories also.
11404d965bbSLois Curfman McInnes 
115f6dfbefdSBarry Smith */
116*ceaaa498SBarry Smith 
117*ceaaa498SBarry Smith // PetscClangLinter pragma disable: -fdoc-sowing-chars
1185e42470aSBarry Smith /*
11904d7464bSBarry Smith   SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
12004d965bbSLois Curfman McInnes   method with a line search.
1215e76c431SBarry Smith 
1222fe279fdSBarry Smith   Input Parameter:
12304d965bbSLois Curfman McInnes . snes - the SNES context
1245e76c431SBarry Smith 
1255e76c431SBarry Smith */
126d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
127d71ae5a4SJacob Faibussowitsch {
128a7cc72afSBarry Smith   PetscInt             maxits, i, lits;
129422a814eSBarry Smith   SNESLineSearchReason lssucceed;
1304279555eSSatish Balay   PetscReal            fnorm, xnorm, ynorm;
131f6c5ca78SBarry Smith   Vec                  Y, X, F;
132f1c6b773SPeter Brune   SNESLineSearch       linesearch;
133c40d0f55SPeter Brune   SNESConvergedReason  reason;
1344279555eSSatish Balay #if defined(PETSC_USE_INFO)
1354279555eSSatish Balay   PetscReal gnorm;
1364279555eSSatish Balay #endif
1375e76c431SBarry Smith 
1383a40ed3dSBarry Smith   PetscFunctionBegin;
1390b121fc5SBarry 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);
140c579b300SPatrick Farrell 
1413d4c4710SBarry Smith   snes->numFailures            = 0;
1423d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
143184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
144184914b5SBarry Smith 
1455e42470aSBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
1465e42470aSBarry Smith   X      = snes->vec_sol;        /* solution vector */
14739e2f89bSBarry Smith   F      = snes->vec_func;       /* residual vector */
148d12b9ff3SPeter Brune   Y      = snes->vec_sol_update; /* newton step */
1495e76c431SBarry Smith 
1509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1514c49b128SBarry Smith   snes->iter = 0;
15275567043SBarry Smith   snes->norm = 0.0;
1539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1549566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
155ce5a860aSPeter Brune 
156ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
157efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1589566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1599566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
160eee3b80eSStefano Zampini     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
161b7281c8aSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1623ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
163b7281c8aSPeter Brune     }
164b7281c8aSPeter Brune 
1659566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(F, NORM_2, &fnorm));
1669566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(F, NORM_2, &fnorm));
167ce5a860aSPeter Brune   } else {
168e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1699566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1701aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
171c1c75074SPeter Brune   }
1721aa26658SKarl Rupp 
1739566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
174422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
1759566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1765e42470aSBarry Smith   snes->norm = fnorm;
1779566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1789566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1793f149594SLisandro Dalcin 
18085385478SLisandro Dalcin   /* test convergence */
181d76a863bSStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
1822d157150SStefano Zampini   PetscCall(SNESMonitor(snes, 0, fnorm));
1833ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
184d034289bSBarry Smith 
1855e76c431SBarry Smith   for (i = 0; i < maxits; i++) {
186329e7e40SMatthew Knepley     /* Call general purpose update function */
187dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
188329e7e40SMatthew Knepley 
189ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
190efd4aadfSBarry Smith     if (snes->npc) {
191efd4aadfSBarry Smith       if (snes->npcside == PC_RIGHT) {
1929566063dSJacob Faibussowitsch         PetscCall(SNESSetInitialFunction(snes->npc, F));
1939566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1949566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1959566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1969566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
197eee3b80eSStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
198c40d0f55SPeter Brune           snes->reason = SNES_DIVERGED_INNER;
1993ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
200c40d0f55SPeter Brune         }
2019566063dSJacob Faibussowitsch         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
202efd4aadfSBarry Smith       } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2039566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, F));
2049566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
205eee3b80eSStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
206b7281c8aSPeter Brune           snes->reason = SNES_DIVERGED_INNER;
2073ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
208b7281c8aSPeter Brune         }
209ce5a860aSPeter Brune       }
210c40d0f55SPeter Brune     }
211c40d0f55SPeter Brune 
212ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2139566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
21407b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
2159566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2169566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, Y));
217422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2189566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
21963a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
22074637425SBarry Smith 
2211baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y));
22274637425SBarry Smith 
2234279555eSSatish Balay #if defined(PETSC_USE_INFO)
2244279555eSSatish Balay     gnorm = fnorm;
2254279555eSSatish Balay #endif
226ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
227d12b9ff3SPeter Brune          X <- X - lambda*Y
228d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
229ea4d3ed3SLois Curfman McInnes     */
2309566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y));
2319566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed));
2329566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
2339566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed));
2347e49c775SBarry Smith     if (snes->reason) break;
235422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
236422a814eSBarry Smith     if (lssucceed) {
237c21ba15cSPeter Brune       if (snes->stol * xnorm > ynorm) {
238c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
2393ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
240c21ba15cSPeter Brune       }
24150ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
242ace3abfcSBarry Smith         PetscBool ismin;
243647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2449566063dSJacob Faibussowitsch         PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin));
2458a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2463add74b1SBarry Smith         if (snes->errorifnotconverged && snes->reason) {
2473add74b1SBarry Smith           PetscViewer monitor;
2489566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
24928b400f6SJacob 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]);
250f7d195e4SLawrence Mitchell           SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]);
2513add74b1SBarry Smith         }
2523505fcc1SBarry Smith         break;
2533505fcc1SBarry Smith       }
25450ffb88aSMatthew Knepley     }
25585385478SLisandro Dalcin     /* Monitor convergence */
2569566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
25785385478SLisandro Dalcin     snes->iter  = i + 1;
25885385478SLisandro Dalcin     snes->norm  = fnorm;
259c1e67a49SFande Kong     snes->ynorm = ynorm;
260c1e67a49SFande Kong     snes->xnorm = xnorm;
2619566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2629566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
2632a997342SPeter Brune     /* Test for convergence */
2642d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
2652d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2663f149594SLisandro Dalcin     if (snes->reason) break;
26716c95cacSBarry Smith   }
2683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2695e76c431SBarry Smith }
270f6dfbefdSBarry Smith 
27104d965bbSLois Curfman McInnes /*
27204d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
27304d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
27404d965bbSLois Curfman McInnes 
27504d965bbSLois Curfman McInnes    Input Parameter:
27604d965bbSLois Curfman McInnes .  snes - the SNES context
27704d965bbSLois Curfman McInnes .  x - the solution vector
27804d965bbSLois Curfman McInnes 
27904d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
28004d965bbSLois Curfman McInnes 
28104d965bbSLois Curfman McInnes  */
282d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
283d71ae5a4SJacob Faibussowitsch {
2843a40ed3dSBarry Smith   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
286efd4aadfSBarry Smith   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
2873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2885e76c431SBarry Smith }
2896b8b9a38SLisandro Dalcin 
290d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONLS(SNES snes)
291d71ae5a4SJacob Faibussowitsch {
2926b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2946b8b9a38SLisandro Dalcin }
2956b8b9a38SLisandro Dalcin 
29604d965bbSLois Curfman McInnes /*
29704d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
29804d7464bSBarry Smith    with SNESCreate_NEWTONLS().
29904d965bbSLois Curfman McInnes 
30004d965bbSLois Curfman McInnes    Input Parameter:
30104d965bbSLois Curfman McInnes .  snes - the SNES context
30204d965bbSLois Curfman McInnes 
303de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30404d965bbSLois Curfman McInnes  */
305d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
306d71ae5a4SJacob Faibussowitsch {
3073a40ed3dSBarry Smith   PetscFunctionBegin;
3089566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONLS(snes));
3099566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3115e76c431SBarry Smith }
31294298653SBarry Smith 
313329e7e40SMatthew Knepley /*
31404d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
31504d965bbSLois Curfman McInnes 
31604d965bbSLois Curfman McInnes    Input Parameters:
31704d965bbSLois Curfman McInnes .  SNES - the SNES context
31804d965bbSLois Curfman McInnes .  viewer - visualization context
31904d965bbSLois Curfman McInnes 
32004d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
32104d965bbSLois Curfman McInnes */
322d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer)
323d71ae5a4SJacob Faibussowitsch {
324ace3abfcSBarry Smith   PetscBool iascii;
325a935fc98SLois Curfman McInnes 
3263a40ed3dSBarry Smith   PetscFunctionBegin;
3279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3289371c9d4SSatish Balay   if (iascii) { }
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33050f6de3fSJed Brown }
33150f6de3fSJed Brown 
33204d965bbSLois Curfman McInnes /*
33304d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
33404d965bbSLois Curfman McInnes 
33504d965bbSLois Curfman McInnes    Input Parameter:
33604d965bbSLois Curfman McInnes .  snes - the SNES context
33704d965bbSLois Curfman McInnes 
33804d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
33904d965bbSLois Curfman McInnes */
340d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject)
341d71ae5a4SJacob Faibussowitsch {
3423a40ed3dSBarry Smith   PetscFunctionBegin;
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3445e42470aSBarry Smith }
3454827ddcaSBarry Smith 
346ebe3b25bSBarry Smith /*MC
34704d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
34804d965bbSLois Curfman McInnes 
3493c7db156SBarry Smith    Options Database Keys:
3505a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
3515a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
352f6dfbefdSBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`)
35382d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
35482d26c24SPeter 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)
3555a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
35682d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
35782d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
358acbee50cSBarry Smith 
359f6dfbefdSBarry Smith     Note:
360f6dfbefdSBarry Smith     This is the default nonlinear solver in `SNES`
36104d965bbSLois Curfman McInnes 
362ee3001cbSBarry Smith    Level: beginner
363ee3001cbSBarry Smith 
364db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()`
365f6dfbefdSBarry Smith           `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()`
366ebe3b25bSBarry Smith M*/
367d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
368d71ae5a4SJacob Faibussowitsch {
36904d7464bSBarry Smith   SNES_NEWTONLS *neP;
370d8d34be6SBarry Smith   SNESLineSearch linesearch;
3715e42470aSBarry Smith 
3723a40ed3dSBarry Smith   PetscFunctionBegin;
37304d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
37404d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
37504d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
37604d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
37704d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
37804d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
3795e42470aSBarry Smith 
380efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
38142f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
382efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
3834fc747eaSLawrence Mitchell 
3849566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
38548a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
386d8d34be6SBarry Smith 
3874fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
3884fc747eaSLawrence Mitchell 
3894dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
3905e42470aSBarry Smith   snes->data = (void *)neP;
3913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3925e42470aSBarry Smith }
393