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