WTFIT
ReebSpace.h
Go to the documentation of this file.
1 
16 #ifndef _REEBSPACE_H
17 #define _REEBSPACE_H
18 
19 // standard includes
20 
21 // base code includes
22 #include <FiberSurface.h>
23 #include <Geometry.h>
24 #include <JacobiSet.h>
25 // to add FiberSurface.h
26 #include <Wrapper.h>
27 
28 #include <map>
29 #include <set>
30 
31 namespace wtfit{
32 
33  class ReebSpace : public Debug{
34 
35  public:
36 
39  rangeArea, // 1
41  };
42 
43  class Sheet0{
44 
45  public:
46 
48  // 0: extrema-sheet, 1: saddle-sheet
49  char type_;
50  vector<int> sheet1List_;
51  vector<int> sheet3List_;
52  };
53 
54  class Sheet1{
55 
56  public:
57 
59  vector<int> edgeList_;
60  // 0: extrema-sheet, 1: saddle-sheet (can be inferred by the number
61  // of 2-sheets attached to it?)
62  vector<int> sheet0List_;
63  // NB: the corresponding 2-sheet should have the same global
64  // indentifier.
65  vector<int> sheet3List_;
66  };
67 
68  class Sheet2{
69 
70  public:
71 
72  bool pruned_;
73  int sheet1Id_;
74  // for each point, x, y, z coordinates
75  // this is how coincident point will be merged afterwards
76  // this list is meant to be temporary. after the merge, we'll use a
77  // global list
78  vector<vector<FiberSurface::Vertex> >
80  vector<vector<FiberSurface::Triangle> >
82  vector<int> sheet3List_;
83  };
84 
85  class Sheet3{
86 
87  public:
88 
90  bool pruned_;
92  vector<int> vertexList_;
93  vector<int> tetList_;
94  vector<int> sheet0List_;
95  vector<int> sheet1List_;
96  vector<int> sheet2List_;
97  vector<int> sheet3List_;
98  vector<int> preMergedSheets_;
99  };
100 
101  ReebSpace();
102 
103  ~ReebSpace();
104 
105  template <class dataTypeU, class dataTypeV>
106  inline int connectivityPreprocessing(
107  const vector<vector<int> > &edgeStarList,
108  vector<vector<pair<int, int> > > &edgeFanLinkEdgeLists,
109  vector<vector<long long int> > &edgeFans,
110  vector<int> &sosOffsets) const;
111 
112  inline bool empty() const{
113  return (currentData_.vertex2sheet0_.size() == 0);
114  }
115 
116  template <class dataTypeU, class dataTypeV> inline int execute();
117 
118  inline const Sheet0* get0sheet(const int &sheetId) const{
119 
120 #ifndef withKamikaze
121  if((sheetId < 0)||(sheetId >= (int) currentData_.sheet0List_.size()))
122  return NULL;
123 #endif
124 
125  return &(currentData_.sheet0List_[sheetId]);
126  }
127 
128  inline const Sheet1* get1sheet(const int &sheetId) const{
129 #ifndef withKamikaze
130  if((sheetId < 0)||(sheetId >= (int) currentData_.sheet1List_.size()))
131  return NULL;
132 #endif
133 
134  return &(currentData_.sheet1List_[sheetId]);
135  }
136 
137  // warning, these are the originals
138  inline const Sheet2* get2sheet(const int &sheetId) const{
139 #ifndef withKamikaze
140  if((sheetId < 0)||(sheetId >= (int) originalData_.sheet2List_.size()))
141  return NULL;
142 #endif
143 
144  return &(originalData_.sheet2List_[sheetId]);
145  }
146 
147  inline const Sheet3* get3sheet(const int &sheetId) const{
148 #ifndef withKamikaze
149  if((sheetId < 0)||(sheetId >= (int) currentData_.sheet3List_.size()))
150  return NULL;
151 #endif
152 
153  return &(currentData_.sheet3List_[sheetId]);
154  }
155 
156  inline const vector<int>* get0sheetSegmentation() const{
158  }
159 
160  const vector<int>* get1sheetSegmentation() const{
161  return &currentData_.edge2sheet1_;
162  }
163 
164  const vector<int>* get3sheetVertexSegmentation() const{
166  }
167 
168  const vector<int>* get3sheetTetSegmentation() const{
169  return &currentData_.tet2sheet3_;
170  }
171 
172  const vector<int>* getEdgeTypes() const{
173  return &currentData_.edgeTypes_;
174  }
175 
176  inline const vector<FiberSurface::Vertex>*
178  return &(fiberSurfaceVertexList_);
179  }
180 
181  inline int getJacobi2Edge(const int &jacobiEdgeId) const{
182  if((jacobiEdgeId < 0)
183  ||(jacobiEdgeId >= (int) jacobi2edges_.size()))
184  return -1;
185  return jacobi2edges_[jacobiEdgeId];
186  }
187 
188  inline const int getNumberOf2sheets() const {
189  return currentData_.sheet2List_.size();
190  }
191 
192 // inline vector<long long int>* getSheetTriangulationCells(){
193 // return &sheet3cells_;
194 // }
195 //
196 // inline vector<float>* getSheetTriangulationPoints(){
197 // return &sheet3points_;
198 // }
199 
200  template <class dataTypeU, class dataTypeV>
201  inline int perturbate(
202  const dataTypeU &uEpsilon = pow(10, -DBL_DIG),
203  const dataTypeV &vEpsilon = pow(10, -DBL_DIG)) const;
204 
205  inline int setEdgeFans(const vector<vector<long long int> > *edgeFans){
206  edgeFans_ = edgeFans;
207  return 0;
208  }
209 
211  const vector<vector<pair<int, int> > > *edgeFanLinkEdgeLists){
212  edgeFanLinkEdgeLists_ = edgeFanLinkEdgeLists;
213  return 0;
214  }
215 
216  inline int setEdgeList(const vector<pair<int, int> > *edgeList){
217  edgeList_ = edgeList;
218  return 0;
219  }
220 
221  inline int setEdgeStars(const vector<vector<int> > *edgeStars){
222  edgeStars_ = edgeStars;
223  return 0;
224  }
225 
226  inline int setExpand3Sheets(const bool &onOff){
227  expand3sheets_ = onOff;
228  return 0;
229  }
230 
231  inline int setInputField(const void *uField, const void *vField){
232 
233  uField_ = uField;
234  vField_ = vField;
235 
236  return 0;
237  }
238 
239  inline int setPointSet(const float *pointSet){
240  pointSet_ = pointSet;
241  return 0;
242  }
243 
244  inline int setSosOffsets(vector<int> *sosOffsets){
245  sosOffsets_ = sosOffsets;
246  return 0;
247  }
248 
249  inline int setTetList(const long long int *tetList){
250  tetList_ = tetList;
251  return 0;
252  }
253 
254  inline int setTetNeighbors(const vector<vector<int> > *tetNeighbors){
255  tetNeighbors_ = tetNeighbors;
256  return 0;
257  }
258 
259  inline int setTetNumber(const int &tetNumber){
260  tetNumber_ = tetNumber;
261  return 0;
262  }
263 
264  inline int setVertexEdgeList(const vector<vector<int> > *vertexEdgeList){
265  vertexEdgeList_ = vertexEdgeList;
266  return 0;
267  }
268 
272  int setVertexNumber(const int &vertexNumber){
273  vertexNumber_ = vertexNumber;
274  return 0;
275  }
276 
277  inline int setVertexStars(const vector<vector<int> > *vertexStars){
278  vertexStars_ = vertexStars;
279  return 0;
280  }
281 
282  template <class dataTypeU, class dataTypeV>
283  inline int simplify(const double &simplificationThreshold,
284  const SimplificationCriterion &criterion = rangeArea);
285 
286 
287  protected:
288 
289  class ReebSpaceData;
290 
291  int compute1sheetsOnly(const vector<pair<int, char> > &jacobiSet,
292  vector<pair<int, int> > &jacobiSetClassification);
293 
294  int compute1sheets(const vector<pair<int, char> > &jacobiSet,
295  vector<pair<int, int> > &jacobiSetClassification);
296 
297  template <class dataTypeU, class dataTypeV>
298  inline int compute2sheets(
299  const vector<pair<int, int> > &jacobiEdges);
300 
301  template <class dataTypeU, class dataTypeV>
302  inline int compute2sheetChambers();
303 
304  int compute3sheet(const int &vertexId,
305  const vector<vector<vector<int> > > &tetTriangles);
306 
307  int compute3sheets(vector<vector<vector<int> > > &tetTriangles);
308 
309  template <class dataTypeU, class dataTypeV>
310  inline int computeGeometricalMeasures(Sheet3 &sheet);
311 
313  ReebSpaceData &data,
314  const int &sheet3Id, const int &sheet0Id);
315 
317  ReebSpaceData &data,
318  const int &sheet3Id, const int &sheet1Id);
319 
321  ReebSpaceData &data,
322  const int &sheet3Id, const int &sheet2Id);
323 
325  ReebSpaceData &data,
326  const int &sheet3Id, const int &otherSheet3Id);
327 
328  int connectSheets();
329 
330  int disconnect1sheetFrom0sheet(ReebSpaceData &data,
331  const int &sheet1Id, const int &sheet0Id, const int &biggerId);
332 
333  int disconnect3sheetFrom0sheet(ReebSpaceData &data,
334  const int &sheet3Id, const int &sheet0Id);
335 
336  int disconnect3sheetFrom1sheet(ReebSpaceData &data,
337  const int &sheet3Id, const int &sheet1Id, const int &biggerId);
338 
339  int disconnect3sheetFrom2sheet(ReebSpaceData &data,
340  const int &sheet3Id, const int &sheet2Id);
341 
342  int disconnect3sheetFrom3sheet(ReebSpaceData &data,
343  const int &sheet3Id, const int &other3SheetId);
344 
345  int flush();
346 
347  int mergeSheets(const int &smallerId, const int &biggerId);
348 
349  int preMergeSheets(const int &sheetId0, const int &sheetId1);
350 
351  int prepareSimplification();
352 
353  int printConnectivity(ostream &stream, const ReebSpaceData &data) const;
354 
355  int simplifySheets(const double &simplificationThreshold,
356  const SimplificationCriterion &simplificationCriterion);
357 
358  int simplifySheet(const int &sheetId,
359  const SimplificationCriterion &simplificationCriterion);
360 
361 // int triangulateTetrahedron(const int &tetId,
362 // const vector<vector<int> > &triangles,
363 // vector<long long int> &outputTets);
364 //
365 // int triangulateThreeSheets();
366 
368 // threeSheetTetNumber_;
370  const void *uField_, *vField_;
371 
372  // input connectivity
373  const vector<pair<int, int> >
375  const vector<vector<int> >
377  const vector<vector<long long int> >
379  const vector<vector<pair<int, int> > >
381  const float *pointSet_;
382  vector<int> *sosOffsets_;
383  const long long int *tetList_;
384  const vector<vector<int> >
386  const vector<vector<int> >
388  const vector<vector<int> >
390 
391  // output segmentation
393 
394  public:
395 
399 
400  vector<int> edge2sheet1_;
401  vector<int> edgeTypes_;
402  vector<int> tet2sheet3_;
403  vector<int> vertex2sheet0_;
404  vector<int> vertex2sheet3_;
405 
406  // structure
407  vector<Sheet0> sheet0List_;
408  vector<Sheet1> sheet1List_;
409  vector<Sheet2> sheet2List_;
410  vector<Sheet3> sheet3List_;
411 
412 // vector<float> sheet3points_;
413 // vector<long long int> sheet3cells_;
414  };
415 
418 
419  // information that does not get simplified
420  vector<pair<int, char> >
422  vector<int> jacobi2edges_;
423 
425  vector<FiberSurface::Vertex>
427  };
428 }
429 
430 // if the package is not a template, comment the following line
431 // #include <ReebSpace.cpp>
432 
433 // template functions
434 template <class dataTypeU, class dataTypeV>
435  inline int ReebSpace::connectivityPreprocessing(
436  const vector<vector<int> > &edgeStarList,
437  vector<vector<pair<int, int> > > &edgeFanLinkEdgeLists,
438  vector<vector<long long int> > &edgeFans,
439  vector<int> &sosOffsets) const{
440 
442 
443  jacobiSet.setWrapper(wrapper_);
444  jacobiSet.setVertexNumber(vertexNumber_);
445  jacobiSet.setTetList(tetList_);
446  jacobiSet.setEdgeList(edgeList_);
447  jacobiSet.connectivityPreprocessing(edgeStarList,
448  edgeFanLinkEdgeLists, edgeFans, sosOffsets);
449 
450  return 0;
451 }
452 
453 template <class dataTypeU, class dataTypeV>
454  inline int ReebSpace::execute(){
455 
456  flush();
457 
458  Timer t;
459 
460  // 1) compute the jacobi set
462 
463  jacobiSet.setWrapper(wrapper_);
464  jacobiSet.setVertexNumber(vertexNumber_);
465  jacobiSet.setInputField(uField_, vField_);
466  jacobiSet.setEdgeList(edgeList_);
468  jacobiSet.setEdgeFans(edgeFans_);
469  jacobiSet.setSosOffsets(sosOffsets_);
470  jacobiSet.execute(jacobiSetEdges_);
471 
472  // 2) compute the list saddle 1-sheets
473  // + list of saddle 0-sheets
474  vector<pair<int, int> > jacobiSetClassification;
476  jacobiSetEdges_, jacobiSetClassification);
477  // at this stage, jacobiSetClassification contains the list of saddle edges
478  // along with their 1-sheet Id.
479 
480  compute2sheets<dataTypeU, dataTypeV>(jacobiSetClassification);
481 // compute2sheetChambers<dataTypeU, dataTypeV>();
482 
483  vector<vector<vector<int> > > tetTriangles;
484  compute3sheets(tetTriangles);
485 
486  {
487  stringstream msg;
488  msg << "[ReebSpace] Data-set (" << vertexNumber_
489  << " points) processed in "
490  << t.getElapsedTime() << " s. TOTAL (" << threadNumber_
491  << " thread(s))."
492  << endl;
493  dMsg(cout, msg.str(), timeMsg);
494  }
495 
496  // post-process for further interaction
497  if((totalArea_ == -1)||(totalVolume_ == -1)||(totalHyperVolume_ == -1)){
498 
499  Timer t;
500 
501 #ifdef withOpenMP
502 #pragma omp parallel for num_threads(threadNumber_)
503 #endif
504  for(int i = 0; i < (int) originalData_.sheet3List_.size(); i++){
505  computeGeometricalMeasures<dataTypeU, dataTypeV>(
507  }
508 
509  for(int i = 0; i < (int) originalData_.sheet3List_.size(); i++){
510  totalArea_ += originalData_.sheet3List_[i].rangeArea_;
511  totalVolume_ += originalData_.sheet3List_[i].domainVolume_;
512  totalHyperVolume_ += originalData_.sheet3List_[i].hyperVolume_;
513  }
514 
515  {
516  stringstream msg;
517  msg << "[ReebSpace] Geometrical measures computed in "
518  << t.getElapsedTime() << " s. (" << threadNumber_
519  << " thread(s))" << endl;
520  dMsg(cout, msg.str(), timeMsg);
521  }
522  }
523 
524  fiberSurface_.finalize<dataTypeU, dataTypeV>();
525 
527 
528  return 0;
529 }
530 
531 template <class dataTypeU, class dataTypeV>
533  const vector<pair<int, int> > &jacobiEdges){
534 
535 #ifndef withKamikaze
536  if(!edgeStars_)
537  return -1;
538  if(!tetList_)
539  return -2;
540  if(!tetNeighbors_)
541  return -3;
542  if(!tetNumber_)
543  return -4;
544  if(!pointSet_)
545  return -5;
546 #endif
547 
548  Timer t;
549 
550  // at this point, they have exactly the same size
552  for(int i = 0; i < (int) originalData_.sheet2List_.size(); i++){
553  originalData_.sheet2List_[i].sheet1Id_ = i;
554  originalData_.sheet2List_[i].pruned_ = false;
555  originalData_.sheet2List_[i].vertexList_.resize(
557  originalData_.sheet2List_[i].sheet1Id_].edgeList_.size());
558  originalData_.sheet2List_[i].triangleList_.resize(
560  originalData_.sheet2List_[i].sheet1Id_].edgeList_.size());
561  for(int j = 0;
562  j < (int) originalData_.sheet2List_[i].vertexList_.size(); j++){
563  originalData_.sheet2List_[i].vertexList_[j].clear();
564  originalData_.sheet2List_[i].triangleList_[j].clear();
565  }
566  }
567 
575  fiberSurface_.setPolygonEdgeNumber(jacobiEdges.size());
576 
577  vector<int> edge2polygonEdgeId(edgeList_->size(), -1);
578  jacobi2edges_.resize(jacobiEdges.size());
579 
580  for(int i = 0; i < (int) jacobiEdges.size(); i++){
581 
582  int edgeId = jacobiEdges[i].first;
583 
584  edge2polygonEdgeId[edgeId] = i;
585  jacobi2edges_[i] = edgeId;
586  }
587 
588 #ifdef withOpenMP
589 #pragma omp parallel for num_threads(threadNumber_)
590 #endif
591  for(int i = 0; i < (int) originalData_.sheet2List_.size(); i++){
592 
593  for(int j = 0;
594  j < (int) originalData_.sheet1List_[
595  originalData_.sheet2List_[i].sheet1Id_].edgeList_.size(); j++){
596 
597  int edgeId = originalData_.sheet1List_[
598  originalData_.sheet2List_[i].sheet1Id_].edgeList_[j];
599 
601  edge2polygonEdgeId[edgeId],
602  &(originalData_.sheet2List_[i].triangleList_[j]));
604  edge2polygonEdgeId[edgeId],
605  &(originalData_.sheet2List_[i].vertexList_[j]));
606  }
607  }
608 
609 #ifdef withOpenMP
610 #pragma omp parallel for num_threads(threadNumber_)
611 #endif
612  for(int i = 0; i < (int) jacobiEdges.size(); i++){
613 
614  int edgeId = jacobiEdges[i].first;
615 
616  pair<double, double> rangePoint0, rangePoint1;
617 
618  rangePoint0.first =
619  ((dataTypeU *) uField_)[(*edgeList_)[edgeId].first];
620  rangePoint0.second =
621  ((dataTypeV *) vField_)[(*edgeList_)[edgeId].first];
622 
623  rangePoint1.first =
624  ((dataTypeU *) uField_)[(*edgeList_)[edgeId].second];
625  rangePoint1.second =
626  ((dataTypeV *) vField_)[(*edgeList_)[edgeId].second];
627 
628  if(originalData_.edgeTypes_[edgeId] == 1){
629  fiberSurface_.computeContour<dataTypeU, dataTypeV>(
630  rangePoint0, rangePoint1,
631  (*edgeStars_)[edgeId],
632  edge2polygonEdgeId[edgeId]);
633  }
634  else{
635  fiberSurface_.computeSurface<dataTypeU, dataTypeV>(
636  rangePoint0, rangePoint1, edge2polygonEdgeId[edgeId]);
637  }
638  }
639 
640 // #ifdef withOpenMP
641 // #pragma omp parallel for num_threads(threadNumber_)
642 // #endif
643 // for(int i = 0; i < (int) originalData_.sheet2List_.size(); i++){
644 //
645 // int sheet1Id = originalData_.sheet2List_[i].sheet1Id_;
646 //
647 // vector<bool> inList(tetNumber_, false);
648 // vector<int> seedList;
649 // vector<pair<pair<double, double>, pair<double, double> > > edgeList;
650 // vector<int> jacobiEdgeIdList;
651 //
652 // for(int j = 0;
653 // j < (int) originalData_.sheet1List_[sheet1Id].edgeList_.size(); j++){
654 // int edgeId = originalData_.sheet1List_[sheet1Id].edgeList_[j];
655 //
656 // if(originalData_.edgeTypes_[edgeId] == 1){
657 // for(int k = 0; k < (int) (*edgeStars_)[edgeId].size(); k++){
658 // int tetId = (*edgeStars_)[edgeId][k];
659 // if(!inList[tetId]){
660 // seedList.push_back(tetId);
661 // }
662 // }
663 // }
664 //
665 // if((originalData_.edgeTypes_[edgeId] != 1)
666 // &&(originalData_.edgeTypes_[edgeId] != -1)){
667 // edgeList.resize(edgeList.size() + 1);
668 //
669 // edgeList.back().first.first =
670 // ((dataTypeU *) uField_)[(*edgeList_)[edgeId].first];
671 // edgeList.back().first.second =
672 // ((dataTypeV *) vField_)[(*edgeList_)[edgeId].first];
673 //
674 // edgeList.back().second.first =
675 // ((dataTypeU *) uField_)[(*edgeList_)[edgeId].second];
676 // edgeList.back().second.second =
677 // ((dataTypeV *) vField_)[(*edgeList_)[edgeId].second];
678 //
679 // jacobiEdgeIdList.push_back(edge2polygonEdgeId[edgeId]);
680 // }
681 // }
682 //
683 // if(seedList.size()){
684 // fiberSurface_.computeContour<dataTypeU, dataTypeV>(
685 // edgeList, seedList,
686 // &jacobiEdgeIdList);
687 // }
688 // }
689 
690  {
691  stringstream msg;
692  msg << "[ReebSpace] Fiber surfaces computed in "
693  << t.getElapsedTime()
694  << " s. overall (" << threadNumber_ << " thread(s))." << endl;
695  dMsg(cout, msg.str(), timeMsg);
696  }
697 
698  return 0;
699 }
700 
701 template <class dataTypeU, class dataTypeV>
703 
704 #ifndef withKamikaze
705  if(!edgeStars_)
706  return -1;
707  if(!tetList_)
708  return -2;
709  if(!tetNeighbors_)
710  return -3;
711  if(!tetNumber_)
712  return -4;
713  if(!pointSet_)
714  return -5;
715 #endif
716 
717  {
718  stringstream msg;
719  msg << "[ReebSpace] Computing chambers' pre-images." << endl;
720  msg << "[ReebSpace] This will take a LONG time." << endl;
721  dMsg(cout, msg.str(), timeMsg);
722  }
723 
724  Timer t;
725 
726  // at this point, they have exactly the same size
728  for(int i = 0; i < (int) originalData_.sheet2List_.size(); i++){
729  originalData_.sheet2List_[i].sheet1Id_ = i;
730  originalData_.sheet2List_[i].pruned_ = false;
731  originalData_.sheet2List_[i].vertexList_.resize(1);
732  originalData_.sheet2List_[i].triangleList_.resize(1);
733  }
734 
743 
744  for(int i = 0; i < (int) originalData_.sheet2List_.size(); i++){
745 
747  &(originalData_.sheet2List_[i].triangleList_[0]));
749  &(originalData_.sheet2List_[i].vertexList_[0]));
750  }
751 
752 #ifdef withOpenMP
753 #pragma omp parallel for num_threads(threadNumber_)
754 #endif
755  for(int i = 0; i < (int) edgeList_->size(); i++){
756 
757  int threadId = omp_get_thread_num();
758 
759  int edgeId = i;
760 
761  pair<double, double> rangePoint0, rangePoint1;
762 
763  rangePoint0.first =
764  ((dataTypeU *) uField_)[(*edgeList_)[edgeId].first];
765  rangePoint0.second =
766  ((dataTypeV *) vField_)[(*edgeList_)[edgeId].first];
767 
768  rangePoint1.first =
769  ((dataTypeU *) uField_)[(*edgeList_)[edgeId].second];
770  rangePoint1.second =
771  ((dataTypeV *) vField_)[(*edgeList_)[edgeId].second];
772 
773  fiberSurface_.computeSurface<dataTypeU, dataTypeV>(
774  rangePoint0, rangePoint1, threadId);
775 
776  // clear the memory now.. otherwise we'll swap and never end
777  originalData_.sheet2List_[threadId].triangleList_[0].clear();
778  originalData_.sheet2List_[threadId].vertexList_[0].clear();
779  }
780 
781  {
782  stringstream msg;
783  msg << "[ReebSpace] Chambers pre-image boundaries computed in "
784  << t.getElapsedTime()
785  << " s. overall (" << threadNumber_ << " thread(s))." << endl;
786  dMsg(cout, msg.str(), timeMsg);
787  }
788 
789  return 0;
790 }
791 
792 template <class dataTypeU, class dataTypeV>
794 
795  sheet.domainVolume_ = 0;
796  sheet.rangeArea_ = 0;
797  sheet.hyperVolume_ = 0;
798 
799  for(int i = 0; i < (int) sheet.tetList_.size(); i++){
800 
801  int tetId = sheet.tetList_[i];
802  vector<pair<double, double> > domainBox, rangeBox;
803  vector<vector<double> > domainPoints(4), rangePoints(4);
804 
805  for(int j = 0; j < 4; j++){
806  domainPoints[j].resize(3);
807  rangePoints[j].resize(2);
808 
809  int vertexId = tetList_[5*tetId + 1 + j];
810  for(int k = 0; k < 3; k++){
811  domainPoints[j][k] = pointSet_[3*vertexId + k];
812  }
813 
814  rangePoints[j][0] = ((dataTypeU *) uField_)[vertexId];
815  rangePoints[j][1] = ((dataTypeV *) vField_)[vertexId];
816  }
817 
818  Geometry::getBoundingBox(domainPoints, domainBox);
819  Geometry::getBoundingBox(rangePoints, rangeBox);
820 
821  sheet.domainVolume_ +=
822  (domainBox[0].second - domainBox[0].first)
823  *(domainBox[1].second - domainBox[1].first)
824  *(domainBox[2].second - domainBox[2].first);
825 
826  sheet.rangeArea_ +=
827  (rangeBox[0].second - rangeBox[0].first)
828  *(rangeBox[1].second - rangeBox[1].first);
829 
830  }
831 
832  if(sheet.domainVolume_){
833  sheet.hyperVolume_ = sheet.rangeArea_/sheet.domainVolume_;
834  }
835  else{
836  sheet.hyperVolume_ = 0;
837  }
838 
839  return 0;
840 }
841 
842 template <class dataTypeU, class dataTypeV>
844  const dataTypeU &uEpsilon, const dataTypeV &vEpsilon) const{
845 
847  jacobiSet.setWrapper(wrapper_);
848  jacobiSet.setInputField(uField_, vField_);
849  jacobiSet.setVertexNumber(vertexNumber_);
850  jacobiSet.perturbate(uEpsilon, vEpsilon);
851 
852  return 0;
853 }
854 
855 template <class dataTypeU, class dataTypeV>
856  inline int ReebSpace::simplify(const double &simplificationThreshold,
857  const SimplificationCriterion &criterion){
858 
859  if((totalArea_ == -1)||(totalVolume_ == -1)||(totalHyperVolume_ == -1)){
860 
861  Timer t;
862 
863 #ifdef withOpenMP
864 #pragma omp parallel for num_threads(threadNumber_)
865 #endif
866  for(int i = 0; i < (int) originalData_.sheet3List_.size(); i++){
867  computeGeometricalMeasures<dataTypeU, dataTypeV>(
869  }
870 
871  for(int i = 0; i < (int) originalData_.sheet3List_.size(); i++){
872  totalArea_ += originalData_.sheet3List_[i].rangeArea_;
873  totalVolume_ += originalData_.sheet3List_[i].domainVolume_;
874  totalHyperVolume_ += originalData_.sheet3List_[i].hyperVolume_;
875  }
876 
877  {
878  stringstream msg;
879  msg << "[ReebSpace] Geometrical measures computed in "
880  << t.getElapsedTime() << " s. (" << threadNumber_
881  << " thread(s))" << endl;
882  dMsg(cout, msg.str(), timeMsg);
883  }
884  }
885 
886  if(!hasConnectedSheets_){
887  connectSheets();
889  }
890 
891  {
892  stringstream msg;
893  msg << "[ReebSpace] Simplifying with criterion ";
894  switch(criterion){
895  case domainVolume:
896  msg << "'Domain Volume'";
897  break;
898  case rangeArea:
899  msg << "'Range Area'";
900  break;
901  case hyperVolume:
902  msg << "'HyperVolume'";
903  break;
904  }
905  msg << " at threshold " << simplificationThreshold << "." << endl;
906  dMsg(cout, msg.str(), timeMsg);
907  }
908 
909  if(!((criterion == currentData_.simplificationCriterion_)
910  &&(simplificationThreshold > currentData_.simplificationThreshold_))){
912  }
913 
914  simplifySheets(simplificationThreshold, criterion);
915 
916  return 0;
917 }
918 
919 #endif // REEBSPACE_H
vector< int > vertex2sheet3_
Definition: ReebSpace.h:404
int setTetList(const long long int *tetList)
Definition: FiberSurface.h:152
const vector< vector< int > > * edgeStars_
Definition: ReebSpace.h:376
int simplifySheets(const double &simplificationThreshold, const SimplificationCriterion &simplificationCriterion)
Definition: ReebSpace.cpp:1427
int setEdgeList(const vector< pair< int, int > > *edgeList)
Definition: JacobiSet.h:57
~ReebSpace()
Definition: ReebSpace.cpp:30
int compute2sheets(const vector< pair< int, int > > &jacobiEdges)
Definition: ReebSpace.h:532
int setVertexEdgeList(const vector< vector< int > > *vertexEdgeList)
Definition: ReebSpace.h:264
vector< int > sheet0List_
Definition: ReebSpace.h:94
int perturbate(const dataTypeU &uEpsilon=pow(10,-DBL_DIG), const dataTypeV &vEpsilon=pow(10,-DBL_DIG)) const
Definition: JacobiSet.cpp:386
int connect3sheetTo3sheet(ReebSpaceData &data, const int &sheet3Id, const int &otherSheet3Id)
Definition: ReebSpace.cpp:785
int preMergeSheets(const int &sheetId0, const int &sheetId1)
Definition: ReebSpace.cpp:1149
Definition: ReebSpace.h:54
double totalVolume_
Definition: ReebSpace.h:369
vector< int > sheet3List_
Definition: ReebSpace.h:97
int finalize(const bool &mergeDuplicatedVertices=false, const bool &removeSmallEdges=false, const bool &edgeFlips=false, const bool &intersectionRemesh=false)
Definition: FiberSurface.h:1589
SimplificationCriterion simplificationCriterion_
Definition: ReebSpace.h:397
bool empty() const
Definition: ReebSpace.h:112
const long long int * tetList_
Definition: ReebSpace.h:383
int setPolygonEdgeNumber(const int &polygonEdgeNumber)
Definition: FiberSurface.h:145
const int msg(const char *msg, const int &debugLevel=infoMsg) const
Definition: Debug.cpp:67
int mergeSheets(const int &smallerId, const int &biggerId)
Definition: ReebSpace.cpp:1046
const Sheet2 * get2sheet(const int &sheetId) const
Definition: ReebSpace.h:138
int Id_
Definition: ReebSpace.h:89
int setEdgeFanLinkEdgeList(const vector< vector< pair< int, int > > > *edgeFanLinkEdgeLists)
Definition: JacobiSet.h:51
int connectSheets()
Definition: ReebSpace.cpp:813
ReebSpaceData currentData_
Definition: ReebSpace.h:417
bool pruned_
Definition: ReebSpace.h:90
int disconnect3sheetFrom3sheet(ReebSpaceData &data, const int &sheet3Id, const int &other3SheetId)
Definition: ReebSpace.cpp:997
int setTetNumber(const int &tetNumber)
Definition: FiberSurface.h:162
const float * pointSet_
Definition: ReebSpace.h:381
bool pruned_
Definition: ReebSpace.h:72
int simplifySheet(const int &sheetId, const SimplificationCriterion &simplificationCriterion)
Definition: ReebSpace.cpp:1304
const int getNumberOf2sheets() const
Definition: ReebSpace.h:188
vector< Sheet3 > sheet3List_
Definition: ReebSpace.h:410
int connectivityPreprocessing(const vector< vector< int > > &edgeStarList, vector< vector< pair< int, int > > > &edgeFanLinkEdgeLists, vector< vector< long long int > > &edgeFans, vector< int > &sosOffsets) const
Definition: JacobiSet.cpp:26
static int getBoundingBox(const vector< vector< double > > &points, vector< pair< double, double >> &bBox)
Definition: Geometry.cpp:387
double totalHyperVolume_
Definition: ReebSpace.h:369
int setEdgeFanLinkEdgeList(const vector< vector< pair< int, int > > > *edgeFanLinkEdgeLists)
Definition: ReebSpace.h:210
const vector< int > * get3sheetTetSegmentation() const
Definition: ReebSpace.h:168
int computeContour(const pair< double, double > &rangePoint0, const pair< double, double > &rangePoint1, const vector< int > &seedTetList, const int &polygonEdgeId=0) const
Definition: FiberSurface.h:1298
int setVertexNumber(const int &vertexNumber)
Definition: ReebSpace.h:272
Definition: ReebSpace.h:68
const vector< vector< int > > * vertexStars_
Definition: ReebSpace.h:389
vector< int > sheet1List_
Definition: ReebSpace.h:95
int setEdgeList(const vector< pair< int, int > > *edgeList)
Definition: ReebSpace.h:216
int printConnectivity(ostream &stream, const ReebSpaceData &data) const
Definition: ReebSpace.cpp:1240
ReebSpace processing package.
Definition: ReebSpace.h:33
double simplificationThreshold_
Definition: ReebSpace.h:398
vector< Sheet1 > sheet1List_
Definition: ReebSpace.h:408
int execute(vector< pair< int, char > > &jacobiSet)
Definition: JacobiSet.cpp:186
int flush()
Definition: ReebSpace.cpp:1014
const vector< vector< pair< int, int > > > * edgeFanLinkEdgeLists_
Definition: ReebSpace.h:380
FiberSurface fiberSurface_
Definition: ReebSpace.h:424
int disconnect3sheetFrom0sheet(ReebSpaceData &data, const int &sheet3Id, const int &sheet0Id)
Definition: ReebSpace.cpp:915
JacobiSet processing package.
Definition: JacobiSet.h:28
ReebSpace()
Definition: ReebSpace.cpp:5
const Sheet0 * get0sheet(const int &sheetId) const
Definition: ReebSpace.h:118
int setEdgeStars(const vector< vector< int > > *edgeStars)
Definition: ReebSpace.h:221
int disconnect3sheetFrom1sheet(ReebSpaceData &data, const int &sheet3Id, const int &sheet1Id, const int &biggerId)
Definition: ReebSpace.cpp:932
virtual const int dMsg(ostream &stream, string msg, const int &debugLevel=infoMsg) const
Definition: Debug.cpp:52
int setVertexNumber(const int &vertexNumber)
Definition: JacobiSet.h:81
int perturbate(const dataTypeU &uEpsilon=pow(10,-DBL_DIG), const dataTypeV &vEpsilon=pow(10,-DBL_DIG)) const
Definition: ReebSpace.h:843
const vector< int > * get0sheetSegmentation() const
Definition: ReebSpace.h:156
SimplificationCriterion
Definition: ReebSpace.h:37
int setInputField(const void *uField, const void *vField)
Definition: ReebSpace.h:231
int getJacobi2Edge(const int &jacobiEdgeId) const
Definition: ReebSpace.h:181
Definition: ReebSpace.h:40
const vector< vector< long long int > > * edgeFans_
Definition: ReebSpace.h:378
int setTetList(const long long int *tetList)
Definition: JacobiSet.h:76
int setEdgeFans(const vector< vector< long long int > > *edgeFans)
Definition: JacobiSet.h:46
vector< int > vertexList_
Definition: ReebSpace.h:92
int compute3sheet(const int &vertexId, const vector< vector< vector< int > > > &tetTriangles)
Definition: ReebSpace.cpp:354
int setTetNeighbors(const vector< vector< int > > *tetNeighbors)
Definition: FiberSurface.h:157
Definition: Debug.h:50
double getElapsedTime()
Definition: Os.h:359
vector< vector< FiberSurface::Triangle > > triangleList_
Definition: ReebSpace.h:81
FiberSurface processing package.
Definition: FiberSurface.h:31
vector< int > edgeTypes_
Definition: ReebSpace.h:401
Definition: Os.h:347
vector< int > sheet0List_
Definition: ReebSpace.h:62
vector< int > tet2sheet3_
Definition: ReebSpace.h:402
bool hasConnectedSheets_
Definition: ReebSpace.h:416
int setTriangleList(const int &polygonEdgeId, vector< Triangle > *triangleList)
Definition: FiberSurface.h:167
int setWrapper(const Wrapper *wrapper)
Definition: Debug.cpp:76
int connect3sheetTo0sheet(ReebSpaceData &data, const int &sheet3Id, const int &sheet0Id)
Definition: ReebSpace.cpp:712
int setTetNumber(const int &tetNumber)
Definition: ReebSpace.h:259
double hyperVolume_
Definition: ReebSpace.h:91
vector< int > * sosOffsets_
Definition: ReebSpace.h:382
vector< int > sheet3List_
Definition: ReebSpace.h:51
vector< int > preMergedSheets_
Definition: ReebSpace.h:98
int compute1sheets(const vector< pair< int, char > > &jacobiSet, vector< pair< int, int > > &jacobiSetClassification)
Definition: ReebSpace.cpp:161
vector< int > sheet2List_
Definition: ReebSpace.h:96
bool pruned_
Definition: ReebSpace.h:58
Definition: ReebSpace.h:39
int compute1sheetsOnly(const vector< pair< int, char > > &jacobiSet, vector< pair< int, int > > &jacobiSetClassification)
Definition: ReebSpace.cpp:34
const vector< int > * get1sheetSegmentation() const
Definition: ReebSpace.h:160
int setGlobalVertexList(vector< Vertex > *globalList)
Definition: FiberSurface.h:106
int setEdgeFans(const vector< vector< long long int > > *edgeFans)
Definition: ReebSpace.h:205
int simplificationId_
Definition: ReebSpace.h:89
vector< int > sheet3List_
Definition: ReebSpace.h:82
int sheet1Id_
Definition: ReebSpace.h:73
int setPointSet(const float *pointSet)
Definition: ReebSpace.h:239
int setInputField(const void *uField, const void *vField)
Definition: FiberSurface.h:111
int setPointSet(const float *pointSet)
Definition: FiberSurface.h:134
int vertexNumber_
Definition: ReebSpace.h:367
vector< int > vertex2sheet0_
Definition: ReebSpace.h:403
int prepareSimplification()
Definition: ReebSpace.cpp:1180
Definition: ReebSpace.h:85
const vector< pair< int, int > > * edgeList_
Definition: ReebSpace.h:374
int pruned_
Definition: ReebSpace.h:47
int compute2sheetChambers()
Definition: ReebSpace.h:702
double totalArea_
Definition: ReebSpace.h:369
vector< int > sheet1List_
Definition: ReebSpace.h:50
bool hasSaddleEdges_
Definition: ReebSpace.h:58
ReebSpaceData originalData_
Definition: ReebSpace.h:417
int connect3sheetTo2sheet(ReebSpaceData &data, const int &sheet3Id, const int &sheet2Id)
Definition: ReebSpace.cpp:761
int preMerger_
Definition: ReebSpace.h:89
vector< int > edge2sheet1_
Definition: ReebSpace.h:400
vector< Sheet2 > sheet2List_
Definition: ReebSpace.h:409
vector< int > tetList_
Definition: ReebSpace.h:93
int vertexId_
Definition: ReebSpace.h:47
Minimalist debugging class.
Definition: Debug.h:39
bool expand3sheets_
Definition: ReebSpace.h:416
const vector< FiberSurface::Vertex > * getFiberSurfaceVertices() const
Definition: ReebSpace.h:177
char type_
Definition: ReebSpace.h:49
int tetNumber_
Definition: ReebSpace.h:367
Definition: ReebSpace.h:43
int setExpand3Sheets(const bool &onOff)
Definition: ReebSpace.h:226
vector< vector< FiberSurface::Vertex > > vertexList_
Definition: ReebSpace.h:79
int disconnect3sheetFrom2sheet(ReebSpaceData &data, const int &sheet3Id, const int &sheet2Id)
Definition: ReebSpace.cpp:980
int setVertexList(const int &polygonEdgeId, vector< Vertex > *vertexList)
Definition: FiberSurface.h:180
const Sheet3 * get3sheet(const int &sheetId) const
Definition: ReebSpace.h:147
int setSosOffsets(vector< int > *sosOffsets)
Definition: JacobiSet.h:69
const void * vField_
Definition: ReebSpace.h:370
int execute()
Definition: ReebSpace.h:454
int setInputField(const void *uField, const void *vField)
Definition: JacobiSet.h:62
vector< int > edgeList_
Definition: ReebSpace.h:59
Definition: CommandLineParser.h:13
Wrapper * wrapper_
Definition: Debug.h:109
const Sheet1 * get1sheet(const int &sheetId) const
Definition: ReebSpace.h:128
const vector< int > * getEdgeTypes() const
Definition: ReebSpace.h:172
int computeSurface(const pair< double, double > &rangePoint0, const pair< double, double > &rangePoint1, const int &polygonEdgeId=0) const
Definition: FiberSurface.h:1443
vector< int > jacobi2edges_
Definition: ReebSpace.h:422
int setTetNeighbors(const vector< vector< int > > *tetNeighbors)
Definition: ReebSpace.h:254
const void * uField_
Definition: ReebSpace.h:370
vector< pair< int, char > > jacobiSetEdges_
Definition: ReebSpace.h:421
const vector< vector< int > > * vertexEdgeList_
Definition: ReebSpace.h:387
vector< int > sheet3List_
Definition: ReebSpace.h:65
double domainVolume_
Definition: ReebSpace.h:91
Definition: ReebSpace.h:392
int setVertexStars(const vector< vector< int > > *vertexStars)
Definition: ReebSpace.h:277
int connectivityPreprocessing(const vector< vector< int > > &edgeStarList, vector< vector< pair< int, int > > > &edgeFanLinkEdgeLists, vector< vector< long long int > > &edgeFans, vector< int > &sosOffsets) const
Definition: ReebSpace.h:435
int setSosOffsets(vector< int > *sosOffsets)
Definition: ReebSpace.h:244
Definition: ReebSpace.h:38
double rangeArea_
Definition: ReebSpace.h:91
const vector< vector< int > > * tetNeighbors_
Definition: ReebSpace.h:385
int computeGeometricalMeasures(Sheet3 &sheet)
Definition: ReebSpace.h:793
int disconnect1sheetFrom0sheet(ReebSpaceData &data, const int &sheet1Id, const int &sheet0Id, const int &biggerId)
Definition: ReebSpace.cpp:892
int compute3sheets(vector< vector< vector< int > > > &tetTriangles)
Definition: ReebSpace.cpp:453
const vector< int > * get3sheetVertexSegmentation() const
Definition: ReebSpace.h:164
int simplify(const double &simplificationThreshold, const SimplificationCriterion &criterion=rangeArea)
Definition: ReebSpace.h:856
int setTetList(const long long int *tetList)
Definition: ReebSpace.h:249
int connect3sheetTo1sheet(ReebSpaceData &data, const int &sheet3Id, const int &sheet1Id)
Definition: ReebSpace.cpp:737
vector< FiberSurface::Vertex > fiberSurfaceVertexList_
Definition: ReebSpace.h:426
int threadNumber_
Definition: Debug.h:108
vector< Sheet0 > sheet0List_
Definition: ReebSpace.h:407