alpaqa 1.0.0a11
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>
10#include <alpaqa/util/print.hpp>
11#include <iostream>
12#include <limits>
13#include <optional>
14#include <stdexcept>
15
16#include <Eigen/Cholesky>
17#include <Eigen/Eigenvalues>
18
19namespace alpaqa {
20
21/// Parameters for the @ref StructuredNewtonDirection class.
22template <Config Conf>
25 /// Minimum eigenvalue of the Hessian, scaled by
26 /// @f$ 1 + |\lambda_\mathrm{max}| @f$, enforced by regularization using
27 /// a multiple of identity.
28 real_t min_eig = std::cbrt(std::numeric_limits<real_t>::epsilon());
29 /// Print the minimum and maximum eigenvalue of the Hessian.
30 bool print_eig = false;
31};
32
33/// Parameters for the @ref StructuredNewtonDirection class.
34template <Config Conf>
37 /// Set this option to a nonzero value to include the Hessian-vector product
38 /// @f$ \nabla^2_{x_\mathcal{J}x_\mathcal{K}}\psi(x) q_\mathcal{K} @f$ from
39 /// equation 12b in @cite pas2022alpaqa, scaled by this parameter.
40 /// Set it to zero to leave out that term.
42};
43
44/// @ingroup grp_DirectionProviders
45template <Config Conf = DefaultConfig>
51
52 struct Params {
55 };
56
59 : reg_params(params.accelerator), direction_params(params.direction) {}
61 const DirectionParams &directionparams = {})
62 : reg_params(params), direction_params(directionparams) {}
63
64 /// @see @ref PANOCDirection::initialize
66 [[maybe_unused]] real_t γ_0, [[maybe_unused]] crvec x_0,
67 [[maybe_unused]] crvec x̂_0, [[maybe_unused]] crvec p_0,
68 [[maybe_unused]] crvec grad_ψx_0) {
70 throw std::invalid_argument(
71 "Structured Newton only supports box-constrained problems");
72 // TODO: support eval_inactive_indices_res_lna
74 throw std::invalid_argument("Structured Newton requires hess_ψ");
75 // Store references to problem and ALM variables
76 this->problem = &problem;
77 this->y.emplace(y);
78 this->Σ.emplace(Σ);
79 // Allocate workspaces
80 const auto n = problem.get_n();
81 JK.resize(n);
82 H_storage.resize(n * n);
83 HJ_storage.resize(n * n);
84 // Store sparsity of H
86 if (nnz_H > 0) {
87 inner_idx_H.resize(nnz_H);
88 outer_ptr_H.resize(n + 1);
89 mvec null{nullptr, 0};
91 throw std::logic_error("Sparse hessians not yet implemented");
92 }
93 }
94
95 /// @see @ref PANOCDirection::has_initial_direction
96 bool has_initial_direction() const { return true; }
97
98 /// @see @ref PANOCDirection::update
99 bool update([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t γₙₑₓₜ,
100 [[maybe_unused]] crvec xₖ, [[maybe_unused]] crvec xₙₑₓₜ,
101 [[maybe_unused]] crvec pₖ, [[maybe_unused]] crvec pₙₑₓₜ,
102 [[maybe_unused]] crvec grad_ψxₖ,
103 [[maybe_unused]] crvec grad_ψxₙₑₓₜ) {
104 return true;
105 }
106
107 /// @see @ref PANOCDirection::apply
108 bool apply(real_t γₖ, crvec xₖ, [[maybe_unused]] crvec x̂ₖ, crvec pₖ,
109 crvec grad_ψxₖ, rvec qₖ) const {
110
111 const auto n = problem->get_n();
112 const auto &C = problem->get_box_C();
113
114 // Find inactive indices J
115 auto nJ = 0;
116 for (index_t i = 0; i < n; ++i) {
117 real_t gd = xₖ(i) - γₖ * grad_ψxₖ(i);
118 if (gd <= C.lowerbound(i)) { // i ∊ J̲ ⊆ K
119 qₖ(i) = pₖ(i); //
120 } else if (C.upperbound(i) <= gd) { // i ∊ J̅ ⊆ K
121 qₖ(i) = pₖ(i); //
122 } else { // i ∊ J
123 JK(nJ++) = i;
124 qₖ(i) = -grad_ψxₖ(i);
125 }
126 }
127
128 // There are no inactive indices J
129 if (nJ == 0) {
130 // No free variables, no Newton step possible
131 return false; // Simply use the projection step
132 }
133
134 // Compute the Hessian
135 mmat H{H_storage.data(), n, n};
137 H_storage);
138
139 // There are no active indices K
140 if (nJ == n) {
141 // If all indices are free, we can factor the entire matrix.
142 // Find the minimum eigenvalue to regularize the Hessian matrix and
143 // make it positive definite.
145 // Find the minimum eigenvalue to regularize the Hessian matrix and
146 // make it positive definite.
147 Eigen::SelfAdjointEigenSolver<mat> eig{H,
148 Eigen::ComputeEigenvectors};
149
150 auto λ_min = eig.eigenvalues().minCoeff(),
151 λ_max = eig.eigenvalues().maxCoeff();
152
154 std::cout << "λ(H): " << float_to_str(λ_min, 3) << ", "
155 << float_to_str(λ_max, 3) << std::endl;
156 // Regularization
157 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
158 // Solve the system
159 qₖ = eig.eigenvectors().transpose() * qₖ;
160 qₖ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qₖ;
161 qₖ = eig.eigenvectors() * qₖ;
162 return true;
163 }
164
165 // There are active indices K
166 crindexvec J = JK.topRows(nJ);
167 rindexvec K = JK.bottomRows(n - nJ);
169
170 // Compute right-hand side of 6.1c
172 qₖ(J).noalias() -=
173 direction_params.hessian_vec_factor * (H(J, K) * qₖ(K));
174
175 // If there are active indices, we need to solve the Newton system with
176 // just the inactive indices.
177 mmat HJ{HJ_storage.data(), nJ, nJ};
178 // We copy the inactive block of the Hessian to a temporary dense matrix.
179 // Since it's symmetric, only the lower part is copied.
180 HJ.template triangularView<Eigen::Lower>() =
181 H(J, J).template triangularView<Eigen::Lower>();
182
184 // Find the minimum eigenvalue to regularize the Hessian matrix and
185 // make it positive definite.
186 Eigen::SelfAdjointEigenSolver<mat> eig{HJ, Eigen::ComputeEigenvectors};
187
188 auto λ_min = eig.eigenvalues().minCoeff(),
189 λ_max = eig.eigenvalues().maxCoeff();
190
192 std::cout << "λ(H_JJ): " << float_to_str(λ_min, 3) << ", "
193 << float_to_str(λ_max, 3) << std::endl;
194 // Regularization
195 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
196 // Solve the system
197 auto qJ = H.col(0).topRows(nJ);
198 qJ = qₖ(J);
199 qJ = eig.eigenvectors().transpose() * qJ;
200 qJ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qJ;
201 qₖ(J) = eig.eigenvectors() * qJ;
202 return true;
203 }
204
205 /// @see @ref PANOCDirection::changed_γ
206 void changed_γ([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t old_γₖ) {
207 }
208
209 /// @see @ref PANOCDirection::reset
210 void reset() {}
211
212 /// @see @ref PANOCDirection::get_name
213 std::string get_name() const {
214 return "StructuredNewtonDirection<" +
215 std::string(config_t::get_name()) + '>';
216 }
217
218 const auto &get_params() const { return direction_params; }
219
220 private:
221 const Problem *problem = nullptr;
222#ifndef _WIN32
223 std::optional<crvec> y = std::nullopt;
224 std::optional<crvec> Σ = std::nullopt;
225#else
226 std::optional<vec> y = std::nullopt;
227 std::optional<vec> Σ = std::nullopt;
228#endif
229
230 mutable indexvec JK;
231 mutable vec H_storage;
234
235 public:
238};
239
240// clang-format off
245// clang-format on
246
247} // namespace alpaqa
248
249#include <alpaqa/inner/panoc.hpp>
250
251namespace alpaqa {
252
253// clang-format off
254ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigd>);
255ALPAQA_IF_FLOAT(ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigf>);)
256ALPAQA_IF_LONGD(ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigl>);)
257ALPAQA_IF_QUADF(ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigq>);)
258// clang-format on
259
260} // namespace alpaqa
void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rindexvec inner_idx, rindexvec outer_ptr, rvec H_values) const
[Optional] Function that evaluates the Hessian of the augmented Lagrangian,
bool provides_get_box_C() const
Returns true if the problem provides an implementation of get_box_C.
length_t get_n() const
[Required] Number of decision variables.
bool provides_eval_hess_ψ() const
Returns true if the problem provides an implementation of eval_hess_ψ.
bool provides_get_box_D() const
Returns true if the problem provides an implementation of get_box_D.
const Box & get_box_C() const
[Optional] Get the rectangular constraint set of the decision variables, .
length_t get_hess_ψ_num_nonzeros() const
[Optional] Function that gets the number of nonzeros of the Hessian of the augmented Lagrangian.
#define USING_ALPAQA_CONFIG(Conf)
Definition: config.hpp:54
#define ALPAQA_IF_QUADF(...)
Definition: config.hpp:171
#define ALPAQA_IF_LONGD(...)
Definition: config.hpp:183
#define ALPAQA_IF_FLOAT(...)
Definition: config.hpp:177
#define ALPAQA_EXPORT_EXTERN_TEMPLATE(...)
Definition: export.hpp:21
typename Conf::mvec mvec
Definition: config.hpp:65
real_t hessian_vec_factor
Set this option to a nonzero value to include the Hessian-vector product from equation 12b in ,...
typename Conf::indexvec indexvec
Definition: config.hpp:76
real_t min_eig
Minimum eigenvalue of the Hessian, scaled by , enforced by regularization using a multiple of identit...
typename Conf::real_t real_t
Definition: config.hpp:63
typename Conf::rindexvec rindexvec
Definition: config.hpp:77
typename Conf::index_t index_t
Definition: config.hpp:75
typename Conf::mmat mmat
Definition: config.hpp:70
typename Conf::length_t length_t
Definition: config.hpp:74
typename Conf::rvec rvec
Definition: config.hpp:67
typename Conf::crvec crvec
Definition: config.hpp:68
typename Conf::vec vec
Definition: config.hpp:64
bool print_eig
Print the minimum and maximum eigenvalue of the Hessian.
std::string float_to_str(F value, int precision)
Definition: print.tpp:67
typename Conf::crindexvec crindexvec
Definition: config.hpp:78
Parameters for the StructuredNewtonDirection class.
Parameters for the StructuredNewtonDirection class.
Double-precision double configuration.
Definition: config.hpp:127
Single-precision float configuration.
Definition: config.hpp:123
long double configuration.
Definition: config.hpp:132
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ₙₑₓₜ)
StructuredNewtonDirection(const Params &params)
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
static void compute_complement(std::span< const index_t > in, std::span< index_t > out)
Definition: index-set.hpp:36