gtsam 4.1.1
gtsam
factorTesting.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
20#pragma once
21
24
25namespace gtsam {
26
38 const Values& values, double delta = 1e-5) {
39 // We will fill a vector of key/Jacobians pairs (a map would sort)
40 std::vector<std::pair<Key, Matrix> > jacobians;
41
42 // Get size
43 const Eigen::VectorXd e = factor.whitenedError(values);
44 const size_t rows = e.size();
45
46 // Loop over all variables
47 const double one_over_2delta = 1.0 / (2.0 * delta);
48 for (Key key : factor) {
49 // Compute central differences using the values struct.
50 VectorValues dX = values.zeroVectors();
51 const size_t cols = dX.dim(key);
52 Matrix J = Matrix::Zero(rows, cols);
53 for (size_t col = 0; col < cols; ++col) {
54 Eigen::VectorXd dx = Eigen::VectorXd::Zero(cols);
55 dx(col) = delta;
56 dX[key] = dx;
57 Values eval_values = values.retract(dX);
58 const Eigen::VectorXd left = factor.whitenedError(eval_values);
59 dx(col) = -delta;
60 dX[key] = dx;
61 eval_values = values.retract(dX);
62 const Eigen::VectorXd right = factor.whitenedError(eval_values);
63 J.col(col) = (left - right) * one_over_2delta;
64 }
65 jacobians.push_back(std::make_pair(key, J));
66 }
67
68 // Next step...return JacobianFactor
69 return JacobianFactor(jacobians, -e);
70}
71
72namespace internal {
73// CPPUnitLite-style test for linearization of a factor
74bool testFactorJacobians(const std::string& name_,
75 const NoiseModelFactor& factor, const gtsam::Values& values, double delta,
76 double tolerance) {
77
78 // Create expected value by numerical differentiation
79 JacobianFactor expected = linearizeNumerically(factor, values, delta);
80
81 // Create actual value by linearize
82 boost::shared_ptr<JacobianFactor> actual = //
83 boost::dynamic_pointer_cast<JacobianFactor>(factor.linearize(values));
84 if (!actual) return false;
85
86 // Check cast result and then equality
87 bool equal = assert_equal(expected, *actual, tolerance);
88
89 // if not equal, test individual jacobians:
90 if (!equal) {
91 for (size_t i = 0; i < actual->size(); i++) {
92 bool i_good = assert_equal((Matrix) (expected.getA(expected.begin() + i)),
93 (Matrix) (actual->getA(actual->begin() + i)), tolerance);
94 if (!i_good) {
95 std::cout << "Mismatch in Jacobian " << i+1 << " (base 1), as shown above" << std::endl;
96 }
97 }
98 }
99
100 return equal;
101}
102}
103
109#define EXPECT_CORRECT_FACTOR_JACOBIANS(factor, values, numerical_derivative_step, tolerance) \
110 { EXPECT(gtsam::internal::testFactorJacobians(name_, factor, values, numerical_derivative_step, tolerance)); }
111
112} // namespace gtsam
Some functions to compute numerical derivatives.
Non-linear factor base classes.
Global functions in a separate testing namespace.
Definition: chartTesting.h:28
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
equals with an tolerance, prints out message if unequal
Definition: Matrix.cpp:42
JacobianFactor linearizeNumerically(const NoiseModelFactor &factor, const Values &values, double delta=1e-5)
Linearize a nonlinear factor using numerical differentiation The benefit of this method is that it do...
Definition: factorTesting.h:37
bool equal(const T &obj1, const T &obj2, double tol)
Call equal on the object.
Definition: Testable.h:84
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:69
A Gaussian factor in the squared-error form.
Definition: JacobianFactor.h:91
This class represents a collection of vector-valued variables associated each with a unique integer i...
Definition: VectorValues.h:74
size_t dim(Key j) const
Return the dimension of variable j.
Definition: VectorValues.h:128
A nonlinear sum-of-squares factor with a zero-mean noise model implementing the density Templated on...
Definition: NonlinearFactor.h:162
Vector whitenedError(const Values &c) const
Vector of errors, whitened This is the raw error, i.e., i.e.
Definition: NonlinearFactor.cpp:98
A non-templated config holding any types of Manifold-group elements.
Definition: Values.h:63
VectorValues zeroVectors() const
Return a VectorValues of zero vectors for each variable in this Values.
Definition: Values.cpp:208
Values retract(const VectorValues &delta) const
Add a delta config to current config and returns a new config.
Definition: Values.cpp:101