testSO3.cpp
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 
19 #include <gtsam/base/Testable.h>
20 #include <gtsam/base/testLie.h>
21 #include <gtsam/geometry/SO3.h>
22 
23 using namespace std::placeholders;
24 using namespace std;
25 using namespace gtsam;
26 
27 //******************************************************************************
28 TEST(SO3, Identity) {
29  const SO3 R;
30  EXPECT_LONGS_EQUAL(3, R.rows());
31  EXPECT_LONGS_EQUAL(3, SO3::dimension);
32  EXPECT_LONGS_EQUAL(3, SO3::Dim());
33  EXPECT_LONGS_EQUAL(3, R.dim());
35 }
36 
37 //******************************************************************************
38 TEST(SO3, Concept) {
40  GTSAM_CONCEPT_ASSERT(IsManifold<SO3>);
42 }
43 
44 //******************************************************************************
45 TEST(SO3, Constructors) {
46  const Vector3 axis(0, 0, 1);
47  const double angle = 2.5;
48  SO3 Q(Eigen::AngleAxisd(angle, axis));
49  SO3 R = SO3::AxisAngle(axis, angle);
51 }
52 
53 /* ************************************************************************* */
54 TEST(SO3, ClosestTo) {
55  Matrix3 M;
56  M << 0.79067393, 0.6051136, -0.0930814, //
57  0.4155925, -0.64214347, -0.64324489, //
58  -0.44948549, 0.47046326, -0.75917576;
59 
60  Matrix expected(3, 3);
61  expected << 0.790687, 0.605096, -0.0931312, //
62  0.415746, -0.642355, -0.643844, //
63  -0.449411, 0.47036, -0.759468;
64 
65  auto actual = SO3::ClosestTo(3 * M);
66  EXPECT(assert_equal(expected, actual.matrix(), 1e-6));
67 }
68 
69 //******************************************************************************
70 namespace {
71 SO3 id;
72 Vector3 z_axis(0, 0, 1), v2(1, 2, 0), v3(1, 2, 3);
73 SO3 R1(Eigen::AngleAxisd(0.1, z_axis));
74 SO3 R2(Eigen::AngleAxisd(0.2, z_axis));
75 } // namespace
76 
77 /* ************************************************************************* */
78 TEST(SO3, ChordalMean) {
79  std::vector<SO3> rotations = {R1, R1.inverse()};
80  EXPECT(assert_equal(SO3(), SO3::ChordalMean(rotations)));
81 }
82 
83 //******************************************************************************
84 // Check that Hat specialization is equal to dynamic version
85 TEST(SO3, Hat) {
86  EXPECT(assert_equal(SO3::Hat(z_axis), SOn::Hat(z_axis)));
87  EXPECT(assert_equal(SO3::Hat(v2), SOn::Hat(v2)));
88  EXPECT(assert_equal(SO3::Hat(v3), SOn::Hat(v3)));
89 }
90 
91 //******************************************************************************
92 // Check that Hat specialization is equal to dynamic version
93 TEST(SO3, Vee) {
94  auto X1 = SOn::Hat(z_axis), X2 = SOn::Hat(v2), X3 = SOn::Hat(v3);
95  EXPECT(assert_equal(SO3::Vee(X1), SOn::Vee(X1)));
96  EXPECT(assert_equal(SO3::Vee(X2), SOn::Vee(X2)));
97  EXPECT(assert_equal(SO3::Vee(X3), SOn::Vee(X3)));
98 }
99 
100 //******************************************************************************
101 TEST(SO3, Local) {
102  Vector3 expected(0, 0, 0.1);
103  Vector3 actual = traits<SO3>::Local(R1, R2);
104  EXPECT(assert_equal((Vector)expected, actual));
105 }
106 
107 //******************************************************************************
108 TEST(SO3, Retract) {
109  Vector3 v(0, 0, 0.1);
110  SO3 actual = traits<SO3>::Retract(R1, v);
111  EXPECT(assert_equal(R2, actual));
112 }
113 
114 //******************************************************************************
115 TEST(SO3, Logmap) {
116  Vector3 expected(0, 0, 0.1);
117  Vector3 actual = SO3::Logmap(R1.between(R2));
118  EXPECT(assert_equal((Vector)expected, actual));
119 }
120 
121 //******************************************************************************
123  Vector3 v(0, 0, 0.1);
124  SO3 actual = R1 * SO3::Expmap(v);
125  EXPECT(assert_equal(R2, actual));
126 }
127 
128 //******************************************************************************
129 TEST(SO3, Invariants) {
130  EXPECT(check_group_invariants(id, id));
131  EXPECT(check_group_invariants(id, R1));
132  EXPECT(check_group_invariants(R2, id));
133  EXPECT(check_group_invariants(R2, R1));
134 
135  EXPECT(check_manifold_invariants(id, id));
136  EXPECT(check_manifold_invariants(id, R1));
137  EXPECT(check_manifold_invariants(R2, id));
138  EXPECT(check_manifold_invariants(R2, R1));
139 }
140 
141 //******************************************************************************
142 TEST(SO3, LieGroupDerivatives) {
147 }
148 
149 //******************************************************************************
150 TEST(SO3, ChartDerivatives) {
151  CHECK_CHART_DERIVATIVES(id, id);
155 }
156 
157 /* ************************************************************************* */
158 TEST(SO3, ExpmapFunctor) {
159  Vector axis = Vector3(0., 1., 0.); // rotation around Y
160  double angle = 3.14 / 4.0;
161  Matrix expected(3,3);
162  expected << 0.707388, 0, 0.706825, 0, 1, 0, -0.706825, 0, 0.707388;
163 
164  // axis angle version
165  so3::ExpmapFunctor f1(axis, angle);
166  SO3 actual1 = f1.expmap();
167  CHECK(assert_equal(expected, actual1.matrix(), 1e-5));
168 
169  // axis angle version, negative angle
170  so3::ExpmapFunctor f2(axis, angle - 2*M_PI);
171  SO3 actual2 = f2.expmap();
172  CHECK(assert_equal(expected, actual2.matrix(), 1e-5));
173 
174  // omega version
175  so3::ExpmapFunctor f3(axis * angle);
176  SO3 actual3 = f3.expmap();
177  CHECK(assert_equal(expected, actual3.matrix(), 1e-5));
178 
179  // omega version, negative angle
180  so3::ExpmapFunctor f4(axis * (angle - 2*M_PI));
181  SO3 actual4 = f4.expmap();
182  CHECK(assert_equal(expected, actual4.matrix(), 1e-5));
183 }
184 
185 /* ************************************************************************* */
186 namespace exmap_derivative {
187 static const Vector3 w(0.1, 0.27, -0.2);
188 }
189 // Left trivialized Derivative of exp(w) wrpt w:
190 // How does exp(w) change when w changes?
191 // We find a y such that: exp(w) exp(y) = exp(w + dw) for dw --> 0
192 // => y = log (exp(-w) * exp(w+dw))
194  using exmap_derivative::w;
195  return SO3::Logmap(SO3::Expmap(-w) * SO3::Expmap(w + dw));
196 }
197 
198 TEST(SO3, ExpmapDerivative) {
199  using exmap_derivative::w;
200  const Matrix actualDexpL = SO3::ExpmapDerivative(w);
201  const Matrix expectedDexpL =
202  numericalDerivative11<Vector3, Vector3>(testDexpL, Vector3::Zero(), 1e-2);
203  EXPECT(assert_equal(expectedDexpL, actualDexpL, 1e-7));
204 
205  const Matrix actualDexpInvL = SO3::LogmapDerivative(w);
206  EXPECT(assert_equal(expectedDexpL.inverse(), actualDexpInvL, 1e-7));
207 }
208 
209 //******************************************************************************
210 TEST(SO3, ExpmapDerivative2) {
211  const Vector3 theta(0.1, 0, 0.1);
212  const Matrix Jexpected = numericalDerivative11<SO3, Vector3>(
213  std::bind(&SO3::Expmap, std::placeholders::_1, nullptr), theta);
214 
215  CHECK(assert_equal(Jexpected, SO3::ExpmapDerivative(theta)));
216  CHECK(assert_equal(Matrix3(Jexpected.transpose()),
217  SO3::ExpmapDerivative(-theta)));
218 }
219 
220 //******************************************************************************
221 TEST(SO3, ExpmapDerivative3) {
222  const Vector3 theta(10, 20, 30);
223  const Matrix Jexpected = numericalDerivative11<SO3, Vector3>(
224  std::bind(&SO3::Expmap, std::placeholders::_1, nullptr), theta);
225 
226  CHECK(assert_equal(Jexpected, SO3::ExpmapDerivative(theta)));
227  CHECK(assert_equal(Matrix3(Jexpected.transpose()),
228  SO3::ExpmapDerivative(-theta)));
229 }
230 
231 //******************************************************************************
232 TEST(SO3, ExpmapDerivative4) {
233  // Iserles05an (Lie-group Methods) says:
234  // scalar is easy: d exp(a(t)) / dt = exp(a(t)) a'(t)
235  // matrix is hard: d exp(A(t)) / dt = exp(A(t)) dexp[-A(t)] A'(t)
236  // where A(t): R -> so(3) is a trajectory in the tangent space of SO(3)
237  // and dexp[A] is a linear map from 3*3 to 3*3 derivatives of se(3)
238  // Hence, the above matrix equation is typed: 3*3 = SO(3) * linear_map(3*3)
239 
240  // In GTSAM, we don't work with the skew-symmetric matrices A directly, but
241  // with 3-vectors.
242  // omega is easy: d Expmap(w(t)) / dt = ExmapDerivative[w(t)] * w'(t)
243 
244  // Let's verify the above formula.
245 
246  auto w = [](double t) { return Vector3(2 * t, sin(t), 4 * t * t); };
247  auto w_dot = [](double t) { return Vector3(2, cos(t), 8 * t); };
248 
249  // We define a function R
250  auto R = [w](double t) { return SO3::Expmap(w(t)); };
251 
252  for (double t = -2.0; t < 2.0; t += 0.3) {
253  const Matrix expected = numericalDerivative11<SO3, double>(R, t);
254  const Matrix actual = SO3::ExpmapDerivative(w(t)) * w_dot(t);
255  CHECK(assert_equal(expected, actual, 1e-7));
256  }
257 }
258 
259 //******************************************************************************
260 TEST(SO3, ExpmapDerivative5) {
261  auto w = [](double t) { return Vector3(2 * t, sin(t), 4 * t * t); };
262  auto w_dot = [](double t) { return Vector3(2, cos(t), 8 * t); };
263 
264  // Now define R as mapping local coordinates to neighborhood around R0.
265  const SO3 R0 = SO3::Expmap(Vector3(0.1, 0.4, 0.2));
266  auto R = [R0, w](double t) { return R0.expmap(w(t)); };
267 
268  for (double t = -2.0; t < 2.0; t += 0.3) {
269  const Matrix expected = numericalDerivative11<SO3, double>(R, t);
270  const Matrix actual = SO3::ExpmapDerivative(w(t)) * w_dot(t);
271  CHECK(assert_equal(expected, actual, 1e-7));
272  }
273 }
274 
275 //******************************************************************************
276 TEST(SO3, ExpmapDerivative6) {
277  const Vector3 thetahat(0.1, 0, 0.1);
278  const Matrix Jexpected = numericalDerivative11<SO3, Vector3>(
279  std::bind(&SO3::Expmap, std::placeholders::_1, nullptr), thetahat);
280  Matrix3 Jactual;
281  SO3::Expmap(thetahat, Jactual);
282  EXPECT(assert_equal(Jexpected, Jactual));
283 }
284 
285 //******************************************************************************
286 TEST(SO3, LogmapDerivative) {
287  const Vector3 thetahat(0.1, 0, 0.1);
288  const SO3 R = SO3::Expmap(thetahat); // some rotation
289  const Matrix Jexpected = numericalDerivative11<Vector, SO3>(
290  std::bind(&SO3::Logmap, std::placeholders::_1, nullptr), R);
291  const Matrix3 Jactual = SO3::LogmapDerivative(thetahat);
292  EXPECT(assert_equal(Jexpected, Jactual));
293 }
294 
295 //******************************************************************************
296 TEST(SO3, JacobianLogmap) {
297  const Vector3 thetahat(0.1, 0, 0.1);
298  const SO3 R = SO3::Expmap(thetahat); // some rotation
299  const Matrix Jexpected = numericalDerivative11<Vector, SO3>(
300  std::bind(&SO3::Logmap, std::placeholders::_1, nullptr), R);
301  Matrix3 Jactual;
302  SO3::Logmap(R, Jactual);
303  EXPECT(assert_equal(Jexpected, Jactual));
304 }
305 
306 //******************************************************************************
307 namespace test_cases {
308 std::vector<Vector3> small{{0, 0, 0}, //
309  {1e-5, 0, 0}, {0, 1e-5, 0}, {0, 0, 1e-5}, //,
310  {1e-4, 0, 0}, {0, 1e-4, 0}, {0, 0, 1e-4}};
311 std::vector<Vector3> large{
312  {0, 0, 0}, {1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {0.1, 0.2, 0.3}};
313 auto omegas = [](bool nearZero) { return nearZero ? small : large; };
314 std::vector<Vector3> vs{{1, 0, 0}, {0, 1, 0}, {0, 0, 1}, {0.4, 0.3, 0.2}};
315 } // namespace test_cases
316 
317 //******************************************************************************
318 TEST(SO3, JacobianInverses) {
319  Matrix HR, HL;
320  for (bool nearZero : {true, false}) {
321  for (const Vector3& omega : test_cases::omegas(nearZero)) {
322  so3::DexpFunctor local(omega, nearZero);
323  EXPECT(assert_equal<Matrix3>(local.rightJacobian().inverse(),
324  local.rightJacobianInverse()));
325  EXPECT(assert_equal<Matrix3>(local.leftJacobian().inverse(),
326  local.leftJacobianInverse()));
327  }
328  }
329 }
330 
331 //******************************************************************************
332 TEST(SO3, CrossB) {
333  Matrix aH1;
334  for (bool nearZero : {true, false}) {
335  std::function<Vector3(const Vector3&, const Vector3&)> f =
336  [=](const Vector3& omega, const Vector3& v) {
337  return so3::DexpFunctor(omega, nearZero).crossB(v);
338  };
339  for (const Vector3& omega : test_cases::omegas(nearZero)) {
340  so3::DexpFunctor local(omega, nearZero);
341  for (const Vector3& v : test_cases::vs) {
342  local.crossB(v, aH1);
344  }
345  }
346  }
347 }
348 
349 //******************************************************************************
350 TEST(SO3, DoubleCrossC) {
351  Matrix aH1;
352  for (bool nearZero : {true, false}) {
353  std::function<Vector3(const Vector3&, const Vector3&)> f =
354  [=](const Vector3& omega, const Vector3& v) {
355  return so3::DexpFunctor(omega, nearZero).doubleCrossC(v);
356  };
357  for (const Vector3& omega : test_cases::omegas(nearZero)) {
358  so3::DexpFunctor local(omega, nearZero);
359  for (const Vector3& v : test_cases::vs) {
360  local.doubleCrossC(v, aH1);
362  }
363  }
364  }
365 }
366 
367 //******************************************************************************
368 TEST(SO3, ApplyDexp) {
369  Matrix aH1, aH2;
370  for (bool nearZero : {true, false}) {
371  std::function<Vector3(const Vector3&, const Vector3&)> f =
372  [=](const Vector3& omega, const Vector3& v) {
373  return so3::DexpFunctor(omega, nearZero).applyDexp(v);
374  };
375  for (const Vector3& omega : test_cases::omegas(nearZero)) {
376  so3::DexpFunctor local(omega, nearZero);
377  for (const Vector3& v : test_cases::vs) {
378  EXPECT(assert_equal(Vector3(local.dexp() * v),
379  local.applyDexp(v, aH1, aH2)));
382  EXPECT(assert_equal(local.dexp(), aH2));
383  }
384  }
385  }
386 }
387 
388 //******************************************************************************
389 TEST(SO3, ApplyInvDexp) {
390  Matrix aH1, aH2;
391  for (bool nearZero : {true, false}) {
392  std::function<Vector3(const Vector3&, const Vector3&)> f =
393  [=](const Vector3& omega, const Vector3& v) {
394  return so3::DexpFunctor(omega, nearZero).applyInvDexp(v);
395  };
396  for (const Vector3& omega : test_cases::omegas(nearZero)) {
397  so3::DexpFunctor local(omega, nearZero);
398  Matrix invJr = local.rightJacobianInverse();
399  for (const Vector3& v : test_cases::vs) {
400  EXPECT(
401  assert_equal(Vector3(invJr * v), local.applyInvDexp(v, aH1, aH2)));
404  EXPECT(assert_equal(invJr, aH2));
405  }
406  }
407  }
408 }
409 
410 //******************************************************************************
411 TEST(SO3, ApplyLeftJacobian) {
412  Matrix aH1, aH2;
413  for (bool nearZero : {true, false}) {
414  std::function<Vector3(const Vector3&, const Vector3&)> f =
415  [=](const Vector3& omega, const Vector3& v) {
416  return so3::DexpFunctor(omega, nearZero).applyLeftJacobian(v);
417  };
418  for (const Vector3& omega : test_cases::omegas(nearZero)) {
419  so3::DexpFunctor local(omega, nearZero);
420  for (const Vector3& v : test_cases::vs) {
421  EXPECT(assert_equal(Vector3(local.leftJacobian() * v),
422  local.applyLeftJacobian(v, aH1, aH2)));
425  EXPECT(assert_equal(local.leftJacobian(), aH2));
426  }
427  }
428  }
429 }
430 
431 //******************************************************************************
432 TEST(SO3, ApplyLeftJacobianInverse) {
433  Matrix aH1, aH2;
434  for (bool nearZero : {true, false}) {
435  std::function<Vector3(const Vector3&, const Vector3&)> f =
436  [=](const Vector3& omega, const Vector3& v) {
437  return so3::DexpFunctor(omega, nearZero).applyLeftJacobianInverse(v);
438  };
439  for (const Vector3& omega : test_cases::omegas(nearZero)) {
440  so3::DexpFunctor local(omega, nearZero);
441  Matrix invJl = local.leftJacobianInverse();
442  for (const Vector3& v : test_cases::vs) {
443  EXPECT(assert_equal(Vector3(invJl * v),
444  local.applyLeftJacobianInverse(v, aH1, aH2)));
447  EXPECT(assert_equal(invJl, aH2));
448  }
449  }
450  }
451 }
452 
453 //******************************************************************************
454 TEST(SO3, vec) {
455  const Vector9 expected = Eigen::Map<const Vector9>(R2.matrix().data());
456  Matrix actualH;
457  const Vector9 actual = R2.vec(actualH);
458  CHECK(assert_equal(expected, actual));
459  std::function<Vector9(const SO3&)> f = [](const SO3& Q) { return Q.vec(); };
460  const Matrix numericalH = numericalDerivative11(f, R2, 1e-5);
461  CHECK(assert_equal(numericalH, actualH));
462 }
463 
464 //******************************************************************************
466  Matrix3 M;
467  M << 1, 2, 3, 4, 5, 6, 7, 8, 9;
468  SO3 R = SO3::Expmap(Vector3(1, 2, 3));
469  const Matrix3 expected = M * R.matrix();
470  Matrix actualH;
471  const Matrix3 actual = so3::compose(M, R, actualH);
472  CHECK(assert_equal(expected, actual));
473  std::function<Matrix3(const Matrix3&)> f = [R](const Matrix3& M) {
474  return so3::compose(M, R);
475  };
476  Matrix numericalH = numericalDerivative11(f, M, 1e-2);
477  CHECK(assert_equal(numericalH, actualH));
478 }
479 
480 //******************************************************************************
481 int main() {
482  TestResult tr;
483  return TestRegistry::runAllTests(tr);
484 }
485 //******************************************************************************
gtsam::SO3
SO< 3 > SO3
Definition: SO3.h:33
TestRegistry::runAllTests
static int runAllTests(TestResult &result)
Definition: TestRegistry.cpp:27
gtsam::LieGroup::between
Class between(const Class &g) const
Definition: Lie.h:52
gtsam::SO::inverse
SO inverse() const
inverse of a rotation = transpose
Definition: SOn.h:196
gtsam::utils.numerical_derivative.X2
X2
Definition: numerical_derivative.py:26
test_constructor::f1
auto f1
Definition: testHybridNonlinearFactor.cpp:56
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
EXPECT_LONGS_EQUAL
#define EXPECT_LONGS_EQUAL(expected, actual)
Definition: Test.h:154
Testable.h
Concept check for values that can be used in unit tests.
ceres::sin
Jet< T, N > sin(const Jet< T, N > &f)
Definition: jet.h:439
EXPECT
#define EXPECT(condition)
Definition: Test.h:150
TestHarness.h
test_cases::omegas
auto omegas
Definition: testPose3.cpp:879
gtsam::numericalDerivative11
internal::FixedSizeMatrix< Y, X >::type numericalDerivative11(std::function< Y(const X &)> h, const X &x, double delta=1e-5)
New-style numerical derivatives using manifold_traits.
Definition: numericalDerivative.h:110
Eigen::AngleAxis
Represents a 3D rotation as a rotation angle around an arbitrary 3D axis.
Definition: ForwardDeclarations.h:290
gtsam::numericalDerivative22
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative22(std::function< Y(const X1 &, const X2 &)> h, const X1 &x1, const X2 &x2, double delta=1e-5)
Definition: numericalDerivative.h:195
f2
double f2(const Vector2 &x)
Definition: testNumericalDerivative.cpp:58
gtsam::Matrix
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
gtsam::Rot2::matrix
Matrix2 matrix() const
Definition: Rot2.cpp:85
test_cases::vs
std::vector< Vector3 > vs
Definition: testPose3.cpp:880
gtsam::Vector3
Eigen::Vector3d Vector3
Definition: Vector.h:44
testLie.h
gtsam::IsGroup
Definition: Group.h:42
gtsam::utils.numerical_derivative.X3
X3
Definition: numerical_derivative.py:27
gtsam::Vector
Eigen::VectorXd Vector
Definition: Vector.h:39
Q
Quaternion Q
Definition: testQuaternion.cpp:27
testDexpL
Vector3 testDexpL(const Vector3 &dw)
Definition: testSO3.cpp:193
ceres::cos
Jet< T, N > cos(const Jet< T, N > &f)
Definition: jet.h:426
exmap_derivative
Definition: testSO3.cpp:186
TEST
TEST(SO3, Identity)
Definition: testSO3.cpp:28
biased_x_rotation::omega
const double omega
Definition: testPreintegratedRotation.cpp:32
so3::R1
SO3 R1
Definition: testShonanFactor.cpp:41
SO3.h
3*3 matrix representation of SO(3)
Expmap
Pose2_ Expmap(const Vector3_ &xi)
Definition: InverseKinematicsExampleExpressions.cpp:47
simple::R0
static Rot3 R0
Definition: testInitializePose3.cpp:48
id
static const Similarity3 id
Definition: testSimilarity3.cpp:44
main
int main()
Definition: testSO3.cpp:481
gtsam::SO::matrix
const MatrixNN & matrix() const
Return matrix.
Definition: SOn.h:158
gtsam::IsLieGroup
Definition: Lie.h:260
cholesky::expected
Matrix expected
Definition: testMatrix.cpp:971
CHECK_LIE_GROUP_DERIVATIVES
#define CHECK_LIE_GROUP_DERIVATIVES(t1, t2)
Definition: testLie.h:82
test_cases::large
std::vector< Vector3 > large
Definition: testPose3.cpp:877
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
TestResult
Definition: TestResult.h:26
tree::f
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
Definition: testExpression.cpp:218
exmap_derivative::w
static const Vector3 w(0.1, 0.27, -0.2)
gtsam::utils.numerical_derivative.X1
X1
Definition: numerical_derivative.py:25
Eigen::Quaternion
The quaternion class used to represent 3D orientations and rotations.
Definition: ForwardDeclarations.h:293
so3::R2
SO3 R2
Definition: testShonanFactor.cpp:43
gtsam
traits
Definition: SFMdata.h:40
gtsam::traits
Definition: Group.h:36
f3
double f3(double x1, double x2)
Definition: testNumericalDerivative.cpp:78
test_cases::small
std::vector< Vector3 > small
Definition: testPose3.cpp:874
gtsam::testing::compose
T compose(const T &t1, const T &t2)
Definition: lieProxies.h:39
CHECK
#define CHECK(condition)
Definition: Test.h:108
gtsam::numericalDerivative21
internal::FixedSizeMatrix< Y, X1 >::type numericalDerivative21(const std::function< Y(const X1 &, const X2 &)> &h, const X1 &x1, const X2 &x2, double delta=1e-5)
Definition: numericalDerivative.h:166
std
Definition: BFloat16.h:88
v2
Vector v2
Definition: testSerializationBase.cpp:39
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
f4
double f4(double x, double y, double z)
Definition: testNumericalDerivative.cpp:107
gtsam::assert_equal
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:41
Eigen::PlainObjectBase::data
EIGEN_DEVICE_FUNC const EIGEN_STRONG_INLINE Scalar * data() const
Definition: PlainObjectBase.h:247
gtsam::SO< 3 >
v3
Vector v3
Definition: testSerializationBase.cpp:40
gtsam::LieGroup::expmap
Class expmap(const TangentVector &v) const
Definition: Lie.h:78
M_PI
#define M_PI
Definition: mconf.h:117
align_3::t
Point2 t(10, 10)
gtsam::SO::vec
VectorN2 vec(OptionalJacobian< internal::NSquaredSO(N), dimension > H={}) const
Definition: SOn-inl.h:88
test_cases
Definition: testPose3.cpp:873
GTSAM_CONCEPT_ASSERT
#define GTSAM_CONCEPT_ASSERT(concept)
Definition: base/concepts.h:22
CHECK_CHART_DERIVATIVES
#define CHECK_CHART_DERIVATIVES(t1, t2)
Definition: testLie.h:85
R
Rot2 R(Rot2::fromAngle(0.1))
HR
#define HR
Definition: mincover.c:26
M
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:51


gtsam
Author(s):
autogenerated on Wed Jan 1 2025 04:07:28