timeSchurFactors.cpp
Go to the documentation of this file.
1 
8 #include "DummyFactor.h"
9 #include <gtsam/base/timing.h>
10 
16 
17 #include <fstream>
18 
19 using namespace std;
20 using namespace gtsam;
21 
22 #define SLOW
23 #define RAW
24 #define HESSIAN
25 #define NUM_ITERATIONS 1000
26 
27 // Create CSV file for results
28 ofstream os("timeSchurFactors.csv");
29 
30 /*************************************************************************************/
31 template<typename CAMERA>
32 void timeAll(size_t m, size_t N) {
33 
34  cout << m << endl;
35 
36  // create F
37  static const int D = CAMERA::dimension;
38  typedef Eigen::Matrix<double, 2, D> Matrix2D;
40  vector <Matrix2D, Eigen::aligned_allocator<Matrix2D>> Fblocks;
41  for (size_t i = 0; i < m; i++) {
42  keys.push_back(i);
43  Fblocks.push_back((i + 1) * Matrix::Ones(2, D));
44  }
45 
46  // create E
47  Matrix E(2 * m, 3);
48  for (size_t i = 0; i < m; i++)
49  E.block < 2, 3 > (2 * i, 0) = Matrix::Ones(2, 3);
50 
51  // Calculate point covariance
52  Matrix P = (E.transpose() * E).inverse();
53 
54  // RHS and sigmas
55  const Vector b = Vector::Constant(2*m,1);
56  const SharedDiagonal model;
57 
58  // parameters for multiplyHessianAdd
59  double alpha = 0.5;
60  VectorValues xvalues, yvalues;
61  for (size_t i = 0; i < m; i++)
62  xvalues.insert(i, Vector::Constant(D,2));
63 
64  // Implicit
65  RegularImplicitSchurFactor<CAMERA> implicitFactor(keys, Fblocks, E, P, b);
66  // JacobianFactor with same error
67  JacobianFactorQ<D, 2> jf(keys, Fblocks, E, P, b, model);
68  // JacobianFactorQR with same error
69  JacobianFactorQR<D, 2> jqr(keys, Fblocks, E, P, b, model);
70  // Hessian
71  HessianFactor hessianFactor(jqr);
72 
73 #define TIME(label,factor,xx,yy) {\
74  for (size_t t = 0; t < N; t++) \
75  factor.multiplyHessianAdd(alpha, xx, yy);\
76  gttic_(label);\
77  for (size_t t = 0; t < N; t++) {\
78  factor.multiplyHessianAdd(alpha, xx, yy);\
79  }\
80  gttoc_(label);\
81  tictoc_getNode(timer, label)\
82  os << timer->secs()/NUM_ITERATIONS << ", ";\
83  }
84 
85 #ifdef SLOW
86  TIME(Implicit, implicitFactor, xvalues, yvalues)
87  TIME(Jacobian, jf, xvalues, yvalues)
88  TIME(JacobianQR, jqr, xvalues, yvalues)
89 #endif
90 
91 #ifdef HESSIAN
92  TIME(Hessian, hessianFactor, xvalues, yvalues)
93 #endif
94 
95 #ifdef OVERHEAD
96  DummyFactor<D> dummy(Fblocks, E, P, b);
97  TIME(Overhead,dummy,xvalues,yvalues)
98 #endif
99 
100 #ifdef RAW
101  { // Raw memory Version
103  for (size_t i = 0; i < m; i++)
104  keys.push_back(i);
105  Vector x = xvalues.vector(keys);
106  double* xdata = x.data();
107 
108  // create a y
109  Vector y = Vector::Zero(m * D);
110  TIME(RawImplicit, implicitFactor, xdata, y.data())
111  TIME(RawJacobianQ, jf, xdata, y.data())
112  TIME(RawJacobianQR, jqr, xdata, y.data())
113  }
114 #endif
115 
116  os << m << endl;
117 
118 } // timeAll
119 
120 /*************************************************************************************/
121 int main(void) {
122 #ifdef SLOW
123  os << "Implicit,";
124  os << "JacobianQ,";
125  os << "JacobianQR,";
126 #endif
127 #ifdef HESSIAN
128  os << "Hessian,";
129 #endif
130 #ifdef OVERHEAD
131  os << "Overhead,";
132 #endif
133 #ifdef RAW
134  os << "RawImplicit,";
135  os << "RawJacobianQ,";
136  os << "RawJacobianQR,";
137 #endif
138  os << "m" << endl;
139  // define images
140  vector < size_t > ms;
141  // ms += 2;
142  // ms += 3, 4, 5, 6, 7, 8, 9, 10;
143  // ms += 11,12,13,14,15,16,17,18,19;
144  // ms += 20, 30, 40, 50;
145  // ms += 20,30,40,50,60,70,80,90,100;
146  for (size_t m = 2; m <= 50; m += 2)
147  ms.push_back(m);
148  //for (size_t m=10;m<=100;m+=10) ms += m;
149  // loop over number of images
150  for(size_t m: ms)
151  timeAll<PinholePose<Cal3Bundler> >(m, NUM_ITERATIONS);
152 }
153 
154 //*************************************************************************************
Matrix3f m
Pinhole camera with known calibration.
Scalar * y
Scalar * b
Definition: benchVecAdd.cpp:17
noiseModel::Diagonal::shared_ptr model
std::vector< T, typename internal::FastDefaultVectorAllocator< T >::type > FastVector
Definition: FastVector.h:34
Eigen::MatrixXd Matrix
Definition: base/Matrix.h:39
Definition: BFloat16.h:88
iterator insert(const std::pair< Key, Vector > &key_value)
static const double ms
#define N
Definition: gksort.c:12
EIGEN_DEVICE_FUNC const InverseReturnType inverse() const
int main(void)
#define NUM_ITERATIONS
Eigen::VectorXd Vector
Definition: Vector.h:38
RealScalar alpha
noiseModel::Diagonal::shared_ptr SharedDiagonal
Definition: NoiseModel.h:743
traits
Definition: chartTesting.h:28
A subclass of GaussianFactor specialized to structureless SFM.
#define TIME(label, factor, xx, yy)
DiscreteKey E(5, 2)
ofstream os("timeSchurFactors.csv")
void timeAll(size_t m, size_t N)
A Gaussian factor using the canonical parameters (information form)
Vector vector() const
const KeyVector keys
FastVector< Key > KeyVector
Define collection type once and for all - also used in wrappers.
Definition: Key.h:86
The matrix class, also used for vectors and row-vectors.
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
Calibration used by Bundler.
Timing utilities.


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:40:03