testInvDepthFactor3.cpp
Go to the documentation of this file.
1 /* ----------------------------------------------------------------------------
2  * GTSAM Copyright 2010, Georgia Tech Research Corporation,
3  * Atlanta, Georgia 30332-0415
4  * All Rights Reserved
5  * Authors: Frank Dellaert, et al. (see THANKS for the full author list)
6  * See LICENSE for the license information
7  * -------------------------------------------------------------------------- */
8 
19 
24 #include <gtsam/inference/Symbol.h>
26 
28 
29 using namespace std;
30 using namespace gtsam;
31 
32 static Cal3_S2::shared_ptr K(new Cal3_S2(1500, 1200, 0, 640, 480));
33 static SharedNoiseModel sigma(noiseModel::Unit::Create(2));
34 
35 // camera pose at (0,0,1) looking straight along the x-axis.
36 Pose3 level_pose = Pose3(Rot3::Ypr(-M_PI/2, 0., -M_PI/2), gtsam::Point3(0,0,1));
38 
41 
42 Matrix factorError(const Pose3& pose, const Vector5& point, double invDepth,
43  const InverseDepthFactor& factor) {
44  return factor.evaluateError(pose, point, invDepth);
45 }
46 
47 /* ************************************************************************* */
48 TEST( InvDepthFactor, optimize) {
49 
50  // landmark 5 meters infront of camera (camera center at (0,0,1))
51  Point3 landmark(5, 0, 1);
52 
53  // get expected projection using pinhole camera
54  Point2 expected_uv = level_camera.project(landmark);
55 
56  InvDepthCamera3<Cal3_S2> inv_camera(level_pose, K);
57  Vector5 inv_landmark((Vector(5) << 0., 0., 1., 0., 0.).finished());
58  // initialize inverse depth with "incorrect" depth of 1/4
59  // in reality this is 1/5, but initial depth is guessed
60  double inv_depth(1./4);
61 
64 
65  graph.emplace_shared<InverseDepthFactor>(expected_uv, sigma,
66  Symbol('x',1), Symbol('l',1), Symbol('d',1), K);
67  graph.emplace_shared<PoseConstraint>(Symbol('x', 1), level_pose);
68  initial.insert(Symbol('x',1), level_pose);
69  initial.insert(Symbol('l',1), inv_landmark);
70  initial.insert(Symbol('d',1), inv_depth);
71 
73  Values result = LevenbergMarquardtOptimizer(graph, initial, lmParams).optimize();
74 
75  // with a single factor the incorrect initialization of 1/4 should not move!
76  EXPECT(assert_equal(initial, result, 1e-9));
77 
79 
80  // add a camera 2 meters to the right
81  Pose3 right_pose = level_pose * Pose3(Rot3(), Point3(2,0,0));
82  PinholeCamera<Cal3_S2> right_camera(right_pose, *K);
83 
84  // projection measurement of landmark into right camera
85  // this measurement disagrees with the depth initialization
86  // and will push it to 1/5
87  Point2 right_uv = right_camera.project(landmark);
88 
90  Symbol('x',2), Symbol('l',1),Symbol('d',1),K));
91  graph.push_back(factor1);
92 
93  graph.emplace_shared<PoseConstraint>(Symbol('x',2),right_pose);
94 
95  initial.insert(Symbol('x',2), right_pose);
96 
97  Values result2 = LevenbergMarquardtOptimizer(graph, initial, lmParams).optimize();
98 
100  result2.at<Vector5>(Symbol('l',1)),
101  result2.at<double>(Symbol('d',1)));
102  EXPECT(assert_equal(landmark, result2_lmk, 1e-9));
103 
104  // TODO: need to add priors to make this work with
105  // Values result2 = optimize<NonlinearFactorGraph>(graph, initial,
106  // NonlinearOptimizationParameters(),MULTIFRONTAL, GN);
107 }
108 
109 
110 /* ************************************************************************* */
111 TEST( InvDepthFactor, Jacobian3D ) {
112 
113  // landmark 5 meters infront of camera (camera center at (0,0,1))
114  Point3 landmark(5, 0, 1);
115 
116  // get expected projection using pinhole camera
117  Point2 expected_uv = level_camera.project(landmark);
118 
119  // get expected landmark representation using backprojection
120  InvDepthCamera3<Cal3_S2> inv_camera(level_pose, K);
121  const auto [inv_landmark, inv_depth] = inv_camera.backproject(expected_uv, 5);
122  Vector5 expected_inv_landmark((Vector(5) << 0., 0., 1., 0., 0.).finished());
123 
124  CHECK(assert_equal(expected_inv_landmark, inv_landmark, 1e-6));
125  CHECK(assert_equal(inv_depth, 1./5, 1e-6));
126 
127  Symbol poseKey('x',1);
128  Symbol pointKey('l',1);
129  Symbol invDepthKey('d',1);
130  InverseDepthFactor factor(expected_uv, sigma, poseKey, pointKey, invDepthKey, K);
131 
132  std::vector<Matrix> actualHs(3);
133  factor.unwhitenedError({{poseKey, genericValue(level_pose)},
134  {pointKey, genericValue(inv_landmark)},
135  {invDepthKey,genericValue(inv_depth)}},
136  actualHs);
137 
138  const Matrix& H1Actual = actualHs.at(0);
139  const Matrix& H2Actual = actualHs.at(1);
140  const Matrix& H3Actual = actualHs.at(2);
141 
142  // Use numerical derivatives to verify the Jacobians
143  Matrix H1Expected, H2Expected, H3Expected;
144 
145  std::function<Matrix(const Pose3 &, const Vector5 &, const double &)>
146  func = std::bind(&factorError, std::placeholders::_1, std::placeholders::_2, std::placeholders::_3, factor);
147  H1Expected = numericalDerivative31(func, level_pose, inv_landmark, inv_depth);
148  H2Expected = numericalDerivative32(func, level_pose, inv_landmark, inv_depth);
149  H3Expected = numericalDerivative33(func, level_pose, inv_landmark, inv_depth);
150 
151  // Verify the Jacobians
152  CHECK(assert_equal(H1Expected, H1Actual, 1e-6))
153  CHECK(assert_equal(H2Expected, H2Actual, 1e-6))
154  CHECK(assert_equal(H3Expected, H3Actual, 1e-6))
155 }
156 
157 /* ************************************************************************* */
158 int main() { TestResult tr; return TestRegistry::runAllTests(tr);}
159 /* ************************************************************************* */
const gtsam::Key poseKey
#define CHECK(condition)
Definition: Test.h:108
Pose3 level_pose
virtual const Values & optimize()
IsDerived< DERIVEDFACTOR > emplace_shared(Args &&... args)
Emplace a shared pointer to factor of given type.
Definition: FactorGraph.h:196
static int runAllTests(TestResult &result)
Factor Graph consisting of non-linear factors.
const ValueType at(Key j) const
Definition: Values-inl.h:204
PinholeCamera< Cal3_S2 > level_camera(level_pose, *K)
Vector2 Point2
Definition: Point2.h:32
IsDerived< DERIVEDFACTOR > push_back(std::shared_ptr< DERIVEDFACTOR > factor)
Add a factor directly using a shared_ptr.
Definition: FactorGraph.h:190
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:40
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
#define M_PI
Definition: main.h:106
Values initial
int main()
static Cal3_S2::shared_ptr K(new Cal3_S2(1500, 1200, 0, 640, 480))
Definition: BFloat16.h:88
Some functions to compute numerical derivatives.
NonlinearFactorGraph graph
Rot3 is a 3D rotation represented as a rotation matrix if the preprocessor symbol GTSAM_USE_QUATERNIO...
Definition: Rot3.h:58
const gtsam::Key pointKey
static SharedNoiseModel sigma(noiseModel::Unit::Create(2))
static Point3 landmark(0, 0, 5)
NonlinearEquality< Pose3 > PoseConstraint
Base class for all pinhole cameras.
Point2 project(const Point3 &pw, OptionalJacobian< 2, 6 > Dpose={}, OptionalJacobian< 2, 3 > Dpoint={}, OptionalJacobian< 2, DimK > Dcal={}) const
project a 3D point from world coordinates into the image
Definition: PinholePose.h:112
internal::FixedSizeMatrix< Y, X3 >::type numericalDerivative33(std::function< Y(const X1 &, const X2 &, const X3 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
Matrix factorError(const Pose3 &pose, const Vector5 &point, double invDepth, const InverseDepthFactor &factor)
std::shared_ptr< Cal3_S2 > shared_ptr
Definition: Cal3_S2.h:39
Values result
std::pair< Vector5, double > backproject(const gtsam::Point2 &pi, const double depth) const
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative32(std::function< Y(const X1 &, const X2 &, const X3 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
#define EXPECT(condition)
Definition: Test.h:150
LevenbergMarquardtParams lmParams
internal::FixedSizeMatrix< Y, X1 >::type numericalDerivative31(std::function< Y(const X1 &, const X2 &, const X3 &)> h, const X1 &x1, const X2 &x2, const X3 &x3, double delta=1e-5)
A nonlinear optimizer that uses the Levenberg-Marquardt trust-region scheme.
Array< double, 1, 3 > e(1./3., 0.5, 2.)
The most common 5DOF 3D->2D calibration.
Definition: Cal3_S2.h:34
Inverse Depth Factor based on Civera09tro, Montiel06rss. Landmarks are initialized from the first cam...
std::shared_ptr< This > shared_ptr
shorthand for a smart pointer to a factor
traits
Definition: chartTesting.h:28
Vector unwhitenedError(const Values &x, OptionalMatrixVecType H=nullptr) const override
GenericValue< T > genericValue(const T &v)
Definition: GenericValue.h:211
TEST(InvDepthFactor, optimize)
static const Pose3 pose(Rot3(Vector3(1, -1, -1).asDiagonal()), Point3(0, 0, 0.5))
int optimize(const SfmData &db, const NonlinearFactorGraph &graph, const Values &initial, bool separateCalibration=false)
Definition: timeSFMBAL.h:64
void insert(Key j, const Value &val)
Definition: Values.cpp:155
Vector3 Point3
Definition: Point3.h:38
Vector evaluateError(const POSE &pose, const Vector5 &point, const INVDEPTH &invDepth, OptionalMatrixType H1, OptionalMatrixType H2, OptionalMatrixType H3) const override
Evaluate error h(x)-z and optionally derivatives.
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
InvDepthFactor3< Pose3, Vector5, double > InverseDepthFactor
noiseModel::Base::shared_ptr SharedNoiseModel
Definition: NoiseModel.h:741


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:38:28