gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
DivideAndConquer.cpp
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 // Triangulation using a divide and conquer algorithm
7 //
8 // The main idea is to be able to merge two Delaunay triangulations
9 // into a single one (see the 'Merge' function). Points are then
10 // separated into two groups, then each group into two, ... until
11 // having 1, 2 or 3 points (the triangulation is then trivial).
12 //
13 // This is used to contruct the initial mesh.
14 //
15 // Warning: point positions must be PERTURBED by a small random
16 // value to avoid 3 aligned points or 4 cocyclical points!
17 
18 #include <stdexcept>
19 #include "GmshMessage.h"
20 #include "DivideAndConquer.h"
21 #include "Numeric.h"
22 #include "fullMatrix.h"
23 #include "robustPredicates.h"
24 #include "BackgroundMeshTools.h"
25 #include "OS.h"
26 #include "GPoint.h"
27 #include "GFace.h"
28 #include "MLine.h"
29 
30 #define Pred(x) ((x)->prev)
31 #define Succ(x) ((x)->next)
32 
34 {
35  DListPeek p = points[a].adjacent;
36 
37  do {
38  if(p == nullptr) return -1;
39  if(p->point_num == b) return Pred(p)->point_num;
40  p = Pred(p);
41  } while(p != points[a].adjacent);
42 
43  return -1;
44 }
45 
47 {
48  DListPeek p = points[a].adjacent;
49 
50  do {
51  if(p == nullptr) return -1;
52  if(p->point_num == b) return Succ(p)->point_num;
53  p = Succ(p);
54  } while(p != points[a].adjacent);
55 
56  return -1;
57 }
58 
60 {
62  if(p == nullptr) return 0;
63 
64  int out = 0;
65  DListPeek copy = p;
66  do {
67  if(p->point_num == f) {
68  points[x].adjacent = p;
69  out = 1;
70  }
71  else
72  p = p->next;
73  } while((p != copy) && !out);
74  return out;
75 }
76 
78 {
79  return (points[x].adjacent)->point_num;
80 }
81 
83 {
84  double pa[2] = {(double)points[x].where.h, (double)points[x].where.v};
85  double pb[2] = {(double)points[y].where.h, (double)points[y].where.v};
86  double pc[2] = {(double)points[check].where.h, (double)points[check].where.v};
87 
88  // we use robust predicates here
89  double result = robustPredicates::orient2d(pa, pb, pc);
90 
91  return result > 0;
92 }
93 
95 {
96  return IsLeftOf(y, x, check);
97 }
98 
100 {
101  PointNumero x, y, z, z1, z2, temp;
102  Segment s;
103 
104  x = vl.end; // vu le tri, c'est le point le + a droite
105  y = vr.begin; // idem, le + a gauche
106  z = First(y);
107  z1 = First(x);
108  z2 = Predecessor(x, z1);
109  for(;;) {
110  if(IsRightOf(x, y, z)) {
111  temp = z;
112  z = Successor(z, y);
113  y = temp;
114  }
115  else if(IsRightOf(x, y, z2)) {
116  temp = z2;
117  z2 = Predecessor(z2, x);
118  x = temp;
119  }
120  else {
121  s.from = x;
122  s.to = y;
123  return s;
124  }
125  }
126 }
127 
129 {
130  PointNumero x, y, z, z1, z2, temp;
131  Segment s;
132 
133  x = vl.end; // vu le tri, c'est le point le + a droite
134  y = vr.begin; // idem, le + a gauche
135  z = First(y);
136  z1 = First(x);
137  z2 = Predecessor(y, z);
138  for(;;) {
139  if(IsLeftOf(x, y, z2)) {
140  temp = z2;
141  z2 = Predecessor(z2, y);
142  y = temp;
143  }
144  else if(IsLeftOf(x, y, z1)) {
145  temp = z1;
146  z1 = Successor(z1, x);
147  x = temp;
148  }
149  else {
150  s.from = x;
151  s.to = y;
152  return s;
153  }
154  }
155 }
156 
157 // return 1 if the point k is NOT in the circumcircle of triangle hij
159 {
160  if((h == i) && (h == j) && (h == k)) {
161  throw std::runtime_error("Identical points in triangulation: "
162  "increase element size or Mesh.RandomFactor");
163  return 0;
164  }
165 
166  double pa[2] = {(double)points[h].where.h, (double)points[h].where.v};
167  double pb[2] = {(double)points[i].where.h, (double)points[i].where.v};
168  double pc[2] = {(double)points[j].where.h, (double)points[j].where.v};
169  double pd[2] = {(double)points[k].where.h, (double)points[k].where.v};
170 
171  // we use robust predicates here
172  double result = robustPredicates::incircle(pa, pb, pc, pd) *
173  robustPredicates::orient2d(pa, pb, pc);
174 
175  return (result < 0) ? 1 : 0;
176 }
177 
178 int DocRecord::Merge(DT vl, DT vr)
179 {
180  Segment bt, ut;
181  int a, b, out;
182  PointNumero r, r1, r2, l, l1, l2;
183 
184  bt = LowerCommonTangent(vl, vr);
185  ut = UpperCommonTangent(vl, vr);
186  l = bt.from; // left endpoint of BT
187  r = bt.to; // right endpoint of BT
188 
189  while((l != ut.from) || (r != ut.to)) {
190  a = b = 0;
191  if(!Insert(l, r)) return 0;
192 
193  r1 = Predecessor(r, l);
194  if(r1 == -1) return 0;
195  if(IsRightOf(l, r, r1))
196  a = 1;
197  else {
198  out = 0;
199  while(!out) {
200  r2 = Predecessor(r, r1);
201  if(r2 == -1) return 0;
202  if(r2 < vr.begin)
203  out = 1;
204  else if(Qtest(l, r, r1, r2))
205  out = 1;
206  else {
207  if(!Delete(r, r1)) return 0;
208  r1 = r2;
209  if(IsRightOf(l, r, r1)) out = a = 1;
210  }
211  }
212  }
213 
214  l1 = Successor(l, r);
215  if(l1 == -1) return 0;
216  if(IsLeftOf(r, l, l1))
217  b = 1;
218  else {
219  out = 0;
220  while(!out) {
221  l2 = Successor(l, l1);
222  if(l2 == -1) return 0;
223  if(l2 > vl.end)
224  out = 1;
225  else if(Qtest(r, l, l1, l2))
226  out = 1;
227  else {
228  if(!Delete(l, l1)) return 0;
229  l1 = l2;
230  if(IsLeftOf(r, l, l1)) out = b = 1;
231  }
232  }
233  }
234 
235  if(a)
236  l = l1;
237  else if(b)
238  r = r1;
239  else {
240  if(Qtest(l, r, r1, l1))
241  r = r1;
242  else
243  l = l1;
244  }
245  }
246  if(!Insert(l, r)) return 0;
247 
248  if(!FixFirst(ut.to, ut.from)) return 0;
249  if(!FixFirst(bt.from, bt.to)) return 0;
250  return 1;
251 }
252 
254 {
255  int n, m;
256  DT dt;
257 
258  dt.begin = left;
259  dt.end = right;
260  n = right - left + 1; // nombre de points a triangulariser
261  switch(n) {
262  case 0:
263  case 1:
264  // 0 ou 1 points -> rien a faire
265  break;
266 
267  case 2: // deux points : cas trivial
268  Insert(left, right);
269  FixFirst(left, right);
270  FixFirst(right, left);
271  break;
272 
273  case 3: // trois points : cas trivial
274  Insert(left, right);
275  Insert(left, left + 1);
276  Insert(left + 1, right);
277  if(IsRightOf(left, right, left + 1)) {
278  FixFirst(left, left + 1);
279  FixFirst(left + 1, right);
280  FixFirst(right, left);
281  }
282  else {
283  FixFirst(left, right);
284  FixFirst(left + 1, left);
285  FixFirst(right, left + 1);
286  }
287  break;
288 
289  default: // plus de trois points : cas recursif
290  m = (left + right) >> 1;
291  if(!Merge(RecurTrig(left, m), RecurTrig(m + 1, right))) break;
292  }
293  return dt;
294 }
295 
296 static int comparePoints(const void *i, const void *j)
297 {
298  double x, y;
299 
300  x = ((PointRecord *)i)->where.h - ((PointRecord *)j)->where.h;
301  if(x == 0.) {
302  y = ((PointRecord *)i)->where.v - ((PointRecord *)j)->where.v;
303  return (y < 0.) ? -1 : 1;
304  }
305  else
306  return (x < 0.) ? -1 : 1;
307 }
308 
309 // this fonction builds the delaunay triangulation for a window
311 {
312  qsort(points, numPoints, sizeof(PointRecord), comparePoints);
313  RecurTrig(0, numPoints - 1);
314  return 1;
315 }
316 
317 // This routine insert the point 'newPoint' in the list dlist,
318 // respecting the clock-wise orientation
320 {
321  DListRecord *p, *newp;
322  double alpha1, alpha, beta, xx, yy;
323  int first;
324 
325  newp = new DListRecord;
326  newp->point_num = newPoint;
327 
328  DListRecord **dlist = &points[centerPoint].adjacent;
329  if(*dlist == nullptr) {
330  *dlist = newp;
331  Pred(*dlist) = newp;
332  Succ(*dlist) = newp;
333  return 1;
334  }
335  if(Succ(*dlist) == *dlist) {
336  Pred(*dlist) = newp;
337  Succ(*dlist) = newp;
338  Pred(newp) = *dlist;
339  Succ(newp) = *dlist;
340  return 1;
341  }
342 
343  // If we are here, the double-linked circular list has 2 or more
344  // elements, so we have to calculate where to put the new one
345 
346  p = *dlist;
347  first = p->point_num;
348 
349  DPoint center = points[centerPoint].where;
350 
351  // first, compute polar coord. of the first point
352  yy = (double)(points[first].where.v - center.v);
353  xx = (double)(points[first].where.h - center.h);
354  alpha1 = atan2(yy, xx);
355 
356  // compute polar coord of the point to insert
357  yy = (double)(points[newPoint].where.v - center.v);
358  xx = (double)(points[newPoint].where.h - center.h);
359  beta = atan2(yy, xx) - alpha1;
360  if(beta <= 0) beta += 2. * M_PI;
361 
362  do {
363  yy = (double)(points[Succ(p)->point_num].where.v - center.v);
364  xx = (double)(points[Succ(p)->point_num].where.h - center.h);
365  alpha = atan2(yy, xx) - alpha1;
366  if(alpha <= -1.e-15 || Succ(p)->point_num == first)
367  alpha += 2. * M_PI;
368  else if(abs(alpha) <= 1e-15 &&
369  IsRightOf(centerPoint, first, Succ(p)->point_num))
370  alpha += 2. * M_PI;
371  if(alpha >= beta + 1e-15) {
372  Succ(newp) = Succ(p);
373  Succ(p) = newp;
374  Pred(newp) = p;
375  Pred(Succ(newp)) = newp;
376  return 1;
377  }
378  else if(alpha >= beta - 1e-15) {
379  // Angles more or less the same. To keep consistency with rest of
380  // algorithm check with IsLeftOf to see which point should be first
381  if(IsRightOf(centerPoint, Succ(p)->point_num, newPoint)) {
382  // New point should become before Succ(p)
383  Succ(newp) = Succ(p);
384  Succ(p) = newp;
385  Pred(newp) = p;
386  Pred(Succ(newp)) = newp;
387  }
388  else {
389  // New point should become after Succ(p)
390  Succ(newp) = Succ(Succ(p));
391  Succ(Succ(p)) = newp;
392  Pred(newp) = Succ(p);
393  Pred(Succ(newp)) = newp;
394  }
395  return 1;
396  }
397  p = Succ(p);
398  } while(p != *dlist);
399 
400  // never here!
401  return 0;
402 }
403 
404 // This routine inserts the point 'a' in the adjency list of 'b' and
405 // the point 'b' in the adjency list of 'a'
407 {
408  int rslt = DListInsert(a, b);
409  rslt &= DListInsert(b, a);
410  return rslt;
411 }
412 
414 {
415  DListPeek p;
416 
417  if(*dlist == nullptr) return 0;
418  if(Succ(*dlist) == *dlist) {
419  if((*dlist)->point_num == oldPoint) {
420  delete *dlist;
421  *dlist = nullptr;
422  return 1;
423  }
424  else
425  return 0;
426  }
427  p = *dlist;
428  do {
429  if(p->point_num == oldPoint) {
430  Succ(Pred(p)) = Succ(p);
431  Pred(Succ(p)) = Pred(p);
432  if(p == *dlist) { *dlist = Succ(p); }
433  delete p;
434  return 1;
435  }
436  p = Succ(p);
437  } while(p != *dlist);
438 
439  return 0;
440 }
441 
442 // This routine removes the point 'a' in the adjency list of 'b' and
443 // the point 'b' in the adjency list of 'a'
445 {
446  int rslt = DListDelete(&points[a].adjacent, b);
447  rslt &= DListDelete(&points[b].adjacent, a);
448  return rslt;
449 }
450 
451 // compte les points sur le polygone convexe
453 {
454  PointNumero p, p2, temp;
455  int i, n = numPoints;
456 
457  if(points[0].adjacent == nullptr) return 0;
458  i = 1;
459  p = 0;
460  p2 = First(0);
461  while((p2 != 0) && (i < n)) {
462  i++;
463  temp = p2;
464  p2 = Successor(p2, p);
465  p = temp;
466  }
467  return (i <= n) ? i : -1;
468 }
469 
470 // compte les points sur le polygone convexe
472 {
473  PointNumero p, p2, temp;
474 
475  if(points[0].adjacent == nullptr) return;
476  int count = 0;
477  p = 0;
478  _hull[count++] = p;
479  p2 = First(0);
480  while((p2 != 0) && (count < numPoints)) {
481  temp = p2;
482  p2 = Successor(p2, p);
483  p = temp;
484  _hull[count++] = p;
485  }
486 }
487 
489 {
490  DListPeek p, temp;
491  int i, max = 0;
492  PointNumero *ptr;
493 
494  p = *dlist;
495  do {
496  max++;
497  p = Pred(p);
498  } while(p != *dlist);
499  ptr = new PointNumero[max + 1];
500  if(ptr == nullptr) return nullptr;
501  p = *dlist;
502  for(i = 0; i < max; i++) {
503  ptr[i] = p->point_num;
504  temp = p;
505  p = Pred(p);
506  delete temp;
507  }
508  ptr[max] = ptr[0];
509  *dlist = nullptr;
510  *n = max;
511  return ptr;
512 }
513 
514 // build ready to use voronoi data
516 {
517  if(_adjacencies) {
518  for(int i = 0; i < numPoints; i++)
519  if(_adjacencies[i].t) delete[] _adjacencies[i].t;
520  delete[] _adjacencies;
521  }
522  if(_hull) delete[] _hull;
524  _hull = new PointNumero[_hullSize];
525  ConvexHull();
526  std::sort(_hull, _hull + _hullSize);
528 
529  for(PointNumero i = 0; i < numPoints; i++)
530  _adjacencies[i].t =
531  ConvertDlistToArray(&points[i].adjacent, &_adjacencies[i].t_length);
532 }
533 
534 void DocRecord::voronoiCell(PointNumero pt, std::vector<SPoint2> &pts) const
535 {
536  if(!_adjacencies) {
537  Msg::Error("No adjacencies were created");
538  return;
539  }
540  const int n = _adjacencies[pt].t_length;
541  for(int j = 0; j < n; j++) {
542  PointNumero a = _adjacencies[pt].t[j];
543  PointNumero b = _adjacencies[pt].t[(j + 1) % n];
544  double pa[2] = {(double)points[a].where.h, (double)points[a].where.v};
545  double pb[2] = {(double)points[b].where.h, (double)points[b].where.v};
546  double pc[2] = {(double)points[pt].where.h, (double)points[pt].where.v};
547  double center[2];
548  circumCenterXY(pa, pb, pc, center);
549  pts.push_back(SPoint2(center[0], center[1]));
550  }
551 }
552 
553 /*
554  Consider N points {X_1, \dots, X_N}. We want to find the point X_P that
555  verifies
556 
557  min max | X_i - X_P | , i=1,\dots,N
558 
559  L2 norm
560 
561  min \int_V || X - X_P||^2 dv
562 
563  => 2 \int_V || X - X_P|| dv = 0 => X_P is the centroid
564 
565  min \int_V || X - X_P||^{2m} dv
566 
567  => 2 \int_V || X - X_P||^{2m-1} dv = 0 => X_P is somewhere ...
568 
569  now in infinite norm, how to find X_P ?
570 */
571 
572 void DocRecord::makePosView(const std::string &fileName, GFace *gf)
573 {
574  FILE *f = Fopen(fileName.c_str(), "w");
575  if(!f) {
576  Msg::Error("Could not open file '%s'", fileName.c_str());
577  return;
578  }
579  if(_adjacencies) {
580  fprintf(f, "View \"voronoi\" {\n");
581  for(PointNumero i = 0; i < numPoints; i++) {
582  std::vector<SPoint2> pts;
583  double pc[2] = {(double)points[i].where.h, (double)points[i].where.v};
584  if(!onHull(i)) {
585  GPoint p0(pc[0], pc[1], 0.0);
586  // if (gf) p0 = gf->point(pc[0], pc[1]);
587  fprintf(f, "SP(%g,%g,%g){%g};\n", p0.x(), p0.y(), p0.z(), (double)i);
588  voronoiCell(i, pts);
589  for(std::size_t j = 0; j < pts.size(); j++) {
590  SPoint2 pp1 = pts[j];
591  SPoint2 pp2 = pts[(j + 1) % pts.size()];
592  GPoint p1(pp1.x(), pp1.y(), 0.0);
593  GPoint p2(pp2.x(), pp2.y(), 0.0);
594  // if (gf) {
595  // p1 = gf->point(p1.x(), p1.y());
596  // p2 = gf->point(p2.x(), p2.y());
597  // }
598  fprintf(f, "SL(%g,%g,%g,%g,%g,%g){%g,%g};\n", p1.x(), p1.y(), p1.z(),
599  p2.x(), p2.y(), p2.z(), (double)i, (double)i);
600  }
601  }
602  else {
603  fprintf(f, "SP(%g,%g,%g){%g};\n", pc[0], pc[1], 0., -(double)i);
604  }
605  }
606  fprintf(f, "};\n");
607  }
608  else {
609  fprintf(f, "View \"scalar\" {\n");
610  for(PointNumero iPoint = 0; iPoint < numPoints; iPoint++) {
611  DListPeek p = points[iPoint].adjacent;
612  if(!p) { continue; }
613  std::vector<PointNumero> adjacentPoints;
614  do {
615  adjacentPoints.push_back(p->point_num);
616  p = Pred(p);
617  } while(p != points[iPoint].adjacent);
618  adjacentPoints.push_back(p->point_num);
619 
620  for(size_t iTriangle = 0; iTriangle < adjacentPoints.size() - 1;
621  iTriangle++) {
622  const PointNumero jPoint = adjacentPoints[iTriangle];
623  const PointNumero kPoint = adjacentPoints[iTriangle + 1];
624  if((jPoint > iPoint) && (kPoint > iPoint) &&
625  (IsRightOf(iPoint, jPoint, kPoint))) {
626  fprintf(f, "ST(%g,%g,%g,%g,%g,%g,%g,%g,%g){%g,%g,%g};\n",
627  points[iPoint].where.h, points[iPoint].where.v, 0.0,
628  points[jPoint].where.h, points[jPoint].where.v, 0.0,
629  points[kPoint].where.h, points[kPoint].where.v, 0.0,
630  (double)iPoint, (double)jPoint, (double)kPoint);
631  }
632  }
633  }
634  fprintf(f, "};\n");
635  }
636  fclose(f);
637 }
638 
639 void DocRecord::printMedialAxis(Octree *_octree, const std::string &fileName,
640  GFace *gf, GEdge *ge)
641 {
642  FILE *f = Fopen(fileName.c_str(), "w");
643  if(!f) {
644  Msg::Error("Could not open file '%s'", fileName.c_str());
645  return;
646  }
647  if(_adjacencies) {
648  fprintf(f, "View \"medial axis\" {\n");
649  for(PointNumero i = 0; i < numPoints; i++) {
650  std::vector<SPoint2> pts;
651  SPoint2 pc((double)points[i].where.h, (double)points[i].where.v);
652  if(!onHull(i)) {
653  GPoint p0(pc[0], pc[1], 0.0);
654  if(gf) p0 = gf->point(pc[0], pc[1]);
655  fprintf(f, "SP(%g,%g,%g){%g};\n", p0.x(), p0.y(), p0.z(), (double)i);
656  voronoiCell(i, pts);
657  for(std::size_t j = 0; j < pts.size(); j++) {
658  SVector3 pp1(pts[j].x(), pts[j].y(), 0.0);
659  SVector3 pp2(pts[(j + 1) % pts.size()].x(),
660  pts[(j + 1) % pts.size()].y(), 0.0);
661  SVector3 v1(pp1.x() - pc.x(), pp1.y() - pc.y(), 0.0);
662  SVector3 v2(pp2.x() - pc.x(), pp2.y() - pc.y(), 0.0);
663  GPoint p1(pp1.x(), pp1.y(), 0.0);
664  GPoint p2(pp2.x(), pp2.y(), 0.0);
665  if(gf) {
666  p1 = gf->point(p1.x(), p1.y());
667  p2 = gf->point(p2.x(), p2.y());
668  }
669  double P1[3] = {p1.x(), p1.y(), p1.z()};
670  double P2[3] = {p2.x(), p2.y(), p2.z()};
671  MElement *m1 = (MElement *)Octree_Search(P1, _octree);
672  MElement *m2 = (MElement *)Octree_Search(P2, _octree);
673  if(m1 && m2) {
674  MVertex *v0 = new MVertex(p1.x(), p1.y(), p1.z());
675  MVertex *v1 = new MVertex(p2.x(), p2.y(), p2.z());
676  ge->lines.push_back(new MLine(v0, v1));
677  ge->mesh_vertices.push_back(v0);
678  ge->mesh_vertices.push_back(v1);
679  fprintf(f, "SL(%g,%g,%g,%g,%g,%g){%g,%g};\n", p1.x(), p1.y(),
680  p1.z(), p2.x(), p2.y(), p2.z(), (double)i, (double)i);
681  }
682  }
683  }
684  }
685  fprintf(f, "};\n");
686  }
687  fclose(f);
688 }
689 
690 void centroidOfOrientedBox(std::vector<SPoint2> &pts, const double &angle,
691  double &xc, double &yc, double &inertia,
692  double &area)
693 {
694  const int N = pts.size();
695 
696  double sina = sin(angle);
697  double cosa = cos(angle);
698 
699  double xmin = cosa * pts[0].x() + sina * pts[0].y();
700  double xmax = cosa * pts[0].x() + sina * pts[0].y();
701  double ymin = -sina * pts[0].x() + cosa * pts[0].y();
702  double ymax = -sina * pts[0].x() + cosa * pts[0].y();
703  for(int j = 1; j < N; j++) {
704  xmin = std::min(cosa * pts[j].x() + sina * pts[j].y(), xmin);
705  ymin = std::min(-sina * pts[j].x() + cosa * pts[j].y(), ymin);
706  xmax = std::max(cosa * pts[j].x() + sina * pts[j].y(), xmax);
707  ymax = std::max(-sina * pts[j].x() + cosa * pts[j].y(), ymax);
708  }
709  double XC = 0.5 * (xmax + xmin);
710  double YC = 0.5 * (ymax + ymin);
711  xc = XC * cosa - YC * sina;
712  yc = XC * sina + YC * cosa;
713  inertia = std::max(xmax - xmin, ymax - ymin);
714  area = (xmax - xmin) * (ymax - ymin);
715 }
716 
717 void centroidOfPolygon(SPoint2 &pc, std::vector<SPoint2> &pts, double &xc,
718  double &yc, double &inertia, double &areaCell,
720 {
721  double area_tot = 0;
722  areaCell = 0.0;
723  SPoint2 center(0, 0);
724  for(std::size_t j = 0; j < pts.size(); j++) {
725  SPoint2 &pa = pts[j];
726  SPoint2 &pb = pts[(j + 1) % pts.size()];
727  const double area = triangle_area2d(pa, pb, pc);
728  const double lc = bgm ? (*bgm)((pa.x() + pb.x() + pc.x()) / 3.0,
729  (pa.y() + pb.y() + pc.y()) / 3.0, 0.0) :
730  1.0;
731  const double fact = 1. / (lc * lc * lc * lc); // rho = 1/lc^4 (emi)
732  areaCell += area;
733  area_tot += area * fact;
734  center += ((pa + pb + pc) * (area * fact / 3.0));
735  }
736  SPoint2 x = center * (1.0 / area_tot);
737  inertia = 0;
738  for(std::size_t j = 0; j < pts.size(); j++) {
739  SPoint2 &pa = pts[j];
740  SPoint2 &pb = pts[(j + 1) % pts.size()];
741  const double area = triangle_area2d(pa, pb, pc);
742 
743  const double b = sqrt((pa.x() - pb.x()) * (pa.x() - pb.x()) +
744  (pa.y() - pb.y()) * (pa.y() - pb.y()));
745  const double h = 2.0 * area / b;
746  const double a = fabs((pb.x() - pa.x()) * (pc.x() - pa.x()) *
747  (pb.y() - pa.y()) * (pc.y() - pa.y())) /
748  b;
749 
750  const double j2 =
751  (h * b * b * b + h * a * b * b + h * a * a * b + b * h * h * h) / 12.0;
752 
753  center = (pa + pb + pc) * (1.0 / 3.0);
754  const SPoint2 dx = x - center;
755  inertia += j2 + area * area * (dx.x() + dx.x() + dx.y() * dx.y());
756  }
757  xc = x.x();
758  yc = x.y();
759 }
760 
761 // Convertir les listes d'adjacence en triangles
763 {
764  // on suppose que n >= 3. points est suppose OK.
765 
766  int n = numPoints, i, j;
767  int count = 0, count2 = 0;
768  PointNumero aa, bb, cc;
769 
770  STriangle *striangle = new STriangle[n];
771 
772  count2 = CountPointsOnHull();
773 
774  // nombre de triangles que l'on doit obtenir
775  count2 = 2 * (n - 1) - count2;
776 
777  // FIXME ::: WHY 2 * ???????????????????
778  triangles = new Triangle[2 * count2];
779 
780  for(i = 0; i < n; i++) {
781  // on cree une liste de points connectes au point i (t) + nombre
782  // de points (t_length)
783  striangle[i].t =
784  ConvertDlistToArray(&points[i].adjacent, &striangle[i].t_length);
785  }
786 
787  // on balaye les noeuds de gauche a droite -> on cree les triangles
788  count = 0;
789  for(i = 0; i < n; i++) {
790  for(j = 0; j < striangle[i].t_length; j++) {
791  if((striangle[i].t[j] > i) && (striangle[i].t[j + 1] > i) &&
792  (IsRightOf(i, striangle[i].t[j], striangle[i].t[j + 1]))) {
793  aa = i;
794  bb = striangle[i].t[j];
795  cc = striangle[i].t[j + 1];
796  triangles[count].a = aa;
797  triangles[count].b = bb;
798  triangles[count].c = cc;
799  count++;
800  }
801  }
802  }
803  numTriangles = count2;
804 
805  for(int i = 0; i < n; i++) delete[] striangle[i].t;
806  delete[] striangle;
807 
808  return 1;
809 }
810 
811 // Cette routine efface toutes les listes d'adjacence de points
813 {
814  int i;
815  DListPeek p, temp;
816 
817  for(i = 0; i < numPoints; i++)
818  if(points[i].adjacent != nullptr) {
819  p = points[i].adjacent;
820  do {
821  temp = p;
822  p = Pred(p);
823  delete temp;
824  } while(p != points[i].adjacent);
825  points[i].adjacent = nullptr;
826  }
827 }
828 
830  : _hullSize(0), _hull(nullptr), _adjacencies(nullptr), numPoints(n),
831  points(nullptr), numTriangles(0), triangles(nullptr)
832 {
833  if(numPoints) points = new PointRecord[numPoints + 3000];
834 }
835 
837 {
838  if(points) delete[] points;
839  if(triangles) delete[] triangles;
840  if(_hull) delete[] _hull;
841  if(_adjacencies) {
842  for(int i = 0; i < numPoints; i++) delete[] _adjacencies[i].t;
843  delete _adjacencies;
844  }
845 }
846 
848 {
849  for(int i = 0; i < numPoints; i++) {
850  if(points[i].adjacent == nullptr) return false;
851  }
852  return true;
853 }
854 
856 {
857  if(numPoints < 3) return;
858  BuildDelaunay();
859 
861  else {
862  Msg::Error("Adjacent nullptrs found");
863  }
864 
865  RemoveAllDList();
866 }
867 
869 {
870  if(numPoints < 3) return;
871  BuildDelaunay();
873 }
874 
876 {
877  if(numPoints != p->size1())
878  throw std::runtime_error("Incompatible number of points");
879  for(int i = 0; i < p->size1(); i++) {
880  x(i) = (*p)(i, 0);
881  y(i) = (*p)(i, 1);
882  data(i) = (*p)(i, 2) < 0 ? (void *)1 : nullptr;
883  }
884 }
885 
887  int iTriangle, std::set<int> &taggedTriangles,
888  std::map<std::pair<void *, void *>, std::pair<int, int> > &edgesToTriangles)
889 {
890  if(taggedTriangles.find(iTriangle) != taggedTriangles.end()) return;
891 
892  taggedTriangles.insert(iTriangle);
893 
894  int a = triangles[iTriangle].a;
895  int b = triangles[iTriangle].b;
896  int c = triangles[iTriangle].c;
897  PointRecord *p[3] = {&points[a], &points[b], &points[c]};
898 
899  for(int j = 0; j < 3; j++) {
900  if(!lookForBoundaryEdge(p[j]->data, p[(j + 1) % 3]->data)) {
901  std::pair<void *, void *> ab =
902  std::make_pair(std::min(p[j]->data, p[(j + 1) % 3]->data),
903  std::max(p[j]->data, p[(j + 1) % 3]->data));
904  std::pair<int, int> &adj = (edgesToTriangles.find(ab))->second;
905  if(iTriangle == adj.first && adj.second != -1)
906  recur_tag_triangles(adj.second, taggedTriangles, edgesToTriangles);
907  else
908  recur_tag_triangles(adj.first, taggedTriangles, edgesToTriangles);
909  }
910  }
911 }
912 
913 std::set<int> DocRecord::tagInterior(double x, double y)
914 {
915  std::map<std::pair<void *, void *>, std::pair<int, int> > edgesToTriangles;
916  int iFirst = 1;
917  for(int i = 0; i < numTriangles; i++) {
918  int a = triangles[i].a;
919  int b = triangles[i].b;
920  int c = triangles[i].c;
921  PointRecord *p[3] = {&points[a], &points[b], &points[c]};
922 
923  // Weisstein, Eric W. "Triangle Interior." From MathWorld--A Wolfram Web
924  // Resource.
925  double x0 = points[a].where.h;
926  double y0 = points[a].where.v;
927  double x1 = points[b].where.h - points[a].where.h;
928  double y1 = points[b].where.v - points[a].where.v;
929  double x2 = points[c].where.h - points[a].where.h;
930  double y2 = points[c].where.v - points[a].where.v;
931  double k1 = ((x * y2 - x2 * y) - (x0 * y2 - x2 * y0)) / (x1 * y2 - x2 * y1);
932  double k2 =
933  -((x * y1 - x1 * y) - (x0 * y1 - x1 * y0)) / (x1 * y2 - x2 * y1);
934  if(k1 > 0.0 && k2 > 0.0 && k1 + k2 < 1.0) iFirst = i;
935  for(int j = 0; j < 3; j++) {
936  std::pair<void *, void *> ab =
937  std::make_pair(std::min(p[j]->data, p[(j + 1) % 3]->data),
938  std::max(p[j]->data, p[(j + 1) % 3]->data));
939  auto it = edgesToTriangles.find(ab);
940  if(it == edgesToTriangles.end()) {
941  edgesToTriangles[ab] = std::make_pair(i, -1);
942  }
943  else {
944  it->second.second = i;
945  }
946  }
947  }
948  std::set<int> taggedTriangles;
949  recur_tag_triangles(iFirst, taggedTriangles, edgesToTriangles);
950  return taggedTriangles;
951 }
952 
953 void DocRecord::concave(double x, double y, GFace *gf)
954 {
955  int index1;
956  int index2;
957  int index3;
958  GEdge *edge;
959  MElement *element;
960  MVertex *vertex1;
961  MVertex *vertex2;
962  std::set<int> set;
963 
964  std::vector<GEdge *> list = gf->edges();
965 
966  for(auto it1 = list.begin(); it1 != list.end(); it1++) {
967  edge = *it1;
968  for(std::size_t i = 0; i < edge->getNumMeshElements(); i++) {
969  element = edge->getMeshElement(i);
970  vertex1 = element->getVertex(0);
971  vertex2 = element->getVertex(1);
972  addBoundaryEdge(vertex1, vertex2);
973  }
974  }
975 
976  for(int i = 0; i < numPoints; i++) { points[i].vicinity.clear(); }
977 
978  try {
980  } catch(const char *err) {
981  Msg::Error("%s", err);
982  }
983 
984  set = tagInterior(x, y);
985  for(auto it2 = set.begin(); it2 != set.end(); it2++) {
986  index1 = triangles[*it2].a;
987  index2 = triangles[*it2].b;
988  index3 = triangles[*it2].c;
989  add(index1, index2);
990  add(index1, index3);
991  add(index2, index1);
992  add(index2, index3);
993  add(index3, index1);
994  add(index3, index2);
995  }
996 }
997 
998 bool DocRecord::contain(int index1, int index2, int index3)
999 {
1000  void *dataA;
1001  void *dataB;
1002  dataA = points[index2].data;
1003  dataB = points[index3].data;
1004  for(std::size_t i = 0; i < points[index1].vicinity.size() - 1; i = i + 2) {
1005  if(points[index1].vicinity[i] == dataA &&
1006  points[index1].vicinity[i + 1] == dataB) {
1007  return 1;
1008  }
1009  else if(points[index1].vicinity[i] == dataB &&
1010  points[index1].vicinity[i + 1] == dataA) {
1011  return 1;
1012  }
1013  }
1014  return 0;
1015 }
1016 
1017 void DocRecord::add(int index1, int index2)
1018 {
1019  void *data;
1020  data = points[index2].data;
1021  points[index1].vicinity.push_back(data);
1022 }
1023 
1025 {
1026  for(int i = 0; i < numPoints; i++) { points[i].flag = 0; }
1027 }
1028 
1030 {
1031  if(points[index].flag == 0) {
1032  points[index].flag = 1;
1033  return 1;
1034  }
1035  else
1036  return 0;
1037 }
1038 
1040 {
1041  int numPoints2 = 0;
1042  for(int i = 0; i < numPoints; i++) {
1043  if(points[i].flag == 0) { numPoints2++; }
1044  }
1045  PointRecord *points2 = new PointRecord[numPoints2];
1046  int index = 0;
1047  for(int i = 0; i < numPoints; i++) {
1048  if(points[i].flag == 0) {
1049  points2[index].where.h = points[i].where.h;
1050  points2[index].where.v = points[i].where.v;
1051  points2[index].data = points[i].data;
1052  points2[index].flag = points[i].flag;
1053  points2[index].identificator = points[i].identificator;
1054  index++;
1055  }
1056  }
1057  delete[] points;
1058  points = points2;
1059  numPoints = numPoints2;
1060 }
1061 
1062 void DocRecord::add_point(double x, double y, GFace *face)
1063 {
1064  PointRecord point;
1065  point.where.h = x;
1066  point.where.v = y;
1067  point.data = new MVertex(x, y, 0.0, (GEntity *)face, 2);
1068  points[numPoints] = point;
1069  numPoints = numPoints + 1;
1070 }
1071 
1073 {
1074  for(int i = 0; i < numPoints; i++) {
1075  MVertex *vertex1 = (MVertex *)points[i].data;
1076  int num = _adjacencies[i].t_length;
1077  for(int j = 0; j < num; j++) {
1078  int index = _adjacencies[i].t[j];
1079  MVertex *vertex2 = (MVertex *)points[index].data;
1080  add_edge(vertex1, vertex2);
1081  }
1082  }
1083 }
1084 
1086 
1088 {
1089  std::vector<GEdge *> const &list = gf->edges();
1090 
1091  for(auto it = list.begin(); it != list.end(); it++) {
1092  GEdge *edge = *it;
1093  for(std::size_t i = 0; i < edge->getNumMeshElements(); i++) {
1094  MElement *element = edge->getMeshElement(i);
1095  MVertex *vertex1 = element->getVertex(0);
1096  MVertex *vertex2 = element->getVertex(1);
1097  bool flag = find_edge(vertex1, vertex2);
1098  if(!flag) return false;
1099  }
1100  }
1101  return 1;
1102 }
DocRecord::FixFirst
int FixFirst(PointNumero x, PointNumero f)
Definition: DivideAndConquer.cpp:59
DocRecord::recur_tag_triangles
void recur_tag_triangles(int, std::set< int > &, std::map< std::pair< void *, void * >, std::pair< int, int > > &)
Definition: DivideAndConquer.cpp:886
DocRecord::Qtest
int Qtest(PointNumero h, PointNumero i, PointNumero j, PointNumero k)
Definition: DivideAndConquer.cpp:158
CDList
Definition: DivideAndConquer.h:40
GFace::edges
virtual std::vector< GEdge * > const & edges() const
Definition: GFace.h:121
GEdge::getMeshElement
MElement * getMeshElement(std::size_t index) const
Definition: GEdge.cpp:173
DocRecord::IsRightOf
int IsRightOf(PointNumero x, PointNumero y, PointNumero check)
Definition: DivideAndConquer.cpp:94
DocRecord::tagInterior
std::set< int > tagInterior(double, double)
Definition: DivideAndConquer.cpp:913
DocRecord::_hullSize
int _hullSize
Definition: DivideAndConquer.h:66
DocRecord::numPoints
int numPoints
Definition: DivideAndConquer.h:93
DocRecord::UpperCommonTangent
Segment UpperCommonTangent(DT vl, DT vr)
Definition: DivideAndConquer.cpp:128
GPoint::y
double y() const
Definition: GPoint.h:22
GFace.h
DocRecord::Delete
int Delete(PointNumero a, PointNumero b)
Definition: DivideAndConquer.cpp:444
circumCenterXY
void circumCenterXY(double *p1, double *p2, double *p3, double *res)
Definition: Numeric.cpp:317
GFace
Definition: GFace.h:33
DocRecord::CountPointsOnHull
int CountPointsOnHull()
Definition: DivideAndConquer.cpp:452
DocRecord::y
double & y(int i)
Definition: DivideAndConquer.h:100
DPoint
Definition: DivideAndConquer.h:22
GEdge::lines
std::vector< MLine * > lines
Definition: GEdge.h:42
angle
double angle(const SVector3 &a, const SVector3 &b)
Definition: SVector3.h:157
SPoint2
Definition: SPoint2.h:12
OS.h
centroidOfPolygon
void centroidOfPolygon(SPoint2 &pc, std::vector< SPoint2 > &pts, double &xc, double &yc, double &inertia, double &areaCell, simpleFunction< double > *bgm)
Definition: DivideAndConquer.cpp:717
DocRecord::clear_edges
void clear_edges()
Definition: DivideAndConquer.cpp:1085
robustPredicates.h
c
static double c(int i, int j, fullMatrix< double > &CA, const std::vector< SPoint3 > &P, const std::vector< SPoint3 > &Q)
Definition: discreteFrechetDistance.cpp:15
MVertex
Definition: MVertex.h:24
DocRecord::delaunay_conformity
bool delaunay_conformity(GFace *)
Definition: DivideAndConquer.cpp:1087
DocRecord::ConvertDListToVoronoiData
void ConvertDListToVoronoiData()
Definition: DivideAndConquer.cpp:515
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
DocRecord::ConvertDlistToArray
PointNumero * ConvertDlistToArray(DListPeek *dlist, int *n)
Definition: DivideAndConquer.cpp:488
DocRecord::mesh_edges
std::set< std::pair< void *, void * > > mesh_edges
Definition: DivideAndConquer.h:145
DocRecord::BuildDelaunay
int BuildDelaunay()
Definition: DivideAndConquer.cpp:310
PointRecord::where
DPoint where
Definition: DivideAndConquer.h:28
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
DocRecord::add
void add(int, int)
Definition: DivideAndConquer.cpp:1017
GFace::point
virtual GPoint point(double par1, double par2) const =0
DocRecord::First
PointNumero First(PointNumero x)
Definition: DivideAndConquer.cpp:77
SVector3
Definition: SVector3.h:16
DocRecord::AdjacentNullptrExists
bool AdjacentNullptrExists()
Definition: DivideAndConquer.cpp:847
DocRecord::lookForBoundaryEdge
bool lookForBoundaryEdge(void *p1, void *p2)
Definition: DivideAndConquer.h:124
DocRecord::ConvertDListToTriangles
int ConvertDListToTriangles()
Definition: DivideAndConquer.cpp:762
DocRecord::points
PointRecord * points
Definition: DivideAndConquer.h:95
DocRecord::DocRecord
DocRecord(int n)
Definition: DivideAndConquer.cpp:829
Segment::to
PointNumero to
Definition: DivideAndConquer.h:57
DivideAndConquer.h
GmshMessage.h
MLine.h
DocRecord::concave
void concave(double, double, GFace *)
Definition: DivideAndConquer.cpp:953
DocRecord::ConvexHull
void ConvexHull()
Definition: DivideAndConquer.cpp:471
DocRecord::initialize
void initialize()
Definition: DivideAndConquer.cpp:1024
GEntity
Definition: GEntity.h:31
DocRecord::IsLeftOf
int IsLeftOf(PointNumero x, PointNumero y, PointNumero check)
Definition: DivideAndConquer.cpp:82
GPoint
Definition: GPoint.h:13
comparePoints
static int comparePoints(const void *i, const void *j)
Definition: DivideAndConquer.cpp:296
CDList::next
DListPeek next
Definition: DivideAndConquer.h:42
DT::begin
PointNumero begin
Definition: DivideAndConquer.h:51
PointRecord::vicinity
std::vector< void * > vicinity
Definition: DivideAndConquer.h:33
PointNumero
int PointNumero
Definition: DivideAndConquer.h:20
DocRecord::_hull
PointNumero * _hull
Definition: DivideAndConquer.h:67
DocRecord::Successor
PointNumero Successor(PointNumero a, PointNumero b)
Definition: DivideAndConquer.cpp:46
Triangle::c
PointNumero c
Definition: DivideAndConquer.h:61
Fopen
FILE * Fopen(const char *f, const char *mode)
Definition: OS.cpp:273
PointRecord::adjacent
DListPeek adjacent
Definition: DivideAndConquer.h:29
MLine
Definition: MLine.h:21
GPoint::z
double z() const
Definition: GPoint.h:23
DocRecord::Insert
int Insert(PointNumero a, PointNumero b)
Definition: DivideAndConquer.cpp:406
fullMatrix< double >
DocRecord::printMedialAxis
void printMedialAxis(Octree *_octree, const std::string &, GFace *gf=nullptr, GEdge *ge=nullptr)
Definition: DivideAndConquer.cpp:639
DocRecord::DListInsert
int DListInsert(PointNumero centerPoint, PointNumero newPoint)
Definition: DivideAndConquer.cpp:319
DocRecord::numTriangles
int numTriangles
Definition: DivideAndConquer.h:96
Triangle::a
PointNumero a
Definition: DivideAndConquer.h:61
Octree_Search
void * Octree_Search(double *pt, Octree *myOctree)
Definition: Octree.cpp:81
DocRecord::setPoints
void setPoints(fullMatrix< double > *p)
Definition: DivideAndConquer.cpp:875
simpleFunction< double >
DocRecord::add_point
void add_point(double, double, GFace *)
Definition: DivideAndConquer.cpp:1062
GEntity::mesh_vertices
std::vector< MVertex * > mesh_vertices
Definition: GEntity.h:56
DocRecord::onHull
int onHull(PointNumero i)
Definition: DivideAndConquer.h:107
PointRecord
Definition: DivideAndConquer.h:27
DocRecord::~DocRecord
~DocRecord()
Definition: DivideAndConquer.cpp:836
Numeric.h
PointRecord::flag
int flag
Definition: DivideAndConquer.h:31
DocRecord::remove_point
bool remove_point(int)
Definition: DivideAndConquer.cpp:1029
triangle_area2d
double triangle_area2d(double p0[2], double p1[2], double p2[2])
Definition: Numeric.cpp:309
DPoint::h
double h
Definition: DivideAndConquer.h:24
DocRecord::contain
bool contain(int, int, int)
Definition: DivideAndConquer.cpp:998
PointRecord::identificator
int identificator
Definition: DivideAndConquer.h:32
DListRecord
struct CDList DListRecord
Definition: DivideAndConquer.h:19
DocRecord::remove_all
void remove_all()
Definition: DivideAndConquer.cpp:1039
MElement
Definition: MElement.h:30
DocRecord::RemoveAllDList
void RemoveAllDList()
Definition: DivideAndConquer.cpp:812
Triangle
Definition: DivideAndConquer.h:60
DocRecord::build_edges
void build_edges()
Definition: DivideAndConquer.cpp:1072
Segment
Definition: DivideAndConquer.h:55
SVector3::x
double x(void) const
Definition: SVector3.h:30
element
Definition: shapeFunctions.h:12
DocRecord::x
double & x(int i)
Definition: DivideAndConquer.h:99
STriangle::t
PointNumero * t
Definition: DivideAndConquer.h:46
PointRecord::data
void * data
Definition: DivideAndConquer.h:30
DocRecord::find_edge
bool find_edge(void *p1, void *p2)
Definition: DivideAndConquer.h:153
Octree
Definition: OctreeInternals.h:51
DocRecord::DListDelete
int DListDelete(DListPeek *dlist, PointNumero oldPoint)
Definition: DivideAndConquer.cpp:413
CDList::point_num
PointNumero point_num
Definition: DivideAndConquer.h:41
Segment::from
PointNumero from
Definition: DivideAndConquer.h:56
DocRecord::makePosView
void makePosView(const std::string &, GFace *gf=nullptr)
Definition: DivideAndConquer.cpp:572
DocRecord::LowerCommonTangent
Segment LowerCommonTangent(DT vl, DT vr)
Definition: DivideAndConquer.cpp:99
DocRecord::add_edge
void add_edge(void *p1, void *p2)
Definition: DivideAndConquer.h:147
DocRecord::Predecessor
PointNumero Predecessor(PointNumero a, PointNumero b)
Definition: DivideAndConquer.cpp:33
DocRecord::triangles
Triangle * triangles
Definition: DivideAndConquer.h:97
DPoint::v
double v
Definition: DivideAndConquer.h:23
SPoint2::x
double x(void) const
Definition: SPoint2.h:86
Pred
#define Pred(x)
Definition: DivideAndConquer.cpp:30
SVector3::y
double y(void) const
Definition: SVector3.h:31
fullMatrix::size1
int size1() const
Definition: fullMatrix.h:274
z
const double z
Definition: GaussQuadratureQuad.cpp:56
STriangle::t_length
int t_length
Definition: DivideAndConquer.h:47
STriangle
Definition: DivideAndConquer.h:45
picojson::copy
void copy(const std::string &s, Iter oi)
Definition: picojson.h:510
GEdge
Definition: GEdge.h:26
BackgroundMeshTools.h
GPoint.h
Succ
#define Succ(x)
Definition: DivideAndConquer.cpp:31
DT::end
PointNumero end
Definition: DivideAndConquer.h:52
GEdge::getNumMeshElements
std::size_t getNumMeshElements() const
Definition: GEdge.h:201
DocRecord::Merge
int Merge(DT vl, DT vr)
Definition: DivideAndConquer.cpp:178
DocRecord::data
void *& data(int i)
Definition: DivideAndConquer.h:101
DocRecord::_adjacencies
STriangle * _adjacencies
Definition: DivideAndConquer.h:92
Triangle::b
PointNumero b
Definition: DivideAndConquer.h:61
GPoint::x
double x() const
Definition: GPoint.h:21
DocRecord::MakeMeshWithPoints
void MakeMeshWithPoints()
Definition: DivideAndConquer.cpp:855
robustPredicates::orient2d
REAL orient2d(REAL *pa, REAL *pb, REAL *pc)
Definition: robustPredicates.cpp:1633
DocRecord::RecurTrig
DT RecurTrig(PointNumero left, PointNumero right)
Definition: DivideAndConquer.cpp:253
fullMatrix.h
DocRecord::addBoundaryEdge
void addBoundaryEdge(void *p1, void *p2)
Definition: DivideAndConquer.h:118
centroidOfOrientedBox
void centroidOfOrientedBox(std::vector< SPoint2 > &pts, const double &angle, double &xc, double &yc, double &inertia, double &area)
Definition: DivideAndConquer.cpp:690
DocRecord::Voronoi
void Voronoi()
Definition: DivideAndConquer.cpp:868
DocRecord::voronoiCell
void voronoiCell(PointNumero pt, std::vector< SPoint2 > &pts) const
Definition: DivideAndConquer.cpp:534
robustPredicates::incircle
REAL incircle(REAL *pa, REAL *pb, REAL *pc, REAL *pd)
Definition: robustPredicates.cpp:3245
DT
Definition: DivideAndConquer.h:50
SPoint2::y
double y(void) const
Definition: SPoint2.h:88