testNoiseModel.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 
22 
24 
25 #include <boost/assign/std/vector.hpp>
26 
27 #include <iostream>
28 #include <limits>
29 
30 using namespace std;
31 using namespace gtsam;
32 using namespace noiseModel;
33 using namespace boost::assign;
34 
35 static const double kSigma = 2, kInverseSigma = 1.0 / kSigma,
37 static const Matrix R = Matrix3::Identity() * kInverseSigma;
38 static const Matrix kCovariance = Matrix3::Identity() * kVariance;
39 static const Vector3 kSigmas(kSigma, kSigma, kSigma);
40 
41 /* ************************************************************************* */
42 TEST(NoiseModel, constructors)
43 {
44  Vector whitened = Vector3(5.0,10.0,15.0);
45  Vector unwhitened = Vector3(10.0,20.0,30.0);
46 
47  // Construct noise models
48  vector<Gaussian::shared_ptr> m;
49  m.push_back(Gaussian::SqrtInformation(R,false));
50  m.push_back(Gaussian::Covariance(kCovariance,false));
51  m.push_back(Gaussian::Information(kCovariance.inverse(),false));
52  m.push_back(Diagonal::Sigmas(kSigmas,false));
53  m.push_back(Diagonal::Variances((Vector3(kVariance, kVariance, kVariance)),false));
54  m.push_back(Diagonal::Precisions(Vector3(prc, prc, prc),false));
55  m.push_back(Isotropic::Sigma(3, kSigma,false));
56  m.push_back(Isotropic::Variance(3, kVariance,false));
57  m.push_back(Isotropic::Precision(3, prc,false));
58 
59  // test kSigmas
60  for(Gaussian::shared_ptr mi: m)
61  EXPECT(assert_equal(kSigmas,mi->sigmas()));
62 
63  // test whiten
64  for(Gaussian::shared_ptr mi: m)
65  EXPECT(assert_equal(whitened,mi->whiten(unwhitened)));
66 
67  // test unwhiten
68  for(Gaussian::shared_ptr mi: m)
69  EXPECT(assert_equal(unwhitened,mi->unwhiten(whitened)));
70 
71  // test squared Mahalanobis distance
72  double distance = 5*5+10*10+15*15;
73  for(Gaussian::shared_ptr mi: m)
74  DOUBLES_EQUAL(distance,mi->squaredMahalanobisDistance(unwhitened),1e-9);
75 
76  // test R matrix
77  for(Gaussian::shared_ptr mi: m)
78  EXPECT(assert_equal(R,mi->R()));
79 
80  // test covariance
81  for(Gaussian::shared_ptr mi: m)
82  EXPECT(assert_equal(kCovariance,mi->covariance()));
83 
84  // test covariance
85  for(Gaussian::shared_ptr mi: m)
86  EXPECT(assert_equal(kCovariance.inverse(),mi->information()));
87 
88  // test Whiten operator
89  Matrix H((Matrix(3, 4) <<
90  0.0, 0.0, 1.0, 1.0,
91  0.0, 1.0, 0.0, 1.0,
92  1.0, 0.0, 0.0, 1.0).finished());
94  for(Gaussian::shared_ptr mi: m)
95  EXPECT(assert_equal(expected,mi->Whiten(H)));
96 
97  // can only test inplace version once :-)
98  m[0]->WhitenInPlace(H);
99  EXPECT(assert_equal(expected,H));
100 }
101 
102 /* ************************************************************************* */
103 TEST(NoiseModel, Unit)
104 {
105  Vector v = Vector3(5.0,10.0,15.0);
106  Gaussian::shared_ptr u(Unit::Create(3));
107  EXPECT(assert_equal(v,u->whiten(v)));
108 }
109 
110 /* ************************************************************************* */
111 TEST(NoiseModel, equals)
112 {
113  Gaussian::shared_ptr g1 = Gaussian::SqrtInformation(R),
114  g2 = Gaussian::SqrtInformation(I_3x3);
115  Diagonal::shared_ptr d1 = Diagonal::Sigmas(Vector3(kSigma, kSigma, kSigma)),
116  d2 = Diagonal::Sigmas(Vector3(0.1, 0.2, 0.3));
117  Isotropic::shared_ptr i1 = Isotropic::Sigma(3, kSigma),
118  i2 = Isotropic::Sigma(3, 0.7);
119 
120  EXPECT(assert_equal(*g1,*g1));
121  EXPECT(assert_inequal(*g1, *g2));
122 
123  EXPECT(assert_equal(*d1,*d1));
124  EXPECT(assert_inequal(*d1,*d2));
125 
126  EXPECT(assert_equal(*i1,*i1));
127  EXPECT(assert_inequal(*i1,*i2));
128 }
129 
130 // TODO enable test once a mechanism for smart constraints exists
132 //TEST(NoiseModel, ConstrainedSmart )
133 //{
134 // Gaussian::shared_ptr nonconstrained = Constrained::MixedSigmas((Vector3(sigma, 0.0, sigma), true);
135 // Diagonal::shared_ptr n1 = boost::dynamic_pointer_cast<Diagonal>(nonconstrained);
136 // Constrained::shared_ptr n2 = boost::dynamic_pointer_cast<Constrained>(nonconstrained);
137 // EXPECT(n1);
138 // EXPECT(!n2);
139 //
140 // Gaussian::shared_ptr constrained = Constrained::MixedSigmas(zero(3), true);
141 // Diagonal::shared_ptr c1 = boost::dynamic_pointer_cast<Diagonal>(constrained);
142 // Constrained::shared_ptr c2 = boost::dynamic_pointer_cast<Constrained>(constrained);
143 // EXPECT(c1);
144 // EXPECT(c2);
145 //}
146 
147 /* ************************************************************************* */
148 TEST(NoiseModel, ConstrainedConstructors )
149 {
151  size_t d = 3;
152  double m = 100.0;
153  Vector3 sigmas(kSigma, 0.0, 0.0);
154  Vector3 mu(200.0, 300.0, 400.0);
155  actual = Constrained::All(d);
156  // TODO: why should this be a thousand ??? Dummy variable?
157  EXPECT(assert_equal(Vector::Constant(d, 1000.0), actual->mu()));
158  EXPECT(assert_equal(Vector::Constant(d, 0), actual->sigmas()));
159  EXPECT(assert_equal(Vector::Constant(d, 0), actual->invsigmas())); // Actually zero as dummy value
160  EXPECT(assert_equal(Vector::Constant(d, 0), actual->precisions())); // Actually zero as dummy value
161 
162  actual = Constrained::All(d, m);
163  EXPECT(assert_equal(Vector::Constant(d, m), actual->mu()));
164 
165  actual = Constrained::All(d, mu);
166  EXPECT(assert_equal(mu, actual->mu()));
167 
168  actual = Constrained::MixedSigmas(mu, sigmas);
169  EXPECT(assert_equal(mu, actual->mu()));
170 
171  actual = Constrained::MixedSigmas(m, sigmas);
172  EXPECT(assert_equal(Vector::Constant(d, m), actual->mu()));
173 }
174 
175 /* ************************************************************************* */
176 TEST(NoiseModel, ConstrainedMixed )
177 {
178  Vector feasible = Vector3(1.0, 0.0, 1.0),
179  infeasible = Vector3(1.0, 1.0, 1.0);
180  Diagonal::shared_ptr d = Constrained::MixedSigmas(Vector3(kSigma, 0.0, kSigma));
181  // NOTE: we catch constrained variables elsewhere, so whitening does nothing
182  EXPECT(assert_equal(Vector3(0.5, 1.0, 0.5),d->whiten(infeasible)));
183  EXPECT(assert_equal(Vector3(0.5, 0.0, 0.5),d->whiten(feasible)));
184 
185  DOUBLES_EQUAL(0.5 * (1000.0 + 0.25 + 0.25),d->loss(d->squaredMahalanobisDistance(infeasible)),1e-9);
186  DOUBLES_EQUAL(0.5, d->squaredMahalanobisDistance(feasible),1e-9);
187  DOUBLES_EQUAL(0.5 * 0.5, d->loss(0.5),1e-9);
188 }
189 
190 /* ************************************************************************* */
191 TEST(NoiseModel, ConstrainedAll )
192 {
193  Vector feasible = Vector3(0.0, 0.0, 0.0),
194  infeasible = Vector3(1.0, 1.0, 1.0);
195 
196  Constrained::shared_ptr i = Constrained::All(3);
197  // NOTE: we catch constrained variables elsewhere, so whitening does nothing
198  EXPECT(assert_equal(Vector3(1.0, 1.0, 1.0),i->whiten(infeasible)));
199  EXPECT(assert_equal(Vector3(0.0, 0.0, 0.0),i->whiten(feasible)));
200 
201  DOUBLES_EQUAL(0.5 * 1000.0 * 3.0,i->loss(i->squaredMahalanobisDistance(infeasible)),1e-9);
202  DOUBLES_EQUAL(0.0, i->squaredMahalanobisDistance(feasible), 1e-9);
203  DOUBLES_EQUAL(0.0, i->loss(0.0),1e-9);
204 }
205 
206 /* ************************************************************************* */
207 namespace exampleQR {
208  // create a matrix to eliminate
209  Matrix Ab = (Matrix(4, 7) <<
210  -1., 0., 1., 0., 0., 0., -0.2,
211  0., -1., 0., 1., 0., 0., 0.3,
212  1., 0., 0., 0., -1., 0., 0.2,
213  0., 1., 0., 0., 0., -1., -0.1).finished();
214  Vector sigmas = (Vector(4) << 0.2, 0.2, 0.1, 0.1).finished();
215 
216  // the matrix AB yields the following factorized version:
217  Matrix Rd = (Matrix(4, 7) <<
218  11.1803, 0.0, -2.23607, 0.0, -8.94427, 0.0, 2.23607,
219  0.0, 11.1803, 0.0, -2.23607, 0.0, -8.94427,-1.56525,
220  0.0, 0.0, 4.47214, 0.0, -4.47214, 0.0, 0.0,
221  0.0, 0.0, 0.0, 4.47214, 0.0, -4.47214, 0.894427).finished();
222 
223  SharedDiagonal diagonal = noiseModel::Diagonal::Sigmas(sigmas);
224 }
225 
226 /* ************************************************************************* */
227 TEST( NoiseModel, QR )
228 {
229  Matrix Ab1 = exampleQR::Ab;
230  Matrix Ab2 = exampleQR::Ab; // otherwise overwritten !
231 
232  // Call Gaussian version
233  SharedDiagonal actual1 = exampleQR::diagonal->QR(Ab1);
234  EXPECT(actual1->isUnit());
235  EXPECT(linear_dependent(exampleQR::Rd,Ab1,1e-4)); // Ab was modified in place !!!
236 
237  // Expected result for constrained version
238  Vector expectedSigmas = (Vector(4) << 0.0894427, 0.0894427, 0.223607, 0.223607).finished();
239  SharedDiagonal expectedModel = noiseModel::Diagonal::Sigmas(expectedSigmas);
240  Matrix expectedRd2 = (Matrix(4, 7) <<
241  1., 0., -0.2, 0., -0.8, 0., 0.2,
242  0., 1., 0.,-0.2, 0., -0.8,-0.14,
243  0., 0., 1., 0., -1., 0., 0.0,
244  0., 0., 0., 1., 0., -1., 0.2).finished();
245 
246  // Call Constrained version
247  SharedDiagonal constrained = noiseModel::Constrained::MixedSigmas(exampleQR::sigmas);
248  SharedDiagonal actual2 = constrained->QR(Ab2);
249  EXPECT(assert_equal(*expectedModel, *actual2, 1e-6));
250  EXPECT(linear_dependent(expectedRd2, Ab2, 1e-6)); // Ab was modified in place !!!
251 }
252 
253 /* ************************************************************************* */
254 TEST(NoiseModel, OverdeterminedQR) {
255  Matrix Ab1(9, 4);
256  Ab1 << 0, 1, 0, 0, //
257  0, 0, 1, 0, //
258  Matrix74::Ones();
259  Matrix Ab2 = Ab1; // otherwise overwritten !
260 
261  // Call Gaussian version
262  Vector9 sigmas = Vector9::Ones() ;
263  SharedDiagonal diagonal = noiseModel::Diagonal::Sigmas(sigmas);
264  SharedDiagonal actual1 = diagonal->QR(Ab1);
265  EXPECT(actual1->isUnit());
266  Matrix expectedRd(9,4);
267  expectedRd << -2.64575131, -2.64575131, -2.64575131, -2.64575131, //
268  0.0, -1, 0, 0, //
269  0.0, 0.0, -1, 0, //
270  Matrix64::Zero();
271  EXPECT(assert_equal(expectedRd, Ab1, 1e-4)); // Ab was modified in place !!!
272 
273  // Expected result for constrained version
274  Vector3 expectedSigmas(0.377964473, 1, 1);
275  SharedDiagonal expectedModel = noiseModel::Diagonal::Sigmas(expectedSigmas);
276 
277  // Call Constrained version
278  SharedDiagonal constrained = noiseModel::Constrained::MixedSigmas(sigmas);
279  SharedDiagonal actual2 = constrained->QR(Ab2);
280  EXPECT(assert_equal(*expectedModel, *actual2, 1e-6));
281  expectedRd.row(0) *= 0.377964473; // not divided by sigma!
282  EXPECT(assert_equal(-expectedRd, Ab2, 1e-6)); // Ab was modified in place !!!
283 }
284 
285 /* ************************************************************************* */
286 TEST( NoiseModel, MixedQR )
287 {
288  // Call Constrained version, with first and third row treated as constraints
289  // Naming the 6 variables u,v,w,x,y,z, we have
290  // u = -z
291  // w = -x
292  // And let's have simple priors on variables
293  Matrix Ab(5,6+1);
294  Ab <<
295  1,0,0,0,0,1, 0, // u+z = 0
296  0,0,0,0,1,0, 0, // y^2
297  0,0,1,1,0,0, 0, // w+x = 0
298  0,1,0,0,0,0, 0, // v^2
299  0,0,0,0,0,1, 0; // z^2
300  Vector mixed_sigmas = (Vector(5) << 0, 1, 0, 1, 1).finished();
301  SharedDiagonal constrained = noiseModel::Constrained::MixedSigmas(mixed_sigmas);
302 
303  // Expected result
304  Vector expectedSigmas = (Vector(5) << 0, 1, 0, 1, 1).finished();
305  SharedDiagonal expectedModel = noiseModel::Diagonal::Sigmas(expectedSigmas);
306  Matrix expectedRd(5, 6+1);
307  expectedRd << 1, 0, 0, 0, 0, 1, 0, //
308  0, 1, 0, 0, 0, 0, 0, //
309  0, 0, 1, 1, 0, 0, 0, //
310  0, 0, 0, 0, 1, 0, 0, //
311  0, 0, 0, 0, 0, 1, 0; //
312 
313  SharedDiagonal actual = constrained->QR(Ab);
314  EXPECT(assert_equal(*expectedModel,*actual,1e-6));
315  EXPECT(linear_dependent(expectedRd,Ab,1e-6)); // Ab was modified in place !!!
316 }
317 
318 /* ************************************************************************* */
319 TEST( NoiseModel, MixedQR2 )
320 {
321  // Let's have three variables x,y,z, but x=z and y=z
322  // Hence, all non-constraints are really measurements on z
323  Matrix Ab(11,3+1);
324  Ab <<
325  1,0,0, 0, //
326  0,1,0, 0, //
327  0,0,1, 0, //
328  -1,0,1, 0, // x=z
329  1,0,0, 0, //
330  0,1,0, 0, //
331  0,0,1, 0, //
332  0,-1,1, 0, // y=z
333  1,0,0, 0, //
334  0,1,0, 0, //
335  0,0,1, 0; //
336 
337  Vector sigmas(11);
338  sigmas.setOnes();
339  sigmas[3] = 0;
340  sigmas[7] = 0;
341  SharedDiagonal constrained = noiseModel::Constrained::MixedSigmas(sigmas);
342 
343  // Expected result
344  Vector3 expectedSigmas(0,0,1.0/3);
345  SharedDiagonal expectedModel = noiseModel::Constrained::MixedSigmas(expectedSigmas);
346  Matrix expectedRd(11, 3+1);
347  expectedRd.setZero();
348  expectedRd.row(0) << -1, 0, 1, 0; // x=z
349  expectedRd.row(1) << 0, -1, 1, 0; // y=z
350  expectedRd.row(2) << 0, 0, 1, 0; // z=0 +/- 1/3
351 
352  SharedDiagonal actual = constrained->QR(Ab);
353  EXPECT(assert_equal(*expectedModel,*actual,1e-6));
354  EXPECT(assert_equal(expectedRd,Ab,1e-6)); // Ab was modified in place !!!
355 }
356 
357 /* ************************************************************************* */
358 TEST( NoiseModel, FullyConstrained )
359 {
360  Matrix Ab(3,7);
361  Ab <<
362  1,0,0,0,0,1, 2, // u+z = 2
363  0,0,1,1,0,0, 4, // w+x = 4
364  0,1,0,1,1,1, 8; // v+x+y+z=8
365  SharedDiagonal constrained = noiseModel::Constrained::All(3);
366 
367  // Expected result
368  SharedDiagonal expectedModel = noiseModel::Diagonal::Sigmas(Vector3 (0,0,0));
369  Matrix expectedRd(3, 7);
370  expectedRd << 1, 0, 0, 0, 0, 1, 2, //
371  0, 1, 0, 1, 1, 1, 8, //
372  0, 0, 1, 1, 0, 0, 4; //
373 
374  SharedDiagonal actual = constrained->QR(Ab);
375  EXPECT(assert_equal(*expectedModel,*actual,1e-6));
376  EXPECT(linear_dependent(expectedRd,Ab,1e-6)); // Ab was modified in place !!!
377 }
378 
379 /* ************************************************************************* */
380 // This matches constraint_eliminate2 in testJacobianFactor
381 TEST(NoiseModel, QRNan )
382 {
383  SharedDiagonal constrained = noiseModel::Constrained::All(2);
384  Matrix Ab = (Matrix25() << 2, 4, 2, 4, 6, 2, 1, 2, 4, 4).finished();
385 
386  SharedDiagonal expected = noiseModel::Constrained::All(2);
387  Matrix expectedAb = (Matrix25() << 1, 2, 1, 2, 3, 0, 1, 0, 0, 2.0/3).finished();
388 
389  SharedDiagonal actual = constrained->QR(Ab);
390  EXPECT(assert_equal(*expected,*actual));
391  EXPECT(linear_dependent(expectedAb,Ab));
392 }
393 
394 /* ************************************************************************* */
395 TEST(NoiseModel, SmartSqrtInformation )
396 {
397  bool smart = true;
398  gtsam::SharedGaussian expected = Unit::Create(3);
399  gtsam::SharedGaussian actual = Gaussian::SqrtInformation(I_3x3, smart);
400  EXPECT(assert_equal(*expected,*actual));
401 }
402 
403 /* ************************************************************************* */
404 TEST(NoiseModel, SmartSqrtInformation2 )
405 {
406  bool smart = true;
407  gtsam::SharedGaussian expected = Unit::Isotropic::Sigma(3,2);
408  gtsam::SharedGaussian actual = Gaussian::SqrtInformation(0.5*I_3x3, smart);
409  EXPECT(assert_equal(*expected,*actual));
410 }
411 
412 /* ************************************************************************* */
413 TEST(NoiseModel, SmartInformation )
414 {
415  bool smart = true;
416  gtsam::SharedGaussian expected = Unit::Isotropic::Variance(3,2);
417  Matrix M = 0.5*I_3x3;
419  gtsam::SharedGaussian actual = Gaussian::Information(M, smart);
420  EXPECT(assert_equal(*expected,*actual));
421 }
422 
423 /* ************************************************************************* */
424 TEST(NoiseModel, SmartCovariance )
425 {
426  bool smart = true;
427  gtsam::SharedGaussian expected = Unit::Create(3);
428  gtsam::SharedGaussian actual = Gaussian::Covariance(I_3x3, smart);
429  EXPECT(assert_equal(*expected,*actual));
430 }
431 
432 /* ************************************************************************* */
433 TEST(NoiseModel, ScalarOrVector )
434 {
435  bool smart = true;
436  SharedGaussian expected = Unit::Create(3);
437  SharedGaussian actual = Gaussian::Covariance(I_3x3, smart);
438  EXPECT(assert_equal(*expected,*actual));
439 }
440 
441 /* ************************************************************************* */
442 TEST(NoiseModel, WhitenInPlace)
443 {
444  Vector sigmas = Vector3(0.1, 0.1, 0.1);
445  SharedDiagonal model = Diagonal::Sigmas(sigmas);
446  Matrix A = I_3x3;
447  model->WhitenInPlace(A);
448  Matrix expected = I_3x3 * 10;
449  EXPECT(assert_equal(expected, A));
450 }
451 
452 /* ************************************************************************* */
453 
454 /*
455  * These tests are responsible for testing the weight functions for the m-estimators in GTSAM.
456  * The weight function is related to the analytic derivative of the loss function. See
457  * https://members.loria.fr/MOBerger/Enseignement/Master2/Documents/ZhangIVC-97-01.pdf
458  * for details. This weight function is required when optimizing cost functions with robust
459  * penalties using iteratively re-weighted least squares.
460  */
461 
462 TEST(NoiseModel, robustFunctionFair)
463 {
464  const double k = 5.0, error1 = 1.0, error2 = 10.0, error3 = -10.0, error4 = -1.0;
465  const mEstimator::Fair::shared_ptr fair = mEstimator::Fair::Create(k);
466  DOUBLES_EQUAL(0.8333333333333333, fair->weight(error1), 1e-8);
467  DOUBLES_EQUAL(0.3333333333333333, fair->weight(error2), 1e-8);
468  // Test negative value to ensure we take absolute value of error.
469  DOUBLES_EQUAL(0.3333333333333333, fair->weight(error3), 1e-8);
470  DOUBLES_EQUAL(0.8333333333333333, fair->weight(error4), 1e-8);
471 
472  DOUBLES_EQUAL(0.441961080151135, fair->loss(error1), 1e-8);
473  DOUBLES_EQUAL(22.534692783297260, fair->loss(error2), 1e-8);
474  DOUBLES_EQUAL(22.534692783297260, fair->loss(error3), 1e-8);
475  DOUBLES_EQUAL(0.441961080151135, fair->loss(error4), 1e-8);
476 }
477 
478 TEST(NoiseModel, robustFunctionHuber)
479 {
480  const double k = 5.0, error1 = 1.0, error2 = 10.0, error3 = -10.0, error4 = -1.0;
481  const mEstimator::Huber::shared_ptr huber = mEstimator::Huber::Create(k);
482  DOUBLES_EQUAL(1.0, huber->weight(error1), 1e-8);
483  DOUBLES_EQUAL(0.5, huber->weight(error2), 1e-8);
484  // Test negative value to ensure we take absolute value of error.
485  DOUBLES_EQUAL(0.5, huber->weight(error3), 1e-8);
486  DOUBLES_EQUAL(1.0, huber->weight(error4), 1e-8);
487 
488  DOUBLES_EQUAL(0.5000, huber->loss(error1), 1e-8);
489  DOUBLES_EQUAL(37.5000, huber->loss(error2), 1e-8);
490  DOUBLES_EQUAL(37.5000, huber->loss(error3), 1e-8);
491  DOUBLES_EQUAL(0.5000, huber->loss(error4), 1e-8);
492 }
493 
494 TEST(NoiseModel, robustFunctionCauchy)
495 {
496  const double k = 5.0, error1 = 1.0, error2 = 10.0, error3 = -10.0, error4 = -1.0;
497  const mEstimator::Cauchy::shared_ptr cauchy = mEstimator::Cauchy::Create(k);
498  DOUBLES_EQUAL(0.961538461538461, cauchy->weight(error1), 1e-8);
499  DOUBLES_EQUAL(0.2000, cauchy->weight(error2), 1e-8);
500  // Test negative value to ensure we take absolute value of error.
501  DOUBLES_EQUAL(0.2000, cauchy->weight(error3), 1e-8);
502  DOUBLES_EQUAL(0.961538461538461, cauchy->weight(error4), 1e-8);
503 
504  DOUBLES_EQUAL(0.490258914416017, cauchy->loss(error1), 1e-8);
505  DOUBLES_EQUAL(20.117973905426254, cauchy->loss(error2), 1e-8);
506  DOUBLES_EQUAL(20.117973905426254, cauchy->loss(error3), 1e-8);
507  DOUBLES_EQUAL(0.490258914416017, cauchy->loss(error4), 1e-8);
508 }
509 
510 TEST(NoiseModel, robustFunctionGemanMcClure)
511 {
512  const double k = 1.0, error1 = 1.0, error2 = 10.0, error3 = -10.0, error4 = -1.0;
513  const mEstimator::GemanMcClure::shared_ptr gmc = mEstimator::GemanMcClure::Create(k);
514  DOUBLES_EQUAL(0.25 , gmc->weight(error1), 1e-8);
515  DOUBLES_EQUAL(9.80296e-5, gmc->weight(error2), 1e-8);
516  DOUBLES_EQUAL(9.80296e-5, gmc->weight(error3), 1e-8);
517  DOUBLES_EQUAL(0.25 , gmc->weight(error4), 1e-8);
518 
519  DOUBLES_EQUAL(0.2500, gmc->loss(error1), 1e-8);
520  DOUBLES_EQUAL(0.495049504950495, gmc->loss(error2), 1e-8);
521  DOUBLES_EQUAL(0.495049504950495, gmc->loss(error3), 1e-8);
522  DOUBLES_EQUAL(0.2500, gmc->loss(error4), 1e-8);
523 }
524 
525 TEST(NoiseModel, robustFunctionWelsch)
526 {
527  const double k = 5.0, error1 = 1.0, error2 = 10.0, error3 = -10.0, error4 = -1.0;
528  const mEstimator::Welsch::shared_ptr welsch = mEstimator::Welsch::Create(k);
529  DOUBLES_EQUAL(0.960789439152323, welsch->weight(error1), 1e-8);
530  DOUBLES_EQUAL(0.018315638888734, welsch->weight(error2), 1e-8);
531  // Test negative value to ensure we take absolute value of error.
532  DOUBLES_EQUAL(0.018315638888734, welsch->weight(error3), 1e-8);
533  DOUBLES_EQUAL(0.960789439152323, welsch->weight(error4), 1e-8);
534 
535  DOUBLES_EQUAL(0.490132010595960, welsch->loss(error1), 1e-8);
536  DOUBLES_EQUAL(12.271054513890823, welsch->loss(error2), 1e-8);
537  DOUBLES_EQUAL(12.271054513890823, welsch->loss(error3), 1e-8);
538  DOUBLES_EQUAL(0.490132010595960, welsch->loss(error4), 1e-8);
539 }
540 
541 TEST(NoiseModel, robustFunctionTukey)
542 {
543  const double k = 5.0, error1 = 1.0, error2 = 10.0, error3 = -10.0, error4 = -1.0;
544  const mEstimator::Tukey::shared_ptr tukey = mEstimator::Tukey::Create(k);
545  DOUBLES_EQUAL(0.9216, tukey->weight(error1), 1e-8);
546  DOUBLES_EQUAL(0.0, tukey->weight(error2), 1e-8);
547  // Test negative value to ensure we take absolute value of error.
548  DOUBLES_EQUAL(0.0, tukey->weight(error3), 1e-8);
549  DOUBLES_EQUAL(0.9216, tukey->weight(error4), 1e-8);
550 
551  DOUBLES_EQUAL(0.480266666666667, tukey->loss(error1), 1e-8);
552  DOUBLES_EQUAL(4.166666666666667, tukey->loss(error2), 1e-8);
553  DOUBLES_EQUAL(4.166666666666667, tukey->loss(error3), 1e-8);
554  DOUBLES_EQUAL(0.480266666666667, tukey->loss(error4), 1e-8);
555 }
556 
557 TEST(NoiseModel, robustFunctionDCS)
558 {
559  const double k = 1.0, error1 = 1.0, error2 = 10.0;
560  const mEstimator::DCS::shared_ptr dcs = mEstimator::DCS::Create(k);
561 
562  DOUBLES_EQUAL(1.0 , dcs->weight(error1), 1e-8);
563  DOUBLES_EQUAL(0.00039211, dcs->weight(error2), 1e-8);
564 
565  DOUBLES_EQUAL(0.5 , dcs->loss(error1), 1e-8);
566  DOUBLES_EQUAL(0.9900990099, dcs->loss(error2), 1e-8);
567 }
568 
569 TEST(NoiseModel, robustFunctionL2WithDeadZone)
570 {
571  const double k = 1.0, e0 = -10.0, e1 = -1.01, e2 = -0.99, e3 = 0.99, e4 = 1.01, e5 = 10.0;
572  const mEstimator::L2WithDeadZone::shared_ptr lsdz = mEstimator::L2WithDeadZone::Create(k);
573 
574  DOUBLES_EQUAL(0.9, lsdz->weight(e0), 1e-8);
575  DOUBLES_EQUAL(0.00990099009, lsdz->weight(e1), 1e-8);
576  DOUBLES_EQUAL(0.0, lsdz->weight(e2), 1e-8);
577  DOUBLES_EQUAL(0.0, lsdz->weight(e3), 1e-8);
578  DOUBLES_EQUAL(0.00990099009, lsdz->weight(e4), 1e-8);
579  DOUBLES_EQUAL(0.9, lsdz->weight(e5), 1e-8);
580 
581  DOUBLES_EQUAL(40.5, lsdz->loss(e0), 1e-8);
582  DOUBLES_EQUAL(0.00005, lsdz->loss(e1), 1e-8);
583  DOUBLES_EQUAL(0.0, lsdz->loss(e2), 1e-8);
584  DOUBLES_EQUAL(0.0, lsdz->loss(e3), 1e-8);
585  DOUBLES_EQUAL(0.00005, lsdz->loss(e4), 1e-8);
586  DOUBLES_EQUAL(40.5, lsdz->loss(e5), 1e-8);
587 }
588 
589 /* ************************************************************************* */
590 TEST(NoiseModel, robustNoiseHuber)
591 {
592  const double k = 10.0, error1 = 1.0, error2 = 100.0;
593  Matrix A = (Matrix(2, 2) << 1.0, 10.0, 100.0, 1000.0).finished();
594  Vector b = Vector2(error1, error2);
595  const Robust::shared_ptr robust = Robust::Create(
596  mEstimator::Huber::Create(k, mEstimator::Huber::Scalar),
597  Unit::Create(2));
598 
599  robust->WhitenSystem(A, b);
600 
601  DOUBLES_EQUAL(error1, b(0), 1e-8);
602  DOUBLES_EQUAL(sqrt(k*error2), b(1), 1e-8);
603 
604  DOUBLES_EQUAL(1.0, A(0,0), 1e-8);
605  DOUBLES_EQUAL(10.0, A(0,1), 1e-8);
606  DOUBLES_EQUAL(sqrt(k*100.0), A(1,0), 1e-8);
607  DOUBLES_EQUAL(sqrt(k/100.0)*1000.0, A(1,1), 1e-8);
608 }
609 
610 TEST(NoiseModel, robustNoiseGemanMcClure)
611 {
612  const double k = 1.0, error1 = 1.0, error2 = 100.0;
613  const double a00 = 1.0, a01 = 10.0, a10 = 100.0, a11 = 1000.0;
614  Matrix A = (Matrix(2, 2) << a00, a01, a10, a11).finished();
615  Vector b = Vector2(error1, error2);
616  const Robust::shared_ptr robust = Robust::Create(
617  mEstimator::GemanMcClure::Create(k, mEstimator::GemanMcClure::Scalar),
618  Unit::Create(2));
619 
620  robust->WhitenSystem(A, b);
621 
622  const double k2 = k*k;
623  const double k4 = k2*k2;
624  const double k2error = k2 + error2*error2;
625 
626  const double sqrt_weight_error1 = sqrt(0.25);
627  const double sqrt_weight_error2 = sqrt(k4/(k2error*k2error));
628 
629  DOUBLES_EQUAL(sqrt_weight_error1*error1, b(0), 1e-8);
630  DOUBLES_EQUAL(sqrt_weight_error2*error2, b(1), 1e-8);
631 
632  DOUBLES_EQUAL(sqrt_weight_error1*a00, A(0,0), 1e-8);
633  DOUBLES_EQUAL(sqrt_weight_error1*a01, A(0,1), 1e-8);
634  DOUBLES_EQUAL(sqrt_weight_error2*a10, A(1,0), 1e-8);
635  DOUBLES_EQUAL(sqrt_weight_error2*a11, A(1,1), 1e-8);
636 }
637 
638 TEST(NoiseModel, robustNoiseDCS)
639 {
640  const double k = 1.0, error1 = 1.0, error2 = 100.0;
641  const double a00 = 1.0, a01 = 10.0, a10 = 100.0, a11 = 1000.0;
642  Matrix A = (Matrix(2, 2) << a00, a01, a10, a11).finished();
643  Vector b = Vector2(error1, error2);
644  const Robust::shared_ptr robust = Robust::Create(
645  mEstimator::DCS::Create(k, mEstimator::DCS::Scalar),
646  Unit::Create(2));
647 
648  robust->WhitenSystem(A, b);
649 
650  const double sqrt_weight = 2.0*k/(k + error2*error2);
651 
652  DOUBLES_EQUAL(error1, b(0), 1e-8);
653  DOUBLES_EQUAL(sqrt_weight*error2, b(1), 1e-8);
654 
655  DOUBLES_EQUAL(a00, A(0,0), 1e-8);
656  DOUBLES_EQUAL(a01, A(0,1), 1e-8);
657  DOUBLES_EQUAL(sqrt_weight*a10, A(1,0), 1e-8);
658  DOUBLES_EQUAL(sqrt_weight*a11, A(1,1), 1e-8);
659 }
660 
661 TEST(NoiseModel, robustNoiseL2WithDeadZone)
662 {
663  double dead_zone_size = 1.0;
664  SharedNoiseModel robust = noiseModel::Robust::Create(
665  noiseModel::mEstimator::L2WithDeadZone::Create(dead_zone_size),
666  Unit::Create(3));
667 
668 /*
669  * TODO(mike): There is currently a bug in GTSAM, where none of the mEstimator classes
670  * implement a loss function, and GTSAM calls the weight function to evaluate the
671  * total penalty, rather than calling the loss function. The weight function should be
672  * used during iteratively reweighted least squares optimization, but should not be used to
673  * evaluate the total penalty. The long-term solution is for all mEstimators to implement
674  * both a weight and a loss function, and for GTSAM to call the loss function when
675  * evaluating the total penalty. This bug causes the test below to fail, so I'm leaving it
676  * commented out until the underlying bug in GTSAM is fixed.
677  *
678  * for (int i = 0; i < 5; i++) {
679  * Vector3 error = Vector3(i, 0, 0);
680  * DOUBLES_EQUAL(0.5*max(0,i-1)*max(0,i-1), robust->distance(error), 1e-8);
681  * }
682  */
683 
684 }
685 
686 TEST(NoiseModel, lossFunctionAtZero)
687 {
688  const double k = 5.0;
689  auto fair = mEstimator::Fair::Create(k);
690  DOUBLES_EQUAL(fair->loss(0), 0, 1e-8);
691  DOUBLES_EQUAL(fair->weight(0), 1, 1e-8);
692  auto huber = mEstimator::Huber::Create(k);
693  DOUBLES_EQUAL(huber->loss(0), 0, 1e-8);
694  DOUBLES_EQUAL(huber->weight(0), 1, 1e-8);
695  auto cauchy = mEstimator::Cauchy::Create(k);
696  DOUBLES_EQUAL(cauchy->loss(0), 0, 1e-8);
697  DOUBLES_EQUAL(cauchy->weight(0), 1, 1e-8);
698  auto gmc = mEstimator::GemanMcClure::Create(k);
699  DOUBLES_EQUAL(gmc->loss(0), 0, 1e-8);
700  DOUBLES_EQUAL(gmc->weight(0), 1, 1e-8);
701  auto welsch = mEstimator::Welsch::Create(k);
702  DOUBLES_EQUAL(welsch->loss(0), 0, 1e-8);
703  DOUBLES_EQUAL(welsch->weight(0), 1, 1e-8);
704  auto tukey = mEstimator::Tukey::Create(k);
705  DOUBLES_EQUAL(tukey->loss(0), 0, 1e-8);
706  DOUBLES_EQUAL(tukey->weight(0), 1, 1e-8);
707  auto dcs = mEstimator::DCS::Create(k);
708  DOUBLES_EQUAL(dcs->loss(0), 0, 1e-8);
709  DOUBLES_EQUAL(dcs->weight(0), 1, 1e-8);
710  // auto lsdz = mEstimator::L2WithDeadZone::Create(k);
711  // DOUBLES_EQUAL(lsdz->loss(0), 0, 1e-8);
712  // DOUBLES_EQUAL(lsdz->weight(0), 1, 1e-8);
713 }
714 
715 
716 /* ************************************************************************* */
717 #define TEST_GAUSSIAN(gaussian)\
718  EQUALITY(info, gaussian->information());\
719  EQUALITY(cov, gaussian->covariance());\
720  EXPECT(assert_equal(white, gaussian->whiten(e)));\
721  EXPECT(assert_equal(e, gaussian->unwhiten(white)));\
722  EXPECT_DOUBLES_EQUAL(251.0, gaussian->squaredMahalanobisDistance(e), 1e-9);\
723  EXPECT_DOUBLES_EQUAL(0.5 * 251.0, gaussian->loss(251.0), 1e-9);\
724  Matrix A = R.inverse(); Vector b = e;\
725  gaussian->WhitenSystem(A, b);\
726  EXPECT(assert_equal(I, A));\
727  EXPECT(assert_equal(white, b));
728 
729 TEST(NoiseModel, NonDiagonalGaussian)
730 {
731  Matrix3 R;
732  R << 6, 5, 4, 0, 3, 2, 0, 0, 1;
733  const Matrix3 info = R.transpose() * R;
734  const Matrix3 cov = info.inverse();
735  const Vector3 e(1, 1, 1), white = R * e;
736  Matrix I = Matrix3::Identity();
737 
738 
739  {
740  SharedGaussian gaussian = Gaussian::SqrtInformation(R);
741  TEST_GAUSSIAN(gaussian);
742  }
743 
744  {
745  SharedGaussian gaussian = Gaussian::Information(info);
746  TEST_GAUSSIAN(gaussian);
747  }
748 
749  {
750  SharedGaussian gaussian = Gaussian::Covariance(cov);
751  TEST_GAUSSIAN(gaussian);
752  }
753 }
754 
755 /* ************************************************************************* */
756 int main() { TestResult tr; return TestRegistry::runAllTests(tr); }
757 /* ************************************************************************* */
Matrix3f m
boost::shared_ptr< Huber > shared_ptr
Provides additional testing facilities for common data structures.
SCALAR Scalar
Definition: bench_gemm.cpp:33
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:38
static const Matrix kCovariance
Scalar * b
Definition: benchVecAdd.cpp:17
static int runAllTests(TestResult &result)
Eigen::Vector3d Vector3
Definition: Vector.h:43
noiseModel::Diagonal::shared_ptr model
Matrix expected
Definition: testMatrix.cpp:974
boost::shared_ptr< Robust > shared_ptr
Definition: NoiseModel.h:655
#define DOUBLES_EQUAL(expected, actual, threshold)
Definition: Test.h:142
double mu
ArrayXcf v
Definition: Cwise_arg.cpp:1
void diagonal(const MatrixType &m)
Definition: diagonal.cpp:12
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:43
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Definition: Half.h:150
TEST(NoiseModel, constructors)
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 set xzeroaxis lt lw set x2zeroaxis lt lw set yzeroaxis lt lw set y2zeroaxis lt lw set tics in set ticslevel set tics set mxtics default set mytics default set mx2tics default set my2tics default set xtics border mirror norotate autofreq set ytics border mirror norotate autofreq set ztics border nomirror norotate autofreq set nox2tics set noy2tics set timestamp bottom norotate set rrange[*:*] noreverse nowriteback set trange[*:*] noreverse nowriteback set urange[*:*] noreverse nowriteback set vrange[*:*] noreverse nowriteback set xlabel matrix size set x2label set timefmt d m y n H
static const double kVariance
static const double prc
boost::shared_ptr< L2WithDeadZone > shared_ptr
static const double kSigma
static const Matrix R
else if n * info
Eigen::VectorXd Vector
Definition: Vector.h:38
boost::shared_ptr< GemanMcClure > shared_ptr
#define EXPECT(condition)
Definition: Test.h:151
static const Vector3 kSigmas(kSigma, kSigma, kSigma)
boost::shared_ptr< Tukey > shared_ptr
Array< double, 1, 3 > e(1./3., 0.5, 2.)
#define TEST_GAUSSIAN(gaussian)
static const Matrix I
Definition: lago.cpp:35
boost::shared_ptr< Cauchy > shared_ptr
boost::shared_ptr< Diagonal > shared_ptr
Definition: NoiseModel.h:301
boost::shared_ptr< Constrained > shared_ptr
Definition: NoiseModel.h:407
SharedDiagonal diagonal
noiseModel::Diagonal::shared_ptr SharedDiagonal
Definition: NoiseModel.h:736
traits
Definition: chartTesting.h:28
bool linear_dependent(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:116
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:42
Eigen::Vector2d Vector2
Definition: Vector.h:42
boost::shared_ptr< DCS > shared_ptr
boost::optional< Vector > checkIfDiagonal(const Matrix M)
Definition: NoiseModel.cpp:50
Pose3 g1(Rot3(), Point3(100.0, 0.0, 300.0))
bool assert_inequal(const Matrix &A, const Matrix &B, double tol)
Definition: Matrix.cpp:62
boost::shared_ptr< Isotropic > shared_ptr
Definition: NoiseModel.h:541
boost::shared_ptr< Gaussian > shared_ptr
Definition: NoiseModel.h:189
static const double kInverseSigma
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > Matrix
boost::shared_ptr< Fair > shared_ptr
Pose3 g2(g1.expmap(h *V1_g1))
boost::shared_ptr< Welsch > shared_ptr
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
int main()
noiseModel::Gaussian::shared_ptr SharedGaussian
Definition: NoiseModel.h:735
noiseModel::Base::shared_ptr SharedNoiseModel
Definition: NoiseModel.h:734


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:48:13