1a99ef635SJonas Heinzmann #include <petsc/private/linesearchimpl.h>
2*76c63389SBarry Smith #include <petsc/private/snesimpl.h>
3a99ef635SJonas Heinzmann #include <petscsnes.h>
4a99ef635SJonas Heinzmann
SNESLineSearchApply_Secant(SNESLineSearch linesearch)5a99ef635SJonas Heinzmann static PetscErrorCode SNESLineSearchApply_Secant(SNESLineSearch linesearch)
6a99ef635SJonas Heinzmann {
7a99ef635SJonas Heinzmann PetscBool changed_y, changed_w;
8a99ef635SJonas Heinzmann Vec X;
9a99ef635SJonas Heinzmann Vec F;
10a99ef635SJonas Heinzmann Vec Y;
11a99ef635SJonas Heinzmann Vec W;
12a99ef635SJonas Heinzmann SNES snes;
13a99ef635SJonas Heinzmann PetscReal gnorm;
14a99ef635SJonas Heinzmann PetscReal ynorm;
15a99ef635SJonas Heinzmann PetscReal xnorm;
16a99ef635SJonas Heinzmann PetscReal minlambda, maxlambda, atol, ltol;
17a99ef635SJonas Heinzmann PetscViewer monitor;
18a99ef635SJonas Heinzmann PetscReal lambda, lambda_old, lambda_mid, lambda_update, delLambda;
19a99ef635SJonas Heinzmann PetscReal fnrm, fnrm_old, fnrm_mid;
20a99ef635SJonas Heinzmann PetscReal delFnrm, delFnrm_old, del2Fnrm;
21a99ef635SJonas Heinzmann PetscInt i, max_it;
22a99ef635SJonas Heinzmann SNESObjectiveFn *objective;
23a99ef635SJonas Heinzmann
24a99ef635SJonas Heinzmann PetscFunctionBegin;
25a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetVecs(linesearch, &X, &F, &Y, &W, NULL));
26a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetNorms(linesearch, &xnorm, &gnorm, &ynorm));
27a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetLambda(linesearch, &lambda));
28a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetSNES(linesearch, &snes));
29a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_SUCCEEDED));
30a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetTolerances(linesearch, &minlambda, &maxlambda, NULL, &atol, <ol, &max_it));
31a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetDefaultMonitor(linesearch, &monitor));
32a99ef635SJonas Heinzmann
33a99ef635SJonas Heinzmann PetscCall(SNESGetObjective(snes, &objective, NULL));
34a99ef635SJonas Heinzmann
35a99ef635SJonas Heinzmann /* precheck */
36a99ef635SJonas Heinzmann PetscCall(SNESLineSearchPreCheck(linesearch, X, Y, &changed_y));
37a99ef635SJonas Heinzmann lambda_old = 0.0;
38a99ef635SJonas Heinzmann if (!objective) {
39a99ef635SJonas Heinzmann fnrm_old = gnorm * gnorm;
40a99ef635SJonas Heinzmann } else {
41a99ef635SJonas Heinzmann PetscCall(SNESComputeObjective(snes, X, &fnrm_old));
42*76c63389SBarry Smith SNESLineSearchCheckObjectiveDomainError(snes, fnrm_old);
43a99ef635SJonas Heinzmann }
44a99ef635SJonas Heinzmann lambda_mid = 0.5 * (lambda + lambda_old);
45a99ef635SJonas Heinzmann
46a99ef635SJonas Heinzmann for (i = 0; i < max_it; i++) {
47*76c63389SBarry Smith /* repeatedly cut lambda until the norm or objective function is not infinity or NaN or lambda is too small */
48a99ef635SJonas Heinzmann while (PETSC_TRUE) {
49a99ef635SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda_mid, Y, X));
50a99ef635SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
51a99ef635SJonas Heinzmann if (!objective) {
52a99ef635SJonas Heinzmann /* compute the norm at the midpoint */
53a99ef635SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
54a99ef635SJonas Heinzmann if (linesearch->ops->vinorm) {
55a99ef635SJonas Heinzmann fnrm_mid = gnorm;
56a99ef635SJonas Heinzmann PetscCall((*linesearch->ops->vinorm)(snes, F, W, &fnrm_mid));
57a99ef635SJonas Heinzmann } else {
58a99ef635SJonas Heinzmann PetscCall(VecNorm(F, NORM_2, &fnrm_mid));
59a99ef635SJonas Heinzmann }
60a99ef635SJonas Heinzmann
61a99ef635SJonas Heinzmann /* compute the norm at the new endpoint */
62a99ef635SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X));
63a99ef635SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
64a99ef635SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, W, F));
65a99ef635SJonas Heinzmann if (linesearch->ops->vinorm) {
66a99ef635SJonas Heinzmann fnrm = gnorm;
67a99ef635SJonas Heinzmann PetscCall((*linesearch->ops->vinorm)(snes, F, W, &fnrm));
68a99ef635SJonas Heinzmann } else {
69a99ef635SJonas Heinzmann PetscCall(VecNorm(F, NORM_2, &fnrm));
70a99ef635SJonas Heinzmann }
71a99ef635SJonas Heinzmann fnrm_mid = fnrm_mid * fnrm_mid;
72a99ef635SJonas Heinzmann fnrm = fnrm * fnrm;
73a99ef635SJonas Heinzmann } else {
74a99ef635SJonas Heinzmann /* compute the objective at the midpoint */
75a99ef635SJonas Heinzmann PetscCall(SNESComputeObjective(snes, W, &fnrm_mid));
76a99ef635SJonas Heinzmann
77a99ef635SJonas Heinzmann /* compute the objective at the new endpoint */
78a99ef635SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X));
79a99ef635SJonas Heinzmann PetscCall(SNESComputeObjective(snes, W, &fnrm));
80a99ef635SJonas Heinzmann }
81a99ef635SJonas Heinzmann
82a99ef635SJonas Heinzmann /* if new endpoint is viable, exit */
83a99ef635SJonas Heinzmann if (!PetscIsInfOrNanReal(fnrm)) break;
84a99ef635SJonas Heinzmann if (monitor) {
85a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
86*76c63389SBarry Smith PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: objective function at lambda = %g is infinity or NaN, cutting lambda\n", (double)lambda));
87a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
88a99ef635SJonas Heinzmann }
89a99ef635SJonas Heinzmann
90a99ef635SJonas Heinzmann /* if smallest allowable lambda gives NaN or Inf, exit line search */
91a99ef635SJonas Heinzmann if (lambda <= minlambda) {
92a99ef635SJonas Heinzmann if (monitor) {
93a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
94*76c63389SBarry Smith PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: objective function at lambda = %g <= lambda_min = %g is infinity or NaN, can not further cut lambda\n", (double)lambda, (double)lambda));
95a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
96a99ef635SJonas Heinzmann }
97a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_REDUCT));
98a99ef635SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
99a99ef635SJonas Heinzmann }
100a99ef635SJonas Heinzmann
101*76c63389SBarry Smith /* forbid the search from ever going back to the "failed" length that generates infinity or NaN */
102a99ef635SJonas Heinzmann maxlambda = .95 * lambda;
103a99ef635SJonas Heinzmann
104a99ef635SJonas Heinzmann /* shrink lambda towards the previous one which was viable */
105a99ef635SJonas Heinzmann lambda = .5 * (lambda + lambda_old);
106a99ef635SJonas Heinzmann lambda_mid = .5 * (lambda + lambda_old);
107a99ef635SJonas Heinzmann }
108a99ef635SJonas Heinzmann
109a99ef635SJonas Heinzmann /* monitoring output */
110a99ef635SJonas Heinzmann if (monitor) {
111a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
112a99ef635SJonas Heinzmann if (!objective) {
113a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: lambdas = [%g, %g, %g], fnorms = [%g, %g, %g]\n", (double)lambda, (double)lambda_mid, (double)lambda_old, (double)PetscSqrtReal(fnrm), (double)PetscSqrtReal(fnrm_mid), (double)PetscSqrtReal(fnrm_old)));
114a99ef635SJonas Heinzmann } else {
115a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: lambdas = [%g, %g, %g], obj = [%g, %g, %g]\n", (double)lambda, (double)lambda_mid, (double)lambda_old, (double)fnrm, (double)fnrm_mid, (double)fnrm_old));
116a99ef635SJonas Heinzmann }
117a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
118a99ef635SJonas Heinzmann }
119a99ef635SJonas Heinzmann
120a99ef635SJonas Heinzmann /* determine change of lambda */
121a99ef635SJonas Heinzmann delLambda = lambda - lambda_old;
122a99ef635SJonas Heinzmann
123a99ef635SJonas Heinzmann /* check change of lambda tolerance */
124a99ef635SJonas Heinzmann if (PetscAbsReal(delLambda) < ltol) {
125a99ef635SJonas Heinzmann if (monitor) {
126a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
127a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(delLambda) = %g < ltol = %g\n", (double)PetscAbsReal(delLambda), (double)ltol));
128a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
129a99ef635SJonas Heinzmann }
130a99ef635SJonas Heinzmann break;
131a99ef635SJonas Heinzmann }
132a99ef635SJonas Heinzmann
133a99ef635SJonas Heinzmann /* compute f'() at the end points using second order one sided differencing */
134a99ef635SJonas Heinzmann delFnrm = (3. * fnrm - 4. * fnrm_mid + 1. * fnrm_old) / delLambda;
135a99ef635SJonas Heinzmann delFnrm_old = (-3. * fnrm_old + 4. * fnrm_mid - 1. * fnrm) / delLambda;
136a99ef635SJonas Heinzmann
137a99ef635SJonas Heinzmann /* compute f''() at the midpoint using centered differencing */
138a99ef635SJonas Heinzmann del2Fnrm = (delFnrm - delFnrm_old) / delLambda;
139a99ef635SJonas Heinzmann
140a99ef635SJonas Heinzmann /* check absolute tolerance */
141a99ef635SJonas Heinzmann if (PetscAbsReal(delFnrm) <= atol) {
142a99ef635SJonas Heinzmann if (monitor) {
143a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
144a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: abs(delFnrm) = %g <= atol = %g\n", (double)PetscAbsReal(delFnrm), (double)atol));
145a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
146a99ef635SJonas Heinzmann }
147a99ef635SJonas Heinzmann break;
148a99ef635SJonas Heinzmann }
149a99ef635SJonas Heinzmann
150a99ef635SJonas Heinzmann /* compute the secant (Newton) update -- always go downhill */
151a99ef635SJonas Heinzmann if (del2Fnrm > 0.) lambda_update = lambda - delFnrm / del2Fnrm;
152a99ef635SJonas Heinzmann else if (del2Fnrm < 0.) lambda_update = lambda + delFnrm / del2Fnrm;
153a99ef635SJonas Heinzmann else {
154a99ef635SJonas Heinzmann if (monitor) {
155a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
156a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: del2Fnrm = 0\n"));
157a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
158a99ef635SJonas Heinzmann }
159a99ef635SJonas Heinzmann break;
160a99ef635SJonas Heinzmann }
161a99ef635SJonas Heinzmann
162a99ef635SJonas Heinzmann /* prevent secant method from stepping out of bounds */
163a99ef635SJonas Heinzmann if (lambda_update < minlambda) lambda_update = 0.5 * (lambda + lambda_old);
164a99ef635SJonas Heinzmann if (lambda_update > maxlambda) {
165a99ef635SJonas Heinzmann lambda_update = maxlambda;
166a99ef635SJonas Heinzmann break;
167a99ef635SJonas Heinzmann }
168a99ef635SJonas Heinzmann
169a99ef635SJonas Heinzmann /* if lambda is NaN or Inf, do not accept update but exit with previous lambda */
170a99ef635SJonas Heinzmann if (PetscIsInfOrNanReal(lambda_update)) {
171a99ef635SJonas Heinzmann if (monitor) {
172a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
173a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: lambda_update is NaN or Inf\n"));
174a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
175a99ef635SJonas Heinzmann }
176a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetReason(linesearch, SNES_LINESEARCH_FAILED_NANORINF));
177a99ef635SJonas Heinzmann break;
178a99ef635SJonas Heinzmann }
179a99ef635SJonas Heinzmann
180a99ef635SJonas Heinzmann /* update the endpoints and the midpoint of the bracketed secant region */
181a99ef635SJonas Heinzmann lambda_old = lambda;
182a99ef635SJonas Heinzmann lambda = lambda_update;
183a99ef635SJonas Heinzmann fnrm_old = fnrm;
184a99ef635SJonas Heinzmann lambda_mid = 0.5 * (lambda + lambda_old);
185a99ef635SJonas Heinzmann
186a99ef635SJonas Heinzmann if ((i == max_it - 1) && monitor) {
187a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
188a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search: reached maximum number of iterations!\n"));
189a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
190a99ef635SJonas Heinzmann }
191a99ef635SJonas Heinzmann }
192a99ef635SJonas Heinzmann
193a99ef635SJonas Heinzmann /* construct the solution */
194a99ef635SJonas Heinzmann PetscCall(VecWAXPY(W, -lambda, Y, X));
195a99ef635SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
196a99ef635SJonas Heinzmann
197a99ef635SJonas Heinzmann /* postcheck */
198a99ef635SJonas Heinzmann PetscCall(SNESLineSearchSetLambda(linesearch, lambda));
199a99ef635SJonas Heinzmann PetscCall(SNESLineSearchPostCheck(linesearch, X, Y, W, &changed_y, &changed_w));
200a99ef635SJonas Heinzmann if (changed_y) {
201a99ef635SJonas Heinzmann if (!changed_w) PetscCall(VecWAXPY(W, -lambda, Y, X));
202a99ef635SJonas Heinzmann if (linesearch->ops->viproject) PetscCall((*linesearch->ops->viproject)(snes, W));
203a99ef635SJonas Heinzmann }
204a99ef635SJonas Heinzmann PetscCall(VecCopy(W, X));
205a99ef635SJonas Heinzmann PetscCall((*linesearch->ops->snesfunc)(snes, X, F));
206a99ef635SJonas Heinzmann
207a99ef635SJonas Heinzmann PetscCall(SNESLineSearchComputeNorms(linesearch));
208a99ef635SJonas Heinzmann PetscCall(SNESLineSearchGetNorms(linesearch, NULL, &gnorm, NULL));
209*76c63389SBarry Smith SNESLineSearchCheckFunctionDomainError(snes, linesearch, gnorm);
210*76c63389SBarry Smith if (monitor) {
211a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIAddTab(monitor, ((PetscObject)linesearch)->tablevel));
212a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIIPrintf(monitor, " Line search terminated: lambda = %g, fnorm = %g\n", (double)lambda, (double)gnorm));
213a99ef635SJonas Heinzmann PetscCall(PetscViewerASCIISubtractTab(monitor, ((PetscObject)linesearch)->tablevel));
214a99ef635SJonas Heinzmann }
215a99ef635SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
216a99ef635SJonas Heinzmann }
217a99ef635SJonas Heinzmann
218a99ef635SJonas Heinzmann /*MC
219a99ef635SJonas Heinzmann SNESLINESEARCHSECANT - Secant search in the L2 norm of the function or the objective function.
220a99ef635SJonas Heinzmann
221a99ef635SJonas Heinzmann Attempts to solve $ \min_{\lambda} f(x_k + \lambda Y_k) $ using the secant method with the initial bracketing of $ \lambda $ between [0,damping].
222a99ef635SJonas Heinzmann $f(x_k + \lambda Y_k)$ is either the squared L2-norm of the function $||F(x_k + \lambda Y_k)||_2^2$,
223a99ef635SJonas Heinzmann or the objective function $G(x_k + \lambda Y_k)$ if it is provided with `SNESSetObjective()`
224a99ef635SJonas Heinzmann Differences of $f()$ are used to approximate the first and second derivative of $f()$ with respect to
225a99ef635SJonas Heinzmann $\lambda$, $f'()$ and $f''()$.
226a99ef635SJonas Heinzmann
227a99ef635SJonas Heinzmann When an objective function is provided $f(w)$ is the objective function otherwise $f(w) = ||F(w)||^2$.
228a99ef635SJonas Heinzmann $x$ is the current step and $y$ is the search direction.
229a99ef635SJonas Heinzmann
230a99ef635SJonas Heinzmann Options Database Keys:
231a99ef635SJonas Heinzmann + -snes_linesearch_max_it <1> - maximum number of iterations within the line search
232a99ef635SJonas Heinzmann . -snes_linesearch_damping <1.0> - initial `lambda` on entry to the line search
233a99ef635SJonas Heinzmann . -snes_linesearch_minlambda <1e\-12> - minimum allowable `lambda`
234a99ef635SJonas Heinzmann . -snes_linesearch_maxlambda <1.0> - maximum `lambda` (scaling of solution update) allowed
235a99ef635SJonas Heinzmann . -snes_linesearch_atol <1e\-15> - absolute tolerance for the secant method $ f'() < atol $
236a99ef635SJonas Heinzmann - -snes_linesearch_ltol <1e\-8> - minimum absolute change in `lambda` allowed
237a99ef635SJonas Heinzmann
238a99ef635SJonas Heinzmann Level: advanced
239a99ef635SJonas Heinzmann
240a99ef635SJonas Heinzmann .seealso: [](ch_snes), `SNESLINESEARCHBT`, `SNESLINESEARCHCP`, `SNESLineSearch`, `SNESLineSearchType`, `SNESLineSearchCreate()`, `SNESLineSearchSetType()`
241a99ef635SJonas Heinzmann M*/
SNESLineSearchCreate_Secant(SNESLineSearch linesearch)242a99ef635SJonas Heinzmann PETSC_EXTERN PetscErrorCode SNESLineSearchCreate_Secant(SNESLineSearch linesearch)
243a99ef635SJonas Heinzmann {
244a99ef635SJonas Heinzmann PetscFunctionBegin;
245a99ef635SJonas Heinzmann linesearch->ops->apply = SNESLineSearchApply_Secant;
246a99ef635SJonas Heinzmann linesearch->ops->destroy = NULL;
247a99ef635SJonas Heinzmann linesearch->ops->setfromoptions = NULL;
248a99ef635SJonas Heinzmann linesearch->ops->reset = NULL;
249a99ef635SJonas Heinzmann linesearch->ops->view = NULL;
250a99ef635SJonas Heinzmann linesearch->ops->setup = NULL;
251a99ef635SJonas Heinzmann
252a99ef635SJonas Heinzmann linesearch->max_it = 1;
253a99ef635SJonas Heinzmann PetscFunctionReturn(PETSC_SUCCESS);
254a99ef635SJonas Heinzmann }
255