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/small_vector.hpp"
14#include "sleipnir/util/symbol_exports.hpp"
15
16namespace slp {
17
27template <int UpLo>
28 requires(UpLo == Eigen::Lower) || (UpLo == (Eigen::Lower | Eigen::Upper))
29class SLEIPNIR_DLLEXPORT Hessian {
30 public:
37 Hessian(Variable variable, Variable wrt) noexcept
38 : Hessian{std::move(variable), VariableMatrix{std::move(wrt)}} {}
39
47 Hessian(Variable variable, SleipnirMatrixLike auto wrt) noexcept
48 : m_variables{detail::AdjointExpressionGraph{variable}
50 m_wrt{wrt} {
51 // Initialize column each expression's adjoint occupies in the Jacobian
52 for (size_t col = 0; col < m_wrt.size(); ++col) {
53 m_wrt[col].expr->col = col;
54 }
55
56 for (auto& variable : m_variables) {
57 m_graphs.emplace_back(variable);
58 }
59
60 // Reset col to -1
61 for (auto& node : m_wrt) {
62 node.expr->col = -1;
63 }
64
65 for (int row = 0; row < m_variables.rows(); ++row) {
66 if (m_variables[row].expr == nullptr) {
67 continue;
68 }
69
70 if (m_variables[row].type() == ExpressionType::LINEAR) {
71 // If the row is linear, compute its gradient once here and cache its
72 // triplets. Constant rows are ignored because their gradients have no
73 // nonzero triplets.
74 m_graphs[row].append_adjoint_triplets(m_cached_triplets, row, m_wrt);
75 } else if (m_variables[row].type() > ExpressionType::LINEAR) {
76 // If the row is quadratic or nonlinear, add it to the list of nonlinear
77 // rows to be recomputed in Value().
78 m_nonlinear_rows.emplace_back(row);
79 }
80 }
81
82 if (m_nonlinear_rows.empty()) {
83 m_H.setFromTriplets(m_cached_triplets.begin(), m_cached_triplets.end());
84 if constexpr (UpLo == Eigen::Lower) {
85 m_H = m_H.triangularView<Eigen::Lower>();
86 }
87 }
88 }
89
99 VariableMatrix result{VariableMatrix::empty, m_variables.rows(),
100 m_wrt.rows()};
101
102 for (int row = 0; row < m_variables.rows(); ++row) {
103 auto grad = m_graphs[row].generate_gradient_tree(m_wrt);
104 for (int col = 0; col < m_wrt.rows(); ++col) {
105 if (grad[col].expr != nullptr) {
106 result[row, col] = std::move(grad[col]);
107 } else {
108 result[row, col] = Variable{0.0};
109 }
110 }
111 }
112
113 return result;
114 }
115
121 const Eigen::SparseMatrix<double>& value() {
122 if (m_nonlinear_rows.empty()) {
123 return m_H;
124 }
125
126 for (auto& graph : m_graphs) {
127 graph.update_values();
128 }
129
130 // Copy the cached triplets so triplets added for the nonlinear rows are
131 // thrown away at the end of the function
132 auto triplets = m_cached_triplets;
133
134 // Compute each nonlinear row of the Hessian
135 for (int row : m_nonlinear_rows) {
136 m_graphs[row].append_adjoint_triplets(triplets, row, m_wrt);
137 }
138
139 if (!triplets.empty()) {
140 m_H.setFromTriplets(triplets.begin(), triplets.end());
141 if constexpr (UpLo == Eigen::Lower) {
142 m_H = m_H.triangularView<Eigen::Lower>();
143 }
144 } else {
145 // setFromTriplets() is a no-op on empty triplets, so explicitly zero out
146 // the storage
147 m_H.setZero();
148 }
149
150 return m_H;
151 }
152
153 private:
154 VariableMatrix m_variables;
155 VariableMatrix m_wrt;
156
157 small_vector<detail::AdjointExpressionGraph> m_graphs;
158
159 Eigen::SparseMatrix<double> m_H{m_variables.rows(), m_wrt.rows()};
160
161 // Cached triplets for gradients of linear rows
162 small_vector<Eigen::Triplet<double>> m_cached_triplets;
163
164 // List of row indices for nonlinear rows whose graients will be computed in
165 // Value()
166 small_vector<int> m_nonlinear_rows;
167};
168
169} // namespace slp
Definition hessian.hpp:29
Hessian(Variable variable, SleipnirMatrixLike auto wrt) noexcept
Definition hessian.hpp:47
Hessian(Variable variable, Variable wrt) noexcept
Definition hessian.hpp:37
VariableMatrix get() const
Definition hessian.hpp:98
const Eigen::SparseMatrix< double > & value()
Definition hessian.hpp:121
Definition variable_matrix.hpp:29
int rows() const
Definition variable_matrix.hpp:951
VariableBlock< VariableMatrix > col(int col)
Definition variable_matrix.hpp:546
static constexpr empty_t empty
Definition variable_matrix.hpp:39
Definition variable.hpp:40
Definition adjoint_expression_graph.hpp:21
VariableMatrix generate_gradient_tree(const VariableMatrix &wrt) const
Definition adjoint_expression_graph.hpp:52
Definition concepts.hpp:30