libigl v2.5.0
Loading...
Searching...
No Matches
WindingNumberTree.h
Go to the documentation of this file.
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
15namespace igl
16{
20 template <
21 typename Point,
22 typename DerivedV,
23 typename DerivedF >
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>
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:
39 std::list<WindingNumberTree * > children;
40 typedef
41 Eigen::Matrix<typename DerivedV::Scalar,Eigen::Dynamic,Eigen::Dynamic>
43 typedef
44 Eigen::Matrix<typename DerivedF::Scalar,Eigen::Dynamic,Eigen::Dynamic>
47 //const Eigen::MatrixXi boundary;
48 // Base mesh vertices
49 DerivedV & V;
50 // Base mesh vertices with duplicates removed
52 // Facets in this bounding volume
54 // Tessellated boundary curve
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(
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 tessllated boundary
100 inline typename DerivedV::Scalar winding_number_boundary(const Point & p) const;
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
138// Implementation
139
140#include "WindingNumberTree.h"
141#include "winding_number.h"
142#include "triangle_fan.h"
143#include "exterior_edges.h"
144
145#include "PI.h"
147
148#include <iostream>
149#include <limits>
150
151//template <typename Point, typename DerivedV, typename DerivedF>
152//WindingNumberMethod WindingNumberTree<Point,DerivedV,DerivedF>::method = EXACT_WINDING_NUMBER_METHOD;
153//template <typename Point, typename DerivedV, typename DerivedF>
154//double WindingNumberTree<Point,DerivedV,DerivedF>::min_max_w = 0;
155template <typename Point, typename DerivedV, typename DerivedF>
156std::map< std::pair<const igl::WindingNumberTree<Point,DerivedV,DerivedF>*,const igl::WindingNumberTree<Point,DerivedV,DerivedF>*>, typename DerivedV::Scalar>
158
159template <typename Point, typename DerivedV, typename DerivedF>
162 parent(NULL),
163 V(dummyV),
164 SV(),
165 F(),
166 cap(),
167 radius(std::numeric_limits<typename DerivedV::Scalar>::infinity()),
168 center(0,0,0)
169{
170}
171
172template <typename Point, typename DerivedV, typename DerivedF>
174 const Eigen::MatrixBase<DerivedV> & _V,
175 const Eigen::MatrixBase<DerivedF> & _F):
177 parent(NULL),
178 V(dummyV),
179 SV(),
180 F(),
181 cap(),
182 radius(std::numeric_limits<typename DerivedV::Scalar>::infinity()),
183 center(0,0,0)
184{
185 set_mesh(_V,_F);
186}
187
188template <typename Point, typename DerivedV, typename DerivedF>
190 const Eigen::MatrixBase<DerivedV> & _V,
191 const Eigen::MatrixBase<DerivedF> & _F)
192{
193 using namespace std;
194 // Remove any exactly duplicate vertices
195 // Q: Can this ever increase the complexity of the boundary?
196 // Q: Would we gain even more by remove almost exactly duplicate vertices?
197 Eigen::Matrix<typename MatrixXF::Scalar,Eigen::Dynamic,1> SVI,SVJ;
198 igl::remove_duplicate_vertices(_V,_F,0.0,SV,SVI,SVJ,F);
200 V = SV;
201}
202
203template <typename Point, typename DerivedV, typename DerivedF>
206 const Eigen::MatrixBase<DerivedF> & _F):
207 method(parent.method),
208 parent(&parent),
209 V(parent.V),
210 SV(),
211 F(_F),
213{
214}
215
216template <typename Point, typename DerivedV, typename DerivedF>
218{
219 delete_children();
220}
221
222template <typename Point, typename DerivedV, typename DerivedF>
224{
225 using namespace std;
226 // Delete children
227 typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::iterator cit = children.begin();
228 while(cit != children.end())
229 {
230 // clear the memory of this item
231 delete (* cit);
232 // erase from list, returns next element in iterator
233 cit = children.erase(cit);
234 }
235}
236
237template <typename Point, typename DerivedV, typename DerivedF>
239{
240 this->method = m;
241 for(auto child : children)
242 {
243 child->set_method(m);
244 }
245}
246
247template <typename Point, typename DerivedV, typename DerivedF>
249{
250 return V;
251}
252
253template <typename Point, typename DerivedV, typename DerivedF>
259
260template <typename Point, typename DerivedV, typename DerivedF>
266
267template <typename Point, typename DerivedV, typename DerivedF>
269{
270 // Don't grow
271 return;
272}
273
274template <typename Point, typename DerivedV, typename DerivedF>
276{
277 return true;
278}
279
280template <typename Point, typename DerivedV, typename DerivedF>
281inline typename DerivedV::Scalar
283{
284 using namespace std;
285 //cout<<"+"<<boundary.rows();
286 // If inside then we need to be careful
287 if(inside(p))
288 {
289 // If not a leaf then recurse
290 if(children.size()>0)
291 {
292 // Recurse on each child and accumulate
293 typename DerivedV::Scalar sum = 0;
294 for(
295 typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::const_iterator cit = children.begin();
296 cit != children.end();
297 cit++)
298 {
299 switch(method)
300 {
302 sum += (*cit)->winding_number(p);
303 break;
306 //if((*cit)->max_simple_abs_winding_number(p) > min_max_w)
307 //{
308 sum += (*cit)->winding_number(p);
309 //}
310 break;
311 default:
312 assert(false);
313 break;
314 }
315 }
316 return sum;
317 }else
318 {
319 return winding_number_all(p);
320 }
321 }else{
322 // Otherwise we can just consider boundary
323 // Q: If we using the "multipole" method should we also subdivide the
324 // boundary case?
325 if((cap.rows() - 2) < F.rows())
326 {
327 switch(method)
328 {
330 return winding_number_boundary(p);
332 {
333 typename DerivedV::Scalar dist = (p-center).norm();
334 // Radius is already an overestimate of inside
335 if(dist>1.0*radius)
336 {
337 return 0;
338 }else
339 {
340 return winding_number_boundary(p);
341 }
342 }
344 {
345 return parent->cached_winding_number(*this,p);
346 }
347 default: assert(false);break;
348 }
349 }else
350 {
351 // doesn't pay off to use boundary
352 return winding_number_all(p);
353 }
354 }
355 return 0;
356}
357
358template <typename Point, typename DerivedV, typename DerivedF>
359inline typename DerivedV::Scalar
364
365template <typename Point, typename DerivedV, typename DerivedF>
366inline typename DerivedV::Scalar
368{
369 using namespace Eigen;
370 using namespace std;
371 return igl::winding_number(V,cap,p);
372}
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
390template <typename Point, typename DerivedV, typename DerivedF>
392{
393 using namespace std;
394 // Print all facets
395 cout<<tab<<"["<<endl<<F<<endl<<"]";
396 // Print children
397 for(
398 typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::iterator cit = children.begin();
399 cit != children.end();
400 cit++)
401 {
402 cout<<","<<endl;
403 (*cit)->print((string(tab)+"").c_str());
404 }
405}
406
407template <typename Point, typename DerivedV, typename DerivedF>
408inline typename DerivedV::Scalar
410{
411 return std::numeric_limits<typename DerivedV::Scalar>::infinity();
412}
413
414template <typename Point, typename DerivedV, typename DerivedF>
415inline typename DerivedV::Scalar
417 const Point & /*p*/) const
418{
419 using namespace std;
420 return numeric_limits<typename DerivedV::Scalar>::infinity();
421}
422
423template <typename Point, typename DerivedV, typename DerivedF>
424inline typename DerivedV::Scalar
427 const Point & p) const
428{
429 using namespace std;
430 // Simple metric for `is_far`
431 //
432 // this that
433 // --------
434 // ----- / | \ .
435 // / r \ / R \ .
436 // | p ! | | ! |
437 // \_____/ \ /
438 // \________/
439 //
440 //
441 // a = angle formed by trapazoid formed by raising sides with lengths r and R
442 // at respective centers.
443 //
444 // a = atan2(R-r,d), where d is the distance between centers
445
446 // That should be bigger (what about parent? what about sister?)
447 bool is_far = this->radius<that.radius;
448 if(is_far)
449 {
450 typename DerivedV::Scalar a = atan2(
451 that.radius - this->radius,
452 (that.center - this->center).norm());
453 assert(a>0);
454 is_far = (a<PI/8.0);
455 }
456
457 if(is_far)
458 {
459 // Not implemented yet
460 pair<const WindingNumberTree*,const WindingNumberTree*> this_that(this,&that);
461 // Need to compute it for first time?
462 if(cached.count(this_that)==0)
463 {
464 cached[this_that] =
465 that.winding_number_boundary(this->center);
466 }
467 return cached[this_that];
468 }else if(children.size() == 0)
469 {
470 // not far and hierarchy ended too soon: can't use cache
471 return that.winding_number_boundary(p);
472 }else
473 {
474 for(
475 typename list<WindingNumberTree<Point,DerivedV,DerivedF>* >::const_iterator cit = children.begin();
476 cit != children.end();
477 cit++)
478 {
479 if((*cit)->inside(p))
480 {
481 return (*cit)->cached_winding_number(that,p);
482 }
483 }
484 // Not inside any children? This can totally happen because bounding boxes
485 // are set to bound contained facets. So sibilings may overlap and their
486 // union may not contain their parent (though, their union is certainly a
487 // subset of their parent).
488 assert(false);
489 }
490 return 0;
491}
492
493// Explicit instantiation of static variable
494template <
495 typename Point,
496 typename DerivedV,
497 typename DerivedF >
499
500#endif
Space partitioning tree for computing winding number hierarchically.
Definition WindingNumberTree.h:25
virtual bool inside(const Point &p) const
Definition WindingNumberTree.h:275
Eigen::Matrix< typename DerivedF::Scalar, Eigen::Dynamic, Eigen::Dynamic > MatrixXF
Definition WindingNumberTree.h:45
virtual DerivedV::Scalar max_simple_abs_winding_number(const Point &p) const
Definition WindingNumberTree.h:416
const WindingNumberTree * parent
Definition WindingNumberTree.h:38
DerivedV & V
Definition WindingNumberTree.h:49
std::list< WindingNumberTree * > children
Definition WindingNumberTree.h:39
WindingNumberTree()
Definition WindingNumberTree.h:160
virtual void set_mesh(const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F)
Definition WindingNumberTree.h:189
const MatrixXF & getcap() const
Definition WindingNumberTree.h:262
WindingNumberMethod method
Definition WindingNumberTree.h:37
virtual DerivedV::Scalar max_abs_winding_number(const Point &p) const
Definition WindingNumberTree.h:409
void print(const char *tab="")
Definition WindingNumberTree.h:391
const MatrixXF & getF() const
Definition WindingNumberTree.h:255
virtual DerivedV::Scalar cached_winding_number(const WindingNumberTree &that, const Point &p) const
Definition WindingNumberTree.h:425
void set_method(const WindingNumberMethod &m)
Definition WindingNumberTree.h:238
MatrixXF cap
Definition WindingNumberTree.h:55
const DerivedV & getV() const
Definition WindingNumberTree.h:248
DerivedV::Scalar winding_number_all(const Point &p) const
Definition WindingNumberTree.h:360
virtual void grow()
Definition WindingNumberTree.h:268
DerivedV::Scalar radius
Definition WindingNumberTree.h:57
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, Eigen::Dynamic > MatrixXS
Definition WindingNumberTree.h:42
static std::map< std::pair< const WindingNumberTree *, const WindingNumberTree * >, typename DerivedV::Scalar > cached
Definition WindingNumberTree.h:32
void delete_children()
Definition WindingNumberTree.h:223
Point center
Definition WindingNumberTree.h:59
virtual ~WindingNumberTree()
Definition WindingNumberTree.h:217
MatrixXF F
Definition WindingNumberTree.h:53
DerivedV::Scalar winding_number(const Point &p) const
Definition WindingNumberTree.h:282
DerivedV::Scalar winding_number_boundary(const Point &p) const
Definition WindingNumberTree.h:367
MatrixXS SV
Definition WindingNumberTree.h:51
static DerivedV dummyV
Definition WindingNumberTree.h:35
Definition AABB.h:17
void exterior_edges(const Eigen::MatrixXi &F, Eigen::MatrixXi &E)
Determines boundary "edges" and also edges with an odd number of occurrences where seeing edge (i,...
void remove_duplicate_vertices(const Eigen::MatrixBase< DerivedV > &V, const double epsilon, Eigen::PlainObjectBase< DerivedSV > &SV, Eigen::PlainObjectBase< DerivedSVI > &SVI, Eigen::PlainObjectBase< DerivedSVJ > &SVJ)
Remove duplicate vertices upto a uniqueness tolerance (epsilon)
void winding_number(const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, const Eigen::MatrixBase< DerivedO > &O, Eigen::PlainObjectBase< DerivedW > &W)
Computes the generalized winding number at each dim-dimensional query point in O with respect to the ...
WindingNumberMethod
Definition WindingNumberMethod.h:13
@ APPROX_SIMPLE_WINDING_NUMBER_METHOD
Definition WindingNumberMethod.h:17
@ APPROX_CACHE_WINDING_NUMBER_METHOD
Definition WindingNumberMethod.h:19
@ EXACT_WINDING_NUMBER_METHOD
Definition WindingNumberMethod.h:15
void triangle_fan(const Eigen::MatrixBase< DerivedE > &E, Eigen::PlainObjectBase< Derivedcap > &cap)
Given a list of faces tessellate all of the "exterior" edges forming another list of.
void sum(const Eigen::SparseMatrix< T > &X, const int dim, Eigen::SparseVector< T > &S)
Sum the columns or rows of a sparse matrix.
constexpr double PI
π
Definition PI.h:18