WindingNumberAABB.h 9.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 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. // # MUTUAL DEPENDENCY ISSUE FOR HEADER ONLY VERSION
  9. // MUST INCLUDE winding_number.h first before guard:
  10. #include "winding_number.h"
  11. #ifndef IGL_WINDINGNUMBERAABB_H
  12. #define IGL_WINDINGNUMBERAABB_H
  13. #include "WindingNumberTree.h"
  14. namespace igl
  15. {
  16. template <typename Point>
  17. class WindingNumberAABB : public WindingNumberTree<Point>
  18. {
  19. protected:
  20. Point min_corner;
  21. Point max_corner;
  22. double total_positive_area;
  23. public:
  24. enum SplitMethod
  25. {
  26. CENTER_ON_LONGEST_AXIS = 0,
  27. MEDIAN_ON_LONGEST_AXIS = 1,
  28. NUM_SPLIT_METHODS = 2
  29. } split_method;
  30. public:
  31. inline WindingNumberAABB():
  32. total_positive_area(std::numeric_limits<double>::infinity()),
  33. split_method(MEDIAN_ON_LONGEST_AXIS)
  34. {}
  35. inline WindingNumberAABB(
  36. const Eigen::MatrixXd & V,
  37. const Eigen::MatrixXi & F);
  38. inline WindingNumberAABB(
  39. const WindingNumberTree<Point> & parent,
  40. const Eigen::MatrixXi & F);
  41. // Initialize some things
  42. inline void set_mesh(
  43. const Eigen::MatrixXd & V,
  44. const Eigen::MatrixXi & F);
  45. inline void init();
  46. inline bool inside(const Point & p) const;
  47. inline virtual void grow();
  48. // Compute min and max corners
  49. inline void compute_min_max_corners();
  50. inline double max_abs_winding_number(const Point & p) const;
  51. inline double max_simple_abs_winding_number(const Point & p) const;
  52. };
  53. }
  54. // Implementation
  55. #include "winding_number.h"
  56. #include "barycenter.h"
  57. #include "median.h"
  58. #include "doublearea.h"
  59. #include "per_face_normals.h"
  60. #include <limits>
  61. #include <vector>
  62. #include <iostream>
  63. // Minimum number of faces in a hierarchy element (this is probably dependent
  64. // on speed of machine and compiler optimization)
  65. #ifndef WindingNumberAABB_MIN_F
  66. # define WindingNumberAABB_MIN_F 100
  67. #endif
  68. template <typename Point>
  69. inline void igl::WindingNumberAABB<Point>::set_mesh(
  70. const Eigen::MatrixXd & V,
  71. const Eigen::MatrixXi & F)
  72. {
  73. igl::WindingNumberTree<Point>::set_mesh(V,F);
  74. init();
  75. }
  76. template <typename Point>
  77. inline void igl::WindingNumberAABB<Point>::init()
  78. {
  79. using namespace Eigen;
  80. assert(max_corner.size() == 3);
  81. assert(min_corner.size() == 3);
  82. compute_min_max_corners();
  83. VectorXd dblA;
  84. doublearea(this->getV(),this->getF(),dblA);
  85. total_positive_area = dblA.sum()/2.0;
  86. }
  87. template <typename Point>
  88. inline igl::WindingNumberAABB<Point>::WindingNumberAABB(
  89. const Eigen::MatrixXd & V,
  90. const Eigen::MatrixXi & F):
  91. WindingNumberTree<Point>(V,F),
  92. min_corner(),
  93. max_corner(),
  94. total_positive_area(std::numeric_limits<double>::infinity()),
  95. split_method(MEDIAN_ON_LONGEST_AXIS)
  96. {
  97. init();
  98. }
  99. template <typename Point>
  100. inline igl::WindingNumberAABB<Point>::WindingNumberAABB(
  101. const WindingNumberTree<Point> & parent,
  102. const Eigen::MatrixXi & F):
  103. WindingNumberTree<Point>(parent,F),
  104. min_corner(),
  105. max_corner(),
  106. total_positive_area(std::numeric_limits<double>::infinity()),
  107. split_method(MEDIAN_ON_LONGEST_AXIS)
  108. {
  109. init();
  110. }
  111. template <typename Point>
  112. inline void igl::WindingNumberAABB<Point>::grow()
  113. {
  114. using namespace std;
  115. using namespace Eigen;
  116. // Clear anything that already exists
  117. this->delete_children();
  118. //cout<<"cap.rows(): "<<this->getcap().rows()<<endl;
  119. //cout<<"F.rows(): "<<this->getF().rows()<<endl;
  120. // Base cases
  121. if(
  122. this->getF().rows() <= (WindingNumberAABB_MIN_F>0?WindingNumberAABB_MIN_F:0) ||
  123. (this->getcap().rows() - 2) >= this->getF().rows())
  124. {
  125. // Don't grow
  126. return;
  127. }
  128. // Compute longest direction
  129. int max_d = -1;
  130. double max_len = -numeric_limits<double>::infinity();
  131. for(int d = 0;d<min_corner.size();d++)
  132. {
  133. if( (max_corner[d] - min_corner[d]) > max_len )
  134. {
  135. max_len = (max_corner[d] - min_corner[d]);
  136. max_d = d;
  137. }
  138. }
  139. // Compute facet barycenters
  140. MatrixXd BC;
  141. barycenter(this->getV(),this->getF(),BC);
  142. // Blerg, why is selecting rows so difficult
  143. double split_value;
  144. // Split in longest direction
  145. switch(split_method)
  146. {
  147. case MEDIAN_ON_LONGEST_AXIS:
  148. // Determine median
  149. median(BC.col(max_d),split_value);
  150. break;
  151. default:
  152. assert(false);
  153. case CENTER_ON_LONGEST_AXIS:
  154. split_value = 0.5*(max_corner[max_d] + min_corner[max_d]);
  155. break;
  156. }
  157. //cout<<"c: "<<0.5*(max_corner[max_d] + min_corner[max_d])<<" "<<
  158. // "m: "<<split_value<<endl;;
  159. vector<int> id( this->getF().rows());
  160. for(int i = 0;i<this->getF().rows();i++)
  161. {
  162. if(BC(i,max_d) <= split_value)
  163. {
  164. id[i] = 0; //left
  165. }else
  166. {
  167. id[i] = 1; //right
  168. }
  169. }
  170. const int lefts = (int) count(id.begin(),id.end(),0);
  171. const int rights = (int) count(id.begin(),id.end(),1);
  172. if(lefts == 0 || rights == 0)
  173. {
  174. // badly balanced base case (could try to recut)
  175. return;
  176. }
  177. assert(lefts+rights == this->getF().rows());
  178. MatrixXi leftF(lefts, this->getF().cols());
  179. MatrixXi rightF(rights,this->getF().cols());
  180. int left_i = 0;
  181. int right_i = 0;
  182. for(int i = 0;i<this->getF().rows();i++)
  183. {
  184. if(id[i] == 0)
  185. {
  186. leftF.row(left_i++) = this->getF().row(i);
  187. }else if(id[i] == 1)
  188. {
  189. rightF.row(right_i++) = this->getF().row(i);
  190. }else
  191. {
  192. assert(false);
  193. }
  194. }
  195. assert(right_i == rightF.rows());
  196. assert(left_i == leftF.rows());
  197. // Finally actually grow children and Recursively grow
  198. WindingNumberAABB<Point> * leftWindingNumberAABB = new WindingNumberAABB<Point>(*this,leftF);
  199. leftWindingNumberAABB->grow();
  200. this->children.push_back(leftWindingNumberAABB);
  201. WindingNumberAABB<Point> * rightWindingNumberAABB = new WindingNumberAABB<Point>(*this,rightF);
  202. rightWindingNumberAABB->grow();
  203. this->children.push_back(rightWindingNumberAABB);
  204. }
  205. template <typename Point>
  206. inline bool igl::WindingNumberAABB<Point>::inside(const Point & p) const
  207. {
  208. assert(p.size() == max_corner.size());
  209. assert(p.size() == min_corner.size());
  210. for(int i = 0;i<p.size();i++)
  211. {
  212. //// Perfect matching is **not** robust
  213. //if( p(i) < min_corner(i) || p(i) >= max_corner(i))
  214. // **MUST** be conservative
  215. if( p(i) < min_corner(i) || p(i) > max_corner(i))
  216. {
  217. return false;
  218. }
  219. }
  220. return true;
  221. }
  222. template <typename Point>
  223. inline void igl::WindingNumberAABB<Point>::compute_min_max_corners()
  224. {
  225. using namespace std;
  226. // initialize corners
  227. for(int d = 0;d<min_corner.size();d++)
  228. {
  229. min_corner[d] = numeric_limits<double>::infinity();
  230. max_corner[d] = -numeric_limits<double>::infinity();
  231. }
  232. this->center = Point(0,0,0);
  233. // Loop over facets
  234. for(int i = 0;i<this->getF().rows();i++)
  235. {
  236. for(int j = 0;j<this->getF().cols();j++)
  237. {
  238. for(int d = 0;d<min_corner.size();d++)
  239. {
  240. min_corner[d] =
  241. this->getV()(this->getF()(i,j),d) < min_corner[d] ?
  242. this->getV()(this->getF()(i,j),d) : min_corner[d];
  243. max_corner[d] =
  244. this->getV()(this->getF()(i,j),d) > max_corner[d] ?
  245. this->getV()(this->getF()(i,j),d) : max_corner[d];
  246. }
  247. // This is biased toward vertices incident on more than one face, but
  248. // perhaps that's good
  249. this->center += this->getV().row(this->getF()(i,j));
  250. }
  251. }
  252. // Average
  253. this->center.array() /= this->getF().size();
  254. //cout<<"min_corner: "<<this->min_corner.transpose()<<endl;
  255. //cout<<"Center: "<<this->center.transpose()<<endl;
  256. //cout<<"max_corner: "<<this->max_corner.transpose()<<endl;
  257. //cout<<"Diag center: "<<((this->max_corner + this->min_corner)*0.5).transpose()<<endl;
  258. //cout<<endl;
  259. this->radius = (max_corner-min_corner).norm()/2.0;
  260. }
  261. template <typename Point>
  262. inline double igl::WindingNumberAABB<Point>::max_abs_winding_number(const Point & p) const
  263. {
  264. using namespace std;
  265. // Only valid if not inside
  266. if(inside(p))
  267. {
  268. return numeric_limits<double>::infinity();
  269. }
  270. // Q: we know the total positive area so what's the most this could project
  271. // to? Remember it could be layered in the same direction.
  272. return numeric_limits<double>::infinity();
  273. }
  274. template <typename Point>
  275. inline double igl::WindingNumberAABB<Point>::max_simple_abs_winding_number(const Point & p) const
  276. {
  277. using namespace std;
  278. using namespace Eigen;
  279. // Only valid if not inside
  280. if(inside(p))
  281. {
  282. return numeric_limits<double>::infinity();
  283. }
  284. // Max simple is the same as sum of positive winding number contributions of
  285. // bounding box
  286. // begin precomputation
  287. //MatrixXd BV((int)pow(2,3),3);
  288. MatrixXd BV((int)(1<<3),3);
  289. BV <<
  290. min_corner[0],min_corner[1],min_corner[2],
  291. min_corner[0],min_corner[1],max_corner[2],
  292. min_corner[0],max_corner[1],min_corner[2],
  293. min_corner[0],max_corner[1],max_corner[2],
  294. max_corner[0],min_corner[1],min_corner[2],
  295. max_corner[0],min_corner[1],max_corner[2],
  296. max_corner[0],max_corner[1],min_corner[2],
  297. max_corner[0],max_corner[1],max_corner[2];
  298. MatrixXi BF(2*2*3,3);
  299. BF <<
  300. 0,6,4,
  301. 0,2,6,
  302. 0,3,2,
  303. 0,1,3,
  304. 2,7,6,
  305. 2,3,7,
  306. 4,6,7,
  307. 4,7,5,
  308. 0,4,5,
  309. 0,5,1,
  310. 1,5,7,
  311. 1,7,3;
  312. MatrixXd BFN;
  313. per_face_normals(BV,BF,BFN);
  314. // end of precomputation
  315. // Only keep those with positive dot products
  316. MatrixXi PBF(BF.rows(),BF.cols());
  317. int pbfi = 0;
  318. Point p2c = 0.5*(min_corner+max_corner)-p;
  319. for(int i = 0;i<BFN.rows();i++)
  320. {
  321. if(p2c.dot(BFN.row(i)) > 0)
  322. {
  323. PBF.row(pbfi++) = BF.row(i);
  324. }
  325. }
  326. PBF.conservativeResize(pbfi,PBF.cols());
  327. double w = numeric_limits<double>::infinity();
  328. igl::winding_number_3(
  329. BV.data(),
  330. BV.rows(),
  331. PBF.data(),
  332. PBF.rows(),
  333. p.data(),
  334. 1,
  335. &w);
  336. return w;
  337. }
  338. //// Explicit instanciation
  339. //template class igl::WindingNumberAABB<Eigen::Vector3d >;
  340. #endif