gtsam 4.2.0
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#include <string>
25#include <vector>
26
27namespace gtsam {
28
40 const Values& values,
41 double delta = 1e-5) {
42 // We will fill a vector of key/Jacobians pairs (a map would sort)
43 std::vector<std::pair<Key, Matrix> > jacobians;
44
45 // Get size
46 const Vector e = factor.whitenedError(values);
47 const size_t rows = e.size();
48
49 // Loop over all variables
50 const double one_over_2delta = 1.0 / (2.0 * delta);
51 for (Key key : factor) {
52 // Compute central differences using the values struct.
53 VectorValues dX = values.zeroVectors();
54 const size_t cols = dX.dim(key);
55 Matrix J = Matrix::Zero(rows, cols);
56 for (size_t col = 0; col < cols; ++col) {
57 Vector dx = Vector::Zero(cols);
58 dx(col) = delta;
59 dX[key] = dx;
60 Values eval_values = values.retract(dX);
61 const Vector left = factor.whitenedError(eval_values);
62 dx(col) = -delta;
63 dX[key] = dx;
64 eval_values = values.retract(dX);
65 const Vector right = factor.whitenedError(eval_values);
66 J.col(col) = (left - right) * one_over_2delta;
67 }
68 jacobians.emplace_back(key, J);
69 }
70
71 // Next step...return JacobianFactor
72 return JacobianFactor(jacobians, -e);
73}
74
75namespace internal {
76// CPPUnitLite-style test for linearization of a factor
77inline bool testFactorJacobians(const std::string& name_,
78 const NoiseModelFactor& factor,
79 const gtsam::Values& values, double delta,
80 double tolerance) {
81 // Create expected value by numerical differentiation
82 JacobianFactor expected = linearizeNumerically(factor, values, delta);
83
84 // Create actual value by linearize
85 auto actual =
86 boost::dynamic_pointer_cast<JacobianFactor>(factor.linearize(values));
87 if (!actual) return false;
88
89 // Check cast result and then equality
90 bool equal = assert_equal(expected, *actual, tolerance);
91
92 // if not equal, test individual jacobians:
93 if (!equal) {
94 for (size_t i = 0; i < actual->size(); i++) {
95 bool i_good =
96 assert_equal((Matrix)(expected.getA(expected.begin() + i)),
97 (Matrix)(actual->getA(actual->begin() + i)), tolerance);
98 if (!i_good) {
99 std::cout << "Mismatch in Jacobian " << i + 1
100 << " (base 1), as shown above" << std::endl;
101 }
102 }
103 }
104
105 return equal;
106}
107} // namespace internal
108
114#define EXPECT_CORRECT_FACTOR_JACOBIANS(factor, values, numerical_derivative_step, tolerance) \
115 { EXPECT(gtsam::internal::testFactorJacobians(name_, factor, values, numerical_derivative_step, tolerance)); }
116
117} // 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:43
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:39
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:100
A Gaussian factor in the squared-error form.
Definition: JacobianFactor.h:91
VectorValues represents a collection of vector-valued variables associated each with a unique integer...
Definition: VectorValues.h:74
size_t dim(Key j) const
Return the dimension of variable j.
Definition: VectorValues.h:130
A nonlinear sum-of-squares factor with a zero-mean noise model implementing the density Templated on...
Definition: NonlinearFactor.h:174
Vector whitenedError(const Values &c) const
Vector of errors, whitened This is the raw error, i.e., i.e.
Definition: NonlinearFactor.cpp:109
A non-templated config holding any types of Manifold-group elements.
Definition: Values.h:65
Values retract(const VectorValues &delta) const
Add a delta config to current config and returns a new config.
Definition: Values.cpp:99
VectorValues zeroVectors() const
Return a VectorValues of zero vectors for each variable in this Values.
Definition: Values.cpp:272
In Gaussian factors, the error function returns either the negative log-likelihood,...
noise model to the factor, and calculates the error by asking the user to implement the method