1af0996ceSBarry Smith #include <petsc/private/linesearchimpl.h>
26a388c36SPeter Brune #include <petscsnes.h>
3bbd5d0b3SPeter Brune
SNESLineSearchApply_CP(SNESLineSearch linesearch)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;
9a99ef635SJonas 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;
12a99ef635SJonas 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));
20a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, &maxlambda, &rtol, &atol, <ol, &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
28a99ef635SJonas 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
45a99ef635SJonas 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
56a99ef635SJonas Heinzmann /* compute change of lambda */
57bbd5d0b3SPeter Brune delLambda = lambda - lambda_old;
58b8ac21a2SPeter Brune
59a99ef635SJonas Heinzmann /* check change of lambda tolerance */
60a99ef635SJonas Heinzmann if (PetscAbsReal(delLambda) < ltol) {
61a99ef635SJonas Heinzmann if (monitor) {
62a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
63a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(dlambda) = %g < ltol = %g\n", (double)PetscAbsReal(delLambda), (double)ltol));
64a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
65a99ef635SJonas Heinzmann }
66a99ef635SJonas Heinzmann break;
67a99ef635SJonas Heinzmann }
68a99ef635SJonas Heinzmann
69a99ef635SJonas Heinzmann /* check relative tolerance */
70a99ef635SJonas Heinzmann if (PetscAbsScalar(fty) / PetscAbsScalar(fty_init) < rtol) {
71a99ef635SJonas Heinzmann if (monitor) {
72a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
73a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(fty/fty_init) = %g <= rtol = %g\n", (double)(PetscAbsScalar(fty) / PetscAbsScalar(fty_init)), (double)rtol));
74a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
75a99ef635SJonas Heinzmann }
76a99ef635SJonas Heinzmann break;
77a99ef635SJonas Heinzmann }
78a99ef635SJonas Heinzmann
79a99ef635SJonas Heinzmann /* check absolute tolerance */
80a99ef635SJonas Heinzmann if (PetscAbsScalar(fty) < atol * ynorm && i > 0) {
81a99ef635SJonas Heinzmann if (monitor) {
82a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
83a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(fty)/||y|| = %g <= atol = %g\n", (double)(PetscAbsScalar(fty) / ynorm), (double)atol));
84a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
85a99ef635SJonas Heinzmann }
86a99ef635SJonas Heinzmann break;
87a99ef635SJonas Heinzmann }
88a99ef635SJonas Heinzmann
89a99ef635SJonas 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
96a99ef635SJonas Heinzmann /* approximate the second derivative */
97b000cd8dSPeter Brune if (linesearch->order == SNES_LINESEARCH_ORDER_LINEAR) {
98a99ef635SJonas Heinzmann /* first order finite difference approximation */
99b8ac21a2SPeter Brune s = (fty - fty_old) / delLambda;
100b000cd8dSPeter Brune } else if (linesearch->order == SNES_LINESEARCH_ORDER_QUADRATIC) {
101a99ef635SJonas 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 }
110a99ef635SJonas Heinzmann
111a99ef635SJonas Heinzmann /* second order finite difference approximation */
112b8ac21a2SPeter Brune s = (3. * fty - 4. * fty_mid1 + fty_old) / delLambda;
113a99ef635SJonas Heinzmann
114b8ac21a2SPeter Brune } else {
115a99ef635SJonas 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 }
124a99ef635SJonas Heinzmann
125a99ef635SJonas 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 }
134a99ef635SJonas Heinzmann
135a99ef635SJonas Heinzmann /* third order finite difference approximation */
136b8ac21a2SPeter Brune s = (2. * fty_mid2 + 3. * fty - 6. * fty_mid1 + fty_old) / (3. * delLambda);
137b8ac21a2SPeter Brune }
138a99ef635SJonas Heinzmann
139a99ef635SJonas Heinzmann /* compute secant update (modifying the search direction if necessary) */
140b8ac21a2SPeter Brune if (PetscRealPart(s) > 0.) s = -s;
141a99ef635SJonas Heinzmann if (s == 0.0) {
142a99ef635SJonas Heinzmann if (monitor) {
143a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
144a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: del2Fnrm = 0\n"));
145a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
146a99ef635SJonas Heinzmann }
147a99ef635SJonas Heinzmann break;
148a99ef635SJonas Heinzmann }
149b8ac21a2SPeter Brune lambda_update = lambda - PetscRealPart(fty / s);
150b8ac21a2SPeter Brune
1511ab86304SZach Atkins /* if step is too small, go the opposite direction */
1521ab86304SZach Atkins if (lambda_update < minlambda) lambda_update = lambda + PetscRealPart(fty / s);
153a99ef635SJonas Heinzmann /* if secant method would step out of bounds, exit with the respective bound */
154a99ef635SJonas Heinzmann if (lambda_update > maxlambda) {
1551ab86304SZach Atkins lambda = maxlambda;
156a99ef635SJonas Heinzmann break;
157a99ef635SJonas Heinzmann }
158b8ac21a2SPeter Brune
159*76c63389SBarry Smith /* if lambda is infinity or NaN, do not accept update but exit with previous lambda */
160a99ef635SJonas Heinzmann if (PetscIsInfOrNanReal(lambda_update)) {
161a99ef635SJonas Heinzmann if (monitor) {
162a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
163*76c63389SBarry Smith PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: lambda_update is infinity or NaN\n"));
164a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
165a99ef635SJonas Heinzmann }
166a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
167a99ef635SJonas Heinzmann break;
168a99ef635SJonas Heinzmann }
169bbd5d0b3SPeter Brune
170bbd5d0b3SPeter Brune /* compute the new state of the line search */
171bbd5d0b3SPeter Brune lambda_old = lambda;
172bbd5d0b3SPeter Brune lambda = lambda_update;
173bbd5d0b3SPeter Brune fty_old = fty;
174bbd5d0b3SPeter Brune }
175a99ef635SJonas Heinzmann
176bbd5d0b3SPeter Brune /* construct the solution */
177ef46b1a6SStefano Zampini PetscCall(VecWAXPY(W, -lambda, Y, X));
1781baa6e33SBarry Smith if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
179bbd5d0b3SPeter Brune /* postcheck */
1806b095a85SStefano Zampini PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
1819566063dSJacob Faibussowitsch PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
1826b095a85SStefano Zampini if (changed_y) {
1836b095a85SStefano Zampini if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
1846b095a85SStefano Zampini if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
185bbd5d0b3SPeter Brune }
1866b095a85SStefano Zampini PetscCall(VecCopy(W, X));
1879566063dSJacob Faibussowitsch PetscCall((*linesearch->ops->snesfunc)(snes, X, F));
188bbd5d0b3SPeter Brune
1899566063dSJacob Faibussowitsch PetscCall(SNESLineSearchComputeNorms(linesearch));
1909566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
191bbd5d0b3SPeter Brune
1926a388c36SPeter Brune if (monitor) {
1939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
1949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(monitor, " Line search terminated: lambda = %g, fnorms = %g\n", (double)lambda, (double)gnorm));
1959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
196bbd5d0b3SPeter Brune }
1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
198bbd5d0b3SPeter Brune }
199bbd5d0b3SPeter Brune
200954494b2SJed Brown /*MC
2011a4f838cSPeter Brune SNESLINESEARCHCP - Critical point line search. This line search assumes that there exists some
202d5def619SJonas Heinzmann artificial $G(x)$ for which the `SNESFunctionFn` $F(x) = grad G(x)$. Therefore, this line search seeks
203a99ef635SJonas 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$.
204cd7522eaSPeter Brune
205cd7522eaSPeter Brune Options Database Keys:
206a99ef635SJonas Heinzmann + -snes_linesearch_minlambda <1e\-12> - the minimum acceptable `lambda` (scaling of solution update)
207a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda <1.0> - the algorithm ensures that `lambda` is never larger than this value
208a99ef635SJonas Heinzmann . -snes_linesearch_damping <1.0> - initial `lambda` on entry to the line search
209a99ef635SJonas Heinzmann . -snes_linesearch_order <1> - order of the approximation in the secant method, must be 1, 2, or 3
210a99ef635SJonas Heinzmann . -snes_linesearch_max_it <1> - the maximum number of secant iterations performed
211a99ef635SJonas Heinzmann . -snes_linesearch_rtol <1e\-8> - relative tolerance for the directional derivative
212a99ef635SJonas Heinzmann . -snes_linesearch_atol <1e\-15> - absolute tolerance for the directional derivative
213a99ef635SJonas Heinzmann - -snes_linesearch_ltol <1e\-8> - minimum absolute change in `lambda` allowed
214cd7522eaSPeter Brune
215420bcc1bSBarry Smith Level: advanced
216420bcc1bSBarry Smith
217cd7522eaSPeter Brune Notes:
218f6dfbefdSBarry Smith This method does NOT use the objective function if it is provided with `SNESSetObjective()`.
219eb23ba39SBarry Smith
220f6dfbefdSBarry Smith This method is the preferred line search for `SNESQN` and `SNESNCG`.
221954494b2SJed Brown
222b9d635d7SJonas Heinzmann .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLINESEARCHBISECTION`
223954494b2SJed Brown M*/
SNESLineSearchCreate_CP(SNESLineSearch linesearch)224d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_CP(SNESLineSearch linesearch)
225d71ae5a4SJacob Faibussowitsch {
226bbd5d0b3SPeter Brune PetscFunctionBegin;
227f1c6b773SPeter Brune linesearch->ops->apply = SNESLineSearchApply_CP;
2280298fd71SBarry Smith linesearch->ops->destroy = NULL;
2290298fd71SBarry Smith linesearch->ops->setfromoptions = NULL;
2300298fd71SBarry Smith linesearch->ops->reset = NULL;
2310298fd71SBarry Smith linesearch->ops->view = NULL;
2320298fd71SBarry Smith linesearch->ops->setup = NULL;
233b000cd8dSPeter Brune linesearch->order = SNES_LINESEARCH_ORDER_LINEAR;
234a491bab8SPeter Brune
235a99ef635SJonas Heinzmann linesearch->max_it = 1;
2363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
237bbd5d0b3SPeter Brune }
238