GTSAM  4.0.2
C++ library for smoothing and mapping (SAM)
MagFactor.h
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 #pragma once
20 
22 #include <gtsam/geometry/Rot2.h>
23 #include <gtsam/geometry/Rot3.h>
24 
25 namespace gtsam {
26 
33 class MagFactor: public NoiseModelFactorN<Rot2> {
34 
35  const Point3 measured_;
36  const Point3 nM_;
37  const Point3 bias_;
38 
39 public:
40 
41  // Provide access to Matrix& version of evaluateError:
42  using NoiseModelFactor1<Rot2>::evaluateError;
43 
53  MagFactor(Key key, const Point3& measured, double scale,
54  const Unit3& direction, const Point3& bias,
55  const SharedNoiseModel& model) :
56  NoiseModelFactorN<Rot2>(model, key), //
57  measured_(measured), nM_(scale * direction), bias_(bias) {
58  }
59 
61  NonlinearFactor::shared_ptr clone() const override {
62  return std::static_pointer_cast<NonlinearFactor>(
63  NonlinearFactor::shared_ptr(new MagFactor(*this)));
64  }
65 
66  static Point3 unrotate(const Rot2& R, const Point3& p,
68  Point3 q = Rot3::Yaw(R.theta()).unrotate(p, HR, {});
69  if (HR) {
70  // assign to temporary first to avoid error in Win-Debug mode
71  Matrix H = HR->col(2);
72  *HR = H;
73  }
74  return q;
75  }
76 
80  Vector evaluateError(const Rot2& nRb, OptionalMatrixType H) const override {
81  // measured bM = nRb� * nM + b
82  Point3 hx = unrotate(nRb, nM_, H) + bias_;
83  return (hx - measured_);
84  }
85 };
86 
92 class MagFactor1: public NoiseModelFactorN<Rot3> {
93 
94  const Point3 measured_;
95  const Point3 nM_;
96  const Point3 bias_;
97 
98 public:
99 
100  // Provide access to Matrix& version of evaluateError:
101  using NoiseModelFactor1<Rot3>::evaluateError;
102 
103 
105  MagFactor1(Key key, const Point3& measured, double scale,
106  const Unit3& direction, const Point3& bias,
107  const SharedNoiseModel& model) :
108  NoiseModelFactorN<Rot3>(model, key), //
109  measured_(measured), nM_(scale * direction), bias_(bias) {
110  }
111 
113  NonlinearFactor::shared_ptr clone() const override {
114  return std::static_pointer_cast<NonlinearFactor>(
115  NonlinearFactor::shared_ptr(new MagFactor1(*this)));
116  }
117 
121  Vector evaluateError(const Rot3& nRb, OptionalMatrixType H) const override {
122  // measured bM = nRb� * nM + b
123  Point3 hx = nRb.unrotate(nM_, H, OptionalNone) + bias_;
124  return (hx - measured_);
125  }
126 };
127 
133 class MagFactor2: public NoiseModelFactorN<Point3, Point3> {
134 
135  const Point3 measured_;
136  const Rot3 bRn_;
137 
138 public:
139 
140  // Provide access to Matrix& version of evaluateError:
141  using NoiseModelFactor2<Point3, Point3>::evaluateError;
142 
143 
145  MagFactor2(Key key1, Key key2, const Point3& measured, const Rot3& nRb,
146  const SharedNoiseModel& model) :
147  NoiseModelFactorN<Point3, Point3>(model, key1, key2), //
148  measured_(measured), bRn_(nRb.inverse()) {
149  }
150 
152  NonlinearFactor::shared_ptr clone() const override {
153  return std::static_pointer_cast<NonlinearFactor>(
154  NonlinearFactor::shared_ptr(new MagFactor2(*this)));
155  }
156 
162  Vector evaluateError(const Point3& nM, const Point3& bias,
163  OptionalMatrixType H1, OptionalMatrixType H2) const override {
164  // measured bM = nRb� * nM + b, where b is unknown bias
165  Point3 hx = bRn_.rotate(nM, OptionalNone, H1) + bias;
166  if (H2)
167  *H2 = I_3x3;
168  return (hx - measured_);
169  }
170 };
171 
177 class MagFactor3: public NoiseModelFactorN<double, Unit3, Point3> {
178 
179  const Point3 measured_;
180  const Rot3 bRn_;
181 
182 public:
183 
184  // Provide access to Matrix& version of evaluateError:
185  using NoiseModelFactor3<double, Unit3, Point3>::evaluateError;
186 
187 
189  MagFactor3(Key key1, Key key2, Key key3, const Point3& measured,
190  const Rot3& nRb, const SharedNoiseModel& model) :
191  NoiseModelFactorN<double, Unit3, Point3>(model, key1, key2, key3), //
192  measured_(measured), bRn_(nRb.inverse()) {
193  }
194 
196  NonlinearFactor::shared_ptr clone() const override {
197  return std::static_pointer_cast<NonlinearFactor>(
198  NonlinearFactor::shared_ptr(new MagFactor3(*this)));
199  }
200 
206  Vector evaluateError(const double& scale, const Unit3& direction,
207  const Point3& bias, OptionalMatrixType H1,
208  OptionalMatrixType H2, OptionalMatrixType H3) const override {
209  // measured bM = nRb� * nM + b, where b is unknown bias
210  Unit3 rotated = bRn_.rotate(direction, OptionalNone, H2);
211  Point3 hx = scale * rotated.point3() + bias;
212  if (H1)
213  *H1 = rotated.point3();
214  if (H2) // H2 is 2*2, but we need 3*2
215  {
216  Matrix H;
217  rotated.point3(H);
218  *H2 = scale * H * (*H2);
219  }
220  if (H3)
221  *H3 = I_3x3;
222  return (hx - measured_);
223  }
224 };
225 
226 }
227 
Definition: NonlinearFactor.h:431
Definition: Rot2.h:35
Vector evaluateError(const double &scale, const Unit3 &direction, const Point3 &bias, OptionalMatrixType H1, OptionalMatrixType H2, OptionalMatrixType H3) const override
vector of errors
Definition: MagFactor.h:206
MagFactor2(Key key1, Key key2, const Point3 &measured, const Rot3 &nRb, const SharedNoiseModel &model)
Definition: MagFactor.h:145
Vector evaluateError(const Point3 &nM, const Point3 &bias, OptionalMatrixType H1, OptionalMatrixType H2) const override
vector of errors
Definition: MagFactor.h:162
Definition: MagFactor.h:92
Key key() const
Definition: NonlinearFactor.h:582
Definition: NonlinearFactor.h:68
Rot3 is a 3D rotation represented as a rotation matrix if the preprocessor symbol GTSAM_USE_QUATERNIO...
Definition: Rot3.h:58
#define OptionalNone
Definition: NonlinearFactor.h:49
Definition: MagFactor.h:133
Represents a 3D point on a unit sphere.
Definition: Unit3.h:42
Matrix * OptionalMatrixType
Definition: NonlinearFactor.h:55
static Rot3 Yaw(double t)
Positive yaw is to right (as in aircraft heading). See ypr.
Definition: Rot3.h:174
Point3 unrotate(const Point3 &p, OptionalJacobian< 3, 3 > H1={}, OptionalJacobian< 3, 3 > H2={}) const
rotate point from world to rotated frame
Vector evaluateError(const Rot2 &nRb, OptionalMatrixType H) const override
vector of errors
Definition: MagFactor.h:80
Vector evaluateError(const Rot3 &nRb, OptionalMatrixType H) const override
vector of errors
Definition: MagFactor.h:121
MagFactor1(Key key, const Point3 &measured, double scale, const Unit3 &direction, const Point3 &bias, const SharedNoiseModel &model)
Definition: MagFactor.h:105
NonlinearFactor::shared_ptr clone() const override
Definition: MagFactor.h:113
Definition: chartTesting.h:28
NonlinearFactor::shared_ptr clone() const override
Definition: MagFactor.h:196
Non-linear factor base classes.
Definition: MagFactor.h:33
NonlinearFactor::shared_ptr clone() const override
Definition: MagFactor.h:61
Point3 point3(OptionalJacobian< 3, 2 > H={}) const
Return unit-norm Point3.
Vector3 Point3
Definition: Point3.h:38
Point3 rotate(const Point3 &p, OptionalJacobian< 3, 3 > H1={}, OptionalJacobian< 3, 3 > H2={}) const
Definition: MagFactor.h:177
std::uint64_t Key
Integer nonlinear key type.
Definition: types.h:102
2D rotation
NonlinearFactor::shared_ptr clone() const override
Definition: MagFactor.h:152
MagFactor3(Key key1, Key key2, Key key3, const Point3 &measured, const Rot3 &nRb, const SharedNoiseModel &model)
Definition: MagFactor.h:189
noiseModel::Base::shared_ptr SharedNoiseModel
Definition: NoiseModel.h:741
double theta() const
Definition: Rot2.h:186
3D rotation represented as a rotation matrix or quaternion
MagFactor(Key key, const Point3 &measured, double scale, const Unit3 &direction, const Point3 &bias, const SharedNoiseModel &model)
Definition: MagFactor.h:53