Graphviz 13.0.0~dev.20241220.2304
Loading...
Searching...
No Matches
delaunay.c
Go to the documentation of this file.
1/*************************************************************************
2 * Copyright (c) 2011 AT&T Intellectual Property
3 * All rights reserved. This program and the accompanying materials
4 * are made available under the terms of the Eclipse Public License v1.0
5 * which accompanies this distribution, and is available at
6 * https://www.eclipse.org/legal/epl-v10.html
7 *
8 * Contributors: Details at https://graphviz.org
9 *************************************************************************/
10
11#include "config.h"
12
13#include <stdbool.h>
14#include <stdio.h>
15#include <stdlib.h>
16#include <string.h>
17#include <math.h>
18#include <cgraph/cgraph.h> /* for agerr() and friends */
19#include <neatogen/delaunay.h>
20#include <util/alloc.h>
21#include <util/sort.h>
22
23#ifdef HAVE_GTS
24#include <gts.h>
25
26static int triangle_is_hole(void *triangle, void *ignored) {
27 GtsTriangle *t = triangle;
28 (void)ignored;
29
30 GtsEdge *e1, *e2, *e3;
31 GtsVertex *v1, *v2, *v3;
32
33 gts_triangle_vertices_edges(t, NULL, &v1, &v2, &v3, &e1, &e2, &e3);
34
35 if ((GTS_IS_CONSTRAINT(e1) && GTS_SEGMENT(e1)->v1 != v1) ||
36 (GTS_IS_CONSTRAINT(e2) && GTS_SEGMENT(e2)->v1 != v2) ||
37 (GTS_IS_CONSTRAINT(e3) && GTS_SEGMENT(e3)->v1 != v3))
38 return TRUE;
39
40 return FALSE;
41}
42
43static unsigned delaunay_remove_holes(GtsSurface *surface) {
44 return gts_surface_foreach_face_remove(surface, triangle_is_hole, NULL);
45}
46
47/* Derived classes for vertices and faces so we can assign integer ids
48 * to make it easier to identify them. In particular, this allows the
49 * segments and faces to refer to vertices using the order in which
50 * they were passed in.
51 */
52typedef struct {
53 GtsVertex v;
54 int idx;
55} GVertex;
56
57typedef struct {
58 GtsVertexClass parent_class;
59} GVertexClass;
60
61static GVertexClass *g_vertex_class(void)
62{
63 static GVertexClass *klass = NULL;
64
65 if (klass == NULL) {
66 GtsObjectClassInfo vertex_info = {
67 "GVertex",
68 sizeof(GVertex),
69 sizeof(GVertexClass),
70 (GtsObjectClassInitFunc) NULL,
71 (GtsObjectInitFunc) NULL,
72 (GtsArgSetFunc) NULL,
73 (GtsArgGetFunc) NULL
74 };
75 klass = gts_object_class_new(GTS_OBJECT_CLASS(gts_vertex_class()),
76 &vertex_info);
77 }
78
79 return klass;
80}
81
82typedef struct {
83 GtsFace v;
84 int idx;
85} GFace;
86
87typedef struct {
88 GtsFaceClass parent_class;
89} GFaceClass;
90
91static GFaceClass *g_face_class(void)
92{
93 static GFaceClass *klass = NULL;
94
95 if (klass == NULL) {
96 GtsObjectClassInfo face_info = {
97 "GFace",
98 sizeof(GFace),
99 sizeof(GFaceClass),
100 (GtsObjectClassInitFunc) NULL,
101 (GtsObjectInitFunc) NULL,
102 (GtsArgSetFunc) NULL,
103 (GtsArgGetFunc) NULL
104 };
105 klass = gts_object_class_new(GTS_OBJECT_CLASS(gts_face_class()),
106 &face_info);
107 }
108
109 return klass;
110}
111
112/* destroy:
113 * Destroy each edge using v, then destroy v
114 */
115static void
116destroy (GtsVertex* v)
117{
118 GSList * i;
119
120 i = v->segments;
121 while (i) {
122 GSList * next = i->next;
123 gts_object_destroy (i->data);
124 i = next;
125 }
126 g_assert (v->segments == NULL);
127 gts_object_destroy(GTS_OBJECT(v));
128}
129
130/* tri:
131 * Main entry point to using GTS for triangulation.
132 * Input is npt points with x and y coordinates stored either separately
133 * in x[] and y[] (sepArr != 0) or consecutively in x[] (sepArr == 0).
134 * Optionally, the input can include nsegs line segments, whose endpoint
135 * indices are supplied in segs[2*i] and segs[2*i+1] yielding a constrained
136 * triangulation.
137 *
138 * The return value is the corresponding gts surface, which can be queries for
139 * the triangles and line segments composing the triangulation.
140 */
141static GtsSurface*
142tri(double *x, double *y, int npt, int *segs, int nsegs, int sepArr)
143{
144 int i;
145 GtsSurface *surface;
146 GVertex **vertices = gv_calloc(npt, sizeof(GVertex *));
147 GtsEdge **edges = gv_calloc(nsegs, sizeof(GtsEdge *));
148 GSList *list = NULL;
149 GtsVertex *v1, *v2, *v3;
150 GtsTriangle *t;
151 GtsVertexClass *vcl = (GtsVertexClass *) g_vertex_class();
152 GtsEdgeClass *ecl = GTS_EDGE_CLASS (gts_constraint_class ());
153
154 if (sepArr) {
155 for (i = 0; i < npt; i++) {
156 GVertex *p = (GVertex *) gts_vertex_new(vcl, x[i], y[i], 0);
157 p->idx = i;
158 vertices[i] = p;
159 }
160 }
161 else {
162 for (i = 0; i < npt; i++) {
163 GVertex *p = (GVertex *) gts_vertex_new(vcl, x[2*i], x[2*i+1], 0);
164 p->idx = i;
165 vertices[i] = p;
166 }
167 }
168
169 /* N.B. Edges need to be created here, presumably before the
170 * the vertices are added to the face. In particular, they cannot
171 * be added created and added vi gts_delaunay_add_constraint() below.
172 */
173 for (i = 0; i < nsegs; i++) {
174 edges[i] = gts_edge_new(ecl,
175 (GtsVertex *)vertices[segs[2 * i]],
176 (GtsVertex *)vertices[segs[2 * i + 1]]);
177 }
178
179 for (i = 0; i < npt; i++)
180 list = g_slist_prepend(list, vertices[i]);
181 t = gts_triangle_enclosing(gts_triangle_class(), list, 100.);
182 g_slist_free(list);
183
184 gts_triangle_vertices(t, &v1, &v2, &v3);
185
186 surface = gts_surface_new(gts_surface_class(),
187 (GtsFaceClass *) g_face_class(),
188 gts_edge_class(),
189 gts_vertex_class());
190 gts_surface_add_face(surface, gts_face_new(gts_face_class(),
191 t->e1, t->e2, t->e3));
192
193 for (i = 0; i < npt; i++) {
194 GtsVertex *v4 = (GtsVertex *)vertices[i];
195 GtsVertex *v = gts_delaunay_add_vertex(surface, v4, NULL);
196
197 /* if v != NULL, it is a previously added pt with the same
198 * coordinates as v4, in which case we replace v4 with v
199 */
200 if (v && v4 != v) {
201 gts_vertex_replace(v4, v);
202 }
203 }
204
205 for (i = 0; i < nsegs; i++) {
206 gts_delaunay_add_constraint(surface,GTS_CONSTRAINT(edges[i]));
207 }
208
209 /* destroy enclosing triangle */
210 gts_allow_floating_vertices = TRUE;
211 gts_allow_floating_edges = TRUE;
212 destroy(v1);
213 destroy(v2);
214 destroy(v3);
215 gts_allow_floating_edges = FALSE;
216 gts_allow_floating_vertices = FALSE;
217
218 if (nsegs)
219 delaunay_remove_holes(surface);
220
221 free (edges);
222 free(vertices);
223 return surface;
224}
225
226typedef struct {
227 int n;
228 v_data *delaunay;
229} estats;
230
231static int cnt_edge(void *edge, void *stats) {
232 GtsSegment *e = edge;
233 estats *sp = stats;
234
235 sp->n++;
236 if (sp->delaunay) {
237 sp->delaunay[((GVertex*)e->v1)->idx].nedges++;
238 sp->delaunay[((GVertex*)e->v2)->idx].nedges++;
239 }
240
241 return 0;
242}
243
244static void
245edgeStats (GtsSurface* s, estats* sp)
246{
247 gts_surface_foreach_edge(s, cnt_edge, sp);
248}
249
250static int add_edge(void *edge, void *data) {
251 GtsSegment *e = edge;
252 v_data *delaunay = data;
253
254 int source = ((GVertex*)e->v1)->idx;
255 int dest = ((GVertex*)e->v2)->idx;
256
257 delaunay[source].edges[delaunay[source].nedges++] = dest;
258 delaunay[dest].edges[delaunay[dest].nedges++] = source;
259
260 return 0;
261}
262
263static v_data *delaunay_triangulation(double *x, double *y, int n) {
264 GtsSurface* s = tri(x, y, n, NULL, 0, 1);
265 int i, nedges;
266 estats stats;
267
268 if (!s) return NULL;
269
270 v_data *delaunay = gv_calloc(n, sizeof(v_data));
271
272 for (i = 0; i < n; i++) {
273 delaunay[i].ewgts = NULL;
274 delaunay[i].nedges = 1;
275 }
276
277 stats.n = 0;
278 stats.delaunay = delaunay;
279 edgeStats (s, &stats);
280 nedges = stats.n;
281 int *edges = gv_calloc(2 * nedges + n, sizeof(int));
282
283 for (i = 0; i < n; i++) {
284 delaunay[i].edges = edges;
285 edges += delaunay[i].nedges;
286 delaunay[i].edges[0] = i;
287 delaunay[i].nedges = 1;
288 }
289 gts_surface_foreach_edge(s, add_edge, delaunay);
290
291 gts_object_destroy (GTS_OBJECT (s));
292
293 return delaunay;
294}
295
296typedef struct {
297 int n;
298 int* edges;
299} estate;
300
301static int addEdge(void *edge, void *state) {
302 GtsSegment *e = edge;
303 estate *es = state;
304
305 int source = ((GVertex*)e->v1)->idx;
306 int dest = ((GVertex*)e->v2)->idx;
307
308 es->edges[2 * es->n] = source;
309 es->edges[2 * es->n + 1] = dest;
310 es->n += 1;
311
312 return 0;
313}
314
315static int vcmp(const void *x, const void *y, void *values) {
316 const int *a = x;
317 const int *b = y;
318 const double *_vals = values;
319 double va = _vals[*a];
320 double vb = _vals[*b];
321
322 if (va < vb) return -1;
323 if (va > vb) return 1;
324 return 0;
325}
326
327/* delaunay_tri:
328 * Given n points whose coordinates are in the x[] and y[]
329 * arrays, compute a Delaunay triangulation of the points.
330 * The number of edges in the triangulation is returned in pnedges.
331 * The return value itself is an array e of 2*(*pnedges) integers,
332 * with edge i having points whose indices are e[2*i] and e[2*i+1].
333 *
334 * If the points are collinear, GTS fails with 0 edges.
335 * In this case, we sort the points by x coordinates (or y coordinates
336 * if the points form a vertical line). We then return a "triangulation"
337 * consisting of the n-1 pairs of adjacent points.
338 */
339int *delaunay_tri(double *x, double *y, int n, int* pnedges)
340{
341 GtsSurface* s = tri(x, y, n, NULL, 0, 1);
342 int nedges;
343 int* edges;
344 estats stats;
345 estate state;
346
347 if (!s) return NULL;
348
349 stats.n = 0;
350 stats.delaunay = NULL;
351 edgeStats (s, &stats);
352 *pnedges = nedges = stats.n;
353
354 if (nedges) {
355 edges = gv_calloc(2 * nedges, sizeof(int));
356 state.n = 0;
357 state.edges = edges;
358 gts_surface_foreach_edge(s, addEdge, &state);
359 }
360 else {
361 int* vs = gv_calloc(n, sizeof(int));
362 int* ip;
363 int i, hd, tl;
364
365 *pnedges = nedges = n-1;
366 ip = edges = gv_calloc(2 * nedges, sizeof(int));
367
368 for (i = 0; i < n; i++)
369 vs[i] = i;
370
371 gv_sort(vs, n, sizeof(int), vcmp, x[0] == x[1] /* vertical line? */ ? y : x);
372
373 tl = vs[0];
374 for (i = 1; i < n; i++) {
375 hd = vs[i];
376 *ip++ = tl;
377 *ip++ = hd;
378 tl = hd;
379 }
380
381 free (vs);
382 }
383
384 gts_object_destroy (GTS_OBJECT (s));
385
386 return edges;
387}
388
389static int cntFace(void *face, void *data) {
390 GFace *fp = face;
391 int *ip = data;
392
393 fp->idx = *ip;
394 *ip += 1;
395
396 return 0;
397}
398
399typedef struct {
400 GtsSurface* s;
401 int* faces;
402 int* neigh;
403} fstate;
404
405typedef struct {
406 int nneigh;
407 int* neigh;
408} ninfo;
409
410static int addNeighbor(void *face, void *ni) {
411 GFace *f = face;
412 ninfo *es = ni;
413
414 es->neigh[es->nneigh] = f->idx;
415 es->nneigh++;
416
417 return 0;
418}
419
420static int addFace(void *face, void *state) {
421 GFace *f = face;
422 fstate *es = state;
423
424 int i, myid = f->idx;
425 int* ip = es->faces + 3*myid;
426 int* neigh = es->neigh + 3*myid;
427 ninfo ni;
428 GtsVertex *v1, *v2, *v3;
429
430 gts_triangle_vertices (&f->v.triangle, &v1, &v2, &v3);
431 *ip++ = ((GVertex*)v1)->idx;
432 *ip++ = ((GVertex*)v2)->idx;
433 *ip++ = ((GVertex*)v3)->idx;
434
435 ni.nneigh = 0;
436 ni.neigh = neigh;
437 gts_face_foreach_neighbor((GtsFace*)f, 0, addNeighbor, &ni);
438 for (i = ni.nneigh; i < 3; i++)
439 neigh[i] = -1;
440
441 return 0;
442}
443
444static int addTri(void *face, void *state) {
445 GFace *f = face;
446 fstate *es = state;
447
448 int myid = f->idx;
449 int* ip = es->faces + 3*myid;
450 GtsVertex *v1, *v2, *v3;
451
452 gts_triangle_vertices (&f->v.triangle, &v1, &v2, &v3);
453 *ip++ = ((GVertex*)v1)->idx;
454 *ip++ = ((GVertex*)v2)->idx;
455 *ip++ = ((GVertex*)v3)->idx;
456
457 return 0;
458}
459
460/* mkSurface:
461 * Given n points whose coordinates are in x[] and y[], and nsegs line
462 * segments whose end point indices are given in segs, return a surface
463 * corresponding the constrained Delaunay triangulation.
464 * The surface records the line segments, the triangles, and the neighboring
465 * triangles.
466 */
467surface_t*
468mkSurface (double *x, double *y, int n, int* segs, int nsegs)
469{
470 GtsSurface* s = tri(x, y, n, segs, nsegs, 1);
471 estats stats;
472 estate state;
473 fstate statf;
474 int nfaces = 0;
475
476 if (!s) return NULL;
477
478 surface_t *sf = gv_alloc(sizeof(surface_t));
479 stats.n = 0;
480 stats.delaunay = NULL;
481 edgeStats (s, &stats);
482 nsegs = stats.n;
483 segs = gv_calloc(2 * nsegs, sizeof(int));
484
485 state.n = 0;
486 state.edges = segs;
487 gts_surface_foreach_edge(s, addEdge, &state);
488
489 gts_surface_foreach_face(s, cntFace, &nfaces);
490
491 int *faces = gv_calloc(3 * nfaces, sizeof(int));
492 int *neigh = gv_calloc(3 * nfaces, sizeof(int));
493
494 statf.faces = faces;
495 statf.neigh = neigh;
496 gts_surface_foreach_face(s, addFace, &statf);
497
498 sf->nedges = nsegs;
499 sf->edges = segs;
500 sf->nfaces = nfaces;
501 sf->faces = faces;
502 sf->neigh = neigh;
503
504 gts_object_destroy (GTS_OBJECT (s));
505
506 return sf;
507}
508
509/* get_triangles:
510 * Given n points whose coordinates are stored as (x[2*i],x[2*i+1]),
511 * compute a Delaunay triangulation of the points.
512 * The number of triangles in the triangulation is returned in tris.
513 * The return value t is an array of 3*(*tris) integers,
514 * with triangle i having points whose indices are t[3*i], t[3*i+1] and t[3*i+2].
515 */
516int*
517get_triangles (double *x, int n, int* tris)
518{
519 GtsSurface* s;
520 int nfaces = 0;
521 fstate statf;
522
523 if (n <= 2) return NULL;
524
525 s = tri(x, NULL, n, NULL, 0, 0);
526 if (!s) return NULL;
527
528 gts_surface_foreach_face(s, cntFace, &nfaces);
529 statf.faces = gv_calloc(3 * nfaces, sizeof(int));
530 gts_surface_foreach_face(s, addTri, &statf);
531
532 gts_object_destroy (GTS_OBJECT (s));
533
534 *tris = nfaces;
535 return statf.faces;
536}
537
538void
540{
541 free (s->edges);
542 free (s->faces);
543 free (s->neigh);
544 free(s);
545}
546#elif defined(HAVE_TRIANGLE)
547#define TRILIBRARY
548#include <triangle.c>
549#include <assert.h>
550#include <sparse/general.h>
551
552int*
553get_triangles (double *x, int n, int* tris)
554{
555 struct triangulateio mid, vorout;
556 int i;
557
558 if (n <= 2) return NULL;
559
560 struct triangulateio in = {.numberofpoints = n};
561 in.pointlist = gv_calloc(in.numberofpoints * 2, sizeof(REAL));
562
563 for (i = 0; i < n; i++){
564 in.pointlist[i*2] = x[i*2];
565 in.pointlist[i*2 + 1] = x[i*2 + 1];
566 }
567 mid.pointlist = NULL; // Not needed if -N switch used.
568 mid.pointattributelist = NULL;
569 mid.pointmarkerlist = NULL; /* Not needed if -N or -B switch used. */
570 mid.trianglelist = NULL; // Not needed if -E switch used.
571 mid.triangleattributelist = NULL;
572 mid.neighborlist = NULL; // Needed only if -n switch used.
573 mid.segmentlist = NULL;
574 mid.segmentmarkerlist = NULL;
575 mid.edgelist = NULL; // Needed only if -e switch used.
576 mid.edgemarkerlist = NULL; // Needed if -e used and -B not used.
577 vorout.pointlist = NULL; // Needed only if -v switch used.
578 vorout.pointattributelist = NULL;
579 vorout.edgelist = NULL; // Needed only if -v switch used.
580 vorout.normlist = NULL; // Needed only if -v switch used.
581
582 /* Triangulate the points. Switches are chosen to read and write a */
583 /* PSLG (p), preserve the convex hull (c), number everything from */
584 /* zero (z), assign a regional attribute to each element (A), and */
585 /* produce an edge list (e), a Voronoi diagram (v), and a triangle */
586 /* neighbor list (n). */
587
588 triangulate("Qenv", &in, &mid, &vorout);
589 assert (mid.numberofcorners == 3);
590
591 *tris = mid.numberoftriangles;
592
593 free(in.pointlist);
594 free(in.pointattributelist);
595 free(in.pointmarkerlist);
596 free(in.regionlist);
597 free(mid.pointlist);
598 free(mid.pointattributelist);
599 free(mid.pointmarkerlist);
600 free(mid.triangleattributelist);
601 free(mid.neighborlist);
602 free(mid.segmentlist);
603 free(mid.segmentmarkerlist);
604 free(mid.edgelist);
605 free(mid.edgemarkerlist);
606 free(vorout.pointlist);
607 free(vorout.pointattributelist);
608 free(vorout.edgelist);
609 free(vorout.normlist);
610
611 return mid.trianglelist;
612}
613
614// maybe it should be replaced by RNG - relative neighborhood graph, or by GG - gabriel graph
615int*
616delaunay_tri (double *x, double *y, int n, int* nedges)
617{
618 int i;
619
620 struct triangulateio in = {
621 .numberofpoints = n,
622 .pointlist = gv_calloc(2 * n, sizeof(REAL))
623 };
624 for (i = 0; i < n; i++) {
625 in.pointlist[2 * i] = x[i];
626 in.pointlist[2 * i + 1] = y[i];
627 }
628
629 struct triangleio out = {.numberofpoints = n};
630
631 triangulate("zQNEeB", &in, &out, NULL);
632
633 *nedges = out.numberofedges;
634 free (in.pointlist);
635 free (in.pointattributelist);
636 free (in.pointmarkerlist);
637 free (in.trianglearealist);
638 free (in.triangleattributelist);
639 free (in.neighborlist);
640 free (in.segmentlist);
641 free (in.segmentmarkerlist);
642 free (in.holelist);
643 free (in.regionlist);
644 free (in.edgemarkerlist);
645 free (in.normlist);
646 free (out.pointattributelist);
647 free (out.pointmarkerlist);
648 free (out.trianglearealist);
649 free (out.triangleattributelist);
650 free (out.neighborlist);
651 free (out.segmentlist);
652 free (out.segmentmarkerlist);
653 free (out.holelist);
654 free (out.regionlist);
655 free (out.edgemarkerlist);
656 free (out.normlist);
657 return out.edgelist;
658}
659
660static v_data *delaunay_triangulation(double *x, double *y, int n) {
661 int nedges;
662 int source, dest;
663 int* edgelist = delaunay_tri (x, y, n, &nedges);
664 int i;
665
666 v_data *delaunay = gv_calloc(n, sizeof(v_data));
667 int *edges = gv_calloc(2 * nedges + n, sizeof(int));
668
669 for (i = 0; i < n; i++) {
670 delaunay[i].ewgts = NULL;
671 delaunay[i].nedges = 1;
672 }
673
674 for (i = 0; i < 2 * nedges; i++)
675 delaunay[edgelist[i]].nedges++;
676
677 for (i = 0; i < n; i++) {
678 delaunay[i].edges = edges;
679 edges += delaunay[i].nedges;
680 delaunay[i].edges[0] = i;
681 delaunay[i].nedges = 1;
682 }
683 for (i = 0; i < nedges; i++) {
684 source = edgelist[2 * i];
685 dest = edgelist[2 * i + 1];
686 delaunay[source].edges[delaunay[source].nedges++] = dest;
687 delaunay[dest].edges[delaunay[dest].nedges++] = source;
688 }
689
690 free(edgelist);
691 return delaunay;
692}
693
694surface_t*
695mkSurface (double *x, double *y, int n, int* segs, int nsegs)
696{
697 agerrorf("mkSurface not yet implemented using Triangle library\n");
698 assert (0);
699 return 0;
700}
701void
703{
704 agerrorf("freeSurface not yet implemented using Triangle library\n");
705 assert (0);
706}
707#else
708static char* err = "Graphviz built without any triangulation library\n";
709int* get_triangles (double *x, int n, int* tris)
710{
711 agerrorf("get_triangles: %s\n", err);
712 return 0;
713}
714static v_data *delaunay_triangulation(double *x, double *y, int n) {
715 agerrorf("delaunay_triangulation: %s\n", err);
716 return 0;
717}
718int *delaunay_tri(double *x, double *y, int n, int* nedges)
719{
720 agerrorf("delaunay_tri: %s\n", err);
721 return 0;
722}
723surface_t*
724mkSurface (double *x, double *y, int n, int* segs, int nsegs)
725{
726 agerrorf("mkSurface: %s\n", err);
727 return 0;
728}
729void
731{
732 agerrorf("freeSurface: %s\n", err);
733}
734#endif
735
736static void remove_edge(v_data * graph, int source, int dest)
737{
738 int i;
739 for (i = 1; i < graph[source].nedges; i++) {
740 if (graph[source].edges[i] == dest) {
741 graph[source].edges[i] = graph[source].edges[--graph[source].nedges];
742 break;
743 }
744 }
745}
746
747v_data *UG_graph(double *x, double *y, int n) {
748 v_data *delaunay;
749 int i;
750 double dist_ij, dist_ik, dist_jk, x_i, y_i, x_j, y_j;
751 int j, k, neighbor_j, neighbor_k;
752
753 if (n == 2) {
754 int *edges = gv_calloc(4, sizeof(int));
755 delaunay = gv_calloc(n, sizeof(v_data));
756 delaunay[0].ewgts = NULL;
757 delaunay[0].edges = edges;
758 delaunay[0].nedges = 2;
759 delaunay[0].edges[0] = 0;
760 delaunay[0].edges[1] = 1;
761 delaunay[1].edges = edges + 2;
762 delaunay[1].ewgts = NULL;
763 delaunay[1].nedges = 2;
764 delaunay[1].edges[0] = 1;
765 delaunay[1].edges[1] = 0;
766 return delaunay;
767 } else if (n == 1) {
768 int *edges = gv_calloc(1, sizeof(int));
769 delaunay = gv_calloc(n, sizeof(v_data));
770 delaunay[0].ewgts = NULL;
771 delaunay[0].edges = edges;
772 delaunay[0].nedges = 1;
773 delaunay[0].edges[0] = 0;
774 return delaunay;
775 }
776
777 delaunay = delaunay_triangulation(x, y, n);
778
779 // remove all edges v-u if there is w, neighbor of u or v, that is closer to both u and v than dist(u,v)
780 for (i = 0; i < n; i++) {
781 x_i = x[i];
782 y_i = y[i];
783 for (j = 1; j < delaunay[i].nedges;) {
784 neighbor_j = delaunay[i].edges[j];
785 x_j = x[neighbor_j];
786 y_j = y[neighbor_j];
787 dist_ij = (x_j - x_i) * (x_j - x_i) + (y_j - y_i) * (y_j - y_i);
788 // now look at i'th neighbors to see whether there is a node in the "forbidden region"
789 // we will also go through neighbor_j's neighbors when we traverse the edge from its other side
790 bool removed = false;
791 for (k = 1; k < delaunay[i].nedges && !removed; k++) {
792 neighbor_k = delaunay[i].edges[k];
793 dist_ik = (x[neighbor_k] - x_i) * (x[neighbor_k] - x_i) +
794 (y[neighbor_k] - y_i) * (y[neighbor_k] - y_i);
795 if (dist_ik < dist_ij) {
796 dist_jk = (x[neighbor_k] - x_j) * (x[neighbor_k] - x_j) +
797 (y[neighbor_k] - y_j) * (y[neighbor_k] - y_j);
798 if (dist_jk < dist_ij) {
799 // remove the edge beteween i and neighbor j
800 delaunay[i].edges[j] = delaunay[i].edges[--delaunay[i].nedges];
801 remove_edge(delaunay, neighbor_j, i);
802 removed = true;
803 }
804 }
805 }
806 if (!removed) {
807 j++;
808 }
809 }
810 }
811 return delaunay;
812}
813
815{
816 if (graph != NULL) {
817 free(graph[0].edges);
818 free(graph[0].ewgts);
819 free(graph);
820 }
821}
822
824{
825 if (graph != NULL) {
826 free(graph[0].edges);
827 free(graph[0].ewgts);
828#ifdef DIGCOLA
829 free(graph[0].edists);
830#endif
831 free(graph);
832 }
833}
834
static void out(agerrlevel_t level, const char *fmt, va_list args)
Report messages using a user-supplied or default write function.
Definition agerror.c:84
Memory allocation wrappers that exit on failure.
static void * gv_calloc(size_t nmemb, size_t size)
Definition alloc.h:26
static void * gv_alloc(size_t size)
Definition alloc.h:47
abstract graph C library, Cgraph API
int * get_triangles(double *x, int n, int *tris)
Definition delaunay.c:709
void freeGraph(v_data *graph)
Definition delaunay.c:814
v_data * UG_graph(double *x, double *y, int n)
Definition delaunay.c:747
void freeGraphData(vtx_data *graph)
Definition delaunay.c:823
static char * err
Definition delaunay.c:708
static v_data * delaunay_triangulation(double *x, double *y, int n)
Definition delaunay.c:714
surface_t * mkSurface(double *x, double *y, int n, int *segs, int nsegs)
Definition delaunay.c:724
static void remove_edge(v_data *graph, int source, int dest)
Definition delaunay.c:736
int * delaunay_tri(double *x, double *y, int n, int *nedges)
Definition delaunay.c:718
void freeSurface(surface_t *s)
Definition delaunay.c:730
void add_edge(edgelist *list, Agedge_t *e)
Definition edgelist.c:57
#define REAL
Definition gmlparse.c:374
void free(void *)
edge
Definition gmlparse.y:240
node NULL
Definition grammar.y:163
void agerrorf(const char *fmt,...)
Definition agerror.c:165
Agraph_t * graph(char *name)
Definition gv.cpp:30
static void addEdge(edge_t *de, edge_t *e)
Definition layout.c:349
struct _tri tri
static tri * addTri(int i, int j, tri *oldp)
static int nedges
total no. of edges used in routing
Definition routespl.c:31
static int triangulate(pointnlink_t **, size_t)
Definition shortest.c:318
static triangles_t tris
Definition shortest.c:55
qsort with carried along context
static void gv_sort(void *base, size_t nmemb, size_t size, int(*compar)(const void *, const void *, void *), void *arg)
qsort with an extra state parameter, ala qsort_r
Definition sort.h:35
Definition legal.c:50
Definition cdt.h:100
int nedges
Definition delaunay.h:16
int * edges
Definition delaunay.h:17
int * faces
Definition delaunay.h:19
int nfaces
Definition delaunay.h:18
int * neigh
Definition delaunay.h:20
int nedges
Definition sparsegraph.h:22
int * edges
Definition sparsegraph.h:23
float * ewgts
Definition sparsegraph.h:24
Definition grammar.c:93