alpaqa pantr
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 true 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, false to leave out that term.
40 bool hessian_vec = true;
41};
42
43/// @ingroup grp_DirectionProviders
44template <Config Conf = DefaultConfig>
50
51 struct Params {
54 };
55
58 : reg_params(params.accelerator), direction_params(params.direction) {}
60 const DirectionParams &directionparams = {})
61 : reg_params(params), direction_params(directionparams) {}
62
63 /// @see @ref PANOCDirection::initialize
65 [[maybe_unused]] real_t γ_0, [[maybe_unused]] crvec x_0,
66 [[maybe_unused]] crvec x̂_0, [[maybe_unused]] crvec p_0,
67 [[maybe_unused]] crvec grad_ψx_0) {
69 throw std::invalid_argument(
70 "Structured Newton only supports box-constrained problems");
71 // TODO: support eval_inactive_indices_res_lna
73 throw std::invalid_argument("Structured Newton requires hess_ψ");
74 // Store references to problem and ALM variables
75 this->problem = &problem;
76 this->y.emplace(y);
77 this->Σ.emplace(Σ);
78 // Allocate workspaces
79 const auto n = problem.get_n();
80 JK.resize(n);
81 H_storage.resize(n * n);
82 HJ_storage.resize(n * n);
83 // Store sparsity of H
85 if (nnz_H > 0) {
86 inner_idx_H.resize(nnz_H);
87 outer_ptr_H.resize(n + 1);
88 mvec null{nullptr, 0};
90 throw std::logic_error("Sparse hessians not yet implemented");
91 }
92 }
93
94 /// @see @ref PANOCDirection::has_initial_direction
95 bool has_initial_direction() const { return true; }
96
97 /// @see @ref PANOCDirection::update
98 bool update([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t γₙₑₓₜ,
99 [[maybe_unused]] crvec xₖ, [[maybe_unused]] crvec xₙₑₓₜ,
100 [[maybe_unused]] crvec pₖ, [[maybe_unused]] crvec pₙₑₓₜ,
101 [[maybe_unused]] crvec grad_ψxₖ,
102 [[maybe_unused]] crvec grad_ψxₙₑₓₜ) {
103 return true;
104 }
105
106 /// @see @ref PANOCDirection::apply
107 bool apply(real_t γₖ, crvec xₖ, [[maybe_unused]] crvec x̂ₖ, crvec pₖ,
108 crvec grad_ψxₖ, rvec qₖ) const {
109
110 const auto n = problem->get_n();
111 const auto &C = problem->get_box_C();
112
113 // Find inactive indices J
114 auto nJ = 0;
115 for (index_t i = 0; i < n; ++i) {
116 real_t gd = xₖ(i) - γₖ * grad_ψxₖ(i);
117 if (gd <= C.lowerbound(i)) { // i ∊ J̲ ⊆ K
118 qₖ(i) = pₖ(i); //
119 } else if (C.upperbound(i) <= gd) { // i ∊ J̅ ⊆ K
120 qₖ(i) = pₖ(i); //
121 } else { // i ∊ J
122 JK(nJ++) = i;
123 qₖ(i) = -grad_ψxₖ(i);
124 }
125 }
126
127 // There are no inactive indices J
128 if (nJ == 0) {
129 // No free variables, no Newton step possible
130 return false; // Simply use the projection step
131 }
132
133 // Compute the Hessian
134 mmat H{H_storage.data(), n, n};
136 H_storage);
137
138 // There are no active indices K
139 if (nJ == n) {
140 // If all indices are free, we can factor the entire matrix.
141 // Find the minimum eigenvalue to regularize the Hessian matrix and
142 // make it positive definite.
144 // Find the minimum eigenvalue to regularize the Hessian matrix and
145 // make it positive definite.
146 Eigen::SelfAdjointEigenSolver<mat> eig{H,
147 Eigen::ComputeEigenvectors};
148
149 auto λ_min = eig.eigenvalues().minCoeff(),
150 λ_max = eig.eigenvalues().maxCoeff();
151
153 std::cout << "λ(H): " << float_to_str(λ_min, 3) << ", "
154 << float_to_str(λ_max, 3) << std::endl;
155 // Regularization
156 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
157 // Solve the system
158 qₖ = eig.eigenvectors().transpose() * qₖ;
159 qₖ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qₖ;
160 qₖ = eig.eigenvectors() * qₖ;
161 return true;
162 }
163
164 // There are active indices K
165 crindexvec J = JK.topRows(nJ);
166 rindexvec K = JK.bottomRows(n - nJ);
168
169 // Compute right-hand side of 6.1c
171 qₖ(J).noalias() -= H(J, K) * qₖ(K);
172
173 // If there are active indices, we need to solve the Newton system with
174 // just the inactive indices.
175 mmat HJ{HJ_storage.data(), nJ, nJ};
176 // We copy the inactive block of the Hessian to a temporary dense matrix.
177 // Since it's symmetric, only the lower part is copied.
178 HJ.template triangularView<Eigen::Lower>() =
179 H(J, J).template triangularView<Eigen::Lower>();
180
182 // Find the minimum eigenvalue to regularize the Hessian matrix and
183 // make it positive definite.
184 Eigen::SelfAdjointEigenSolver<mat> eig{HJ, Eigen::ComputeEigenvectors};
185
186 auto λ_min = eig.eigenvalues().minCoeff(),
187 λ_max = eig.eigenvalues().maxCoeff();
188
190 std::cout << "λ(H_JJ): " << float_to_str(λ_min, 3) << ", "
191 << float_to_str(λ_max, 3) << std::endl;
192 // Regularization
193 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
194 // Solve the system
195 auto qJ = H.col(0).topRows(nJ);
196 qJ = qₖ(J);
197 qJ = eig.eigenvectors().transpose() * qJ;
198 qJ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qJ;
199 qₖ(J) = eig.eigenvectors() * qJ;
200 return true;
201 }
202
203 /// @see @ref PANOCDirection::changed_γ
204 void changed_γ([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t old_γₖ) {
205 }
206
207 /// @see @ref PANOCDirection::reset
208 void reset() {}
209
210 /// @see @ref PANOCDirection::get_name
211 std::string get_name() const {
212 return "StructuredNewtonDirection<" +
213 std::string(config_t::get_name()) + '>';
214 }
215
216 const auto &get_params() const { return direction_params; }
217
218 private:
219 const Problem *problem = nullptr;
220#ifndef _WIN32
221 std::optional<crvec> y = std::nullopt;
222 std::optional<crvec> Σ = std::nullopt;
223#else
224 std::optional<vec> y = std::nullopt;
225 std::optional<vec> Σ = std::nullopt;
226#endif
227
228 mutable indexvec JK;
229 mutable vec H_storage;
232
233 public:
236};
237
242#ifdef ALPAQA_WITH_QUAD_PRECISION
244#endif
245
246} // namespace alpaqa
247
248#include <alpaqa/inner/panoc.hpp>
249
250namespace alpaqa {
251
252// clang-format off
253ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<DefaultConfig>);
254ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigf>);
255ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigd>);
256ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigl>);
257#ifdef ALPAQA_WITH_QUAD_PRECISION
258ALPAQA_EXPORT_EXTERN_TEMPLATE(class, PANOCSolver, StructuredNewtonDirection<EigenConfigq>);
259#endif
260// clang-format on
261
262} // 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:42
#define ALPAQA_EXPORT_EXTERN_TEMPLATE(...)
Definition: export.hpp:21
typename Conf::mvec mvec
Definition: config.hpp:53
bool hessian_vec
Set this option to true to include the Hessian-vector product from equation 12b in ,...
typename Conf::indexvec indexvec
Definition: config.hpp:64
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:51
typename Conf::rindexvec rindexvec
Definition: config.hpp:65
typename Conf::index_t index_t
Definition: config.hpp:63
typename Conf::mmat mmat
Definition: config.hpp:58
typename Conf::length_t length_t
Definition: config.hpp:62
typename Conf::rvec rvec
Definition: config.hpp:55
typename Conf::crvec crvec
Definition: config.hpp:56
typename Conf::vec vec
Definition: config.hpp:52
bool print_eig
Print the minimum and maximum eigenvalue of the Hessian.
std::string float_to_str(F value, int precision)
Definition: print.tpp:66
typename Conf::crindexvec crindexvec
Definition: config.hpp:66
Parameters for the StructuredNewtonDirection class.
Parameters for the StructuredNewtonDirection class.
Double-precision double configuration.
Definition: config.hpp:115
Single-precision float configuration.
Definition: config.hpp:111
long double configuration.
Definition: config.hpp:120
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