Skip to content

Commit a3ac0f9

Browse files
add the bata translation factor, with unit tests
1 parent 4f66a49 commit a3ac0f9

6 files changed

+506
-12
lines changed

gtsam/sfm/TranslationFactor.h

+72-4
Original file line numberDiff line numberDiff line change
@@ -18,6 +18,7 @@
1818
* @brief Binary factor for a relative translation direction measurement.
1919
*/
2020

21+
#include <gtsam/base/Vector.h>
2122
#include <gtsam/geometry/Point3.h>
2223
#include <gtsam/geometry/Unit3.h>
2324
#include <gtsam/linear/NoiseModel.h>
@@ -36,8 +37,6 @@ namespace gtsam {
3637
* normalized(Tb - Ta) - w_aZb.point3()
3738
*
3839
* @ingroup sfm
39-
*
40-
*
4140
*/
4241
class TranslationFactor : public NoiseModelFactorN<Point3, Point3> {
4342
private:
@@ -56,7 +55,7 @@ class TranslationFactor : public NoiseModelFactorN<Point3, Point3> {
5655
* @brief Caclulate error: (norm(Tb - Ta) - measurement)
5756
* where Tb and Ta are Point3 translations and measurement is
5857
* the Unit3 translation direction from a to b.
59-
*
58+
*
6059
* @param Ta translation for key a
6160
* @param Tb translation for key b
6261
* @param H1 optional jacobian in Ta
@@ -69,7 +68,8 @@ class TranslationFactor : public NoiseModelFactorN<Point3, Point3> {
6968
boost::optional<Matrix&> H2 = boost::none) const override {
7069
const Point3 dir = Tb - Ta;
7170
Matrix33 H_predicted_dir;
72-
const Point3 predicted = normalize(dir, H1 || H2 ? &H_predicted_dir : nullptr);
71+
const Point3 predicted =
72+
normalize(dir, H1 || H2 ? &H_predicted_dir : nullptr);
7373
if (H1) *H1 = -H_predicted_dir;
7474
if (H2) *H2 = H_predicted_dir;
7575
return predicted - measured_w_aZb_;
@@ -83,4 +83,72 @@ class TranslationFactor : public NoiseModelFactorN<Point3, Point3> {
8383
"Base", boost::serialization::base_object<Base>(*this));
8484
}
8585
}; // \ TranslationFactor
86+
87+
/**
88+
* Binary factor for a relative translation direction measurement
89+
* w_aZb. The measurement equation is
90+
* w_aZb = scale * (Tb - Ta)
91+
* i.e., w_aZb is the translation direction from frame A to B, in world
92+
* coordinates.
93+
*
94+
* Instead of normalizing the Tb - Ta vector, we multiply it by a scalar
95+
* which is also jointly optimized. This is inspired by the work on
96+
* BATA (Zhuang et al, CVPR 2018).
97+
*
98+
* @ingroup sfm
99+
*/
100+
class BilinearAngleTranslationFactor
101+
: public NoiseModelFactorN<Point3, Point3, Vector1> {
102+
private:
103+
typedef NoiseModelFactorN<Point3, Point3, Vector1> Base;
104+
Point3 measured_w_aZb_;
105+
106+
public:
107+
/// default constructor
108+
BilinearAngleTranslationFactor() {}
109+
110+
BilinearAngleTranslationFactor(Key a, Key b, Key scale_key, const Unit3& w_aZb,
111+
const SharedNoiseModel& noiseModel)
112+
: Base(noiseModel, a, b, scale_key),
113+
measured_w_aZb_(w_aZb.point3()) {}
114+
115+
/**
116+
* @brief Caclulate error: (scale * (Tb - Ta) - measurement)
117+
* where Tb and Ta are Point3 translations and measurement is
118+
* the Unit3 translation direction from a to b.
119+
*
120+
* @param Ta translation for key a
121+
* @param Tb translation for key b
122+
* @param H1 optional jacobian in Ta
123+
* @param H2 optional jacobian in Tb
124+
* @return * Vector
125+
*/
126+
Vector evaluateError(
127+
const Point3& Ta, const Point3& Tb, const Vector1& scale,
128+
boost::optional<Matrix&> H1 = boost::none,
129+
boost::optional<Matrix&> H2 = boost::none,
130+
boost::optional<Matrix&> H3 = boost::none) const override {
131+
// Ideally we should use a positive real valued scalar datatype for scale.
132+
const double abs_scale = abs(scale[0]);
133+
const Point3 predicted = (Tb - Ta) * abs_scale;
134+
if (H1) {
135+
*H1 = -Matrix3::Identity() * abs_scale;
136+
}
137+
if (H2) {
138+
*H2 = Matrix3::Identity() * abs_scale;
139+
}
140+
if (H3) {
141+
*H3 = scale[0] >= 0 ? (Tb - Ta) : (Ta - Tb);
142+
}
143+
return predicted - measured_w_aZb_;
144+
}
145+
146+
private:
147+
friend class boost::serialization::access;
148+
template <class ARCHIVE>
149+
void serialize(ARCHIVE& ar, const unsigned int /*version*/) {
150+
ar& boost::serialization::make_nvp(
151+
"Base", boost::serialization::base_object<Base>(*this));
152+
}
153+
}; // \ BilinearAngleTranslationFactor
86154
} // namespace gtsam

gtsam/sfm/TranslationRecovery.cpp

+16-2
Original file line numberDiff line numberDiff line change
@@ -101,9 +101,17 @@ NonlinearFactorGraph TranslationRecovery::buildGraph(
101101
NonlinearFactorGraph graph;
102102

103103
// Add translation factors for input translation directions.
104+
uint i = 0;
104105
for (auto edge : relativeTranslations) {
105-
graph.emplace_shared<TranslationFactor>(edge.key1(), edge.key2(),
106-
edge.measured(), edge.noiseModel());
106+
if (use_bilinear_translation_factor_) {
107+
graph.emplace_shared<BilinearAngleTranslationFactor>(
108+
edge.key1(), edge.key2(), Symbol('S', i), edge.measured(),
109+
edge.noiseModel());
110+
} else {
111+
graph.emplace_shared<TranslationFactor>(
112+
edge.key1(), edge.key2(), edge.measured(), edge.noiseModel());
113+
}
114+
i++;
107115
}
108116
return graph;
109117
}
@@ -163,6 +171,12 @@ Values TranslationRecovery::initializeRandomly(
163171
insert(edge.key1());
164172
insert(edge.key2());
165173
}
174+
175+
if (use_bilinear_translation_factor_) {
176+
for (uint i = 0; i < relativeTranslations.size(); i++) {
177+
initial.insert<Vector1>(Symbol('S', i), Vector1(1.0));
178+
}
179+
}
166180
return initial;
167181
}
168182

