xref: /petsc/src/snes/linesearch/impls/nleqerr/linesearchnleqerr.c (revision dcf2fd19addf68f961057cf0cf7f6217b6b499d5)
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 
22d4c6564cSPatrick Farrell #undef __FUNCT__
237cbffc34SPatrick Farrell #define __FUNCT__ "SNESLineSearchReset_NLEQERR"
247cbffc34SPatrick Farrell static PetscErrorCode SNESLineSearchReset_NLEQERR(SNESLineSearch linesearch)
257cbffc34SPatrick Farrell {
267cbffc34SPatrick Farrell   SNESLineSearch_NLEQERR *nleqerr;
277cbffc34SPatrick Farrell 
287cbffc34SPatrick Farrell   nleqerr   = (SNESLineSearch_NLEQERR*)linesearch->data;
297cbffc34SPatrick Farrell   nleqerr->mu_curr    = 0.0;
307cbffc34SPatrick Farrell   nleqerr->norm_delta_x_prev = -1.0;
317cbffc34SPatrick Farrell   nleqerr->norm_bar_delta_x_prev = -1.0;
327cbffc34SPatrick Farrell   PetscFunctionReturn(0);
337cbffc34SPatrick Farrell }
347cbffc34SPatrick Farrell 
357cbffc34SPatrick Farrell #undef __FUNCT__
36d4c6564cSPatrick Farrell #define __FUNCT__ "SNESLineSearchApply_NLEQERR"
37d4c6564cSPatrick Farrell static PetscErrorCode  SNESLineSearchApply_NLEQERR(SNESLineSearch linesearch)
38d4c6564cSPatrick Farrell {
39d4c6564cSPatrick Farrell   PetscBool         changed_y,changed_w;
40d4c6564cSPatrick Farrell   PetscErrorCode    ierr;
41d4c6564cSPatrick Farrell   Vec               X,F,Y,W,G;
42d4c6564cSPatrick Farrell   SNES              snes;
43d4c6564cSPatrick Farrell   PetscReal         fnorm, xnorm, ynorm, gnorm, wnorm;
44d4c6564cSPatrick Farrell   PetscReal         lambda, minlambda, stol;
45d4c6564cSPatrick Farrell   PetscViewer       monitor;
467cbffc34SPatrick Farrell   PetscInt          max_its, count, snes_iteration;
47d4c6564cSPatrick Farrell   PetscReal         theta, mudash, lambdadash;
48d4c6564cSPatrick Farrell   SNESLineSearch_NLEQERR *nleqerr;
49d4c6564cSPatrick Farrell   KSPConvergedReason kspreason;
50d4c6564cSPatrick Farrell 
51d4c6564cSPatrick Farrell   PetscFunctionBegin;
52d4c6564cSPatrick Farrell   ierr = PetscCitationsRegister(NLEQERR_citation, &NLEQERR_cited);CHKERRQ(ierr);
53d4c6564cSPatrick Farrell 
54d4c6564cSPatrick Farrell   ierr = SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G);CHKERRQ(ierr);
55d4c6564cSPatrick Farrell   ierr = SNESLineSearchGetNorms(linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
56d4c6564cSPatrick Farrell   ierr = SNESLineSearchGetLambda(linesearch, &lambda);CHKERRQ(ierr);
57d4c6564cSPatrick Farrell   ierr = SNESLineSearchGetSNES(linesearch, &snes);CHKERRQ(ierr);
58*dcf2fd19SBarry Smith   ierr = SNESLineSearchGetDefaultMonitor(linesearch, &monitor);CHKERRQ(ierr);
59d4c6564cSPatrick Farrell   ierr = SNESLineSearchGetTolerances(linesearch,&minlambda,NULL,NULL,NULL,NULL,&max_its);CHKERRQ(ierr);
60d4c6564cSPatrick Farrell   ierr = SNESGetTolerances(snes,NULL,NULL,&stol,NULL,NULL);CHKERRQ(ierr);
61d4c6564cSPatrick Farrell   nleqerr   = (SNESLineSearch_NLEQERR*)linesearch->data;
62d4c6564cSPatrick Farrell 
637cbffc34SPatrick Farrell   /* reset the state of the Lipschitz estimates */
647cbffc34SPatrick Farrell   ierr = SNESGetIterationNumber(snes, &snes_iteration);CHKERRQ(ierr);
657cbffc34SPatrick Farrell   if (snes_iteration == 0) {
667cbffc34SPatrick Farrell     ierr = SNESLineSearchReset_NLEQERR(linesearch);CHKERRQ(ierr);
677cbffc34SPatrick Farrell   }
687cbffc34SPatrick Farrell 
69d4c6564cSPatrick Farrell   /* precheck */
70d4c6564cSPatrick Farrell   ierr = SNESLineSearchPreCheck(linesearch,X,Y,&changed_y);CHKERRQ(ierr);
71422a814eSBarry Smith   ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
72d4c6564cSPatrick Farrell 
73d4c6564cSPatrick Farrell   ierr = VecNormBegin(Y, NORM_2, &ynorm);CHKERRQ(ierr);
74d4c6564cSPatrick Farrell   ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr);
75d4c6564cSPatrick Farrell   ierr = VecNormEnd(Y, NORM_2, &ynorm);CHKERRQ(ierr);
76d4c6564cSPatrick Farrell   ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr);
77d4c6564cSPatrick Farrell 
78d4c6564cSPatrick Farrell   /* Note: Y is *minus* the Newton step. For whatever reason PETSc doesn't solve with the minus on
79d4c6564cSPatrick Farrell      the RHS. */
80d4c6564cSPatrick Farrell 
81d4c6564cSPatrick Farrell   if (ynorm == 0.0) {
82d4c6564cSPatrick Farrell     if (monitor) {
83d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
84d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Initial direction and size is 0\n");CHKERRQ(ierr);
85d4c6564cSPatrick Farrell       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
86d4c6564cSPatrick Farrell     }
87d4c6564cSPatrick Farrell     ierr = VecCopy(X,W);CHKERRQ(ierr);
88d4c6564cSPatrick Farrell     ierr = VecCopy(F,G);CHKERRQ(ierr);
89d4c6564cSPatrick Farrell     ierr = SNESLineSearchSetNorms(linesearch,xnorm,fnorm,ynorm);CHKERRQ(ierr);
90e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
91d4c6564cSPatrick Farrell     PetscFunctionReturn(0);
92d4c6564cSPatrick Farrell   }
93d4c6564cSPatrick Farrell 
94d4c6564cSPatrick Farrell   /* At this point, we've solved the Newton system for delta_x, and we assume that
95d4c6564cSPatrick Farrell      its norm is greater than the solution tolerance (otherwise we wouldn't be in
96d4c6564cSPatrick Farrell      here). So let's go ahead and estimate the Lipschitz constant.
97d4c6564cSPatrick Farrell 
98d4c6564cSPatrick Farrell      W contains bar_delta_x_prev at this point. */
99d4c6564cSPatrick Farrell 
100d4c6564cSPatrick Farrell   if (monitor) {
101d4c6564cSPatrick Farrell     ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
102d4c6564cSPatrick Farrell     ierr = PetscViewerASCIIPrintf(monitor,"    Line search: norm of Newton step: %14.12e\n", (double) ynorm);CHKERRQ(ierr);
103d4c6564cSPatrick Farrell     ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
104d4c6564cSPatrick Farrell   }
105d4c6564cSPatrick Farrell 
106d4c6564cSPatrick Farrell   /* this needs information from a previous iteration, so can't do it on the first one */
107d4c6564cSPatrick Farrell   if (nleqerr->norm_delta_x_prev > 0 && nleqerr->norm_bar_delta_x_prev > 0) {
108d4c6564cSPatrick Farrell     ierr = VecWAXPY(G, +1.0, Y, W);CHKERRQ(ierr); /* bar_delta_x - delta_x; +1 because Y is -delta_x */
109d4c6564cSPatrick Farrell     ierr = VecNormBegin(G, NORM_2, &gnorm);CHKERRQ(ierr);
110d4c6564cSPatrick Farrell     ierr = VecNormEnd(G, NORM_2, &gnorm);CHKERRQ(ierr);
111d4c6564cSPatrick Farrell 
112d4c6564cSPatrick Farrell     nleqerr->mu_curr = nleqerr->lambda_prev * (nleqerr->norm_delta_x_prev * nleqerr->norm_bar_delta_x_prev) / (gnorm * ynorm);
113d4c6564cSPatrick Farrell     lambda = PetscMin(1.0, nleqerr->mu_curr);
114d4c6564cSPatrick Farrell 
115d4c6564cSPatrick Farrell     if (monitor) {
116d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
117d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: Lipschitz estimate: %14.12e; lambda: %14.12e\n", (double) nleqerr->mu_curr, (double) lambda);CHKERRQ(ierr);
118d4c6564cSPatrick Farrell       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
119d4c6564cSPatrick Farrell     }
120bf8ddb9bSPatrick Farrell   } else {
121d4c6564cSPatrick Farrell     lambda = linesearch->damping;
122d4c6564cSPatrick Farrell   }
123d4c6564cSPatrick Farrell 
124d4c6564cSPatrick Farrell   /* The main while loop of the algorithm.
125d4c6564cSPatrick Farrell      At the end of this while loop, G should have the accepted new X in it. */
126d4c6564cSPatrick Farrell 
127d4c6564cSPatrick Farrell   count = 0;
12859ea5d13SPatrick Farrell   while (PETSC_TRUE) {
129d4c6564cSPatrick Farrell     if (monitor) {
130d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
131d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: entering iteration with lambda: %14.12e\n", lambda);CHKERRQ(ierr);
132d4c6564cSPatrick Farrell       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
133d4c6564cSPatrick Farrell     }
134d4c6564cSPatrick Farrell 
135d4c6564cSPatrick Farrell     /* Check that we haven't performed too many iterations */
136d4c6564cSPatrick Farrell     count += 1;
137bf8ddb9bSPatrick Farrell     if (count >= max_its) {
138d4c6564cSPatrick Farrell       if (monitor) {
139d4c6564cSPatrick Farrell         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
140d4c6564cSPatrick Farrell         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: maximum iterations reached\n");CHKERRQ(ierr);
141d4c6564cSPatrick Farrell         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
142d4c6564cSPatrick Farrell       }
143e9b602ebSSatish Balay       ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT);CHKERRQ(ierr);
144d4c6564cSPatrick Farrell       PetscFunctionReturn(0);
145d4c6564cSPatrick Farrell     }
146d4c6564cSPatrick Farrell 
147d4c6564cSPatrick Farrell     /* Now comes the Regularity Test. */
148d4c6564cSPatrick Farrell     if (lambda <= minlambda) {
149d4c6564cSPatrick Farrell       /* This isn't what is suggested by Deuflhard, but it works better in my experience */
150d4c6564cSPatrick Farrell       if (monitor) {
151d4c6564cSPatrick Farrell         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
152d4c6564cSPatrick Farrell         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: lambda has reached lambdamin, taking full Newton step\n");CHKERRQ(ierr);
153d4c6564cSPatrick Farrell         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
154d4c6564cSPatrick Farrell       }
155d4c6564cSPatrick Farrell       lambda = 1.0;
156d4c6564cSPatrick Farrell       ierr = VecWAXPY(G, -lambda, Y, X);CHKERRQ(ierr);
1577cbffc34SPatrick Farrell 
1587cbffc34SPatrick Farrell       /* and clean up the state for next time */
1597cbffc34SPatrick Farrell       ierr = SNESLineSearchReset_NLEQERR(linesearch);CHKERRQ(ierr);
160d4c6564cSPatrick Farrell       break;
161d4c6564cSPatrick Farrell     }
162d4c6564cSPatrick Farrell 
163d4c6564cSPatrick Farrell     /* Compute new trial iterate */
164d4c6564cSPatrick Farrell     ierr = VecWAXPY(W, -lambda, Y, X);CHKERRQ(ierr);
165d4c6564cSPatrick Farrell     ierr = SNESComputeFunction(snes, W, G);CHKERRQ(ierr);
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 */
168d4c6564cSPatrick Farrell     ierr = KSPSolve(snes->ksp, G, W);CHKERRQ(ierr);
169d4c6564cSPatrick Farrell     ierr = KSPGetConvergedReason(snes->ksp, &kspreason);CHKERRQ(ierr);
170d4c6564cSPatrick Farrell     if (kspreason < 0) {
171d4c6564cSPatrick Farrell       ierr = PetscInfo(snes,"Solution for \\bar{delta x}^{k+1} failed.");CHKERRQ(ierr);
172d4c6564cSPatrick Farrell     }
173d4c6564cSPatrick Farrell 
174d4c6564cSPatrick Farrell     /* W now contains -bar_delta_x_curr. */
175d4c6564cSPatrick Farrell 
176d4c6564cSPatrick Farrell     ierr = VecNorm(W, NORM_2, &wnorm);CHKERRQ(ierr);
177d4c6564cSPatrick Farrell     if (monitor) {
178d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
179d4c6564cSPatrick Farrell       ierr = PetscViewerASCIIPrintf(monitor,"    Line search: norm of simplified Newton update: %14.12e\n", (double) wnorm);CHKERRQ(ierr);
180d4c6564cSPatrick Farrell       ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
181d4c6564cSPatrick Farrell     }
182d4c6564cSPatrick Farrell 
183d4c6564cSPatrick Farrell     /* compute the monitoring quantities theta and mudash. */
184d4c6564cSPatrick Farrell 
185d4c6564cSPatrick Farrell     theta = wnorm / ynorm;
186d4c6564cSPatrick Farrell 
187d4c6564cSPatrick Farrell     ierr = VecWAXPY(G, -(1.0 - lambda), Y, W);CHKERRQ(ierr);
188d4c6564cSPatrick Farrell     ierr = VecNorm(G, NORM_2, &gnorm);CHKERRQ(ierr);
189d4c6564cSPatrick Farrell 
190d4c6564cSPatrick Farrell     mudash = (0.5 * ynorm * lambda * lambda) / gnorm;
191d4c6564cSPatrick Farrell 
192d4c6564cSPatrick Farrell     /* Check for termination of the linesearch */
193d4c6564cSPatrick Farrell     if (theta >= 1.0) {
194d4c6564cSPatrick Farrell       /* need to go around again with smaller lambda */
195d4c6564cSPatrick Farrell       if (monitor) {
196d4c6564cSPatrick Farrell         ierr = PetscViewerASCIIAddTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
197d4c6564cSPatrick Farrell         ierr = PetscViewerASCIIPrintf(monitor,"    Line search: monotonicity check failed, ratio: %14.12e\n", (double) theta);CHKERRQ(ierr);
198d4c6564cSPatrick Farrell         ierr = PetscViewerASCIISubtractTab(monitor,((PetscObject)linesearch)->tablevel);CHKERRQ(ierr);
199d4c6564cSPatrick Farrell       }
200d4c6564cSPatrick Farrell       lambda = PetscMin(mudash, 0.5 * lambda);
201d4c6564cSPatrick Farrell       lambda = PetscMax(lambda, minlambda);
202d4c6564cSPatrick Farrell       /* continue through the loop, i.e. go back to regularity test */
203bf8ddb9bSPatrick Farrell     } else {
204d4c6564cSPatrick Farrell       /* linesearch terminated */
205d4c6564cSPatrick Farrell       lambdadash = PetscMin(1.0, mudash);
206d4c6564cSPatrick Farrell 
207d4c6564cSPatrick Farrell       if (lambdadash == 1.0 && lambda == 1.0 && wnorm <= stol) {
208d4c6564cSPatrick Farrell         /* store the updated state, X - Y - W, in G:
209d4c6564cSPatrick Farrell            I need to keep W for the next linesearch */
210bf8ddb9bSPatrick Farrell         ierr = VecCopy(X, G);CHKERRQ(ierr);
211d4c6564cSPatrick Farrell         ierr = VecAXPY(G, -1.0, Y);CHKERRQ(ierr);
212d4c6564cSPatrick Farrell         ierr = VecAXPY(G, -1.0, W);CHKERRQ(ierr);
213d4c6564cSPatrick Farrell         break;
214d4c6564cSPatrick Farrell       }
215d4c6564cSPatrick Farrell 
216d4c6564cSPatrick Farrell       /* Deuflhard suggests to add the following:
217d4c6564cSPatrick Farrell       else if (lambdadash >= 4.0 * lambda) {
218d4c6564cSPatrick Farrell         lambda = lambdadash;
219d4c6564cSPatrick Farrell       }
220d4c6564cSPatrick Farrell       to continue through the loop, i.e. go back to regularity test.
221d4c6564cSPatrick Farrell       I deliberately exclude this, as I have practical experience of this
222d4c6564cSPatrick Farrell       getting stuck in infinite loops (on e.g. an Allen--Cahn problem). */
223d4c6564cSPatrick Farrell 
224d4c6564cSPatrick Farrell       else {
225d4c6564cSPatrick Farrell         /* accept iterate without adding on, i.e. don't use bar_delta_x;
226d4c6564cSPatrick Farrell            again, I need to keep W for the next linesearch */
227d4c6564cSPatrick Farrell         ierr = VecWAXPY(G, -lambda, Y, X);CHKERRQ(ierr);
228d4c6564cSPatrick Farrell         break;
229d4c6564cSPatrick Farrell       }
230d4c6564cSPatrick Farrell     }
231d4c6564cSPatrick Farrell   }
232d4c6564cSPatrick Farrell 
233d4c6564cSPatrick Farrell   if (linesearch->ops->viproject) {
234d4c6564cSPatrick Farrell     ierr = (*linesearch->ops->viproject)(snes, G);CHKERRQ(ierr);
235d4c6564cSPatrick Farrell   }
236d4c6564cSPatrick Farrell 
237d4c6564cSPatrick Farrell   /* W currently contains -bar_delta_u. Scale it so that it contains bar_delta_u. */
238d4c6564cSPatrick Farrell   ierr = VecScale(W, -1.0);CHKERRQ(ierr);
239d4c6564cSPatrick Farrell 
240d4c6564cSPatrick Farrell   /* postcheck */
241d4c6564cSPatrick Farrell   ierr = SNESLineSearchPostCheck(linesearch,X,Y,G,&changed_y,&changed_w);CHKERRQ(ierr);
242d4c6564cSPatrick Farrell   if (changed_y || changed_w) {
243e9b602ebSSatish Balay     ierr = SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_USER);CHKERRQ(ierr);
244d4c6564cSPatrick Farrell     ierr = PetscInfo(snes,"Changing the search direction here doesn't make sense.\n");CHKERRQ(ierr);
245d4c6564cSPatrick Farrell     PetscFunctionReturn(0);
246d4c6564cSPatrick Farrell   }
247d4c6564cSPatrick Farrell 
248d4c6564cSPatrick Farrell   /* copy the solution and information from this iteration over */
249d4c6564cSPatrick Farrell   nleqerr->norm_delta_x_prev = ynorm;
250d4c6564cSPatrick Farrell   nleqerr->norm_bar_delta_x_prev = wnorm;
251d4c6564cSPatrick Farrell   nleqerr->lambda_prev = lambda;
252d4c6564cSPatrick Farrell 
253d4c6564cSPatrick Farrell   ierr = VecCopy(G, X);CHKERRQ(ierr);
254d4c6564cSPatrick Farrell   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
255d4c6564cSPatrick Farrell   ierr = VecNorm(X, NORM_2, &xnorm);CHKERRQ(ierr);
256d4c6564cSPatrick Farrell   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
257d4c6564cSPatrick Farrell   ierr = SNESLineSearchSetLambda(linesearch, lambda);CHKERRQ(ierr);
258d4c6564cSPatrick Farrell   ierr = SNESLineSearchSetNorms(linesearch, xnorm, fnorm, ynorm);CHKERRQ(ierr);
259d4c6564cSPatrick Farrell   PetscFunctionReturn(0);
260d4c6564cSPatrick Farrell }
261d4c6564cSPatrick Farrell 
262d4c6564cSPatrick Farrell #undef __FUNCT__
263d4c6564cSPatrick Farrell #define __FUNCT__ "SNESLineSearchView_NLEQERR"
264d4c6564cSPatrick Farrell PetscErrorCode SNESLineSearchView_NLEQERR(SNESLineSearch linesearch, PetscViewer viewer)
265d4c6564cSPatrick Farrell {
266d4c6564cSPatrick Farrell   PetscErrorCode          ierr;
267d4c6564cSPatrick Farrell   PetscBool               iascii;
268d4c6564cSPatrick Farrell   SNESLineSearch_NLEQERR *nleqerr;
269d4c6564cSPatrick Farrell 
270d4c6564cSPatrick Farrell   PetscFunctionBegin;
271d4c6564cSPatrick Farrell   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
272d4c6564cSPatrick Farrell   nleqerr   = (SNESLineSearch_NLEQERR*)linesearch->data;
273d4c6564cSPatrick Farrell   if (iascii) {
274d4c6564cSPatrick Farrell     ierr = PetscViewerASCIIPrintf(viewer, "  NLEQ-ERR affine-covariant linesearch");CHKERRQ(ierr);
275d4c6564cSPatrick Farrell     ierr = PetscViewerASCIIPrintf(viewer, "  current local Lipschitz estimate omega=%e\n", (double)nleqerr->mu_curr);CHKERRQ(ierr);
276d4c6564cSPatrick Farrell   }
277d4c6564cSPatrick Farrell   PetscFunctionReturn(0);
278d4c6564cSPatrick Farrell }
279d4c6564cSPatrick Farrell 
280d4c6564cSPatrick Farrell #undef __FUNCT__
281d4c6564cSPatrick Farrell #define __FUNCT__ "SNESLineSearchDestroy_NLEQERR"
282d4c6564cSPatrick Farrell static PetscErrorCode SNESLineSearchDestroy_NLEQERR(SNESLineSearch linesearch)
283d4c6564cSPatrick Farrell {
284d4c6564cSPatrick Farrell   PetscErrorCode ierr;
285d4c6564cSPatrick Farrell 
286d4c6564cSPatrick Farrell   PetscFunctionBegin;
287d4c6564cSPatrick Farrell   ierr = PetscFree(linesearch->data);CHKERRQ(ierr);
288d4c6564cSPatrick Farrell   PetscFunctionReturn(0);
289d4c6564cSPatrick Farrell }
290d4c6564cSPatrick Farrell 
291d4c6564cSPatrick Farrell #undef __FUNCT__
292d4c6564cSPatrick Farrell #define __FUNCT__ "SNESLineSearchCreate_NLEQERR"
293d4c6564cSPatrick Farrell /*MC
294d4c6564cSPatrick Farrell    SNESLINESEARCHNLEQERR - Error-oriented affine-covariant globalised Newton algorithm of Deuflhard (2011).
295d4c6564cSPatrick Farrell 
296d4c6564cSPatrick Farrell    This linesearch is intended for Newton-type methods which are affine covariant. Affine covariance
297d4c6564cSPatrick Farrell    means that Newton's method will give the same iterations for F(x) = 0 and AF(x) = 0 for a nonsingular
298d4c6564cSPatrick Farrell    matrix A. This is a fundamental property; the philosophy of this linesearch is that globalisations
299d4c6564cSPatrick Farrell    of Newton's method should carefully preserve it.
300d4c6564cSPatrick Farrell 
301d4c6564cSPatrick Farrell    For a discussion of the theory behind this algorithm, see
302d4c6564cSPatrick Farrell 
303d4c6564cSPatrick Farrell    @book{deuflhard2011,
304d4c6564cSPatrick Farrell      title={Newton Methods for Nonlinear Problems},
305d4c6564cSPatrick Farrell      author={Deuflhard, P.},
306d4c6564cSPatrick Farrell      volume={35},
307d4c6564cSPatrick Farrell      year={2011},
308d4c6564cSPatrick Farrell      publisher={Springer-Verlag},
309d4c6564cSPatrick Farrell      address={Berlin, Heidelberg}
310d4c6564cSPatrick Farrell    }
311d4c6564cSPatrick Farrell 
312d4c6564cSPatrick Farrell    Pseudocode is given on page 148.
313d4c6564cSPatrick Farrell 
314d4c6564cSPatrick Farrell    Options Database Keys:
315d4c6564cSPatrick Farrell +  -snes_linesearch_damping<1.0> - initial step length
316d4c6564cSPatrick Farrell -  -snes_linesearch_minlambda<1e-12> - minimum step length allowed
317d4c6564cSPatrick Farrell 
318d4c6564cSPatrick Farrell    Contributed by Patrick Farrell <patrick.farrell@maths.ox.ac.uk>
319d4c6564cSPatrick Farrell 
320d4c6564cSPatrick Farrell    Level: advanced
321d4c6564cSPatrick Farrell 
322d4c6564cSPatrick Farrell .keywords: SNES, SNESLineSearch, damping
323d4c6564cSPatrick Farrell 
324d4c6564cSPatrick Farrell .seealso: SNESLineSearchCreate(), SNESLineSearchSetType()
325d4c6564cSPatrick Farrell M*/
326d4c6564cSPatrick Farrell PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_NLEQERR(SNESLineSearch linesearch)
327d4c6564cSPatrick Farrell {
328d4c6564cSPatrick Farrell 
329d4c6564cSPatrick Farrell   SNESLineSearch_NLEQERR *nleqerr;
330d4c6564cSPatrick Farrell   PetscErrorCode    ierr;
331d4c6564cSPatrick Farrell 
332d4c6564cSPatrick Farrell   PetscFunctionBegin;
333d4c6564cSPatrick Farrell   linesearch->ops->apply          = SNESLineSearchApply_NLEQERR;
334d4c6564cSPatrick Farrell   linesearch->ops->destroy        = SNESLineSearchDestroy_NLEQERR;
335d4c6564cSPatrick Farrell   linesearch->ops->setfromoptions = NULL;
336d4c6564cSPatrick Farrell   linesearch->ops->reset          = SNESLineSearchReset_NLEQERR;
337d4c6564cSPatrick Farrell   linesearch->ops->view           = SNESLineSearchView_NLEQERR;
338d4c6564cSPatrick Farrell   linesearch->ops->setup          = NULL;
339d4c6564cSPatrick Farrell 
340d4c6564cSPatrick Farrell   ierr = PetscNewLog(linesearch,&nleqerr);CHKERRQ(ierr);
341d4c6564cSPatrick Farrell 
342d4c6564cSPatrick Farrell   linesearch->data    = (void*)nleqerr;
343d4c6564cSPatrick Farrell   linesearch->max_its = 40;
3447cbffc34SPatrick Farrell   SNESLineSearchReset_NLEQERR(linesearch);
345d4c6564cSPatrick Farrell   PetscFunctionReturn(0);
346d4c6564cSPatrick Farrell }
347