WTFIT
FiberSurface.h
Go to the documentation of this file.
1 
15 #ifndef _FIBERSURFACE_H
16 #define _FIBERSURFACE_H
17 
18 // standard includes
19 #include <parallel/algorithm>
20 #include <queue>
21 
22 // base code includes
23 #include <Geometry.h>
24 #ifdef withrangeDrivenOctree
25 #include <RangeDrivenOctree.h>
26 #endif
27 #include <Wrapper.h>
28 
29 namespace wtfit{
30 
31  class FiberSurface : public Debug{
32 
33  public:
34 
35  class Vertex{
36 
37  public:
40  // TODO also encode the vertex ids of the triangle of the input mesh
41  // where this point has been computed (for constrained triangulation)
42  pair<int, int> meshEdge_;
43  double p_[3], t_;
44  pair<double, double>
45  uv_;
46  };
47 
48  class Triangle{
49 
50  public:
51 
53  };
54 
55  FiberSurface();
56 
57  ~FiberSurface();
58 
59 #ifdef withrangeDrivenOctree
60  template <class dataTypeU, class dataTypeV>
61  inline int buildOctree();
62 #endif
63 
64  template <class dataTypeU, class dataTypeV>
65  inline int computeContour(
66  const pair<double, double> &rangePoint0,
67  const pair<double, double> &rangePoint1,
68  const vector<int> &seedTetList,
69  const int &polygonEdgeId = 0) const;
70 
71  template <class dataTypeU, class dataTypeV>
72  inline int computeContour(
73  const vector<pair<pair<double, double>, pair<double, double> > >
74  &edgeList,
75  const vector<int> &seedTetList,
76  const vector<int> *edgeIdList = NULL) const;
77 
78  template <class dataTypeU, class dataTypeV>
79  inline int computeSurface(
80  const pair<double, double> &rangePoint0,
81  const pair<double, double> &rangePoint1,
82  const int &polygonEdgeId = 0) const;
83 
84  template <class dataTypeU, class dataTypeV>
85  inline int computeSurface();
86 
87  template <class dataTypeU, class dataTypeV>
88  inline int finalize(
89  const bool &mergeDuplicatedVertices = false,
90  const bool &removeSmallEdges = false,
91  const bool &edgeFlips = false,
92  const bool &intersectionRemesh = false);
93 
94 #ifdef withrangeDrivenOctree
95  inline int flushOctree(){
96  return octree_.flush();
97  }
98 #endif
99 
100  template <class dataTypeU, class dataTypeV>
101  inline int processTetrahedron(const int &tetId,
102  const pair<double, double> &rangePoint0,
103  const pair<double, double> &rangePoint1,
104  const int &polygonEdgeId = 0) const;
105 
106  inline int setGlobalVertexList(vector<Vertex> *globalList){
107  globalVertexList_ = globalList;
108  return 0;
109  }
110 
111  inline int setInputField(const void *uField, const void *vField){
112 
113  uField_ = uField;
114  vField_ = vField;
115 
116  return 0;
117  }
118 
119  inline int setPointMerging(const bool &onOff){
120  pointSnapping_ = onOff;
121  return 0;
122  }
123 
124  inline int setPointMergingThreshold(const double &threshold){
125  pointSnappingThreshold_ = threshold;
126  return 0;
127  }
128 
129  inline int setPointNumber(const int &number){
130  pointNumber_ = number;
131  return 0;
132  }
133 
134  inline int setPointSet(const float *pointSet){
135  pointSet_ = pointSet;
136  return 0;
137  }
138 
139  inline int setPolygon(const vector<pair<pair<double, double>,
140  pair<double, double> > > *polygon){
141  polygon_ = polygon;
142  return 0;
143  }
144 
145  inline int setPolygonEdgeNumber(const int &polygonEdgeNumber){
146  polygonEdgeNumber_ = polygonEdgeNumber;
147  polygonEdgeVertexLists_.resize(polygonEdgeNumber, NULL);
148  polygonEdgeTriangleLists_.resize(polygonEdgeNumber, NULL);
149  return 0;
150  }
151 
152  inline int setTetList(const long long int *tetList){
153  tetList_ = tetList;
154  return 0;
155  }
156 
157  inline int setTetNeighbors(const vector<vector<int> > *tetNeighbors){
158  tetNeighbors_ = tetNeighbors;
159  return 0;
160  }
161 
162  inline int setTetNumber(const int &tetNumber){
163  tetNumber_ = tetNumber;
164  return 0;
165  }
166 
167  inline int setTriangleList(const int &polygonEdgeId,
168  vector<Triangle> *triangleList){
169 
170 #ifndef withKamikaze
171  if((polygonEdgeId >= 0)
172  &&(polygonEdgeId < (int) polygonEdgeTriangleLists_.size()))
173 #endif
174  polygonEdgeTriangleLists_[polygonEdgeId] = triangleList;
175 
176 
177  return 0;
178  }
179 
180  inline int setVertexList(const int &polygonEdgeId,
181  vector<Vertex> *vertexList){
182 
183 #ifndef withKamikaze
184  if((polygonEdgeId >= 0)
185  &&(polygonEdgeId < (int) polygonEdgeVertexLists_.size()))
186 #endif
187  polygonEdgeVertexLists_[polygonEdgeId] = vertexList;
188 
189  return 0;
190  }
191 
192  protected:
193 
194  typedef struct _intersectionTriangle{
195  int caseId_;
196  // use negative values for new triangles
199  // use negative values for new vertices
200  int vertexIds_[3];
201  pair<double, double> uv_[3];
202  double t_[3];
203  double p_[3][3];
204  pair<double, double> intersection_;
206 
207  template <class dataTypeU, class dataTypeV>
208  inline int computeBaseTriangle(
209  const int &tetId,
210  const int &localEdgeId0, const double &t0,
211  const double &u0, const double &v0,
212  const int &localEdgeId1, const double &t1,
213  const double &u1, const double &v1,
214  const int &localEdgeId2, const double &t2,
215  const double &u2, const double &v2,
216  vector<vector<double> > &basePoints,
217  vector<pair<double, double> > &basePointProections,
218  vector<double> &basePointParameterization,
219  vector<pair<int, int> > &baseEdges) const;
220 
221  template <class dataTypeU, class dataTYpeV>
222  inline int computeCase0(const int &polygonEdgeId, const int &tetId,
223  const int &localEdgeId0, const double &t0,
224  const double &u0, const double &v0,
225  const int &localEdgeId1, const double &t1,
226  const double &u1, const double &v1,
227  const int &localEdgeId2, const double &t2,
228  const double &u2, const double &v2) const;
229 
230  template <class dataTypeU, class dataTYpeV>
231  inline int computeCase1(const int &polygonEdgeId, const int &tetId,
232  const int &localEdgeId0, const double &t0,
233  const double &u0, const double &v0,
234  const int &localEdgeId1, const double &t1,
235  const double &u1, const double &v1,
236  const int &localEdgeId2, const double &t2,
237  const double &u2, const double &v2) const;
238 
239  template <class dataTypeU, class dataTYpeV>
240  inline int computeCase2(const int &polygonEdgeId, const int &tetId,
241  const int &localEdgeId0, const double &t0,
242  const double &u0, const double &v0,
243  const int &localEdgeId1, const double &t1,
244  const double &u1, const double &v1,
245  const int &localEdgeId2, const double &t2,
246  const double &u2, const double &v2) const;
247 
248  template <class dataTypeU, class dataTYpeV>
249  inline int computeCase3(const int &polygonEdgeId, const int &tetId,
250  const int &localEdgeId0, const double &t0,
251  const double &u0, const double &v0,
252  const int &localEdgeId1, const double &t1,
253  const double &u1, const double &v1,
254  const int &localEdgeId2, const double &t2,
255  const double &u2, const double &v2) const;
256 
257  template <class dataTypeU, class dataTYpeV>
258  inline int computeCase4(const int &polygonEdgeId, const int &tetId,
259  const int &localEdgeId0, const double &t0,
260  const double &u0, const double &v0,
261  const int &localEdgeId1, const double &t1,
262  const double &u1, const double &v1,
263  const int &localEdgeId2, const double &t2,
264  const double &u2, const double &v2) const;
265 
266 #ifdef withrangeDrivenOctree
267  template <class dataTypeU, class dataTypeV>
268  inline int computeSurfaceWithOctree(
269  const pair<double, double> &rangePoint0,
270  const pair<double, double> &rangePoint1,
271  const int &polygonEdgeId) const;
272 #endif
273 
274  int computeTriangleFiber(const int &tetId, const int &triangleId,
275  const pair<double, double> &intersection,
276  const vector<vector<IntersectionTriangle> > &tetIntersections,
277  vector<double> &pA, vector<double> &pB, int &pivotVertexId,
278  bool &edgeFiber) const;
279 
281  const int &tetId, const int &triangleId0, const int &triangleId1,
282  const int &polygonEdgeId0, const int &polygonEdgeId1,
283  const pair<double, double> &intersection,
284  int &newVertexNumber, int &newTriangleNumber,
285  vector<vector<IntersectionTriangle> > &tetIntersections,
286  vector<vector<Vertex> > &tetNewVertices) const;
287 
289  const int &tetId, const int &triangleId, const int &polygonEdgeId,
290  const pair<double, double> &intersection,
291  const vector<double> &pA, const vector<double> &pB,
292  const int &pivotVertexId,
293  int &newVertexNumber, int &newTriangleNumber,
294  vector<vector<IntersectionTriangle> > &tetIntersections,
295  vector<vector<Vertex> > &tetNewVertices) const;
296 
298  const int &tetId, const int &triangleId,
299  const int &vertexId0, const int &vertexId1, const int &vertexId2,
300  const vector<vector<Vertex> > &tetNewVertices,
301  int &newTriangleNumber,
302  vector<vector<IntersectionTriangle> > &tetIntersections,
303  const pair<double, double> *intersection = NULL) const;
304 
305  int flipEdges() const;
306 
307  int flipEdges(vector<pair<int, int> > &triangles) const;
308 
309  int getNumberOfCommonVertices(const int &tetId,
310  const int &triangleId0, const int &triangleId1,
311  const vector<vector<IntersectionTriangle> > &tetIntersections) const;
312 
313  int getTriangleRangeExtremities(const int &tetId, const int &triangleId,
314  const vector<vector<IntersectionTriangle> > &tetIntersections,
315  pair<double, double> &extremity0,
316  pair<double, double> &extremity1) const;
317 
319  const double *p0, const double *p1, const double *p2) const;
320 
322  const vector<double> &p0,
323  const pair<double, double> &uv0,
324  const double &t0,
325  const vector<double> &p1,
326  const pair<double, double> &uv1,
327  const double &t1,
328  const double &t,
329  Vertex &v) const;
330 
332  const int &source, const int &destination,
333  const int &pivotVertexId,
334  const vector<pair<int, int> > &starNeighbors) const;
335 
336  bool isEdgeFlippable(
337  const int &edgeVertexId0, const int &edgeVertexId1,
338  const int &otherVertexId0, const int &otherVertexId1) const;
339 
341  const int &tetId, const int &triangleId,
342  const vector<vector<IntersectionTriangle> > &tetIntersections,
343  const vector<vector<Vertex> > &tetNewVertices,
344  const int &vertexId0, const int &vertexId1, const int &vertexId2) const{
345 
346  vector<vector<double> > points(3);
347  for(int i = 0; i < 3; i++){
348  int vertexId = vertexId0;
349  if(i == 1) vertexId = vertexId1;
350  if(i == 2) vertexId = vertexId2;
351 
352  points[i].resize(3);
353  if(vertexId >= 0){
354  for(int j = 0; j < 3; j++){
355  points[i][j] =
356  tetIntersections[tetId][triangleId].p_[vertexId][j];
357  }
358  }
359  else{
360  for(int j = 0; j < 3; j++){
361  points[i][j] = tetNewVertices[tetId][(-vertexId) - 1].p_[j];
362  }
363  }
364  }
365 
367  points[0].data(), points[1].data(), points[2].data());
368  }
369 
370  int mergeEdges(const double &distanceThreshold) const;
371 
372  int mergeVertices(const double &distanceThreshold) const;
373 
374  template <class dataTypeU, class dataTypeV>
375  inline int remeshIntersections() const;
376 
377  int snapToBasePoint(const vector<vector<double> > &basePoints,
378  const vector<pair<double, double> > &uv,
379  const vector<double> &t,
380  Vertex &v) const;
381 
382  int snapVertexBarycentrics(const double &distanceThreshold) const;
383 
384  int snapVertexBarycentrics(const int &tetId,
385  const vector<pair<int, int> > &triangles,
386  const double &distanceThreshold) const;
387 
388 
390 
392  const void *uField_, *vField_;
393  const float *pointSet_;
394  const long long int *tetList_;
395  const vector<vector<int> >
398 
400 
401  const vector<pair<pair<double, double>, pair<double, double> > >
403 
404  vector<Vertex> *globalVertexList_;
405  vector<vector<Vertex> *>
407  vector<vector<Triangle> *>
409 
410 #ifdef withrangeDrivenOctree
411  RangeDrivenOctree octree_;
412 #endif
413  };
414 }
415 
416 #ifdef withrangeDrivenOctree
417 template <class dataTypeU, class dataTypeV>
418  inline int FiberSurface::buildOctree(){
419 
420  if(octree_.empty()){
421 
422  octree_.setDebugLevel(debugLevel_);
423  octree_.setThreadNumber(threadNumber_);
424  octree_.setCellList(tetList_);
425  octree_.setCellNumber(tetNumber_);
426  octree_.setPointList(pointSet_);
427  octree_.setVertexNumber(pointNumber_);
428  octree_.setRange(uField_, vField_);
429 
430  octree_.build<dataTypeU, dataTypeV>();
431  }
432 
433  return 0;
434 }
435 #endif
436 
437 template <class dataTypeU, class dataTypeV>
438  inline int FiberSurface::computeBaseTriangle(const int &tetId,
439  const int &localEdgeId0, const double &t0,
440  const double &u0, const double &v0,
441  const int &localEdgeId1, const double &t1,
442  const double &u1, const double &v1,
443  const int &localEdgeId2, const double &t2,
444  const double &u2, const double &v2,
445  vector<vector<double> > &basePoints,
446  vector<pair<double, double> > &basePointProjections,
447  vector<double> &basePointParameterization,
448  vector<pair<int, int> > &baseEdges) const{
449 
450  basePoints.resize(3);
451  basePointProjections.resize(3);
452  basePointParameterization.resize(3);
453  baseEdges.resize(3);
454 
455  for(int i = 0; i < 3; i++){
456 
457  int vertexId0, vertexId1;
458 
459  switch(i){
460 
461  case 0:
462  vertexId0 = tetList_[5*tetId + 1
463  + edgeImplicitEncoding_[2*localEdgeId0]];
464  vertexId1 = tetList_[5*tetId + 1
465  + edgeImplicitEncoding_[2*localEdgeId0 + 1]];
466  basePointProjections[i].first = u0;
467  basePointProjections[i].second = v0;
468  basePointParameterization[i] = t0;
469  break;
470 
471  case 1:
472  vertexId0 = tetList_[5*tetId + 1
473  + edgeImplicitEncoding_[2*localEdgeId1]];
474  vertexId1 = tetList_[5*tetId + 1
475  + edgeImplicitEncoding_[2*localEdgeId1 + 1]];
476  basePointProjections[i].first = u1;
477  basePointProjections[i].second = v1;
478  basePointParameterization[i] = t1;
479  break;
480 
481  case 2:
482  vertexId0 = tetList_[5*tetId + 1
483  + edgeImplicitEncoding_[2*localEdgeId2]];
484  vertexId1 = tetList_[5*tetId + 1
485  + edgeImplicitEncoding_[2*localEdgeId2 + 1]];
486  basePointProjections[i].first = u2;
487  basePointProjections[i].second = v2;
488  basePointParameterization[i] = t2;
489  break;
490  }
491 
492  vector<double> baryCentrics;
493  vector<double> p0(2), p1(2), p(2);
494  p0[0] = ((dataTypeU *) uField_)[vertexId0];
495  p0[1] = ((dataTypeV *) vField_)[vertexId0];
496  p1[0] = ((dataTypeU *) uField_)[vertexId1];
497  p1[1] = ((dataTypeV *) vField_)[vertexId1];
498  p[0] = basePointProjections[i].first;
499  p[1] = basePointProjections[i].second;
501  p0.data(), p1.data(), p.data(), baryCentrics, 2);
502 
503  basePoints[i].resize(3);
504  for(int j = 0; j < 3; j++){
505  basePoints[i][j] =
506  baryCentrics[0]*pointSet_[3*vertexId0 + j]
507  + baryCentrics[1]*pointSet_[3*vertexId1 + j];
508  }
509 
510  if(vertexId0 < vertexId1){
511  baseEdges[i] = pair<int, int>(vertexId0, vertexId1);
512  }
513  else{
514  baseEdges[i] = pair<int, int>(vertexId1, vertexId0);
515  }
516  }
517 
518  return 0;
519 }
520 
521 template<class dataTypeU, class dataTypeV>
523  const int &polygonEdgeId, const int &tetId,
524  const int &localEdgeId0, const double &t0,
525  const double &u0, const double &v0,
526  const int &localEdgeId1, const double &t1,
527  const double &u1, const double &v1,
528  const int &localEdgeId2, const double &t2,
529  const double &u2, const double &v2) const{
530 
531  // that one's easy, make just one triangle
532  int vertexId = (*polygonEdgeVertexLists_[polygonEdgeId]).size();
533 
534  // alloc 1 more triangle
535  (*polygonEdgeTriangleLists_[polygonEdgeId]).resize(
536  (*polygonEdgeTriangleLists_[polygonEdgeId]).size() + 1);
537  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().tetId_ = tetId;
538  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().caseId_ = 0;
539  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().polygonEdgeId_
540  = polygonEdgeId;
541 
542  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().vertexIds_[0] = vertexId;
543  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().vertexIds_[1] =
544  vertexId + 1;
545  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().vertexIds_[2] =
546  vertexId + 2;
547 
548  // alloc 3 more vertices
549  (*polygonEdgeVertexLists_[polygonEdgeId]).resize(vertexId + 3);
550  for(int i = 0; i < 3; i++){
551  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].isBasePoint_ = true;
552  (*polygonEdgeVertexLists_[polygonEdgeId])[
553  vertexId + i].isIntersectionPoint_ = false;
554  }
555 
556  // get the vertex coordinates
557  for(int i = 0; i < 3; i++){
558 
559  int vertexId0, vertexId1;
560 
561  switch(i){
562  case 0:
563  vertexId0 = tetList_[5*tetId + 1
564  + edgeImplicitEncoding_[2*localEdgeId0]];
565  vertexId1 = tetList_[5*tetId + 1
566  + edgeImplicitEncoding_[2*localEdgeId0 + 1]];
567  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.first = u0;
568  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.second = v0;
569  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t0;
570  break;
571 
572  case 1:
573  vertexId0 = tetList_[5*tetId + 1
574  + edgeImplicitEncoding_[2*localEdgeId1]];
575  vertexId1 = tetList_[5*tetId + 1
576  + edgeImplicitEncoding_[2*localEdgeId1 + 1]];
577  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.first = u1;
578  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.second = v1;
579  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t1;
580  break;
581 
582  case 2:
583  vertexId0 = tetList_[5*tetId + 1
584  + edgeImplicitEncoding_[2*localEdgeId2]];
585  vertexId1 = tetList_[5*tetId + 1
586  + edgeImplicitEncoding_[2*localEdgeId2 + 1]];
587  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.first = u2;
588  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.second = v2;
589  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t2;
590  break;
591  }
592 
593  vector<double> baryCentrics;
594  vector<double> p0(2), p1(2), p(2);
595  p0[0] = ((dataTypeU *) uField_)[vertexId0];
596  p0[1] = ((dataTypeV *) vField_)[vertexId0];
597  p1[0] = ((dataTypeU *) uField_)[vertexId1];
598  p1[1] = ((dataTypeV *) vField_)[vertexId1];
599  p[0] = (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.first;
600  p[1] = (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_.second;
602  p0.data(), p1.data(), p.data(), baryCentrics, 2);
603 
604  for(int j = 0; j < 3; j++){
605  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].p_[j] =
606  baryCentrics[0]*pointSet_[3*vertexId0 + j]
607  + baryCentrics[1]*pointSet_[3*vertexId1 + j];
608  }
609 
610  if(vertexId0 < vertexId1)
611  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_
612  = pair<int, int>(vertexId0, vertexId1);
613  else
614  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_
615  = pair<int, int>(vertexId1, vertexId0);
616  }
617 
618  // return the number of created vertices
619  return 3;
620 }
621 
622 template <class dataTypeU, class dataTypeV>
624  const int &polygonEdgeId, const int &tetId,
625  const int &localEdgeId0, const double &t0,
626  const double &u0, const double &v0,
627  const int &localEdgeId1, const double &t1,
628  const double &u1, const double &v1,
629  const int &localEdgeId2, const double &t2,
630  const double &u2, const double &v2) const{
631 
632  int vertexId = (*polygonEdgeVertexLists_[polygonEdgeId]).size();
633 
634  // alloc 5 more vertices
635  (*polygonEdgeVertexLists_[polygonEdgeId]).resize(vertexId + 5);
636  for(int i = 0; i < 5; i++){
637  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].isBasePoint_ = true;
638  (*polygonEdgeVertexLists_[polygonEdgeId])[
639  vertexId + i].isIntersectionPoint_ = false;
640  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_
641  = pair<int, int>(-1, -1);
642  }
643 
644  // alloc 3 more triangles
645  int triangleId = (*polygonEdgeTriangleLists_[polygonEdgeId]).size();
646  (*polygonEdgeTriangleLists_[polygonEdgeId]).resize(triangleId + 3);
647 
648  for(int i = 0; i < 3; i++){
649 
650  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].tetId_ = tetId;
651  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].caseId_ = 1;
652  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].polygonEdgeId_
653  = polygonEdgeId;
654 
655  switch(i){
656  case 0:
657  (*polygonEdgeTriangleLists_[polygonEdgeId])[
658  triangleId + i].vertexIds_[0] = vertexId;
659  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId
660  + i].vertexIds_[1] = vertexId + 1;
661  (*polygonEdgeTriangleLists_[polygonEdgeId])[
662  triangleId + i].vertexIds_[2] = vertexId + 2;
663  break;
664  case 1:
665  (*polygonEdgeTriangleLists_[polygonEdgeId])[
666  triangleId + i].vertexIds_[0] = vertexId + 1;
667  (*polygonEdgeTriangleLists_[polygonEdgeId])[
668  triangleId + i].vertexIds_[1] = vertexId + 2;
669  (*polygonEdgeTriangleLists_[polygonEdgeId])[
670  triangleId + i].vertexIds_[2] = vertexId + 3;
671  break;
672  case 2:
673  (*polygonEdgeTriangleLists_[polygonEdgeId])[
674  triangleId + i].vertexIds_[0] = vertexId + 2;
675  (*polygonEdgeTriangleLists_[polygonEdgeId])[
676  triangleId + i].vertexIds_[1] = vertexId + 3;
677  (*polygonEdgeTriangleLists_[polygonEdgeId])[
678  triangleId + i].vertexIds_[2] = vertexId + 4;
679  break;
680  }
681  }
682 
683  // compute the base triangle vertices like in case 1
684  vector<vector<double> > basePoints(3);
685  vector<pair<double, double> > basePointProjections(3);
686  vector<double> basePointParameterization(3);
687  vector<pair<int, int> > baseEdges(3);
688 
689  computeBaseTriangle<dataTypeU, dataTypeV>(tetId,
690  localEdgeId0, t0, u0, v0,
691  localEdgeId1, t1, u1, v1,
692  localEdgeId2, t2, u2, v2,
693  basePoints, basePointProjections, basePointParameterization, baseEdges);
694 
695  // find the pivot vertex for this case
696  int pivotVertexId = -1;
697 
698  if((t0 >= 0)&&(t0 <= 1)){
699  pivotVertexId = 0;
700  }
701  if((t1 >= 0)&&(t1 <= 1)){
702  pivotVertexId = 1;
703  }
704  if((t2 >= 0)&&(t2 <= 1)){
705  pivotVertexId = 2;
706  }
707 
708  // now get the vertex coordinates
709  for(int i = 0; i < 5; i++){
710 
711  int vertexId0, vertexId1;
712  double t;
713 
714  if(!i){
715  // just take the pivot vertex
716  for(int j = 0; j < 3; j++){
717  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].p_[j] =
718  basePoints[pivotVertexId][j];
719  }
720 
721  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].t_ =
722  basePointParameterization[pivotVertexId];
723  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].uv_ =
724  basePointProjections[pivotVertexId];
725  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].meshEdge_ =
726  baseEdges[pivotVertexId];
727  }
728  else{
729 
730  switch(i){
731 
732  case 1:
733  // interpolation between pivotVertexId and (pivotVertexId-1)%3
734  // if pivot is positive: interpolate at 1
735  // if not interpolate at 0
736  vertexId0 = pivotVertexId;
737  vertexId1 = (pivotVertexId + 2)%3;
738  if(basePointParameterization[(pivotVertexId + 2)%3] > 1){
739  t = 1;
740  }
741  else{
742  t = 0;
743  }
744  break;
745 
746  case 2:
747  // interpolation between pivotVertexId and (pivotVertexId+1)%3
748  // if pivot is positive: interpolate at 1
749  // if not interpolate at 0
750  vertexId0 = pivotVertexId;
751  vertexId1 = (pivotVertexId + 1)%3;
752  if(basePointParameterization[(pivotVertexId + 1)%3] > 1){
753  t = 1;
754  }
755  else{
756  t = 0;
757  }
758  break;
759 
760  case 3:
761  vertexId0 = (pivotVertexId + 2)%3;
762  vertexId1 = (pivotVertexId + 1)%3;
763  if(basePointParameterization[(pivotVertexId + 2)%3] < 0){
764  t = 0;
765  }
766  else{
767  t = 1;
768  }
769  break;
770 
771  case 4:
772  vertexId0 = (pivotVertexId + 2)%3;
773  vertexId1 = (pivotVertexId + 1)%3;
774  if(basePointParameterization[(pivotVertexId + 2)%3] < 0){
775  t = 1;
776  }
777  else{
778  t = 0;
779  }
780  break;
781  }
782 
783  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t;
784 
786  basePoints[vertexId0],
787  basePointProjections[vertexId0],
788  basePointParameterization[vertexId0],
789  basePoints[vertexId1],
790  basePointProjections[vertexId1],
791  basePointParameterization[vertexId1],
792  t,
793  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
794 // snapToBasePoint(
795 // basePoints, basePointProjections, basePointParameterization,
796 // (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
797 
798  }
799  }
800 
801  // return the number of created vertices
802  return 5;
803 }
804 
805 template <class dataTypeU, class dataTypeV>
807  const int &polygonEdgeId, const int &tetId,
808  const int &localEdgeId0, const double &t0,
809  const double &u0, const double &v0,
810  const int &localEdgeId1, const double &t1,
811  const double &u1, const double &v1,
812  const int &localEdgeId2, const double &t2,
813  const double &u2, const double &v2) const{
814 
815  int vertexId = (*polygonEdgeVertexLists_[polygonEdgeId]).size();
816 
817  // alloc 4 more vertices
818  (*polygonEdgeVertexLists_[polygonEdgeId]).resize(vertexId + 4);
819  for(int i = 0; i < 4; i++){
820  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].isBasePoint_ = true;
821  (*polygonEdgeVertexLists_[polygonEdgeId])[
822  vertexId + i].isIntersectionPoint_ = false;
823  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_
824  = pair<int, int>(-1, -1);
825  }
826 
827  // alloc 2 more triangles
828  int triangleId = (*polygonEdgeTriangleLists_[polygonEdgeId]).size();
829  (*polygonEdgeTriangleLists_[polygonEdgeId]).resize(triangleId + 2);
830 
831  for(int i = 0; i < 2; i++){
832 
833  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].tetId_ = tetId;
834  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].caseId_ = 2;
835  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].polygonEdgeId_
836  = polygonEdgeId;
837 
838  if(!i){
839  (*polygonEdgeTriangleLists_[polygonEdgeId])[
840  triangleId + i].vertexIds_[0] = vertexId;
841  (*polygonEdgeTriangleLists_[polygonEdgeId])[
842  triangleId + i].vertexIds_[1] = vertexId + 1;
843  (*polygonEdgeTriangleLists_[polygonEdgeId])[
844  triangleId + i].vertexIds_[2] = vertexId + 2;
845  }
846  else{
847  (*polygonEdgeTriangleLists_[polygonEdgeId])[
848  triangleId + i].vertexIds_[0] = vertexId + 1;
849  (*polygonEdgeTriangleLists_[polygonEdgeId])[
850  triangleId + i].vertexIds_[1] = vertexId + 3;
851  (*polygonEdgeTriangleLists_[polygonEdgeId])[
852  triangleId + i].vertexIds_[2] = vertexId + 2;
853  }
854  }
855 
856  // compute the base triangle vertices like in case 1
857  vector<vector<double> > basePoints(3);
858  vector<pair<double, double> > basePointProjections(3);
859  vector<double> basePointParameterization(3);
860  vector<pair<int, int> > baseEdges(3);
861 
862  computeBaseTriangle<dataTypeU, dataTypeV>(tetId,
863  localEdgeId0, t0, u0, v0,
864  localEdgeId1, t1, u1, v1,
865  localEdgeId2, t2, u2, v2,
866  basePoints, basePointProjections, basePointParameterization, baseEdges);
867 
868  // find the pivot for this case
869  bool isPivotPositive = false;
870  int pivotVertexId = -1;
871 
872  if(((t0 < 0)&&((t1 < 0)||(t2 < 0)))
873  ||((t1 < 0)&&((t0 < 0)||(t2 < 0)))
874  ||((t2 < 0)&&((t1 < 0)||(t0 < 0)))){
875  isPivotPositive = true;
876  }
877  if(isPivotPositive){
878  if(t0 >= 1)
879  pivotVertexId = 0;
880  }
881  else{
882  if(t0 <= 0)
883  pivotVertexId = 0;
884  }
885  if(isPivotPositive){
886  if(t1 >= 1)
887  pivotVertexId = 1;
888  }
889  else{
890  if(t1 <= 0)
891  pivotVertexId = 1;
892  }
893  if(isPivotPositive){
894  if(t2 >= 1)
895  pivotVertexId = 2;
896  }
897  else{
898  if(t2 <= 0)
899  pivotVertexId = 2;
900  }
901 
902  // now get the vertex coordinates
903  for(int i = 0; i < 4; i++){
904 
905  int vertexId0, vertexId1;
906  double t;
907 
908  switch(i){
909 
910  case 0:
911  // interpolation between pivotVertexId and (pivotVertexId-1)%3
912  // if pivot is positive: interpolate at 1
913  // if not interpolate at 0
914  vertexId0 = pivotVertexId;
915  vertexId1 = (pivotVertexId + 2)%3;
916  if(isPivotPositive)
917  t = 1;
918  else
919  t = 0;
920  break;
921 
922  case 1:
923  // interpolation between pivotVertexId and (pivotVertexId+1)%3
924  // if pivot is positive: interpolate at 1
925  // if not interpolate at 0
926  vertexId0 = pivotVertexId;
927  vertexId1 = (pivotVertexId + 1)%3;
928  if(isPivotPositive)
929  t = 1;
930  else
931  t = 0;
932  break;
933 
934  case 2:
935  // interpolation between pivotVertexId and (pivotVertexId-1)%3
936  // if pivot is positive: interpolate at 0
937  // if not interpolate at 1
938  vertexId0 = pivotVertexId;
939  vertexId1 = (pivotVertexId + 2)%3;
940  if(isPivotPositive)
941  t = 0;
942  else
943  t = 1;
944  break;
945 
946  case 3:
947  // interpolation between pivotVertexId and (pivotVertexId+1)%3
948  // if pivot is positive: interpolate at 0
949  // if not interpolate at 1
950  vertexId0 = pivotVertexId;
951  vertexId1 = (pivotVertexId + 1)%3;
952  if(isPivotPositive)
953  t = 0;
954  else
955  t = 1;
956  break;
957  }
958 
959  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t;
960 
962  basePoints[vertexId0],
963  basePointProjections[vertexId0],
964  basePointParameterization[vertexId0],
965  basePoints[vertexId1],
966  basePointProjections[vertexId1],
967  basePointParameterization[vertexId1],
968  t,
969  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
970 // snapToBasePoint(
971 // basePoints, basePointProjections, basePointParameterization,
972 // (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
973 
974  }
975 
976  // return the number of created vertices
977  return 4;
978 }
979 
980 template <class dataTypeU, class dataTypeV>
982  const int &polygonEdgeId, const int &tetId,
983  const int &localEdgeId0, const double &t0,
984  const double &u0, const double &v0,
985  const int &localEdgeId1, const double &t1,
986  const double &u1, const double &v1,
987  const int &localEdgeId2, const double &t2,
988  const double &u2, const double &v2) const{
989 
990  int vertexId = (*polygonEdgeVertexLists_[polygonEdgeId]).size();
991 
992  // alloc 3 more vertices
993  (*polygonEdgeVertexLists_[polygonEdgeId]).resize(vertexId + 3);
994  for(int i = 0; i < 3; i++){
995  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].isBasePoint_ = true;
996  (*polygonEdgeVertexLists_[polygonEdgeId])[
997  vertexId + i].isIntersectionPoint_ = false;
998  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_
999  = pair<int, int>(-1, -1);
1000  }
1001 
1002  // alloc 1 more triangle
1003  int triangleId = (*polygonEdgeTriangleLists_[polygonEdgeId]).size();
1004  (*polygonEdgeTriangleLists_[polygonEdgeId]).resize(triangleId + 1);
1005 
1006 
1007  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId].tetId_ = tetId;
1008  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId].caseId_ = 3;
1009  (*polygonEdgeTriangleLists_[polygonEdgeId]).back().polygonEdgeId_
1010  = polygonEdgeId;
1011 
1012  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId].vertexIds_[0] =
1013  vertexId;
1014  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId].vertexIds_[1] =
1015  vertexId + 1;
1016  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId].vertexIds_[2] =
1017  vertexId + 2;
1018 
1019  // compute the base triangle vertices like in case 1
1020  vector<vector<double> > basePoints(3);
1021  vector<pair<double, double> > basePointProjections(3);
1022  vector<double> basePointParameterization(3);
1023  vector<pair<int, int> > baseEdges(3);
1024 
1025  computeBaseTriangle<dataTypeU, dataTypeV>(tetId,
1026  localEdgeId0, t0, u0, v0,
1027  localEdgeId1, t1, u1, v1,
1028  localEdgeId2, t2, u2, v2,
1029  basePoints, basePointProjections, basePointParameterization, baseEdges);
1030 
1031  // now find the pivot
1032  bool isPivotPositive = false;
1033  int pivotVertexId = -1;
1034 
1035  if((t0 <= 1)&&(t0 >= 0)){
1036  pivotVertexId = 0;
1037  }
1038  else{
1039  if(t0 < 0)
1040  isPivotPositive = true;
1041  }
1042  if((t1 <= 1)&&(t1 >= 0)){
1043  pivotVertexId = 1;
1044  }
1045  else{
1046  if(t1 < 0)
1047  isPivotPositive = true;
1048  }
1049  if((t2 <= 1)&&(t2 >= 0)){
1050  pivotVertexId = 2;
1051  }
1052  else{
1053  if(t2 < 0)
1054  isPivotPositive = true;
1055  }
1056 
1057  // now get the vertex coordinates
1058  for(int i = 0; i < 3; i++){
1059 
1060  int vertexId0, vertexId1;
1061  double t;
1062 
1063  if(!i){
1064  // special case of the pivot vertex
1065  for(int j = 0; j < 3; j++){
1066  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].p_[j] =
1067  basePoints[pivotVertexId][j];
1068  }
1069 
1070  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].t_ =
1071  basePointParameterization[pivotVertexId];
1072  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].uv_ =
1073  basePointProjections[pivotVertexId];
1074  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId].meshEdge_ =
1075  baseEdges[pivotVertexId];
1076  }
1077  else{
1078  if(i == 1){
1079  // interpolation between pivotVertexId and pivotVertexId+1
1080  // if pivot is positive, interpolate at value 0
1081  // else interpolate at value 1
1082  vertexId0 = pivotVertexId;
1083  vertexId1 = (pivotVertexId + 1)%3;
1084  if(isPivotPositive)
1085  t = 0;
1086  else
1087  t = 1;
1088  }
1089  if(i == 2){
1090  // interpolation between pivotVertexId and pivotVertexId-1
1091  // if pivot is positive, interpolate at value 0
1092  // else interpolate at value 1
1093  vertexId0 = pivotVertexId;
1094  vertexId1 = (pivotVertexId + 2)%3;
1095  if(isPivotPositive)
1096  t = 0;
1097  else
1098  t = 1;
1099  }
1100 
1101  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t;
1102 
1104  basePoints[vertexId0],
1105  basePointProjections[vertexId0],
1106  basePointParameterization[vertexId0],
1107  basePoints[vertexId1],
1108  basePointProjections[vertexId1],
1109  basePointParameterization[vertexId1],
1110  t,
1111  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
1112 // snapToBasePoint(
1113 // basePoints, basePointProjections, basePointParameterization,
1114 // (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
1115 
1116  }
1117  }
1118 
1119  // return the number of created vertices
1120  return 3;
1121 }
1122 
1123 template <class dataTypeU, class dataTypeV>
1125  const int &polygonEdgeId, const int &tetId,
1126  const int &localEdgeId0, const double &t0,
1127  const double &u0, const double &v0,
1128  const int &localEdgeId1, const double &t1,
1129  const double &u1, const double &v1,
1130  const int &localEdgeId2, const double &t2,
1131  const double &u2, const double &v2) const{
1132 
1133  int vertexId = (*polygonEdgeVertexLists_[polygonEdgeId]).size();
1134 
1135  // alloc 4 more vertices
1136  (*polygonEdgeVertexLists_[polygonEdgeId]).resize(vertexId + 4);
1137  for(int i = 0; i < 4; i++){
1138  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].isBasePoint_ = true;
1139  (*polygonEdgeVertexLists_[polygonEdgeId])[
1140  vertexId + i].isIntersectionPoint_ = false;
1141  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_
1142  = pair<int, int>(-1, -1);
1143  }
1144 
1145  // alloc 2 more triangles
1146  int triangleId = (*polygonEdgeTriangleLists_[polygonEdgeId]).size();
1147  (*polygonEdgeTriangleLists_[polygonEdgeId]).resize(triangleId + 2);
1148 
1149  for(int i = 0; i < 2; i++){
1150 
1151  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].tetId_ = tetId;
1152  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].caseId_ = 4;
1153  (*polygonEdgeTriangleLists_[polygonEdgeId])[triangleId + i].polygonEdgeId_
1154  = polygonEdgeId;
1155 
1156  if(!i){
1157  (*polygonEdgeTriangleLists_[polygonEdgeId])[
1158  triangleId + i].vertexIds_[0] = vertexId;
1159  (*polygonEdgeTriangleLists_[polygonEdgeId])[
1160  triangleId + i].vertexIds_[1] = vertexId + 1;
1161  (*polygonEdgeTriangleLists_[polygonEdgeId])[
1162  triangleId + i].vertexIds_[2] = vertexId + 2;
1163  }
1164  else{
1165  (*polygonEdgeTriangleLists_[polygonEdgeId])[
1166  triangleId + i].vertexIds_[0] = vertexId + 1;
1167  (*polygonEdgeTriangleLists_[polygonEdgeId])[
1168  triangleId + i].vertexIds_[1] = vertexId + 3;
1169  (*polygonEdgeTriangleLists_[polygonEdgeId])[
1170  triangleId + i].vertexIds_[2] = vertexId + 2;
1171  }
1172  }
1173 
1174  // compute the base triangle vertices like in case 1
1175  vector<vector<double> > basePoints(3);
1176  vector<pair<double, double> > basePointProjections(3);
1177  vector<double> basePointParameterization(3);
1178  vector<pair<int, int> > baseEdges(3);
1179 
1180  computeBaseTriangle<dataTypeU, dataTypeV>(tetId,
1181  localEdgeId0, t0, u0, v0,
1182  localEdgeId1, t1, u1, v1,
1183  localEdgeId2, t2, u2, v2,
1184  basePoints, basePointProjections, basePointParameterization, baseEdges);
1185 
1186  // find the pivot vertex for this case
1187  bool isPivotPositive = false;
1188  int pivotVertexId = -1;
1189 
1190  if(t0 > 1){
1191  pivotVertexId = 0;
1192  isPivotPositive = true;
1193  }
1194  else if(t0 < 0){
1195  pivotVertexId = 0;
1196  isPivotPositive = false;
1197  }
1198 
1199  if(t1 > 1){
1200  pivotVertexId = 1;
1201  isPivotPositive = true;
1202  }
1203  else if(t1 < 0){
1204  pivotVertexId = 1;
1205  isPivotPositive = false;
1206  }
1207  if(t2 > 1){
1208  pivotVertexId = 2;
1209  isPivotPositive = true;
1210  }
1211  else if(t2 < 0){
1212  pivotVertexId = 2;
1213  isPivotPositive = false;
1214  }
1215 
1216  // now get the vertex coordinates depending on the case
1217  for(int i = 0; i < 4; i++){
1218 
1219  int vertexId0, vertexId1;
1220  double t;
1221 
1222  if(i < 2){
1223  if(!i){
1224  // interpolation between pivotVertexId and (pivotVertexId-1)%3
1225  // if pivot is positive: interpolate at 1
1226  // if not interpolate at 0
1227  vertexId0 = pivotVertexId;
1228  vertexId1 = (pivotVertexId + 2)%3;
1229  if(isPivotPositive)
1230  t = 1;
1231  else
1232  t = 0;
1233  }
1234  else{
1235  // interpolation between pivotVertexId and (pivotVertexId+1)%3
1236  // if pivot is positive: interpolate at 1
1237  // if not interpolate at 0
1238  vertexId0 = pivotVertexId;
1239  vertexId1 = (pivotVertexId + 1)%3;
1240  if(isPivotPositive)
1241  t = 1;
1242  else
1243  t = 0;
1244  }
1245 
1246  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ = t;
1247 
1249  basePoints[vertexId0],
1250  basePointProjections[vertexId0],
1251  basePointParameterization[vertexId0],
1252  basePoints[vertexId1],
1253  basePointProjections[vertexId1],
1254  basePointParameterization[vertexId1],
1255  t,
1256  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
1257 // snapToBasePoint(
1258 // basePoints, basePointProjections, basePointParameterization,
1259 // (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i]);
1260 
1261  }
1262  else{
1263  if(i == 2){
1264  // take (pivotVertexId-1)%3
1265  for(int j = 0; j < 3; j++){
1266  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].p_[j] =
1267  basePoints[(pivotVertexId + 2)%3][j];
1268  }
1269 
1270  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ =
1271  basePointParameterization[(pivotVertexId + 2)%3];
1272  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_ =
1273  basePointProjections[(pivotVertexId + 2)%3];
1274  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_ =
1275  baseEdges[(pivotVertexId + 2)%3];
1276  }
1277  else{
1278  // take (pivtoVertexId+1)%3
1279  for(int j = 0; j < 3; j++){
1280  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].p_[j] =
1281  basePoints[(pivotVertexId + 1)%3][j];
1282  }
1283  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].t_ =
1284  basePointParameterization[(pivotVertexId + 1)%3];
1285  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].uv_ =
1286  basePointProjections[(pivotVertexId + 1)%3];
1287  (*polygonEdgeVertexLists_[polygonEdgeId])[vertexId + i].meshEdge_ =
1288  baseEdges[(pivotVertexId + 1)%3];
1289  }
1290  }
1291  }
1292 
1293  // return the number of created vertices
1294  return 4;
1295 }
1296 
1297 template <class dataTypeU, class dataTypeV>
1299  const pair<double, double> &rangePoint0,
1300  const pair<double, double> &rangePoint1,
1301  const vector<int> &seedTetList,
1302  const int &polygonEdgeId) const{
1303 
1304 #ifndef withKamikaze
1305  if(!tetNeighbors_)
1306  return -1;
1307  if(!tetNumber_)
1308  return -2;
1309  if(!tetList_)
1310  return -3;
1311  if(!uField_)
1312  return -4;
1313  if(!vField_)
1314  return -5;
1315  if(!pointSet_)
1316  return -6;
1317  if(!polygonEdgeNumber_)
1318  return -7;
1319  if(!globalVertexList_)
1320  return -8;
1321 #endif
1322 
1323  vector<bool> visitedTets(tetNumber_, false);
1324  queue<int> tetQueue;
1325 
1326  // init the queue
1327  for(int i = 0; i < (int) seedTetList.size(); i++){
1328  tetQueue.push(seedTetList[i]);
1329  }
1330 
1331  int createdVertices = 0;
1332 
1333  do{
1334 
1335  int tetId = tetQueue.front();
1336  tetQueue.pop();
1337 
1338  if(!visitedTets[tetId]){
1339 
1340  createdVertices = processTetrahedron<dataTypeU, dataTypeV>(
1341  tetId, rangePoint0, rangePoint1, polygonEdgeId);
1342 
1343  if(createdVertices){
1344  // only propagate if we created a triangle
1345  for(int i = 0; i < (int) (*tetNeighbors_)[tetId].size(); i++){
1346  if(!visitedTets[(*tetNeighbors_)[tetId][i]])
1347  tetQueue.push((*tetNeighbors_)[tetId][i]);
1348  }
1349  }
1350 
1351  visitedTets[tetId] = true;
1352  }
1353 
1354  }while(tetQueue.size());
1355 
1356  return 0;
1357 }
1358 
1359 template <class dataTypeU, class dataTypeV> int FiberSurface::computeContour(
1360  const vector<pair<pair<double, double>, pair<double, double> > > &edgeList,
1361  const vector<int> &seedTetList,
1362  const vector<int> *edgeIdList) const{
1363 
1364 #ifndef withKamikaze
1365  if(!tetNeighbors_)
1366  return -1;
1367  if(!tetNumber_)
1368  return -2;
1369  if(!tetList_)
1370  return -3;
1371  if(!uField_)
1372  return -4;
1373  if(!vField_)
1374  return -5;
1375  if(!pointSet_)
1376  return -6;
1377  if(!polygonEdgeNumber_)
1378  return -7;
1379  if(!globalVertexList_)
1380  return -8;
1381 #endif
1382 
1383  vector<bool> visitedTets(tetNumber_, false);
1384  queue<int> tetQueue;
1385 
1386  // init the queue
1387  for(int i = 0; i < (int) seedTetList.size(); i++){
1388  tetQueue.push(seedTetList[i]);
1389  }
1390 
1391  int createdVertices = 0;
1392 
1393  do{
1394 
1395  int tetId = tetQueue.front();
1396  tetQueue.pop();
1397 
1398  if(!visitedTets[tetId]){
1399 
1400  vector<vector<int> > threadedTetQueue(edgeList.size());
1401 #ifdef withOpenMP
1402 #pragma omp parallel for num_threads(threadNumber_)
1403 #endif
1404  for(int i = 0; i < (int) edgeList.size(); i++){
1405 
1406  int polygonEdgeId = 0;
1407 
1408  if(edgeIdList){
1409  polygonEdgeId = (*edgeIdList)[i];
1410  }
1411 
1412  createdVertices = processTetrahedron<dataTypeU, dataTypeV>(
1413  tetId,
1414  edgeList[i].first,
1415  edgeList[i].second,
1416  polygonEdgeId);
1417 
1418  if(createdVertices){
1419  // only propagate if we created a triangle
1420  for(int j = 0; j < (int) (*tetNeighbors_)[tetId].size(); j++){
1421  if(!visitedTets[(*tetNeighbors_)[tetId][j]]){
1422  threadedTetQueue[i].push_back((*tetNeighbors_)[tetId][j]);
1423  }
1424  }
1425  }
1426  }
1427 
1428  visitedTets[tetId] = true;
1429 
1430  for(int i = 0; i < (int) threadedTetQueue.size(); i++){
1431  for(int j = 0; j < (int) threadedTetQueue[i].size(); j++){
1432  tetQueue.push(threadedTetQueue[i][j]);
1433  }
1434  }
1435  }
1436 
1437  }while(tetQueue.size());
1438 
1439  return 0;
1440 }
1441 
1442 template <class dataTypeU, class dataTypeV>
1444  const pair<double, double> &rangePoint0,
1445  const pair<double, double> &rangePoint1,
1446  const int &polygonEdgeId) const {
1447 
1448 #ifndef withKamikaze
1449  if(!tetNumber_)
1450  return -1;
1451  if(!tetList_)
1452  return -2;
1453  if(!uField_)
1454  return -3;
1455  if(!vField_)
1456  return -4;
1457  if(!pointSet_)
1458  return -5;
1459  if(!polygonEdgeNumber_)
1460  return -6;
1461  if(!globalVertexList_)
1462  return -7;
1463 #endif
1464 
1465 #ifdef withOpenMP
1466 #pragma omp parallel for num_threads(threadNumber_)
1467 #endif
1468  for(int i = 0; i < tetNumber_; i++){
1469 
1470  processTetrahedron<dataTypeU, dataTypeV>(
1471  i, rangePoint0, rangePoint1, polygonEdgeId);
1472  }
1473 
1474  return 0;
1475 }
1476 
1477 template <class dataTypeU, class dataTypeV>
1479 
1480 #ifndef withKamikaze
1481  if(!tetNumber_)
1482  return -1;
1483  if(!tetList_)
1484  return -2;
1485  if(!uField_)
1486  return -3;
1487  if(!vField_)
1488  return -4;
1489  if(!pointSet_)
1490  return -5;
1491  if(!polygon_)
1492  return -6;
1493  if(polygonEdgeNumber_ != (int) polygon_->size())
1494  return -7;
1495  if(!globalVertexList_)
1496  return -8;
1497 #endif
1498 
1499  Timer t;
1500 
1501 #ifdef withrangeDrivenOctree
1502  if(!octree_.empty()){
1503 
1504 #ifdef withOpenMP
1505 #pragma omp parallel for num_threads(threadNumber_)
1506 #endif
1507  for(int i = 0; i < polygonEdgeNumber_; i++){
1508 
1509  computeSurfaceWithOctree<dataTypeU, dataTypeV>(
1510  (*polygon_)[i].first, (*polygon_)[i].second, i);
1511  }
1512  }
1513  else{
1514  // regular extraction (the octree has not been computed)
1515 #ifdef withOpenMP
1516 #pragma omp parallel for num_threads(threadNumber_)
1517 #endif
1518  for(int i = 0; i < polygonEdgeNumber_; i++){
1519 
1520  computeSurface<dataTypeU, dataTypeV>(
1521  (*polygon_)[i].first, (*polygon_)[i].second, i);
1522  }
1523  }
1524 
1525 #else
1526 #ifdef withOpenMP
1527 #pragma omp parallel for num_threads(threadNumber_)
1528 #endif
1529  for(int i = 0; i < polygonEdgeNumber_; i++){
1530 
1531  computeSurface<dataTypeU, dataTypeV>(
1532  (*polygon_)[i].first, (*polygon_)[i].second, i);
1533  }
1534 #endif
1535 
1536  finalize<dataTypeU, dataTypeV>(pointSnapping_, NULL, NULL, NULL);
1537 
1538  {
1539  stringstream msg;
1540  msg << "[FiberSurface] FiberSurface extracted in "
1541  << t.getElapsedTime() << " s. (" << globalVertexList_->size()
1542  << " vertices, " << threadNumber_
1543  << " thread(s))" << endl;
1544  dMsg(cout, msg.str(), timeMsg);
1545  }
1546 
1547  return 0;
1548 }
1549 
1550 #ifdef withrangeDrivenOctree
1551 template <class dataTypeU, class dataTypeV>
1552  inline int FiberSurface::computeSurfaceWithOctree(
1553  const pair<double, double> &rangePoint0,
1554  const pair<double, double> &rangePoint1,
1555  const int &polygonEdgeId) const {
1556 
1557 #ifndef withKamikaze
1558  if(!tetNumber_)
1559  return -1;
1560  if(!tetList_)
1561  return -2;
1562  if(!uField_)
1563  return -3;
1564  if(!vField_)
1565  return -4;
1566  if(!pointSet_)
1567  return -5;
1568  if(!polygonEdgeNumber_)
1569  return -6;
1570  if(!globalVertexList_)
1571  return -7;
1572 #endif
1573 
1574  vector<int> tetList;
1575  octree_.rangeSegmentQuery(rangePoint0, rangePoint1, tetList);
1576 
1577 #ifdef withOpenMP
1578 #pragma omp parallel for num_threads(threadNumber_)
1579 #endif
1580  for(int i = 0; i < (int) tetList.size(); i++){
1581  processTetrahedron<dataTypeU, dataTypeV>(
1582  tetList[i], rangePoint0, rangePoint1, polygonEdgeId);
1583  }
1584 
1585  return 0;
1586 }
1587 #endif
1588 
1589 template <class dataTypeU, class dataTypeV> int FiberSurface::finalize(
1590  const bool &mergeDuplicatedVertices,
1591  const bool &removeSmallEdges,
1592  const bool &edgeFlips,
1593  const bool &intersectionRemesh){
1594 
1595  // make only one vertex list
1596  int fiberSurfaceVertexNumber = 0;
1597  for(int i = 0; i < (int) polygonEdgeVertexLists_.size(); i++){
1598  fiberSurfaceVertexNumber += (*polygonEdgeVertexLists_[i]).size();
1599  }
1600 
1601  (*globalVertexList_).resize(fiberSurfaceVertexNumber);
1602  fiberSurfaceVertexNumber = 0;
1603  for(int i = 0; i < (int) polygonEdgeVertexLists_.size(); i++){
1604  for(int j = 0; j < (int) polygonEdgeVertexLists_[i]->size(); j++){
1605  (*polygonEdgeVertexLists_[i])[j].polygonEdgeId_ = i;
1606  (*polygonEdgeVertexLists_[i])[j].localId_ = j;
1607  (*polygonEdgeVertexLists_[i])[j].globalId_ = fiberSurfaceVertexNumber;
1608  (*globalVertexList_)[fiberSurfaceVertexNumber] =
1609  (*polygonEdgeVertexLists_[i])[j];
1610  fiberSurfaceVertexNumber++;
1611  }
1612  }
1613  for(int i = 0; i < (int) polygonEdgeTriangleLists_.size(); i++){
1614  for(int j = 0; j < (int) polygonEdgeTriangleLists_[i]->size(); j++){
1615  for(int k = 0; k < 3; k++){
1616  (*polygonEdgeTriangleLists_[i])[j].vertexIds_[k] =
1618  (*polygonEdgeTriangleLists_[i])[j].vertexIds_[k]].globalId_;
1619 
1620  }
1621  }
1622  }
1623 
1624  // NOTE:
1625  // 1) in a first step, make everything work in a POST process.
1626  // this is what really matters for tet gen.
1627  // 2) probably come back and apply in a pre-process (more degenerate
1628  // intersection cases).
1629  if(intersectionRemesh){
1630  remeshIntersections<dataTypeU, dataTypeV>();
1631  }
1632 
1633  if((mergeDuplicatedVertices)||(removeSmallEdges)){
1634 // we need to have a complex to perform the edge collapses
1636  }
1637 
1638  if(edgeFlips)
1639  flipEdges();
1640 
1641  if(removeSmallEdges)
1643 
1644  // now we can release the memory for the threaded vertices
1645  for(int i = 0; i < (int) polygonEdgeVertexLists_.size(); i++){
1646  polygonEdgeVertexLists_[i]->clear();
1647  }
1648 
1649  return 0;
1650 }
1651 
1652 template <class dataTypeU, class dataTypeV>
1653  inline int FiberSurface::processTetrahedron(const int &tetId,
1654  const pair<double, double> &rangePoint0,
1655  const pair<double, double> &rangePoint1,
1656  const int &polygonEdgeId) const{
1657 
1658  double rangeEdge[2];
1659  rangeEdge[0] = rangePoint0.first - rangePoint1.first;
1660  rangeEdge[1] = rangePoint0.second - rangePoint1.second;
1661 
1662  double rangeNormal[2];
1663  rangeNormal[0] = -rangeEdge[1];
1664  rangeNormal[1] = rangeEdge[0];
1665 
1666  // 1. compute the distance to the range line carrying the saddleEdge
1667  int upperNumber = 0;
1668  int lowerNumber = 0;
1669  int equalVertexLocalId = -1;
1670  double d[4];
1671  for(int i = 0; i < 4; i++){
1672 
1673  int vertexId = tetList_[5*tetId + 1 + i];
1674 
1675  double projectedVertex[2];
1676  projectedVertex[0] = ((dataTypeU *) uField_)[vertexId];
1677  projectedVertex[1] = ((dataTypeV *) vField_)[vertexId];
1678 
1679  double vertexRangeEdge[2];
1680  vertexRangeEdge[0] = projectedVertex[0] - rangePoint0.first;
1681  vertexRangeEdge[1] = projectedVertex[1] - rangePoint0.second;
1682 
1683  d[i] = vertexRangeEdge[0]*rangeNormal[0]
1684  + vertexRangeEdge[1]*rangeNormal[1];
1685 
1686  if(fabs(d[i]) < pow(10, -DBL_DIG))
1687  d[i] = 0;
1688 
1689  if(d[i] > 0)
1690  upperNumber++;
1691  if(d[i] < 0)
1692  lowerNumber++;
1693 
1694  if(d[i] == 0){
1695  equalVertexLocalId = i;
1696  }
1697  }
1698 
1699  // 2. compute the base triangle(s)
1700  if(!((upperNumber == 0)||(lowerNumber == 0))){
1701 
1702  // the fiber surface is passing through this tetrahedron.
1703  vector<bool> lonelyVertex(4, false);
1704  vector<int> triangleEdgeNumbers(2, 0);
1705  vector<vector<int> > triangleEdges(2);
1706  triangleEdges[0].resize(3, -1);
1707  triangleEdges[1].resize(3, -1);
1708 
1709  // implicit edge encoding
1710  // 0: O-1
1711  // 1: 0-2
1712  // 2: 0-3
1713  // 3: 3-1 [order!]
1714  // 4: 2-1 [order!]
1715  // 5: 2-3
1716  int edgeCounter = 0;
1717  for(int i = 0; i < 4; i++){
1718 
1719  int jStart = i + 1;
1720  int jEnd = 4;
1721  int jStep = 1;
1722 
1723  if(i == 1){
1724  // special ordering here
1725  // (to facilitate the creation of valid base triangles)
1726  // any two consecutive edges shall share a vertex
1727  jStart = 3;
1728  jEnd = i;
1729  jStep = -1;
1730  }
1731 
1732  for(int j = jStart; j != jEnd; j += jStep){
1733 
1734  if(((d[i] > 0)&&(d[j] < 0))||((d[i] < 0)&&(d[j] > 0))){
1735 
1736  // the edge is crossed by a base triangle
1737  if(triangleEdgeNumbers[0] == 3){
1738  triangleEdges[1][triangleEdgeNumbers[1]] = edgeCounter;
1739  triangleEdgeNumbers[1]++;
1740  }
1741  else{
1742  triangleEdges[0][triangleEdgeNumbers[0]] = edgeCounter;
1743  triangleEdgeNumbers[0]++;
1744  }
1745  }
1746 
1747  if((d[i] == 0)&&(d[j] == 0)){
1748  // special case of a seed tet containing the jacobi edge.
1749  // the entire edge is on the fiber surface.
1750  // let's put this edge twice.
1751  // NOTE: in such a case, we're producing only one triangle.
1752  // so we just need to add that to the first triangle.
1753  triangleEdges[0][triangleEdgeNumbers[0]] = edgeCounter;
1754  triangleEdgeNumbers[0]++;
1755 
1756  triangleEdges[0][triangleEdgeNumbers[0]] = edgeCounter;
1757  triangleEdgeNumbers[0]++;
1758  }
1759 
1760  edgeCounter++;
1761  }
1762  }
1763 
1764  // post-process in the case of a second base triangle
1765  if(triangleEdges[1][0] != -1){
1766  if(triangleEdges[1][1] == -1){
1767  // we need to consider the edges of the first triangle which share
1768  // a vertex with the second triangle's edge.
1769  // given the ordering, the following edge pairs are forbidden:
1770  // (opposite edges of the tetrahedron)
1771  // 0/5
1772  // 1/3
1773  // 2/4
1774  int forbiddenEdge = -1;
1775  switch(triangleEdges[1][0]){
1776  case 0:
1777  forbiddenEdge = 5;
1778  break;
1779  case 1:
1780  forbiddenEdge = 3;
1781  break;
1782  case 2:
1783  forbiddenEdge = 4;
1784  break;
1785  case 3:
1786  forbiddenEdge = 1;
1787  break;
1788  case 4:
1789  forbiddenEdge = 2;
1790  break;
1791  case 5:
1792  forbiddenEdge = 0;
1793  break;
1794  }
1795  for(int i = 0; i < (int) triangleEdges[0].size(); i++){
1796  if(triangleEdges[0][i] != forbiddenEdge){
1797  if(triangleEdges[1][1] != -1){
1798  triangleEdges[1][2] = triangleEdges[0][i];
1799  break;
1800  }
1801  else{
1802  triangleEdges[1][1] = triangleEdges[0][i];
1803  }
1804  }
1805  }
1806  }
1807  }
1808 
1809  // post-process in case of exactly one vertex on the jacobi edge
1810  for(int i = 0; i < 2; i++){
1811  if(triangleEdges[i][0] != -1){
1812  // the jacobi vertex has to be the last one
1813  if(triangleEdges[i][2] == -1){
1814  // add whatever edge connected to equalVertexLocalId
1815  switch(equalVertexLocalId){
1816  case 0:
1817  case 1:
1818  triangleEdges[i][2] = 0;
1819  break;
1820  case 2:
1821  case 3:
1822  triangleEdges[i][2] = 5;
1823  break;
1824  }
1825  }
1826  }
1827  }
1828 
1829  // 3. crop the resulting triangles to the saddleEdge
1830  // for each edge recorded previously, get the u,v coordinates of the
1831  // intersection point.
1832  // take its barycentric coordinates on the saddle edge.
1833  // see if it lies in it (in between 0 and 1)
1834  // figure 7 of the paper
1835  double d0, d1;
1836  pair<double, double> uv0, uv1;
1837  vector<pair<double, double> > uv(3);
1838  vector<double> t(3);
1839 
1840  int createdVertices = 0;
1841 
1842  for(int i = 0; i < 2; i++){
1843  if(triangleEdges[i][0] != -1){
1844  // this is a valid triangle, let's go ahead.
1845 
1846  int lowerVertexNumber = 0;
1847  int upperVertexNumber = 0;
1848  int greyVertexNumber = 0;
1849  // iterate over the edges and compute the edge intersection (range)
1850  for(int j = 0; j < (int) triangleEdges[i].size(); j++){
1851 
1852  if((j < (int) triangleEdges[i].size() - 1)
1853  &&(triangleEdges[i][j] == triangleEdges[i][j + 1])){
1854  // special case of a jacobi edge
1855  uv[j].first = ((dataTypeU *) uField_)[tetList_[
1856  5*tetId + 1 + edgeImplicitEncoding_[2*triangleEdges[i][j]]]];
1857  uv[j].second = ((dataTypeV *) vField_)[tetList_[
1858  5*tetId + 1 + edgeImplicitEncoding_[2*triangleEdges[i][j]]]];
1859 
1860  uv[j + 1].first = ((dataTypeU *) uField_)[tetList_[
1861  5*tetId + 1
1862  + edgeImplicitEncoding_[2*triangleEdges[i][j] + 1]]];
1863  uv[j + 1].second = ((dataTypeV *) vField_)[tetList_[
1864  5*tetId + 1
1865  + edgeImplicitEncoding_[2*triangleEdges[i][j] + 1]]];
1866 
1867  }
1868  else if((!j)||
1869  ((j)&&(triangleEdges[i][j] != triangleEdges[i][j - 1]))){
1870 
1871  // regular intersection case
1872  d0 = d[edgeImplicitEncoding_[2*triangleEdges[i][j]]];
1873  uv0.first = ((dataTypeU *) uField_)[tetList_[
1874  5*tetId + 1 + edgeImplicitEncoding_[2*triangleEdges[i][j]]]];
1875  uv0.second = ((dataTypeV *) vField_)[tetList_[
1876  5*tetId + 1 + edgeImplicitEncoding_[2*triangleEdges[i][j]]]];
1877 
1878  d1 = d[edgeImplicitEncoding_[2*triangleEdges[i][j] + 1]];
1879  uv1.first = ((dataTypeU *) uField_)[tetList_[
1880  5*tetId + 1
1881  + edgeImplicitEncoding_[2*triangleEdges[i][j] + 1]]];
1882  uv1.second = ((dataTypeV *) vField_)[tetList_[
1883  5*tetId + 1
1884  + edgeImplicitEncoding_[2*triangleEdges[i][j] + 1]]];
1885 
1886  uv[j].first = uv0.first
1887  + (d0/(d0 - d1))*(uv1.first - uv0.first);
1888  uv[j].second = uv0.second
1889  + (d0/(d0 - d1))*(uv1.second - uv0.second);
1890  }
1891 
1892  // now determine the line parameterization of this intersection
1893  // point on the saddle edge
1894  if(fabs(rangePoint1.first - rangePoint0.first)
1895  > fabs(rangePoint1.second - rangePoint0.second)){
1896  t[j] = (uv[j].first - rangePoint0.first)
1897  / (rangePoint1.first - rangePoint0.first);
1898  }
1899  else{
1900  t[j] = (uv[j].second - rangePoint0.second)
1901  / (rangePoint1.second - rangePoint0.second);
1902  }
1903 
1904  if((t[j] <= 1)&&(t[j] >= 0))
1905  greyVertexNumber++;
1906  else if(t[j] < 0)
1907  lowerVertexNumber++;
1908  else
1909  upperVertexNumber++;
1910  }
1911  // at this point, we know the uv coordinates (and the edge param)
1912  // of each vertex of the base triangle.
1913  // we can proceed with the cropping
1914 
1915 
1916  // 4. triangulate the result
1917  if(greyVertexNumber == 3){
1918  createdVertices += computeCase0<dataTypeU, dataTypeV>(
1919  polygonEdgeId, tetId,
1920  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1921  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1922  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1923  }
1924  else if(lowerVertexNumber == 3){
1925  // well do nothing (empty triangle)
1926  }
1927  else if(upperVertexNumber == 3){
1928  // well do nothing (empty triangle)
1929  }
1930  else if((lowerVertexNumber == 1)
1931  &&(upperVertexNumber == 1)
1932  &&(greyVertexNumber == 1)){
1933  createdVertices += computeCase1<dataTypeU, dataTypeV>(
1934  polygonEdgeId, tetId,
1935  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1936  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1937  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1938  }
1939  else if((lowerVertexNumber == 2)&&(upperVertexNumber == 1)){
1940  createdVertices += computeCase2<dataTypeU, dataTypeV>(
1941  polygonEdgeId, tetId,
1942  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1943  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1944  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1945  }
1946  else if((lowerVertexNumber == 1)&&(upperVertexNumber == 2)){
1947  createdVertices += computeCase2<dataTypeU, dataTypeV>(
1948  polygonEdgeId, tetId,
1949  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1950  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1951  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1952  }
1953  else if((greyVertexNumber == 1)
1954  &&(lowerVertexNumber == 2)){
1955  createdVertices += computeCase3<dataTypeU, dataTypeV>(
1956  polygonEdgeId, tetId,
1957  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1958  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1959  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1960  }
1961  else if((greyVertexNumber == 1)
1962  &&(upperVertexNumber == 2)){
1963  createdVertices += computeCase3<dataTypeU, dataTypeV>(
1964  polygonEdgeId, tetId,
1965  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1966  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1967  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1968  }
1969  else if((greyVertexNumber == 2)
1970  &&(lowerVertexNumber == 1)){
1971  createdVertices += computeCase4<dataTypeU, dataTypeV>(
1972  polygonEdgeId, tetId,
1973  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1974  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1975  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1976  }
1977  else if((greyVertexNumber == 2)
1978  &&(upperVertexNumber == 1)){
1979  createdVertices += computeCase4<dataTypeU, dataTypeV>(
1980  polygonEdgeId, tetId,
1981  triangleEdges[i][0], t[0], uv[0].first, uv[0].second,
1982  triangleEdges[i][1], t[1], uv[1].first, uv[1].second,
1983  triangleEdges[i][2], t[2], uv[2].first, uv[2].second);
1984  }
1985  }
1986  }
1987 
1988  // in case of 2 triangles, remesh locally
1989  // NOTE: here we just snap vertices together if they are colinear
1990  // the mergeFiberSurfaces() function will take care of removing any
1991  // duplicate
1992  if((triangleEdges[1][0] != -1)&&(createdVertices > 3)){
1993 
1994  vector<int> createdVertexList(createdVertices);
1995  for(int i = 0; i < (int) createdVertices; i++){
1996  createdVertexList[i] =
1997  polygonEdgeVertexLists_[polygonEdgeId]->size() - 1 - i;
1998  }
1999 
2000  vector<bool> snappedVertices(createdVertices, false);
2001 
2002  for(int i = 0; i < createdVertices; i++){
2003 
2004  vector<int> colinearVertices;
2005  if(!snappedVertices[i]){
2006  colinearVertices.push_back(i);
2007  for(int j = 0; j < createdVertices; j++){
2008  if((i != j)&&(!snappedVertices[j])){
2009  // not the same vertex
2010  // not snapped already
2011 
2012  if((*polygonEdgeVertexLists_[polygonEdgeId])[
2013  createdVertexList[i]].t_ ==
2014  (*polygonEdgeVertexLists_[polygonEdgeId])[
2015  createdVertexList[j]].t_){
2016  colinearVertices.push_back(j);
2017  }
2018  }
2019  }
2020  }
2021 
2022  if((colinearVertices.size() == 4)||(colinearVertices.size() == 3)){
2023  // 3 co-linear vertices with 1 duplicate
2024 
2025  // we just need to find the pair of duplicates and snap both of
2026  // them to another vertex
2027  pair<int, int> minPair;
2028  double minDistance = -1;
2029  for(int j = 0; j < (int) colinearVertices.size(); j++){
2030  for(int k = 0; k < (int) colinearVertices.size(); k++){
2031  if(j != k){
2032 
2033  double distance = Geometry::distance(
2034  (*polygonEdgeVertexLists_[polygonEdgeId])[
2035  createdVertexList[colinearVertices[j]]].p_,
2036  (*polygonEdgeVertexLists_[polygonEdgeId])[
2037  createdVertexList[colinearVertices[k]]].p_);
2038 
2039 // bool basePointSnap = true;
2040 // for(int l = 0; l < 3; l++){
2041 // if((*polygonEdgeVertexLists_[polygonEdgeId])[
2042 // createdVertexList[colinearVertices[j]]].p_[l] !=
2043 // (*polygonEdgeVertexLists_[polygonEdgeId])[
2044 // createdVertexList[colinearVertices[k]]].p_[l]){
2045 // basePointSnap = false;
2046 // break;
2047 // }
2048 // }
2049 
2050 // if(!basePointSnap){
2051  if((minDistance < 0)||(distance < minDistance)){
2052  minDistance = distance;
2053  minPair.first = j;
2054  minPair.second = k;
2055  }
2056 // }
2057  }
2058  }
2059  }
2060  if((minDistance != -1)&&(minDistance < pow10(-DBL_DIG))){
2061 // if((minDistance != -1)&&(minDistance < pointSnappingThreshold_)){
2062  // snap them to another colinear vertex
2063  for(int j = 0; j < (int) colinearVertices.size(); j++){
2064  if((j != minPair.first)&&(j != minPair.second)){
2065  // snap minPair.first to j
2066 
2067  for(int k = 0; k < 3; k++){
2068  (*polygonEdgeVertexLists_[polygonEdgeId])[
2069  createdVertexList[colinearVertices[minPair.first]]].p_[k] =
2070  (*polygonEdgeVertexLists_[polygonEdgeId])[
2071  createdVertexList[colinearVertices[j]]].p_[k];
2072  }
2073  (*polygonEdgeVertexLists_[polygonEdgeId])[
2074  createdVertexList[colinearVertices[minPair.first]]].uv_ =
2075  (*polygonEdgeVertexLists_[polygonEdgeId])[
2076  createdVertexList[colinearVertices[j]]].uv_;
2077 
2078  (*polygonEdgeVertexLists_[polygonEdgeId])[
2079  createdVertexList[colinearVertices[minPair.first]]].t_ =
2080  (*polygonEdgeVertexLists_[polygonEdgeId])[
2081  createdVertexList[colinearVertices[j]]].t_;
2082 
2083  (*polygonEdgeVertexLists_[polygonEdgeId])[
2084  createdVertexList[
2085  colinearVertices[minPair.first]]].isBasePoint_ =
2086  (*polygonEdgeVertexLists_[polygonEdgeId])[
2087  createdVertexList[colinearVertices[j]]].isBasePoint_;
2088 
2089  (*polygonEdgeVertexLists_[polygonEdgeId])[
2090  createdVertexList[
2091  colinearVertices[minPair.first]]].isIntersectionPoint_ =
2092  (*polygonEdgeVertexLists_[polygonEdgeId])[
2093  createdVertexList[
2094  colinearVertices[j]]].isIntersectionPoint_;
2095  if((*polygonEdgeVertexLists_[polygonEdgeId])[
2096  createdVertexList[
2097  colinearVertices[j]]].meshEdge_.first != -1){
2098  (*polygonEdgeVertexLists_[polygonEdgeId])[
2099  createdVertexList[
2100  colinearVertices[minPair.first]]].meshEdge_ =
2101  (*polygonEdgeVertexLists_[polygonEdgeId])[
2102  createdVertexList[
2103  colinearVertices[j]]].meshEdge_;
2104  }
2105 
2106 
2107  // snap minPair.second to j
2108  for(int k = 0; k < 3; k++){
2109  (*polygonEdgeVertexLists_[polygonEdgeId])[
2110  createdVertexList[colinearVertices[minPair.second]]].p_[k]
2111  = (*polygonEdgeVertexLists_[polygonEdgeId])[
2112  createdVertexList[colinearVertices[j]]].p_[k];
2113  }
2114  (*polygonEdgeVertexLists_[polygonEdgeId])[
2115  createdVertexList[colinearVertices[minPair.second]]].uv_ =
2116  (*polygonEdgeVertexLists_[polygonEdgeId])[
2117  createdVertexList[colinearVertices[j]]].uv_;
2118 
2119  (*polygonEdgeVertexLists_[polygonEdgeId])[
2120  createdVertexList[colinearVertices[minPair.second]]].t_ =
2121  (*polygonEdgeVertexLists_[polygonEdgeId])[
2122  createdVertexList[colinearVertices[j]]].t_;
2123 
2124  (*polygonEdgeVertexLists_[polygonEdgeId])[
2125  createdVertexList[
2126  colinearVertices[minPair.second]]].isBasePoint_ =
2127  (*polygonEdgeVertexLists_[polygonEdgeId])[
2128  createdVertexList[colinearVertices[j]]].isBasePoint_;
2129 
2130  (*polygonEdgeVertexLists_[polygonEdgeId])[
2131  createdVertexList[
2132  colinearVertices[minPair.second]]].isIntersectionPoint_ =
2133  (*polygonEdgeVertexLists_[polygonEdgeId])[
2134  createdVertexList[
2135  colinearVertices[j]]].isIntersectionPoint_;
2136  if((*polygonEdgeVertexLists_[polygonEdgeId])[
2137  createdVertexList[
2138  colinearVertices[j]]].meshEdge_.first != -1){
2139  (*polygonEdgeVertexLists_[polygonEdgeId])[
2140  createdVertexList[
2141  colinearVertices[minPair.second]]].meshEdge_ =
2142  (*polygonEdgeVertexLists_[polygonEdgeId])[
2143  createdVertexList[
2144  colinearVertices[j]]].meshEdge_;
2145  }
2146 
2147  snappedVertices[colinearVertices[minPair.first]] = true;
2148  snappedVertices[colinearVertices[minPair.second]] = true;
2149 
2150  // we're done
2151  break;
2152  }
2153  }
2154  }
2155  }
2156  }
2157  }
2158 
2159  return createdVertices;
2160  }
2161 
2162  return 0;
2163 }
2164 
2165 template <class dataTypeU, class dataTypeV>
2167 
2168 #ifndef withKamikaze
2169  if(!tetNumber_)
2170  return -1;
2171 #endif
2172 
2173  // Algorithm
2174  // 1) loop over each triangle and mark the containing tetrahedron
2175  // 2) for each tet, in parallel, check for pairwise triangle intersections
2176  // (based on the triangles' range projection)
2177  // Given a pair of intersecting triangles:
2178  // 3) given the point of intersection, take one of its coordinates (u or v)
2179  // and compute an iso-contour I on both triangles
2180  // 4) express the barycentric coordinates of I in both triangles, 2 cases:
2181  // a) the intersection is the entire fiber (old code)
2182  // b) the intersection is a segment of fiber
2183 
2184  // NOTE:
2185  // the topological aspect of the code is OK.
2186  // if any bug, it's very likely to be geometry accuracy related.
2187 
2188  vector<vector<IntersectionTriangle> > tetIntersections(tetNumber_);
2189  vector<vector<Vertex> > tetNewVertices(tetNumber_);
2190 
2191  // fill the information prior to the parallel pass
2192  for(int i = 0; i < (int) polygonEdgeTriangleLists_.size(); i++){
2193  for(int j = 0; j < (int) polygonEdgeTriangleLists_[i]->size(); j++){
2194 
2195  int tetId = (*polygonEdgeTriangleLists_[i])[j].tetId_;
2196 
2197  tetIntersections[tetId].resize(tetIntersections[tetId].size() + 1);
2198 
2199  tetIntersections[tetId].back().caseId_ =
2200  (*polygonEdgeTriangleLists_[i])[j].caseId_;
2201  tetIntersections[tetId].back().polygonEdgeId_ = i;
2202  tetIntersections[tetId].back().triangleId_ = j;
2203  for(int k = 0; k < 3; k++){
2204  tetIntersections[tetId].back().vertexIds_[k] =
2205  (*polygonEdgeTriangleLists_[i])[j].vertexIds_[k];
2206  tetIntersections[tetId].back().uv_[k] =
2207  (*globalVertexList_)[
2208  tetIntersections[tetId].back().vertexIds_[k]].uv_;
2209  tetIntersections[tetId].back().t_[k] =
2210  (*globalVertexList_)[
2211  tetIntersections[tetId].back().vertexIds_[k]].t_;
2212  for(int l = 0; l < 3; l++){
2213  tetIntersections[tetId].back().p_[k][l] =
2214  (*globalVertexList_)[
2215  tetIntersections[tetId].back().vertexIds_[k]].p_[l];
2216  }
2217  }
2218  tetIntersections[tetId].back().intersection_.first = -DBL_MAX;
2219  tetIntersections[tetId].back().intersection_.second = -DBL_MAX;
2220  }
2221  }
2222 
2223  vector<int> tetList;
2224  for(int i = 0; i < (int) tetIntersections.size(); i++){
2225  if(tetIntersections[i].size() > 1)
2226  tetList.push_back(i);
2227  }
2228 
2229 #ifdef withOpenMP
2230 #pragma omp parallel for num_threads(threadNumber_)
2231 #endif
2232  for(int i = 0; i < (int) tetList.size(); i++){
2233 
2234  int tetId = tetList[i];
2235 
2236  // pre-process by merging nearby vertices...?
2237 
2238  int newTriangleNumber = 1;
2239  int newVertexNumber = 1;
2240 
2241  for(int j = 0; j < (int) tetIntersections[tetId].size(); j++){
2242 
2243  if(j > 1000){
2244  stringstream msg;
2245  msg << "[FiberSurface] Preventing an infinite loop!" << endl;
2246  msg << "[FiberSurface] More than 1000 re-meshed triangles in tet #"
2247  << tetId << " :(" << endl;
2248  msg << "[FiberSurface] Extra-thin triangles keep on intersecting?!"
2249  << endl;
2250  dMsg(cerr, msg.str(), infoMsg);
2251  break;
2252  }
2253 
2254  // if we re-mesh, we add new triangles at the end of the list.
2255  // there's no need to check intersections with those.
2256  int originalTriangleNumber = (int) tetIntersections[tetId].size();
2257 
2258  for(int k = 0; k < originalTriangleNumber; k++){
2259 
2260  int polygonEdgeId0 = tetIntersections[tetId][j].polygonEdgeId_;
2261  int polygonEdgeId1 = tetIntersections[tetId][k].polygonEdgeId_;
2262 
2263  if((j != k)&&(polygonEdgeId0 != polygonEdgeId1)){
2264  // cases 3, 4 and 6 of the fiber surface table (multiple triangle
2265  // per tet given a single edge). we don't need to re-mesh that.
2266 
2267  // grab the range projection for the triangle j
2268  pair<double, double> edge0point0, edge0point1;
2269 
2270  getTriangleRangeExtremities(tetId, j, tetIntersections,
2271  edge0point0, edge0point1);
2272 
2273  // now do the same thing for the triangle k
2274  pair<double, double> edge1point0, edge1point1;
2275 
2276  getTriangleRangeExtremities(tetId, k, tetIntersections,
2277  edge1point0, edge1point1);
2278 
2279  // compute the intersection
2280  pair<double, double> intersection;
2281  bool hasIntersection =
2283  edge0point0.first, edge0point0.second,
2284  edge0point1.first, edge0point1.second,
2285  edge1point0.first, edge1point0.second,
2286  edge1point1.first, edge1point1.second,
2287  intersection.first, intersection.second);
2288 
2289  if((hasIntersection)
2290  // check if that intersection has been registered before
2291  // in the end, only one intersection per triangle, no matter what
2292  /*&&(((fabs(tetIntersections[tetId][j].intersection_.first
2293  - intersection.first) > pow(10, -FLT_DIG))
2294  ||(fabs(tetIntersections[tetId][j].intersection_.second
2295  - intersection.second) > pow(10, -FLT_DIG)))
2296  &&((fabs(tetIntersections[tetId][k].intersection_.first
2297  - intersection.first) > pow(10, -FLT_DIG))
2298  ||(fabs(tetIntersections[tetId][k].intersection_.second
2299  - intersection.second) > pow(10, -FLT_DIG))))*/){
2300 
2301  computeTriangleIntersection(tetId, j, k,
2302  polygonEdgeId0, polygonEdgeId1, intersection,
2303  newVertexNumber, newTriangleNumber,
2304  tetIntersections, tetNewVertices);
2305  }
2306  }
2307  }
2308  }
2309  }
2310 
2311  // now copy the new vertices
2312  for(int i = 0; i < (int) tetNewVertices.size(); i++){
2313  for(int j = 0; j < (int) tetNewVertices[i].size(); j++){
2314 
2315  int localId = (*globalVertexList_).size();
2316  tetNewVertices[i][j].localId_ = localId;
2317  (*globalVertexList_).push_back(tetNewVertices[i][j]);
2318  }
2319  }
2320 
2321  for(int i = 0; i < (int) tetIntersections.size(); i++){
2322 
2323  for(int j = 0; j < (int) tetIntersections[i].size(); j++){
2324 
2325  if(((tetIntersections[i][j].intersection_.first != -DBL_MAX)
2326  &&(tetIntersections[i][j].intersection_.second != -DBL_MAX))
2327  ||(tetIntersections[i][j].triangleId_ < 0)){
2328 
2329  int triangleId = tetIntersections[i][j].triangleId_;
2330 
2331  if(triangleId < 0){
2332 
2333  // this is a new triangle
2334  triangleId =
2336  tetIntersections[i][j].polygonEdgeId_]).size();
2338  tetIntersections[i][j].polygonEdgeId_]).resize(triangleId + 1);
2340  tetIntersections[i][j].polygonEdgeId_]).back().tetId_ = i;
2342  tetIntersections[i][j].polygonEdgeId_]).back().caseId_ =
2343  tetIntersections[i][j].caseId_;
2344  }
2345 
2346  for(int k = 0; k < 3; k++){
2347 
2348  int vertexId = tetIntersections[i][j].vertexIds_[k];
2349 
2350  if(vertexId < 0){
2351  // newly created vertex
2352  vertexId = tetNewVertices[i][-(vertexId + 1)].localId_;
2353  }
2355  tetIntersections[i][j].polygonEdgeId_])[triangleId].vertexIds_[k]
2356  = vertexId;
2357  }
2358  }
2359  }
2360  }
2361 
2362  return 0;
2363 }
2364 
2365 #endif // FIBERSURFACE_H
int setTetList(const long long int *tetList)
Definition: FiberSurface.h:152
int polygonEdgeId_
Definition: FiberSurface.h:198
int computeSurface()
Definition: FiberSurface.h:1478
int polygonEdgeId_
Definition: FiberSurface.h:52
bool pointSnapping_
Definition: FiberSurface.h:389
double p_[3]
Definition: FiberSurface.h:43
int computeCase4(const int &polygonEdgeId, const int &tetId, const int &localEdgeId0, const double &t0, const double &u0, const double &v0, const int &localEdgeId1, const double &t1, const double &u1, const double &v1, const int &localEdgeId2, const double &t2, const double &u2, const double &v2) const
Definition: FiberSurface.h:1124
int tetNumber_
Definition: FiberSurface.h:391
int setPointMergingThreshold(const double &threshold)
Definition: FiberSurface.h:124
pair< double, double > intersection_
Definition: FiberSurface.h:204
const void * uField_
Definition: FiberSurface.h:392
int finalize(const bool &mergeDuplicatedVertices=false, const bool &removeSmallEdges=false, const bool &edgeFlips=false, const bool &intersectionRemesh=false)
Definition: FiberSurface.h:1589
static int computeBarycentricCoordinates(const double *p0, const double *p1, const double *p, vector< double > &baryCentrics, const int &dimension=3)
Definition: Geometry.cpp:121
int edgeImplicitEncoding_[12]
Definition: FiberSurface.h:397
pair< double, double > uv_[3]
Definition: FiberSurface.h:201
int polygonEdgeId_
Definition: FiberSurface.h:39
bool isEdgeAngleCollapsible(const int &source, const int &destination, const int &pivotVertexId, const vector< pair< int, int > > &starNeighbors) const
Definition: FiberSurface.cpp:1080
int setPolygonEdgeNumber(const int &polygonEdgeNumber)
Definition: FiberSurface.h:145
const int msg(const char *msg, const int &debugLevel=infoMsg) const
Definition: Debug.cpp:67
static double distance(const double *p0, const double *p1, const int &dimension=3)
Definition: Geometry.cpp:347
int snapVertexBarycentrics(const double &distanceThreshold) const
Definition: FiberSurface.cpp:1636
Definition: FiberSurface.h:194
double p_[3][3]
Definition: FiberSurface.h:203
int vertexIds_[3]
Definition: FiberSurface.h:52
int tetId_
Definition: FiberSurface.h:52
int computeCase3(const int &polygonEdgeId, const int &tetId, const int &localEdgeId0, const double &t0, const double &u0, const double &v0, const int &localEdgeId1, const double &t1, const double &u1, const double &v1, const int &localEdgeId2, const double &t2, const double &u2, const double &v2) const
Definition: FiberSurface.h:981
int setTetNumber(const int &tetNumber)
Definition: FiberSurface.h:162
int vertexIds_[3]
Definition: FiberSurface.h:200
FiberSurface()
Definition: FiberSurface.cpp:88
const long long int * tetList_
Definition: FiberSurface.h:394
int setPointMerging(const bool &onOff)
Definition: FiberSurface.h:119
int setPolygon(const vector< pair< pair< double, double >, pair< double, double > > > *polygon)
Definition: FiberSurface.h:139
bool isBasePoint_
Definition: FiberSurface.h:38
int getNumberOfCommonVertices(const int &tetId, const int &triangleId0, const int &triangleId1, const vector< vector< IntersectionTriangle > > &tetIntersections) const
Definition: FiberSurface.cpp:129
const float * pointSet_
Definition: FiberSurface.h:393
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 mergeVertices(const double &distanceThreshold) const
Definition: FiberSurface.cpp:1374
struct wtfit::FiberSurface::_intersectionTriangle IntersectionTriangle
static bool isTriangleColinear(const double *p0, const double *p1, const double *p2, const double *tolerance=NULL)
Definition: Geometry.cpp:489
Definition: Debug.h:52
int pointNumber_
Definition: FiberSurface.h:391
virtual const int dMsg(ostream &stream, string msg, const int &debugLevel=infoMsg) const
Definition: Debug.cpp:52
int caseId_
Definition: FiberSurface.h:195
int computeCase1(const int &polygonEdgeId, const int &tetId, const int &localEdgeId0, const double &t0, const double &u0, const double &v0, const int &localEdgeId1, const double &t1, const double &u1, const double &v1, const int &localEdgeId2, const double &t2, const double &u2, const double &v2) const
Definition: FiberSurface.h:623
pair< int, int > meshEdge_
Definition: FiberSurface.h:42
int computeCase0(const int &polygonEdgeId, const int &tetId, const int &localEdgeId0, const double &t0, const double &u0, const double &v0, const int &localEdgeId1, const double &t1, const double &u1, const double &v1, const int &localEdgeId2, const double &t2, const double &u2, const double &v2) const
Definition: FiberSurface.h:522
int processTetrahedron(const int &tetId, const pair< double, double > &rangePoint0, const pair< double, double > &rangePoint1, const int &polygonEdgeId=0) const
Definition: FiberSurface.h:1653
double pointSnappingThreshold_
Definition: FiberSurface.h:399
int setTetNeighbors(const vector< vector< int > > *tetNeighbors)
Definition: FiberSurface.h:157
Definition: Debug.h:50
vector< Vertex > * globalVertexList_
Definition: FiberSurface.h:404
double getElapsedTime()
Definition: Os.h:359
FiberSurface processing package.
Definition: FiberSurface.h:31
Definition: Os.h:347
pair< double, double > uv_
Definition: FiberSurface.h:45
int setTriangleList(const int &polygonEdgeId, vector< Triangle > *triangleList)
Definition: FiberSurface.h:167
int getTriangleRangeExtremities(const int &tetId, const int &triangleId, const vector< vector< IntersectionTriangle > > &tetIntersections, pair< double, double > &extremity0, pair< double, double > &extremity1) const
Definition: FiberSurface.cpp:972
const vector< vector< int > > * tetNeighbors_
Definition: FiberSurface.h:396
bool isIntersectionTriangleColinear(const int &tetId, const int &triangleId, const vector< vector< IntersectionTriangle > > &tetIntersections, const vector< vector< Vertex > > &tetNewVertices, const int &vertexId0, const int &vertexId1, const int &vertexId2) const
Definition: FiberSurface.h:340
int triangleId_
Definition: FiberSurface.h:197
int flipEdges() const
Definition: FiberSurface.cpp:712
bool isIntersectionPoint_
Definition: FiberSurface.h:38
int mergeEdges(const double &distanceThreshold) const
Definition: FiberSurface.cpp:1166
int computeTriangleIntersection(const int &tetId, const int &triangleId0, const int &triangleId1, const int &polygonEdgeId0, const int &polygonEdgeId1, const pair< double, double > &intersection, int &newVertexNumber, int &newTriangleNumber, vector< vector< IntersectionTriangle > > &tetIntersections, vector< vector< Vertex > > &tetNewVertices) const
Definition: FiberSurface.cpp:296
int setGlobalVertexList(vector< Vertex > *globalList)
Definition: FiberSurface.h:106
int computeBaseTriangle(const int &tetId, const int &localEdgeId0, const double &t0, const double &u0, const double &v0, const int &localEdgeId1, const double &t1, const double &u1, const double &v1, const int &localEdgeId2, const double &t2, const double &u2, const double &v2, vector< vector< double > > &basePoints, vector< pair< double, double > > &basePointProections, vector< double > &basePointParameterization, vector< pair< int, int > > &baseEdges) const
Definition: FiberSurface.h:438
int setInputField(const void *uField, const void *vField)
Definition: FiberSurface.h:111
double edgeCollapseThreshold_
Definition: FiberSurface.h:399
int setPointSet(const float *pointSet)
Definition: FiberSurface.h:134
int setPointNumber(const int &number)
Definition: FiberSurface.h:129
Minimalist debugging class.
Definition: Debug.h:39
static bool computeSegmentIntersection(const double &xA, const double &yA, const double &xB, const double &yB, const double &xC, const double &yC, const double &xD, const double &yD, double &x, double &y)
Definition: Geometry.cpp:251
int polygonEdgeNumber_
Definition: FiberSurface.h:391
int computeCase2(const int &polygonEdgeId, const int &tetId, const int &localEdgeId0, const double &t0, const double &u0, const double &v0, const int &localEdgeId1, const double &t1, const double &u1, const double &v1, const int &localEdgeId2, const double &t2, const double &u2, const double &v2) const
Definition: FiberSurface.h:806
bool isEdgeFlippable(const int &edgeVertexId0, const int &edgeVertexId1, const int &otherVertexId0, const int &otherVertexId1) const
Definition: FiberSurface.cpp:1127
int caseId_
Definition: FiberSurface.h:52
vector< vector< Vertex > * > polygonEdgeVertexLists_
Definition: FiberSurface.h:406
int setVertexList(const int &polygonEdgeId, vector< Vertex > *vertexList)
Definition: FiberSurface.h:180
~FiberSurface()
Definition: FiberSurface.cpp:125
int createNewIntersectionTriangle(const int &tetId, const int &triangleId, const int &vertexId0, const int &vertexId1, const int &vertexId2, const vector< vector< Vertex > > &tetNewVertices, int &newTriangleNumber, vector< vector< IntersectionTriangle > > &tetIntersections, const pair< double, double > *intersection=NULL) const
Definition: FiberSurface.cpp:641
bool hasDuplicatedVertices(const double *p0, const double *p1, const double *p2) const
Definition: FiberSurface.cpp:1043
double t_[3]
Definition: FiberSurface.h:202
int computeTriangleFiber(const int &tetId, const int &triangleId, const pair< double, double > &intersection, const vector< vector< IntersectionTriangle > > &tetIntersections, vector< double > &pA, vector< double > &pB, int &pivotVertexId, bool &edgeFiber) const
Definition: FiberSurface.cpp:166
const vector< pair< pair< double, double >, pair< double, double > > > * polygon_
Definition: FiberSurface.h:402
Definition: CommandLineParser.h:13
const void * vField_
Definition: FiberSurface.h:392
int interpolateBasePoints(const vector< double > &p0, const pair< double, double > &uv0, const double &t0, const vector< double > &p1, const pair< double, double > &uv1, const double &t1, const double &t, Vertex &v) const
Definition: FiberSurface.cpp:1058
int remeshIntersections() const
Definition: FiberSurface.h:2166
int localId_
Definition: FiberSurface.h:39
double t_
Definition: FiberSurface.h:43
int globalId_
Definition: FiberSurface.h:39
Definition: FiberSurface.h:35
Definition: FiberSurface.h:48
vector< vector< Triangle > * > polygonEdgeTriangleLists_
Definition: FiberSurface.h:408
int snapToBasePoint(const vector< vector< double > > &basePoints, const vector< pair< double, double > > &uv, const vector< double > &t, Vertex &v) const
Definition: FiberSurface.cpp:1605
int threadNumber_
Definition: Debug.h:108
int debugLevel_
Definition: Debug.h:108