alpaqa 1.0.0a16
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
62 : reg_params(params.accelerator), direction_params(params.direction) {}
66
67 /// @see @ref PANOCDirection::initialize
73 throw std::invalid_argument(
74 "Structured Newton only supports box-constrained problems");
75 // TODO: support eval_inactive_indices_res_lna
77 throw std::invalid_argument("Structured Newton requires hess_ψ");
78 // Store references to problem and ALM variables
79 this->problem = &problem;
80 this->y.emplace(y);
81 this->Σ.emplace(Σ);
82 // Allocate workspaces
83 const auto n = problem.get_n();
84 JK.resize(n);
85 H.resize(n, n);
86 HJ_storage.resize(n * n);
87 if (!is_dense(problem.get_hess_ψ_sparsity()))
88 throw std::logic_error("Sparse hessians not yet implemented");
89 }
90
91 /// @see @ref PANOCDirection::has_initial_direction
92 bool has_initial_direction() const { return true; }
93
94 /// @see @ref PANOCDirection::update
102
103 /// @see @ref PANOCDirection::apply
105 crvec grad_ψxₖ, rvec qₖ) const {
106
107 const auto n = problem->get_n();
108 const auto &C = problem->get_box_C();
109
110 // Find inactive indices J
111 auto nJ = 0;
112 for (index_t i = 0; i < n; ++i) {
113 real_t gd = xₖ(i) - γₖ * grad_ψxₖ(i);
114 if (gd <= C.lowerbound(i)) { // i ∊ J̲ ⊆ K
115 qₖ(i) = pₖ(i); //
116 } else if (C.upperbound(i) <= gd) { // i ∊ J̅ ⊆ K
117 qₖ(i) = pₖ(i); //
118 } else { // i ∊ J
119 JK(nJ++) = i;
120 qₖ(i) = -grad_ψxₖ(i);
121 }
122 }
123
124 // There are no inactive indices J
125 if (nJ == 0) {
126 // No free variables, no Newton step possible
127 return false; // Simply use the projection step
128 }
129
130 // Compute the Hessian
131 problem->eval_hess_ψ(xₖ, *y, *Σ, 1, H.reshaped());
132
133 // There are no active indices K
134 if (nJ == n) {
135 // If all indices are free, we can factor the entire matrix.
136 // Find the minimum eigenvalue to regularize the Hessian matrix and
137 // make it positive definite.
139 // Find the minimum eigenvalue to regularize the Hessian matrix and
140 // make it positive definite.
141 Eigen::SelfAdjointEigenSolver<mat> eig{H,
142 Eigen::ComputeEigenvectors};
143
144 auto λ_min = eig.eigenvalues().minCoeff(),
145 λ_max = eig.eigenvalues().maxCoeff();
146
148 std::cout << "λ(H): " << float_to_str(λ_min, 3) << ", "
149 << float_to_str(λ_max, 3) << std::endl;
150 // Regularization
151 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
152 // Solve the system
153 qₖ = eig.eigenvectors().transpose() * qₖ;
154 qₖ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qₖ;
155 qₖ = eig.eigenvectors() * qₖ;
156 return true;
157 }
158
159 // There are active indices K
160 crindexvec J = JK.topRows(nJ);
161 rindexvec K = JK.bottomRows(n - nJ);
163
164 // Compute right-hand side of 6.1c
166 qₖ(J).noalias() -=
168
169 // If there are active indices, we need to solve the Newton system with
170 // just the inactive indices.
171 mmat HJ{HJ_storage.data(), nJ, nJ};
172 // We copy the inactive block of the Hessian to a temporary dense matrix.
173 // Since it's symmetric, only the lower part is copied.
175 H(J, J).template triangularView<Eigen::Lower>();
176
178 // Find the minimum eigenvalue to regularize the Hessian matrix and
179 // make it positive definite.
180 Eigen::SelfAdjointEigenSolver<mat> eig{HJ, Eigen::ComputeEigenvectors};
181
182 auto λ_min = eig.eigenvalues().minCoeff(),
183 λ_max = eig.eigenvalues().maxCoeff();
184
186 std::cout << "λ(H_JJ): " << float_to_str(λ_min, 3) << ", "
187 << float_to_str(λ_max, 3) << std::endl;
188 // Regularization
189 real_t ε = reg_params.min_eig * (1 + std::abs(λ_max)); // TODO
190 // Solve the system
191 auto qJ = H.col(0).topRows(nJ);
192 qJ = qₖ(J);
193 qJ = eig.eigenvectors().transpose() * qJ;
194 qJ = eig.eigenvalues().cwiseMax(ε).asDiagonal().inverse() * qJ;
195 qₖ(J) = eig.eigenvectors() * qJ;
196 return true;
197 }
198
199 /// @see @ref PANOCDirection::changed_γ
202
203 /// @see @ref PANOCDirection::reset
204 void reset() {}
205
206 /// @see @ref PANOCDirection::get_name
207 std::string get_name() const {
208 return "StructuredNewtonDirection<" +
209 std::string(config_t::get_name()) + '>';
210 }
211
212 const auto &get_params() const { return direction_params; }
213
214 private:
215 const Problem *problem = nullptr;
216#ifndef _WIN32
217 std::optional<crvec> y = std::nullopt;
218 std::optional<crvec> Σ = std::nullopt;
219#else
220 std::optional<vec> y = std::nullopt;
221 std::optional<vec> Σ = std::nullopt;
222#endif
223
224 mutable indexvec JK;
225 mutable mat H;
227
228 public:
231};
232
233// clang-format off
238// clang-format on
239
240} // namespace alpaqa
241
242#include <alpaqa/inner/panoc.hpp>
243
244namespace alpaqa {
245
246// clang-format off
251// clang-format on
252
253} // 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:63
#define ALPAQA_IF_QUADF(...)
Definition config.hpp:207
#define ALPAQA_IF_LONGD(...)
Definition config.hpp:219
#define ALPAQA_IF_FLOAT(...)
Definition config.hpp:213
#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 ,...
real_t min_eig
Minimum eigenvalue of the Hessian, scaled by , enforced by regularization using a multiple of identit...
bool print_eig
Print the minimum and maximum eigenvalue of the Hessian.
Parameters for the StructuredNewtonDirection class.
Parameters for the StructuredNewtonDirection class.
typename Conf::mat mat
Definition config.hpp:79
typename Conf::indexvec indexvec
Definition config.hpp:91
typename Conf::real_t real_t
Definition config.hpp:72
typename Conf::rindexvec rindexvec
Definition config.hpp:92
typename Conf::index_t index_t
Definition config.hpp:90
typename Conf::mmat mmat
Definition config.hpp:80
constexpr const auto inf
Definition config.hpp:98
typename Conf::rvec rvec
Definition config.hpp:77
typename Conf::crvec crvec
Definition config.hpp:78
typename Conf::vec vec
Definition config.hpp:74
std::string float_to_str(F value, int precision)
Definition print.tpp:67
typename Conf::crindexvec crindexvec
Definition config.hpp:93
Double-precision double configuration.
Definition config.hpp:160
Single-precision float configuration.
Definition config.hpp:156
long double configuration.
Definition config.hpp:165
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