xref: /petsc/src/snes/linesearch/impls/cp/linesearchcp.c (revision a99ef635ef3fb31080710bc8e0dbbde292199c98)
1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
26a388c36SPeter Brune #include <petscsnes.h>
3bbd5d0b3SPeter Brune 
4d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESLineSearchApply_CP(SNESLineSearch linesearch)
5d71ae5a4SJacob Faibussowitsch {
6422a814eSBarry Smith   PetscBool   changed_y, changed_w;
76a388c36SPeter Brune   Vec         X, Y, F, W;
86a388c36SPeter Brune   SNES        snes;
9*a99ef635SJonas Heinzmann   PetscReal   xnorm, ynorm, gnorm, minlambda, maxlambda, rtol, atol, ltol;
10bbd5d0b3SPeter Brune   PetscReal   lambda, lambda_old, lambda_update, delLambda;
11b8ac21a2SPeter Brune   PetscScalar fty, fty_init, fty_old, fty_mid1, fty_mid2, s;
12*a99ef635SJonas Heinzmann   PetscInt    i, max_it;
136a388c36SPeter Brune   PetscViewer monitor;
14bbd5d0b3SPeter Brune 
15bbd5d0b3SPeter Brune   PetscFunctionBegin;
169566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
179566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
189566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
199566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
20*a99ef635SJonas Heinzmann   PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, &maxlambda, &rtol, &atol, &ltol, &max_it));
219566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
229566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
236a388c36SPeter Brune 
24bbd5d0b3SPeter Brune   /* precheck */
259566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
26bbd5d0b3SPeter Brune   lambda_old = 0.0;
2762893cf2SPeter Brune 
28*a99ef635SJonas Heinzmann   /* evaluate initial point */
29d5def619SJonas Heinzmann   if (linesearch->ops->vidirderiv) {
30d5def619SJonas Heinzmann     PetscCall((*linesearch->ops->vidirderiv)(snes, F, X, Y, &fty_old));
31d5def619SJonas Heinzmann   } else {
329566063dSJacob Faibussowitsch     PetscCall(VecDot(F, Y, &fty_old));
33d5def619SJonas Heinzmann   }
3496b7d5c2SJed Brown   if (PetscAbsScalar(fty_old) < atol * ynorm) {
359e0a2ef9SBarry Smith     if (monitor) {
369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
3763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search terminated at initial point because dot(F,Y) = %g < atol*||y|| = %g\n", (double)PetscAbsScalar(fty_old), (double)(atol * ynorm)));
389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
399e0a2ef9SBarry Smith     }
409566063dSJacob Faibussowitsch     PetscCall(SNESSetConvergedReason(linesearch->snes, SNES_CONVERGED_FNORM_ABS));
413ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
429e0a2ef9SBarry Smith   }
4362893cf2SPeter Brune   fty_init = fty_old;
44c0b2db79SJed Brown 
45*a99ef635SJonas Heinzmann   for (i = 0; i < max_it; i++) {
46bbd5d0b3SPeter Brune     /* compute the norm at lambda */
47ef46b1a6SStefano Zampini     PetscCall(VecWAXPY(W, -lambda, Y, X));
481baa6e33SBarry Smith     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
499566063dSJacob Faibussowitsch     PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
50d5def619SJonas Heinzmann     if (linesearch->ops->vidirderiv) {
51d5def619SJonas Heinzmann       PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty));
52d5def619SJonas Heinzmann     } else {
539566063dSJacob Faibussowitsch       PetscCall(VecDot(F, Y, &fty));
54d5def619SJonas Heinzmann     }
55bbd5d0b3SPeter Brune 
56*a99ef635SJonas Heinzmann     /* compute change of lambda */
57bbd5d0b3SPeter Brune     delLambda = lambda - lambda_old;
58b8ac21a2SPeter Brune 
59*a99ef635SJonas Heinzmann     /* check change of lambda tolerance */
60*a99ef635SJonas Heinzmann     if (PetscAbsReal(delLambda) < ltol) {
61*a99ef635SJonas Heinzmann       if (monitor) {
62*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
63*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: abs(dlambda) = %g < ltol = %g\n", (double)PetscAbsReal(delLambda), (double)ltol));
64*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
65*a99ef635SJonas Heinzmann       }
66*a99ef635SJonas Heinzmann       break;
67*a99ef635SJonas Heinzmann     }
68*a99ef635SJonas Heinzmann 
69*a99ef635SJonas Heinzmann     /* check relative tolerance */
70*a99ef635SJonas Heinzmann     if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) {
71*a99ef635SJonas Heinzmann       if (monitor) {
72*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
73*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: abs(fty/fty_init) = %g <= rtol  = %g\n", (double)(PetscAbsScalar(fty) / PetscAbsScalar(fty_init)), (double)rtol));
74*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
75*a99ef635SJonas Heinzmann       }
76*a99ef635SJonas Heinzmann       break;
77*a99ef635SJonas Heinzmann     }
78*a99ef635SJonas Heinzmann 
79*a99ef635SJonas Heinzmann     /* check absolute tolerance */
80*a99ef635SJonas Heinzmann     if (PetscAbsScalar(fty) < atol * ynorm && i > 0) {
81*a99ef635SJonas Heinzmann       if (monitor) {
82*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
83*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: abs(fty)/||y|| = %g <= atol = %g\n", (double)(PetscAbsScalar(fty) / ynorm), (double)atol));
84*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
85*a99ef635SJonas Heinzmann       }
86*a99ef635SJonas Heinzmann       break;
87*a99ef635SJonas Heinzmann     }
88*a99ef635SJonas Heinzmann 
89*a99ef635SJonas Heinzmann     /* print iteration information */
906a388c36SPeter Brune     if (monitor) {
919566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
929566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: lambdas = [%g, %g], ftys = [%g, %g]\n", (double)lambda, (double)lambda_old, (double)PetscRealPart(fty), (double)PetscRealPart(fty_old)));
939566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
94bbd5d0b3SPeter Brune     }
95bbd5d0b3SPeter Brune 
96*a99ef635SJonas Heinzmann     /* approximate the second derivative */
97b000cd8dSPeter Brune     if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) {
98*a99ef635SJonas Heinzmann       /* first order finite difference approximation */
99b8ac21a2SPeter Brune       s = (fty - fty_old) / delLambda;
100b000cd8dSPeter Brune     } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
101*a99ef635SJonas Heinzmann       /* evaluate function at midpoint 0.5 * (lambda + lambda_old) */
102ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X));
1031baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
1049566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
105d5def619SJonas Heinzmann       if (linesearch->ops->vidirderiv) {
106d5def619SJonas Heinzmann         PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty_mid1));
107d5def619SJonas Heinzmann       } else {
1089566063dSJacob Faibussowitsch         PetscCall(VecDot(F, Y, &fty_mid1));
109d5def619SJonas Heinzmann       }
110*a99ef635SJonas Heinzmann 
111*a99ef635SJonas Heinzmann       /* second order finite difference approximation */
112b8ac21a2SPeter Brune       s = (3. * fty - 4. * fty_mid1 + fty_old) / delLambda;
113*a99ef635SJonas Heinzmann 
114b8ac21a2SPeter Brune     } else {
115*a99ef635SJonas Heinzmann       /* evaluate function at midpoint 0.5 * (lambda + lambda_old) */
116ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -0.5 * (lambda + lambda_old), Y, X));
1171baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
1189566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
119d5def619SJonas Heinzmann       if (linesearch->ops->vidirderiv) {
120d5def619SJonas Heinzmann         PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty_mid1));
121d5def619SJonas Heinzmann       } else {
1229566063dSJacob Faibussowitsch         PetscCall(VecDot(F, Y, &fty_mid1));
123d5def619SJonas Heinzmann       }
124*a99ef635SJonas Heinzmann 
125*a99ef635SJonas Heinzmann       /* evaluate function at lambda + 0.5 * (lambda - lambda_old) */
126ef46b1a6SStefano Zampini       PetscCall(VecWAXPY(W, -(lambda + 0.5 * (lambda - lambda_old)), Y, X));
1271baa6e33SBarry Smith       if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
1289566063dSJacob Faibussowitsch       PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
129d5def619SJonas Heinzmann       if (linesearch->ops->vidirderiv) {
130d5def619SJonas Heinzmann         PetscCall((*linesearch->ops->vidirderiv)(snes, F, W, Y, &fty_mid2));
131d5def619SJonas Heinzmann       } else {
1329566063dSJacob Faibussowitsch         PetscCall(VecDot(F, Y, &fty_mid2));
133d5def619SJonas Heinzmann       }
134*a99ef635SJonas Heinzmann 
135*a99ef635SJonas Heinzmann       /* third order finite difference approximation */
136b8ac21a2SPeter Brune       s = (2. * fty_mid2 + 3. * fty - 6. * fty_mid1 + fty_old) / (3. * delLambda);
137b8ac21a2SPeter Brune     }
138*a99ef635SJonas Heinzmann 
139*a99ef635SJonas Heinzmann     /* compute secant update (modifying the search direction if necessary) */
140b8ac21a2SPeter Brune     if (PetscRealPart(s) > 0.) s = -s;
141*a99ef635SJonas Heinzmann     if (s == 0.0) {
142*a99ef635SJonas Heinzmann       if (monitor) {
143*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
144*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: del2Fnrm = 0\n"));
145*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
146*a99ef635SJonas Heinzmann       }
147*a99ef635SJonas Heinzmann       break;
148*a99ef635SJonas Heinzmann     }
149b8ac21a2SPeter Brune     lambda_update = lambda - PetscRealPart(fty / s);
150b8ac21a2SPeter Brune 
151*a99ef635SJonas Heinzmann     /* if secant method would step out of bounds, exit with the respective bound */
152*a99ef635SJonas Heinzmann     if (lambda_update < minlambda) {
153*a99ef635SJonas Heinzmann       lambda_update = minlambda;
154*a99ef635SJonas Heinzmann       break;
155*a99ef635SJonas Heinzmann     }
156*a99ef635SJonas Heinzmann     if (lambda_update > maxlambda) {
157*a99ef635SJonas Heinzmann       lambda_update = maxlambda;
158*a99ef635SJonas Heinzmann       break;
159*a99ef635SJonas Heinzmann     }
160b8ac21a2SPeter Brune 
161*a99ef635SJonas Heinzmann     /* if lambda is NaN or Inf, do not accept update but exit with previous lambda */
162*a99ef635SJonas Heinzmann     if (PetscIsInfOrNanReal(lambda_update)) {
163*a99ef635SJonas Heinzmann       if (monitor) {
164*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
165*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search: lambda_update is NaN or Inf\n"));
166*a99ef635SJonas Heinzmann         PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
167*a99ef635SJonas Heinzmann       }
168*a99ef635SJonas Heinzmann       PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
169*a99ef635SJonas Heinzmann       break;
170*a99ef635SJonas Heinzmann     }
171bbd5d0b3SPeter Brune 
172bbd5d0b3SPeter Brune     /* compute the new state of the line search */
173bbd5d0b3SPeter Brune     lambda_old = lambda;
174bbd5d0b3SPeter Brune     lambda     = lambda_update;
175bbd5d0b3SPeter Brune     fty_old    = fty;
176bbd5d0b3SPeter Brune   }
177*a99ef635SJonas Heinzmann 
178bbd5d0b3SPeter Brune   /* construct the solution */
179ef46b1a6SStefano Zampini   PetscCall(VecWAXPY(W, -lambda, Y, X));
1801baa6e33SBarry Smith   if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
181bbd5d0b3SPeter Brune   /* postcheck */
1826b095a85SStefano Zampini   PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
1839566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
1846b095a85SStefano Zampini   if (changed_y) {
1856b095a85SStefano Zampini     if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
1866b095a85SStefano Zampini     if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
187bbd5d0b3SPeter Brune   }
1886b095a85SStefano Zampini   PetscCall(VecCopy(W, X));
1899566063dSJacob Faibussowitsch   PetscCall((*linesearch->ops->snesfunc)(snes, X, F));
190bbd5d0b3SPeter Brune 
1919566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchComputeNorms(linesearch));
1929566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
193bbd5d0b3SPeter Brune 
1946a388c36SPeter Brune   if (monitor) {
1959566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1969566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(monitor, "    Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm));
1979566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
198bbd5d0b3SPeter Brune   }
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200bbd5d0b3SPeter Brune }
201bbd5d0b3SPeter Brune 
202954494b2SJed Brown /*MC
2031a4f838cSPeter Brune    SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
204d5def619SJonas Heinzmann    artificial $G(x)$ for which the `SNESFunctionFn` $F(x) = grad G(x)$. Therefore, this line search seeks
205*a99ef635SJonas Heinzmann    to find roots of the directional derivative via a secant method, that is $F(x_k - \lambda Y_k) \cdot Y_k / ||Y_k|| = 0$.
206cd7522eaSPeter Brune 
207cd7522eaSPeter Brune    Options Database Keys:
208*a99ef635SJonas Heinzmann +  -snes_linesearch_minlambda <1e\-12> - the minimum acceptable `lambda` (scaling of solution update)
209*a99ef635SJonas Heinzmann .  -snes_linesearch_maxlambda <1.0>    - the algorithm ensures that `lambda` is never larger than this value
210*a99ef635SJonas Heinzmann .  -snes_linesearch_damping <1.0>      - initial `lambda` on entry to the line search
211*a99ef635SJonas Heinzmann .  -snes_linesearch_order <1>          - order of the approximation in the secant method, must be 1, 2, or 3
212*a99ef635SJonas Heinzmann .  -snes_linesearch_max_it <1>         - the maximum number of secant iterations performed
213*a99ef635SJonas Heinzmann .  -snes_linesearch_rtol <1e\-8>       - relative tolerance for the directional derivative
214*a99ef635SJonas Heinzmann .  -snes_linesearch_atol <1e\-15>      - absolute tolerance for the directional derivative
215*a99ef635SJonas Heinzmann -  -snes_linesearch_ltol <1e\-8>       - minimum absolute change in `lambda` allowed
216cd7522eaSPeter Brune 
217420bcc1bSBarry Smith    Level: advanced
218420bcc1bSBarry Smith 
219cd7522eaSPeter Brune    Notes:
220f6dfbefdSBarry Smith    This method does NOT use the objective function if it is provided with `SNESSetObjective()`.
221eb23ba39SBarry Smith 
222f6dfbefdSBarry Smith    This method is the preferred line search for `SNESQN` and `SNESNCG`.
223954494b2SJed Brown 
224b9d635d7SJonas Heinzmann .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLINESEARCHBISECTION`
225954494b2SJed Brown M*/
226d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
227d71ae5a4SJacob Faibussowitsch {
228bbd5d0b3SPeter Brune   PetscFunctionBegin;
229f1c6b773SPeter Brune   linesearch->ops->apply          = SNESLineSearchApply_CP;
2300298fd71SBarry Smith   linesearch->ops->destroy        = NULL;
2310298fd71SBarry Smith   linesearch->ops->setfromoptions = NULL;
2320298fd71SBarry Smith   linesearch->ops->reset          = NULL;
2330298fd71SBarry Smith   linesearch->ops->view           = NULL;
2340298fd71SBarry Smith   linesearch->ops->setup          = NULL;
235b000cd8dSPeter Brune   linesearch->order               = SNES_LINESEARCH_ORDER_LINEAR;
236a491bab8SPeter Brune 
237*a99ef635SJonas Heinzmann   linesearch->max_it = 1;
2383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
239bbd5d0b3SPeter Brune }
240