Nonconvex constrained optimization
Loading...
Searching...
No Matches
convex-newton.hpp
Go to the documentation of this file.
1#pragma once
2
4#include <alpaqa/export.hpp>
12#include <iostream>
13#include <optional>
14#include <stdexcept>
15
16#include <Eigen/Cholesky>
17#include <Eigen/Eigenvalues>
18#include <Eigen/src/Cholesky/LDLT.h>
19
20namespace alpaqa {
21
22/// Parameters for the @ref ConvexNewtonDirection class.
23/// @ingroup grp_Parameters
24template <Config Conf>
27 real_t ζ = real_t(1e-8);
29 bool ldlt = false;
30};
31
32/// Parameters for the @ref ConvexNewtonDirection class.
33/// @ingroup grp_Parameters
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 bool quadratic = false;
43};
44
45/// @ingroup grp_DirectionProviders
46template <Config Conf = DefaultConfig>
52
57
62 const DirectionParams &directionparams = {})
63 : reg_params(params), direction_params(directionparams) {}
64
65 /// @see @ref PANOCDirection::initialize
66 void initialize(const Problem &problem, [[maybe_unused]] crvec y,
67 [[maybe_unused]] crvec Σ, [[maybe_unused]] real_t γ_0,
68 [[maybe_unused]] crvec x_0, [[maybe_unused]] crvec x̂_0,
69 [[maybe_unused]] crvec p_0,
70 [[maybe_unused]] crvec grad_ψx_0) {
71 if (problem.get_num_constraints() != 0)
72 throw std::invalid_argument("Convex Newton direction does not "
73 "support general constraints");
74 if (!problem.provides_eval_inactive_indices_res_lna())
75 throw std::invalid_argument("Convex Newton direction requires "
76 "eval_inactive_indices_res_lna");
77 if (!problem.provides_eval_lagrangian_hessian())
78 throw std::invalid_argument("Structured Newton requires "
79 "eval_lagrangian_hessian");
80 // Store reference to problem
81 this->problem = &problem;
82 // Allocate workspaces
83 const auto n = problem.get_num_variables();
84 JK.resize(n);
85 JK_old.resize(n);
86 nJ_old = -1;
87 HJ_storage.resize(n * n);
88 work.resize(n);
89 auto sparsity = problem.get_lagrangian_hessian_sparsity();
90 if (!is_dense(sparsity))
91 std::cerr << "Sparse hessians not yet implemented, converting to "
92 "dense matrix (may be very slow)\n";
93 H_eval.emplace(sparsity);
94 have_hess = false;
95 }
96
97 /// @see @ref PANOCDirection::has_initial_direction
98 bool has_initial_direction() const { return true; }
99
100 /// @see @ref PANOCDirection::update
101 bool update([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t γₙₑₓₜ,
102 [[maybe_unused]] crvec xₖ, [[maybe_unused]] crvec xₙₑₓₜ,
103 [[maybe_unused]] crvec pₖ, [[maybe_unused]] crvec pₙₑₓₜ,
104 [[maybe_unused]] crvec grad_ψxₖ,
105 [[maybe_unused]] crvec grad_ψxₙₑₓₜ) {
106 return true;
107 }
108
109 template <class Solver>
110 void solve(rmat H, rvec q) const {
112 Solver ll{H};
113 if (ll.info() != Eigen::Success)
114 throw std::runtime_error("Cholesky factorization failed. "
115 "Is the problem convex?");
116 ll.solveInPlace(q);
117 }
118
119 /// @see @ref PANOCDirection::apply
120 bool apply(real_t γₖ, crvec xₖ, [[maybe_unused]] crvec x̂ₖ, crvec pₖ,
121 crvec grad_ψxₖ, rvec qₖ) const {
122 length_t n = xₖ.size();
123 // Evaluate the Hessian
124 if (!have_hess) {
125 assert(H_eval);
126 const auto &y = null_vec<config_t>;
127 H.emplace(H_eval->eval([&](rvec v) {
128 problem->eval_lagrangian_hessian(xₖ, y, 1, v);
129 }));
130 have_hess = direction_params.quadratic;
131 }
132 assert(H);
133 // Find inactive indices J
134 auto nJ = problem->eval_inactive_indices_res_lna(γₖ, xₖ, grad_ψxₖ, JK);
135 auto J = JK.topRows(nJ);
136 auto HJ = HJ_storage.topRows(nJ * nJ).reshaped(nJ, nJ);
137 HJ = (*H)(J, J);
138 // Regularize the Hessian
139 real_t res_sq = pₖ.squaredNorm() / (γₖ * γₖ);
140 real_t reg = reg_params.ζ * std::pow(res_sq, reg_params.ν / 2);
141 HJ += reg * mat::Identity(nJ, nJ);
142 // Compute the right-hand side
143 qₖ = pₖ;
144 rvec w = work.topRows(nJ);
145 if (direction_params.hessian_vec_factor != 0) {
146 rindexvec K = JK.bottomRows(n - nJ);
148 w = (real_t(1) / γₖ) * pₖ(J) -
149 direction_params.hessian_vec_factor * ((*H)(J, K) * qₖ(K));
150 } else {
151 w = (real_t(1) / γₖ) * pₖ(J);
152 }
153 // Solve the system
154 auto symmetry = H_eval->converter.converter.get_sparsity().symmetry;
155 if (symmetry == sparsity::Symmetry::Upper)
158 else
161 qₖ(J) = w;
162 return true;
163 }
164
165 /// @see @ref PANOCDirection::changed_γ
166 void changed_γ([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t old_γₖ) {
167 }
168
169 /// @see @ref PANOCDirection::reset
170 void reset() {}
171
172 /// @see @ref PANOCDirection::get_name
173 std::string get_name() const {
174 return "ConvexNewtonDirection<" + std::string(config_t::get_name()) +
175 '>';
176 }
177
178 const auto &get_params() const { return direction_params; }
179
180 private:
181 const Problem *problem = nullptr;
182
184 mutable index_t nJ_old = -1;
185 mutable std::optional<rmat> H;
186 mutable std::optional<sparsity::DenseEvaluator<config_t>> H_eval;
188 mutable bool have_hess = false;
189
190 public:
193};
194
195// clang-format off
200// clang-format on
201
202} // 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 ConvexNewtonDirection class.
Parameters for the ConvexNewtonDirection class.
typename Conf::rmat rmat
Definition config.hpp:96
typename Conf::indexvec indexvec
Definition config.hpp:105
typename Conf::real_t real_t
Definition config.hpp:86
const rvec< Conf > null_vec
Global empty vector for convenience.
Definition config.hpp:193
typename Conf::rindexvec rindexvec
Definition config.hpp:106
typename Conf::index_t index_t
Definition config.hpp:104
typename Conf::length_t length_t
Definition config.hpp:103
typename Conf::rvec rvec
Definition config.hpp:91
typename Conf::crvec crvec
Definition config.hpp:92
typename Conf::vec vec
Definition config.hpp:88
void changed_γ(real_t γₖ, real_t old_γₖ)
ConvexNewtonRegularizationParams< config_t > AcceleratorParams
ConvexNewtonDirection(const Params &params)
const auto & get_params() const
std::optional< sparsity::DenseEvaluator< config_t > > H_eval
bool update(real_t γₖ, real_t γₙₑₓₜ, crvec xₖ, crvec xₙₑₓₜ, crvec pₖ, crvec pₙₑₓₜ, crvec grad_ψxₖ, crvec grad_ψxₙₑₓₜ)
ConvexNewtonDirection(const AcceleratorParams &params, const DirectionParams &directionparams={})
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
void solve(rmat H, rvec q) const
TypeErasedProblem< config_t > Problem
ConvexNewtonDirectionParams< config_t > DirectionParams
Double-precision double configuration.
Definition config.hpp:176
Single-precision float configuration.
Definition config.hpp:172
long double configuration.
Definition config.hpp:181
static void compute_complement(std::span< const index_t > in, std::span< index_t > out)
Definition index-set.hpp:37