Skip to content

Commit 5e419e1

Browse files
Merge pull request #1828 from borglab/ta-bata2
Add the BATA translation averaging factor with unit tests
2 parents e3dd4e1 + 14d8870 commit 5e419e1

6 files changed

+510
-16
lines changed

gtsam/sfm/TranslationFactor.h

+76-9
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,16 +37,13 @@ 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:
4443
typedef NoiseModelFactorN<Point3, Point3> Base;
4544
Point3 measured_w_aZb_;
4645

4746
public:
48-
4947
// Provide access to the Matrix& version of evaluateError:
5048
using NoiseModelFactor2<Point3, Point3>::evaluateError;
5149

@@ -60,20 +58,20 @@ class TranslationFactor : public NoiseModelFactorN<Point3, Point3> {
6058
* @brief Caclulate error: (norm(Tb - Ta) - measurement)
6159
* where Tb and Ta are Point3 translations and measurement is
6260
* the Unit3 translation direction from a to b.
63-
*
61+
*
6462
* @param Ta translation for key a
6563
* @param Tb translation for key b
6664
* @param H1 optional jacobian in Ta
6765
* @param H2 optional jacobian in Tb
6866
* @return * Vector
6967
*/
70-
Vector evaluateError(
71-
const Point3& Ta, const Point3& Tb,
72-
OptionalMatrixType H1,
73-
OptionalMatrixType H2) const override {
68+
Vector evaluateError(const Point3& Ta, const Point3& Tb,
69+
OptionalMatrixType H1,
70+
OptionalMatrixType H2) const override {
7471
const Point3 dir = Tb - Ta;
7572
Matrix33 H_predicted_dir;
76-
const Point3 predicted = normalize(dir, H1 || H2 ? &H_predicted_dir : nullptr);
73+
const Point3 predicted =
74+
normalize(dir, H1 || H2 ? &H_predicted_dir : nullptr);
7775
if (H1) *H1 = -H_predicted_dir;
7876
if (H2) *H2 = H_predicted_dir;
7977
return predicted - measured_w_aZb_;
@@ -89,4 +87,73 @@ class TranslationFactor : public NoiseModelFactorN<Point3, Point3> {
8987
}
9088
#endif
9189
}; // \ TranslationFactor
90+
91+
/**
92+
* Binary factor for a relative translation direction measurement
93+
* w_aZb. The measurement equation is
94+
* w_aZb = scale * (Tb - Ta)
95+
* i.e., w_aZb is the translation direction from frame A to B, in world
96+
* coordinates.
97+
*
98+
* Instead of normalizing the Tb - Ta vector, we multiply it by a scalar
99+
* which is also jointly optimized. This is inspired by the work on
100+
* BATA (Zhuang et al, CVPR 2018).
101+
*
102+
* @ingroup sfm
103+
*/
104+
class BilinearAngleTranslationFactor
105+
: public NoiseModelFactorN<Point3, Point3, Vector1> {
106+
private:
107+
typedef NoiseModelFactorN<Point3, Point3, Vector1> Base;
108+
Point3 measured_w_aZb_;
109+
110+
public:
111+
/// default constructor
112+
BilinearAngleTranslationFactor() {}
113+
114+
BilinearAngleTranslationFactor(Key a, Key b, Key scale_key,
115+
const Unit3& w_aZb,
116+
const SharedNoiseModel& noiseModel)
117+
: Base(noiseModel, a, b, scale_key), measured_w_aZb_(w_aZb.point3()) {}
118+
119+
// Provide access to the Matrix& version of evaluateError:
120+
using NoiseModelFactor2<Point3, Point3, Vector1>::evaluateError;
121+
122+
/**
123+
* @brief Caclulate error: (scale * (Tb - Ta) - measurement)
124+
* where Tb and Ta are Point3 translations and measurement is
125+
* the Unit3 translation direction from a to b.
126+
*
127+
* @param Ta translation for key a
128+
* @param Tb translation for key b
129+
* @param H1 optional jacobian in Ta
130+
* @param H2 optional jacobian in Tb
131+
* @return * Vector
132+
*/
133+
Vector evaluateError(const Point3& Ta, const Point3& Tb, const Vector1& scale,
134+
OptionalMatrixType H1, OptionalMatrixType H2,
135+
OptionalMatrixType H3) const override {
136+
// Ideally we should use a positive real valued scalar datatype for scale.
137+
const double abs_scale = abs(scale[0]);
138+
const Point3 predicted = (Tb - Ta) * abs_scale;
139+
if (H1) {
140+
*H1 = -Matrix3::Identity() * abs_scale;
141+
}
142+
if (H2) {
143+
*H2 = Matrix3::Identity() * abs_scale;
144+
}
145+
if (H3) {
146+
*H3 = scale[0] >= 0 ? (Tb - Ta) : (Ta - Tb);
147+
}
148+
return predicted - measured_w_aZb_;
149+
}
150+
151+
private:
152+
friend class boost::serialization::access;
153+
template <class ARCHIVE>
154+
void serialize(ARCHIVE& ar, const unsigned int /*version*/) {
155+
ar& boost::serialization::make_nvp(
156+
"Base", boost::serialization::base_object<Base>(*this));
157+
}
158+
}; // \ BilinearAngleTranslationFactor
92159
} // 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+
uint64_t 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 (uint64_t 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 GTSAM_EXPORT 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 = false)
73+
: lmParams_(lmParams),
74+
use_bilinear_translation_factor_(use_bilinear_translation_factor) {}
7175

7276
/**
7377
* @brief Default constructor.

gtsam/sfm/sfm.i

+4-2
Original file line numberDiff line numberDiff line change
@@ -23,7 +23,7 @@ virtual class SfmTrack : gtsam::SfmTrack2d {
2323
SfmTrack();
2424
SfmTrack(const gtsam::Point3& pt);
2525
const Point3& point3() const;
26-
26+
2727
Point3 p;
2828

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

40-
#include <gtsam/nonlinear/Values.h>
4140
#include <gtsam/nonlinear/NonlinearFactorGraph.h>
41+
#include <gtsam/nonlinear/Values.h>
4242
#include <gtsam/sfm/SfmData.h>
4343
class SfmData {
4444
SfmData();
@@ -268,6 +268,8 @@ class MFAS {
268268
#include <gtsam/sfm/TranslationRecovery.h>
269269

270270
class TranslationRecovery {
271+
TranslationRecovery(const gtsam::LevenbergMarquardtParams& lmParams,
272+
const bool use_bilinear_translation_factor);
271273
TranslationRecovery(const gtsam::LevenbergMarquardtParams& lmParams);
272274
TranslationRecovery(); // default params.
273275
void addPrior(const gtsam::BinaryMeasurementsUnit3& relativeTranslations,

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)