Actual source code: owarmijo.h
1: #pragma once
3: /* Context for an Armijo (nonmonotone) linesearch for orthant wise unconstrained
4: minimization.
6: Given a function f, the current iterate x, and a descent direction d:
7: Find the smallest i in 0, 1, 2, ..., such that:
9: f(x + (beta**i)d) <= f(x) + (sigma*beta**i)<grad f(x),d>
11: The nonmonotone modification of this linesearch replaces the f(x) term
12: with a reference value, R, and seeks to find the smallest i such that:
14: f(x + (beta**i)d) <= R + (sigma*beta**i)<grad f(x),d>
16: This modification does effect neither the convergence nor rate of
17: convergence of an algorithm when R is chosen appropriately. Essentially,
18: R must decrease on average in some sense. The benefit of a nonmonotone
19: linesearch is that local minimizers can be avoided (by allowing increase
20: in function value), and typically, fewer iterations are performed in
21: the main code.
23: The reference value is chosen based upon some historical information
24: consisting of function values for previous iterates. The amount of
25: historical information used is determined by the memory size where the
26: memory is used to store the previous function values. The memory is
27: initialized to alpha*f(x^0) for some alpha >= 1, with alpha=1 signifying
28: that we always force decrease from the initial point.
30: The reference value can be the maximum value in the memory or can be
31: chosen to provide some mean descent. Elements are removed from the
32: memory with a replacement policy that either removes the oldest
33: value in the memory (FIFO), or the largest value in the memory (MRU).
35: Additionally, we can add a watchdog strategy to the search, which
36: essentially accepts small directions and only checks the nonmonotonic
37: descent criteria every m-steps. This strategy is NOT implemented in
38: the code.
40: Finally, care must be taken when steepest descent directions are used.
41: For example, when the Newton direction does not satisfy a sufficient
42: descent criteria. The code will apply the same test regardless of
43: the direction. This type of search may not be appropriate for all
44: algorithms. For example, when a gradient direction is used, we may
45: want to revert to the best point found and reset the memory so that
46: we stay in an appropriate level set after using a gradient steps.
47: This type of search is currently NOT supported by the code.
49: References:
50: + * - Armijo, "Minimization of Functions Having Lipschitz Continuous
51: First-Partial Derivatives," Pacific Journal of Mathematics, volume 16,
52: pages 1-3, 1966.
53: . * - Ferris and Lucidi, "Nonmonotone Stabilization Methods for Nonlinear
54: Equations," Journal of Optimization Theory and Applications, volume 81,
55: pages 53-71, 1994.
56: . * - Grippo, Lampariello, and Lucidi, "A Nonmonotone Line Search Technique
57: for Newton's Method," SIAM Journal on Numerical Analysis, volume 23,
58: pages 707-716, 1986.
59: - * - Grippo, Lampariello, and Lucidi, "A Class of Nonmonotone Stabilization
60: Methods in Unconstrained Optimization," Numerische Mathematik, volume 59,
61: pages 779-805, 1991. */
62: #include <petsc/private/taolinesearchimpl.h>
63: typedef struct {
64: PetscReal *memory;
66: PetscReal alpha; /* Initial reference factor >= 1 */
67: PetscReal beta; /* Steplength determination < 1 */
68: PetscReal beta_inf; /* Steplength determination < 1 */
69: PetscReal sigma; /* Acceptance criteria < 1) */
70: PetscReal minimumStep; /* Minimum step size */
71: PetscReal lastReference; /* Reference value of last iteration */
73: PetscInt memorySize; /* Number of functions kept in memory */
74: PetscInt current; /* Current element for FIFO */
75: PetscInt referencePolicy; /* Integer for reference calculation rule */
76: PetscInt replacementPolicy; /* Policy for replacing values in memory */
78: PetscBool nondescending;
79: PetscBool memorySetup;
81: Vec x; /* Maintain reference to variable vector to check for changes */
82: Vec work;
83: } TaoLineSearch_OWARMIJO;
85: static PetscErrorCode ProjWork_OWLQN(Vec w, Vec x, Vec gv, PetscReal *gdx);