timeSchurFactors.cpp 3.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157
  1. /**
  2. * @file timeSchurFactors.cpp
  3. * @brief Time various Schur-complement Jacobian factors
  4. * @author Frank Dellaert
  5. * @date Oct 27, 2013
  6. */
  7. #include "DummyFactor.h"
  8. #include <gtsam/base/timing.h>
  9. #include <gtsam/slam/JacobianFactorQ.h>
  10. #include "gtsam/slam/JacobianFactorQR.h"
  11. #include <gtsam/slam/RegularImplicitSchurFactor.h>
  12. #include <gtsam/geometry/Cal3Bundler.h>
  13. #include <gtsam/geometry/PinholePose.h>
  14. #include <boost/assign/list_of.hpp>
  15. #include <boost/assign/std/vector.hpp>
  16. #include <fstream>
  17. using namespace std;
  18. using namespace boost::assign;
  19. using namespace gtsam;
  20. #define SLOW
  21. #define RAW
  22. #define HESSIAN
  23. #define NUM_ITERATIONS 1000
  24. // Create CSV file for results
  25. ofstream os("timeSchurFactors.csv");
  26. /*************************************************************************************/
  27. template<typename CAMERA>
  28. void timeAll(size_t m, size_t N) {
  29. cout << m << endl;
  30. // create F
  31. static const int D = CAMERA::dimension;
  32. typedef Eigen::Matrix<double, 2, D> Matrix2D;
  33. KeyVector keys;
  34. vector <Matrix2D, Eigen::aligned_allocator<Matrix2D>> Fblocks;
  35. for (size_t i = 0; i < m; i++) {
  36. keys.push_back(i);
  37. Fblocks.push_back((i + 1) * Matrix::Ones(2, D));
  38. }
  39. // create E
  40. Matrix E(2 * m, 3);
  41. for (size_t i = 0; i < m; i++)
  42. E.block < 2, 3 > (2 * i, 0) = Matrix::Ones(2, 3);
  43. // Calculate point covariance
  44. Matrix P = (E.transpose() * E).inverse();
  45. // RHS and sigmas
  46. const Vector b = Vector::Constant(2*m,1);
  47. const SharedDiagonal model;
  48. // parameters for multiplyHessianAdd
  49. double alpha = 0.5;
  50. VectorValues xvalues, yvalues;
  51. for (size_t i = 0; i < m; i++)
  52. xvalues.insert(i, Vector::Constant(D,2));
  53. // Implicit
  54. RegularImplicitSchurFactor<CAMERA> implicitFactor(keys, Fblocks, E, P, b);
  55. // JacobianFactor with same error
  56. JacobianFactorQ<D, 2> jf(keys, Fblocks, E, P, b, model);
  57. // JacobianFactorQR with same error
  58. JacobianFactorQR<D, 2> jqr(keys, Fblocks, E, P, b, model);
  59. // Hessian
  60. HessianFactor hessianFactor(jqr);
  61. #define TIME(label,factor,xx,yy) {\
  62. for (size_t t = 0; t < N; t++) \
  63. factor.multiplyHessianAdd(alpha, xx, yy);\
  64. gttic_(label);\
  65. for (size_t t = 0; t < N; t++) {\
  66. factor.multiplyHessianAdd(alpha, xx, yy);\
  67. }\
  68. gttoc_(label);\
  69. tictoc_getNode(timer, label)\
  70. os << timer->secs()/NUM_ITERATIONS << ", ";\
  71. }
  72. #ifdef SLOW
  73. TIME(Implicit, implicitFactor, xvalues, yvalues)
  74. TIME(Jacobian, jf, xvalues, yvalues)
  75. TIME(JacobianQR, jqr, xvalues, yvalues)
  76. #endif
  77. #ifdef HESSIAN
  78. TIME(Hessian, hessianFactor, xvalues, yvalues)
  79. #endif
  80. #ifdef OVERHEAD
  81. DummyFactor<D> dummy(Fblocks, E, P, b);
  82. TIME(Overhead,dummy,xvalues,yvalues)
  83. #endif
  84. #ifdef RAW
  85. { // Raw memory Version
  86. FastVector < Key > keys;
  87. for (size_t i = 0; i < m; i++)
  88. keys += i;
  89. Vector x = xvalues.vector(keys);
  90. double* xdata = x.data();
  91. // create a y
  92. Vector y = Vector::Zero(m * D);
  93. TIME(RawImplicit, implicitFactor, xdata, y.data())
  94. TIME(RawJacobianQ, jf, xdata, y.data())
  95. TIME(RawJacobianQR, jqr, xdata, y.data())
  96. }
  97. #endif
  98. os << m << endl;
  99. } // timeAll
  100. /*************************************************************************************/
  101. int main(void) {
  102. #ifdef SLOW
  103. os << "Implicit,";
  104. os << "JacobianQ,";
  105. os << "JacobianQR,";
  106. #endif
  107. #ifdef HESSIAN
  108. os << "Hessian,";
  109. #endif
  110. #ifdef OVERHEAD
  111. os << "Overhead,";
  112. #endif
  113. #ifdef RAW
  114. os << "RawImplicit,";
  115. os << "RawJacobianQ,";
  116. os << "RawJacobianQR,";
  117. #endif
  118. os << "m" << endl;
  119. // define images
  120. vector < size_t > ms;
  121. // ms += 2;
  122. // ms += 3, 4, 5, 6, 7, 8, 9, 10;
  123. // ms += 11,12,13,14,15,16,17,18,19;
  124. // ms += 20, 30, 40, 50;
  125. // ms += 20,30,40,50,60,70,80,90,100;
  126. for (size_t m = 2; m <= 50; m += 2)
  127. ms += m;
  128. //for (size_t m=10;m<=100;m+=10) ms += m;
  129. // loop over number of images
  130. for(size_t m: ms)
  131. timeAll<PinholePose<Cal3Bundler> >(m, NUM_ITERATIONS);
  132. }
  133. //*************************************************************************************