alpaqa 1.0.0a10
Nonconvex constrained optimization
Loading...
Searching...
No Matches
newton-tr.hpp
Go to the documentation of this file.
1#pragma once
2
7#include <cmath>
8#include <limits>
9#include <optional>
10#include <stdexcept>
11
12namespace alpaqa {
13
14/// Parameters for the @ref NewtonTRDirection class.
15template <Config Conf>
20 bool finite_diff = false;
22 std::sqrt(std::numeric_limits<real_t>::epsilon());
23};
24
25/// @ingroup grp_DirectionProviders
26template <Config Conf>
29
33
34 struct Params {
37 };
38
39 NewtonTRDirection() = default;
41 : steihaug(params.accelerator), direction_params(params.direction) {}
43 const DirectionParams &directionparams = {})
44 : steihaug(params), direction_params(directionparams) {}
45
46 /// @see @ref PANTRDirection::initialize
47 void initialize(const Problem &problem, [[maybe_unused]] crvec y,
48 [[maybe_unused]] crvec Σ, [[maybe_unused]] real_t γ_0,
49 [[maybe_unused]] crvec x_0, [[maybe_unused]] crvec x̂_0,
50 [[maybe_unused]] crvec p_0,
51 [[maybe_unused]] crvec grad_ψx_0) {
54 throw std::invalid_argument("NewtonTR without finite differences "
55 "requires Problem::eval_hess_ψ_prod()");
57 throw std::invalid_argument(
58 "NewtonTR requires "
59 "Problem::eval_inactive_indices_res_lna()");
60 // Store references to problem and ALM variables
61 this->problem = &problem;
62 this->y.emplace(y);
63 this->Σ.emplace(Σ);
64 // Resize workspaces
65 const auto n = problem.get_n(), m = problem.get_m();
66 JK_sto.resize(n);
67 rJ_sto.resize(n);
68 qJ_sto.resize(n);
69 work.resize(n);
70 work_2.resize(n);
73 work_n_fd.resize(n);
74 work_m_fd.resize(m);
75 }
76 }
77
78 /// @see @ref PANTRDirection::has_initial_direction
79 bool has_initial_direction() const { return true; }
80
81 /// @see @ref PANTRDirection::update
82 bool update([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t γₙₑₓₜ,
83 [[maybe_unused]] crvec xₖ, [[maybe_unused]] crvec xₙₑₓₜ,
84 [[maybe_unused]] crvec pₖ, [[maybe_unused]] crvec pₙₑₓₜ,
85 [[maybe_unused]] crvec grad_ψxₖ,
86 [[maybe_unused]] crvec grad_ψxₙₑₓₜ) {
87 return true;
88 }
89
90 /// @see @ref PANTRDirection::apply
91 real_t apply([[maybe_unused]] real_t γₖ, [[maybe_unused]] crvec xₖ,
92 [[maybe_unused]] crvec x̂ₖ, crvec pₖ,
93 [[maybe_unused]] crvec grad_ψxₖ, real_t radius,
94 rvec qₖ) const {
95
96 if (!std::isfinite(radius))
97 throw std::logic_error("Invalid trust radius");
98 if (radius < std::numeric_limits<real_t>::epsilon())
99 throw std::logic_error("Trust radius too small");
100
101 // Newton with exact Hessian
102
103 // Find inactive and active constraints
104 const auto n = problem->get_n();
105 index_t nJ =
106 problem->eval_inactive_indices_res_lna(γₖ, xₖ, grad_ψxₖ, JK_sto);
107 crindexvec J = JK_sto.topRows(nJ);
108 rindexvec K = JK_sto.bottomRows(n - nJ);
110 auto rJ = rJ_sto.topRows(nJ);
111 auto qJ = qJ_sto.topRows(nJ);
112 rJ = (-real_t(1) / γₖ) * pₖ(J);
113 qₖ(K) = pₖ(K);
114 qₖ(J).setZero();
115 real_t norm_qK_sq = pₖ(K).squaredNorm();
116
117 // Hessian-vector term
120 real_t ε = (1 + grad_ψxₖ.norm()) *
122 work = xₖ + ε * qₖ;
124 work_m_fd);
125 rJ.noalias() += (work_2 - grad_ψxₖ)(J) *
127 } else {
128 problem->eval_hess_ψ_prod(xₖ, *y, *Σ, 1, qₖ, work);
129 rJ.noalias() += work(J) * direction_params.hessian_vec_factor;
130 }
131 }
132
133 // Hessian-vector product on subset J
134 auto hess_vec_mult = [&](crvec p, rvec Bp) {
136 real_t ε = (1 + grad_ψxₖ.norm()) *
138 work = xₖ;
139 work(J) += ε * p;
141 work_m_fd);
142 Bp.topRows(nJ) = (work_2 - grad_ψxₖ)(J) / ε;
143 } else {
144 work.setZero();
145 work(J) = p;
146 problem->eval_hess_ψ_prod(xₖ, *y, *Σ, 1, work, work_2);
147 Bp.topRows(nJ) = work_2(J);
148 }
149 };
150
151 // Steihaug conjugate gradients
152 real_t qJ_model = steihaug.solve(rJ, hess_vec_mult, radius, qJ);
153 qₖ(J) = qJ;
154 return qJ_model - norm_qK_sq / (2 * γₖ);
155 }
156
157 /// @see @ref PANTRDirection::changed_γ
158 void changed_γ([[maybe_unused]] real_t γₖ, [[maybe_unused]] real_t old_γₖ) {
160 throw std::invalid_argument("NewtonTRDirection does not support "
161 "rescale_on_step_size_changes");
162 }
163
164 /// @see @ref PANTRDirection::reset
165 void reset() {}
166
167 /// @see @ref PANTRDirection::get_name
168 std::string get_name() const {
169 return "NewtonTRDirection<" + std::string(config_t::get_name()) + '>';
170 }
171
172 auto get_params() const {
173 return std::tie(steihaug.params, direction_params);
174 }
175
178 const Problem *problem = nullptr;
179#ifndef _WIN32
180 std::optional<crvec> y = std::nullopt;
181 std::optional<crvec> Σ = std::nullopt;
182#else
183 std::optional<vec> y = std::nullopt;
184 std::optional<vec> Σ = std::nullopt;
185#endif
187 mutable vec rJ_sto;
188 mutable vec qJ_sto;
190};
191
192} // namespace alpaqa
bool provides_eval_hess_ψ_prod() const
Returns true if the problem provides an implementation of eval_hess_ψ_prod.
length_t get_n() const
[Required] Number of decision variables.
length_t get_m() const
[Required] Number of constraints.
bool provides_eval_inactive_indices_res_lna() const
Returns true if the problem provides an implementation of eval_inactive_indices_res_lna.
index_t eval_inactive_indices_res_lna(real_t γ, crvec x, crvec grad_ψ, rindexvec J) const
[Optional] Function that computes the inactive indices for the evaluation of the linear Newton appro...
void eval_grad_ψ(crvec x, crvec y, crvec Σ, rvec grad_ψ, rvec work_n, rvec work_m) const
[Optional] Calculate the gradient ∇ψ(x).
void eval_hess_ψ_prod(crvec x, crvec y, crvec Σ, real_t scale, crvec v, rvec Hv) const
[Optional] Function that evaluates the Hessian of the augmented Lagrangian multiplied by a vector,
#define USING_ALPAQA_CONFIG(Conf)
Definition: config.hpp:54
typename Conf::indexvec indexvec
Definition: config.hpp:76
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::rvec rvec
Definition: config.hpp:67
typename Conf::crvec crvec
Definition: config.hpp:68
typename Conf::vec vec
Definition: config.hpp:64
typename Conf::crindexvec crindexvec
Definition: config.hpp:78
Parameters for the NewtonTRDirection class.
Definition: newton-tr.hpp:16
std::string get_name() const
Definition: newton-tr.hpp:168
void changed_γ(real_t γₖ, real_t old_γₖ)
Definition: newton-tr.hpp:158
SteihaugCG< config_t > steihaug
Definition: newton-tr.hpp:176
real_t apply(real_t γₖ, crvec xₖ, crvec x̂ₖ, crvec pₖ, crvec grad_ψxₖ, real_t radius, rvec qₖ) const
Definition: newton-tr.hpp:91
NewtonTRDirection(const Params &params)
Definition: newton-tr.hpp:40
NewtonTRDirection(const AcceleratorParams &params, const DirectionParams &directionparams={})
Definition: newton-tr.hpp:42
DirectionParams direction_params
Definition: newton-tr.hpp:177
bool update(real_t γₖ, real_t γₙₑₓₜ, crvec xₖ, crvec xₙₑₓₜ, crvec pₖ, crvec pₙₑₓₜ, crvec grad_ψxₖ, crvec grad_ψxₙₑₓₜ)
Definition: newton-tr.hpp:82
std::optional< crvec > y
Definition: newton-tr.hpp:180
void initialize(const Problem &problem, crvec y, crvec Σ, real_t γ_0, crvec x_0, crvec x̂_0, crvec p_0, crvec grad_ψx_0)
Definition: newton-tr.hpp:47
const Problem * problem
Definition: newton-tr.hpp:178
bool has_initial_direction() const
Definition: newton-tr.hpp:79
std::optional< crvec > Σ
Definition: newton-tr.hpp:181
Steihaug conjugate gradients procedure based on https://github.com/scipy/scipy/blob/583e70a50573169fc...
Definition: steihaugcg.hpp:20
void resize(length_t n)
Definition: steihaugcg.hpp:31
real_t solve(const auto &grad, const HessFun &hess_prod, real_t trust_radius, rvec step) const
Definition: steihaugcg.hpp:40
static void compute_complement(std::span< const index_t > in, std::span< index_t > out)
Definition: index-set.hpp:36