fast_winding_number.cpp 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371
  1. #include "fast_winding_number.h"
  2. #include "build_octree.h"
  3. #include "knn_octree.h"
  4. #include "parallel_for.h"
  5. #include "copyleft/cgal/point_areas_and_normals.h"
  6. #include <vector>
  7. namespace igl {
  8. template <typename DerivedP, typename DerivedA, typename DerivedN,
  9. typename Index, typename DerivedCM, typename DerivedR, typename DerivedEC>
  10. IGL_INLINE void fast_winding_number(const Eigen::MatrixBase<DerivedP>& P,
  11. const Eigen::MatrixBase<DerivedN>& N,
  12. const Eigen::MatrixBase<DerivedA>& A,
  13. const std::vector<std::vector<Index> > & point_indices,
  14. const std::vector<Eigen::Matrix<Index,8,1>,
  15. Eigen::aligned_allocator<Eigen::Matrix<Index,8,1> > > & children,
  16. const int expansion_order,
  17. Eigen::PlainObjectBase<DerivedCM>& CM,
  18. Eigen::PlainObjectBase<DerivedR>& R,
  19. Eigen::PlainObjectBase<DerivedEC>& EC
  20. ){
  21. typedef typename DerivedP::Scalar real_p;
  22. typedef typename DerivedN::Scalar real_n;
  23. typedef typename DerivedA::Scalar real_a;
  24. typedef typename DerivedCM::Scalar real_cm;
  25. typedef typename DerivedR::Scalar real_r;
  26. typedef typename DerivedEC::Scalar real_ec;
  27. typedef Eigen::Matrix<real_p,1,3> RowVec3p;
  28. int m = children.size();
  29. int num_terms;
  30. if(expansion_order == 0){
  31. num_terms = 3;
  32. } else if(expansion_order ==1){
  33. num_terms = 3 + 9;
  34. } else if(expansion_order == 2){
  35. num_terms = 3 + 9 + 27;
  36. } else {
  37. assert(false);
  38. }
  39. R.resize(m);
  40. CM.resize(m,3);
  41. EC.resize(m,num_terms);
  42. EC.setZero(m,num_terms);
  43. std::function< void(const int) > helper;
  44. helper = [&helper,
  45. &P,&N,&A,&expansion_order,&point_indices,&children,&EC,&R,&CM]
  46. (const int index)-> void
  47. {
  48. Eigen::Matrix<real_cm,1,3> masscenter;
  49. masscenter << 0,0,0;
  50. Eigen::Matrix<real_ec,1,3> zeroth_expansion;
  51. masscenter << 0,0,0;
  52. real_p areatotal = 0.0;
  53. for(int j = 0; j < point_indices.at(index).size(); j++){
  54. int curr_point_index = point_indices.at(index).at(j);
  55. areatotal += A(curr_point_index);
  56. masscenter += A(curr_point_index)*P.row(curr_point_index);
  57. zeroth_expansion += A(curr_point_index)*N.row(curr_point_index);
  58. }
  59. masscenter = masscenter/areatotal;
  60. CM.row(index) = masscenter;
  61. EC.block(index,0,1,3) = zeroth_expansion;
  62. real_r max_norm = 0;
  63. real_r curr_norm;
  64. for(int i = 0; i < point_indices.at(index).size(); i++){
  65. //Get max distance from center of mass:
  66. int curr_point_index = point_indices.at(index).at(i);
  67. Eigen::Matrix<real_r,1,3> point =
  68. P.row(curr_point_index)-masscenter;
  69. curr_norm = point.norm();
  70. if(curr_norm > max_norm){
  71. max_norm = curr_norm;
  72. }
  73. //Calculate higher order terms if necessary
  74. Eigen::Matrix<real_ec,3,3> TempCoeffs;
  75. if(EC.cols() >= (3+9)){
  76. TempCoeffs = A(curr_point_index)*point.transpose()*
  77. N.row(curr_point_index);
  78. EC.block(index,3,1,9) +=
  79. Eigen::Map<Eigen::Matrix<real_ec,1,9> >(TempCoeffs.data(),
  80. TempCoeffs.size());
  81. }
  82. if(EC.cols() == (3+9+27)){
  83. for(int k = 0; k < 3; k++){
  84. TempCoeffs = 0.5 * point(k) * (A(curr_point_index)*
  85. point.transpose()*N.row(curr_point_index));
  86. EC.block(index,12+9*k,1,9) += Eigen::Map<
  87. Eigen::Matrix<real_ec,1,9> >(TempCoeffs.data(),
  88. TempCoeffs.size());
  89. }
  90. }
  91. }
  92. R(index) = max_norm;
  93. if(children.at(index)(0) != -1)
  94. {
  95. for(int i = 0; i < 8; i++){
  96. int child = children.at(index)(i);
  97. helper(child);
  98. }
  99. }
  100. };
  101. helper(0);
  102. }
  103. template <typename DerivedP, typename DerivedA, typename DerivedN,
  104. typename Index, typename DerivedCM, typename DerivedR, typename DerivedEC,
  105. typename DerivedQ, typename BetaType, typename DerivedWN>
  106. IGL_INLINE void fast_winding_number(const Eigen::MatrixBase<DerivedP>& P,
  107. const Eigen::MatrixBase<DerivedN>& N,
  108. const Eigen::MatrixBase<DerivedA>& A,
  109. const std::vector<std::vector<Index> > & point_indices,
  110. const std::vector<Eigen::Matrix<Index,8,1>,
  111. Eigen::aligned_allocator<Eigen::Matrix<Index,8,1> > > & children,
  112. const Eigen::MatrixBase<DerivedCM>& CM,
  113. const Eigen::MatrixBase<DerivedR>& R,
  114. const Eigen::MatrixBase<DerivedEC>& EC,
  115. const Eigen::MatrixBase<DerivedQ>& Q,
  116. const BetaType beta,
  117. Eigen::PlainObjectBase<DerivedWN>& WN
  118. ){
  119. typedef typename DerivedP::Scalar real_p;
  120. typedef typename DerivedN::Scalar real_n;
  121. typedef typename DerivedA::Scalar real_a;
  122. typedef typename DerivedCM::Scalar real_cm;
  123. typedef typename DerivedR::Scalar real_r;
  124. typedef typename DerivedEC::Scalar real_ec;
  125. typedef typename DerivedQ::Scalar real_q;
  126. typedef typename DerivedWN::Scalar real_wn;
  127. typedef Eigen::Matrix<real_q,1,3> RowVec;
  128. typedef Eigen::Matrix<real_ec,3,3> EC_3by3;
  129. auto direct_eval = [](const RowVec & loc,
  130. const Eigen::Matrix<real_ec,1,3> & anorm){
  131. real_wn wn = (loc(0)*anorm(0)+loc(1)*anorm(1)+loc(2)*anorm(2))
  132. /(4.0*M_PI*std::pow(loc.norm(),3));
  133. if(std::isnan(wn)){
  134. return 0.5;
  135. }else{
  136. return wn;
  137. }
  138. };
  139. auto expansion_eval = [&direct_eval](const RowVec & loc,
  140. const Eigen::RowVectorXd & EC){
  141. real_wn wn = direct_eval(loc,EC.head<3>());
  142. double r = loc.norm();
  143. if(EC.size()>3){
  144. Eigen::Matrix<real_ec,3,3> SecondDerivative =
  145. Eigen::Matrix<real_ec,3,3>::Identity()/(4.0*M_PI*std::pow(r,3));
  146. SecondDerivative += -3.0*loc.transpose()*loc/(4.0*M_PI*std::pow(r,5));
  147. Eigen::Matrix<real_ec,1,9> derivative_vector =
  148. Eigen::Map<Eigen::Matrix<real_ec,1,9> >(SecondDerivative.data(),
  149. SecondDerivative.size());
  150. wn += derivative_vector.cwiseProduct(EC.segment<9>(3)).sum();
  151. }
  152. if(EC.size()>3+9){
  153. Eigen::Matrix<real_ec,3,3> ThirdDerivative;
  154. for(int i = 0; i < 3; i++){
  155. ThirdDerivative =
  156. 15.0*loc(i)*loc.transpose()*loc/(4.0*M_PI*std::pow(r,7));
  157. Eigen::Matrix<real_ec,3,3> Diagonal;
  158. Diagonal << loc(i), 0, 0,
  159. 0, loc(i), 0,
  160. 0, 0, loc(i);
  161. Eigen::Matrix<real_ec,3,3> RowCol;
  162. RowCol.setZero(3,3);
  163. RowCol.row(i) = loc;
  164. RowCol = RowCol + RowCol.transpose();
  165. ThirdDerivative +=
  166. -3.0/(4.0*M_PI*std::pow(r,5))*(RowCol+Diagonal);
  167. Eigen::Matrix<real_ec,1,9> derivative_vector =
  168. Eigen::Map<Eigen::Matrix<real_ec,1,9> >(ThirdDerivative.data(),
  169. ThirdDerivative.size());
  170. wn += derivative_vector.cwiseProduct(
  171. EC.segment<9>(12 + i*9)).sum();
  172. }
  173. }
  174. return wn;
  175. };
  176. int m = Q.rows();
  177. WN.resize(m,1);
  178. std::function< real_wn(const RowVec, const std::vector<int>) > helper;
  179. helper = [&helper,
  180. &P,&N,&A,
  181. &point_indices,&children,
  182. &CM,&R,&EC,&beta,
  183. &direct_eval,&expansion_eval]
  184. (const RowVec query, const std::vector<int> near_indices)-> real_wn
  185. {
  186. std::vector<int> new_near_indices;
  187. real_wn wn = 0;
  188. for(int i = 0; i < near_indices.size(); i++){
  189. int index = near_indices.at(i);
  190. //Leaf Case, Brute force
  191. if(children.at(index)(0) == -1){
  192. for(int j = 0; j < point_indices.at(index).size(); j++){
  193. int curr_row = point_indices.at(index).at(j);
  194. wn += direct_eval(P.row(curr_row)-query,
  195. N.row(curr_row)*A(curr_row));
  196. }
  197. }
  198. //Non-Leaf Case
  199. else {
  200. for(int child = 0; child < 8; child++){
  201. int child_index = children.at(index)(child);
  202. if(point_indices.at(child_index).size() > 0){
  203. if((CM.row(child_index)-query).norm() > beta*R(child_index)){
  204. if(children.at(child_index)(0) == -1){
  205. for(int j=0;j<point_indices.at(child_index).size();j++){
  206. int curr_row = point_indices.at(child_index).at(j);
  207. wn += direct_eval(P.row(curr_row)-query,
  208. N.row(curr_row)*A(curr_row));
  209. }
  210. }else{
  211. wn += expansion_eval(CM.row(child_index)-query,
  212. EC.row(child_index));
  213. }
  214. }else {
  215. new_near_indices.emplace_back(child_index);
  216. }
  217. }
  218. }
  219. }
  220. }
  221. if(new_near_indices.size() > 0){
  222. wn += helper(query,new_near_indices);
  223. }
  224. return wn;
  225. };
  226. if(beta > 0){
  227. std::vector<int> near_indices_start = {0};
  228. igl::parallel_for(m,[&](int iter){
  229. WN(iter) = helper(Q.row(iter),near_indices_start);
  230. },1000);
  231. } else {
  232. igl::parallel_for(m,[&](int iter){
  233. double wn = 0;
  234. for(int j = 0; j <P.rows(); j++){
  235. wn += direct_eval(P.row(j)-Q.row(iter),N.row(j)*A(j));
  236. }
  237. WN(iter) = wn;
  238. },1000);
  239. }
  240. }
  241. template <typename DerivedP, typename DerivedA, typename DerivedN,
  242. typename DerivedQ, typename BetaType, typename DerivedWN>
  243. IGL_INLINE void fast_winding_number(const Eigen::MatrixBase<DerivedP>& P,
  244. const Eigen::MatrixBase<DerivedN>& N,
  245. const Eigen::MatrixBase<DerivedA>& A,
  246. const Eigen::MatrixBase<DerivedQ>& Q,
  247. const int expansion_order,
  248. const BetaType beta,
  249. Eigen::PlainObjectBase<DerivedWN>& WN
  250. ){
  251. typedef typename DerivedWN::Scalar real;
  252. typedef typename Eigen::Matrix<real,1,3> RowVec;
  253. std::vector<std::vector<int> > point_indices;
  254. std::vector<Eigen::Matrix<int,8,1>,
  255. Eigen::aligned_allocator<Eigen::Matrix<int,8,1> > > children;
  256. std::vector<RowVec, Eigen::aligned_allocator<RowVec> > centers;
  257. std::vector<real> widths;
  258. Eigen::MatrixXi I;
  259. Eigen::MatrixXd NotUsed;
  260. build_octree(P,point_indices,children,centers,widths);
  261. Eigen::MatrixXd EC, CM;
  262. Eigen::VectorXd R;
  263. fast_winding_number(P,N,A,point_indices,children,
  264. expansion_order,CM,R,EC);
  265. fast_winding_number(P,N,A,point_indices,children,CM,R,EC,Q,beta,WN);
  266. }
  267. template <typename DerivedP, typename DerivedA, typename DerivedN,
  268. typename DerivedQ, typename DerivedWN>
  269. IGL_INLINE void fast_winding_number(const Eigen::MatrixBase<DerivedP>& P,
  270. const Eigen::MatrixBase<DerivedN>& N,
  271. const Eigen::MatrixBase<DerivedA>& A,
  272. const Eigen::MatrixBase<DerivedQ>& Q,
  273. Eigen::PlainObjectBase<DerivedWN>& WN
  274. ){
  275. fast_winding_number(P,N,A,Q,2,2.0,WN);
  276. }
  277. template <typename DerivedP, typename DerivedN, typename DerivedQ,
  278. typename BetaType, typename DerivedWN>
  279. IGL_INLINE void fast_winding_number(const Eigen::MatrixBase<DerivedP>& P,
  280. const Eigen::MatrixBase<DerivedN>& N,
  281. const Eigen::MatrixBase<DerivedQ>& Q,
  282. const int expansion_order,
  283. const BetaType beta,
  284. Eigen::PlainObjectBase<DerivedWN>& WN
  285. ){
  286. typedef typename DerivedWN::Scalar real;
  287. typedef typename Eigen::Matrix<real,1,3> RowVec;
  288. std::vector<std::vector<int> > point_indices;
  289. std::vector<Eigen::Matrix<int,8,1>,
  290. Eigen::aligned_allocator<Eigen::Matrix<int,8,1> > > children;
  291. std::vector<RowVec, Eigen::aligned_allocator<RowVec> > centers;
  292. std::vector<real> widths;
  293. Eigen::MatrixXi I;
  294. Eigen::Matrix<real,Eigen::Dynamic,Eigen::Dynamic> NotUsed;
  295. Eigen::Matrix<real,Eigen::Dynamic,1> A;
  296. build_octree(P,point_indices,children,centers,widths);
  297. knn_octree(P,21,point_indices,children,centers,widths,I);
  298. copyleft::cgal::point_areas_and_normals(P,I,N,A,NotUsed);
  299. Eigen::MatrixXd EC, CM;
  300. Eigen::VectorXd R;
  301. fast_winding_number(P,N,A,point_indices,children,expansion_order,CM,R,EC);
  302. fast_winding_number(P,N,A,point_indices,children,CM,R,EC,Q,beta,WN);
  303. }
  304. template <typename DerivedP, typename DerivedN,
  305. typename DerivedQ, typename DerivedWN>
  306. IGL_INLINE void fast_winding_number(
  307. const Eigen::MatrixBase<DerivedP>& P,
  308. const Eigen::MatrixBase<DerivedN>& N,
  309. const Eigen::MatrixBase<DerivedQ>& Q,
  310. Eigen::PlainObjectBase<DerivedWN>& WN
  311. ){
  312. fast_winding_number(P,N,Q,2,2.0,WN);
  313. }
  314. }