gtsam 4.1.1
gtsam
SparseEigen.h
Go to the documentation of this file.
1/* ----------------------------------------------------------------------------
2
3 * GTSAM Copyright 2010, Georgia Tech Research Corporation,
4 * Atlanta, Georgia 30332-0415
5 * All Rights Reserved
6 * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
7
8 * See LICENSE for the license information
9
10 * -------------------------------------------------------------------------- */
11
24#pragma once
25
28
29#include <Eigen/Sparse>
30
31namespace gtsam {
32
35typedef Eigen::SparseMatrix<double, Eigen::ColMajor, int> SparseEigen;
36
39 const GaussianFactorGraph &gfg, const Ordering &ordering) {
40 gttic_(SparseEigen_sparseJacobianEigen);
41 // intermediate `entries` vector is kind of unavoidable due to how expensive
42 // factor->rows() is, which prevents us from populating SparseEigen directly.
43 size_t nrows, ncols;
44 auto entries = gfg.sparseJacobian(ordering, nrows, ncols);
45 // declare sparse matrix
46 SparseEigen Ab(nrows, ncols);
47 // See Eigen::set_from_triplets. This is about 5% faster.
48 // pass 1: count the nnz per inner-vector
49 std::vector<int> nnz(ncols, 0);
50 for (const auto &entry : entries) nnz[std::get<1>(entry)]++;
51 Ab.reserve(nnz);
52 // pass 2: insert the elements
53 for (const auto &entry : entries)
54 Ab.insert(std::get<0>(entry), std::get<1>(entry)) = std::get<2>(entry);
55 return Ab;
56}
57
58SparseEigen sparseJacobianEigen(const GaussianFactorGraph &gfg) {
59 gttic_(SparseEigen_sparseJacobianEigen_defaultOrdering);
60 return sparseJacobianEigen(gfg, Ordering(gfg.keys()));
61}
62
63} // namespace gtsam
Linear Factor Graph where all factors are Gaussians.
Factor Graph Values.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
Eigen::SparseMatrix< double, Eigen::ColMajor, int > SparseEigen
Eigen-format sparse matrix.
Definition: SparseEigen.h:35
SparseEigen sparseJacobianEigen(const GaussianFactorGraph &gfg, const Ordering &ordering)
Constructs an Eigen-format SparseMatrix of a GaussianFactorGraph.
Definition: SparseEigen.h:38
Definition: Ordering.h:34
A Linear Factor Graph is a factor graph where all factors are Gaussian, i.e.
Definition: GaussianFactorGraph.h:69
std::vector< std::tuple< int, int, double > > sparseJacobian(const Ordering &ordering, size_t &nrows, size_t &ncols) const
Returns a sparse augmented Jacbian matrix as a vector of i, j, and s, where i(k) and j(k) are the bas...
Definition: GaussianFactorGraph.cpp:103