testBetweenFactor.cpp
Go to the documentation of this file.
1 
9 #include <gtsam/geometry/Rot3.h>
10 #include <gtsam/inference/Symbol.h>
13 
14 using namespace gtsam;
15 using namespace gtsam::symbol_shorthand;
16 using namespace gtsam::noiseModel;
17 
22  Rot3 R1 = Rot3::Rodrigues(0.1, 0.2, 0.3);
23  Rot3 R2 = Rot3::Rodrigues(0.4, 0.5, 0.6);
24  Rot3 noise = Rot3(); // Rot3::Rodrigues(0.01, 0.01, 0.01); // Uncomment to make unit test fail
25  Rot3 measured = R1.between(R2)*noise ;
26 
27  BetweenFactor<Rot3> factor(R(1), R(2), measured, Isotropic::Sigma(3, 0.05));
28  Matrix actualH1, actualH2;
29  Vector actual = factor.evaluateError(R1, R2, actualH1, actualH2);
30 
31  Vector expected = Rot3::Logmap(measured.inverse() * R1.between(R2));
32  EXPECT(assert_equal(expected,actual/*, 1e-100*/)); // Uncomment to make unit test fail
33 
34  Matrix numericalH1 = numericalDerivative21<Vector3,Rot3,Rot3>(
35  boost::function<Vector(const Rot3&, const Rot3&)>(boost::bind(
36  &BetweenFactor<Rot3>::evaluateError, factor, _1, _2, boost::none,
37  boost::none)), R1, R2, 1e-5);
38  EXPECT(assert_equal(numericalH1,actualH1, 1E-5));
39 
40  Matrix numericalH2 = numericalDerivative22<Vector3,Rot3,Rot3>(
41  boost::function<Vector(const Rot3&, const Rot3&)>(boost::bind(
42  &BetweenFactor<Rot3>::evaluateError, factor, _1, _2, boost::none,
43  boost::none)), R1, R2, 1e-5);
44  EXPECT(assert_equal(numericalH2,actualH2, 1E-5));
45 }
46 
47 /* ************************************************************************* */
48 /*
49 // Constructor scalar
50 TEST(BetweenFactor, ConstructorScalar) {
51  SharedNoiseModel model;
52  double measured_value = 0.0;
53  BetweenFactor<double> factor(1, 2, measured_value, model);
54 }
55 
56 // Constructor vector3
57 TEST(BetweenFactor, ConstructorVector3) {
58  SharedNoiseModel model = noiseModel::Isotropic::Sigma(3, 1.0);
59  Vector3 measured_value(1, 2, 3);
60  BetweenFactor<Vector3> factor(1, 2, measured_value, model);
61 }
62 
63 // Constructor dynamic sized vector
64 TEST(BetweenFactor, ConstructorDynamicSizeVector) {
65  SharedNoiseModel model = noiseModel::Isotropic::Sigma(5, 1.0);
66  Vector measured_value(5); measured_value << 1, 2, 3, 4, 5;
67  BetweenFactor<Vector> factor(1, 2, measured_value, model);
68 }
69 */
70 
71 /* ************************************************************************* */
72 int main() {
73  TestResult tr;
74  return TestRegistry::runAllTests(tr);
75 }
76 /* ************************************************************************* */
Key E(std::uint64_t j)
static int runAllTests(TestResult &result)
Matrix expected
Definition: testMatrix.cpp:974
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:43
Some functions to compute numerical derivatives.
Rot3 inverse() const
inverse of a rotation
Definition: Rot3.h:311
static Rot3 Rodrigues(const Vector3 &w)
Definition: Rot3.h:243
Eigen::VectorXd Vector
Definition: Vector.h:38
#define EXPECT(condition)
Definition: Test.h:151
Array< double, 1, 3 > e(1./3., 0.5, 2.)
Q R1(Eigen::AngleAxisd(1, Q_z_axis))
traits
Definition: chartTesting.h:28
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:42
int main()
static Vector3 Logmap(const Rot3 &R, OptionalJacobian< 3, 3 > H=boost::none)
Definition: Rot3M.cpp:158
Key R(std::uint64_t j)
Class between(const Class &g) const
Definition: Lie.h:51
TEST(LPInitSolver, InfiniteLoopSingleVar)
Point3 measured
Vector evaluateError(const T &p1, const T &p2, boost::optional< Matrix & > H1=boost::none, boost::optional< Matrix & > H2=boost::none) const override
evaluate error, returns vector of errors size of tangent space
Q R2(Eigen::AngleAxisd(2, Vector3(0, 1, 0)))
3D rotation represented as a rotation matrix or quaternion
All noise models live in the noiseModel namespace.


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:46:18