smooth
A C++ library for Lie theory
Loading...
Searching...
No Matches
lie_sparse.hpp
Go to the documentation of this file.
1// Copyright (C) 2022 Petter Nilsson. MIT License.
2
3#pragma once
4
5#include <Eigen/Sparse>
6
7#include "detail/utils.hpp"
8#include "lie_groups.hpp"
9
27SMOOTH_BEGIN_NAMESPACE
28
29namespace traits {
30
34template<typename G>
36{};
37
38} // namespace traits
39
43
47template<LieGroup G>
48inline std::array<Eigen::SparseMatrix<Scalar<G>>, Dof<G>> generators_sparse =
49 []() -> std::array<Eigen::SparseMatrix<Scalar<G>>, Dof<G>> {
50 std::array<Eigen::SparseMatrix<Scalar<G>>, Dof<G>> ret;
51 for (auto i = 0u; i < Dof<G>; ++i) {
52 ret[i] = ad<G>(Tangent<G>::Unit(i)).sparseView();
53 ret[i].makeCompressed();
54 }
55 return ret;
56}();
57
61template<LieGroup G>
62inline Eigen::SparseMatrix<Scalar<G>> ad_sparse_pattern = [] {
63 Eigen::SparseMatrix<Scalar<G>> ret(Dof<G>, Dof<G>);
64 if constexpr (!IsCommutative<G>) {
65 for (const auto & gen : generators_sparse<G>) { ret += gen; }
66 }
67 ret.makeCompressed();
68 ret.coeffs().setZero();
69 return ret;
70}();
71
82template<LieGroup G>
83void ad_sparse(Eigen::SparseMatrix<Scalar<G>> & sp, const Tangent<G> & a);
84
88template<LieGroup G>
89inline Eigen::SparseMatrix<Scalar<G>> d_exp_sparse_pattern = [] {
90 Eigen::SparseMatrix<Scalar<G>> ret(Dof<G>, Dof<G>);
91 if constexpr (IsCommutative<G>) {
92 // identity matrix
93 for (auto i = 0u; i < Dof<G>; ++i) { ret.insert(i, i) = Scalar<G>(1); }
94 } else if constexpr (requires { traits::lie_sparse<G>::d_exp_sparse_pattern; }) {
95 // use specialized method if one exists
97 } else {
98 // fall back on dense pattern
99 for (auto i = 0u; i < Dof<G>; ++i) {
100 for (auto j = 0u; j < Dof<G>; ++j) { ret.insert(i, j) = (i == j) ? Scalar<G>(1) : Scalar<G>(0); }
101 }
102 }
103 ret.makeCompressed();
104 return ret;
105}();
106
121template<LieGroup G, bool Inv = false>
122void dr_exp_sparse(Eigen::SparseMatrix<Scalar<G>> & sp, const Tangent<G> & a, Eigen::Index i0 = 0);
123
137template<LieGroup G>
138void dr_expinv_sparse(Eigen::SparseMatrix<Scalar<G>> & sp, const Tangent<G> & a, Eigen::Index i0 = 0);
139
143template<LieGroup G>
144inline Eigen::SparseMatrix<Scalar<G>> d2_exp_sparse_pattern = [] {
145 Eigen::SparseMatrix<Scalar<G>> ret(Dof<G>, Dof<G> * Dof<G>);
146 if constexpr (IsCommutative<G>) {
147 // zero matrix--do nothing
148 } else if constexpr (requires { traits::lie_sparse<G>::d2_exp_sparse_pattern; }) {
149 // use spacialized method
151 } else {
152 // dense pattern
153 ret = Eigen::Matrix<Scalar<G>, Dof<G>, Dof<G> * Dof<G>>::Ones().sparseView();
154 }
155 ret.makeCompressed();
156 return ret;
157}();
158
173template<LieGroup G, bool Inv = false>
174void d2r_exp_sparse(Eigen::SparseMatrix<Scalar<G>> & sp, const Tangent<G> & a, Eigen::Index i0 = 0);
175
189template<LieGroup G>
190inline void d2r_expinv_sparse(Eigen::SparseMatrix<Scalar<G>> & sp, const Tangent<G> & a, Eigen::Index i0 = 0);
191
192SMOOTH_END_NAMESPACE
193
194#include "detail/lie_group_sparse_impl.hpp"
Meta header to include all LieGroup concept model specifications.
void dr_exp_sparse(Eigen::SparseMatrix< Scalar< G > > &sp, const Tangent< G > &a, Eigen::Index i0=0)
Sparse dr_exp.
Eigen::SparseMatrix< Scalar< G > > d_exp_sparse_pattern
Sparsity pattern of dr_exp_sparse, dr_expinv_sparse (inline variable).
void d2r_expinv_sparse(Eigen::SparseMatrix< Scalar< G > > &sp, const Tangent< G > &a, Eigen::Index i0=0)
Sparse d2r_expinv.
void ad_sparse(Eigen::SparseMatrix< Scalar< G > > &sp, const Tangent< G > &a)
Sparse ad.
void dr_expinv_sparse(Eigen::SparseMatrix< Scalar< G > > &sp, const Tangent< G > &a, Eigen::Index i0=0)
Sparse dr_expinv.
std::array< Eigen::SparseMatrix< Scalar< G > >, Dof< G > > generators_sparse
Generators of lie algebra (inline variable).
Eigen::SparseMatrix< Scalar< G > > ad_sparse_pattern
Sparsity pattern of ad_sparse (inline variable).
void d2r_exp_sparse(Eigen::SparseMatrix< Scalar< G > > &sp, const Tangent< G > &a, Eigen::Index i0=0)
Sparse d2r_exp.
Eigen::SparseMatrix< Scalar< G > > d2_exp_sparse_pattern
Sparsity pattern of d2r_exp_sparse(), d2r_expinv_sparse() (inline variable).
Eigen::Vector< Scalar< M >, Dof< M > > Tangent
Tangent as a Dof-length vector.
Definition manifold.hpp:106
typename traits::man< M >::Scalar Scalar
Manifold scalar type.
Definition manifold.hpp:88
Traits that defines sparsity patterns for various groups.