xref: /petsc/src/snes/impls/ls/ls.c (revision 5f3c5e7ab1713b2b36ec2007ece43899b4f0dcb3)
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 */
1165e42470aSBarry Smith /*
11704d7464bSBarry Smith   SNESSolve_NEWTONLS - Solves a nonlinear system with a truncated Newton
11804d965bbSLois Curfman McInnes   method with a line search.
1195e76c431SBarry Smith 
1202fe279fdSBarry Smith   Input Parameter:
12104d965bbSLois Curfman McInnes . snes - the SNES context
1225e76c431SBarry Smith 
12304d965bbSLois Curfman McInnes   Output Parameter:
124c2a78f1aSLois Curfman McInnes .  outits - number of iterations until termination
12504d965bbSLois Curfman McInnes 
12604d965bbSLois Curfman McInnes    Application Interface Routine: SNESSolve()
1275e76c431SBarry Smith 
1285e76c431SBarry Smith */
129d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSolve_NEWTONLS(SNES snes)
130d71ae5a4SJacob Faibussowitsch {
131a7cc72afSBarry Smith   PetscInt             maxits, i, lits;
132422a814eSBarry Smith   SNESLineSearchReason lssucceed;
1334279555eSSatish Balay   PetscReal            fnorm, xnorm, ynorm;
134f6c5ca78SBarry Smith   Vec                  Y, X, F;
135f1c6b773SPeter Brune   SNESLineSearch       linesearch;
136c40d0f55SPeter Brune   SNESConvergedReason  reason;
1374279555eSSatish Balay #if defined(PETSC_USE_INFO)
1384279555eSSatish Balay   PetscReal gnorm;
1394279555eSSatish Balay #endif
1405e76c431SBarry Smith 
1413a40ed3dSBarry Smith   PetscFunctionBegin;
1420b121fc5SBarry 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);
143c579b300SPatrick Farrell 
1443d4c4710SBarry Smith   snes->numFailures            = 0;
1453d4c4710SBarry Smith   snes->numLinearSolveFailures = 0;
146184914b5SBarry Smith   snes->reason                 = SNES_CONVERGED_ITERATING;
147184914b5SBarry Smith 
1485e42470aSBarry Smith   maxits = snes->max_its;        /* maximum number of iterations */
1495e42470aSBarry Smith   X      = snes->vec_sol;        /* solution vector */
15039e2f89bSBarry Smith   F      = snes->vec_func;       /* residual vector */
151d12b9ff3SPeter Brune   Y      = snes->vec_sol_update; /* newton step */
1525e76c431SBarry Smith 
1539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1544c49b128SBarry Smith   snes->iter = 0;
15575567043SBarry Smith   snes->norm = 0.0;
1569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1579566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
158ce5a860aSPeter Brune 
159ce5a860aSPeter Brune   /* compute the preconditioned function first in the case of left preconditioning with preconditioned function */
160efd4aadfSBarry Smith   if (snes->npc && snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_PRECONDITIONED) {
1619566063dSJacob Faibussowitsch     PetscCall(SNESApplyNPC(snes, X, NULL, F));
1629566063dSJacob Faibussowitsch     PetscCall(SNESGetConvergedReason(snes->npc, &reason));
163eee3b80eSStefano Zampini     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
164*5f3c5e7aSBarry Smith       PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
1653ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
166b7281c8aSPeter Brune     }
167b7281c8aSPeter Brune 
1689566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(F, NORM_2, &fnorm));
1699566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(F, NORM_2, &fnorm));
170ce5a860aSPeter Brune   } else {
171e4ed7901SPeter Brune     if (!snes->vec_func_init_set) {
1729566063dSJacob Faibussowitsch       PetscCall(SNESComputeFunction(snes, X, F));
1731aa26658SKarl Rupp     } else snes->vec_func_init_set = PETSC_FALSE;
174c1c75074SPeter Brune   }
1751aa26658SKarl Rupp 
1769566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
177422a814eSBarry Smith   SNESCheckFunctionNorm(snes, fnorm);
1789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
1795e42470aSBarry Smith   snes->norm = fnorm;
1809566063dSJacob Faibussowitsch   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
1819566063dSJacob Faibussowitsch   PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0));
1823f149594SLisandro Dalcin 
18385385478SLisandro Dalcin   /* test convergence */
184d76a863bSStefano Zampini   PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm));
1852d157150SStefano Zampini   PetscCall(SNESMonitor(snes, 0, fnorm));
1863ba16761SJacob Faibussowitsch   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
187d034289bSBarry Smith 
1885e76c431SBarry Smith   for (i = 0; i < maxits; i++) {
189329e7e40SMatthew Knepley     /* Call general purpose update function */
190dbbe0bcdSBarry Smith     PetscTryTypeMethod(snes, update, snes->iter);
191329e7e40SMatthew Knepley 
192ce5a860aSPeter Brune     /* apply the nonlinear preconditioner */
193efd4aadfSBarry Smith     if (snes->npc) {
194efd4aadfSBarry Smith       if (snes->npcside == PC_RIGHT) {
1959566063dSJacob Faibussowitsch         PetscCall(SNESSetInitialFunction(snes->npc, F));
1969566063dSJacob Faibussowitsch         PetscCall(PetscLogEventBegin(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1979566063dSJacob Faibussowitsch         PetscCall(SNESSolve(snes->npc, snes->vec_rhs, X));
1989566063dSJacob Faibussowitsch         PetscCall(PetscLogEventEnd(SNES_NPCSolve, snes->npc, X, snes->vec_rhs, 0));
1999566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
200eee3b80eSStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
201*5f3c5e7aSBarry Smith           PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
2023ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
203c40d0f55SPeter Brune         }
2049566063dSJacob Faibussowitsch         PetscCall(SNESGetNPCFunction(snes, F, &fnorm));
205efd4aadfSBarry Smith       } else if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_UNPRECONDITIONED) {
2069566063dSJacob Faibussowitsch         PetscCall(SNESApplyNPC(snes, X, F, F));
2079566063dSJacob Faibussowitsch         PetscCall(SNESGetConvergedReason(snes->npc, &reason));
208eee3b80eSStefano Zampini         if (reason < 0 && reason != SNES_DIVERGED_MAX_IT && reason != SNES_DIVERGED_TR_DELTA) {
209*5f3c5e7aSBarry Smith           PetscCall(SNESSetConvergedReason(snes, SNES_DIVERGED_INNER));
2103ba16761SJacob Faibussowitsch           PetscFunctionReturn(PETSC_SUCCESS);
211b7281c8aSPeter Brune         }
212ce5a860aSPeter Brune       }
213c40d0f55SPeter Brune     }
214c40d0f55SPeter Brune 
215ea4d3ed3SLois Curfman McInnes     /* Solve J Y = F, where J is Jacobian matrix */
2169566063dSJacob Faibussowitsch     PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre));
21707b62357SFande Kong     SNESCheckJacobianDomainerror(snes);
2189566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
2199566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, F, Y));
220422a814eSBarry Smith     SNESCheckKSPSolve(snes);
2219566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(snes->ksp, &lits));
22263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "iter=%" PetscInt_FMT ", linear solve iterations=%" PetscInt_FMT "\n", snes->iter, lits));
22374637425SBarry Smith 
2241baa6e33SBarry Smith     if (PetscLogPrintInfo) PetscCall(SNESNEWTONLSCheckResidual_Private(snes, snes->jacobian, F, Y));
22574637425SBarry Smith 
2264279555eSSatish Balay #if defined(PETSC_USE_INFO)
2274279555eSSatish Balay     gnorm = fnorm;
2284279555eSSatish Balay #endif
229ea4d3ed3SLois Curfman McInnes     /* Compute a (scaled) negative update in the line search routine:
230d12b9ff3SPeter Brune          X <- X - lambda*Y
231d12b9ff3SPeter Brune        and evaluate F = function(X) (depends on the line search).
232ea4d3ed3SLois Curfman McInnes     */
2339566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchApply(linesearch, X, F, &fnorm, Y));
2349566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetReason(linesearch, &lssucceed));
2359566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
2369566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "fnorm=%18.16e, gnorm=%18.16e, ynorm=%18.16e, lssucceed=%d\n", (double)gnorm, (double)fnorm, (double)ynorm, (int)lssucceed));
2377e49c775SBarry Smith     if (snes->reason) break;
238422a814eSBarry Smith     SNESCheckFunctionNorm(snes, fnorm);
239422a814eSBarry Smith     if (lssucceed) {
240c21ba15cSPeter Brune       if (snes->stol * xnorm > ynorm) {
241c21ba15cSPeter Brune         snes->reason = SNES_CONVERGED_SNORM_RELATIVE;
2423ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
243c21ba15cSPeter Brune       }
24450ffb88aSMatthew Knepley       if (++snes->numFailures >= snes->maxFailures) {
245ace3abfcSBarry Smith         PetscBool ismin;
246647a2e1fSBarry Smith         snes->reason = SNES_DIVERGED_LINE_SEARCH;
2479566063dSJacob Faibussowitsch         PetscCall(SNESNEWTONLSCheckLocalMin_Private(snes, snes->jacobian, F, fnorm, &ismin));
2488a5d9ee4SBarry Smith         if (ismin) snes->reason = SNES_DIVERGED_LOCAL_MIN;
2493add74b1SBarry Smith         if (snes->errorifnotconverged && snes->reason) {
2503add74b1SBarry Smith           PetscViewer monitor;
2519566063dSJacob Faibussowitsch           PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
25228b400f6SJacob 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]);
253f7d195e4SLawrence Mitchell           SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_NOT_CONVERGED, "SNESSolve has not converged due %s.", SNESConvergedReasons[snes->reason]);
2543add74b1SBarry Smith         }
2553505fcc1SBarry Smith         break;
2563505fcc1SBarry Smith       }
25750ffb88aSMatthew Knepley     }
25885385478SLisandro Dalcin     /* Monitor convergence */
2599566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
26085385478SLisandro Dalcin     snes->iter  = i + 1;
26185385478SLisandro Dalcin     snes->norm  = fnorm;
262c1e67a49SFande Kong     snes->ynorm = ynorm;
263c1e67a49SFande Kong     snes->xnorm = xnorm;
2649566063dSJacob Faibussowitsch     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
2659566063dSJacob Faibussowitsch     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, lits));
2662a997342SPeter Brune     /* Test for convergence */
2672d157150SStefano Zampini     PetscCall(SNESConverged(snes, snes->iter, xnorm, ynorm, fnorm));
2682d157150SStefano Zampini     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
2693f149594SLisandro Dalcin     if (snes->reason) break;
27016c95cacSBarry Smith   }
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2725e76c431SBarry Smith }
273f6dfbefdSBarry Smith 
27404d965bbSLois Curfman McInnes /*
27504d7464bSBarry Smith    SNESSetUp_NEWTONLS - Sets up the internal data structures for the later use
27604d7464bSBarry Smith    of the SNESNEWTONLS nonlinear solver.
27704d965bbSLois Curfman McInnes 
27804d965bbSLois Curfman McInnes    Input Parameter:
27904d965bbSLois Curfman McInnes .  snes - the SNES context
28004d965bbSLois Curfman McInnes .  x - the solution vector
28104d965bbSLois Curfman McInnes 
28204d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetUp()
28304d965bbSLois Curfman McInnes 
28404d965bbSLois Curfman McInnes  */
285d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetUp_NEWTONLS(SNES snes)
286d71ae5a4SJacob Faibussowitsch {
2873a40ed3dSBarry Smith   PetscFunctionBegin;
2889566063dSJacob Faibussowitsch   PetscCall(SNESSetUpMatrices(snes));
289efd4aadfSBarry Smith   if (snes->npcside == PC_LEFT && snes->functype == SNES_FUNCTION_DEFAULT) snes->functype = SNES_FUNCTION_PRECONDITIONED;
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2915e76c431SBarry Smith }
2926b8b9a38SLisandro Dalcin 
293d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESReset_NEWTONLS(SNES snes)
294d71ae5a4SJacob Faibussowitsch {
2956b8b9a38SLisandro Dalcin   PetscFunctionBegin;
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2976b8b9a38SLisandro Dalcin }
2986b8b9a38SLisandro Dalcin 
29904d965bbSLois Curfman McInnes /*
30004d7464bSBarry Smith    SNESDestroy_NEWTONLS - Destroys the private SNES_NEWTONLS context that was created
30104d7464bSBarry Smith    with SNESCreate_NEWTONLS().
30204d965bbSLois Curfman McInnes 
30304d965bbSLois Curfman McInnes    Input Parameter:
30404d965bbSLois Curfman McInnes .  snes - the SNES context
30504d965bbSLois Curfman McInnes 
306de49b36dSLois Curfman McInnes    Application Interface Routine: SNESDestroy()
30704d965bbSLois Curfman McInnes  */
308d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESDestroy_NEWTONLS(SNES snes)
309d71ae5a4SJacob Faibussowitsch {
3103a40ed3dSBarry Smith   PetscFunctionBegin;
3119566063dSJacob Faibussowitsch   PetscCall(SNESReset_NEWTONLS(snes));
3129566063dSJacob Faibussowitsch   PetscCall(PetscFree(snes->data));
3133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3145e76c431SBarry Smith }
31594298653SBarry Smith 
316329e7e40SMatthew Knepley /*
31704d7464bSBarry Smith    SNESView_NEWTONLS - Prints info from the SNESNEWTONLS data structure.
31804d965bbSLois Curfman McInnes 
31904d965bbSLois Curfman McInnes    Input Parameters:
32004d965bbSLois Curfman McInnes .  SNES - the SNES context
32104d965bbSLois Curfman McInnes .  viewer - visualization context
32204d965bbSLois Curfman McInnes 
32304d965bbSLois Curfman McInnes    Application Interface Routine: SNESView()
32404d965bbSLois Curfman McInnes */
325d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NEWTONLS(SNES snes, PetscViewer viewer)
326d71ae5a4SJacob Faibussowitsch {
327ace3abfcSBarry Smith   PetscBool iascii;
328a935fc98SLois Curfman McInnes 
3293a40ed3dSBarry Smith   PetscFunctionBegin;
3309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3319371c9d4SSatish Balay   if (iascii) { }
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33350f6de3fSJed Brown }
33450f6de3fSJed Brown 
33504d965bbSLois Curfman McInnes /*
33604d7464bSBarry Smith    SNESSetFromOptions_NEWTONLS - Sets various parameters for the SNESNEWTONLS method.
33704d965bbSLois Curfman McInnes 
33804d965bbSLois Curfman McInnes    Input Parameter:
33904d965bbSLois Curfman McInnes .  snes - the SNES context
34004d965bbSLois Curfman McInnes 
34104d965bbSLois Curfman McInnes    Application Interface Routine: SNESSetFromOptions()
34204d965bbSLois Curfman McInnes */
343d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NEWTONLS(SNES snes, PetscOptionItems *PetscOptionsObject)
344d71ae5a4SJacob Faibussowitsch {
3453a40ed3dSBarry Smith   PetscFunctionBegin;
3463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3475e42470aSBarry Smith }
3484827ddcaSBarry Smith 
349ebe3b25bSBarry Smith /*MC
35004d7464bSBarry Smith       SNESNEWTONLS - Newton based nonlinear solver that uses a line search
35104d965bbSLois Curfman McInnes 
3523c7db156SBarry Smith    Options Database Keys:
3535a9b6599SPeter Brune +   -snes_linesearch_type <bt> - bt,basic.  Select line search type
3545a9b6599SPeter Brune .   -snes_linesearch_order <3> - 2, 3. Selects the order of the line search for bt
355f6dfbefdSBarry Smith .   -snes_linesearch_norms <true> - Turns on/off computation of the norms for basic linesearch (`SNESLineSearchSetComputeNorms()`)
35682d26c24SPeter Brune .   -snes_linesearch_alpha <alpha> - Sets alpha used in determining if reduction in function norm is sufficient
35782d26c24SPeter 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)
3585a9b6599SPeter Brune .   -snes_linesearch_minlambda <minlambda>  - Sets the minimum lambda the line search will tolerate
35982d26c24SPeter Brune .   -snes_linesearch_monitor - print information about progress of line searches
36082d26c24SPeter Brune -   -snes_linesearch_damping - damping factor used for basic line search
361acbee50cSBarry Smith 
362f6dfbefdSBarry Smith     Note:
363f6dfbefdSBarry Smith     This is the default nonlinear solver in `SNES`
36404d965bbSLois Curfman McInnes 
365ee3001cbSBarry Smith    Level: beginner
366ee3001cbSBarry Smith 
367db781477SPatrick Sanan .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESNEWTONTR`, `SNESQN`, `SNESLineSearchSetType()`, `SNESLineSearchSetOrder()`
368f6dfbefdSBarry Smith           `SNESLineSearchSetPostCheck()`, `SNESLineSearchSetPreCheck()` `SNESLineSearchSetComputeNorms()`, `SNESGetLineSearch()`
369ebe3b25bSBarry Smith M*/
370d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NEWTONLS(SNES snes)
371d71ae5a4SJacob Faibussowitsch {
37204d7464bSBarry Smith   SNES_NEWTONLS *neP;
373d8d34be6SBarry Smith   SNESLineSearch linesearch;
3745e42470aSBarry Smith 
3753a40ed3dSBarry Smith   PetscFunctionBegin;
37604d7464bSBarry Smith   snes->ops->setup          = SNESSetUp_NEWTONLS;
37704d7464bSBarry Smith   snes->ops->solve          = SNESSolve_NEWTONLS;
37804d7464bSBarry Smith   snes->ops->destroy        = SNESDestroy_NEWTONLS;
37904d7464bSBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_NEWTONLS;
38004d7464bSBarry Smith   snes->ops->view           = SNESView_NEWTONLS;
38104d7464bSBarry Smith   snes->ops->reset          = SNESReset_NEWTONLS;
3825e42470aSBarry Smith 
383efd4aadfSBarry Smith   snes->npcside = PC_RIGHT;
38442f4f86dSBarry Smith   snes->usesksp = PETSC_TRUE;
385efd4aadfSBarry Smith   snes->usesnpc = PETSC_TRUE;
3864fc747eaSLawrence Mitchell 
3879566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes, &linesearch));
38848a46eb9SPierre Jolivet   if (!((PetscObject)linesearch)->type_name) PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHBT));
389d8d34be6SBarry Smith 
3904fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
3914fc747eaSLawrence Mitchell 
3924dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&neP));
3935e42470aSBarry Smith   snes->data = (void *)neP;
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3955e42470aSBarry Smith }
396