testRegularJacobianFactor.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 
24 
26 
27 #include <boost/assign/std/vector.hpp>
28 #include <boost/assign/list_of.hpp>
29 #include <boost/range/iterator_range.hpp>
30 #include <boost/range/adaptor/map.hpp>
31 
32 using namespace std;
33 using namespace gtsam;
34 using namespace boost::assign;
35 
36 static const size_t fixedDim = 3;
37 static const size_t nrKeys = 3;
38 
39 // Keys are assumed to be from 0 to n
40 namespace {
41  namespace simple {
42  // Terms we'll use
43  const vector<pair<Key, Matrix> > terms = list_of<pair<Key,Matrix> >
44  (make_pair(0, Matrix3::Identity()))
45  (make_pair(1, 2*Matrix3::Identity()))
46  (make_pair(2, 3*Matrix3::Identity()));
47 
48  // RHS and sigmas
49  const Vector b = (Vector(3) << 1., 2., 3.).finished();
50  const SharedDiagonal noise = noiseModel::Diagonal::Sigmas((Vector(3) << 0.5,0.5,0.5).finished());
51  }
52 
53  namespace simple2 {
54  // Terms
55  const vector<pair<Key, Matrix> > terms2 = list_of<pair<Key,Matrix> >
56  (make_pair(0, 2*Matrix3::Identity()))
57  (make_pair(1, 4*Matrix3::Identity()))
58  (make_pair(2, 6*Matrix3::Identity()));
59 
60  // RHS
61  const Vector b2 = (Vector(3) << 2., 4., 6.).finished();
62  }
63 }
64 
65 /* ************************************************************************* */
66 // Convert from double* to VectorValues
67 VectorValues double2vv(const double* x,
68  const size_t nrKeys, const size_t dim) {
69  // create map with dimensions
70  std::map<gtsam::Key, size_t> dims;
71  for (size_t i = 0; i < nrKeys; i++)
72  dims.insert(make_pair(i, dim));
73 
74  size_t n = nrKeys*dim;
75  Vector xVec(n);
76  for (size_t i = 0; i < n; i++){
77  xVec(i) = x[i];
78  }
79  return VectorValues(xVec, dims);
80 }
81 
82 /* ************************************************************************* */
83 void vv2double(const VectorValues& vv, double* y,
84  const size_t nrKeys, const size_t dim) {
85  // create map with dimensions
86  std::map<gtsam::Key, size_t> dims;
87  for (size_t i = 0; i < nrKeys; i++)
88  dims.insert(make_pair(i, dim));
89 
90  Vector yvector = vv.vector(dims);
91  size_t n = nrKeys*dim;
92  for (size_t j = 0; j < n; j++)
93  y[j] = yvector[j];
94 }
95 
96 /* ************************************************************************* */
97 TEST(RegularJacobianFactor, constructorNway)
98 {
99  using namespace simple;
100  JacobianFactor factor(terms[0].first, terms[0].second,
101  terms[1].first, terms[1].second, terms[2].first, terms[2].second, b, noise);
102  RegularJacobianFactor<fixedDim> regularFactor(terms, b, noise);
103 
104  LONGS_EQUAL((long)terms[2].first, (long)regularFactor.keys().back());
105  EXPECT(assert_equal(terms[2].second, regularFactor.getA(regularFactor.end() - 1)));
106  EXPECT(assert_equal(b, factor.getb()));
107  EXPECT(assert_equal(b, regularFactor.getb()));
108  EXPECT(noise == factor.get_model());
109  EXPECT(noise == regularFactor.get_model());
110 }
111 
112 /* ************************************************************************* */
113 TEST(RegularJacobianFactor, hessianDiagonal)
114 {
115  using namespace simple;
116  JacobianFactor factor(terms[0].first, terms[0].second,
117  terms[1].first, terms[1].second, terms[2].first, terms[2].second, b, noise);
118  RegularJacobianFactor<fixedDim> regularFactor(terms, b, noise);
119 
120  // we compute hessian diagonal from the standard Jacobian
121  VectorValues expectedHessianDiagonal = factor.hessianDiagonal();
122 
123  // we compare against the Raw memory access implementation of hessianDiagonal
124  double actualValue[9]={0};
125  regularFactor.hessianDiagonal(actualValue);
126  VectorValues actualHessianDiagonalRaw = double2vv(actualValue,nrKeys,fixedDim);
127  EXPECT(assert_equal(expectedHessianDiagonal, actualHessianDiagonalRaw));
128 }
129 
130 /* ************************************************************************* */
131 TEST(RegularJacobian, gradientAtZero)
132 {
133  using namespace simple;
134  JacobianFactor factor(terms[0].first, terms[0].second,
135  terms[1].first, terms[1].second, terms[2].first, terms[2].second, b, noise);
136  RegularJacobianFactor<fixedDim> regularFactor(terms, b, noise);
137 
138  // we compute gradient at zero from the standard Jacobian
139  VectorValues expectedGradientAtZero = factor.gradientAtZero();
140 
141  //EXPECT(assert_equal(expectedGradientAtZero, regularFactor.gradientAtZero()));
142 
143  // we compare against the Raw memory access implementation of gradientAtZero
144  double actualValue[9]={0};
145  regularFactor.gradientAtZero(actualValue);
146  VectorValues actualGradientAtZeroRaw = double2vv(actualValue,nrKeys,fixedDim);
147  EXPECT(assert_equal(expectedGradientAtZero, actualGradientAtZeroRaw));
148 }
149 
150 /* ************************************************************************* */
151 TEST(RegularJacobian, gradientAtZero_multiFactors)
152 {
153  using namespace simple;
154  JacobianFactor factor(terms[0].first, terms[0].second,
155  terms[1].first, terms[1].second, terms[2].first, terms[2].second, b, noise);
156  RegularJacobianFactor<fixedDim> regularFactor(terms, b, noise);
157 
158  // we compute gradient at zero from the standard Jacobian
159  VectorValues expectedGradientAtZero = factor.gradientAtZero();
160 
161  // we compare against the Raw memory access implementation of gradientAtZero
162  double actualValue[9]={0};
163  regularFactor.gradientAtZero(actualValue);
164  VectorValues actualGradientAtZeroRaw = double2vv(actualValue,nrKeys,fixedDim);
165  EXPECT(assert_equal(expectedGradientAtZero, actualGradientAtZeroRaw));
166 
167  // One more factor
168  using namespace simple2;
169  JacobianFactor factor2(terms2[0].first, terms2[0].second,
170  terms2[1].first, terms2[1].second, terms2[2].first, terms2[2].second, b2, noise);
171  RegularJacobianFactor<fixedDim> regularFactor2(terms2, b2, noise);
172 
173  // we accumulate computed gradient at zero from the standard Jacobian
174  VectorValues expectedGradientAtZero2 = expectedGradientAtZero.add(factor2.gradientAtZero());
175 
176  // we compare against the Raw memory access implementation of gradientAtZero
177  regularFactor2.gradientAtZero(actualValue);
178  VectorValues actualGradientAtZeroRaw2 = double2vv(actualValue,nrKeys,fixedDim);
179  EXPECT(assert_equal(expectedGradientAtZero2, actualGradientAtZeroRaw2));
180 
181 }
182 
183 /* ************************************************************************* */
184 TEST(RegularJacobian, multiplyHessianAdd)
185 {
186  using namespace simple;
187  JacobianFactor factor(terms[0].first, terms[0].second,
188  terms[1].first, terms[1].second, terms[2].first, terms[2].second, b, noise);
189  RegularJacobianFactor<fixedDim> regularFactor(terms, b, noise);
190 
191  // arbitrary vector X
192  VectorValues X;
193  X.insert(0, (Vector(3) << 10.,20.,30.).finished());
194  X.insert(1, (Vector(3) << 10.,20.,30.).finished());
195  X.insert(2, (Vector(3) << 10.,20.,30.).finished());
196 
197  // arbitrary vector Y
198  VectorValues Y;
199  Y.insert(0, (Vector(3) << 10.,10.,10.).finished());
200  Y.insert(1, (Vector(3) << 20.,20.,20.).finished());
201  Y.insert(2, (Vector(3) << 30.,30.,30.).finished());
202 
203  // multiplyHessianAdd Y += alpha*A'A*X
204  double alpha = 2.0;
205  VectorValues expectedMHA = Y;
206  factor.multiplyHessianAdd(alpha, X, expectedMHA);
207 
208  // create data for raw memory access
209  double XRaw[9];
210  vv2double(X, XRaw, nrKeys, fixedDim);
211 
212  // test 1st version: multiplyHessianAdd(double alpha, const double* x, double* y)
213  double actualMHARaw[9];
214  vv2double(Y, actualMHARaw, nrKeys, fixedDim);
215  regularFactor.multiplyHessianAdd(alpha, XRaw, actualMHARaw);
216  VectorValues actualMHARawVV = double2vv(actualMHARaw,nrKeys,fixedDim);
217  EXPECT(assert_equal(expectedMHA,actualMHARawVV));
218 
219  // test 2nd version: multiplyHessianAdd(double alpha, const double* x, double* y, std::vector<size_t> keys)
220  double actualMHARaw2[9];
221  vv2double(Y, actualMHARaw2, nrKeys, fixedDim);
222  vector<size_t> dims;
223  size_t accumulatedDim = 0;
224  for (size_t i = 0; i < nrKeys+1; i++){
225  dims.push_back(accumulatedDim);
226  accumulatedDim += fixedDim;
227  }
228  regularFactor.multiplyHessianAdd(alpha, XRaw, actualMHARaw2, dims);
229  VectorValues actualMHARawVV2 = double2vv(actualMHARaw2,nrKeys,fixedDim);
230  EXPECT(assert_equal(expectedMHA,actualMHARawVV2));
231 }
232 
233 /* ************************************************************************* */
234 int main() { TestResult tr; return TestRegistry::runAllTests(tr);}
235 /* ************************************************************************* */
Provides additional testing facilities for common data structures.
static const size_t fixedDim
Scalar * y
Scalar * b
Definition: benchVecAdd.cpp:17
static int runAllTests(TestResult &result)
int n
const_iterator end() const
Definition: Factor.h:130
void multiplyHessianAdd(double alpha, const VectorValues &x, VectorValues &y) const override
Definition: Half.h:150
TEST(RegularJacobianFactor, constructorNway)
static const size_t nrKeys
void multiplyHessianAdd(double alpha, const VectorValues &x, VectorValues &y) const override
Factor Graph Values.
VectorValues gradientAtZero() const override
Expose base class gradientAtZero.
constexpr int first(int i)
Implementation details for constexpr functions.
Eigen::VectorXd Vector
Definition: Vector.h:38
const constBVector getb() const
#define EXPECT(condition)
Definition: Test.h:151
Vector2 b2(4,-5)
VectorValues add(const VectorValues &c) const
RealScalar alpha
JacobianFactor class with fixed sized blcoks.
Conditional Gaussian Base class.
Linear Factor Graph where all factors are Gaussians.
constABlock getA(const_iterator variable) const
const SharedDiagonal & get_model() const
#define LONGS_EQUAL(expected, actual)
Definition: Test.h:135
const mpreal dim(const mpreal &a, const mpreal &b, mp_rnd_t r=mpreal::get_default_rnd())
Definition: mpreal.h:2201
noiseModel::Diagonal::shared_ptr SharedDiagonal
Definition: NoiseModel.h:736
traits
Definition: chartTesting.h:28
bool assert_equal(const Matrix &expected, const Matrix &actual, double tol)
Definition: Matrix.cpp:42
void vv2double(const VectorValues &vv, double *y, const size_t nrKeys, const size_t dim)
const KeyVector & keys() const
Access the factor&#39;s involved variable keys.
Definition: Factor.h:124
void hessianDiagonal(double *d) const override
Raw memory access version of hessianDiagonal.
JacobianFactor factor2(keyX, A21, keyY, A22, b2, noiseModel::Isotropic::Sigma(2, sigma2))
VectorValues gradientAtZero() const override
A&#39;*b for Jacobian.
#define X
Definition: icosphere.cpp:20
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
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
iterator insert(const std::pair< Key, Vector > &key_value)
std::ptrdiff_t j
VectorValues double2vv(const double *x, const size_t nrKeys, const size_t dim)
void hessianDiagonal(double *d) const override
Raw memory access version of hessianDiagonal.
Vector vector() const


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:49:17