fit_rotations.cpp 6.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <alecjacobson@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "fit_rotations.h"
  9. #include "polar_svd3x3.h"
  10. #include <igl/repmat.h>
  11. #include <igl/verbose.h>
  12. #include <igl/polar_dec.h>
  13. #include <igl/polar_svd.h>
  14. #include <igl/matlab_format.h>
  15. #include <igl/C_STR.h>
  16. #include <iostream>
  17. template <typename DerivedS, typename DerivedD>
  18. IGL_INLINE void igl::fit_rotations(
  19. const Eigen::PlainObjectBase<DerivedS> & S,
  20. const bool single_precision,
  21. Eigen::PlainObjectBase<DerivedD> & R)
  22. {
  23. using namespace std;
  24. const int dim = S.cols();
  25. const int nr = S.rows()/dim;
  26. assert(nr * dim == S.rows());
  27. assert(dim == 3);
  28. // resize output
  29. R.resize(dim,dim*nr); // hopefully no op (should be already allocated)
  30. //std::cout<<"S=["<<std::endl<<S<<std::endl<<"];"<<std::endl;
  31. //MatrixXd si(dim,dim);
  32. Eigen::Matrix<typename DerivedS::Scalar,3,3> si;// = Eigen::Matrix3d::Identity();
  33. // loop over number of rotations we're computing
  34. for(int r = 0;r<nr;r++)
  35. {
  36. // build this covariance matrix
  37. for(int i = 0;i<dim;i++)
  38. {
  39. for(int j = 0;j<dim;j++)
  40. {
  41. si(i,j) = S(i*nr+r,j);
  42. }
  43. }
  44. typedef Eigen::Matrix<typename DerivedD::Scalar,3,3> Mat3;
  45. typedef Eigen::Matrix<typename DerivedD::Scalar,3,1> Vec3;
  46. Mat3 ri;
  47. if(single_precision)
  48. {
  49. polar_svd3x3(si, ri);
  50. }else
  51. {
  52. Mat3 ti,ui,vi;
  53. Vec3 _;
  54. igl::polar_svd(si,ri,ti,ui,_,vi);
  55. // Check for reflection
  56. if(ri.determinant() < 0)
  57. {
  58. vi.col(1) *= -1.;
  59. ri = ui * vi.transpose();
  60. }
  61. }
  62. assert(ri.determinant() >= 0);
  63. // Not sure why polar_dec computes transpose...
  64. R.block(0,r*dim,dim,dim) = ri.block(0,0,dim,dim).transpose();
  65. //cout<<matlab_format(si,C_STR("si_"<<r))<<endl;
  66. //cout<<matlab_format(ri.transpose().eval(),C_STR("ri_"<<r))<<endl;
  67. }
  68. }
  69. template <typename DerivedS, typename DerivedD>
  70. IGL_INLINE void igl::fit_rotations_planar(
  71. const Eigen::PlainObjectBase<DerivedS> & S,
  72. Eigen::PlainObjectBase<DerivedD> & R)
  73. {
  74. using namespace std;
  75. const int dim = S.cols();
  76. const int nr = S.rows()/dim;
  77. assert(dim == 2 && "_planar input should be 2D");
  78. assert(nr * dim == S.rows());
  79. // resize output
  80. R.resize(dim,dim*nr); // hopefully no op (should be already allocated)
  81. Eigen::Matrix<typename DerivedS::Scalar,2,2> si;
  82. // loop over number of rotations we're computing
  83. for(int r = 0;r<nr;r++)
  84. {
  85. // build this covariance matrix
  86. for(int i = 0;i<2;i++)
  87. {
  88. for(int j = 0;j<2;j++)
  89. {
  90. si(i,j) = S(i*nr+r,j);
  91. }
  92. }
  93. typedef Eigen::Matrix<typename DerivedD::Scalar,2,2> Mat2;
  94. typedef Eigen::Matrix<typename DerivedD::Scalar,2,1> Vec2;
  95. Mat2 ri,ti,ui,vi;
  96. Vec2 _;
  97. igl::polar_svd(si,ri,ti,ui,_,vi);
  98. #ifndef FIT_ROTATIONS_ALLOW_FLIPS
  99. // Check for reflection
  100. if(ri.determinant() < 0)
  101. {
  102. vi.col(1) *= -1.;
  103. ri = ui * vi.transpose();
  104. }
  105. assert(ri.determinant() >= 0);
  106. #endif
  107. // Not sure why polar_dec computes transpose...
  108. R.block(0,r*dim,2,2) = ri.transpose();
  109. }
  110. }
  111. #ifdef __SSE__
  112. IGL_INLINE void igl::fit_rotations_SSE(
  113. const Eigen::MatrixXf & S,
  114. Eigen::MatrixXf & R)
  115. {
  116. const int cStep = 4;
  117. assert(S.cols() == 3);
  118. const int dim = 3; //S.cols();
  119. const int nr = S.rows()/dim;
  120. assert(nr * dim == S.rows());
  121. // resize output
  122. R.resize(dim,dim*nr); // hopefully no op (should be already allocated)
  123. Eigen::Matrix<float, 3*cStep, 3> siBig;
  124. // using SSE decompose cStep matrices at a time:
  125. int r = 0;
  126. for( ; r<nr; r+=cStep)
  127. {
  128. int numMats = cStep;
  129. if (r + cStep >= nr) numMats = nr - r;
  130. // build siBig:
  131. for (int k=0; k<numMats; k++)
  132. {
  133. for(int i = 0;i<dim;i++)
  134. {
  135. for(int j = 0;j<dim;j++)
  136. {
  137. siBig(i + 3*k, j) = S(i*nr + r + k, j);
  138. }
  139. }
  140. }
  141. Eigen::Matrix<float, 3*cStep, 3> ri;
  142. polar_svd3x3_sse(siBig, ri);
  143. for (int k=0; k<cStep; k++)
  144. assert(ri.block(3*k, 0, 3, 3).determinant() >= 0);
  145. // Not sure why polar_dec computes transpose...
  146. for (int k=0; k<numMats; k++)
  147. {
  148. R.block(0, (r + k)*dim, dim, dim) = ri.block(3*k, 0, dim, dim).transpose();
  149. }
  150. }
  151. }
  152. IGL_INLINE void igl::fit_rotations_SSE(
  153. const Eigen::MatrixXd & S,
  154. Eigen::MatrixXd & R)
  155. {
  156. const Eigen::MatrixXf Sf = S.cast<float>();
  157. Eigen::MatrixXf Rf;
  158. fit_rotations_SSE(Sf,Rf);
  159. R = Rf.cast<double>();
  160. }
  161. #endif
  162. #ifdef __AVX__
  163. IGL_INLINE void igl::fit_rotations_AVX(
  164. const Eigen::MatrixXf & S,
  165. Eigen::MatrixXf & R)
  166. {
  167. const int cStep = 8;
  168. assert(S.cols() == 3);
  169. const int dim = 3; //S.cols();
  170. const int nr = S.rows()/dim;
  171. assert(nr * dim == S.rows());
  172. // resize output
  173. R.resize(dim,dim*nr); // hopefully no op (should be already allocated)
  174. Eigen::Matrix<float, 3*cStep, 3> siBig;
  175. // using SSE decompose cStep matrices at a time:
  176. int r = 0;
  177. for( ; r<nr; r+=cStep)
  178. {
  179. int numMats = cStep;
  180. if (r + cStep >= nr) numMats = nr - r;
  181. // build siBig:
  182. for (int k=0; k<numMats; k++)
  183. {
  184. for(int i = 0;i<dim;i++)
  185. {
  186. for(int j = 0;j<dim;j++)
  187. {
  188. siBig(i + 3*k, j) = S(i*nr + r + k, j);
  189. }
  190. }
  191. }
  192. Eigen::Matrix<float, 3*cStep, 3> ri;
  193. polar_svd3x3_avx(siBig, ri);
  194. for (int k=0; k<cStep; k++)
  195. assert(ri.block(3*k, 0, 3, 3).determinant() >= 0);
  196. // Not sure why polar_dec computes transpose...
  197. for (int k=0; k<numMats; k++)
  198. {
  199. R.block(0, (r + k)*dim, dim, dim) = ri.block(3*k, 0, dim, dim).transpose();
  200. }
  201. }
  202. }
  203. #endif
  204. #ifndef IGL_HEADER_ONLY
  205. // Explicit template instanciation
  206. template void igl::fit_rotations<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  207. template void igl::fit_rotations_planar<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  208. template void igl::fit_rotations_planar<Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, -1, 0, -1, -1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&);
  209. #endif