WindingNumberTree.h 14 KB

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