7#include <Eigen/SparseCore>
26template <
typename Scalar,
int UpLo>
27 requires(UpLo == Eigen::Lower) || (UpLo == (Eigen::Lower | Eigen::Upper))
44 detail::GradientExpressionGraph<Scalar>{variable}.generate_tree(
50 for (
size_t col = 0; col < m_wrt.size(); ++col) {
51 m_wrt[col].expr->col = col;
54 for (
auto& variable : m_variables) {
55 m_graphs.emplace_back(variable);
59 for (
auto& node : m_wrt) {
63 for (
int row = 0; row < m_variables.rows(); ++row) {
64 if (m_variables[row].expr ==
nullptr) {
68 if (m_variables[row].
type() == ExpressionType::LINEAR) {
72 m_graphs[row].append_triplets(m_cached_triplets, row, m_wrt);
73 }
else if (m_variables[row].
type() > ExpressionType::LINEAR) {
76 m_nonlinear_rows.emplace_back(row);
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>();
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]);
104 result(row, col) =
Variable{Scalar(0)};
115 const Eigen::SparseMatrix<Scalar>&
value() {
116 if (m_nonlinear_rows.empty()) {
120 for (
auto& graph : m_graphs) {
121 graph.update_values();
126 auto triplets = m_cached_triplets;
129 for (
int row : m_nonlinear_rows) {
130 m_graphs[row].append_triplets(triplets, row, m_wrt);
133 m_H.setFromTriplets(triplets.begin(), triplets.end());
134 if constexpr (UpLo == Eigen::Lower) {
135 m_H = m_H.template triangularView<Eigen::Lower>();
147 Eigen::SparseMatrix<Scalar> m_H{m_variables.
rows(), m_wrt.
rows()};
159Hessian<double, Eigen::Lower | Eigen::Upper>;
#define EXPORT_TEMPLATE_DECLARE(export)
Definition SymbolExports.hpp:94
#define slp_assert(condition)
Aborts in C++.
Definition assert.hpp:25
VariableMatrix< Scalar > get() const
Returns the Hessian as a VariableMatrix.
Definition hessian.hpp:94
const Eigen::SparseMatrix< Scalar > & value()
Evaluates the Hessian at wrt's value.
Definition hessian.hpp:115
Hessian(Variable< Scalar > variable, SleipnirMatrixLike< Scalar > auto wrt)
Constructs a Hessian object.
Definition hessian.hpp:42
Hessian(Variable< Scalar > variable, Variable< Scalar > wrt)
Constructs a Hessian object.
Definition hessian.hpp:34
An autodiff variable pointing to an expression node.
Definition variable.hpp:47
A matrix of autodiff variables.
Definition variable_matrix.hpp:33
int rows() const
Returns the number of rows in the matrix.
Definition variable_matrix.hpp:972
Definition concepts.hpp:33
Converts a string literal into a format string that will be parsed at compile time and converted into...
Definition printf.h:50
type
Definition base.h:985
wpi::util::SmallVector< T > small_vector
Definition small_vector.hpp:10
static constexpr empty_t empty
Designates an uninitialized VariableMatrix.
Definition empty.hpp:11
Definition to_underlying.hpp:7
Definition StringMap.hpp:773
#define SLEIPNIR_DLLEXPORT
Definition symbol_exports.hpp:34