Sleipnir C++ API
Loading...
Searching...
No Matches
newton.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <chrono>
6#include <cmath>
7#include <functional>
8#include <limits>
9#include <span>
10
11#include <Eigen/Core>
12#include <Eigen/SparseCore>
13#include <gch/small_vector.hpp>
14
15#include "sleipnir/optimization/solver/exit_status.hpp"
16#include "sleipnir/optimization/solver/iteration_info.hpp"
17#include "sleipnir/optimization/solver/newton_matrix_callbacks.hpp"
18#include "sleipnir/optimization/solver/options.hpp"
19#include "sleipnir/optimization/solver/util/error_estimate.hpp"
20#include "sleipnir/optimization/solver/util/filter.hpp"
21#include "sleipnir/optimization/solver/util/kkt_error.hpp"
22#include "sleipnir/optimization/solver/util/regularized_ldlt.hpp"
23#include "sleipnir/util/assert.hpp"
24#include "sleipnir/util/print_diagnostics.hpp"
25#include "sleipnir/util/scope_exit.hpp"
26#include "sleipnir/util/scoped_profiler.hpp"
27#include "sleipnir/util/solve_profiler.hpp"
28#include "sleipnir/util/symbol_exports.hpp"
29
30// See docs/algorithms.md#Works_cited for citation definitions.
31
32namespace slp {
33
52template <typename Scalar>
53ExitStatus newton(
54 const NewtonMatrixCallbacks<Scalar>& matrix_callbacks,
55 std::span<std::function<bool(const IterationInfo<Scalar>& info)>>
56 iteration_callbacks,
57 const Options& options, Eigen::Vector<Scalar, Eigen::Dynamic>& x) {
58 using DenseVector = Eigen::Vector<Scalar, Eigen::Dynamic>;
59 using SparseMatrix = Eigen::SparseMatrix<Scalar>;
60 using SparseVector = Eigen::SparseVector<Scalar>;
61
62 using std::isfinite;
63
64 const auto solve_start_time = std::chrono::steady_clock::now();
65
66 gch::small_vector<SolveProfiler> solve_profilers;
67 solve_profilers.emplace_back("solver");
68 solve_profilers.emplace_back(" ↳ setup");
69 solve_profilers.emplace_back(" ↳ iteration");
70 solve_profilers.emplace_back(" ↳ feasibility ✓");
71 solve_profilers.emplace_back(" ↳ iter callbacks");
72 solve_profilers.emplace_back(" ↳ KKT matrix decomp");
73 solve_profilers.emplace_back(" ↳ KKT system solve");
74 solve_profilers.emplace_back(" ↳ line search");
75 solve_profilers.emplace_back(" ↳ next iter prep");
76 solve_profilers.emplace_back(" ↳ f(x)");
77 solve_profilers.emplace_back(" ↳ ∇f(x)");
78 solve_profilers.emplace_back(" ↳ ∇²ₓₓL");
79
80 auto& solver_prof = solve_profilers[0];
81 auto& setup_prof = solve_profilers[1];
82 auto& inner_iter_prof = solve_profilers[2];
83 auto& feasibility_check_prof = solve_profilers[3];
84 auto& iter_callbacks_prof = solve_profilers[4];
85 auto& kkt_matrix_decomp_prof = solve_profilers[5];
86 auto& kkt_system_solve_prof = solve_profilers[6];
87 auto& line_search_prof = solve_profilers[7];
88 auto& next_iter_prep_prof = solve_profilers[8];
89
90 // Set up profiled matrix callbacks
91#ifndef SLEIPNIR_DISABLE_DIAGNOSTICS
92 auto& f_prof = solve_profilers[9];
93 auto& g_prof = solve_profilers[10];
94 auto& H_prof = solve_profilers[11];
95
96 NewtonMatrixCallbacks<Scalar> matrices{
97 [&](const DenseVector& x) -> Scalar {
98 ScopedProfiler prof{f_prof};
99 return matrix_callbacks.f(x);
100 },
101 [&](const DenseVector& x) -> SparseVector {
102 ScopedProfiler prof{g_prof};
103 return matrix_callbacks.g(x);
104 },
105 [&](const DenseVector& x) -> SparseMatrix {
106 ScopedProfiler prof{H_prof};
107 return matrix_callbacks.H(x);
108 }};
109#else
110 const auto& matrices = matrix_callbacks;
111#endif
112
113 solver_prof.start();
114 setup_prof.start();
115
116 Scalar f = matrices.f(x);
117
118 int num_decision_variables = x.rows();
119
120 SparseVector g = matrices.g(x);
121 SparseMatrix H = matrices.H(x);
122
123 // Ensure matrix callback dimensions are consistent
124 slp_assert(g.rows() == num_decision_variables);
125 slp_assert(H.rows() == num_decision_variables);
126 slp_assert(H.cols() == num_decision_variables);
127
128 // Check whether initial guess has finite f(xₖ)
129 if (!isfinite(f)) {
130 return ExitStatus::NONFINITE_INITIAL_COST_OR_CONSTRAINTS;
131 }
132
133 int iterations = 0;
134
135 Filter<Scalar> filter;
136
137 RegularizedLDLT<Scalar> solver{num_decision_variables, 0};
138
139 // Variables for determining when a step is acceptable
140 constexpr Scalar α_reduction_factor(0.5);
141 constexpr Scalar α_min(1e-20);
142
143 // Error estimate
144 Scalar E_0 = std::numeric_limits<Scalar>::infinity();
145
146 setup_prof.stop();
147
148 // Prints final solver diagnostics when the solver exits
149 scope_exit exit{[&] {
150 if (options.diagnostics) {
151 solver_prof.stop();
152 if (iterations > 0) {
153 print_bottom_iteration_diagnostics();
154 }
155 print_solver_diagnostics(solve_profilers);
156 }
157 }};
158
159 while (E_0 > Scalar(options.tolerance)) {
160 ScopedProfiler inner_iter_profiler{inner_iter_prof};
161 ScopedProfiler feasibility_check_profiler{feasibility_check_prof};
162
163 // Check for diverging iterates
164 if (x.template lpNorm<Eigen::Infinity>() > Scalar(1e10) || !x.allFinite()) {
165 return ExitStatus::DIVERGING_ITERATES;
166 }
167
168 feasibility_check_profiler.stop();
169 ScopedProfiler iter_callbacks_profiler{iter_callbacks_prof};
170
171 // Call iteration callbacks
172 for (const auto& callback : iteration_callbacks) {
173 if (callback({iterations, x, g, H, SparseMatrix{}, SparseMatrix{}})) {
174 return ExitStatus::CALLBACK_REQUESTED_STOP;
175 }
176 }
177
178 iter_callbacks_profiler.stop();
179 ScopedProfiler kkt_matrix_decomp_profiler{kkt_matrix_decomp_prof};
180
181 // Solve the Newton-KKT system
182 //
183 // Hpˣ = −∇f
184 solver.compute(H);
185
186 kkt_matrix_decomp_profiler.stop();
187 ScopedProfiler kkt_system_solve_profiler{kkt_system_solve_prof};
188
189 DenseVector p_x = solver.solve(-g);
190
191 kkt_system_solve_profiler.stop();
192 ScopedProfiler line_search_profiler{line_search_prof};
193
194 constexpr Scalar α_max(1);
195 Scalar α = α_max;
196
197 // Loop until a step is accepted. If a step becomes acceptable, the loop
198 // will exit early.
199 while (1) {
200 DenseVector trial_x = x + α * p_x;
201
202 Scalar trial_f = matrices.f(trial_x);
203
204 // If f(xₖ + αpₖˣ) isn't finite, reduce step size immediately
205 if (!isfinite(trial_f)) {
206 // Reduce step size
207 α *= α_reduction_factor;
208
209 if (α < α_min) {
210 return ExitStatus::LINE_SEARCH_FAILED;
211 }
212 continue;
213 }
214
215 // Check whether filter accepts trial iterate
216 if (filter.try_add(FilterEntry{trial_f}, α)) {
217 // Accept step
218 break;
219 }
220
221 // Reduce step size
222 α *= α_reduction_factor;
223
224 // If step size hit a minimum, check if the KKT error was reduced. If it
225 // wasn't, report bad line search.
226 if (α < α_min) {
227 Scalar current_kkt_error = kkt_error<Scalar>(g);
228
229 DenseVector trial_x = x + α_max * p_x;
230
231 Scalar next_kkt_error = kkt_error<Scalar>(matrices.g(trial_x));
232
233 // If the step using αᵐᵃˣ reduced the KKT error, accept it anyway
234 if (next_kkt_error <= Scalar(0.999) * current_kkt_error) {
235 α = α_max;
236
237 // Accept step
238 break;
239 }
240
241 return ExitStatus::LINE_SEARCH_FAILED;
242 }
243 }
244
245 line_search_profiler.stop();
246
247 // xₖ₊₁ = xₖ + αₖpₖˣ
248 x += α * p_x;
249
250 // Update autodiff for Hessian
251 f = matrices.f(x);
252 g = matrices.g(x);
253 H = matrices.H(x);
254
255 ScopedProfiler next_iter_prep_profiler{next_iter_prep_prof};
256
257 // Update the error estimate
258 E_0 = error_estimate<Scalar>(g);
259
260 next_iter_prep_profiler.stop();
261 inner_iter_profiler.stop();
262
263 if (options.diagnostics) {
264 print_iteration_diagnostics(iterations, IterationType::NORMAL,
265 inner_iter_profiler.current_duration(), E_0,
266 f, Scalar(0), Scalar(0), Scalar(0),
267 solver.hessian_regularization(), α, α_max,
268 α_reduction_factor, Scalar(1));
269 }
270
271 ++iterations;
272
273 // Check for max iterations
274 if (iterations >= options.max_iterations) {
275 return ExitStatus::MAX_ITERATIONS_EXCEEDED;
276 }
277
278 // Check for max wall clock time
279 if (std::chrono::steady_clock::now() - solve_start_time > options.timeout) {
280 return ExitStatus::TIMEOUT;
281 }
282 }
283
284 return ExitStatus::SUCCESS;
285}
286
287extern template SLEIPNIR_DLLEXPORT ExitStatus
288newton(const NewtonMatrixCallbacks<double>& matrix_callbacks,
289 std::span<std::function<bool(const IterationInfo<double>& info)>>
290 iteration_callbacks,
291 const Options& options, Eigen::Vector<double, Eigen::Dynamic>& x);
292
293} // namespace slp