gtsam/sfm/TranslationRecovery.h

+6-2
Original file line numberDiff line numberDiff line change
@@ -60,14 +60,18 @@ class TranslationRecovery {
6060
// Parameters.
6161
LevenbergMarquardtParams lmParams_;
6262

63+
const bool use_bilinear_translation_factor_ = false;
64+
6365
public:
6466
/**
6567
* @brief Construct a new Translation Recovery object
6668
*
6769
* @param lmParams parameters for optimization.
6870
*/
69-
TranslationRecovery(const LevenbergMarquardtParams &lmParams)
70-
: lmParams_(lmParams) {}
71+
TranslationRecovery(const LevenbergMarquardtParams &lmParams,
72+
bool use_bilinear_translation_factor)
73+
: lmParams_(lmParams),
74+
use_bilinear_translation_factor_(use_bilinear_translation_factor) {}
7175

7276
/**
7377
* @brief Default constructor.

gtsam/sfm/sfm.i

+4-3
Original file line numberDiff line numberDiff line change
@@ -24,7 +24,7 @@ virtual class SfmTrack : gtsam::SfmTrack2d {
2424
SfmTrack();
2525
SfmTrack(const gtsam::Point3& pt);
2626
const Point3& point3() const;
27-
27+
2828
Point3 p;
2929

3030
double r;
@@ -38,8 +38,8 @@ virtual class SfmTrack : gtsam::SfmTrack2d {
3838
bool equals(const gtsam::SfmTrack& expected, double tol) const;
3939
};
4040

41-
#include <gtsam/nonlinear/Values.h>
4241
#include <gtsam/nonlinear/NonlinearFactorGraph.h>
42+
#include <gtsam/nonlinear/Values.h>
4343
#include <gtsam/sfm/SfmData.h>
4444
class SfmData {
4545
SfmData();
@@ -288,7 +288,8 @@ class MFAS {
288288
#include <gtsam/sfm/TranslationRecovery.h>
289289

290290
class TranslationRecovery {
291-
TranslationRecovery(const gtsam::LevenbergMarquardtParams& lmParams);
291+
TranslationRecovery(const gtsam::LevenbergMarquardtParams& lmParams,
292+
const bool use_bilinear_translation_factor);
292293
TranslationRecovery(); // default params.
293294
void addPrior(const gtsam::BinaryMeasurementsUnit3& relativeTranslations,
294295
const double scale,

gtsam/sfm/tests/testTranslationFactor.cpp

+74-1
Original file line numberDiff line numberDiff line change
@@ -30,7 +30,7 @@ using namespace gtsam;
3030
static SharedNoiseModel model(noiseModel::Isotropic::Sigma(3, 0.05));
3131

3232
// Keys are deliberately *not* in sorted order to test that case.
33-
static const Key kKey1(2), kKey2(1);
33+
static const Key kKey1(2), kKey2(1), kEdgeKey(3);
3434
static const Unit3 kMeasured(1, 0, 0);
3535

3636
/* ************************************************************************* */
@@ -99,6 +99,79 @@ TEST(TranslationFactor, Jacobian) {
9999
EXPECT(assert_equal(H2Expected, H2Actual, 1e-9));
100100
}
101101

102+
103+
/* ************************************************************************* */
104+
TEST(BilinearAngleTranslationFactor, Constructor) {
105+
BilinearAngleTranslationFactor factor(kKey1, kKey2, kEdgeKey, kMeasured, model);
106+
}
107+
108+
/* ************************************************************************* */
109+
TEST(BilinearAngleTranslationFactor, ZeroError) {
110+
// Create a factor
111+
BilinearAngleTranslationFactor factor(kKey1, kKey2, kEdgeKey, kMeasured, model);
112+
113+
// Set the linearization
114+
Point3 T1(1, 0, 0), T2(2, 0, 0);
115+
Vector1 scale(1.0);
116+
117+
// Use the factor to calculate the error
118+
Vector actualError(factor.evaluateError(T1, T2, scale));
119+
120+
// Verify we get the expected error
121+
Vector expected = (Vector3() << 0, 0, 0).finished();
122+
EXPECT(assert_equal(expected, actualError, 1e-9));
123+
}
124+
125+
/* ************************************************************************* */
126+
TEST(BilinearAngleTranslationFactor, NonZeroError) {
127+
// create a factor
128+
BilinearAngleTranslationFactor factor(kKey1, kKey2, kEdgeKey, kMeasured, model);
129+
130+
// set the linearization
131+
Point3 T1(0, 1, 1), T2(0, 2, 2);
132+
Vector1 scale(1.0 / sqrt(2));
133+
134+
// use the factor to calculate the error
135+
Vector actualError(factor.evaluateError(T1, T2, scale));
136+
137+
// verify we get the expected error
138+
Vector expected = (Vector3() << -1, 1 / sqrt(2), 1 / sqrt(2)).finished();
139+
EXPECT(assert_equal(expected, actualError, 1e-9));
140+
}
141+
142+
/* ************************************************************************* */
143+
Vector bilinearAngleFactorError(const Point3 &T1, const Point3 &T2, const Vector1 &scale,
144+
const BilinearAngleTranslationFactor &factor) {
145+
return factor.evaluateError(T1, T2, scale);
146+
}
147+
148+
TEST(BilinearAngleTranslationFactor, Jacobian) {
149+
// Create a factor
150+
BilinearAngleTranslationFactor factor(kKey1, kKey2, kEdgeKey, kMeasured, model);
151+
152+
// Set the linearization
153+
Point3 T1(1, 0, 0), T2(2, 0, 0);
154+
Vector1 scale(1.0);
155+
156+
// Use the factor to calculate the Jacobians
157+
Matrix H1Actual, H2Actual, H3Actual;
158+
factor.evaluateError(T1, T2, scale, H1Actual, H2Actual, H3Actual);
159+
160+
// Use numerical derivatives to calculate the Jacobians
161+
Matrix H1Expected, H2Expected, H3Expected;
162+
H1Expected = numericalDerivative11<Vector, Point3>(
163+
std::bind(&bilinearAngleFactorError, std::placeholders::_1, T2, scale, factor), T1);
164+
H2Expected = numericalDerivative11<Vector, Point3>(
165+
std::bind(&bilinearAngleFactorError, T1, std::placeholders::_1, scale, factor), T2);
166+
H3Expected = numericalDerivative11<Vector, Vector1>(
167+
std::bind(&bilinearAngleFactorError, T1, T2, std::placeholders::_1, factor), scale);
168+
169+
// Verify the Jacobians are correct
170+
EXPECT(assert_equal(H1Expected, H1Actual, 1e-9));
171+
EXPECT(assert_equal(H2Expected, H2Actual, 1e-9));
172+
EXPECT(assert_equal(H3Expected, H3Actual, 1e-9));
173+
}
174+
102175
/* ************************************************************************* */
103176
int main() {
104177
TestResult tr;

0 commit comments

Comments
 (0)