manifold.cc 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316
  1. #include "ceres/manifold.h"
  2. #include <algorithm>
  3. #include <cmath>
  4. #include "ceres/internal/eigen.h"
  5. #include "ceres/internal/fixed_array.h"
  6. #include "glog/logging.h"
  7. namespace ceres {
  8. namespace {
  9. struct CeresQuaternionOrder {
  10. static constexpr int kW = 0;
  11. static constexpr int kX = 1;
  12. static constexpr int kY = 2;
  13. static constexpr int kZ = 3;
  14. };
  15. struct EigenQuaternionOrder {
  16. static constexpr int kW = 3;
  17. static constexpr int kX = 0;
  18. static constexpr int kY = 1;
  19. static constexpr int kZ = 2;
  20. };
  21. template <typename Order>
  22. inline void QuaternionPlusImpl(const double* x,
  23. const double* delta,
  24. double* x_plus_delta) {
  25. // x_plus_delta = QuaternionProduct(q_delta, x), where q_delta is the
  26. // quaternion constructed from delta.
  27. const double norm_delta = std::hypot(delta[0], delta[1], delta[2]);
  28. if (std::fpclassify(norm_delta) == FP_ZERO) {
  29. // No change in rotation: return the quaternion as is.
  30. std::copy_n(x, 4, x_plus_delta);
  31. return;
  32. }
  33. const double sin_delta_by_delta = (std::sin(norm_delta) / norm_delta);
  34. double q_delta[4];
  35. q_delta[Order::kW] = std::cos(norm_delta);
  36. q_delta[Order::kX] = sin_delta_by_delta * delta[0];
  37. q_delta[Order::kY] = sin_delta_by_delta * delta[1];
  38. q_delta[Order::kZ] = sin_delta_by_delta * delta[2];
  39. x_plus_delta[Order::kW] =
  40. q_delta[Order::kW] * x[Order::kW] - q_delta[Order::kX] * x[Order::kX] -
  41. q_delta[Order::kY] * x[Order::kY] - q_delta[Order::kZ] * x[Order::kZ];
  42. x_plus_delta[Order::kX] =
  43. q_delta[Order::kW] * x[Order::kX] + q_delta[Order::kX] * x[Order::kW] +
  44. q_delta[Order::kY] * x[Order::kZ] - q_delta[Order::kZ] * x[Order::kY];
  45. x_plus_delta[Order::kY] =
  46. q_delta[Order::kW] * x[Order::kY] - q_delta[Order::kX] * x[Order::kZ] +
  47. q_delta[Order::kY] * x[Order::kW] + q_delta[Order::kZ] * x[Order::kX];
  48. x_plus_delta[Order::kZ] =
  49. q_delta[Order::kW] * x[Order::kZ] + q_delta[Order::kX] * x[Order::kY] -
  50. q_delta[Order::kY] * x[Order::kX] + q_delta[Order::kZ] * x[Order::kW];
  51. }
  52. template <typename Order>
  53. inline void QuaternionPlusJacobianImpl(const double* x, double* jacobian_ptr) {
  54. Eigen::Map<Eigen::Matrix<double, 4, 3, Eigen::RowMajor>> jacobian(
  55. jacobian_ptr);
  56. jacobian(Order::kW, 0) = -x[Order::kX];
  57. jacobian(Order::kW, 1) = -x[Order::kY];
  58. jacobian(Order::kW, 2) = -x[Order::kZ];
  59. jacobian(Order::kX, 0) = x[Order::kW];
  60. jacobian(Order::kX, 1) = x[Order::kZ];
  61. jacobian(Order::kX, 2) = -x[Order::kY];
  62. jacobian(Order::kY, 0) = -x[Order::kZ];
  63. jacobian(Order::kY, 1) = x[Order::kW];
  64. jacobian(Order::kY, 2) = x[Order::kX];
  65. jacobian(Order::kZ, 0) = x[Order::kY];
  66. jacobian(Order::kZ, 1) = -x[Order::kX];
  67. jacobian(Order::kZ, 2) = x[Order::kW];
  68. }
  69. template <typename Order>
  70. inline void QuaternionMinusImpl(const double* y,
  71. const double* x,
  72. double* y_minus_x) {
  73. // ambient_y_minus_x = QuaternionProduct(y, -x) where -x is the conjugate of
  74. // x.
  75. double ambient_y_minus_x[4];
  76. ambient_y_minus_x[Order::kW] =
  77. y[Order::kW] * x[Order::kW] + y[Order::kX] * x[Order::kX] +
  78. y[Order::kY] * x[Order::kY] + y[Order::kZ] * x[Order::kZ];
  79. ambient_y_minus_x[Order::kX] =
  80. -y[Order::kW] * x[Order::kX] + y[Order::kX] * x[Order::kW] -
  81. y[Order::kY] * x[Order::kZ] + y[Order::kZ] * x[Order::kY];
  82. ambient_y_minus_x[Order::kY] =
  83. -y[Order::kW] * x[Order::kY] + y[Order::kX] * x[Order::kZ] +
  84. y[Order::kY] * x[Order::kW] - y[Order::kZ] * x[Order::kX];
  85. ambient_y_minus_x[Order::kZ] =
  86. -y[Order::kW] * x[Order::kZ] - y[Order::kX] * x[Order::kY] +
  87. y[Order::kY] * x[Order::kX] + y[Order::kZ] * x[Order::kW];
  88. const double u_norm = std::hypot(ambient_y_minus_x[Order::kX],
  89. ambient_y_minus_x[Order::kY],
  90. ambient_y_minus_x[Order::kZ]);
  91. if (std::fpclassify(u_norm) != FP_ZERO) {
  92. const double theta = std::atan2(u_norm, ambient_y_minus_x[Order::kW]);
  93. y_minus_x[0] = theta * ambient_y_minus_x[Order::kX] / u_norm;
  94. y_minus_x[1] = theta * ambient_y_minus_x[Order::kY] / u_norm;
  95. y_minus_x[2] = theta * ambient_y_minus_x[Order::kZ] / u_norm;
  96. } else {
  97. std::fill_n(y_minus_x, 3, 0.0);
  98. }
  99. }
  100. template <typename Order>
  101. inline void QuaternionMinusJacobianImpl(const double* x, double* jacobian_ptr) {
  102. Eigen::Map<Eigen::Matrix<double, 3, 4, Eigen::RowMajor>> jacobian(
  103. jacobian_ptr);
  104. jacobian(0, Order::kW) = -x[Order::kX];
  105. jacobian(0, Order::kX) = x[Order::kW];
  106. jacobian(0, Order::kY) = -x[Order::kZ];
  107. jacobian(0, Order::kZ) = x[Order::kY];
  108. jacobian(1, Order::kW) = -x[Order::kY];
  109. jacobian(1, Order::kX) = x[Order::kZ];
  110. jacobian(1, Order::kY) = x[Order::kW];
  111. jacobian(1, Order::kZ) = -x[Order::kX];
  112. jacobian(2, Order::kW) = -x[Order::kZ];
  113. jacobian(2, Order::kX) = -x[Order::kY];
  114. jacobian(2, Order::kY) = x[Order::kX];
  115. jacobian(2, Order::kZ) = x[Order::kW];
  116. }
  117. } // namespace
  118. Manifold::~Manifold() = default;
  119. bool Manifold::RightMultiplyByPlusJacobian(const double* x,
  120. const int num_rows,
  121. const double* ambient_matrix,
  122. double* tangent_matrix) const {
  123. const int tangent_size = TangentSize();
  124. if (tangent_size == 0) {
  125. return true;
  126. }
  127. const int ambient_size = AmbientSize();
  128. Matrix plus_jacobian(ambient_size, tangent_size);
  129. if (!PlusJacobian(x, plus_jacobian.data())) {
  130. return false;
  131. }
  132. MatrixRef(tangent_matrix, num_rows, tangent_size) =
  133. ConstMatrixRef(ambient_matrix, num_rows, ambient_size) * plus_jacobian;
  134. return true;
  135. }
  136. SubsetManifold::SubsetManifold(const int size,
  137. const std::vector<int>& constant_parameters)
  138. : tangent_size_(size - constant_parameters.size()),
  139. constancy_mask_(size, false) {
  140. if (constant_parameters.empty()) {
  141. return;
  142. }
  143. std::vector<int> constant = constant_parameters;
  144. std::sort(constant.begin(), constant.end());
  145. CHECK_GE(constant.front(), 0) << "Indices indicating constant parameter must "
  146. "be greater than equal to zero.";
  147. CHECK_LT(constant.back(), size)
  148. << "Indices indicating constant parameter must be less than the size "
  149. << "of the parameter block.";
  150. CHECK(std::adjacent_find(constant.begin(), constant.end()) == constant.end())
  151. << "The set of constant parameters cannot contain duplicates";
  152. for (auto index : constant_parameters) {
  153. constancy_mask_[index] = true;
  154. }
  155. }
  156. int SubsetManifold::AmbientSize() const { return constancy_mask_.size(); }
  157. int SubsetManifold::TangentSize() const { return tangent_size_; }
  158. bool SubsetManifold::Plus(const double* x,
  159. const double* delta,
  160. double* x_plus_delta) const {
  161. const int ambient_size = AmbientSize();
  162. for (int i = 0, j = 0; i < ambient_size; ++i) {
  163. if (constancy_mask_[i]) {
  164. x_plus_delta[i] = x[i];
  165. } else {
  166. x_plus_delta[i] = x[i] + delta[j++];
  167. }
  168. }
  169. return true;
  170. }
  171. bool SubsetManifold::PlusJacobian(const double* /*x*/,
  172. double* plus_jacobian) const {
  173. if (tangent_size_ == 0) {
  174. return true;
  175. }
  176. const int ambient_size = AmbientSize();
  177. MatrixRef m(plus_jacobian, ambient_size, tangent_size_);
  178. m.setZero();
  179. for (int r = 0, c = 0; r < ambient_size; ++r) {
  180. if (!constancy_mask_[r]) {
  181. m(r, c++) = 1.0;
  182. }
  183. }
  184. return true;
  185. }
  186. bool SubsetManifold::RightMultiplyByPlusJacobian(const double* /*x*/,
  187. const int num_rows,
  188. const double* ambient_matrix,
  189. double* tangent_matrix) const {
  190. if (tangent_size_ == 0) {
  191. return true;
  192. }
  193. const int ambient_size = AmbientSize();
  194. for (int r = 0; r < num_rows; ++r) {
  195. for (int idx = 0, c = 0; idx < ambient_size; ++idx) {
  196. if (!constancy_mask_[idx]) {
  197. tangent_matrix[r * tangent_size_ + c++] =
  198. ambient_matrix[r * ambient_size + idx];
  199. }
  200. }
  201. }
  202. return true;
  203. }
  204. bool SubsetManifold::Minus(const double* y,
  205. const double* x,
  206. double* y_minus_x) const {
  207. if (tangent_size_ == 0) {
  208. return true;
  209. }
  210. const int ambient_size = AmbientSize();
  211. for (int i = 0, j = 0; i < ambient_size; ++i) {
  212. if (!constancy_mask_[i]) {
  213. y_minus_x[j++] = y[i] - x[i];
  214. }
  215. }
  216. return true;
  217. }
  218. bool SubsetManifold::MinusJacobian(const double* /*x*/,
  219. double* minus_jacobian) const {
  220. const int ambient_size = AmbientSize();
  221. MatrixRef m(minus_jacobian, tangent_size_, ambient_size);
  222. m.setZero();
  223. for (int c = 0, r = 0; c < ambient_size; ++c) {
  224. if (!constancy_mask_[c]) {
  225. m(r++, c) = 1.0;
  226. }
  227. }
  228. return true;
  229. }
  230. bool QuaternionManifold::Plus(const double* x,
  231. const double* delta,
  232. double* x_plus_delta) const {
  233. QuaternionPlusImpl<CeresQuaternionOrder>(x, delta, x_plus_delta);
  234. return true;
  235. }
  236. bool QuaternionManifold::PlusJacobian(const double* x, double* jacobian) const {
  237. QuaternionPlusJacobianImpl<CeresQuaternionOrder>(x, jacobian);
  238. return true;
  239. }
  240. bool QuaternionManifold::Minus(const double* y,
  241. const double* x,
  242. double* y_minus_x) const {
  243. QuaternionMinusImpl<CeresQuaternionOrder>(y, x, y_minus_x);
  244. return true;
  245. }
  246. bool QuaternionManifold::MinusJacobian(const double* x,
  247. double* jacobian) const {
  248. QuaternionMinusJacobianImpl<CeresQuaternionOrder>(x, jacobian);
  249. return true;
  250. }
  251. bool EigenQuaternionManifold::Plus(const double* x,
  252. const double* delta,
  253. double* x_plus_delta) const {
  254. QuaternionPlusImpl<EigenQuaternionOrder>(x, delta, x_plus_delta);
  255. return true;
  256. }
  257. bool EigenQuaternionManifold::PlusJacobian(const double* x,
  258. double* jacobian) const {
  259. QuaternionPlusJacobianImpl<EigenQuaternionOrder>(x, jacobian);
  260. return true;
  261. }
  262. bool EigenQuaternionManifold::Minus(const double* y,
  263. const double* x,
  264. double* y_minus_x) const {
  265. QuaternionMinusImpl<EigenQuaternionOrder>(y, x, y_minus_x);
  266. return true;
  267. }
  268. bool EigenQuaternionManifold::MinusJacobian(const double* x,
  269. double* jacobian) const {
  270. QuaternionMinusJacobianImpl<EigenQuaternionOrder>(x, jacobian);
  271. return true;
  272. }
  273. } // namespace ceres