gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
meshPolyMesh.h
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #ifndef MESH_POLYMESH_H
7 #define MESH_POLYMESH_H
8 
9 #include <vector>
10 #include <algorithm>
11 #include <stack>
12 #include <stdio.h>
13 #include "SVector3.h"
14 
15 class PolyMesh {
16 public:
17  class HalfEdge;
18  class Face;
19  class Vertex;
20 
21  class Vertex {
22  public:
23  Vertex(double x, double y, double z, int _d = -1)
24  : position(x, y, z), he(NULL), data(_d)
25  {
26  }
28  PolyMesh::HalfEdge *he; // one incident half edge
29  int data;
30  };
31 
32  class HalfEdge {
33  public:
35  : v(vv), f(NULL), prev(NULL), next(NULL), opposite(NULL), data(-1)
36  {
37  }
38  Vertex *v; // origin
39  Face *f; // incident face
40  HalfEdge *prev; // previous half edge on the face
41  HalfEdge *next; // next half edge on the face
42  HalfEdge *opposite; // opposite half edge (twin)
43  int data;
44  SVector3 d() const
45  {
46  SVector3 t = next->v->position - v->position;
47  t.normalize();
48  return t;
49  }
50  };
51 
52  class Face {
53  public:
54  Face(HalfEdge *e) : he(e), data(-1) {}
55  HalfEdge *he; // one half edge of the face
56  int data;
57  };
58 
59  std::vector<Vertex *> vertices;
60  std::vector<HalfEdge *> hedges;
61  std::vector<Face *> faces;
62  std::vector<SVector3> high_order_nodes;
63 
64  void reset()
65  {
66  for(auto it : vertices) delete it;
67  for(auto it : hedges) delete it;
68  for(auto it : faces) delete it;
69  }
70 
71  ~PolyMesh() { reset(); }
72 
73  void print4debug(const int debugTag)
74  {
75  char name[256];
76  sprintf(name, "polyMesh%d.pos", debugTag);
77  FILE *f = fopen(name, "w");
78  fprintf(f, "View \" %s \"{\n", name);
79  for(auto it : faces) {
80  HalfEdge *he0 = it->he;
81  HalfEdge *he1 = it->he->next;
82  HalfEdge *he2 = it->he->next->next;
83  fprintf(f, "ST(%g,%g,0,%g,%g,0,%g,%g,0){%d,%d,%d};\n",
84  he0->v->position.x(), he0->v->position.y(), he1->v->position.x(),
85  he1->v->position.y(), he2->v->position.x(), he2->v->position.y(),
86  it->data, it->data, it->data);
87  }
88  for(auto it : hedges) {
89  HalfEdge *he = it;
90  if(he->data >= 0) {
91  fprintf(f, "SL(%g,%g,0,%g,%g,0){%d,%d};\n", he->v->position.x(),
92  he->v->position.y(), he->opposite->v->position.x(),
93  he->opposite->v->position.y(), he->data, he->data);
94  }
95  }
96 
97  fprintf(f, "};\n");
98  fclose(f);
99  }
100 
101  // compute the degree of a given vertex v
102  inline int degree(const Vertex *v) const
103  {
104  HalfEdge *he = v->he;
105  size_t count = 0;
106  do {
107  he = he->opposite;
108  if(he == NULL) return -1;
109  he = he->next;
110  count++;
111  } while(he != v->he);
112  return count;
113  }
114 
115  inline int num_sides(const HalfEdge *he) const
116  {
117  size_t count = 0;
118  const HalfEdge *start = he;
119  do {
120  count++;
121  he = he->next;
122  } while(he != start);
123  return count;
124  }
125 
126  // compute the normal of an internal vertex v
127  inline SVector3 normal(const Vertex *v) const
128  {
129  SVector3 n(0, 0, 0);
130  HalfEdge *he = v->he;
131  do {
132  SVector3 n1 = he->d();
133  he = he->opposite;
134  if(he == NULL) return -1;
135  he = he->next;
136  n += crossprod(n1, he->d());
137  } while(he != v->he);
138  n.normalize();
139  return n;
140  }
141 
142  inline HalfEdge *getEdge(Vertex *v0, Vertex *v1)
143  {
144  HalfEdge *he = v0->he;
145  do {
146  if(he->next->v == v1) return he;
147  he = he->opposite;
148  if(he == NULL) return NULL;
149  he = he->next;
150  } while(he != v0->he);
151  return NULL;
152  }
153 
154  inline void createFace(Face *f, Vertex *v0, Vertex *v1, Vertex *v2,
155  HalfEdge *he0, HalfEdge *he1, HalfEdge *he2)
156  {
157  he0->v = v0;
158  he1->v = v1;
159  he2->v = v2;
160  v0->he = he0;
161  v1->he = he1;
162  v2->he = he2;
163 
164  he0->next = he1;
165  he1->prev = he0;
166  he1->next = he2;
167  he2->prev = he1;
168  he2->next = he0;
169  he0->prev = he2;
170  he0->f = he1->f = he2->f = f;
171  f->he = he0;
172  }
173 
174  // swap without asking questions
175  //
176  // he1
177  // v2 ------->------ v3
178  // | \ |
179  // | \ he0 | he2
180  // heo2| \ |
181  // | heo0 \ |
182  // | \ |
183  // v1 -----<------- v0
184  // heo1
185  //
186  // he1
187  // --------------
188  // | / |
189  // | he0 / | he2
190  // heo2| / |
191  // | /heo0 |
192  // |/ |
193  // --------------
194  // heo1
195  //
196 
197  inline int swap_edge(HalfEdge *he0)
198  {
199  HalfEdge *heo0 = he0->opposite;
200  if(heo0 == NULL) return -1;
201 
202  HalfEdge *he1 = he0->next;
203  HalfEdge *he2 = he1->next;
204  HalfEdge *heo1 = heo0->next;
205  HalfEdge *heo2 = heo1->next;
206 
207  Vertex *v0 = heo1->v;
208  Vertex *v1 = heo2->v;
209  Vertex *v2 = heo0->v;
210  Vertex *v3 = he2->v;
211 
212  createFace(he0->f, v0, v1, v3, heo1, heo0, he2);
213  createFace(heo2->f, v1, v2, v3, heo2, he1, he0);
214  return 0;
215  }
216 
217  inline int merge_faces(HalfEdge *he)
218  {
219  PolyMesh::HalfEdge *heo = he->opposite;
220 
221  if(heo == nullptr) return -1;
222 
223  PolyMesh::Face *to_delete = heo->f;
224 
225  do {
226  heo->f = he->f;
227  heo = heo->next;
228  } while(heo != he->opposite);
229 
230  he->next->prev = heo->prev;
231  heo->prev->next = he->next;
232  he->prev->next = heo->next;
233  heo->next->prev = he->prev;
234 
235  he->f->he = he->next;
236  he->v->he = heo->next;
237  heo->v->he = he->next;
238 
239  // remove afterwards...
240  he->v = nullptr;
241  heo->v = nullptr;
242  to_delete->he = nullptr;
243  return 0;
244  }
245 
246  void cleanv()
247  {
248  std::vector<Vertex *> uv;
249  for(auto v : vertices) {
250  if(v->he)
251  uv.push_back(v);
252  else
253  delete v;
254  }
255  vertices = uv;
256  }
257 
258  void cleanh()
259  {
260  std::vector<HalfEdge *> uh;
261  for(auto h : hedges) {
262  if(h->f)
263  uh.push_back(h);
264  else
265  delete h;
266  }
267  hedges = uh;
268  }
269 
270  void cleanf()
271  {
272  std::vector<Face *> uf;
273  for(auto f : faces) {
274  if(f->he)
275  uf.push_back(f);
276  else
277  delete f;
278  }
279  faces = uf;
280  }
281 
282  void clean()
283  {
284  cleanv();
285  cleanh();
286  cleanf();
287  }
288 
289  inline int split_edge(HalfEdge *he0m, const SVector3 &position, int data)
290  {
291  HalfEdge *he1m = he0m->opposite;
292  if(he1m == nullptr) return -1;
293 
294  Vertex *mid = new Vertex(position.x(), position.y(), position.z(), data);
295  vertices.push_back(mid);
296 
297  HalfEdge *he12 = he0m->next;
298  HalfEdge *he20 = he0m->next->next;
299  HalfEdge *he03 = he0m->opposite->next;
300  HalfEdge *he31 = he0m->opposite->next->next;
301 
302  // if(he03->v != he0m->v) Msg::Error("error 1");
303  // if(he1m->v != he12->v) Msg::Error("error 2");
304 
305  Vertex *v0 = he03->v;
306  Vertex *v1 = he12->v;
307  Vertex *v2 = he20->v;
308  Vertex *v3 = he31->v;
309 
310  HalfEdge *hem0 = new HalfEdge(mid);
311  HalfEdge *hem1 = new HalfEdge(mid);
312  HalfEdge *hem2 = new HalfEdge(mid);
313  HalfEdge *hem3 = new HalfEdge(mid);
314 
315  HalfEdge *he2m = new HalfEdge(v2);
316  HalfEdge *he3m = new HalfEdge(v3);
317 
318  he0m->opposite = hem0;
319  hem0->opposite = he0m;
320  he1m->opposite = hem1;
321  hem1->opposite = he1m;
322  he2m->opposite = hem2;
323  hem2->opposite = he2m;
324  he3m->opposite = hem3;
325  hem3->opposite = he3m;
326 
327  hedges.push_back(hem0);
328  hedges.push_back(hem1);
329  hedges.push_back(hem2);
330  hedges.push_back(hem3);
331  hedges.push_back(he2m);
332  hedges.push_back(he3m);
333 
334  Face *f0m2 = he0m->f;
335  Face *f1m3 = he1m->f;
336  Face *f2m1 = new Face(he2m);
337  Face *f3m0 = new Face(he3m);
338  faces.push_back(f2m1);
339  faces.push_back(f3m0);
340 
341  createFace(f0m2, v0, mid, v2, he0m, hem2, he20);
342  createFace(f1m3, v1, mid, v3, he1m, hem3, he31);
343  createFace(f2m1, v2, mid, v1, he2m, hem1, he12);
344  createFace(f3m0, v3, mid, v0, he3m, hem0, he03);
345  return 0;
346  }
347 
348  //
349  // v0 he0
350  // ------------------>------ v1
351  // | /
352  // | /
353  // | v /
354  // |he2 /
355  // | / he1
356  // | /
357  // | /
358  // |/
359  // v2
360 
361  inline void initialize_rectangle(double xmin, double xmax, double ymin,
362  double ymax)
363  {
364  reset();
365  Vertex *v_mm = new PolyMesh::Vertex(xmin, ymin, 0);
366  vertices.push_back(v_mm);
367  Vertex *v_mM = new PolyMesh::Vertex(xmin, ymax, 0);
368  vertices.push_back(v_mM);
369  Vertex *v_MM = new PolyMesh::Vertex(xmax, ymax, 0);
370  vertices.push_back(v_MM);
371  Vertex *v_Mm = new PolyMesh::Vertex(xmax, ymin, 0);
372  vertices.push_back(v_Mm);
373  HalfEdge *mm_MM = new HalfEdge(v_mm);
374  HalfEdge *MM_Mm = new HalfEdge(v_MM);
375  HalfEdge *Mm_mm = new HalfEdge(v_Mm);
376  hedges.push_back(mm_MM);
377  hedges.push_back(MM_Mm);
378  hedges.push_back(Mm_mm);
379  Face *f0 = new Face(mm_MM);
380  faces.push_back(f0);
381  createFace(f0, v_mm, v_MM, v_Mm, mm_MM, MM_Mm, Mm_mm);
382 
383  HalfEdge *MM_mm = new HalfEdge(v_MM);
384  HalfEdge *mm_mM = new HalfEdge(v_mm);
385  HalfEdge *mM_MM = new HalfEdge(v_mM);
386  hedges.push_back(MM_mm);
387  hedges.push_back(mm_mM);
388  hedges.push_back(mM_MM);
389  Face *f1 = new Face(MM_mm);
390  faces.push_back(f1);
391  createFace(f1, v_MM, v_mm, v_mM, MM_mm, mm_mM, mM_MM);
392 
393  MM_mm->opposite = mm_MM;
394  mm_MM->opposite = MM_mm;
395  }
396 
397  inline int split_triangle(int index, double x, double y, double z, Face *f,
398  int (*doSwap)(PolyMesh::HalfEdge *, void *) = NULL,
399  void *data = NULL,
400  std::vector<HalfEdge *> *_t = NULL)
401  {
402  Vertex *v = new PolyMesh::Vertex(x, y, z); // one more vertex
403  v->data = -1;
404 
405  vertices.push_back(v);
406 
407  HalfEdge *he0 = f->he;
408  HalfEdge *he1 = he0->next;
409  HalfEdge *he2 = he1->next;
410 
411  Vertex *v0 = he0->v;
412  Vertex *v1 = he1->v;
413  Vertex *v2 = he2->v;
414  HalfEdge *hev0 = new PolyMesh::HalfEdge(v);
415  HalfEdge *hev1 = new PolyMesh::HalfEdge(v);
416  HalfEdge *hev2 = new PolyMesh::HalfEdge(v);
417 
418  HalfEdge *he0v = new PolyMesh::HalfEdge(v0);
419  HalfEdge *he1v = new PolyMesh::HalfEdge(v1);
420  HalfEdge *he2v = new PolyMesh::HalfEdge(v2);
421 
422  hedges.push_back(hev0);
423  hedges.push_back(hev1);
424  hedges.push_back(hev2);
425  hedges.push_back(he0v);
426  hedges.push_back(he1v);
427  hedges.push_back(he2v);
428 
429  hev0->opposite = he0v;
430  he0v->opposite = hev0;
431  hev1->opposite = he1v;
432  he1v->opposite = hev1;
433  hev2->opposite = he2v;
434  he2v->opposite = hev2;
435 
436  Face *f0 = f;
437  f->he = hev0;
438  Face *f1 = new Face(hev1);
439  Face *f2 = new Face(hev2);
440  f1->data = f2->data = f0->data;
441 
442  faces.push_back(f1);
443  faces.push_back(f2);
444 
445  createFace(f0, v0, v1, v, he0, he1v, hev0);
446  createFace(f1, v1, v2, v, he1, he2v, hev1);
447  createFace(f2, v2, v0, v, he2, he0v, hev2);
448 
449  if(doSwap) {
450  std::stack<HalfEdge *> _stack;
451  _stack.push(he0);
452  _stack.push(he1);
453  _stack.push(he2);
454  std::vector<HalfEdge *> _touched;
455  while(!_stack.empty()) {
456  HalfEdge *he = _stack.top();
457  _touched.push_back(he);
458  _stack.pop();
459  // printf("do we swap %g %g --> %g %g ?\n",
460  // he->v->position.x(),he->v->position.y(),
461  // he->next->v->position.x(),he->next->v->position.y());
462  if(doSwap(he, data) == 1) {
463  // printf("YES\n");
464  swap_edge(he);
465 
466  HalfEdge *H[2] = {he, he->opposite};
467 
468  for(int k = 0; k < 2; k++) {
469  if(H[k] == NULL) continue;
470  HalfEdge *heb = H[k]->next;
471  HalfEdge *hebo = heb->opposite;
472 
473  if(std::find(_touched.begin(), _touched.end(), heb) ==
474  _touched.end() &&
475  std::find(_touched.begin(), _touched.end(), hebo) ==
476  _touched.end()) {
477  _stack.push(heb);
478  }
479 
480  HalfEdge *hec = heb->next;
481  HalfEdge *heco = hec->opposite;
482 
483  if(std::find(_touched.begin(), _touched.end(), hec) ==
484  _touched.end() &&
485  std::find(_touched.begin(), _touched.end(), heco) ==
486  _touched.end()) {
487  _stack.push(hec);
488  }
489  }
490  }
491  }
492  if(_t) *_t = _touched;
493  }
494  return 0;
495  }
496 };
497 
500  {
501  PolyMesh::Vertex *l10 = std::min(l1->v, l1->next->v);
502  PolyMesh::Vertex *l11 = std::max(l1->v, l1->next->v);
503  PolyMesh::Vertex *l20 = std::min(l2->v, l2->next->v);
504  PolyMesh::Vertex *l21 = std::max(l2->v, l2->next->v);
505  if(l10 < l20) return true;
506  if(l10 > l20) return false;
507  if(l11 > l21) return true;
508  return false;
509  }
510 };
511 
514  {
515  PolyMesh::Vertex *l10 = std::min(l1->v, l1->next->v);
516  PolyMesh::Vertex *l11 = std::max(l1->v, l1->next->v);
517  PolyMesh::Vertex *l20 = std::min(l2->v, l2->next->v);
518  PolyMesh::Vertex *l21 = std::max(l2->v, l2->next->v);
519  if(l10 == l20 && l11 == l21) return true;
520  return false;
521  }
522 };
523 
524 #endif
crossprod
SVector3 crossprod(const SVector3 &a, const SVector3 &b)
Definition: SVector3.h:150
PolyMesh::normal
SVector3 normal(const Vertex *v) const
Definition: meshPolyMesh.h:127
PolyMesh::num_sides
int num_sides(const HalfEdge *he) const
Definition: meshPolyMesh.h:115
PolyMesh::Face::data
int data
Definition: meshPolyMesh.h:56
PolyMesh::swap_edge
int swap_edge(HalfEdge *he0)
Definition: meshPolyMesh.h:197
PolyMesh::reset
void reset()
Definition: meshPolyMesh.h:64
PolyMesh::HalfEdge::v
Vertex * v
Definition: meshPolyMesh.h:38
PolyMesh::Face
Definition: meshPolyMesh.h:52
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
PolyMesh::HalfEdge::d
SVector3 d() const
Definition: meshPolyMesh.h:44
SVector3
Definition: SVector3.h:16
PolyMesh::clean
void clean()
Definition: meshPolyMesh.h:282
SVector3.h
HalfEdgePtrEqual
Definition: meshPolyMesh.h:512
PolyMesh::HalfEdge::next
HalfEdge * next
Definition: meshPolyMesh.h:41
PolyMesh::Vertex::Vertex
Vertex(double x, double y, double z, int _d=-1)
Definition: meshPolyMesh.h:23
PolyMesh::cleanv
void cleanv()
Definition: meshPolyMesh.h:246
PolyMesh::Vertex
Definition: meshPolyMesh.h:21
HalfEdgePtrLessThan
Definition: meshPolyMesh.h:498
PolyMesh::Face::Face
Face(HalfEdge *e)
Definition: meshPolyMesh.h:54
PolyMesh::faces
std::vector< Face * > faces
Definition: meshPolyMesh.h:61
PolyMesh::getEdge
HalfEdge * getEdge(Vertex *v0, Vertex *v1)
Definition: meshPolyMesh.h:142
PolyMesh::cleanh
void cleanh()
Definition: meshPolyMesh.h:258
PolyMesh::initialize_rectangle
void initialize_rectangle(double xmin, double xmax, double ymin, double ymax)
Definition: meshPolyMesh.h:361
SVector3::x
double x(void) const
Definition: SVector3.h:30
PolyMesh::Vertex::data
int data
Definition: meshPolyMesh.h:29
PolyMesh::HalfEdge::prev
HalfEdge * prev
Definition: meshPolyMesh.h:40
PolyMesh::Face::he
HalfEdge * he
Definition: meshPolyMesh.h:55
PolyMesh::hedges
std::vector< HalfEdge * > hedges
Definition: meshPolyMesh.h:60
HalfEdgePtrLessThan::operator()
bool operator()(PolyMesh::HalfEdge *l1, PolyMesh::HalfEdge *l2) const
Definition: meshPolyMesh.h:499
PolyMesh::~PolyMesh
~PolyMesh()
Definition: meshPolyMesh.h:71
SVector3::y
double y(void) const
Definition: SVector3.h:31
z
const double z
Definition: GaussQuadratureQuad.cpp:56
PolyMesh::HalfEdge::HalfEdge
HalfEdge(Vertex *vv)
Definition: meshPolyMesh.h:34
HalfEdgePtrEqual::operator()
bool operator()(PolyMesh::HalfEdge *l1, PolyMesh::HalfEdge *l2) const
Definition: meshPolyMesh.h:513
PolyMesh::print4debug
void print4debug(const int debugTag)
Definition: meshPolyMesh.h:73
SVector3::z
double z(void) const
Definition: SVector3.h:32
PolyMesh::HalfEdge::data
int data
Definition: meshPolyMesh.h:43
PolyMesh::HalfEdge::opposite
HalfEdge * opposite
Definition: meshPolyMesh.h:42
PolyMesh::split_triangle
int split_triangle(int index, double x, double y, double z, Face *f, int(*doSwap)(PolyMesh::HalfEdge *, void *)=NULL, void *data=NULL, std::vector< HalfEdge * > *_t=NULL)
Definition: meshPolyMesh.h:397
PolyMesh::Vertex::he
PolyMesh::HalfEdge * he
Definition: meshPolyMesh.h:28
PolyMesh::high_order_nodes
std::vector< SVector3 > high_order_nodes
Definition: meshPolyMesh.h:62
PolyMesh::createFace
void createFace(Face *f, Vertex *v0, Vertex *v1, Vertex *v2, HalfEdge *he0, HalfEdge *he1, HalfEdge *he2)
Definition: meshPolyMesh.h:154
PolyMesh::cleanf
void cleanf()
Definition: meshPolyMesh.h:270
PolyMesh::merge_faces
int merge_faces(HalfEdge *he)
Definition: meshPolyMesh.h:217
PolyMesh
Definition: meshPolyMesh.h:15
PolyMesh::degree
int degree(const Vertex *v) const
Definition: meshPolyMesh.h:102
PolyMesh::split_edge
int split_edge(HalfEdge *he0m, const SVector3 &position, int data)
Definition: meshPolyMesh.h:289
PolyMesh::Vertex::position
SVector3 position
Definition: meshPolyMesh.h:27
PolyMesh::vertices
std::vector< Vertex * > vertices
Definition: meshPolyMesh.h:59
PolyMesh::HalfEdge
Definition: meshPolyMesh.h:32
SVector3::normalize
double normalize()
Definition: SVector3.h:38
PolyMesh::HalfEdge::f
Face * f
Definition: meshPolyMesh.h:39