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#include <gch/small_vector.hpp>
9
10#include "sleipnir/autodiff/gradient_expression_graph.hpp"
11#include "sleipnir/autodiff/variable.hpp"
12#include "sleipnir/autodiff/variable_matrix.hpp"
13#include "sleipnir/util/assert.hpp"
14#include "sleipnir/util/concepts.hpp"
15
16namespace slp {
17
26template <typename Scalar, int UpLo>
27 requires(UpLo == Eigen::Lower) || (UpLo == (Eigen::Lower | Eigen::Upper))
28class Hessian {
29 public:
36
43 : m_variables{
44 detail::GradientExpressionGraph<Scalar>{variable}.generate_tree(
45 wrt)},
46 m_wrt{wrt} {
47 slp_assert(m_wrt.cols() == 1);
48
49 // Initialize column each expression's adjoint occupies in the Jacobian
50 for (size_t col = 0; col < m_wrt.size(); ++col) {
51 m_wrt[col].expr->col = col;
52 }
53
54 for (auto& variable : m_variables) {
55 m_graphs.emplace_back(variable);
56 }
57
58 // Reset col to -1
59 for (auto& node : m_wrt) {
60 node.expr->col = -1;
61 }
62
63 for (int row = 0; row < m_variables.rows(); ++row) {
64 if (m_variables[row].expr == nullptr) {
65 continue;
66 }
67
68 if (m_variables[row].type() == ExpressionType::LINEAR) {
69 // If the row is linear, compute its gradient once here and cache its
70 // triplets. Constant rows are ignored because their gradients have no
71 // nonzero triplets.
72 m_graphs[row].append_triplets(m_cached_triplets, row, m_wrt);
73 } else if (m_variables[row].type() > ExpressionType::LINEAR) {
74 // If the row is quadratic or nonlinear, add it to the list of nonlinear
75 // rows to be recomputed in value().
76 m_nonlinear_rows.emplace_back(row);
77 }
78 }
79
80 if (m_nonlinear_rows.empty()) {
81 m_H.setFromTriplets(m_cached_triplets.begin(), m_cached_triplets.end());
82 if constexpr (UpLo == Eigen::Lower) {
83 m_H = m_H.template triangularView<Eigen::Lower>();
84 }
85 }
86 }
87
95 VariableMatrix<Scalar> result{detail::empty, m_variables.rows(),
96 m_wrt.rows()};
97
98 for (int row = 0; row < m_variables.rows(); ++row) {
99 auto grad = m_graphs[row].generate_tree(m_wrt);
100 for (int col = 0; col < m_wrt.rows(); ++col) {
101 if (grad[col].expr != nullptr) {
102 result[row, col] = std::move(grad[col]);
103 } else {
104 result[row, col] = Variable{Scalar(0)};
105 }
106 }
107 }
108
109 return result;
110 }
111
115 const Eigen::SparseMatrix<Scalar>& value() {
116 if (m_nonlinear_rows.empty()) {
117 return m_H;
118 }
119
120 for (auto& graph : m_graphs) {
121 graph.update_values();
122 }
123
124 // Copy the cached triplets so triplets added for the nonlinear rows are
125 // thrown away at the end of the function
126 auto triplets = m_cached_triplets;
127
128 // Compute each nonlinear row of the Hessian
129 for (int row : m_nonlinear_rows) {
130 m_graphs[row].append_triplets(triplets, row, m_wrt);
131 }
132
133 m_H.setFromTriplets(triplets.begin(), triplets.end());
134 if constexpr (UpLo == Eigen::Lower) {
135 m_H = m_H.template triangularView<Eigen::Lower>();
136 }
137
138 return m_H;
139 }
140
141 private:
142 VariableMatrix<Scalar> m_variables;
144
145 gch::small_vector<detail::GradientExpressionGraph<Scalar>> m_graphs;
146
147 Eigen::SparseMatrix<Scalar> m_H{m_variables.rows(), m_wrt.rows()};
148
149 // Cached triplets for gradients of linear rows
150 gch::small_vector<Eigen::Triplet<Scalar>> m_cached_triplets;
151
152 // List of row indices for nonlinear rows whose graients will be computed in
153 // value()
154 gch::small_vector<int> m_nonlinear_rows;
155};
156
157// @cond Suppress Doxygen
158extern template class EXPORT_TEMPLATE_DECLARE(SLEIPNIR_DLLEXPORT)
159Hessian<double, Eigen::Lower | Eigen::Upper>;
160// @endcond
161
162} // namespace slp
Definition hessian.hpp:28
VariableMatrix< Scalar > get() const
Definition hessian.hpp:94
const Eigen::SparseMatrix< Scalar > & value()
Definition hessian.hpp:115
Hessian(Variable< Scalar > variable, SleipnirMatrixLike< Scalar > auto wrt)
Definition hessian.hpp:42
Hessian(Variable< Scalar > variable, Variable< Scalar > wrt)
Definition hessian.hpp:34
Definition intrusive_shared_ptr.hpp:27
Definition variable_matrix.hpp:33
Definition variable.hpp:47
Definition concepts.hpp:33