41 double delta = 1e-5) {
43 std::vector<std::pair<Key, Matrix> > jacobians;
47 const size_t rows = e.size();
50 const double one_over_2delta = 1.0 / (2.0 * delta);
51 for (
Key key : factor) {
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);
64 eval_values = values.
retract(dX);
66 J.col(col) = (left - right) * one_over_2delta;
68 jacobians.emplace_back(key, J);
77inline bool testFactorJacobians(
const std::string& name_,
86 boost::dynamic_pointer_cast<JacobianFactor>(factor.linearize(values));
87 if (!actual)
return false;
94 for (
size_t i = 0; i < actual->size(); i++) {
96 assert_equal((Matrix)(expected.getA(expected.begin() + i)),
97 (Matrix)(actual->getA(actual->begin() + i)), tolerance);
99 std::cout <<
"Mismatch in Jacobian " << i + 1
100 <<
" (base 1), as shown above" << std::endl;
114#define EXPECT_CORRECT_FACTOR_JACOBIANS(factor, values, numerical_derivative_step, tolerance) \
115 { EXPECT(gtsam::internal::testFactorJacobians(name_, factor, values, numerical_derivative_step, tolerance)); }
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