Go to the documentation of this file.
28 using namespace std::placeholders;
30 using namespace gtsam;
43 static const Vector9
kZeroXi = Vector9::Zero();
45 static const
Point3 V(3, 0.4, -2.2);
46 static const
Point3 P(0.2, 0.7, -2);
48 static const
Point3 V2(-6.5, 3.5, 6.2);
49 static const
Point3 P2(3.5, -8.2, 4.2);
58 std::bind(&NavState::Create, std::placeholders::_1, std::placeholders::_2,
59 std::placeholders::_3,
nullptr,
nullptr,
nullptr);
78 std::bind(&NavState::FromPoseVelocity, std::placeholders::_1,
79 std::placeholders::_2,
nullptr,
nullptr);
93 eH = numericalDerivative11<Rot3, NavState>(
103 eH = numericalDerivative11<Point3, NavState>(
114 eH = numericalDerivative11<Velocity3, NavState>(
125 eH = numericalDerivative11<Velocity3, NavState>(
141 xi << 0.1, 0.1, 0.1, 0.2, 0.3, 0.4, -0.1, -0.2, -0.3;
160 std::placeholders::_2,
nullptr,
nullptr);
165 const Vector9 xi2 = -3.0*
xi;
172 std::bind(&NavState::localCoordinates, std::placeholders::_1,
173 std::placeholders::_2,
nullptr,
nullptr);
199 NavState nav_state_a(Rot3::Identity(), {0.0, 1.0, 2.0}, {1.0, -1.0, 1.0});
200 NavState nav_state_b(Rot3::Rx(
M_PI_4), {0.0, 1.0, 3.0}, {1.0, -1.0, 2.0});
201 NavState nav_state_c(Rot3::Ry(
M_PI / 180.0), {1.0, 1.0, 2.0},
204 auto ab_c = (nav_state_a * nav_state_b) * nav_state_c;
205 auto a_bc = nav_state_a * (nav_state_b * nav_state_c);
213 Matrix actualDcompose1, actualDcompose2;
235 Matrix actualDcompose1, actualDcompose2;
250 NavState nav_state_a(Rot3::Identity(), {0.0, 1.0, 2.0}, {1.0, -1.0, 1.0});
251 NavState nav_state_b(Rot3::Rx(
M_PI_4), {0.0, 1.0, 3.0}, {1.0, -1.0, 2.0});
252 NavState nav_state_c(Rot3::Ry(
M_PI / 180.0), {1.0, 1.0, 2.0},
255 auto a_inv = nav_state_a.
inverse();
256 auto a_a_inv = nav_state_a * a_inv;
259 auto b_inv = nav_state_b.inverse();
260 auto b_b_inv = nav_state_b * b_inv;
264 Matrix actual =
T.inverse(actualDinverse).matrix();
275 Rot3 R = Rot3::Rodrigues(0.3, 0.4, -0.5);
282 T.inverse(actualDinverse);
306 Matrix actualDBetween1, actualDBetween2;
307 actual =
T2.
between(
T3, actualDBetween1, actualDBetween2);
326 static const double dt = 2.0;
329 std::placeholders::_2,
nullptr);
364 double dt = 2.0, dt2 =
dt *
dt;
365 Vector3 n_omega(0.0, 0.0, 1.0), n_t(1.0, 0.0, 0.0), n_v(0.0, 1.0, 0.0);
366 Vector3 n_aCorr1 = -2.0 * n_omega.cross(n_v),
367 n_aCorr2 = -n_omega.cross(n_omega.cross(n_t));
368 Rot3 nRb =
Rot3(-1.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 1.0, 0.0),
372 Vector3 n_dP1e = 0.5 * dt2 * n_aCorr1,
373 n_dP2e = 0.5 * dt2 * (n_aCorr1 + n_aCorr2),
374 n_dV1e =
dt * n_aCorr1,
375 n_dV2e =
dt * (n_aCorr1 + n_aCorr2);
379 b_dP1e = bRn * n_dP1e, b_dP2e = bRn * n_dP2e,
380 b_dV1e = bRn * n_dV1e, b_dV2e = bRn * n_dV2e;
384 Vector9 dXi1a = kState2.coriolis(
dt, n_omega,
false),
385 dXi2a = kState2.coriolis(
dt, n_omega,
true);
386 Vector3 dRa = NavState::dR(dXi1a),
387 b_dP1a = NavState::dP(dXi1a), b_dV1a = NavState::dV(dXi1a),
388 b_dP2a = NavState::dP(dXi2a), b_dV2a = NavState::dV(dXi2a);
401 xi << 0.1, 0.1, 0.1, 0.2, 0.3, 0.4, -0.1, -0.2, -0.3;
404 std::function<Vector9(
const NavState&,
const Vector9&)> correctPIM =
405 std::bind(&NavState::correctPIM, std::placeholders::_1,
418 std::ostringstream
os;
422 expected =
"R: [\n\t1, 0, 0;\n\t0, 1, 0;\n\t0, 0, 1\n]\np: 0 0 0\nv: 0 0 0";
432 std::string
R =
"R: [\n\t1, 0, 0;\n\t0, 1, 0;\n\t0, 0, 1\n]\n";
433 std::string
p =
"p: 1 2 3\n";
434 std::string
v =
"v: 1 2 3\n";
441 #ifndef GTSAM_POSE3_EXPMAP
490 (
Vector(9) << 0.0, 0.0, 0.1, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0).finished());
529 (
Vector(9) << 0.1, 0.1, 0.1, 0.4, 0.2, 0.8, 0.4, 0.2, 0.8).finished();
539 d << 1, 2, 3, 4, 5, 6, 7, 8, 9;
541 const Rot3 R = Rot3::Retract(
d.head<3>());
548 d12 << 1, 2, 3, 4, 5, 6, 7, 8, 9;
555 Vector d12 = Vector9::Constant(0.1);
590 (
Vector(9) << 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9).finished();
611 expected << wh, Z_3x3, Z_3x3, vh, wh, Z_3x3, rh, Z_3x3, wh;
619 w << 0.1, 0.2, 0.3, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0;
622 std::function<
NavState(
const Vector9&)>
f = [](
const Vector9&
w) {
626 numericalDerivative21<NavState, Vector9, OptionalJacobian<9, 9> >(
635 w << 0.1, 0.2, 0.3, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0;
640 return NavState::Logmap(
p);
643 numericalDerivative21<Vector9, NavState, OptionalJacobian<9, 9> >(
644 &NavState::Logmap,
p, {});
652 EXPECT(check_group_invariants(
id,
id));
653 EXPECT(check_group_invariants(
id,
T3));
654 EXPECT(check_group_invariants(
T2,
id));
657 EXPECT(check_manifold_invariants(
id,
id));
658 EXPECT(check_manifold_invariants(
id,
T3));
659 EXPECT(check_manifold_invariants(
T2,
id));
static int runAllTests(TestResult &result)
ADT create(const Signature &signature)
Class between(const Class &g) const
static const Pose3 kPose(kAttitude, kPosition)
static const Pose3 pose(Rot3(Vector3(1, -1, -1).asDiagonal()), Point3(0, 0, 0.5))
const EIGEN_DEVICE_FUNC InverseReturnType inverse() const
Point3 expectedP(0.29552, 0.0446635, 1)
static const NavState T(R, P2, V2)
Array< double, 1, 3 > e(1./3., 0.5, 2.)
static const double d[K][N]
#define GTSAM_CONCEPT_TESTABLE_INST(T)
Jet< T, N > sin(const Jet< T, N > &f)
#define EXPECT(condition)
static const Velocity3 kVelocity(0.4, 0.5, 0.6)
Vector3 Velocity3
Velocity is currently typedef'd to Vector3.
bool assert_print_equal(const std::string &expected, const V &actual, const std::string &s="")
Class compose(const Class &g) const
static const Vector3 kOmegaCoriolis(0.02, 0.03, 0.04)
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
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.
static const Rot3 kAttitude
static const Point3 kPosition(1.0, 2.0, 3.0)
internal::FixedSizeMatrix< Y, X2 >::type numericalDerivative22(std::function< Y(const X1 &, const X2 &)> h, const X1 &x1, const X2 &x2, double delta=1e-5)
Matrix9 AdjointMap() const
def retract(a, np.ndarray xi)
Velocity3 bodyVelocity(OptionalJacobian< 3, 9 > H={}) const
NavState expected(expectedR, expectedV, expectedP)
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
ofstream os("timeSchurFactors.csv")
Matrix3 skewSymmetric(double wx, double wy, double wz)
int BLASFUNC() drot(int *, double *, int *, double *, int *, double *, double *)
Point3 unrotate(const Point3 &p, OptionalJacobian< 3, 3 > H1={}, OptionalJacobian< 3, 3 > H2={}) const
rotate point from world to rotated frame
Jet< T, N > cos(const Jet< T, N > &f)
Provides additional testing facilities for common data structures.
static const NavState kState1(kAttitude, kPosition, kVelocity)
Vector9 correctPIM(const Vector9 &pim, double dt, const Vector3 &n_gravity, const std::optional< Vector3 > &omegaCoriolis, bool use2ndOrderCoriolis=false, OptionalJacobian< 9, 9 > H1={}, OptionalJacobian< 9, 9 > H2={}) const
static const Point3 P2(3.5, -8.2, 4.2)
Rot3 attitude(const NavState &X, OptionalJacobian< 3, 9 > H)
Pose2_ Expmap(const Vector3_ &xi)
static const Similarity3 id
Some functions to compute numerical derivatives.
const Rot3 & attitude(OptionalJacobian< 3, 9 > H={}) const
Vector9 Adjoint(const Vector9 &xi_b, OptionalJacobian< 9, 9 > H_this={}, OptionalJacobian< 9, 9 > H_xib={}) const
Rot3 is a 3D rotation represented as a rotation matrix if the preprocessor symbol GTSAM_USE_QUATERNIO...
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)
static const NavState T2(Rot3::Rodrigues(0.3, 0.2, 0.1), P2, V2)
#define GTSAM_CONCEPT_LIE_INST(T)
Navigation state composing of attitude, position, and velocity.
static const Vector3 kGravity(0, 0, 9.81)
T interpolate(const T &X, const T &Y, double t, typename MakeOptionalJacobian< T, T >::type Hx={}, typename MakeOptionalJacobian< T, T >::type Hy={})
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)
Velocity3 velocity(const NavState &X, OptionalJacobian< 3, 9 > H)
TEST(NavState, Constructor)
#define CHECK_LIE_GROUP_DERIVATIVES(t1, t2)
Vector3 bodyVelocity(const Pose3 &w_t_b, const Vector3 &vec_w, OptionalJacobian< 3, 6 > Hpose={}, OptionalJacobian< 3, 3 > Hvel={})
static const Similarity3 T5(R, P, 10)
static const Vector9 kZeroXi
static const NavState state2(x2, v2)
Point3 expectedV(0.29552, 0.0446635, 1)
static const NavState kIdentity
TangentVector logmap(const Class &g) const
Rot3 inverse() const
inverse of a rotation
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
static const Point3 V2(-6.5, 3.5, 6.2)
static const Point3 P(0.2, 0.7, -2)
bool equals(const NavState &other, double tol=1e-8) const
equals
const Point3 & position(OptionalJacobian< 3, 9 > H={}) const
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)
Rot3 expectedR(c, -s, 0, s, c, 0, 0, 0, 1)
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy y set format x g set format y g set format x2 g set format y2 g set format z g set angles radians set nogrid set key title set key left top Right noreverse box linetype linewidth samplen spacing width set nolabel set noarrow set nologscale set logscale x set set pointsize set encoding default set nopolar set noparametric set set set set surface set nocontour set clabel set mapping cartesian set nohidden3d set cntrparam order set cntrparam linear set cntrparam levels auto set cntrparam points set size set origin
std::function< Vector9(const NavState &, const bool &)> coriolis
static const Point3 V(3, 0.4, -2.2)
void Print(const CONTAINER &keys, const string &s, const KeyFormatter &keyFormatter)
Array< int, Dynamic, 1 > v
Point3 position(const NavState &X, OptionalJacobian< 3, 9 > H)
Matrix4 matrix() const
Calculate 4*4 matrix group equivalent.
const Velocity3 & velocity(OptionalJacobian< 3, 9 > H={}) const
T between(const T &t1, const T &t2)
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
NavState inverse() const
inverse transformation with derivatives
NavState retract(const Vector9 &v, OptionalJacobian< 9, 9 > H1={}, OptionalJacobian< 9, 9 > H2={}) const
retract with optional derivatives
Class expmap(const TangentVector &v) const
static const Similarity3 T1(R, Point3(3.5, -8.2, 4.2), 1)
Vector9 localCoordinates(const NavState &g, OptionalJacobian< 9, 9 > H1={}, OptionalJacobian< 9, 9 > H2={}) const
localCoordinates with optional derivatives
Vector9 coriolis(double dt, const Vector3 &omega, bool secondOrder=false, OptionalJacobian< 9, 9 > H={}) const
Compute tangent space contribution due to Coriolis forces.
Provides convenient mappings of common member functions for testing.
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
#define CHECK_CHART_DERIVATIVES(t1, t2)
#define ROT3_DEFAULT_COORDINATES_MODE
static const NavState T3(Rot3::Rodrigues(-90, 0, 0), Point3(5, 6, 7), Point3(1, 2, 3))
gtsam
Author(s):
autogenerated on Thu Dec 19 2024 04:06:48