WindingNumberTree.h 16 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523
  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. #ifndef IGL_WINDINGNUMBERTREE_H
  9. #define IGL_WINDINGNUMBERTREE_H
  10. #include <list>
  11. #include <map>
  12. #include <Eigen/Dense>
  13. #include "WindingNumberMethod.h"
  14. namespace igl
  15. {
  16. // Space partitioning tree for computing winding number hierarchically.
  17. //
  18. // Templates:
  19. // Point type for points in space, e.g. Eigen::Vector3d
  20. template <
  21. typename Point,
  22. typename DerivedV,
  23. typename DerivedF >
  24. class WindingNumberTree
  25. {
  26. public:
  27. // Method to use (see enum above)
  28. //static double min_max_w;
  29. static std::map<
  30. std::pair<const WindingNumberTree*,const WindingNumberTree*>,
  31. typename DerivedV::Scalar>
  32. cached;
  33. // This is only need to fill in references, it should never actually be touched
  34. // and shouldn't cause race conditions. (This is a hack, but I think it's "safe")
  35. static DerivedV dummyV;
  36. protected:
  37. WindingNumberMethod method;
  38. const WindingNumberTree * parent;
  39. std::list<WindingNumberTree * > children;
  40. typedef
  41. Eigen::Matrix<typename DerivedV::Scalar,Eigen::Dynamic,Eigen::Dynamic>
  42. MatrixXS;
  43. typedef
  44. Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,Eigen::Dynamic>
  45. MatrixXF;
  46. //// List of boundary edges (recall edges are vertices in 2d)
  47. //const Eigen::MatrixXi boundary;
  48. // Base mesh vertices
  49. DerivedV & V;
  50. // Base mesh vertices with duplicates removed
  51. MatrixXS SV;
  52. // Facets in this bounding volume
  53. MatrixXF F;
  54. // Tesselated boundary curve
  55. MatrixXF cap;
  56. // Upper Bound on radius of enclosing ball
  57. typename DerivedV::Scalar radius;
  58. // (Approximate) center (of mass)
  59. Point center;
  60. public:
  61. inline WindingNumberTree();
  62. // For root
  63. inline WindingNumberTree(
  64. const Eigen::MatrixBase<DerivedV> & V,
  65. const Eigen::MatrixBase<DerivedF> & F);
  66. // For chilluns
  67. inline WindingNumberTree(
  68. const WindingNumberTree<Point,DerivedV,DerivedF> & parent,
  69. const Eigen::MatrixBase<DerivedF> & F);
  70. inline virtual ~WindingNumberTree();
  71. inline void delete_children();
  72. inline virtual void set_mesh(
  73. const Eigen::MatrixBase<DerivedV> & V,
  74. const Eigen::MatrixBase<DerivedF> & F);
  75. // Set method
  76. inline void set_method( const WindingNumberMethod & m);
  77. public:
  78. inline const DerivedV & getV() const;
  79. inline const MatrixXF & getF() const;
  80. inline const MatrixXF & getcap() const;
  81. // Grow the Tree recursively
  82. inline virtual void grow();
  83. // Determine whether a given point is inside the bounding
  84. //
  85. // Inputs:
  86. // p query point
  87. // Returns true if the point p is inside this bounding volume
  88. inline virtual bool inside(const Point & p) const;
  89. // Compute the (partial) winding number of a given point p
  90. // According to method
  91. //
  92. // Inputs:
  93. // p query point
  94. // Returns winding number
  95. inline typename DerivedV::Scalar winding_number(const Point & p) const;
  96. // Same as above, but always computes winding number using exact method
  97. // (sum over every facet)
  98. inline typename DerivedV::Scalar winding_number_all(const Point & p) const;
  99. // Same as above, but always computes using sum over tesslated boundary
  100. inline typename DerivedV::Scalar winding_number_boundary(const Point & p) const;
  101. //// Same as winding_number above, but if max_simple_abs_winding_number is
  102. //// less than some threshold min_max_w just return 0 (colloquially the "fast
  103. //// multipole method)
  104. ////
  105. ////
  106. //// Inputs:
  107. //// p query point
  108. //// min_max_w minimum max simple w to be processed
  109. //// Returns approximate winding number
  110. //double winding_number_approx_simple(
  111. // const Point & p,
  112. // const double min_max_w);
  113. // Print contents of Tree
  114. //
  115. // Optional input:
  116. // tab tab to show depth
  117. inline void print(const char * tab="");
  118. // Determine max absolute winding number
  119. //
  120. // Inputs:
  121. // p query point
  122. // Returns max winding number of
  123. inline virtual typename DerivedV::Scalar max_abs_winding_number(const Point & p) const;
  124. // Same as above, but stronger assumptions on (V,F). Assumes (V,F) is a
  125. // simple polyhedron
  126. inline virtual typename DerivedV::Scalar max_simple_abs_winding_number(const Point & p) const;
  127. // Compute or read cached winding number for point p with respect to mesh
  128. // in bounding box, recursing according to approximation criteria
  129. //
  130. // Inputs:
  131. // p query point
  132. // that WindingNumberTree containing mesh w.r.t. which we're computing w.n.
  133. // Returns cached winding number
  134. inline virtual typename DerivedV::Scalar cached_winding_number(const WindingNumberTree & that, const Point & p) const;
  135. };
  136. }
  137. // Implementation
  138. #include "WindingNumberTree.h"
  139. #include "winding_number.h"
  140. #include "triangle_fan.h"
  141. #include "exterior_edges.h"
  142. #include <igl/PI.h>
  143. #include <igl/remove_duplicate_vertices.h>
  144. #include <iostream>
  145. #include <limits>
  146. //template <typename Point, typename DerivedV, typename DerivedF>
  147. //WindingNumberMethod WindingNumberTree<Point,DerivedV,DerivedF>::method = EXACT_WINDING_NUMBER_METHOD;
  148. //template <typename Point, typename DerivedV, typename DerivedF>
  149. //double WindingNumberTree<Point,DerivedV,DerivedF>::min_max_w = 0;
  150. template <typename Point, typename DerivedV, typename DerivedF>
  151. std::map< std::pair<const igl::WindingNumberTree<Point,DerivedV,DerivedF>*,const igl::WindingNumberTree<Point,DerivedV,DerivedF>*>, typename DerivedV::Scalar>
  152. igl::WindingNumberTree<Point,DerivedV,DerivedF>::cached;
  153. template <typename Point, typename DerivedV, typename DerivedF>
  154. inline igl::WindingNumberTree<Point,DerivedV,DerivedF>::WindingNumberTree():
  155. method(EXACT_WINDING_NUMBER_METHOD),
  156. parent(NULL),
  157. V(dummyV),
  158. SV(),
  159. F(),
  160. //boundary(igl::boundary_facets<Eigen::MatrixXi,Eigen::MatrixXi>(F))
  161. cap(),
  162. radius(std::numeric_limits<typename DerivedV::Scalar>::infinity()),
  163. center(0,0,0)
  164. {
  165. }
  166. template <typename Point, typename DerivedV, typename DerivedF>
  167. inline igl::WindingNumberTree<Point,DerivedV,DerivedF>::WindingNumberTree(
  168. const Eigen::MatrixBase<DerivedV> & _V,
  169. const Eigen::MatrixBase<DerivedF> & _F):
  170. method(EXACT_WINDING_NUMBER_METHOD),
  171. parent(NULL),
  172. V(dummyV),
  173. SV(),
  174. F(),
  175. //boundary(igl::boundary_facets<Eigen::MatrixXi,Eigen::MatrixXi>(F))
  176. cap(),
  177. radius(std::numeric_limits<typename DerivedV::Scalar>::infinity()),
  178. center(0,0,0)
  179. {
  180. set_mesh(_V,_F);
  181. }
  182. template <typename Point, typename DerivedV, typename DerivedF>
  183. inline void igl::WindingNumberTree<Point,DerivedV,DerivedF>::set_mesh(
  184. const Eigen::MatrixBase<DerivedV> & _V,
  185. const Eigen::MatrixBase<DerivedF> & _F)
  186. {
  187. using namespace std;
  188. // Remove any exactly duplicate vertices
  189. // Q: Can this ever increase the complexity of the boundary?
  190. // Q: Would we gain even more by remove almost exactly duplicate vertices?
  191. MatrixXF SF,SVI,SVJ;
  192. igl::remove_duplicate_vertices(_V,_F,0.0,SV,SVI,SVJ,F);
  193. triangle_fan(igl::exterior_edges(F),cap);
  194. V = SV;
  195. }
  196. template <typename Point, typename DerivedV, typename DerivedF>
  197. inline igl::WindingNumberTree<Point,DerivedV,DerivedF>::WindingNumberTree(
  198. const igl::WindingNumberTree<Point,DerivedV,DerivedF> & parent,
  199. const Eigen::MatrixBase<DerivedF> & _F):
  200. method(parent.method),
  201. parent(&parent),
  202. V(parent.V),
  203. SV(),
  204. F(_F),
  205. cap(triangle_fan(igl::exterior_edges(_F)))
  206. {
  207. }
  208. template <typename Point, typename DerivedV, typename DerivedF>
  209. inline igl::WindingNumberTree<Point,DerivedV,DerivedF>::~WindingNumberTree()
  210. {
  211. delete_children();
  212. }
  213. template <typename Point, typename DerivedV, typename DerivedF>
  214. inline void igl::WindingNumberTree<Point,DerivedV,DerivedF>::delete_children()
  215. {
  216. using namespace std;
  217. // Delete children
  218. typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::iterator cit = children.begin();
  219. while(cit != children.end())
  220. {
  221. // clear the memory of this item
  222. delete (* cit);
  223. // erase from list, returns next element in iterator
  224. cit = children.erase(cit);
  225. }
  226. }
  227. template <typename Point, typename DerivedV, typename DerivedF>
  228. inline void igl::WindingNumberTree<Point,DerivedV,DerivedF>::set_method(const WindingNumberMethod & m)
  229. {
  230. this->method = m;
  231. for(auto child : children)
  232. {
  233. child->set_method(m);
  234. }
  235. }
  236. template <typename Point, typename DerivedV, typename DerivedF>
  237. inline const DerivedV & igl::WindingNumberTree<Point,DerivedV,DerivedF>::getV() const
  238. {
  239. return V;
  240. }
  241. template <typename Point, typename DerivedV, typename DerivedF>
  242. inline const typename igl::WindingNumberTree<Point,DerivedV,DerivedF>::MatrixXF&
  243. igl::WindingNumberTree<Point,DerivedV,DerivedF>::getF() const
  244. {
  245. return F;
  246. }
  247. template <typename Point, typename DerivedV, typename DerivedF>
  248. inline const typename igl::WindingNumberTree<Point,DerivedV,DerivedF>::MatrixXF&
  249. igl::WindingNumberTree<Point,DerivedV,DerivedF>::getcap() const
  250. {
  251. return cap;
  252. }
  253. template <typename Point, typename DerivedV, typename DerivedF>
  254. inline void igl::WindingNumberTree<Point,DerivedV,DerivedF>::grow()
  255. {
  256. // Don't grow
  257. return;
  258. }
  259. template <typename Point, typename DerivedV, typename DerivedF>
  260. inline bool igl::WindingNumberTree<Point,DerivedV,DerivedF>::inside(const Point & /*p*/) const
  261. {
  262. return true;
  263. }
  264. template <typename Point, typename DerivedV, typename DerivedF>
  265. inline typename DerivedV::Scalar
  266. igl::WindingNumberTree<Point,DerivedV,DerivedF>::winding_number(const Point & p) const
  267. {
  268. using namespace std;
  269. //cout<<"+"<<boundary.rows();
  270. // If inside then we need to be careful
  271. if(inside(p))
  272. {
  273. // If not a leaf then recurse
  274. if(children.size()>0)
  275. {
  276. // Recurse on each child and accumulate
  277. typename DerivedV::Scalar sum = 0;
  278. for(
  279. typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::const_iterator cit = children.begin();
  280. cit != children.end();
  281. cit++)
  282. {
  283. switch(method)
  284. {
  285. case EXACT_WINDING_NUMBER_METHOD:
  286. sum += (*cit)->winding_number(p);
  287. break;
  288. case APPROX_SIMPLE_WINDING_NUMBER_METHOD:
  289. case APPROX_CACHE_WINDING_NUMBER_METHOD:
  290. //if((*cit)->max_simple_abs_winding_number(p) > min_max_w)
  291. //{
  292. sum += (*cit)->winding_number(p);
  293. //}
  294. break;
  295. default:
  296. assert(false);
  297. break;
  298. }
  299. }
  300. return sum;
  301. }else
  302. {
  303. return winding_number_all(p);
  304. }
  305. }else{
  306. // Otherwise we can just consider boundary
  307. // Q: If we using the "multipole" method should we also subdivide the
  308. // boundary case?
  309. if((cap.rows() - 2) < F.rows())
  310. {
  311. switch(method)
  312. {
  313. case EXACT_WINDING_NUMBER_METHOD:
  314. return winding_number_boundary(p);
  315. case APPROX_SIMPLE_WINDING_NUMBER_METHOD:
  316. {
  317. typename DerivedV::Scalar dist = (p-center).norm();
  318. // Radius is already an overestimate of inside
  319. if(dist>1.0*radius)
  320. {
  321. return 0;
  322. }else
  323. {
  324. return winding_number_boundary(p);
  325. }
  326. }
  327. case APPROX_CACHE_WINDING_NUMBER_METHOD:
  328. {
  329. return parent->cached_winding_number(*this,p);
  330. }
  331. default: assert(false);break;
  332. }
  333. }else
  334. {
  335. // doesn't pay off to use boundary
  336. return winding_number_all(p);
  337. }
  338. }
  339. return 0;
  340. }
  341. template <typename Point, typename DerivedV, typename DerivedF>
  342. inline typename DerivedV::Scalar
  343. igl::WindingNumberTree<Point,DerivedV,DerivedF>::winding_number_all(const Point & p) const
  344. {
  345. typename DerivedV::Scalar w = 0;
  346. igl::winding_number_3(
  347. V.data(),
  348. V.rows(),
  349. F.data(),
  350. F.rows(),
  351. p.data(),
  352. 1,
  353. &w);
  354. return w;
  355. }
  356. template <typename Point, typename DerivedV, typename DerivedF>
  357. inline typename DerivedV::Scalar
  358. igl::WindingNumberTree<Point,DerivedV,DerivedF>::winding_number_boundary(const Point & p) const
  359. {
  360. using namespace Eigen;
  361. using namespace std;
  362. typename DerivedV::Scalar w = 0;
  363. // `cap` is already flipped inside out, so we don't need to flip sign of w
  364. igl::winding_number_3(
  365. V.data(),
  366. V.rows(),
  367. cap.data(),
  368. cap.rows(),
  369. &p[0],
  370. 1,
  371. &w);
  372. return w;
  373. }
  374. //template <typename Point, typename DerivedV, typename DerivedF>
  375. //inline double igl::WindingNumberTree<Point,DerivedV,DerivedF>::winding_number_approx_simple(
  376. // const Point & p,
  377. // const double min_max_w)
  378. //{
  379. // using namespace std;
  380. // if(max_simple_abs_winding_number(p) > min_max_w)
  381. // {
  382. // return winding_number(p);
  383. // }else
  384. // {
  385. // cout<<"Skipped! "<<max_simple_abs_winding_number(p)<<"<"<<min_max_w<<endl;
  386. // return 0;
  387. // }
  388. //}
  389. template <typename Point, typename DerivedV, typename DerivedF>
  390. inline void igl::WindingNumberTree<Point,DerivedV,DerivedF>::print(const char * tab)
  391. {
  392. using namespace std;
  393. // Print all facets
  394. cout<<tab<<"["<<endl<<F<<endl<<"]";
  395. // Print children
  396. for(
  397. typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::iterator cit = children.begin();
  398. cit != children.end();
  399. cit++)
  400. {
  401. cout<<","<<endl;
  402. (*cit)->print((string(tab)+"").c_str());
  403. }
  404. }
  405. template <typename Point, typename DerivedV, typename DerivedF>
  406. inline typename DerivedV::Scalar
  407. igl::WindingNumberTree<Point,DerivedV,DerivedF>::max_abs_winding_number(const Point & /*p*/) const
  408. {
  409. return std::numeric_limits<typename DerivedV::Scalar>::infinity();
  410. }
  411. template <typename Point, typename DerivedV, typename DerivedF>
  412. inline typename DerivedV::Scalar
  413. igl::WindingNumberTree<Point,DerivedV,DerivedF>::max_simple_abs_winding_number(
  414. const Point & /*p*/) const
  415. {
  416. using namespace std;
  417. return numeric_limits<typename DerivedV::Scalar>::infinity();
  418. }
  419. template <typename Point, typename DerivedV, typename DerivedF>
  420. inline typename DerivedV::Scalar
  421. igl::WindingNumberTree<Point,DerivedV,DerivedF>::cached_winding_number(
  422. const igl::WindingNumberTree<Point,DerivedV,DerivedF> & that,
  423. const Point & p) const
  424. {
  425. using namespace std;
  426. // Simple metric for `is_far`
  427. //
  428. // this that
  429. // --------
  430. // ----- / | \ .
  431. // / r \ / R \ .
  432. // | p ! | | ! |
  433. // \_____/ \ /
  434. // \________/
  435. //
  436. //
  437. // a = angle formed by trapazoid formed by raising sides with lengths r and R
  438. // at respective centers.
  439. //
  440. // a = atan2(R-r,d), where d is the distance between centers
  441. // That should be bigger (what about parent? what about sister?)
  442. bool is_far = this->radius<that.radius;
  443. if(is_far)
  444. {
  445. typename DerivedV::Scalar a = atan2(
  446. that.radius - this->radius,
  447. (that.center - this->center).norm());
  448. assert(a>0);
  449. is_far = (a<PI/8.0);
  450. }
  451. if(is_far)
  452. {
  453. // Not implemented yet
  454. pair<const WindingNumberTree*,const WindingNumberTree*> this_that(this,&that);
  455. // Need to compute it for first time?
  456. if(cached.count(this_that)==0)
  457. {
  458. cached[this_that] =
  459. that.winding_number_boundary(this->center);
  460. }
  461. return cached[this_that];
  462. }else if(children.size() == 0)
  463. {
  464. // not far and hierarchy ended too soon: can't use cache
  465. return that.winding_number_boundary(p);
  466. }else
  467. {
  468. for(
  469. typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::const_iterator cit = children.begin();
  470. cit != children.end();
  471. cit++)
  472. {
  473. if((*cit)->inside(p))
  474. {
  475. return (*cit)->cached_winding_number(that,p);
  476. }
  477. }
  478. // Not inside any children? This can totally happen because bounding boxes
  479. // are set to bound contained facets. So sibilings may overlap and their
  480. // union may not contain their parent (though, their union is certainly a
  481. // subset of their parent).
  482. assert(false);
  483. }
  484. return 0;
  485. }
  486. // Explicit instanciation of static variable
  487. template <
  488. typename Point,
  489. typename DerivedV,
  490. typename DerivedF >
  491. DerivedV igl::WindingNumberTree<Point,DerivedV,DerivedF>::dummyV;
  492. #endif