Nonconvex constrained optimization
Loading...
Searching...
No Matches
structured-newton.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <alpaqa/export.hpp>
11#include <alpaqa/util/print.hpp>
12#include <iostream>
13#include <limits>
14#include <optional>
15#include <stdexcept>
16
17#include <Eigen/Cholesky>
18#include <Eigen/Eigenvalues>
19
20namespace alpaqa {
21
22/// Parameters for the @ref StructuredNewtonDirection class.
23/// @ingroup grp_Parameters
24template <Config Conf>
27 /// Minimum eigenvalue of the Hessian, scaled by
28 /// @f$ 1 + |\lambda_\mathrm{max}| @f$, enforced by regularization using
29 /// a multiple of identity.
30 real_t min_eig = std::cbrt(std::numeric_limits<real_t>::epsilon());
31 /// Print the minimum and maximum eigenvalue of the Hessian.
32 bool print_eig = false;
33};
34
35/// Parameters for the @ref StructuredNewtonDirection class.
36/// @ingroup grp_Parameters
37template <Config Conf>
40 /// Set this option to a nonzero value to include the Hessian-vector product
41 /// @f$ \nabla^2_{x_\mathcal{J}x_\mathcal{K}}\psi(x) q_\mathcal{K} @f$ from
42 /// equation 12b in @cite pas2022alpaqa, scaled by this parameter.
43 /// Set it to zero to leave out that term.
45};
46
47/// @ingroup grp_DirectionProviders
48template <Config Conf = DefaultConfig>
54
59
64 const DirectionParams &directionparams = {})
65 : reg_params(params), direction_params(directionparams) {}
66
67 /// @see @ref PANOCDirection::initialize
69 [[maybe_unused]] real_t γ_0, [[maybe_unused]] crvec x_0,
70 [[maybe_unused]] crvec x̂_0, [[maybe_unused]] crvec p_0,
71 [[maybe_unused]] crvec grad_ψx_0) {
72 if (!(problem.provides_get_variable_bounds() &&
73 problem.provides_get_general_bounds()))
74 throw std::invalid_argument(
75 "Structured Newton only supports box-constrained problems");
76 // TODO: support eval_inactive_indices_res_lna
77 if (!problem.supports_eval_augmented_lagrangian_hessian())
78 throw std::invalid_argument("Structured Newton requires hess_ψ");
79 // Store references to problem and ALM variables
80 this->problem = &problem;
81 this->y.emplace(y);
82 this->Σ.emplace(Σ);
83 // Allocate workspaces
84 const auto n = problem.get_num_variables();
85 JK.resize(n);
86 H.resize(n, n);
87 HJ_storage.resize(n * n);
88 if (!is_dense(problem.get_augmented_lagrangian_hessian_sparsity()))
89 throw std::logic_error("Sparse hessians not yet implemented");
90 }
91
92 /// @see @ref PANOCDirection::has_initial_direction
93 bool has_initial_direction() const { return true; }
94
95 /// @see @ref PANOCDirection::update
96 bool update([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t γₙₑₓₜ,
97 [[maybe_unused]] crvec xₖ, [[maybe_unused]] crvec xₙₑₓₜ,
98 [[maybe_unused]] crvec pₖ, [[maybe_unused]] crvec pₙₑₓₜ,
99 [[maybe_unused]] crvec grad_ψxₖ,
100 [[maybe_unused]] crvec grad_ψxₙₑₓₜ) {
101 return true;
102 }
103
104 /// @see @ref PANOCDirection::apply
105 bool apply(real_t γₖ, crvec xₖ, [[maybe_unused]] crvec x̂ₖ, crvec pₖ,
106 crvec grad_ψxₖ, rvec qₖ) const {
107
108 const auto n = problem->get_num_variables();
109 const auto &C = problem->get_variable_bounds();
110
111 // Find inactive indices J
112 auto nJ = 0;
113 for (index_t i = 0; i < n; ++i) {
114 real_t gd = xₖ(i) - γₖ * grad_ψxₖ(i);
115 if (gd <= C.lower(i)) { // i ∊ J̲ ⊆ K
116 qₖ(i) = pₖ(i); //
117 } else if (C.upper(i) <= gd) { // i ∊ J̅ ⊆ K
118 qₖ(i) = pₖ(i); //
119 } else { // i ∊ J
120 JK(nJ++) = i;
121 qₖ(i) = -grad_ψxₖ(i);
122 }
123 }
124
125 // There are no inactive indices J
126 if (nJ == 0) {
127 // No free variables, no Newton step possible
128 return false; // Simply use the projection step
129 }
130
131 // Compute the Hessian
132 problem->eval_augmented_lagrangian_hessian(xₖ, *y, *Σ, 1, H.reshaped());
133
134 // There are no active indices K
135 if (nJ == n) {
136 // If all indices are free, we can factor the entire matrix.
137 // Find the minimum eigenvalue to regularize the Hessian matrix and
138 // make it positive definite.
140 // Find the minimum eigenvalue to regularize the Hessian matrix and
141 // make it positive definite.
142 Eigen::SelfAdjointEigenSolver<mat> eig{H,
143 Eigen::ComputeEigenvectors};
144
145 auto λ_min = eig.eigenvalues().minCoeff(),
146 λ_max = eig.eigenvalues().maxCoeff();
147
148 if (reg_params.print_eig)
149 std::cout << "λ(H): " << float_to_str(λ_min, 3) << ", "
150 << float_to_str(λ_max, 3) << std::endl;
151 // Regularization
152 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
153 // Solve the system
154 qₖ = eig.eigenvectors().transpose() * qₖ;
155 qₖ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qₖ;
156 qₖ = eig.eigenvectors() * qₖ;
157 return true;
158 }
159
160 // There are active indices K
161 crindexvec J = JK.topRows(nJ);
162 rindexvec K = JK.bottomRows(n - nJ);
164
165 // Compute right-hand side of 6.1c
166 if (direction_params.hessian_vec_factor != 0)
167 qₖ(J).noalias() -=
168 direction_params.hessian_vec_factor * (H(J, K) * qₖ(K));
169
170 // If there are active indices, we need to solve the Newton system with
171 // just the inactive indices.
172 mmat HJ{HJ_storage.data(), nJ, nJ};
173 // We copy the inactive block of the Hessian to a temporary dense matrix.
174 // Since it's symmetric, only the lower part is copied.
175 HJ.template triangularView<Eigen::Lower>() =
176 H(J, J).template triangularView<Eigen::Lower>();
177
179 // Find the minimum eigenvalue to regularize the Hessian matrix and
180 // make it positive definite.
181 Eigen::SelfAdjointEigenSolver<mat> eig{HJ, Eigen::ComputeEigenvectors};
182
183 auto λ_min = eig.eigenvalues().minCoeff(),
184 λ_max = eig.eigenvalues().maxCoeff();
185
186 if (reg_params.print_eig)
187 std::cout << "λ(H_JJ): " << float_to_str(λ_min, 3) << ", "
188 << float_to_str(λ_max, 3) << std::endl;
189 // Regularization
190 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
191 // Solve the system
192 auto qJ = H.col(0).topRows(nJ);
193 qJ = qₖ(J);
194 qJ = eig.eigenvectors().transpose() * qJ;
195 qJ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qJ;
196 qₖ(J) = eig.eigenvectors() * qJ;
197 return true;
198 }
199
200 /// @see @ref PANOCDirection::changed_γ
201 void changed_γ([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t old_γₖ) {
202 }
203
204 /// @see @ref PANOCDirection::reset
205 void reset() {}
206
207 /// @see @ref PANOCDirection::get_name
208 std::string get_name() const {
209 return "StructuredNewtonDirection<" +
210 std::string(config_t::get_name()) + '>';
211 }
212
213 const auto &get_params() const { return direction_params; }
214
215 private:
216 const Problem *problem = nullptr;
217#ifndef _WIN32
218 std::optional<crvec> y = std::nullopt;
219 std::optional<crvec> Σ = std::nullopt;
220#else
221 std::optional<vec> y = std::nullopt;
222 std::optional<vec> Σ = std::nullopt;
223#endif
224
225 mutable indexvec JK;
226 mutable mat H;
228
229 public:
232};
233
234// clang-format off
239// clang-format on
240
241} // namespace alpaqa
The main polymorphic minimization problem interface.
#define USING_ALPAQA_CONFIG(Conf)
Definition config.hpp:77
#define ALPAQA_IF_QUADF(...)
Definition config.hpp:223
#define ALPAQA_IF_LONGD(...)
Definition config.hpp:235
#define ALPAQA_IF_FLOAT(...)
Definition config.hpp:229
#define ALPAQA_EXPORT_EXTERN_TEMPLATE(...)
Definition export.hpp:23
Parameters for the StructuredNewtonDirection class.
Parameters for the StructuredNewtonDirection class.
typename Conf::mat mat
Definition config.hpp:93
typename Conf::indexvec indexvec
Definition config.hpp:105
typename Conf::real_t real_t
Definition config.hpp:86
typename Conf::rindexvec rindexvec
Definition config.hpp:106
typename Conf::index_t index_t
Definition config.hpp:104
typename Conf::mmat mmat
Definition config.hpp:94
typename Conf::rvec rvec
Definition config.hpp:91
typename Conf::crvec crvec
Definition config.hpp:92
typename Conf::vec vec
Definition config.hpp:88
typename Conf::crindexvec crindexvec
Definition config.hpp:107
Double-precision double configuration.
Definition config.hpp:176
Single-precision float configuration.
Definition config.hpp:172
long double configuration.
Definition config.hpp:181
void changed_γ(real_t γₖ, real_t old_γₖ)
StructuredNewtonDirection(const AcceleratorParams &params, const DirectionParams &directionparams={})
bool update(real_t γₖ, real_t γₙₑₓₜ, crvec xₖ, crvec xₙₑₓₜ, crvec pₖ, crvec pₙₑₓₜ, crvec grad_ψxₖ, crvec grad_ψxₙₑₓₜ)
StructuredNewtonDirectionParams< config_t > DirectionParams
StructuredNewtonDirection(const Params &params)
StructuredNewtonRegularizationParams< config_t > AcceleratorParams
void initialize(const Problem &problem, crvec y, crvec Σ, real_t γ_0, crvec x_0, crvec x̂_0, crvec p_0, crvec grad_ψx_0)
bool apply(real_t γₖ, crvec xₖ, crvec x̂ₖ, crvec pₖ, crvec grad_ψxₖ, rvec qₖ) const
TypeErasedProblem< config_t > Problem
static void compute_complement(std::span< const index_t > in, std::span< index_t > out)
Definition index-set.hpp:37