alpaqa 1.0.0a12
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.
23template <Config Conf>
26 /// Minimum eigenvalue of the Hessian, scaled by
27 /// @f$ 1 + |\lambda_\mathrm{max}| @f$, enforced by regularization using
28 /// a multiple of identity.
29 real_t min_eig = std::cbrt(std::numeric_limits<real_t>::epsilon());
30 /// Print the minimum and maximum eigenvalue of the Hessian.
31 bool print_eig = false;
32};
33
34/// Parameters for the @ref StructuredNewtonDirection class.
35template <Config Conf>
38 /// Set this option to a nonzero value to include the Hessian-vector product
39 /// @f$ \nabla^2_{x_\mathcal{J}x_\mathcal{K}}\psi(x) q_\mathcal{K} @f$ from
40 /// equation 12b in @cite pas2022alpaqa, scaled by this parameter.
41 /// Set it to zero to leave out that term.
43};
44
45/// @ingroup grp_DirectionProviders
46template <Config Conf = DefaultConfig>
52
57
60 : reg_params(params.accelerator), direction_params(params.direction) {}
64
65 /// @see @ref PANOCDirection::initialize
71 throw std::invalid_argument(
72 "Structured Newton only supports box-constrained problems");
73 // TODO: support eval_inactive_indices_res_lna
75 throw std::invalid_argument("Structured Newton requires hess_ψ");
76 // Store references to problem and ALM variables
77 this->problem = &problem;
78 this->y.emplace(y);
79 this->Σ.emplace(Σ);
80 // Allocate workspaces
81 const auto n = problem.get_n();
82 JK.resize(n);
83 H.resize(n, n);
84 HJ_storage.resize(n * n);
85 if (!is_dense(problem.get_hess_ψ_sparsity()))
86 throw std::logic_error("Sparse hessians not yet implemented");
87 }
88
89 /// @see @ref PANOCDirection::has_initial_direction
90 bool has_initial_direction() const { return true; }
91
92 /// @see @ref PANOCDirection::update
100
101 /// @see @ref PANOCDirection::apply
103 crvec grad_ψxₖ, rvec qₖ) const {
104
105 const auto n = problem->get_n();
106 const auto &C = problem->get_box_C();
107
108 // Find inactive indices J
109 auto nJ = 0;
110 for (index_t i = 0; i < n; ++i) {
111 real_t gd = xₖ(i) - γₖ * grad_ψxₖ(i);
112 if (gd <= C.lowerbound(i)) { // i ∊ J̲ ⊆ K
113 qₖ(i) = pₖ(i); //
114 } else if (C.upperbound(i) <= gd) { // i ∊ J̅ ⊆ K
115 qₖ(i) = pₖ(i); //
116 } else { // i ∊ J
117 JK(nJ++) = i;
118 qₖ(i) = -grad_ψxₖ(i);
119 }
120 }
121
122 // There are no inactive indices J
123 if (nJ == 0) {
124 // No free variables, no Newton step possible
125 return false; // Simply use the projection step
126 }
127
128 // Compute the Hessian
129 problem->eval_hess_ψ(xₖ, *y, *Σ, 1, H.reshaped());
130
131 // There are no active indices K
132 if (nJ == n) {
133 // If all indices are free, we can factor the entire matrix.
134 // Find the minimum eigenvalue to regularize the Hessian matrix and
135 // make it positive definite.
137 // Find the minimum eigenvalue to regularize the Hessian matrix and
138 // make it positive definite.
139 Eigen::SelfAdjointEigenSolver<mat> eig{H,
140 Eigen::ComputeEigenvectors};
141
142 auto λ_min = eig.eigenvalues().minCoeff(),
143 λ_max = eig.eigenvalues().maxCoeff();
144
146 std::cout << "λ(H): " << float_to_str(λ_min, 3) << ", "
147 << float_to_str(λ_max, 3) << std::endl;
148 // Regularization
149 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
150 // Solve the system
151 qₖ = eig.eigenvectors().transpose() * qₖ;
152 qₖ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qₖ;
153 qₖ = eig.eigenvectors() * qₖ;
154 return true;
155 }
156
157 // There are active indices K
158 crindexvec J = JK.topRows(nJ);
159 rindexvec K = JK.bottomRows(n - nJ);
161
162 // Compute right-hand side of 6.1c
164 qₖ(J).noalias() -=
166
167 // If there are active indices, we need to solve the Newton system with
168 // just the inactive indices.
169 mmat HJ{HJ_storage.data(), nJ, nJ};
170 // We copy the inactive block of the Hessian to a temporary dense matrix.
171 // Since it's symmetric, only the lower part is copied.
173 H(J, J).template triangularView<Eigen::Lower>();
174
176 // Find the minimum eigenvalue to regularize the Hessian matrix and
177 // make it positive definite.
178 Eigen::SelfAdjointEigenSolver<mat> eig{HJ, Eigen::ComputeEigenvectors};
179
180 auto λ_min = eig.eigenvalues().minCoeff(),
181 λ_max = eig.eigenvalues().maxCoeff();
182
184 std::cout << "λ(H_JJ): " << float_to_str(λ_min, 3) << ", "
185 << float_to_str(λ_max, 3) << std::endl;
186 // Regularization
187 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
188 // Solve the system
189 auto qJ = H.col(0).topRows(nJ);
190 qJ = qₖ(J);
191 qJ = eig.eigenvectors().transpose() * qJ;
192 qJ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qJ;
193 qₖ(J) = eig.eigenvectors() * qJ;
194 return true;
195 }
196
197 /// @see @ref PANOCDirection::changed_γ
200
201 /// @see @ref PANOCDirection::reset
202 void reset() {}
203
204 /// @see @ref PANOCDirection::get_name
205 std::string get_name() const {
206 return "StructuredNewtonDirection<" +
207 std::string(config_t::get_name()) + '>';
208 }
209
210 const auto &get_params() const { return direction_params; }
211
212 private:
213 const Problem *problem = nullptr;
214#ifndef _WIN32
215 std::optional<crvec> y = std::nullopt;
216 std::optional<crvec> Σ = std::nullopt;
217#else
218 std::optional<vec> y = std::nullopt;
219 std::optional<vec> Σ = std::nullopt;
220#endif
221
222 mutable indexvec JK;
223 mutable mat H;
225
226 public:
229};
230
231// clang-format off
236// clang-format on
237
238} // namespace alpaqa
239
240#include <alpaqa/inner/panoc.hpp>
241
242namespace alpaqa {
243
244// clang-format off
249// clang-format on
250
251} // namespace alpaqa
bool provides_get_box_C() const
Returns true if the problem provides an implementation of get_box_C.
Sparsity get_hess_ψ_sparsity() const
[Optional] Function that returns (a view of) the sparsity pattern of the Hessian of the augmented Lag...
length_t get_n() const
[Required] Number of decision variables.
void eval_hess_ψ(crvec x, crvec y, crvec Σ, real_t scale, rvec H_values) const
[Optional] Function that evaluates the nonzero values of the Hessian of the augmented Lagrangian,
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, .
#define USING_ALPAQA_CONFIG(Conf)
Definition config.hpp:56
#define ALPAQA_IF_QUADF(...)
Definition config.hpp:182
#define ALPAQA_IF_LONGD(...)
Definition config.hpp:194
#define ALPAQA_IF_FLOAT(...)
Definition config.hpp:188
#define ALPAQA_EXPORT_EXTERN_TEMPLATE(...)
Definition export.hpp:21
real_t hessian_vec_factor
Set this option to a nonzero value to include the Hessian-vector product from equation 12b in ,...
typename Conf::mat mat
Definition config.hpp:71
typename Conf::indexvec indexvec
Definition config.hpp:78
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:65
typename Conf::rindexvec rindexvec
Definition config.hpp:79
typename Conf::index_t index_t
Definition config.hpp:77
typename Conf::mmat mmat
Definition config.hpp:72
constexpr const auto inf
Definition config.hpp:85
typename Conf::rvec rvec
Definition config.hpp:69
typename Conf::crvec crvec
Definition config.hpp:70
typename Conf::vec vec
Definition config.hpp:66
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:80
Parameters for the StructuredNewtonDirection class.
Parameters for the StructuredNewtonDirection class.
Double-precision double configuration.
Definition config.hpp:135
Single-precision float configuration.
Definition config.hpp:131
long double configuration.
Definition config.hpp:140
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