xref: /petsc/src/snes/linesearch/impls/nleqerr/linesearchnleqerr.c (revision bd89dbf26d8a5efecb980364933175da61864cd7)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h> /*I  "petscsnes.h"  I*/
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
3d4c6564cSPatrick Farrell 
4d4c6564cSPatrick Farrell typedef struct {
5d4c6564cSPatrick Farrell   PetscReal norm_delta_x_prev;     /* norm of previous update */
6d4c6564cSPatrick Farrell   PetscReal norm_bar_delta_x_prev; /* norm of previous bar update */
7d4c6564cSPatrick Farrell   PetscReal mu_curr;               /* current local Lipschitz estimate */
8d4c6564cSPatrick Farrell   PetscReal lambda_prev;           /* previous step length: for some reason SNESLineSearchGetLambda returns 1 instead of the previous step length */
9d4c6564cSPatrick Farrell } SNESLineSearch_NLEQERR;
10d4c6564cSPatrick Farrell 
1157a6bf86SPatrick Farrell static PetscBool  NLEQERR_cited      = PETSC_FALSE;
1257a6bf86SPatrick Farrell static const char NLEQERR_citation[] = "@book{deuflhard2011,\n"
13d4c6564cSPatrick Farrell                                        "  title = {Newton Methods for Nonlinear Problems},\n"
14d4c6564cSPatrick Farrell                                        "  author = {Peter Deuflhard},\n"
15d4c6564cSPatrick Farrell                                        "  volume = 35,\n"
16d4c6564cSPatrick Farrell                                        "  year = 2011,\n"
17d4c6564cSPatrick Farrell                                        "  isbn = {978-3-642-23898-7},\n"
18d4c6564cSPatrick Farrell                                        "  doi  = {10.1007/978-3-642-23899-4},\n"
19d4c6564cSPatrick Farrell                                        "  publisher = {Springer-Verlag},\n"
20d4c6564cSPatrick Farrell                                        "  address = {Berlin, Heidelberg}\n}\n";
21d4c6564cSPatrick Farrell 
SNESLineSearchReset_NLEQERR(SNESLineSearch linesearch)22d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchReset_NLEQERR(SNESLineSearch linesearch)
23d71ae5a4SJacob Faibussowitsch {
2470d8d27cSBarry Smith   SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR *)linesearch->data;
257cbffc34SPatrick Farrell 
2670d8d27cSBarry Smith   PetscFunctionBegin;
277cbffc34SPatrick Farrell   nleqerr->mu_curr               = 0.0;
287cbffc34SPatrick Farrell   nleqerr->norm_delta_x_prev     = -1.0;
297cbffc34SPatrick Farrell   nleqerr->norm_bar_delta_x_prev = -1.0;
303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
317cbffc34SPatrick Farrell }
327cbffc34SPatrick Farrell 
SNESLineSearchApply_NLEQERR(SNESLineSearch linesearch)33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_NLEQERR(SNESLineSearch linesearch)
34d71ae5a4SJacob Faibussowitsch {
35d4c6564cSPatrick Farrell   PetscBool               changed_y, changed_w;
36d4c6564cSPatrick Farrell   Vec                     X, F, Y, W, G;
37d4c6564cSPatrick Farrell   SNES                    snes;
38d4c6564cSPatrick Farrell   PetscReal               fnorm, xnorm, ynorm, gnorm, wnorm;
39d4c6564cSPatrick Farrell   PetscReal               lambda, minlambda, stol;
40d4c6564cSPatrick Farrell   PetscViewer             monitor;
41a99ef635SJonas Heinzmann   PetscInt                max_it, count, snes_iteration;
42d4c6564cSPatrick Farrell   PetscReal               theta, mudash, lambdadash;
4370d8d27cSBarry Smith   SNESLineSearch_NLEQERR *nleqerr = (SNESLineSearch_NLEQERR *)linesearch->data;
44d4c6564cSPatrick Farrell   KSPConvergedReason      kspreason;
45d4c6564cSPatrick Farrell 
46d4c6564cSPatrick Farrell   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(NLEQERR_citation, &NLEQERR_cited));
48d4c6564cSPatrick Farrell 
499566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G));
509566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm));
519566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
529566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
539566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
54a99ef635SJonas Heinzmann   PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, NULL, NULL, NULL, NULL, &max_it));
559566063dSJacob Faibussowitsch   PetscCall(SNESGetTolerances(snes, NULL, NULL, &stol, NULL, NULL));
56d4c6564cSPatrick Farrell 
577cbffc34SPatrick Farrell   /* reset the state of the Lipschitz estimates */
589566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &snes_iteration));
5948a46eb9SPierre Jolivet   if (!snes_iteration) PetscCall(SNESLineSearchReset_NLEQERR(linesearch));
607cbffc34SPatrick Farrell 
61d4c6564cSPatrick Farrell   /* precheck */
629566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
639566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
64d4c6564cSPatrick Farrell 
659566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(Y, NORM_2, &ynorm));
669566063dSJacob Faibussowitsch   PetscCall(VecNormBegin(X, NORM_2, &xnorm));
679566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(Y, NORM_2, &ynorm));
689566063dSJacob Faibussowitsch   PetscCall(VecNormEnd(X, NORM_2, &xnorm));
69d4c6564cSPatrick Farrell 
7070d8d27cSBarry Smith   /* Note: Y is *minus* the Newton step. For whatever reason PETSc doesn't solve with the minus on  the RHS. */
71d4c6564cSPatrick Farrell 
72d4c6564cSPatrick Farrell   if (ynorm == 0.0) {
73d4c6564cSPatrick Farrell     if (monitor) {
749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Initial direction and size is 0\n"));
769566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
77d4c6564cSPatrick Farrell     }
789566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, W));
799566063dSJacob Faibussowitsch     PetscCall(VecCopy(F, G));
809566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm));
819566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
823ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
83d4c6564cSPatrick Farrell   }
84d4c6564cSPatrick Farrell 
85d4c6564cSPatrick Farrell   /* At this point, we've solved the Newton system for delta_x, and we assume that
86d4c6564cSPatrick Farrell      its norm is greater than the solution tolerance (otherwise we wouldn't be in
87d4c6564cSPatrick Farrell      here). So let's go ahead and estimate the Lipschitz constant.
88d4c6564cSPatrick Farrell 
89d4c6564cSPatrick Farrell      W contains bar_delta_x_prev at this point. */
90d4c6564cSPatrick Farrell 
91d4c6564cSPatrick Farrell   if (monitor) {
929566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
939566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: norm of Newton step: %14.12e\n", (double)ynorm));
949566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
95d4c6564cSPatrick Farrell   }
96d4c6564cSPatrick Farrell 
97d4c6564cSPatrick Farrell   /* this needs information from a previous iteration, so can't do it on the first one */
98d4c6564cSPatrick Farrell   if (nleqerr->norm_delta_x_prev > 0 && nleqerr->norm_bar_delta_x_prev > 0) {
999566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(G, +1.0, Y, W)); /* bar_delta_x - delta_x; +1 because Y is -delta_x */
1009566063dSJacob Faibussowitsch     PetscCall(VecNormBegin(G, NORM_2, &gnorm));
1019566063dSJacob Faibussowitsch     PetscCall(VecNormEnd(G, NORM_2, &gnorm));
102d4c6564cSPatrick Farrell 
103d4c6564cSPatrick Farrell     nleqerr->mu_curr = nleqerr->lambda_prev * (nleqerr->norm_delta_x_prev * nleqerr->norm_bar_delta_x_prev) / (gnorm * ynorm);
104d4c6564cSPatrick Farrell     lambda           = PetscMin(1.0, nleqerr->mu_curr);
105d4c6564cSPatrick Farrell 
106d4c6564cSPatrick Farrell     if (monitor) {
1079566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1089566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: Lipschitz estimate: %14.12e; lambda: %14.12e\n", (double)nleqerr->mu_curr, (double)lambda));
1099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
110d4c6564cSPatrick Farrell     }
111bf8ddb9bSPatrick Farrell   } else {
112d4c6564cSPatrick Farrell     lambda = linesearch->damping;
113d4c6564cSPatrick Farrell   }
114d4c6564cSPatrick Farrell 
115d4c6564cSPatrick Farrell   /* The main while loop of the algorithm.
116d4c6564cSPatrick Farrell      At the end of this while loop, G should have the accepted new X in it. */
117d4c6564cSPatrick Farrell 
118d4c6564cSPatrick Farrell   count = 0;
11959ea5d13SPatrick Farrell   while (PETSC_TRUE) {
120d4c6564cSPatrick Farrell     if (monitor) {
1219566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
12263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: entering iteration with lambda: %14.12e\n", (double)lambda));
1239566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
124d4c6564cSPatrick Farrell     }
125d4c6564cSPatrick Farrell 
126d4c6564cSPatrick Farrell     /* Check that we haven't performed too many iterations */
127d4c6564cSPatrick Farrell     count += 1;
128a99ef635SJonas Heinzmann     if (count >= max_it) {
129d4c6564cSPatrick Farrell       if (monitor) {
1309566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1319566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: maximum iterations reached\n"));
1329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
133d4c6564cSPatrick Farrell       }
1349566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
1353ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
136d4c6564cSPatrick Farrell     }
137d4c6564cSPatrick Farrell 
138d4c6564cSPatrick Farrell     /* Now comes the Regularity Test. */
139d4c6564cSPatrick Farrell     if (lambda <= minlambda) {
140d4c6564cSPatrick Farrell       /* This isn't what is suggested by Deuflhard, but it works better in my experience */
141d4c6564cSPatrick Farrell       if (monitor) {
1429566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
143a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: lambda has reached minlambda, taking full Newton step\n"));
1449566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
145d4c6564cSPatrick Farrell       }
146d4c6564cSPatrick Farrell       lambda = 1.0;
1479566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(G, -lambda, Y, X));
1487cbffc34SPatrick Farrell 
1497cbffc34SPatrick Farrell       /* and clean up the state for next time */
1509566063dSJacob Faibussowitsch       PetscCall(SNESLineSearchReset_NLEQERR(linesearch));
15170d8d27cSBarry Smith       /*
15270d8d27cSBarry Smith          The clang static analyzer detected a problem here; once the loop is broken the values
15370d8d27cSBarry Smith          nleqerr->norm_delta_x_prev     = ynorm;
15470d8d27cSBarry Smith          nleqerr->norm_bar_delta_x_prev = wnorm;
15570d8d27cSBarry Smith          are set, but wnorm has not even been computed.
15670d8d27cSBarry Smith          I don't know if this is the correct fix but by setting ynorm and wnorm to -1.0 at
15770d8d27cSBarry Smith          least the linesearch object is kept in the state set by the SNESLineSearchReset_NLEQERR() call above
15870d8d27cSBarry Smith       */
15970d8d27cSBarry Smith       ynorm = wnorm = -1.0;
160d4c6564cSPatrick Farrell       break;
161d4c6564cSPatrick Farrell     }
162d4c6564cSPatrick Farrell 
163d4c6564cSPatrick Farrell     /* Compute new trial iterate */
1649566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(W, -lambda, Y, X));
1659566063dSJacob Faibussowitsch     PetscCall(SNESComputeFunction(snes, W, G));
166d4c6564cSPatrick Farrell 
167d4c6564cSPatrick Farrell     /* Solve linear system for bar_delta_x_curr: old Jacobian, new RHS. Note absence of minus sign, compared to Deuflhard, in keeping with PETSc convention */
1689566063dSJacob Faibussowitsch     PetscCall(KSPSolve(snes->ksp, G, W));
1699566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(snes->ksp, &kspreason));
1709d3446b2SPierre Jolivet     if (kspreason < 0) PetscCall(PetscInfo(snes, "Solution for \\bar{delta x}^{k+1} failed.\n"));
171d4c6564cSPatrick Farrell 
172d4c6564cSPatrick Farrell     /* W now contains -bar_delta_x_curr. */
173d4c6564cSPatrick Farrell 
1749566063dSJacob Faibussowitsch     PetscCall(VecNorm(W, NORM_2, &wnorm));
175d4c6564cSPatrick Farrell     if (monitor) {
1769566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1779566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: norm of simplified Newton update: %14.12e\n", (double)wnorm));
1789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
179d4c6564cSPatrick Farrell     }
180d4c6564cSPatrick Farrell 
181d4c6564cSPatrick Farrell     /* compute the monitoring quantities theta and mudash. */
182d4c6564cSPatrick Farrell 
183d4c6564cSPatrick Farrell     theta = wnorm / ynorm;
184d4c6564cSPatrick Farrell 
1859566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(G, -(1.0 - lambda), Y, W));
1869566063dSJacob Faibussowitsch     PetscCall(VecNorm(G, NORM_2, &gnorm));
187d4c6564cSPatrick Farrell 
188d4c6564cSPatrick Farrell     mudash = (0.5 * ynorm * lambda * lambda) / gnorm;
189d4c6564cSPatrick Farrell 
190d4c6564cSPatrick Farrell     /* Check for termination of the linesearch */
191d4c6564cSPatrick Farrell     if (theta >= 1.0) {
192d4c6564cSPatrick Farrell       /* need to go around again with smaller lambda */
193d4c6564cSPatrick Farrell       if (monitor) {
1949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: monotonicity check failed, ratio: %14.12e\n", (double)theta));
1969566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
197d4c6564cSPatrick Farrell       }
198d4c6564cSPatrick Farrell       lambda = PetscMin(mudash, 0.5 * lambda);
199d4c6564cSPatrick Farrell       lambda = PetscMax(lambda, minlambda);
200d4c6564cSPatrick Farrell       /* continue through the loop, i.e. go back to regularity test */
201bf8ddb9bSPatrick Farrell     } else {
202d4c6564cSPatrick Farrell       /* linesearch terminated */
203d4c6564cSPatrick Farrell       lambdadash = PetscMin(1.0, mudash);
204d4c6564cSPatrick Farrell 
205d4c6564cSPatrick Farrell       if (lambdadash == 1.0 && lambda == 1.0 && wnorm <= stol) {
206d4c6564cSPatrick Farrell         /* store the updated state, X - Y - W, in G:
207d4c6564cSPatrick Farrell            I need to keep W for the next linesearch */
208ef46b1a6SStefano Zampini         PetscCall(VecWAXPY(G, -1.0, Y, X));
2099566063dSJacob Faibussowitsch         PetscCall(VecAXPY(G, -1.0, W));
210d4c6564cSPatrick Farrell         break;
211d4c6564cSPatrick Farrell       }
212d4c6564cSPatrick Farrell 
213d4c6564cSPatrick Farrell       /* Deuflhard suggests to add the following:
214d4c6564cSPatrick Farrell       else if (lambdadash >= 4.0 * lambda) {
215d4c6564cSPatrick Farrell         lambda = lambdadash;
216d4c6564cSPatrick Farrell       }
217d4c6564cSPatrick Farrell       to continue through the loop, i.e. go back to regularity test.
218d4c6564cSPatrick Farrell       I deliberately exclude this, as I have practical experience of this
219d4c6564cSPatrick Farrell       getting stuck in infinite loops (on e.g. an Allen--Cahn problem). */
220d4c6564cSPatrick Farrell 
221d4c6564cSPatrick Farrell       else {
222d4c6564cSPatrick Farrell         /* accept iterate without adding on, i.e. don't use bar_delta_x;
223d4c6564cSPatrick Farrell            again, I need to keep W for the next linesearch */
2249566063dSJacob Faibussowitsch         PetscCall(VecWAXPY(G, -lambda, Y, X));
225d4c6564cSPatrick Farrell         break;
226d4c6564cSPatrick Farrell       }
227d4c6564cSPatrick Farrell     }
228d4c6564cSPatrick Farrell   }
229d4c6564cSPatrick Farrell 
2301baa6e33SBarry Smith   if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, G));
231d4c6564cSPatrick Farrell 
232d4c6564cSPatrick Farrell   /* W currently contains -bar_delta_u. Scale it so that it contains bar_delta_u. */
2339566063dSJacob Faibussowitsch   PetscCall(VecScale(W, -1.0));
234d4c6564cSPatrick Farrell 
235d4c6564cSPatrick Farrell   /* postcheck */
2369566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, G, &changed_y, &changed_w));
237d4c6564cSPatrick Farrell   if (changed_y || changed_w) {
2389566063dSJacob Faibussowitsch     PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_USER));
2399566063dSJacob Faibussowitsch     PetscCall(PetscInfo(snes, "Changing the search direction here doesn't make sense.\n"));
2403ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
241d4c6564cSPatrick Farrell   }
242d4c6564cSPatrick Farrell 
243d4c6564cSPatrick Farrell   /* copy the solution and information from this iteration over */
244d4c6564cSPatrick Farrell   nleqerr->norm_delta_x_prev     = ynorm;
245d4c6564cSPatrick Farrell   nleqerr->norm_bar_delta_x_prev = wnorm;
246d4c6564cSPatrick Farrell   nleqerr->lambda_prev           = lambda;
247d4c6564cSPatrick Farrell 
2489566063dSJacob Faibussowitsch   PetscCall(VecCopy(G, X));
2499566063dSJacob Faibussowitsch   PetscCall(SNESComputeFunction(snes, X, F));
2509566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &xnorm));
2519566063dSJacob Faibussowitsch   PetscCall(VecNorm(F, NORM_2, &fnorm));
252*76c63389SBarry Smith   SNESLineSearchCheckFunctionDomainError(snes, linesearch, fnorm);
2539566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
25457508eceSPierre Jolivet   PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm < 0 ? PETSC_INFINITY : ynorm));
2553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
256d4c6564cSPatrick Farrell }
257d4c6564cSPatrick Farrell 
SNESLineSearchView_NLEQERR(SNESLineSearch linesearch,PetscViewer viewer)25866976f2fSJacob Faibussowitsch static PetscErrorCode SNESLineSearchView_NLEQERR(SNESLineSearch linesearch, PetscViewer viewer)
259d71ae5a4SJacob Faibussowitsch {
2609f196a02SMartin Diehl   PetscBool               isascii;
261d4c6564cSPatrick Farrell   SNESLineSearch_NLEQERR *nleqerr;
262d4c6564cSPatrick Farrell 
263d4c6564cSPatrick Farrell   PetscFunctionBegin;
2649f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
265d4c6564cSPatrick Farrell   nleqerr = (SNESLineSearch_NLEQERR *)linesearch->data;
2669f196a02SMartin Diehl   if (isascii) {
2679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  NLEQ-ERR affine-covariant linesearch"));
2689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  current local Lipschitz estimate omega=%e\n", (double)nleqerr->mu_curr));
269d4c6564cSPatrick Farrell   }
2703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
271d4c6564cSPatrick Farrell }
272d4c6564cSPatrick Farrell 
SNESLineSearchDestroy_NLEQERR(SNESLineSearch linesearch)273d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchDestroy_NLEQERR(SNESLineSearch linesearch)
274d71ae5a4SJacob Faibussowitsch {
275d4c6564cSPatrick Farrell   PetscFunctionBegin;
2769566063dSJacob Faibussowitsch   PetscCall(PetscFree(linesearch->data));
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
278d4c6564cSPatrick Farrell }
279d4c6564cSPatrick Farrell 
280d4c6564cSPatrick Farrell /*MC
281420bcc1bSBarry Smith    SNESLINESEARCHNLEQERR - Error-oriented affine-covariant globalised Newton algorithm of Deuflhard {cite}`deuflhard2011`
282d4c6564cSPatrick Farrell 
283d4c6564cSPatrick Farrell    This linesearch is intended for Newton-type methods which are affine covariant. Affine covariance
284420bcc1bSBarry Smith    means that Newton's method will give the same iterations for F(x) = 0 and AF(x) = 0 for any nonsingular
285d4c6564cSPatrick Farrell    matrix A. This is a fundamental property; the philosophy of this linesearch is that globalisations
286d4c6564cSPatrick Farrell    of Newton's method should carefully preserve it.
287d4c6564cSPatrick Farrell 
288d4c6564cSPatrick Farrell    Options Database Keys:
289a99ef635SJonas Heinzmann +  -snes_linesearch_damping <1.0>      - initial `lambda`
290a99ef635SJonas Heinzmann .  -snes_linesearch_max_it <40>        - maximum number of iterations for the line search
291a99ef635SJonas Heinzmann -  -snes_linesearch_minlambda <1e\-12> - minimum `lambda` allowed
292d4c6564cSPatrick Farrell 
293d4c6564cSPatrick Farrell    Level: advanced
294d4c6564cSPatrick Farrell 
295f6dfbefdSBarry Smith    Note:
296f6dfbefdSBarry Smith    Contributed by Patrick Farrell <patrick.farrell@maths.ox.ac.uk>
297f6dfbefdSBarry Smith 
298420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNES`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
299d4c6564cSPatrick Farrell M*/
SNESLineSearchCreate_NLEQERR(SNESLineSearch linesearch)300d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NLEQERR(SNESLineSearch linesearch)
301d71ae5a4SJacob Faibussowitsch {
302d4c6564cSPatrick Farrell   SNESLineSearch_NLEQERR *nleqerr;
303d4c6564cSPatrick Farrell 
304d4c6564cSPatrick Farrell   PetscFunctionBegin;
305d4c6564cSPatrick Farrell   linesearch->ops->apply          = SNESLineSearchApply_NLEQERR;
306d4c6564cSPatrick Farrell   linesearch->ops->destroy        = SNESLineSearchDestroy_NLEQERR;
307d4c6564cSPatrick Farrell   linesearch->ops->setfromoptions = NULL;
308d4c6564cSPatrick Farrell   linesearch->ops->reset          = SNESLineSearchReset_NLEQERR;
309d4c6564cSPatrick Farrell   linesearch->ops->view           = SNESLineSearchView_NLEQERR;
310d4c6564cSPatrick Farrell   linesearch->ops->setup          = NULL;
311d4c6564cSPatrick Farrell 
3124dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&nleqerr));
313d4c6564cSPatrick Farrell 
314d4c6564cSPatrick Farrell   linesearch->data   = (void *)nleqerr;
315a99ef635SJonas Heinzmann   linesearch->max_it = 40;
3163ba16761SJacob Faibussowitsch   PetscCall(SNESLineSearchReset_NLEQERR(linesearch));
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318d4c6564cSPatrick Farrell }
319