Sleipnir C++ API
Loading...
Searching...
No Matches
hessian.hpp
1// Copyright (c) Sleipnir contributors
2
3#pragma once
4
5#include <utility>
6
7#include <Eigen/SparseCore>
8
9#include "sleipnir/autodiff/adjoint_expression_graph.hpp"
10#include "sleipnir/autodiff/variable.hpp"
11#include "sleipnir/autodiff/variable_matrix.hpp"
12#include "sleipnir/util/concepts.hpp"
13#include "sleipnir/util/scoped_profiler.hpp"
14#include "sleipnir/util/small_vector.hpp"
15#include "sleipnir/util/solve_profiler.hpp"
16#include "sleipnir/util/symbol_exports.hpp"
17
18namespace slp {
19
29template <int UpLo>
30 requires(UpLo == Eigen::Lower) || (UpLo == (Eigen::Lower | Eigen::Upper))
31class SLEIPNIR_DLLEXPORT Hessian {
32 public:
39 Hessian(Variable variable, Variable wrt) noexcept
40 : Hessian{std::move(variable), VariableMatrix{std::move(wrt)}} {}
41
49 Hessian(Variable variable, SleipnirMatrixLike auto wrt) noexcept
50 : m_variables{detail::AdjointExpressionGraph{variable}
52 m_wrt{wrt} {
53 // Initialize column each expression's adjoint occupies in the Jacobian
54 for (size_t col = 0; col < m_wrt.size(); ++col) {
55 m_wrt[col].expr->col = col;
56 }
57
58 for (auto& variable : m_variables) {
59 m_graphs.emplace_back(variable);
60 }
61
62 // Reset col to -1
63 for (auto& node : m_wrt) {
64 node.expr->col = -1;
65 }
66
67 for (int row = 0; row < m_variables.rows(); ++row) {
68 if (m_variables[row].expr == nullptr) {
69 continue;
70 }
71
72 if (m_variables[row].type() == ExpressionType::LINEAR) {
73 // If the row is linear, compute its gradient once here and cache its
74 // triplets. Constant rows are ignored because their gradients have no
75 // nonzero triplets.
76 m_graphs[row].append_adjoint_triplets(m_cached_triplets, row, m_wrt);
77 } else if (m_variables[row].type() > ExpressionType::LINEAR) {
78 // If the row is quadratic or nonlinear, add it to the list of nonlinear
79 // rows to be recomputed in Value().
80 m_nonlinear_rows.emplace_back(row);
81 }
82 }
83
84 if (m_nonlinear_rows.empty()) {
85 m_H.setFromTriplets(m_cached_triplets.begin(), m_cached_triplets.end());
86 if constexpr (UpLo == Eigen::Lower) {
87 m_H = m_H.triangularView<Eigen::Lower>();
88 }
89 }
90
91 m_profilers.emplace_back("");
92 m_profilers.emplace_back(" ↳ graph update");
93 m_profilers.emplace_back(" ↳ adjoints");
94 m_profilers.emplace_back(" ↳ matrix build");
95 }
96
106 VariableMatrix result{VariableMatrix::empty, m_variables.rows(),
107 m_wrt.rows()};
108
109 for (int row = 0; row < m_variables.rows(); ++row) {
110 auto grad = m_graphs[row].generate_gradient_tree(m_wrt);
111 for (int col = 0; col < m_wrt.rows(); ++col) {
112 if (grad[col].expr != nullptr) {
113 result(row, col) = std::move(grad[col]);
114 } else {
115 result(row, col) = Variable{0.0};
116 }
117 }
118 }
119
120 return result;
121 }
122
128 const Eigen::SparseMatrix<double>& value() {
129 ScopedProfiler value_profiler{m_profilers[0]};
130
131 if (m_nonlinear_rows.empty()) {
132 return m_H;
133 }
134
135 ScopedProfiler graph_update_profiler{m_profilers[1]};
136
137 for (auto& graph : m_graphs) {
138 graph.update_values();
139 }
140
141 graph_update_profiler.stop();
142 ScopedProfiler adjoints_profiler{m_profilers[2]};
143
144 // Copy the cached triplets so triplets added for the nonlinear rows are
145 // thrown away at the end of the function
146 auto triplets = m_cached_triplets;
147
148 // Compute each nonlinear row of the Hessian
149 for (int row : m_nonlinear_rows) {
150 m_graphs[row].append_adjoint_triplets(triplets, row, m_wrt);
151 }
152
153 adjoints_profiler.stop();
154 ScopedProfiler matrix_build_profiler{m_profilers[3]};
155
156 if (!triplets.empty()) {
157 m_H.setFromTriplets(triplets.begin(), triplets.end());
158 if constexpr (UpLo == Eigen::Lower) {
159 m_H = m_H.triangularView<Eigen::Lower>();
160 }
161 } else {
162 // setFromTriplets() is a no-op on empty triplets, so explicitly zero out
163 // the storage
164 m_H.setZero();
165 }
166
167 return m_H;
168 }
169
175 const small_vector<SolveProfiler>& get_profilers() const {
176 return m_profilers;
177 }
178
179 private:
180 VariableMatrix m_variables;
181 VariableMatrix m_wrt;
182
183 small_vector<detail::AdjointExpressionGraph> m_graphs;
184
185 Eigen::SparseMatrix<double> m_H{m_variables.rows(), m_wrt.rows()};
186
187 // Cached triplets for gradients of linear rows
188 small_vector<Eigen::Triplet<double>> m_cached_triplets;
189
190 // List of row indices for nonlinear rows whose graients will be computed in
191 // Value()
192 small_vector<int> m_nonlinear_rows;
193
194 small_vector<SolveProfiler> m_profilers;
195};
196
197} // namespace slp
Definition hessian.hpp:31
const small_vector< SolveProfiler > & get_profilers() const
Definition hessian.hpp:175
Hessian(Variable variable, SleipnirMatrixLike auto wrt) noexcept
Definition hessian.hpp:49
Hessian(Variable variable, Variable wrt) noexcept
Definition hessian.hpp:39
VariableMatrix get() const
Definition hessian.hpp:105
const Eigen::SparseMatrix< double > & value()
Definition hessian.hpp:128
Definition scoped_profiler.hpp:19
void stop()
Definition scoped_profiler.hpp:57
Definition variable_matrix.hpp:29
int rows() const
Definition variable_matrix.hpp:907
VariableBlock< VariableMatrix > col(int col)
Definition variable_matrix.hpp:546
static constexpr empty_t empty
Definition variable_matrix.hpp:39
Definition variable.hpp:41
Definition adjoint_expression_graph.hpp:21
VariableMatrix generate_gradient_tree(const VariableMatrix &wrt) const
Definition adjoint_expression_graph.hpp:52
Definition concepts.hpp:23