1b9d635d7SJonas Heinzmann #include <petsc/private/linesearchimpl.h>
2b9d635d7SJonas Heinzmann #include <petsc/private/snesimpl.h>
3b9d635d7SJonas Heinzmann
SNESLineSearchApply_Bisection(SNESLineSearch linesearch)4b9d635d7SJonas Heinzmann static PetscErrorCode SNESLineSearchApply_Bisection(SNESLineSearch linesearch)
5b9d635d7SJonas Heinzmann {
6b9d635d7SJonas Heinzmann PetscBool changed_y, changed_w;
7b9d635d7SJonas Heinzmann Vec X, F, Y, W, G;
8b9d635d7SJonas Heinzmann SNES snes;
9b9d635d7SJonas Heinzmann PetscReal ynorm;
1076c63389SBarry Smith PetscReal lambda_left, lambda, lambda_right, lambda_old, fnorm;
11d5def619SJonas Heinzmann PetscScalar fty_left, fty, fty_initial;
12b9d635d7SJonas Heinzmann PetscViewer monitor;
13b9d635d7SJonas Heinzmann PetscReal rtol, atol, ltol;
14a99ef635SJonas Heinzmann PetscInt it, max_it;
15b9d635d7SJonas Heinzmann
16b9d635d7SJonas Heinzmann PetscFunctionBegin;
17b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, &G));
18b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
19b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
20a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetTolerances(linesearch, NULL, NULL, &rtol, &atol, <ol, &max_it));
21b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
22*1e1ea152SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
23b9d635d7SJonas Heinzmann
24b9d635d7SJonas Heinzmann /* pre-check */
25b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
26b9d635d7SJonas Heinzmann
27b9d635d7SJonas Heinzmann /* compute ynorm to normalize search direction */
28b9d635d7SJonas Heinzmann PetscCall(VecNorm(Y, NORM_2, &ynorm));
29b9d635d7SJonas Heinzmann
30b9d635d7SJonas Heinzmann /* initialize interval for bisection */
31b9d635d7SJonas Heinzmann lambda_left = 0.0;
32b9d635d7SJonas Heinzmann lambda_right = lambda;
33b9d635d7SJonas Heinzmann
34b9d635d7SJonas Heinzmann /* compute fty at left end of interval */
35d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) {
36d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, F, X, Y, &fty_left));
37d5def619SJonas Heinzmann } else {
38d5def619SJonas Heinzmann PetscCall(VecDot(F, Y, &fty_left));
39d5def619SJonas Heinzmann }
40b9d635d7SJonas Heinzmann fty_initial = fty_left;
41b9d635d7SJonas Heinzmann
42b9d635d7SJonas Heinzmann /* compute fty at right end of interval (initial lambda) */
43b9d635d7SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X));
44b9d635d7SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
45b9d635d7SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
46b9d635d7SJonas Heinzmann if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
47b9d635d7SJonas Heinzmann PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations during line search!\n"));
48b9d635d7SJonas Heinzmann snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
49b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
50b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
51b9d635d7SJonas Heinzmann }
52d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) {
53d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, G, W, Y, &fty));
54d5def619SJonas Heinzmann } else {
55d5def619SJonas Heinzmann PetscCall(VecDot(G, Y, &fty));
56d5def619SJonas Heinzmann }
57b9d635d7SJonas Heinzmann /* check whether sign changes in interval */
58d5def619SJonas Heinzmann if (!PetscIsInfOrNanScalar(fty) && (PetscRealPart(fty_left * fty) > 0.0)) {
59b9d635d7SJonas Heinzmann /* no change of sign: accept full step */
60b9d635d7SJonas Heinzmann if (monitor) {
61b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
62d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: sign of fty does not change in step interval, accepting full step\n"));
63b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
64b9d635d7SJonas Heinzmann }
65b9d635d7SJonas Heinzmann } else {
66b9d635d7SJonas Heinzmann /* change of sign: iteratively bisect interval */
67b9d635d7SJonas Heinzmann lambda_old = 0.0;
68b9d635d7SJonas Heinzmann it = 0;
69b9d635d7SJonas Heinzmann
70b9d635d7SJonas Heinzmann while (PETSC_TRUE) {
7176c63389SBarry Smith /* check for infinity or NaN */
72d5def619SJonas Heinzmann if (PetscIsInfOrNanScalar(fty)) {
73b9d635d7SJonas Heinzmann if (monitor) {
74b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
7576c63389SBarry Smith PetscCall(PetscViewerASCIIPrintf(monitor, " Line search fty is infinity or NaN!\n"));
76b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
77b9d635d7SJonas Heinzmann }
7876c63389SBarry Smith PetscCheck(!snes->errorifnotconverged, PetscObjectComm((PetscObject)snes), PETSC_ERR_CONV_FAILED, "infinity or NaN in function evaluation");
7976c63389SBarry Smith if (snes->functiondomainerror) {
8076c63389SBarry Smith PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION_DOMAIN));
8176c63389SBarry Smith snes->functiondomainerror = PETSC_FALSE;
8276c63389SBarry Smith } else PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
83b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
84b9d635d7SJonas Heinzmann break;
85b9d635d7SJonas Heinzmann }
86b9d635d7SJonas Heinzmann
87b9d635d7SJonas Heinzmann /* check absolute tolerance */
88d5def619SJonas Heinzmann if (PetscAbsScalar(fty) <= atol * ynorm) {
89b9d635d7SJonas Heinzmann if (monitor) {
90b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
91d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(fty)/||y|| = %g <= atol = %g\n", (double)(PetscAbsScalar(fty) / ynorm), (double)atol));
92b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
93b9d635d7SJonas Heinzmann }
94b9d635d7SJonas Heinzmann break;
95b9d635d7SJonas Heinzmann }
96b9d635d7SJonas Heinzmann
97b9d635d7SJonas Heinzmann /* check relative tolerance */
98d5def619SJonas Heinzmann if (PetscAbsScalar(fty) / PetscAbsScalar(fty_initial) <= rtol) {
99b9d635d7SJonas Heinzmann if (monitor) {
100b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
101d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(fty/fty_initial) = %g <= rtol = %g\n", (double)(PetscAbsScalar(fty) / PetscAbsScalar(fty_initial)), (double)rtol));
102b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
103b9d635d7SJonas Heinzmann }
104b9d635d7SJonas Heinzmann break;
105b9d635d7SJonas Heinzmann }
106b9d635d7SJonas Heinzmann
107b9d635d7SJonas Heinzmann /* check maximum number of iterations */
108a99ef635SJonas Heinzmann if (it > max_it) {
109b9d635d7SJonas Heinzmann if (monitor) {
110b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
111b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: maximum iterations reached\n"));
112b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
113b9d635d7SJonas Heinzmann }
114b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
115b9d635d7SJonas Heinzmann break;
116b9d635d7SJonas Heinzmann }
117b9d635d7SJonas Heinzmann
118b9d635d7SJonas Heinzmann /* check change of lambda tolerance */
119b9d635d7SJonas Heinzmann if (PetscAbsReal(lambda - lambda_old) < ltol) {
120b9d635d7SJonas Heinzmann if (monitor) {
121b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
122b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(dlambda) = %g < ltol = %g\n", (double)PetscAbsReal(lambda - lambda_old), (double)ltol));
123b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
124b9d635d7SJonas Heinzmann }
125b9d635d7SJonas Heinzmann break;
126b9d635d7SJonas Heinzmann }
127b9d635d7SJonas Heinzmann
128b9d635d7SJonas Heinzmann /* determine direction of bisection (not necessary for 0th iteration) */
129b9d635d7SJonas Heinzmann if (it > 0) {
130d5def619SJonas Heinzmann if (PetscRealPart(fty * fty_left) <= 0.0) {
131b9d635d7SJonas Heinzmann lambda_right = lambda;
132b9d635d7SJonas Heinzmann } else {
133b9d635d7SJonas Heinzmann lambda_left = lambda;
134b9d635d7SJonas Heinzmann /* also update fty_left for direction check in next iteration */
135b9d635d7SJonas Heinzmann fty_left = fty;
136b9d635d7SJonas Heinzmann }
137b9d635d7SJonas Heinzmann }
138b9d635d7SJonas Heinzmann
139b9d635d7SJonas Heinzmann /* bisect interval */
140b9d635d7SJonas Heinzmann lambda_old = lambda;
141b9d635d7SJonas Heinzmann lambda = 0.5 * (lambda_left + lambda_right);
142b9d635d7SJonas Heinzmann
143b9d635d7SJonas Heinzmann /* compute fty at new lambda */
144b9d635d7SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X));
145b9d635d7SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
146b9d635d7SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, W, G));
147b9d635d7SJonas Heinzmann if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
148b9d635d7SJonas Heinzmann PetscCall(PetscInfo(snes, "Exceeded maximum function evaluations during line search!\n"));
149b9d635d7SJonas Heinzmann snes->reason = SNES_DIVERGED_FUNCTION_COUNT;
150b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_FUNCTION));
151b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
152b9d635d7SJonas Heinzmann }
153d5def619SJonas Heinzmann if (linesearch->ops->vidirderiv) {
154d5def619SJonas Heinzmann PetscCall((*linesearch->ops->vidirderiv)(snes, G, W, Y, &fty));
155d5def619SJonas Heinzmann } else {
156d5def619SJonas Heinzmann PetscCall(VecDot(G, Y, &fty));
157d5def619SJonas Heinzmann }
158b9d635d7SJonas Heinzmann
159b9d635d7SJonas Heinzmann /* print iteration information */
160b9d635d7SJonas Heinzmann if (monitor) {
161b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
162d5def619SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " %3" PetscInt_FMT " Line search: fty/||y|| = %g, lambda = %g\n", it, (double)(PetscRealPart(fty) / ynorm), (double)lambda));
163b9d635d7SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
164b9d635d7SJonas Heinzmann }
165b9d635d7SJonas Heinzmann
166b9d635d7SJonas Heinzmann /* count up iteration */
167b9d635d7SJonas Heinzmann it++;
168b9d635d7SJonas Heinzmann }
169b9d635d7SJonas Heinzmann }
170b9d635d7SJonas Heinzmann
171b9d635d7SJonas Heinzmann /* post-check */
172b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
173b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
174b9d635d7SJonas Heinzmann if (changed_y) {
175b9d635d7SJonas Heinzmann if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
176b9d635d7SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
177b9d635d7SJonas Heinzmann }
178b9d635d7SJonas Heinzmann
179b9d635d7SJonas Heinzmann /* update solution*/
180b9d635d7SJonas Heinzmann PetscCall(VecCopy(W, X));
181b9d635d7SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, X, F));
182b9d635d7SJonas Heinzmann PetscCall(SNESLineSearchComputeNorms(linesearch));
18376c63389SBarry Smith PetscCall(SNESLineSearchGetNorms(linesearch, NULL, &fnorm, NULL));
18476c63389SBarry Smith SNESLineSearchCheckFunctionDomainError(snes, linesearch, fnorm);
185b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
186b9d635d7SJonas Heinzmann }
187b9d635d7SJonas Heinzmann
188b9d635d7SJonas Heinzmann /*MC
189b9d635d7SJonas Heinzmann SNESLINESEARCHBISECTION - Bisection line search.
190d5def619SJonas Heinzmann Similar to the critical point line search, `SNESLINESEARCHCP`, the bisection line search assumes that there exists some $G(x)$ for which the `SNESFunctionFn` $F(x) = grad G(x)$.
191a99ef635SJonas Heinzmann This line search seeks to find the root of the directional derivative, that is $F(x_k - \lambda Y_k) \cdot Y_k / ||Y_k|| = 0$, along the search direction $Y_k$ through bisection.
192b9d635d7SJonas Heinzmann
193b9d635d7SJonas Heinzmann Options Database Keys:
194a99ef635SJonas Heinzmann + -snes_linesearch_max_it <50> - maximum number of bisection iterations for the line search
195a99ef635SJonas Heinzmann . -snes_linesearch_damping <1.0> - initial `lambda` on entry to the line search
196d5def619SJonas Heinzmann . -snes_linesearch_rtol <1e\-8> - relative tolerance for the directional derivative
197d5def619SJonas Heinzmann . -snes_linesearch_atol <1e\-6> - absolute tolerance for the directional derivative
198*1e1ea152SJonas Heinzmann - -snes_linesearch_ltol <1e\-6> - minimum absolute change in `lambda` allowed (this is an alternative to setting a maximum number of iterations)
199b9d635d7SJonas Heinzmann
200b9d635d7SJonas Heinzmann Level: intermediate
201b9d635d7SJonas Heinzmann
202a99ef635SJonas Heinzmann Notes:
203a99ef635SJonas Heinzmann `lambda` is the scaling of the search direction (vector) that is computed by this algorithm.
204a99ef635SJonas Heinzmann If there is no change of sign in the directional derivative from $\lambda=0$ to the initial `lambda` (the damping), then the initial `lambda` will be used.
205a99ef635SJonas Heinzmann Hence, this line search will always give a `lambda` in the interval $[0, damping]$.
206b9d635d7SJonas Heinzmann This method does NOT use the objective function if it is provided with `SNESSetObjective()`.
207b9d635d7SJonas Heinzmann
208b9d635d7SJonas Heinzmann .seealso: [](ch_snes), `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`, `SNESLINESEARCHCP`
209b9d635d7SJonas Heinzmann M*/
SNESLineSearchCreate_Bisection(SNESLineSearch linesearch)210b9d635d7SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Bisection(SNESLineSearch linesearch)
211b9d635d7SJonas Heinzmann {
212b9d635d7SJonas Heinzmann PetscFunctionBegin;
213b9d635d7SJonas Heinzmann linesearch->ops->apply = SNESLineSearchApply_Bisection;
214b9d635d7SJonas Heinzmann linesearch->ops->destroy = NULL;
215b9d635d7SJonas Heinzmann linesearch->ops->setfromoptions = NULL;
216b9d635d7SJonas Heinzmann linesearch->ops->reset = NULL;
217b9d635d7SJonas Heinzmann linesearch->ops->view = NULL;
218b9d635d7SJonas Heinzmann linesearch->ops->setup = NULL;
219b9d635d7SJonas Heinzmann
220b9d635d7SJonas Heinzmann /* set default option values */
221a99ef635SJonas Heinzmann linesearch->max_it = 50;
222b9d635d7SJonas Heinzmann linesearch->rtol = 1e-8;
223b9d635d7SJonas Heinzmann linesearch->atol = 1e-6;
224b9d635d7SJonas Heinzmann linesearch->ltol = 1e-6;
225b9d635d7SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
226b9d635d7SJonas Heinzmann }
227