284#define FILENAMESIZE 2048
289#define INPUTLINESIZE 1024
295#define TRIPERBLOCK 4092
296#define SUBSEGPERBLOCK 508
297#define VERTEXPERBLOCK 4092
298#define VIRUSPERBLOCK 1020
300#define BADSUBSEGPERBLOCK 252
302#define BADTRIPERBLOCK 4092
304#define FLIPSTACKERPERBLOCK 252
306#define SPLAYNODEPERBLOCK 508
313#define SEGMENTVERTEX 1
315#define DEADVERTEX -32768
316#define UNDEADVERTEX -32767
336#define SAMPLEFACTOR 11
346#define PI 3.141592653589793238462643383279502884197169399375105820974944592308
350#define SQUAREROOTTWO 1.4142135623730950488016887242096980785696718753769480732
354#define ONETHIRD 0.333333333333333333333333333333333333333333333333333333333333
366#if defined(_MSC_EXTENSIONS)
367 #define DELTA_EPOCH_IN_MICROSECS 11644473600000000Ui64
369 #define DELTA_EPOCH_IN_MICROSECS 11644473600000000ULL
380int gettimeofday(
struct timeval *tv,
struct timezone *tz)
386 unsigned __int64 tmpres = 0;
387 static int tzflag = 0;
391 GetSystemTimeAsFileTime(&ft);
396 tmpres |= ft.dwHighDateTime;
398 tmpres |= ft.dwLowDateTime;
405 tmpres -= DELTA_EPOCH_IN_MICROSECS;
409 tv->tv_sec = (long)(tmpres / 1000000UL);
410 tv->tv_usec = (long)(tmpres % 1000000UL);
422 tz->tz_minuteswest = _timezone / 60;
423 tz->tz_dsttime = _daylight;
435#include <fpu_control.h>
452enum locateresult {INTRIANGLE, ONEDGE, ONVERTEX, OUTSIDE};
460enum insertvertexresult {SUCCESSFULVERTEX, ENCROACHINGVERTEX, VIOLATINGVERTEX,
468enum finddirectionresult {WITHIN, LEFTCOLLINEAR, RIGHTCOLLINEAR};
587typedef REAL **triangle;
604typedef REAL **subseg;
628 vertex subsegorg, subsegdest;
638 vertex triangorg, triangdest, triangapex;
708 VOID **firstblock, **nowblock;
717 long items, maxitems;
718 int unallocateditems;
728REAL ccwerrboundA, ccwerrboundB, ccwerrboundC;
729REAL iccerrboundA, iccerrboundB, iccerrboundC;
730REAL o3derrboundA, o3derrboundB, o3derrboundC;
734unsigned long randomseed;
760 int nextnonemptyq[4096];
769 REAL xmin, xmax, ymin, ymax;
794 long counterclockcount;
797 long circumcentercount;
802 vertex infvertex1, infvertex2, infvertex3;
807 triangle *dummytribase;
814 subseg *dummysubbase;
819 struct otri recenttri;
862 int poly, refine, quality, vararea, fixedarea, usertest;
863 int regionattrib, convex, weighted, jettison;
865 int edgesout, voronoi, neighbors, geomview;
866 int nobound, nopolywritten, nonodewritten, noelewritten, noiterationnum;
867 int noholes, noexact, conformdel;
868 int incremental, sweepline, dwyer;
876 REAL minangle, goodangle, offconstant;
882 char innodefilename[FILENAMESIZE];
883 char inelefilename[FILENAMESIZE];
884 char inpolyfilename[FILENAMESIZE];
885 char areafilename[FILENAMESIZE];
886 char outnodefilename[FILENAMESIZE];
887 char outelefilename[FILENAMESIZE];
888 char outpolyfilename[FILENAMESIZE];
889 char edgefilename[FILENAMESIZE];
890 char vnodefilename[FILENAMESIZE];
891 char vedgefilename[FILENAMESIZE];
892 char neighborfilename[FILENAMESIZE];
893 char offfilename[FILENAMESIZE];
1011int plus1mod3[3] = {1, 2, 0};
1012int minus1mod3[3] = {2, 0, 1};
1021#define decode(ptr, otri) \
1022 (otri).orient = (int) ((unsigned long) (ptr) & (unsigned long) 3l); \
1023 (otri).tri = (triangle *) \
1024 ((unsigned long) (ptr) ^ (unsigned long) (otri).orient)
1030#define encode(otri) \
1031 (triangle) ((unsigned long) (otri).tri | (unsigned long) (otri).orient)
1041#define sym(otri1, otri2) \
1042 ptr = (otri1).tri[(otri1).orient]; \
1045#define symself(otri) \
1046 ptr = (otri).tri[(otri).orient]; \
1051#define lnext(otri1, otri2) \
1052 (otri2).tri = (otri1).tri; \
1053 (otri2).orient = plus1mod3[(otri1).orient]
1055#define lnextself(otri) \
1056 (otri).orient = plus1mod3[(otri).orient]
1060#define lprev(otri1, otri2) \
1061 (otri2).tri = (otri1).tri; \
1062 (otri2).orient = minus1mod3[(otri1).orient]
1064#define lprevself(otri) \
1065 (otri).orient = minus1mod3[(otri).orient]
1071#define onext(otri1, otri2) \
1072 lprev(otri1, otri2); \
1075#define onextself(otri) \
1083#define oprev(otri1, otri2) \
1084 sym(otri1, otri2); \
1087#define oprevself(otri) \
1095#define dnext(otri1, otri2) \
1096 sym(otri1, otri2); \
1099#define dnextself(otri) \
1107#define dprev(otri1, otri2) \
1108 lnext(otri1, otri2); \
1111#define dprevself(otri) \
1119#define rnext(otri1, otri2) \
1120 sym(otri1, otri2); \
1124#define rnextself(otri) \
1133#define rprev(otri1, otri2) \
1134 sym(otri1, otri2); \
1138#define rprevself(otri) \
1146#define org(otri, vertexptr) \
1147 vertexptr = (vertex) (otri).tri[plus1mod3[(otri).orient] + 3]
1149#define dest(otri, vertexptr) \
1150 vertexptr = (vertex) (otri).tri[minus1mod3[(otri).orient] + 3]
1152#define apex(otri, vertexptr) \
1153 vertexptr = (vertex) (otri).tri[(otri).orient + 3]
1155#define setorg(otri, vertexptr) \
1156 (otri).tri[plus1mod3[(otri).orient] + 3] = (triangle) vertexptr
1158#define setdest(otri, vertexptr) \
1159 (otri).tri[minus1mod3[(otri).orient] + 3] = (triangle) vertexptr
1161#define setapex(otri, vertexptr) \
1162 (otri).tri[(otri).orient + 3] = (triangle) vertexptr
1166#define bond(otri1, otri2) \
1167 (otri1).tri[(otri1).orient] = encode(otri2); \
1168 (otri2).tri[(otri2).orient] = encode(otri1)
1175#define dissolve(otri) \
1176 (otri).tri[(otri).orient] = (triangle) m->dummytri
1180#define otricopy(otri1, otri2) \
1181 (otri2).tri = (otri1).tri; \
1182 (otri2).orient = (otri1).orient
1186#define otriequal(otri1, otri2) \
1187 (((otri1).tri == (otri2).tri) && \
1188 ((otri1).orient == (otri2).orient))
1193#define infect(otri) \
1194 (otri).tri[6] = (triangle) \
1195 ((unsigned long) (otri).tri[6] | (unsigned long) 2l)
1197#define uninfect(otri) \
1198 (otri).tri[6] = (triangle) \
1199 ((unsigned long) (otri).tri[6] & ~ (unsigned long) 2l)
1203#define infected(otri) \
1204 (((unsigned long) (otri).tri[6] & (unsigned long) 2l) != 0l)
1208#define elemattribute(otri, attnum) \
1209 ((REAL *) (otri).tri)[m->elemattribindex + (attnum)]
1211#define setelemattribute(otri, attnum, value) \
1212 ((REAL *) (otri).tri)[m->elemattribindex + (attnum)] = value
1216#define areabound(otri) ((REAL *) (otri).tri)[m->areaboundindex]
1218#define setareabound(otri, value) \
1219 ((REAL *) (otri).tri)[m->areaboundindex] = value
1226#define deadtri(tria) ((tria)[1] == (triangle) NULL)
1228#define killtri(tria) \
1229 (tria)[1] = (triangle) NULL; \
1230 (tria)[3] = (triangle) NULL
1241#define sdecode(sptr, osub) \
1242 (osub).ssorient = (int) ((unsigned long) (sptr) & (unsigned long) 1l); \
1243 (osub).ss = (subseg *) \
1244 ((unsigned long) (sptr) & ~ (unsigned long) 3l)
1250#define sencode(osub) \
1251 (subseg) ((unsigned long) (osub).ss | (unsigned long) (osub).ssorient)
1255#define ssym(osub1, osub2) \
1256 (osub2).ss = (osub1).ss; \
1257 (osub2).ssorient = 1 - (osub1).ssorient
1259#define ssymself(osub) \
1260 (osub).ssorient = 1 - (osub).ssorient
1265#define spivot(osub1, osub2) \
1266 sptr = (osub1).ss[(osub1).ssorient]; \
1267 sdecode(sptr, osub2)
1269#define spivotself(osub) \
1270 sptr = (osub).ss[(osub).ssorient]; \
1276#define snext(osub1, osub2) \
1277 sptr = (osub1).ss[1 - (osub1).ssorient]; \
1278 sdecode(sptr, osub2)
1280#define snextself(osub) \
1281 sptr = (osub).ss[1 - (osub).ssorient]; \
1287#define sorg(osub, vertexptr) \
1288 vertexptr = (vertex) (osub).ss[2 + (osub).ssorient]
1290#define sdest(osub, vertexptr) \
1291 vertexptr = (vertex) (osub).ss[3 - (osub).ssorient]
1293#define setsorg(osub, vertexptr) \
1294 (osub).ss[2 + (osub).ssorient] = (subseg) vertexptr
1296#define setsdest(osub, vertexptr) \
1297 (osub).ss[3 - (osub).ssorient] = (subseg) vertexptr
1299#define segorg(osub, vertexptr) \
1300 vertexptr = (vertex) (osub).ss[4 + (osub).ssorient]
1302#define segdest(osub, vertexptr) \
1303 vertexptr = (vertex) (osub).ss[5 - (osub).ssorient]
1305#define setsegorg(osub, vertexptr) \
1306 (osub).ss[4 + (osub).ssorient] = (subseg) vertexptr
1308#define setsegdest(osub, vertexptr) \
1309 (osub).ss[5 - (osub).ssorient] = (subseg) vertexptr
1315#define mark(osub) (* (int *) ((osub).ss + 8))
1317#define setmark(osub, value) \
1318 * (int *) ((osub).ss + 8) = value
1322#define sbond(osub1, osub2) \
1323 (osub1).ss[(osub1).ssorient] = sencode(osub2); \
1324 (osub2).ss[(osub2).ssorient] = sencode(osub1)
1329#define sdissolve(osub) \
1330 (osub).ss[(osub).ssorient] = (subseg) m->dummysub
1334#define subsegcopy(osub1, osub2) \
1335 (osub2).ss = (osub1).ss; \
1336 (osub2).ssorient = (osub1).ssorient
1340#define subsegequal(osub1, osub2) \
1341 (((osub1).ss == (osub2).ss) && \
1342 ((osub1).ssorient == (osub2).ssorient))
1349#define deadsubseg(sub) ((sub)[1] == (subseg) NULL)
1351#define killsubseg(sub) \
1352 (sub)[1] = (subseg) NULL; \
1353 (sub)[2] = (subseg) NULL
1361#define tspivot(otri, osub) \
1362 sptr = (subseg) (otri).tri[6 + (otri).orient]; \
1368#define stpivot(osub, otri) \
1369 ptr = (triangle) (osub).ss[6 + (osub).ssorient]; \
1374#define tsbond(otri, osub) \
1375 (otri).tri[6 + (otri).orient] = (triangle) sencode(osub); \
1376 (osub).ss[6 + (osub).ssorient] = (subseg) encode(otri)
1380#define tsdissolve(otri) \
1381 (otri).tri[6 + (otri).orient] = (triangle) m->dummysub
1385#define stdissolve(osub) \
1386 (osub).ss[6 + (osub).ssorient] = (subseg) m->dummytri
1392#define vertexmark(vx) ((int *) (vx))[m->vertexmarkindex]
1394#define setvertexmark(vx, value) \
1395 ((int *) (vx))[m->vertexmarkindex] = value
1397#define vertextype(vx) ((int *) (vx))[m->vertexmarkindex + 1]
1399#define setvertextype(vx, value) \
1400 ((int *) (vx))[m->vertexmarkindex + 1] = value
1402#define vertex2tri(vx) ((triangle *) (vx))[m->vertex2triindex]
1404#define setvertex2tri(vx, value) \
1405 ((triangle *) (vx))[m->vertex2triindex] = value
1444#ifdef ANSI_DECLARATORS
1445int triunsuitable(vertex triorg, vertex tridest, vertex triapex, REAL area)
1447int triunsuitable(triorg, tridest, triapex, area)
1455 REAL dxoa, dxda, dxod;
1456 REAL dyoa, dyda, dyod;
1457 REAL oalen, dalen, odlen;
1460 dxoa = triorg[0] - triapex[0];
1461 dyoa = triorg[1] - triapex[1];
1462 dxda = tridest[0] - triapex[0];
1463 dyda = tridest[1] - triapex[1];
1464 dxod = triorg[0] - tridest[0];
1465 dyod = triorg[1] - tridest[1];
1467 oalen = dxoa * dxoa + dyoa * dyoa;
1468 dalen = dxda * dxda + dyda * dyda;
1469 odlen = dxod * dxod + dyod * dyod;
1471 maxlen = (dalen > oalen) ? dalen : oalen;
1472 maxlen = (odlen > maxlen) ? odlen : maxlen;
1474 if (maxlen > 0.05 * (triorg[0] * triorg[0] + triorg[1] * triorg[1]) + 0.02) {
1491#ifdef ANSI_DECLARATORS
1492void triexit(
int status)
1502#ifdef ANSI_DECLARATORS
1503VOID *trimalloc(
int size)
1505VOID *trimalloc(size)
1512 memptr = (VOID *) malloc((
unsigned int) size);
1513 if (memptr == (VOID *) NULL) {
1514 printf(
"Error: Out of memory.\n");
1520#ifdef ANSI_DECLARATORS
1521void trifree(VOID *memptr)
1551 printf(
"triangle [-pAcjevngBPNEIOXzo_lQVh] input_file\n");
1553 printf(
"triangle [-pAcjevngBPNEIOXzo_iFlCQVh] input_file\n");
1557 printf(
"triangle [-prq__a__uAcDjevngBPNEIOXzo_YS__lQVh] input_file\n");
1559 printf(
"triangle [-prq__a__uAcDjevngBPNEIOXzo_YS__iFlsCQVh] input_file\n");
1563 printf(
" -p Triangulates a Planar Straight Line Graph (.poly file).\n");
1565 printf(
" -r Refines a previously generated mesh.\n");
1567 " -q Quality mesh generation. A minimum angle may be specified.\n");
1568 printf(
" -a Applies a maximum triangle area constraint.\n");
1569 printf(
" -u Applies a user-defined triangle constraint.\n");
1572 " -A Applies attributes to identify triangles in certain regions.\n");
1573 printf(
" -c Encloses the convex hull with segments.\n");
1575 printf(
" -D Conforming Delaunay: all triangles are truly Delaunay.\n");
1581 printf(
" -j Jettison unused vertices from output .node file.\n");
1582 printf(
" -e Generates an edge list.\n");
1583 printf(
" -v Generates a Voronoi diagram.\n");
1584 printf(
" -n Generates a list of triangle neighbors.\n");
1585 printf(
" -g Generates an .off file for Geomview.\n");
1586 printf(
" -B Suppresses output of boundary information.\n");
1587 printf(
" -P Suppresses output of .poly file.\n");
1588 printf(
" -N Suppresses output of .node file.\n");
1589 printf(
" -E Suppresses output of .ele file.\n");
1590 printf(
" -I Suppresses mesh iteration numbers.\n");
1591 printf(
" -O Ignores holes in .poly file.\n");
1592 printf(
" -X Suppresses use of exact arithmetic.\n");
1593 printf(
" -z Numbers all items starting from zero (rather than one).\n");
1594 printf(
" -o2 Generates second-order subparametric elements.\n");
1596 printf(
" -Y Suppresses boundary segment splitting.\n");
1597 printf(
" -S Specifies maximum number of added Steiner points.\n");
1600 printf(
" -i Uses incremental method, rather than divide-and-conquer.\n");
1601 printf(
" -F Uses Fortune's sweepline algorithm, rather than d-and-c.\n");
1603 printf(
" -l Uses vertical cuts only, rather than alternating cuts.\n");
1607 " -s Force segments into mesh by splitting (instead of using CDT).\n");
1609 printf(
" -C Check consistency of final mesh.\n");
1611 printf(
" -Q Quiet: No terminal output except errors.\n");
1612 printf(
" -V Verbose: Detailed information on what I'm doing.\n");
1613 printf(
" -h Help: Detailed instructions for Triangle.\n");
1629 printf(
"Triangle\n");
1631"A Two-Dimensional Quality Mesh Generator and Delaunay Triangulator.\n");
1632 printf(
"Version 1.6\n\n");
1634"Copyright 1993, 1995, 1997, 1998, 2002, 2005 Jonathan Richard Shewchuk\n");
1635 printf(
"2360 Woolsey #H / Berkeley, California 94705-1927\n");
1636 printf(
"Bugs/comments to jrs@cs.berkeley.edu\n");
1638"Created as part of the Quake project (tools for earthquake simulation).\n");
1640"Supported in part by NSF Grant CMS-9318163 and an NSERC 1967 Scholarship.\n");
1641 printf(
"There is no warranty whatsoever. Use at your own risk.\n");
1643 printf(
"This executable is compiled for single precision arithmetic.\n\n\n");
1645 printf(
"This executable is compiled for double precision arithmetic.\n\n\n");
1648"Triangle generates exact Delaunay triangulations, constrained Delaunay\n");
1650"triangulations, conforming Delaunay triangulations, Voronoi diagrams, and\n");
1652"high-quality triangular meshes. The latter can be generated with no small\n"
1655"or large angles, and are thus suitable for finite element analysis. If no\n"
1658"command line switch is specified, your .node input file is read, and the\n");
1660"Delaunay triangulation is returned in .node and .ele output files. The\n");
1661 printf(
"command syntax is:\n\n");
1662 printf(
"triangle [-prq__a__uAcDjevngBPNEIOXzo_YS__iFlsCQVh] input_file\n\n");
1664"Underscores indicate that numbers may optionally follow certain switches.\n");
1666"Do not leave any space between a switch and its numeric parameter.\n");
1668"input_file must be a file with extension .node, or extension .poly if the\n");
1670"-p switch is used. If -r is used, you must supply .node and .ele files,\n");
1672"and possibly a .poly file and an .area file as well. The formats of these\n"
1674 printf(
"files are described below.\n\n");
1675 printf(
"Command Line Switches:\n\n");
1677" -p Reads a Planar Straight Line Graph (.poly file), which can specify\n"
1680" vertices, segments, holes, regional attributes, and regional area\n");
1682" constraints. Generates a constrained Delaunay triangulation (CDT)\n"
1685" fitting the input; or, if -s, -q, -a, or -u is used, a conforming\n");
1687" constrained Delaunay triangulation (CCDT). If you want a truly\n");
1689" Delaunay (not just constrained Delaunay) triangulation, use -D as\n");
1691" well. When -p is not used, Triangle reads a .node file by default.\n"
1694" -r Refines a previously generated mesh. The mesh is read from a .node\n"
1697" file and an .ele file. If -p is also used, a .poly file is read\n");
1699" and used to constrain segments in the mesh. If -a is also used\n");
1701" (with no number following), an .area file is read and used to\n");
1703" impose area constraints on the mesh. Further details on refinement\n"
1705 printf(
" appear below.\n");
1707" -q Quality mesh generation by Delaunay refinement (a hybrid of Paul\n");
1709" Chew's and Jim Ruppert's algorithms). Adds vertices to the mesh to\n"
1712" ensure that all angles are between 20 and 140 degrees. An\n");
1714" alternative bound on the minimum angle, replacing 20 degrees, may\n");
1716" be specified after the `q'. The specified angle may include a\n");
1718" decimal point, but not exponential notation. Note that a bound of\n"
1721" theta degrees on the smallest angle also implies a bound of\n");
1723" (180 - 2 theta) on the largest angle. If the minimum angle is 28.6\n"
1726" degrees or smaller, Triangle is mathematically guaranteed to\n");
1728" terminate (assuming infinite precision arithmetic--Triangle may\n");
1730" fail to terminate if you run out of precision). In practice,\n");
1732" Triangle often succeeds for minimum angles up to 34 degrees. For\n");
1734" some meshes, however, you might need to reduce the minimum angle to\n"
1737" avoid problems associated with insufficient floating-point\n");
1738 printf(
" precision.\n");
1740" -a Imposes a maximum triangle area. If a number follows the `a', no\n");
1742" triangle is generated whose area is larger than that number. If no\n"
1745" number is specified, an .area file (if -r is used) or .poly file\n");
1747" (if -r is not used) specifies a set of maximum area constraints.\n");
1749" An .area file contains a separate area constraint for each\n");
1751" triangle, and is useful for refining a finite element mesh based on\n"
1754" a posteriori error estimates. A .poly file can optionally contain\n"
1757" an area constraint for each segment-bounded region, thereby\n");
1759" controlling triangle densities in a first triangulation of a PSLG.\n"
1762" You can impose both a fixed area constraint and a varying area\n");
1764" constraint by invoking the -a switch twice, once with and once\n");
1766" without a number following. Each area specified may include a\n");
1767 printf(
" decimal point.\n");
1769" -u Imposes a user-defined constraint on triangle size. There are two\n"
1772" ways to use this feature. One is to edit the triunsuitable()\n");
1774" procedure in triangle.c to encode any constraint you like, then\n");
1776" recompile Triangle. The other is to compile triangle.c with the\n");
1778" EXTERNAL_TEST symbol set (compiler switch -DEXTERNAL_TEST), then\n");
1780" link Triangle with a separate object file that implements\n");
1782" triunsuitable(). In either case, the -u switch causes the user-\n");
1783 printf(
" defined test to be applied to every triangle.\n");
1785" -A Assigns an additional floating-point attribute to each triangle\n");
1787" that identifies what segment-bounded region each triangle belongs\n");
1789" to. Attributes are assigned to regions by the .poly file. If a\n");
1791" region is not explicitly marked by the .poly file, triangles in\n");
1793" that region are assigned an attribute of zero. The -A switch has\n");
1795" an effect only when the -p switch is used and the -r switch is not.\n"
1798" -c Creates segments on the convex hull of the triangulation. If you\n");
1800" are triangulating a vertex set, this switch causes a .poly file to\n"
1803" be written, containing all edges of the convex hull. If you are\n");
1805" triangulating a PSLG, this switch specifies that the whole convex\n");
1807" hull of the PSLG should be triangulated, regardless of what\n");
1809" segments the PSLG has. If you do not use this switch when\n");
1811" triangulating a PSLG, Triangle assumes that you have identified the\n"
1814" region to be triangulated by surrounding it with segments of the\n");
1816" input PSLG. Beware: if you are not careful, this switch can cause\n"
1819" the introduction of an extremely thin angle between a PSLG segment\n"
1822" and a convex hull segment, which can cause overrefinement (and\n");
1824" possibly failure if Triangle runs out of precision). If you are\n");
1826" refining a mesh, the -c switch works differently: it causes a\n");
1828" .poly file to be written containing the boundary edges of the mesh\n"
1830 printf(
" (useful if no .poly file was read).\n");
1832" -D Conforming Delaunay triangulation: use this switch if you want to\n"
1835" ensure that all the triangles in the mesh are Delaunay, and not\n");
1837" merely constrained Delaunay; or if you want to ensure that all the\n"
1840" Voronoi vertices lie within the triangulation. (Some finite volume\n"
1843" methods have this requirement.) This switch invokes Ruppert's\n");
1845" original algorithm, which splits every subsegment whose diametral\n");
1847" circle is encroached. It usually increases the number of vertices\n"
1849 printf(
" and triangles.\n");
1851" -j Jettisons vertices that are not part of the final triangulation\n");
1853" from the output .node file. By default, Triangle copies all\n");
1855" vertices in the input .node file to the output .node file, in the\n");
1857" same order, so their indices do not change. The -j switch prevents\n"
1860" duplicated input vertices, or vertices `eaten' by holes, from\n");
1862" appearing in the output .node file. Thus, if two input vertices\n");
1864" have exactly the same coordinates, only the first appears in the\n");
1866" output. If any vertices are jettisoned, the vertex numbering in\n");
1868" the output .node file differs from that of the input .node file.\n");
1870" -e Outputs (to an .edge file) a list of edges of the triangulation.\n");
1872" -v Outputs the Voronoi diagram associated with the triangulation.\n");
1874" Does not attempt to detect degeneracies, so some Voronoi vertices\n");
1876" may be duplicated. See the discussion of Voronoi diagrams below.\n");
1878" -n Outputs (to a .neigh file) a list of triangles neighboring each\n");
1879 printf(
" triangle.\n");
1881" -g Outputs the mesh to an Object File Format (.off) file, suitable for\n"
1883 printf(
" viewing with the Geometry Center's Geomview package.\n");
1885" -B No boundary markers in the output .node, .poly, and .edge output\n");
1887" files. See the detailed discussion of boundary markers below.\n");
1889" -P No output .poly file. Saves disk space, but you lose the ability\n");
1891" to maintain constraining segments on later refinements of the mesh.\n"
1893 printf(
" -N No output .node file.\n");
1894 printf(
" -E No output .ele file.\n");
1896" -I No iteration numbers. Suppresses the output of .node and .poly\n");
1898" files, so your input files won't be overwritten. (If your input is\n"
1901" a .poly file only, a .node file is written.) Cannot be used with\n");
1903" the -r switch, because that would overwrite your input .ele file.\n");
1905" Shouldn't be used with the -q, -a, -u, or -s switch if you are\n");
1907" using a .node file for input, because no .node file is written, so\n"
1909 printf(
" there is no record of any added Steiner points.\n");
1910 printf(
" -O No holes. Ignores the holes in the .poly file.\n");
1912" -X No exact arithmetic. Normally, Triangle uses exact floating-point\n"
1915" arithmetic for certain tests if it thinks the inexact tests are not\n"
1918" accurate enough. Exact arithmetic ensures the robustness of the\n");
1920" triangulation algorithms, despite floating-point roundoff error.\n");
1922" Disabling exact arithmetic with the -X switch causes a small\n");
1924" improvement in speed and creates the possibility that Triangle will\n"
1926 printf(
" fail to produce a valid mesh. Not recommended.\n");
1928" -z Numbers all items starting from zero (rather than one). Note that\n"
1931" this switch is normally overridden by the value used to number the\n"
1934" first vertex of the input .node or .poly file. However, this\n");
1936" switch is useful when calling Triangle from another program.\n");
1938" -o2 Generates second-order subparametric elements with six nodes each.\n"
1941" -Y No new vertices on the boundary. This switch is useful when the\n");
1943" mesh boundary must be preserved so that it conforms to some\n");
1945" adjacent mesh. Be forewarned that you will probably sacrifice much\n"
1948" of the quality of the mesh; Triangle will try, but the resulting\n");
1950" mesh may contain poorly shaped triangles. Works well if all the\n");
1952" boundary vertices are closely spaced. Specify this switch twice\n");
1954" (`-YY') to prevent all segment splitting, including internal\n");
1955 printf(
" boundaries.\n");
1957" -S Specifies the maximum number of Steiner points (vertices that are\n");
1959" not in the input, but are added to meet the constraints on minimum\n"
1962" angle and maximum area). The default is to allow an unlimited\n");
1964" number. If you specify this switch with no number after it,\n");
1966" the limit is set to zero. Triangle always adds vertices at segment\n"
1969" intersections, even if it needs to use more vertices than the limit\n"
1972" you set. When Triangle inserts segments by splitting (-s), it\n");
1974" always adds enough vertices to ensure that all the segments of the\n"
1976 printf(
" PLSG are recovered, ignoring the limit if necessary.\n");
1978" -i Uses an incremental rather than a divide-and-conquer algorithm to\n");
1980" construct a Delaunay triangulation. Try it if the divide-and-\n");
1981 printf(
" conquer algorithm fails.\n");
1983" -F Uses Steven Fortune's sweepline algorithm to construct a Delaunay\n");
1985" triangulation. Warning: does not use exact arithmetic for all\n");
1986 printf(
" calculations. An exact result is not guaranteed.\n");
1988" -l Uses only vertical cuts in the divide-and-conquer algorithm. By\n");
1990" default, Triangle alternates between vertical and horizontal cuts,\n"
1993" which usually improve the speed except with vertex sets that are\n");
1995" small or short and wide. This switch is primarily of theoretical\n");
1996 printf(
" interest.\n");
1998" -s Specifies that segments should be forced into the triangulation by\n"
2001" recursively splitting them at their midpoints, rather than by\n");
2003" generating a constrained Delaunay triangulation. Segment splitting\n"
2006" is true to Ruppert's original algorithm, but can create needlessly\n"
2009" small triangles. This switch is primarily of theoretical interest.\n"
2012" -C Check the consistency of the final mesh. Uses exact arithmetic for\n"
2015" checking, even if the -X switch is used. Useful if you suspect\n");
2016 printf(
" Triangle is buggy.\n");
2018" -Q Quiet: Suppresses all explanation of what Triangle is doing,\n");
2019 printf(
" unless an error occurs.\n");
2021" -V Verbose: Gives detailed information about what Triangle is doing.\n"
2024" Add more `V's for increasing amount of detail. `-V' is most\n");
2026" useful; itgives information on algorithmic progress and much more\n");
2028" detailed statistics. `-VV' gives vertex-by-vertex details, and\n");
2030" prints so much that Triangle runs much more slowly. `-VVVV' gives\n"
2032 printf(
" information only a debugger could love.\n");
2033 printf(
" -h Help: Displays these instructions.\n");
2035 printf(
"Definitions:\n");
2038" A Delaunay triangulation of a vertex set is a triangulation whose\n");
2040" vertices are the vertex set, that covers the convex hull of the vertex\n");
2042" set. A Delaunay triangulation has the property that no vertex lies\n");
2044" inside the circumscribing circle (circle that passes through all three\n");
2045 printf(
" vertices) of any triangle in the triangulation.\n\n");
2047" A Voronoi diagram of a vertex set is a subdivision of the plane into\n");
2049" polygonal cells (some of which may be unbounded, meaning infinitely\n");
2051" large), where each cell is the set of points in the plane that are closer\n"
2054" to some input vertex than to any other input vertex. The Voronoi diagram\n"
2056 printf(
" is a geometric dual of the Delaunay triangulation.\n\n");
2058" A Planar Straight Line Graph (PSLG) is a set of vertices and segments.\n");
2060" Segments are simply edges, whose endpoints are all vertices in the PSLG.\n"
2063" Segments may intersect each other only at their endpoints. The file\n");
2064 printf(
" format for PSLGs (.poly files) is described below.\n\n");
2066" A constrained Delaunay triangulation (CDT) of a PSLG is similar to a\n");
2068" Delaunay triangulation, but each PSLG segment is present as a single edge\n"
2071" of the CDT. (A constrained Delaunay triangulation is not truly a\n");
2073" Delaunay triangulation, because some of its triangles might not be\n");
2075" Delaunay.) By definition, a CDT does not have any vertices other than\n");
2077" those specified in the input PSLG. Depending on context, a CDT might\n");
2079" cover the convex hull of the PSLG, or it might cover only a segment-\n");
2080 printf(
" bounded region (e.g. a polygon).\n\n");
2082" A conforming Delaunay triangulation of a PSLG is a triangulation in which\n"
2085" each triangle is truly Delaunay, and each PSLG segment is represented by\n"
2088" a linear contiguous sequence of edges of the triangulation. New vertices\n"
2091" (not part of the PSLG) may appear, and each input segment may have been\n");
2093" subdivided into shorter edges (subsegments) by these additional vertices.\n"
2096" The new vertices are frequently necessary to maintain the Delaunay\n");
2097 printf(
" property while ensuring that every segment is represented.\n\n");
2099" A conforming constrained Delaunay triangulation (CCDT) of a PSLG is a\n");
2101" triangulation of a PSLG whose triangles are constrained Delaunay. New\n");
2102 printf(
" vertices may appear, and input segments may be subdivided into\n");
2104" subsegments, but not to guarantee that segments are respected; rather, to\n"
2107" improve the quality of the triangles. The high-quality meshes produced\n");
2109" by the -q switch are usually CCDTs, but can be made conforming Delaunay\n");
2110 printf(
" with the -D switch.\n\n");
2111 printf(
"File Formats:\n\n");
2113" All files may contain comments prefixed by the character '#'. Vertices,\n"
2116" triangles, edges, holes, and maximum area constraints must be numbered\n");
2118" consecutively, starting from either 1 or 0. Whichever you choose, all\n");
2120" input files must be consistent; if the vertices are numbered from 1, so\n");
2122" must be all other objects. Triangle automatically detects your choice\n");
2124" while reading the .node (or .poly) file. (When calling Triangle from\n");
2126" another program, use the -z switch if you wish to number objects from\n");
2127 printf(
" zero.) Examples of these file formats are given below.\n\n");
2128 printf(
" .node files:\n");
2130" First line: <# of vertices> <dimension (must be 2)> <# of attributes>\n"
2133" <# of boundary markers (0 or 1)>\n"
2136" Remaining lines: <vertex #> <x> <y> [attributes] [boundary marker]\n");
2139" The attributes, which are typically floating-point values of physical\n");
2141" quantities (such as mass or conductivity) associated with the nodes of\n"
2144" a finite element mesh, are copied unchanged to the output mesh. If -q,\n"
2147" -a, -u, -D, or -s is selected, each new Steiner point added to the mesh\n"
2149 printf(
" has attributes assigned to it by linear interpolation.\n\n");
2151" If the fourth entry of the first line is `1', the last column of the\n");
2153" remainder of the file is assumed to contain boundary markers. Boundary\n"
2156" markers are used to identify boundary vertices and vertices resting on\n"
2159" PSLG segments; a complete description appears in a section below. The\n"
2162" .node file produced by Triangle contains boundary markers in the last\n");
2163 printf(
" column unless they are suppressed by the -B switch.\n\n");
2164 printf(
" .ele files:\n");
2166" First line: <# of triangles> <nodes per triangle> <# of attributes>\n");
2168" Remaining lines: <triangle #> <node> <node> <node> ... [attributes]\n");
2171" Nodes are indices into the corresponding .node file. The first three\n");
2173" nodes are the corner vertices, and are listed in counterclockwise order\n"
2176" around each triangle. (The remaining nodes, if any, depend on the type\n"
2178 printf(
" of finite element used.)\n\n");
2180" The attributes are just like those of .node files. Because there is no\n"
2183" simple mapping from input to output triangles, Triangle attempts to\n");
2185" interpolate attributes, and may cause a lot of diffusion of attributes\n"
2188" among nearby triangles as the triangulation is refined. Attributes do\n"
2190 printf(
" not diffuse across segments, so attributes used to identify\n");
2191 printf(
" segment-bounded regions remain intact.\n\n");
2193" In .ele files produced by Triangle, each triangular element has three\n");
2195" nodes (vertices) unless the -o2 switch is used, in which case\n");
2197" subparametric quadratic elements with six nodes each are generated.\n");
2199" The first three nodes are the corners in counterclockwise order, and\n");
2201" the fourth, fifth, and sixth nodes lie on the midpoints of the edges\n");
2203" opposite the first, second, and third vertices, respectively.\n");
2205 printf(
" .poly files:\n");
2207" First line: <# of vertices> <dimension (must be 2)> <# of attributes>\n"
2210" <# of boundary markers (0 or 1)>\n"
2213" Following lines: <vertex #> <x> <y> [attributes] [boundary marker]\n");
2214 printf(
" One line: <# of segments> <# of boundary markers (0 or 1)>\n");
2216" Following lines: <segment #> <endpoint> <endpoint> [boundary marker]\n");
2217 printf(
" One line: <# of holes>\n");
2218 printf(
" Following lines: <hole #> <x> <y>\n");
2220" Optional line: <# of regional attributes and/or area constraints>\n");
2222" Optional following lines: <region #> <x> <y> <attribute> <max area>\n");
2225" A .poly file represents a PSLG, as well as some additional information.\n"
2228" The first section lists all the vertices, and is identical to the\n");
2230" format of .node files. <# of vertices> may be set to zero to indicate\n"
2233" that the vertices are listed in a separate .node file; .poly files\n");
2235" produced by Triangle always have this format. A vertex set represented\n"
2238" this way has the advantage that it may easily be triangulated with or\n");
2240" without segments (depending on whether the -p switch is invoked).\n");
2243" The second section lists the segments. Segments are edges whose\n");
2245" presence in the triangulation is enforced. (Depending on the choice of\n"
2248" switches, segment might be subdivided into smaller edges). Each\n");
2250" segment is specified by listing the indices of its two endpoints. This\n"
2253" means that you must include its endpoints in the vertex list. Each\n");
2254 printf(
" segment, like each point, may have a boundary marker.\n\n");
2256" If -q, -a, -u, and -s are not selected, Triangle produces a constrained\n"
2259" Delaunay triangulation (CDT), in which each segment appears as a single\n"
2262" edge in the triangulation. If -q, -a, -u, or -s is selected, Triangle\n"
2265" produces a conforming constrained Delaunay triangulation (CCDT), in\n");
2267" which segments may be subdivided into smaller edges. If -D is\n");
2269" selected, Triangle produces a conforming Delaunay triangulation, so\n");
2271" that every triangle is Delaunay, and not just constrained Delaunay.\n");
2274" The third section lists holes (and concavities, if -c is selected) in\n");
2276" the triangulation. Holes are specified by identifying a point inside\n");
2278" each hole. After the triangulation is formed, Triangle creates holes\n");
2280" by eating triangles, spreading out from each hole point until its\n");
2282" progress is blocked by segments in the PSLG. You must be careful to\n");
2284" enclose each hole in segments, or your whole triangulation might be\n");
2286" eaten away. If the two triangles abutting a segment are eaten, the\n");
2288" segment itself is also eaten. Do not place a hole directly on a\n");
2289 printf(
" segment; if you do, Triangle chooses one side of the segment\n");
2290 printf(
" arbitrarily.\n\n");
2292" The optional fourth section lists regional attributes (to be assigned\n");
2294" to all triangles in a region) and regional constraints on the maximum\n");
2296" triangle area. Triangle reads this section only if the -A switch is\n");
2298" used or the -a switch is used without a number following it, and the -r\n"
2301" switch is not used. Regional attributes and area constraints are\n");
2303" propagated in the same manner as holes: you specify a point for each\n");
2305" attribute and/or constraint, and the attribute and/or constraint\n");
2307" affects the whole region (bounded by segments) containing the point.\n");
2309" If two values are written on a line after the x and y coordinate, the\n");
2311" first such value is assumed to be a regional attribute (but is only\n");
2313" applied if the -A switch is selected), and the second value is assumed\n"
2316" to be a regional area constraint (but is only applied if the -a switch\n"
2319" is selected). You may specify just one value after the coordinates,\n");
2321" which can serve as both an attribute and an area constraint, depending\n"
2324" on the choice of switches. If you are using the -A and -a switches\n");
2326" simultaneously and wish to assign an attribute to some region without\n");
2327 printf(
" imposing an area constraint, use a negative maximum area.\n\n");
2329" When a triangulation is created from a .poly file, you must either\n");
2331" enclose the entire region to be triangulated in PSLG segments, or\n");
2333" use the -c switch, which automatically creates extra segments that\n");
2335" enclose the convex hull of the PSLG. If you do not use the -c switch,\n"
2338" Triangle eats all triangles that are not enclosed by segments; if you\n");
2340" are not careful, your whole triangulation may be eaten away. If you do\n"
2343" use the -c switch, you can still produce concavities by the appropriate\n"
2346" placement of holes just inside the boundary of the convex hull.\n");
2349" An ideal PSLG has no intersecting segments, nor any vertices that lie\n");
2351" upon segments (except, of course, the endpoints of each segment). You\n"
2354" aren't required to make your .poly files ideal, but you should be aware\n"
2357" of what can go wrong. Segment intersections are relatively safe--\n");
2359" Triangle calculates the intersection points for you and adds them to\n");
2361" the triangulation--as long as your machine's floating-point precision\n");
2363" doesn't become a problem. You are tempting the fates if you have three\n"
2366" segments that cross at the same location, and expect Triangle to figure\n"
2369" out where the intersection point is. Thanks to floating-point roundoff\n"
2372" error, Triangle will probably decide that the three segments intersect\n"
2375" at three different points, and you will find a minuscule triangle in\n");
2377" your output--unless Triangle tries to refine the tiny triangle, uses\n");
2379" up the last bit of machine precision, and fails to terminate at all.\n");
2381" You're better off putting the intersection point in the input files,\n");
2383" and manually breaking up each segment into two. Similarly, if you\n");
2385" place a vertex at the middle of a segment, and hope that Triangle will\n"
2388" break up the segment at that vertex, you might get lucky. On the other\n"
2391" hand, Triangle might decide that the vertex doesn't lie precisely on\n");
2393" the segment, and you'll have a needle-sharp triangle in your output--or\n"
2395 printf(
" a lot of tiny triangles if you're generating a quality mesh.\n");
2398" When Triangle reads a .poly file, it also writes a .poly file, which\n");
2400" includes all the subsegments--the edges that are parts of input\n");
2402" segments. If the -c switch is used, the output .poly file also\n");
2404" includes all of the edges on the convex hull. Hence, the output .poly\n"
2407" file is useful for finding edges associated with input segments and for\n"
2410" setting boundary conditions in finite element simulations. Moreover,\n");
2412" you will need the output .poly file if you plan to refine the output\n");
2414" mesh, and don't want segments to be missing in later triangulations.\n");
2416 printf(
" .area files:\n");
2417 printf(
" First line: <# of triangles>\n");
2418 printf(
" Following lines: <triangle #> <maximum area>\n");
2421" An .area file associates with each triangle a maximum area that is used\n"
2424" for mesh refinement. As with other file formats, every triangle must\n");
2426" be represented, and the triangles must be numbered consecutively. A\n");
2428" triangle may be left unconstrained by assigning it a negative maximum\n");
2429 printf(
" area.\n\n");
2430 printf(
" .edge files:\n");
2431 printf(
" First line: <# of edges> <# of boundary markers (0 or 1)>\n");
2433" Following lines: <edge #> <endpoint> <endpoint> [boundary marker]\n");
2436" Endpoints are indices into the corresponding .node file. Triangle can\n"
2439" produce .edge files (use the -e switch), but cannot read them. The\n");
2441" optional column of boundary markers is suppressed by the -B switch.\n");
2444" In Voronoi diagrams, one also finds a special kind of edge that is an\n");
2446" infinite ray with only one endpoint. For these edges, a different\n");
2447 printf(
" format is used:\n\n");
2448 printf(
" <edge #> <endpoint> -1 <direction x> <direction y>\n\n");
2450" The `direction' is a floating-point vector that indicates the direction\n"
2452 printf(
" of the infinite ray.\n\n");
2453 printf(
" .neigh files:\n");
2455" First line: <# of triangles> <# of neighbors per triangle (always 3)>\n"
2458" Following lines: <triangle #> <neighbor> <neighbor> <neighbor>\n");
2461" Neighbors are indices into the corresponding .ele file. An index of -1\n"
2464" indicates no neighbor (because the triangle is on an exterior\n");
2466" boundary). The first neighbor of triangle i is opposite the first\n");
2467 printf(
" corner of triangle i, and so on.\n\n");
2469" Triangle can produce .neigh files (use the -n switch), but cannot read\n"
2471 printf(
" them.\n\n");
2472 printf(
"Boundary Markers:\n\n");
2474" Boundary markers are tags used mainly to identify which output vertices\n");
2476" and edges are associated with which PSLG segment, and to identify which\n");
2478" vertices and edges occur on a boundary of the triangulation. A common\n");
2480" use is to determine where boundary conditions should be applied to a\n");
2482" finite element mesh. You can prevent boundary markers from being written\n"
2484 printf(
" into files produced by Triangle by using the -B switch.\n\n");
2486" The boundary marker associated with each segment in an output .poly file\n"
2488 printf(
" and each edge in an output .edge file is chosen as follows:\n");
2490" - If an output edge is part or all of a PSLG segment with a nonzero\n");
2492" boundary marker, then the edge is assigned the same marker.\n");
2494" - Otherwise, if the edge lies on a boundary of the triangulation\n");
2496" (even the boundary of a hole), then the edge is assigned the marker\n");
2497 printf(
" one (1).\n");
2498 printf(
" - Otherwise, the edge is assigned the marker zero (0).\n");
2500" The boundary marker associated with each vertex in an output .node file\n");
2501 printf(
" is chosen as follows:\n");
2503" - If a vertex is assigned a nonzero boundary marker in the input file,\n"
2506" then it is assigned the same marker in the output .node file.\n");
2508" - Otherwise, if the vertex lies on a PSLG segment (even if it is an\n");
2510" endpoint of the segment) with a nonzero boundary marker, then the\n");
2512" vertex is assigned the same marker. If the vertex lies on several\n");
2513 printf(
" such segments, one of the markers is chosen arbitrarily.\n");
2515" - Otherwise, if the vertex occurs on a boundary of the triangulation,\n");
2516 printf(
" then the vertex is assigned the marker one (1).\n");
2517 printf(
" - Otherwise, the vertex is assigned the marker zero (0).\n");
2520" If you want Triangle to determine for you which vertices and edges are on\n"
2523" the boundary, assign them the boundary marker zero (or use no markers at\n"
2526" all) in your input files. In the output files, all boundary vertices,\n");
2527 printf(
" edges, and segments will be assigned the value one.\n\n");
2528 printf(
"Triangulation Iteration Numbers:\n\n");
2530" Because Triangle can read and refine its own triangulations, input\n");
2532" and output files have iteration numbers. For instance, Triangle might\n");
2534" read the files mesh.3.node, mesh.3.ele, and mesh.3.poly, refine the\n");
2536" triangulation, and output the files mesh.4.node, mesh.4.ele, and\n");
2537 printf(
" mesh.4.poly. Files with no iteration number are treated as if\n");
2539" their iteration number is zero; hence, Triangle might read the file\n");
2541" points.node, triangulate it, and produce the files points.1.node and\n");
2542 printf(
" points.1.ele.\n\n");
2544" Iteration numbers allow you to create a sequence of successively finer\n");
2546" meshes suitable for multigrid methods. They also allow you to produce a\n"
2549" sequence of meshes using error estimate-driven mesh refinement.\n");
2552" If you're not using refinement or quality meshing, and you don't like\n");
2554" iteration numbers, use the -I switch to disable them. This switch also\n");
2556" disables output of .node and .poly files to prevent your input files from\n"
2559" being overwritten. (If the input is a .poly file that contains its own\n");
2561" points, a .node file is written. This can be quite convenient for\n");
2562 printf(
" computing CDTs or quality meshes.)\n\n");
2563 printf(
"Examples of How to Use Triangle:\n\n");
2565" `triangle dots' reads vertices from dots.node, and writes their Delaunay\n"
2568" triangulation to dots.1.node and dots.1.ele. (dots.1.node is identical\n");
2570" to dots.node.) `triangle -I dots' writes the triangulation to dots.ele\n");
2572" instead. (No additional .node file is needed, so none is written.)\n");
2575" `triangle -pe object.1' reads a PSLG from object.1.poly (and possibly\n");
2577" object.1.node, if the vertices are omitted from object.1.poly) and writes\n"
2580" its constrained Delaunay triangulation to object.2.node and object.2.ele.\n"
2583" The segments are copied to object.2.poly, and all edges are written to\n");
2584 printf(
" object.2.edge.\n\n");
2586" `triangle -pq31.5a.1 object' reads a PSLG from object.poly (and possibly\n"
2589" object.node), generates a mesh whose angles are all between 31.5 and 117\n"
2592" degrees and whose triangles all have areas of 0.1 or less, and writes the\n"
2595" mesh to object.1.node and object.1.ele. Each segment may be broken up\n");
2596 printf(
" into multiple subsegments; these are written to object.1.poly.\n");
2599" Here is a sample file `box.poly' describing a square with a square hole:\n"
2603" # A box with eight vertices in 2D, no attributes, one boundary marker.\n"
2605 printf(
" 8 2 0 1\n");
2606 printf(
" # Outer box has these vertices:\n");
2607 printf(
" 1 0 0 0\n");
2608 printf(
" 2 0 3 0\n");
2609 printf(
" 3 3 0 0\n");
2610 printf(
" 4 3 3 33 # A special marker for this vertex.\n");
2611 printf(
" # Inner square has these vertices:\n");
2612 printf(
" 5 1 1 0\n");
2613 printf(
" 6 1 2 0\n");
2614 printf(
" 7 2 1 0\n");
2615 printf(
" 8 2 2 0\n");
2616 printf(
" # Five segments with boundary markers.\n");
2618 printf(
" 1 1 2 5 # Left side of outer box.\n");
2619 printf(
" # Square hole has these segments:\n");
2620 printf(
" 2 5 7 0\n");
2621 printf(
" 3 7 8 0\n");
2622 printf(
" 4 8 6 10\n");
2623 printf(
" 5 6 5 0\n");
2624 printf(
" # One hole in the middle of the inner square.\n");
2626 printf(
" 1 1.5 1.5\n");
2629" Note that some segments are missing from the outer square, so you must\n");
2631" use the `-c' switch. After `triangle -pqc box.poly', here is the output\n"
2634" file `box.1.node', with twelve vertices. The last four vertices were\n");
2636" added to meet the angle constraint. Vertices 1, 2, and 9 have markers\n");
2638" from segment 1. Vertices 6 and 8 have markers from segment 4. All the\n");
2640" other vertices but 4 have been marked to indicate that they lie on a\n");
2641 printf(
" boundary.\n\n");
2642 printf(
" 12 2 0 1\n");
2643 printf(
" 1 0 0 5\n");
2644 printf(
" 2 0 3 5\n");
2645 printf(
" 3 3 0 1\n");
2646 printf(
" 4 3 3 33\n");
2647 printf(
" 5 1 1 1\n");
2648 printf(
" 6 1 2 10\n");
2649 printf(
" 7 2 1 1\n");
2650 printf(
" 8 2 2 10\n");
2651 printf(
" 9 0 1.5 5\n");
2652 printf(
" 10 1.5 0 1\n");
2653 printf(
" 11 3 1.5 1\n");
2654 printf(
" 12 1.5 3 1\n");
2655 printf(
" # Generated by triangle -pqc box.poly\n");
2657 printf(
" Here is the output file `box.1.ele', with twelve triangles.\n");
2659 printf(
" 12 3 0\n");
2660 printf(
" 1 5 6 9\n");
2661 printf(
" 2 10 3 7\n");
2662 printf(
" 3 6 8 12\n");
2663 printf(
" 4 9 1 5\n");
2664 printf(
" 5 6 2 9\n");
2665 printf(
" 6 7 3 11\n");
2666 printf(
" 7 11 4 8\n");
2667 printf(
" 8 7 5 10\n");
2668 printf(
" 9 12 2 6\n");
2669 printf(
" 10 8 7 11\n");
2670 printf(
" 11 5 1 10\n");
2671 printf(
" 12 8 4 12\n");
2672 printf(
" # Generated by triangle -pqc box.poly\n\n");
2674" Here is the output file `box.1.poly'. Note that segments have been added\n"
2677" to represent the convex hull, and some segments have been subdivided by\n");
2679" newly added vertices. Note also that <# of vertices> is set to zero to\n");
2680 printf(
" indicate that the vertices should be read from the .node file.\n");
2682 printf(
" 0 2 0 1\n");
2684 printf(
" 1 1 9 5\n");
2685 printf(
" 2 5 7 1\n");
2686 printf(
" 3 8 7 1\n");
2687 printf(
" 4 6 8 10\n");
2688 printf(
" 5 5 6 1\n");
2689 printf(
" 6 3 10 1\n");
2690 printf(
" 7 4 11 1\n");
2691 printf(
" 8 2 12 1\n");
2692 printf(
" 9 9 2 5\n");
2693 printf(
" 10 10 1 1\n");
2694 printf(
" 11 11 3 1\n");
2695 printf(
" 12 12 4 1\n");
2697 printf(
" 1 1.5 1.5\n");
2698 printf(
" # Generated by triangle -pqc box.poly\n");
2700 printf(
"Refinement and Area Constraints:\n");
2703" The -r switch causes a mesh (.node and .ele files) to be read and\n");
2705" refined. If the -p switch is also used, a .poly file is read and used to\n"
2708" specify edges that are constrained and cannot be eliminated (although\n");
2710" they can be subdivided into smaller edges) by the refinement process.\n");
2713" When you refine a mesh, you generally want to impose tighter constraints.\n"
2716" One way to accomplish this is to use -q with a larger angle, or -a\n");
2718" followed by a smaller area than you used to generate the mesh you are\n");
2720" refining. Another way to do this is to create an .area file, which\n");
2722" specifies a maximum area for each triangle, and use the -a switch\n");
2724" (without a number following). Each triangle's area constraint is applied\n"
2727" to that triangle. Area constraints tend to diffuse as the mesh is\n");
2729" refined, so if there are large variations in area constraint between\n");
2731" adjacent triangles, you may not get the results you want. In that case,\n"
2734" consider instead using the -u switch and writing a C procedure that\n");
2735 printf(
" determines which triangles are too large.\n\n");
2737" If you are refining a mesh composed of linear (three-node) elements, the\n"
2740" output mesh contains all the nodes present in the input mesh, in the same\n"
2743" order, with new nodes added at the end of the .node file. However, the\n");
2745" refinement is not hierarchical: there is no guarantee that each output\n");
2747" element is contained in a single input element. Often, an output element\n"
2750" can overlap two or three input elements, and some input edges are not\n");
2752" present in the output mesh. Hence, a sequence of refined meshes forms a\n"
2755" hierarchy of nodes, but not a hierarchy of elements. If you refine a\n");
2757" mesh of higher-order elements, the hierarchical property applies only to\n"
2760" the nodes at the corners of an element; the midpoint nodes on each edge\n");
2761 printf(
" are discarded before the mesh is refined.\n\n");
2763" Maximum area constraints in .poly files operate differently from those in\n"
2766" .area files. A maximum area in a .poly file applies to the whole\n");
2768" (segment-bounded) region in which a point falls, whereas a maximum area\n");
2770" in an .area file applies to only one triangle. Area constraints in .poly\n"
2773" files are used only when a mesh is first generated, whereas area\n");
2775" constraints in .area files are used only to refine an existing mesh, and\n"
2778" are typically based on a posteriori error estimates resulting from a\n");
2779 printf(
" finite element simulation on that mesh.\n\n");
2781" `triangle -rq25 object.1' reads object.1.node and object.1.ele, then\n");
2783" refines the triangulation to enforce a 25 degree minimum angle, and then\n"
2786" writes the refined triangulation to object.2.node and object.2.ele.\n");
2789" `triangle -rpaa6.2 z.3' reads z.3.node, z.3.ele, z.3.poly, and z.3.area.\n"
2792" After reconstructing the mesh and its subsegments, Triangle refines the\n");
2794" mesh so that no triangle has area greater than 6.2, and furthermore the\n");
2796" triangles satisfy the maximum area constraints in z.3.area. No angle\n");
2798" bound is imposed at all. The output is written to z.4.node, z.4.ele, and\n"
2800 printf(
" z.4.poly.\n\n");
2802" The sequence `triangle -qa1 x', `triangle -rqa.3 x.1', `triangle -rqa.1\n");
2804" x.2' creates a sequence of successively finer meshes x.1, x.2, and x.3,\n");
2805 printf(
" suitable for multigrid.\n\n");
2806 printf(
"Convex Hulls and Mesh Boundaries:\n\n");
2808" If the input is a vertex set (not a PSLG), Triangle produces its convex\n");
2810" hull as a by-product in the output .poly file if you use the -c switch.\n");
2812" There are faster algorithms for finding a two-dimensional convex hull\n");
2813 printf(
" than triangulation, of course, but this one comes for free.\n\n");
2815" If the input is an unconstrained mesh (you are using the -r switch but\n");
2817" not the -p switch), Triangle produces a list of its boundary edges\n");
2819" (including hole boundaries) as a by-product when you use the -c switch.\n");
2821" If you also use the -p switch, the output .poly file contains all the\n");
2822 printf(
" segments from the input .poly file as well.\n\n");
2823 printf(
"Voronoi Diagrams:\n\n");
2825" The -v switch produces a Voronoi diagram, in files suffixed .v.node and\n");
2827" .v.edge. For example, `triangle -v points' reads points.node, produces\n");
2829" its Delaunay triangulation in points.1.node and points.1.ele, and\n");
2831" produces its Voronoi diagram in points.1.v.node and points.1.v.edge. The\n"
2834" .v.node file contains a list of all Voronoi vertices, and the .v.edge\n");
2836" file contains a list of all Voronoi edges, some of which may be infinite\n"
2839" rays. (The choice of filenames makes it easy to run the set of Voronoi\n");
2840 printf(
" vertices through Triangle, if so desired.)\n\n");
2842" This implementation does not use exact arithmetic to compute the Voronoi\n"
2845" vertices, and does not check whether neighboring vertices are identical.\n"
2848" Be forewarned that if the Delaunay triangulation is degenerate or\n");
2850" near-degenerate, the Voronoi diagram may have duplicate vertices or\n");
2851 printf(
" crossing edges.\n\n");
2853" The result is a valid Voronoi diagram only if Triangle's output is a true\n"
2856" Delaunay triangulation. The Voronoi output is usually meaningless (and\n");
2858" may contain crossing edges and other pathology) if the output is a CDT or\n"
2861" CCDT, or if it has holes or concavities. If the triangulated domain is\n");
2863" convex and has no holes, you can use -D switch to force Triangle to\n");
2865" construct a conforming Delaunay triangulation instead of a CCDT, so the\n");
2866 printf(
" Voronoi diagram will be valid.\n\n");
2867 printf(
"Mesh Topology:\n\n");
2869" You may wish to know which triangles are adjacent to a certain Delaunay\n");
2871" edge in an .edge file, which Voronoi cells are adjacent to a certain\n");
2873" Voronoi edge in a .v.edge file, or which Voronoi cells are adjacent to\n");
2875" each other. All of this information can be found by cross-referencing\n");
2877" output files with the recollection that the Delaunay triangulation and\n");
2878 printf(
" the Voronoi diagram are planar duals.\n\n");
2880" Specifically, edge i of an .edge file is the dual of Voronoi edge i of\n");
2882" the corresponding .v.edge file, and is rotated 90 degrees counterclock-\n");
2884" wise from the Voronoi edge. Triangle j of an .ele file is the dual of\n");
2886" vertex j of the corresponding .v.node file. Voronoi cell k is the dual\n");
2887 printf(
" of vertex k of the corresponding .node file.\n\n");
2889" Hence, to find the triangles adjacent to a Delaunay edge, look at the\n");
2891" vertices of the corresponding Voronoi edge. If the endpoints of a\n");
2893" Voronoi edge are Voronoi vertices 2 and 6 respectively, then triangles 2\n"
2896" and 6 adjoin the left and right sides of the corresponding Delaunay edge,\n"
2899" respectively. To find the Voronoi cells adjacent to a Voronoi edge, look\n"
2902" at the endpoints of the corresponding Delaunay edge. If the endpoints of\n"
2905" a Delaunay edge are input vertices 7 and 12, then Voronoi cells 7 and 12\n"
2908" adjoin the right and left sides of the corresponding Voronoi edge,\n");
2910" respectively. To find which Voronoi cells are adjacent to each other,\n");
2911 printf(
" just read the list of Delaunay edges.\n\n");
2913" Triangle does not write a list of the edges adjoining each Voronoi cell,\n"
2916" but you can reconstructed it straightforwardly. For instance, to find\n");
2918" all the edges of Voronoi cell 1, search the output .edge file for every\n");
2920" edge that has input vertex 1 as an endpoint. The corresponding dual\n");
2922" edges in the output .v.edge file form the boundary of Voronoi cell 1.\n");
2925" For each Voronoi vertex, the .neigh file gives a list of the three\n");
2927" Voronoi vertices attached to it. You might find this more convenient\n");
2928 printf(
" than the .v.edge file.\n\n");
2929 printf(
"Quadratic Elements:\n\n");
2931" Triangle generates meshes with subparametric quadratic elements if the\n");
2933" -o2 switch is specified. Quadratic elements have six nodes per element,\n"
2936" rather than three. `Subparametric' means that the edges of the triangles\n"
2939" are always straight, so that subparametric quadratic elements are\n");
2941" geometrically identical to linear elements, even though they can be used\n"
2944" with quadratic interpolating functions. The three extra nodes of an\n");
2946" element fall at the midpoints of the three edges, with the fourth, fifth,\n"
2949" and sixth nodes appearing opposite the first, second, and third corners\n");
2950 printf(
" respectively.\n\n");
2951 printf(
"Domains with Small Angles:\n\n");
2953" If two input segments adjoin each other at a small angle, clearly the -q\n"
2956" switch cannot remove the small angle. Moreover, Triangle may have no\n");
2958" choice but to generate additional triangles whose smallest angles are\n");
2960" smaller than the specified bound. However, these triangles only appear\n");
2962" between input segments separated by small angles. Moreover, if you\n");
2964" request a minimum angle of theta degrees, Triangle will generally produce\n"
2967" no angle larger than 180 - 2 theta, even if it is forced to compromise on\n"
2969 printf(
" the minimum angle.\n\n");
2970 printf(
"Statistics:\n\n");
2972" After generating a mesh, Triangle prints a count of entities in the\n");
2974" output mesh, including the number of vertices, triangles, edges, exterior\n"
2977" boundary edges (i.e. subsegments on the boundary of the triangulation,\n");
2979" including hole boundaries), interior boundary edges (i.e. subsegments of\n"
2982" input segments not on the boundary), and total subsegments. If you've\n");
2984" forgotten the statistics for an existing mesh, run Triangle on that mesh\n"
2987" with the -rNEP switches to read the mesh and print the statistics without\n"
2990" writing any files. Use -rpNEP if you've got a .poly file for the mesh.\n");
2993" The -V switch produces extended statistics, including a rough estimate\n");
2995" of memory use, the number of calls to geometric predicates, and\n");
2997" histograms of the angles and the aspect ratios of the triangles in the\n");
2998 printf(
" mesh.\n\n");
2999 printf(
"Exact Arithmetic:\n\n");
3001" Triangle uses adaptive exact arithmetic to perform what computational\n");
3003" geometers call the `orientation' and `incircle' tests. If the floating-\n"
3006" point arithmetic of your machine conforms to the IEEE 754 standard (as\n");
3008" most workstations do), and does not use extended precision internal\n");
3010" floating-point registers, then your output is guaranteed to be an\n");
3012" absolutely true Delaunay or constrained Delaunay triangulation, roundoff\n"
3015" error notwithstanding. The word `adaptive' implies that these arithmetic\n"
3018" routines compute the result only to the precision necessary to guarantee\n"
3021" correctness, so they are usually nearly as fast as their approximate\n");
3022 printf(
" counterparts.\n\n");
3024" May CPUs, including Intel x86 processors, have extended precision\n");
3026" floating-point registers. These must be reconfigured so their precision\n"
3029" is reduced to memory precision. Triangle does this if it is compiled\n");
3030 printf(
" correctly. See the makefile for details.\n\n");
3032" The exact tests can be disabled with the -X switch. On most inputs, this\n"
3035" switch reduces the computation time by about eight percent--it's not\n");
3037" worth the risk. There are rare difficult inputs (having many collinear\n");
3039" and cocircular vertices), however, for which the difference in speed\n");
3041" could be a factor of two. Be forewarned that these are precisely the\n");
3043" inputs most likely to cause errors if you use the -X switch. Hence, the\n"
3045 printf(
" -X switch is not recommended.\n\n");
3047" Unfortunately, the exact tests don't solve every numerical problem.\n");
3049" Exact arithmetic is not used to compute the positions of new vertices,\n");
3051" because the bit complexity of vertex coordinates would grow without\n");
3053" bound. Hence, segment intersections aren't computed exactly; in very\n");
3055" unusual cases, roundoff error in computing an intersection point might\n");
3057" actually lead to an inverted triangle and an invalid triangulation.\n");
3059" (This is one reason to specify your own intersection points in your .poly\n"
3062" files.) Similarly, exact arithmetic is not used to compute the vertices\n"
3064 printf(
" of the Voronoi diagram.\n\n");
3066" Another pair of problems not solved by the exact arithmetic routines is\n");
3068" underflow and overflow. If Triangle is compiled for double precision\n");
3070" arithmetic, I believe that Triangle's geometric predicates work correctly\n"
3073" if the exponent of every input coordinate falls in the range [-148, 201].\n"
3076" Underflow can silently prevent the orientation and incircle tests from\n");
3078" being performed exactly, while overflow typically causes a floating\n");
3079 printf(
" exception.\n\n");
3080 printf(
"Calling Triangle from Another Program:\n\n");
3081 printf(
" Read the file triangle.h for details.\n\n");
3082 printf(
"Troubleshooting:\n\n");
3083 printf(
" Please read this section before mailing me bugs.\n\n");
3084 printf(
" `My output mesh has no triangles!'\n\n");
3086" If you're using a PSLG, you've probably failed to specify a proper set\n"
3089" of bounding segments, or forgotten to use the -c switch. Or you may\n");
3091" have placed a hole badly, thereby eating all your triangles. To test\n");
3092 printf(
" these possibilities, try again with the -c and -O switches.\n");
3094" Alternatively, all your input vertices may be collinear, in which case\n"
3096 printf(
" you can hardly expect to triangulate them.\n\n");
3097 printf(
" `Triangle doesn't terminate, or just crashes.'\n\n");
3099" Bad things can happen when triangles get so small that the distance\n");
3101" between their vertices isn't much larger than the precision of your\n");
3103" machine's arithmetic. If you've compiled Triangle for single-precision\n"
3106" arithmetic, you might do better by recompiling it for double-precision.\n"
3109" Then again, you might just have to settle for more lenient constraints\n"
3112" on the minimum angle and the maximum area than you had planned.\n");
3115" You can minimize precision problems by ensuring that the origin lies\n");
3117" inside your vertex set, or even inside the densest part of your\n");
3119" mesh. If you're triangulating an object whose x-coordinates all fall\n");
3121" between 6247133 and 6247134, you're not leaving much floating-point\n");
3122 printf(
" precision for Triangle to work with.\n\n");
3124" Precision problems can occur covertly if the input PSLG contains two\n");
3126" segments that meet (or intersect) at an extremely small angle, or if\n");
3128" such an angle is introduced by the -c switch. If you don't realize\n");
3130" that a tiny angle is being formed, you might never discover why\n");
3132" Triangle is crashing. To check for this possibility, use the -S switch\n"
3135" (with an appropriate limit on the number of Steiner points, found by\n");
3137" trial-and-error) to stop Triangle early, and view the output .poly file\n"
3140" with Show Me (described below). Look carefully for regions where dense\n"
3143" clusters of vertices are forming and for small angles between segments.\n"
3146" Zoom in closely, as such segments might look like a single segment from\n"
3148 printf(
" a distance.\n\n");
3150" If some of the input values are too large, Triangle may suffer a\n");
3152" floating exception due to overflow when attempting to perform an\n");
3154" orientation or incircle test. (Read the section on exact arithmetic\n");
3156" above.) Again, I recommend compiling Triangle for double (rather\n");
3157 printf(
" than single) precision arithmetic.\n\n");
3159" Unexpected problems can arise if you use quality meshing (-q, -a, or\n");
3161" -u) with an input that is not segment-bounded--that is, if your input\n");
3163" is a vertex set, or you're using the -c switch. If the convex hull of\n"
3166" your input vertices has collinear vertices on its boundary, an input\n");
3168" vertex that you think lies on the convex hull might actually lie just\n");
3170" inside the convex hull. If so, the vertex and the nearby convex hull\n");
3172" edge form an extremely thin triangle. When Triangle tries to refine\n");
3174" the mesh to enforce angle and area constraints, Triangle might generate\n"
3177" extremely tiny triangles, or it might fail because of insufficient\n");
3178 printf(
" floating-point precision.\n\n");
3180" `The numbering of the output vertices doesn't match the input vertices.'\n"
3184" You may have had duplicate input vertices, or you may have eaten some\n");
3186" of your input vertices with a hole, or by placing them outside the area\n"
3189" enclosed by segments. In any case, you can solve the problem by not\n");
3190 printf(
" using the -j switch.\n\n");
3192" `Triangle executes without incident, but when I look at the resulting\n");
3194" mesh, it has overlapping triangles or other geometric inconsistencies.'\n");
3197" If you select the -X switch, Triangle occasionally makes mistakes due\n");
3199" to floating-point roundoff error. Although these errors are rare,\n");
3201" don't use the -X switch. If you still have problems, please report the\n"
3203 printf(
" bug.\n\n");
3205" `Triangle executes without incident, but when I look at the resulting\n");
3206 printf(
" Voronoi diagram, it has overlapping edges or other geometric\n");
3207 printf(
" inconsistencies.'\n");
3210" If your input is a PSLG (-p), you can only expect a meaningful Voronoi\n"
3213" diagram if the domain you are triangulating is convex and free of\n");
3215" holes, and you use the -D switch to construct a conforming Delaunay\n");
3216 printf(
" triangulation (instead of a CDT or CCDT).\n\n");
3218" Strange things can happen if you've taken liberties with your PSLG. Do\n");
3220" you have a vertex lying in the middle of a segment? Triangle sometimes\n");
3222" copes poorly with that sort of thing. Do you want to lay out a collinear\n"
3225" row of evenly spaced, segment-connected vertices? Have you simply\n");
3227" defined one long segment connecting the leftmost vertex to the rightmost\n"
3230" vertex, and a bunch of vertices lying along it? This method occasionally\n"
3233" works, especially with horizontal and vertical lines, but often it\n");
3235" doesn't, and you'll have to connect each adjacent pair of vertices with a\n"
3237 printf(
" separate segment. If you don't like it, tough.\n\n");
3239" Furthermore, if you have segments that intersect other than at their\n");
3241" endpoints, try not to let the intersections fall extremely close to PSLG\n"
3243 printf(
" vertices or each other.\n\n");
3245" If you have problems refining a triangulation not produced by Triangle:\n");
3247" Are you sure the triangulation is geometrically valid? Is it formatted\n");
3249" correctly for Triangle? Are the triangles all listed so the first three\n"
3252" vertices are their corners in counterclockwise order? Are all of the\n");
3254" triangles constrained Delaunay? Triangle's Delaunay refinement algorithm\n"
3256 printf(
" assumes that it starts with a CDT.\n\n");
3257 printf(
"Show Me:\n\n");
3259" Triangle comes with a separate program named `Show Me', whose primary\n");
3261" purpose is to draw meshes on your screen or in PostScript. Its secondary\n"
3264" purpose is to check the validity of your input files, and do so more\n");
3266" thoroughly than Triangle does. Unlike Triangle, Show Me requires that\n");
3268" you have the X Windows system. Sorry, Microsoft Windows users.\n");
3270 printf(
"Triangle on the Web:\n");
3272 printf(
" To see an illustrated version of these instructions, check out\n");
3274 printf(
" http://www.cs.cmu.edu/~quake/triangle.html\n");
3276 printf(
"A Brief Plea:\n");
3279" If you use Triangle, and especially if you use it to accomplish real\n");
3281" work, I would like very much to hear from you. A short letter or email\n");
3283" (to jrs@cs.berkeley.edu) describing how you use Triangle will mean a lot\n"
3286" to me. The more people I know are using this program, the more easily I\n"
3289" can justify spending time on improvements, which in turn will benefit\n");
3291" you. Also, I can put you on a list to receive email whenever a new\n");
3292 printf(
" version of Triangle is available.\n\n");
3294" If you use a mesh generated by Triangle in a publication, please include\n"
3297" an acknowledgment as well. And please spell Triangle with a capital `T'!\n"
3300" If you want to include a citation, use `Jonathan Richard Shewchuk,\n");
3302" ``Triangle: Engineering a 2D Quality Mesh Generator and Delaunay\n");
3304" Triangulator,'' in Applied Computational Geometry: Towards Geometric\n");
3306" Engineering (Ming C. Lin and Dinesh Manocha, editors), volume 1148 of\n");
3308" Lecture Notes in Computer Science, pages 203-222, Springer-Verlag,\n");
3310" Berlin, May 1996. (From the First ACM Workshop on Applied Computational\n"
3312 printf(
" Geometry.)'\n\n");
3313 printf(
"Research credit:\n\n");
3315" Of course, I can take credit for only a fraction of the ideas that made\n");
3317" this mesh generator possible. Triangle owes its existence to the efforts\n"
3320" of many fine computational geometers and other researchers, including\n");
3322" Marshall Bern, L. Paul Chew, Kenneth L. Clarkson, Boris Delaunay, Rex A.\n"
3325" Dwyer, David Eppstein, Steven Fortune, Leonidas J. Guibas, Donald E.\n");
3327" Knuth, Charles L. Lawson, Der-Tsai Lee, Gary L. Miller, Ernst P. Mucke,\n");
3329" Steven E. Pav, Douglas M. Priest, Jim Ruppert, Isaac Saias, Bruce J.\n");
3331" Schachter, Micha Sharir, Peter W. Shor, Daniel D. Sleator, Jorge Stolfi,\n"
3333 printf(
" Robert E. Tarjan, Alper Ungor, Christopher J. Van Wyk, Noel J.\n");
3335" Walkington, and Binhai Zhu. See the comments at the beginning of the\n");
3336 printf(
" source code for references.\n\n");
3350 printf(
" Please report this bug to jrs@cs.berkeley.edu\n");
3351 printf(
" Include the message above, your input data set, and the exact\n");
3352 printf(
" command line you used to run Triangle.\n");
3363#ifdef ANSI_DECLARATORS
3364void parsecommandline(
int argc,
char **argv,
struct behavior *b)
3366void parsecommandline(argc, argv, b)
3381 char workstring[FILENAMESIZE];
3383 b->poly = b->refine = b->quality = 0;
3384 b->vararea = b->fixedarea = b->usertest = 0;
3385 b->regionattrib = b->convex = b->weighted = b->jettison = 0;
3387 b->edgesout = b->voronoi = b->neighbors = b->geomview = 0;
3388 b->nobound = b->nopolywritten = b->nonodewritten = b->noelewritten = 0;
3389 b->noiterationnum = 0;
3390 b->noholes = b->noexact = 0;
3391 b->incremental = b->sweepline = 0;
3401 b->quiet = b->verbose = 0;
3403 b->innodefilename[0] =
'\0';
3406 for (i = STARTINDEX; i < argc; i++) {
3408 if (argv[i][0] ==
'-') {
3410 for (j = STARTINDEX; argv[i][j] !=
'\0'; j++) {
3411 if (argv[i][j] ==
'p') {
3415 if (argv[i][j] ==
'r') {
3418 if (argv[i][j] ==
'q') {
3420 if (((argv[i][j + 1] >=
'0') && (argv[i][j + 1] <=
'9')) ||
3421 (argv[i][j + 1] ==
'.')) {
3423 while (((argv[i][j + 1] >=
'0') && (argv[i][j + 1] <=
'9')) ||
3424 (argv[i][j + 1] ==
'.')) {
3426 workstring[k] = argv[i][j];
3429 workstring[k] =
'\0';
3430 b->minangle = (REAL) strtod(workstring, (
char **) NULL);
3435 if (argv[i][j] ==
'a') {
3437 if (((argv[i][j + 1] >=
'0') && (argv[i][j + 1] <=
'9')) ||
3438 (argv[i][j + 1] ==
'.')) {
3441 while (((argv[i][j + 1] >=
'0') && (argv[i][j + 1] <=
'9')) ||
3442 (argv[i][j + 1] ==
'.')) {
3444 workstring[k] = argv[i][j];
3447 workstring[k] =
'\0';
3448 b->maxarea = (REAL) strtod(workstring, (
char **) NULL);
3449 if (b->maxarea <= 0.0) {
3450 printf(
"Error: Maximum area must be greater than zero.\n");
3457 if (argv[i][j] ==
'u') {
3462 if (argv[i][j] ==
'A') {
3463 b->regionattrib = 1;
3465 if (argv[i][j] ==
'c') {
3468 if (argv[i][j] ==
'w') {
3471 if (argv[i][j] ==
'W') {
3474 if (argv[i][j] ==
'j') {
3477 if (argv[i][j] ==
'z') {
3480 if (argv[i][j] ==
'e') {
3483 if (argv[i][j] ==
'v') {
3486 if (argv[i][j] ==
'n') {
3489 if (argv[i][j] ==
'g') {
3492 if (argv[i][j] ==
'B') {
3495 if (argv[i][j] ==
'P') {
3496 b->nopolywritten = 1;
3498 if (argv[i][j] ==
'N') {
3499 b->nonodewritten = 1;
3501 if (argv[i][j] ==
'E') {
3502 b->noelewritten = 1;
3505 if (argv[i][j] ==
'I') {
3506 b->noiterationnum = 1;
3509 if (argv[i][j] ==
'O') {
3512 if (argv[i][j] ==
'X') {
3515 if (argv[i][j] ==
'o') {
3516 if (argv[i][j + 1] ==
'2') {
3522 if (argv[i][j] ==
'Y') {
3525 if (argv[i][j] ==
'S') {
3527 while ((argv[i][j + 1] >=
'0') && (argv[i][j + 1] <=
'9')) {
3529 b->steiner = b->steiner * 10 + (int) (argv[i][j] -
'0');
3534 if (argv[i][j] ==
'i') {
3537 if (argv[i][j] ==
'F') {
3541 if (argv[i][j] ==
'l') {
3546 if (argv[i][j] ==
's') {
3549 if ((argv[i][j] ==
'D') || (argv[i][j] ==
'L')) {
3554 if (argv[i][j] ==
'C') {
3558 if (argv[i][j] ==
'Q') {
3561 if (argv[i][j] ==
'V') {
3565 if ((argv[i][j] ==
'h') || (argv[i][j] ==
'H') ||
3566 (argv[i][j] ==
'?')) {
3573 strncpy(b->innodefilename, argv[i], FILENAMESIZE - 1);
3574 b->innodefilename[FILENAMESIZE - 1] =
'\0';
3579 if (b->innodefilename[0] ==
'\0') {
3582 if (!strcmp(&b->innodefilename[strlen(b->innodefilename) - 5],
".node")) {
3583 b->innodefilename[strlen(b->innodefilename) - 5] =
'\0';
3585 if (!strcmp(&b->innodefilename[strlen(b->innodefilename) - 5],
".poly")) {
3586 b->innodefilename[strlen(b->innodefilename) - 5] =
'\0';
3590 if (!strcmp(&b->innodefilename[strlen(b->innodefilename) - 4],
".ele")) {
3591 b->innodefilename[strlen(b->innodefilename) - 4] =
'\0';
3594 if (!strcmp(&b->innodefilename[strlen(b->innodefilename) - 5],
".area")) {
3595 b->innodefilename[strlen(b->innodefilename) - 5] =
'\0';
3602 b->usesegments = b->poly || b->refine || b->quality || b->convex;
3603 b->goodangle = cos(b->minangle * PI / 180.0);
3604 if (b->goodangle == 1.0) {
3605 b->offconstant = 0.0;
3607 b->offconstant = 0.475 * sqrt((1.0 + b->goodangle) / (1.0 - b->goodangle));
3609 b->goodangle *= b->goodangle;
3610 if (b->refine && b->noiterationnum) {
3612 "Error: You cannot use the -I switch when refining a triangulation.\n");
3617 if (!b->refine && !b->poly) {
3622 if (b->refine || !b->poly) {
3623 b->regionattrib = 0;
3627 if (b->weighted && (b->poly || b->quality)) {
3630 printf(
"Warning: weighted triangulations (-w, -W) are incompatible\n");
3631 printf(
" with PSLGs (-p) and meshing (-q, -a, -u). Weights ignored.\n"
3635 if (b->jettison && b->nonodewritten && !b->quiet) {
3636 printf(
"Warning: -j and -N switches are somewhat incompatible.\n");
3637 printf(
" If any vertices are jettisoned, you will need the output\n");
3638 printf(
" .node file to reconstruct the new node indices.");
3642 strcpy(b->inpolyfilename, b->innodefilename);
3643 strcpy(b->inelefilename, b->innodefilename);
3644 strcpy(b->areafilename, b->innodefilename);
3646 strcpy(workstring, b->innodefilename);
3648 while (workstring[j] !=
'\0') {
3649 if ((workstring[j] ==
'.') && (workstring[j + 1] !=
'\0')) {
3655 if (increment > 0) {
3658 if ((workstring[j] >=
'0') && (workstring[j] <=
'9')) {
3659 meshnumber = meshnumber * 10 + (int) (workstring[j] -
'0');
3664 }
while (workstring[j] !=
'\0');
3666 if (b->noiterationnum) {
3667 strcpy(b->outnodefilename, b->innodefilename);
3668 strcpy(b->outelefilename, b->innodefilename);
3669 strcpy(b->edgefilename, b->innodefilename);
3670 strcpy(b->vnodefilename, b->innodefilename);
3671 strcpy(b->vedgefilename, b->innodefilename);
3672 strcpy(b->neighborfilename, b->innodefilename);
3673 strcpy(b->offfilename, b->innodefilename);
3674 strcat(b->outnodefilename,
".node");
3675 strcat(b->outelefilename,
".ele");
3676 strcat(b->edgefilename,
".edge");
3677 strcat(b->vnodefilename,
".v.node");
3678 strcat(b->vedgefilename,
".v.edge");
3679 strcat(b->neighborfilename,
".neigh");
3680 strcat(b->offfilename,
".off");
3681 }
else if (increment == 0) {
3682 strcpy(b->outnodefilename, b->innodefilename);
3683 strcpy(b->outpolyfilename, b->innodefilename);
3684 strcpy(b->outelefilename, b->innodefilename);
3685 strcpy(b->edgefilename, b->innodefilename);
3686 strcpy(b->vnodefilename, b->innodefilename);
3687 strcpy(b->vedgefilename, b->innodefilename);
3688 strcpy(b->neighborfilename, b->innodefilename);
3689 strcpy(b->offfilename, b->innodefilename);
3690 strcat(b->outnodefilename,
".1.node");
3691 strcat(b->outpolyfilename,
".1.poly");
3692 strcat(b->outelefilename,
".1.ele");
3693 strcat(b->edgefilename,
".1.edge");
3694 strcat(b->vnodefilename,
".1.v.node");
3695 strcat(b->vedgefilename,
".1.v.edge");
3696 strcat(b->neighborfilename,
".1.neigh");
3697 strcat(b->offfilename,
".1.off");
3699 workstring[increment] =
'%';
3700 workstring[increment + 1] =
'd';
3701 workstring[increment + 2] =
'\0';
3702 sprintf(b->outnodefilename, workstring, meshnumber + 1);
3703 strcpy(b->outpolyfilename, b->outnodefilename);
3704 strcpy(b->outelefilename, b->outnodefilename);
3705 strcpy(b->edgefilename, b->outnodefilename);
3706 strcpy(b->vnodefilename, b->outnodefilename);
3707 strcpy(b->vedgefilename, b->outnodefilename);
3708 strcpy(b->neighborfilename, b->outnodefilename);
3709 strcpy(b->offfilename, b->outnodefilename);
3710 strcat(b->outnodefilename,
".node");
3711 strcat(b->outpolyfilename,
".poly");
3712 strcat(b->outelefilename,
".ele");
3713 strcat(b->edgefilename,
".edge");
3714 strcat(b->vnodefilename,
".v.node");
3715 strcat(b->vedgefilename,
".v.edge");
3716 strcat(b->neighborfilename,
".neigh");
3717 strcat(b->offfilename,
".off");
3719 strcat(b->innodefilename,
".node");
3720 strcat(b->inpolyfilename,
".poly");
3721 strcat(b->inelefilename,
".ele");
3722 strcat(b->areafilename,
".area");
3745#ifdef ANSI_DECLARATORS
3748void printtriangle(m, b, t)
3755 struct otri printtri;
3756 struct osub printsh;
3759 printf(
"triangle x%lx with orientation %d:\n", (
unsigned long) t->tri,
3761 decode(t->tri[0], printtri);
3762 if (printtri.tri == m->dummytri) {
3763 printf(
" [0] = Outer space\n");
3765 printf(
" [0] = x%lx %d\n", (
unsigned long) printtri.tri,
3768 decode(t->tri[1], printtri);
3769 if (printtri.tri == m->dummytri) {
3770 printf(
" [1] = Outer space\n");
3772 printf(
" [1] = x%lx %d\n", (
unsigned long) printtri.tri,
3775 decode(t->tri[2], printtri);
3776 if (printtri.tri == m->dummytri) {
3777 printf(
" [2] = Outer space\n");
3779 printf(
" [2] = x%lx %d\n", (
unsigned long) printtri.tri,
3783 org(*t, printvertex);
3784 if (printvertex == (vertex) NULL)
3785 printf(
" Origin[%d] = NULL\n", (t->orient + 1) % 3 + 3);
3787 printf(
" Origin[%d] = x%lx (%.12g, %.12g)\n",
3788 (t->orient + 1) % 3 + 3, (
unsigned long) printvertex,
3789 printvertex[0], printvertex[1]);
3790 dest(*t, printvertex);
3791 if (printvertex == (vertex) NULL)
3792 printf(
" Dest [%d] = NULL\n", (t->orient + 2) % 3 + 3);
3794 printf(
" Dest [%d] = x%lx (%.12g, %.12g)\n",
3795 (t->orient + 2) % 3 + 3, (
unsigned long) printvertex,
3796 printvertex[0], printvertex[1]);
3797 apex(*t, printvertex);
3798 if (printvertex == (vertex) NULL)
3799 printf(
" Apex [%d] = NULL\n", t->orient + 3);
3801 printf(
" Apex [%d] = x%lx (%.12g, %.12g)\n",
3802 t->orient + 3, (
unsigned long) printvertex,
3803 printvertex[0], printvertex[1]);
3805 if (b->usesegments) {
3806 sdecode(t->tri[6], printsh);
3807 if (printsh.ss != m->dummysub) {
3808 printf(
" [6] = x%lx %d\n", (
unsigned long) printsh.ss,
3811 sdecode(t->tri[7], printsh);
3812 if (printsh.ss != m->dummysub) {
3813 printf(
" [7] = x%lx %d\n", (
unsigned long) printsh.ss,
3816 sdecode(t->tri[8], printsh);
3817 if (printsh.ss != m->dummysub) {
3818 printf(
" [8] = x%lx %d\n", (
unsigned long) printsh.ss,
3824 printf(
" Area constraint: %.4g\n", areabound(*t));
3839#ifdef ANSI_DECLARATORS
3842void printsubseg(m, b, s)
3849 struct osub printsh;
3850 struct otri printtri;
3853 printf(
"subsegment x%lx with orientation %d and mark %d:\n",
3854 (
unsigned long) s->ss, s->ssorient, mark(*s));
3855 sdecode(s->ss[0], printsh);
3856 if (printsh.ss == m->dummysub) {
3857 printf(
" [0] = No subsegment\n");
3859 printf(
" [0] = x%lx %d\n", (
unsigned long) printsh.ss,
3862 sdecode(s->ss[1], printsh);
3863 if (printsh.ss == m->dummysub) {
3864 printf(
" [1] = No subsegment\n");
3866 printf(
" [1] = x%lx %d\n", (
unsigned long) printsh.ss,
3870 sorg(*s, printvertex);
3871 if (printvertex == (vertex) NULL)
3872 printf(
" Origin[%d] = NULL\n", 2 + s->ssorient);
3874 printf(
" Origin[%d] = x%lx (%.12g, %.12g)\n",
3875 2 + s->ssorient, (
unsigned long) printvertex,
3876 printvertex[0], printvertex[1]);
3877 sdest(*s, printvertex);
3878 if (printvertex == (vertex) NULL)
3879 printf(
" Dest [%d] = NULL\n", 3 - s->ssorient);
3881 printf(
" Dest [%d] = x%lx (%.12g, %.12g)\n",
3882 3 - s->ssorient, (
unsigned long) printvertex,
3883 printvertex[0], printvertex[1]);
3885 decode(s->ss[6], printtri);
3886 if (printtri.tri == m->dummytri) {
3887 printf(
" [6] = Outer space\n");
3889 printf(
" [6] = x%lx %d\n", (
unsigned long) printtri.tri,
3892 decode(s->ss[7], printtri);
3893 if (printtri.tri == m->dummytri) {
3894 printf(
" [7] = Outer space\n");
3896 printf(
" [7] = x%lx %d\n", (
unsigned long) printtri.tri,
3900 segorg(*s, printvertex);
3901 if (printvertex == (vertex) NULL)
3902 printf(
" Segment origin[%d] = NULL\n", 4 + s->ssorient);
3904 printf(
" Segment origin[%d] = x%lx (%.12g, %.12g)\n",
3905 4 + s->ssorient, (
unsigned long) printvertex,
3906 printvertex[0], printvertex[1]);
3907 segdest(*s, printvertex);
3908 if (printvertex == (vertex) NULL)
3909 printf(
" Segment dest [%d] = NULL\n", 5 - s->ssorient);
3911 printf(
" Segment dest [%d] = x%lx (%.12g, %.12g)\n",
3912 5 - s->ssorient, (
unsigned long) printvertex,
3913 printvertex[0], printvertex[1]);
3933#ifdef ANSI_DECLARATORS
3941 pool->firstblock = (VOID **) NULL;
3942 pool->nowblock = (VOID **) NULL;
3943 pool->nextitem = (VOID *) NULL;
3944 pool->deaditemstack = (VOID *) NULL;
3945 pool->pathblock = (VOID **) NULL;
3946 pool->pathitem = (VOID *) NULL;
3947 pool->alignbytes = 0;
3948 pool->itembytes = 0;
3949 pool->itemsperblock = 0;
3950 pool->itemsfirstblock = 0;
3953 pool->unallocateditems = 0;
3954 pool->pathitemsleft = 0;
3967#ifdef ANSI_DECLARATORS
3970void poolrestart(pool)
3975 unsigned long alignptr;
3981 pool->nowblock = pool->firstblock;
3983 alignptr = (
unsigned long) (pool->nowblock + 1);
3985 pool->nextitem = (VOID *)
3986 (alignptr + (
unsigned long) pool->alignbytes -
3987 (alignptr % (
unsigned long) pool->alignbytes));
3989 pool->unallocateditems = pool->itemsfirstblock;
3991 pool->deaditemstack = (VOID *) NULL;
4013#ifdef ANSI_DECLARATORS
4014void poolinit(
struct memorypool* pool,
int bytecount,
int itemcount,
4015 int firstitemcount,
unsigned alignment)
4017void poolinit(pool, bytecount, itemcount, firstitemcount, alignment)
4030 if (alignment >
sizeof(VOID *)) {
4031 pool->alignbytes = alignment;
4033 pool->alignbytes =
sizeof(VOID *);
4035 pool->itembytes = ((bytecount - 1) / pool->alignbytes + 1) *
4037 pool->itemsperblock = itemcount;
4038 if (firstitemcount == 0) {
4039 pool->itemsfirstblock = itemcount;
4041 pool->itemsfirstblock = firstitemcount;
4047 pool->firstblock = (VOID **)
4048 trimalloc(pool->itemsfirstblock * pool->itembytes + (
int)
sizeof(VOID *) +
4051 *(pool->firstblock) = (VOID *) NULL;
4061#ifdef ANSI_DECLARATORS
4064void pooldeinit(pool)
4069 while (pool->firstblock != (VOID **) NULL) {
4070 pool->nowblock = (VOID **) *(pool->firstblock);
4071 trifree((VOID *) pool->firstblock);
4072 pool->firstblock = pool->nowblock;
4082#ifdef ANSI_DECLARATORS
4085VOID* poolalloc(pool)
4092 unsigned long alignptr;
4096 if (pool->deaditemstack != (VOID *) NULL) {
4097 newitem = pool->deaditemstack;
4098 pool->deaditemstack = * (VOID **) pool->deaditemstack;
4101 if (pool->unallocateditems == 0) {
4103 if (*(pool->nowblock) == (VOID *) NULL) {
4105 newblock = (VOID **) trimalloc(pool->itemsperblock * pool->itembytes +
4106 (
int)
sizeof(VOID *) +
4108 *(pool->nowblock) = (VOID *) newblock;
4110 *newblock = (VOID *) NULL;
4114 pool->nowblock = (VOID **) *(pool->nowblock);
4117 alignptr = (
unsigned long) (pool->nowblock + 1);
4119 pool->nextitem = (VOID *)
4120 (alignptr + (
unsigned long) pool->alignbytes -
4121 (alignptr % (
unsigned long) pool->alignbytes));
4123 pool->unallocateditems = pool->itemsperblock;
4127 newitem = pool->nextitem;
4129 pool->nextitem = (VOID *) ((
char *) pool->nextitem + pool->itembytes);
4130 pool->unallocateditems--;
4145#ifdef ANSI_DECLARATORS
4146void pooldealloc(
struct memorypool* pool, VOID *dyingitem)
4148void pooldealloc(pool, dyingitem)
4155 *((VOID **) dyingitem) = pool->deaditemstack;
4156 pool->deaditemstack = dyingitem;
4168#ifdef ANSI_DECLARATORS
4171void traversalinit(pool)
4176 unsigned long alignptr;
4179 pool->pathblock = pool->firstblock;
4181 alignptr = (
unsigned long) (pool->pathblock + 1);
4183 pool->pathitem = (VOID *)
4184 (alignptr + (
unsigned long) pool->alignbytes -
4185 (alignptr % (
unsigned long) pool->alignbytes));
4187 pool->pathitemsleft = pool->itemsfirstblock;
4204#ifdef ANSI_DECLARATORS
4213 unsigned long alignptr;
4216 if (pool->pathitem == pool->nextitem) {
4217 return (VOID *) NULL;
4221 if (pool->pathitemsleft == 0) {
4223 pool->pathblock = (VOID **) *(pool->pathblock);
4225 alignptr = (
unsigned long) (pool->pathblock + 1);
4227 pool->pathitem = (VOID *)
4228 (alignptr + (
unsigned long) pool->alignbytes -
4229 (alignptr % (
unsigned long) pool->alignbytes));
4231 pool->pathitemsleft = pool->itemsperblock;
4234 newitem = pool->pathitem;
4236 pool->pathitem = (VOID *) ((
char *) pool->pathitem + pool->itembytes);
4237 pool->pathitemsleft--;
4269#ifdef ANSI_DECLARATORS
4270void dummyinit(
struct mesh *m,
struct behavior *b,
int trianglebytes,
4273void dummyinit(m, b, trianglebytes, subsegbytes)
4281 unsigned long alignptr;
4284 m->dummytribase = (triangle *) trimalloc(trianglebytes +
4285 m->triangles.alignbytes);
4287 alignptr = (
unsigned long) m->dummytribase;
4288 m->dummytri = (triangle *)
4289 (alignptr + (
unsigned long) m->triangles.alignbytes -
4290 (alignptr % (
unsigned long) m->triangles.alignbytes));
4295 m->dummytri[0] = (triangle) m->dummytri;
4296 m->dummytri[1] = (triangle) m->dummytri;
4297 m->dummytri[2] = (triangle) m->dummytri;
4299 m->dummytri[3] = (triangle) NULL;
4300 m->dummytri[4] = (triangle) NULL;
4301 m->dummytri[5] = (triangle) NULL;
4303 if (b->usesegments) {
4307 m->dummysubbase = (subseg *) trimalloc(subsegbytes +
4308 m->subsegs.alignbytes);
4310 alignptr = (
unsigned long) m->dummysubbase;
4311 m->dummysub = (subseg *)
4312 (alignptr + (
unsigned long) m->subsegs.alignbytes -
4313 (alignptr % (
unsigned long) m->subsegs.alignbytes));
4318 m->dummysub[0] = (subseg) m->dummysub;
4319 m->dummysub[1] = (subseg) m->dummysub;
4321 m->dummysub[2] = (subseg) NULL;
4322 m->dummysub[3] = (subseg) NULL;
4323 m->dummysub[4] = (subseg) NULL;
4324 m->dummysub[5] = (subseg) NULL;
4326 m->dummysub[6] = (subseg) m->dummytri;
4327 m->dummysub[7] = (subseg) m->dummytri;
4329 * (
int *) (m->dummysub + 8) = 0;
4333 m->dummytri[6] = (triangle) m->dummysub;
4334 m->dummytri[7] = (triangle) m->dummysub;
4335 m->dummytri[8] = (triangle) m->dummysub;
4349#ifdef ANSI_DECLARATORS
4350void initializevertexpool(
struct mesh *m,
struct behavior *b)
4352void initializevertexpool(m, b)
4363 m->vertexmarkindex = ((m->mesh_dim + m->nextras) *
sizeof(REAL) +
4366 vertexsize = (m->vertexmarkindex + 2) *
sizeof(
int);
4370 m->vertex2triindex = (vertexsize +
sizeof(triangle) - 1) /
4372 vertexsize = (m->vertex2triindex + 1) *
sizeof(triangle);
4376 poolinit(&m->vertices, vertexsize, VERTEXPERBLOCK,
4377 m->invertices > VERTEXPERBLOCK ? m->invertices : VERTEXPERBLOCK,
4392#ifdef ANSI_DECLARATORS
4393void initializetrisubpools(
struct mesh *m,
struct behavior *b)
4395void initializetrisubpools(m, b)
4407 m->highorderindex = 6 + (b->usesegments * 3);
4409 trisize = ((b->order + 1) * (b->order + 2) / 2 + (m->highorderindex - 3)) *
4413 m->elemattribindex = (trisize +
sizeof(REAL) - 1) /
sizeof(REAL);
4417 m->areaboundindex = m->elemattribindex + m->eextras + b->regionattrib;
4421 trisize = (m->areaboundindex + 1) *
sizeof(REAL);
4422 }
else if (m->eextras + b->regionattrib > 0) {
4423 trisize = m->areaboundindex *
sizeof(REAL);
4429 if ((b->voronoi || b->neighbors) &&
4430 (trisize < 6 *
sizeof(triangle) +
sizeof(
int))) {
4431 trisize = 6 *
sizeof(triangle) +
sizeof(
int);
4435 poolinit(&m->triangles, trisize, TRIPERBLOCK,
4436 (2 * m->invertices - 2) > TRIPERBLOCK ? (2 * m->invertices - 2) :
4439 if (b->usesegments) {
4442 poolinit(&m->subsegs, 8 *
sizeof(triangle) +
sizeof(
int),
4443 SUBSEGPERBLOCK, SUBSEGPERBLOCK, 4);
4446 dummyinit(m, b, m->triangles.itembytes, m->subsegs.itembytes);
4449 dummyinit(m, b, m->triangles.itembytes, 0);
4459#ifdef ANSI_DECLARATORS
4460void triangledealloc(
struct mesh *m, triangle *dyingtriangle)
4462void triangledealloc(m, dyingtriangle)
4464triangle *dyingtriangle;
4470 killtri(dyingtriangle);
4471 pooldealloc(&m->triangles, (VOID *) dyingtriangle);
4480#ifdef ANSI_DECLARATORS
4481triangle *triangletraverse(
struct mesh *m)
4483triangle *triangletraverse(m)
4488 triangle *newtriangle;
4491 newtriangle = (triangle *) traverse(&m->triangles);
4492 if (newtriangle == (triangle *) NULL) {
4493 return (triangle *) NULL;
4495 }
while (deadtri(newtriangle));
4505#ifdef ANSI_DECLARATORS
4506void subsegdealloc(
struct mesh *m, subseg *dyingsubseg)
4508void subsegdealloc(m, dyingsubseg)
4516 killsubseg(dyingsubseg);
4517 pooldealloc(&m->subsegs, (VOID *) dyingsubseg);
4526#ifdef ANSI_DECLARATORS
4527subseg *subsegtraverse(
struct mesh *m)
4529subseg *subsegtraverse(m)
4537 newsubseg = (subseg *) traverse(&m->subsegs);
4538 if (newsubseg == (subseg *) NULL) {
4539 return (subseg *) NULL;
4541 }
while (deadsubseg(newsubseg));
4551#ifdef ANSI_DECLARATORS
4552void vertexdealloc(
struct mesh *m, vertex dyingvertex)
4554void vertexdealloc(m, dyingvertex)
4562 setvertextype(dyingvertex, DEADVERTEX);
4563 pooldealloc(&m->vertices, (VOID *) dyingvertex);
4572#ifdef ANSI_DECLARATORS
4573vertex vertextraverse(
struct mesh *m)
4575vertex vertextraverse(m)
4583 newvertex = (vertex) traverse(&m->vertices);
4584 if (newvertex == (vertex) NULL) {
4585 return (vertex) NULL;
4587 }
while (vertextype(newvertex) == DEADVERTEX);
4600#ifdef ANSI_DECLARATORS
4601void badsubsegdealloc(
struct mesh *m,
struct badsubseg *dyingseg)
4603void badsubsegdealloc(m, dyingseg)
4611 dyingseg->subsegorg = (vertex) NULL;
4612 pooldealloc(&m->badsubsegs, (VOID *) dyingseg);
4625#ifdef ANSI_DECLARATORS
4636 newseg = (
struct badsubseg *) traverse(&m->badsubsegs);
4637 if (newseg == (
struct badsubseg *) NULL) {
4640 }
while (newseg->subsegorg == (vertex) NULL);
4658#ifdef ANSI_DECLARATORS
4659vertex getvertex(
struct mesh *m,
struct behavior *b,
int number)
4661vertex getvertex(m, b, number)
4670 unsigned long alignptr;
4673 getblock = m->vertices.firstblock;
4674 current = b->firstnumber;
4677 if (current + m->vertices.itemsfirstblock <= number) {
4678 getblock = (VOID **) *getblock;
4679 current += m->vertices.itemsfirstblock;
4680 while (current + m->vertices.itemsperblock <= number) {
4681 getblock = (VOID **) *getblock;
4682 current += m->vertices.itemsperblock;
4687 alignptr = (
unsigned long) (getblock + 1);
4688 foundvertex = (
char *) (alignptr + (
unsigned long) m->vertices.alignbytes -
4689 (alignptr % (
unsigned long) m->vertices.alignbytes));
4690 return (vertex) (foundvertex + m->vertices.itembytes * (number - current));
4699#ifdef ANSI_DECLARATORS
4702void triangledeinit(m, b)
4708 pooldeinit(&m->triangles);
4709 trifree((VOID *) m->dummytribase);
4710 if (b->usesegments) {
4711 pooldeinit(&m->subsegs);
4712 trifree((VOID *) m->dummysubbase);
4714 pooldeinit(&m->vertices);
4717 pooldeinit(&m->badsubsegs);
4718 if ((b->minangle > 0.0) || b->vararea || b->fixedarea || b->usertest) {
4719 pooldeinit(&m->badtriangles);
4720 pooldeinit(&m->flipstackers);
4740#ifdef ANSI_DECLARATORS
4743void maketriangle(m, b, newotri)
4746struct otri *newotri;
4752 newotri->tri = (triangle *) poolalloc(&m->triangles);
4754 newotri->tri[0] = (triangle) m->dummytri;
4755 newotri->tri[1] = (triangle) m->dummytri;
4756 newotri->tri[2] = (triangle) m->dummytri;
4758 newotri->tri[3] = (triangle) NULL;
4759 newotri->tri[4] = (triangle) NULL;
4760 newotri->tri[5] = (triangle) NULL;
4761 if (b->usesegments) {
4764 newotri->tri[6] = (triangle) m->dummysub;
4765 newotri->tri[7] = (triangle) m->dummysub;
4766 newotri->tri[8] = (triangle) m->dummysub;
4768 for (i = 0; i < m->eextras; i++) {
4769 setelemattribute(*newotri, i, 0.0);
4772 setareabound(*newotri, -1.0);
4775 newotri->orient = 0;
4784#ifdef ANSI_DECLARATORS
4785void makesubseg(
struct mesh *m,
struct osub *newsubseg)
4787void makesubseg(m, newsubseg)
4789struct osub *newsubseg;
4793 newsubseg->ss = (subseg *) poolalloc(&m->subsegs);
4796 newsubseg->ss[0] = (subseg) m->dummysub;
4797 newsubseg->ss[1] = (subseg) m->dummysub;
4799 newsubseg->ss[2] = (subseg) NULL;
4800 newsubseg->ss[3] = (subseg) NULL;
4801 newsubseg->ss[4] = (subseg) NULL;
4802 newsubseg->ss[5] = (subseg) NULL;
4804 newsubseg->ss[6] = (subseg) m->dummytri;
4805 newsubseg->ss[7] = (subseg) m->dummytri;
4807 setmark(*newsubseg, 0);
4809 newsubseg->ssorient = 0;
4833#define Absolute(a) ((a) >= 0.0 ? (a) : -(a))
4849#define Fast_Two_Sum_Tail(a, b, x, y) \
4853#define Fast_Two_Sum(a, b, x, y) \
4854 x = (REAL) (a + b); \
4855 Fast_Two_Sum_Tail(a, b, x, y)
4857#define Two_Sum_Tail(a, b, x, y) \
4858 bvirt = (REAL) (x - a); \
4859 avirt = x - bvirt; \
4860 bround = b - bvirt; \
4861 around = a - avirt; \
4864#define Two_Sum(a, b, x, y) \
4865 x = (REAL) (a + b); \
4866 Two_Sum_Tail(a, b, x, y)
4868#define Two_Diff_Tail(a, b, x, y) \
4869 bvirt = (REAL) (a - x); \
4870 avirt = x + bvirt; \
4871 bround = bvirt - b; \
4872 around = a - avirt; \
4875#define Two_Diff(a, b, x, y) \
4876 x = (REAL) (a - b); \
4877 Two_Diff_Tail(a, b, x, y)
4879#define Split(a, ahi, alo) \
4880 c = (REAL) (splitter * a); \
4881 abig = (REAL) (c - a); \
4885#define Two_Product_Tail(a, b, x, y) \
4886 Split(a, ahi, alo); \
4887 Split(b, bhi, blo); \
4888 err1 = x - (ahi * bhi); \
4889 err2 = err1 - (alo * bhi); \
4890 err3 = err2 - (ahi * blo); \
4891 y = (alo * blo) - err3
4893#define Two_Product(a, b, x, y) \
4894 x = (REAL) (a * b); \
4895 Two_Product_Tail(a, b, x, y)
4900#define Two_Product_Presplit(a, b, bhi, blo, x, y) \
4901 x = (REAL) (a * b); \
4902 Split(a, ahi, alo); \
4903 err1 = x - (ahi * bhi); \
4904 err2 = err1 - (alo * bhi); \
4905 err3 = err2 - (ahi * blo); \
4906 y = (alo * blo) - err3
4910#define Square_Tail(a, x, y) \
4911 Split(a, ahi, alo); \
4912 err1 = x - (ahi * ahi); \
4913 err3 = err1 - ((ahi + ahi) * alo); \
4914 y = (alo * alo) - err3
4916#define Square(a, x, y) \
4917 x = (REAL) (a * a); \
4918 Square_Tail(a, x, y)
4923#define Two_One_Sum(a1, a0, b, x2, x1, x0) \
4924 Two_Sum(a0, b , _i, x0); \
4925 Two_Sum(a1, _i, x2, x1)
4927#define Two_One_Diff(a1, a0, b, x2, x1, x0) \
4928 Two_Diff(a0, b , _i, x0); \
4929 Two_Sum( a1, _i, x2, x1)
4931#define Two_Two_Sum(a1, a0, b1, b0, x3, x2, x1, x0) \
4932 Two_One_Sum(a1, a0, b0, _j, _0, x0); \
4933 Two_One_Sum(_j, _0, b1, x3, x2, x1)
4935#define Two_Two_Diff(a1, a0, b1, b0, x3, x2, x1, x0) \
4936 Two_One_Diff(a1, a0, b0, _j, _0, x0); \
4937 Two_One_Diff(_j, _0, b1, x3, x2, x1)
4941#define Two_One_Product(a1, a0, b, x3, x2, x1, x0) \
4942 Split(b, bhi, blo); \
4943 Two_Product_Presplit(a0, b, bhi, blo, _i, x0); \
4944 Two_Product_Presplit(a1, b, bhi, blo, _j, _0); \
4945 Two_Sum(_i, _0, _k, x1); \
4946 Fast_Two_Sum(_j, _k, x3, x2)
4970 REAL check, lastcheck;
4978 _control87(_PC_24, _MCW_PC);
4980 _control87(_PC_53, _MCW_PC);
5009 every_other = !every_other;
5010 check = 1.0 + epsilon;
5011 }
while ((check != 1.0) && (check != lastcheck));
5014 resulterrbound = (3.0 + 8.0 * epsilon) * epsilon;
5015 ccwerrboundA = (3.0 + 16.0 * epsilon) * epsilon;
5016 ccwerrboundB = (2.0 + 12.0 * epsilon) * epsilon;
5017 ccwerrboundC = (9.0 + 64.0 * epsilon) * epsilon * epsilon;
5018 iccerrboundA = (10.0 + 96.0 * epsilon) * epsilon;
5019 iccerrboundB = (4.0 + 48.0 * epsilon) * epsilon;
5020 iccerrboundC = (44.0 + 576.0 * epsilon) * epsilon * epsilon;
5021 o3derrboundA = (7.0 + 56.0 * epsilon) * epsilon;
5022 o3derrboundB = (3.0 + 28.0 * epsilon) * epsilon;
5023 o3derrboundC = (26.0 + 288.0 * epsilon) * epsilon * epsilon;
5040#ifdef ANSI_DECLARATORS
5041int fast_expansion_sum_zeroelim(
int elen, REAL *e,
int flen, REAL *f, REAL *h)
5043int fast_expansion_sum_zeroelim(elen, e, flen, f, h)
5056 REAL avirt, bround, around;
5057 int eindex, findex, hindex;
5062 eindex = findex = 0;
5063 if ((fnow > enow) == (fnow > -enow)) {
5071 if ((eindex < elen) && (findex < flen)) {
5072 if ((fnow > enow) == (fnow > -enow)) {
5073 Fast_Two_Sum(enow, Q, Qnew, hh);
5076 Fast_Two_Sum(fnow, Q, Qnew, hh);
5083 while ((eindex < elen) && (findex < flen)) {
5084 if ((fnow > enow) == (fnow > -enow)) {
5085 Two_Sum(Q, enow, Qnew, hh);
5088 Two_Sum(Q, fnow, Qnew, hh);
5097 while (eindex < elen) {
5098 Two_Sum(Q, enow, Qnew, hh);
5105 while (findex < flen) {
5106 Two_Sum(Q, fnow, Qnew, hh);
5113 if ((Q != 0.0) || (hindex == 0)) {
5134#ifdef ANSI_DECLARATORS
5135int scale_expansion_zeroelim(
int elen, REAL *e, REAL b, REAL *h)
5137int scale_expansion_zeroelim(elen, e, b, h)
5145 INEXACT REAL Q, sum;
5147 INEXACT REAL product1;
5152 REAL avirt, bround, around;
5155 REAL ahi, alo, bhi, blo;
5156 REAL err1, err2, err3;
5159 Two_Product_Presplit(e[0], b, bhi, blo, Q, hh);
5164 for (eindex = 1; eindex < elen; eindex++) {
5166 Two_Product_Presplit(enow, b, bhi, blo, product1, product0);
5167 Two_Sum(Q, product0, sum, hh);
5171 Fast_Two_Sum(product1, sum, Q, hh);
5176 if ((Q != 0.0) || (hindex == 0)) {
5190#ifdef ANSI_DECLARATORS
5191REAL estimate(
int elen, REAL *e)
5193REAL estimate(elen, e)
5203 for (eindex = 1; eindex < elen; eindex++) {
5229#ifdef ANSI_DECLARATORS
5230REAL counterclockwiseadapt(vertex pa, vertex pb, vertex pc, REAL detsum)
5232REAL counterclockwiseadapt(pa, pb, pc, detsum)
5240 INEXACT REAL acx, acy, bcx, bcy;
5241 REAL acxtail, acytail, bcxtail, bcytail;
5242 INEXACT REAL detleft, detright;
5243 REAL detlefttail, detrighttail;
5245 REAL B[4], C1[8], C2[12], D[16];
5247 int C1length, C2length, Dlength;
5250 INEXACT REAL s1, t1;
5254 REAL avirt, bround, around;
5257 REAL ahi, alo, bhi, blo;
5258 REAL err1, err2, err3;
5259 INEXACT REAL _i, _j;
5262 acx = (REAL) (pa[0] - pc[0]);
5263 bcx = (REAL) (pb[0] - pc[0]);
5264 acy = (REAL) (pa[1] - pc[1]);
5265 bcy = (REAL) (pb[1] - pc[1]);
5267 Two_Product(acx, bcy, detleft, detlefttail);
5268 Two_Product(acy, bcx, detright, detrighttail);
5270 Two_Two_Diff(detleft, detlefttail, detright, detrighttail,
5271 B3, B[2], B[1], B[0]);
5274 det = estimate(4, B);
5275 errbound = ccwerrboundB * detsum;
5276 if ((det >= errbound) || (-det >= errbound)) {
5280 Two_Diff_Tail(pa[0], pc[0], acx, acxtail);
5281 Two_Diff_Tail(pb[0], pc[0], bcx, bcxtail);
5282 Two_Diff_Tail(pa[1], pc[1], acy, acytail);
5283 Two_Diff_Tail(pb[1], pc[1], bcy, bcytail);
5285 if ((acxtail == 0.0) && (acytail == 0.0)
5286 && (bcxtail == 0.0) && (bcytail == 0.0)) {
5290 errbound = ccwerrboundC * detsum + resulterrbound * Absolute(det);
5291 det += (acx * bcytail + bcy * acxtail)
5292 - (acy * bcxtail + bcx * acytail);
5293 if ((det >= errbound) || (-det >= errbound)) {
5297 Two_Product(acxtail, bcy, s1, s0);
5298 Two_Product(acytail, bcx, t1, t0);
5299 Two_Two_Diff(s1, s0, t1, t0, u3, u[2], u[1], u[0]);
5301 C1length = fast_expansion_sum_zeroelim(4, B, 4, u, C1);
5303 Two_Product(acx, bcytail, s1, s0);
5304 Two_Product(acy, bcxtail, t1, t0);
5305 Two_Two_Diff(s1, s0, t1, t0, u3, u[2], u[1], u[0]);
5307 C2length = fast_expansion_sum_zeroelim(C1length, C1, 4, u, C2);
5309 Two_Product(acxtail, bcytail, s1, s0);
5310 Two_Product(acytail, bcxtail, t1, t0);
5311 Two_Two_Diff(s1, s0, t1, t0, u3, u[2], u[1], u[0]);
5313 Dlength = fast_expansion_sum_zeroelim(C2length, C2, 4, u, D);
5315 return(D[Dlength - 1]);
5318#ifdef ANSI_DECLARATORS
5319REAL counterclockwise(
struct mesh *m,
struct behavior *b,
5320 vertex pa, vertex pb, vertex pc)
5322REAL counterclockwise(m, b, pa, pb, pc)
5331 REAL detleft, detright, det;
5332 REAL detsum, errbound;
5334 m->counterclockcount++;
5336 detleft = (pa[0] - pc[0]) * (pb[1] - pc[1]);
5337 detright = (pa[1] - pc[1]) * (pb[0] - pc[0]);
5338 det = detleft - detright;
5344 if (detleft > 0.0) {
5345 if (detright <= 0.0) {
5348 detsum = detleft + detright;
5350 }
else if (detleft < 0.0) {
5351 if (detright >= 0.0) {
5354 detsum = -detleft - detright;
5360 errbound = ccwerrboundA * detsum;
5361 if ((det >= errbound) || (-det >= errbound)) {
5365 return counterclockwiseadapt(pa, pb, pc, detsum);
5387#ifdef ANSI_DECLARATORS
5388REAL incircleadapt(vertex pa, vertex pb, vertex pc, vertex pd, REAL permanent)
5390REAL incircleadapt(pa, pb, pc, pd, permanent)
5399 INEXACT REAL adx, bdx, cdx, ady, bdy, cdy;
5402 INEXACT REAL bdxcdy1, cdxbdy1, cdxady1, adxcdy1, adxbdy1, bdxady1;
5403 REAL bdxcdy0, cdxbdy0, cdxady0, adxcdy0, adxbdy0, bdxady0;
5404 REAL bc[4], ca[4], ab[4];
5405 INEXACT REAL bc3, ca3, ab3;
5406 REAL axbc[8], axxbc[16], aybc[8], ayybc[16], adet[32];
5407 int axbclen, axxbclen, aybclen, ayybclen, alen;
5408 REAL bxca[8], bxxca[16], byca[8], byyca[16], bdet[32];
5409 int bxcalen, bxxcalen, bycalen, byycalen, blen;
5410 REAL cxab[8], cxxab[16], cyab[8], cyyab[16], cdet[32];
5411 int cxablen, cxxablen, cyablen, cyyablen, clen;
5414 REAL fin1[1152], fin2[1152];
5415 REAL *finnow, *finother, *finswap;
5418 REAL adxtail, bdxtail, cdxtail, adytail, bdytail, cdytail;
5419 INEXACT REAL adxadx1, adyady1, bdxbdx1, bdybdy1, cdxcdx1, cdycdy1;
5420 REAL adxadx0, adyady0, bdxbdx0, bdybdy0, cdxcdx0, cdycdy0;
5421 REAL aa[4], bb[4], cc[4];
5422 INEXACT REAL aa3, bb3, cc3;
5423 INEXACT REAL ti1, tj1;
5426 INEXACT REAL u3, v3;
5427 REAL temp8[8], temp16a[16], temp16b[16], temp16c[16];
5428 REAL temp32a[32], temp32b[32], temp48[48], temp64[64];
5429 int temp8len, temp16alen, temp16blen, temp16clen;
5430 int temp32alen, temp32blen, temp48len, temp64len;
5431 REAL axtbb[8], axtcc[8], aytbb[8], aytcc[8];
5432 int axtbblen, axtcclen, aytbblen, aytcclen;
5433 REAL bxtaa[8], bxtcc[8], bytaa[8], bytcc[8];
5434 int bxtaalen, bxtcclen, bytaalen, bytcclen;
5435 REAL cxtaa[8], cxtbb[8], cytaa[8], cytbb[8];
5436 int cxtaalen, cxtbblen, cytaalen, cytbblen;
5437 REAL axtbc[8], aytbc[8], bxtca[8], bytca[8], cxtab[8], cytab[8];
5438 int axtbclen=0, aytbclen=0, bxtcalen=0, bytcalen=0, cxtablen=0, cytablen=0;
5439 REAL axtbct[16], aytbct[16], bxtcat[16], bytcat[16], cxtabt[16], cytabt[16];
5440 int axtbctlen, aytbctlen, bxtcatlen, bytcatlen, cxtabtlen, cytabtlen;
5441 REAL axtbctt[8], aytbctt[8], bxtcatt[8];
5442 REAL bytcatt[8], cxtabtt[8], cytabtt[8];
5443 int axtbcttlen, aytbcttlen, bxtcattlen, bytcattlen, cxtabttlen, cytabttlen;
5444 REAL abt[8], bct[8], cat[8];
5445 int abtlen, bctlen, catlen;
5446 REAL abtt[4], bctt[4], catt[4];
5447 int abttlen, bcttlen, cattlen;
5448 INEXACT REAL abtt3, bctt3, catt3;
5452 REAL avirt, bround, around;
5455 REAL ahi, alo, bhi, blo;
5456 REAL err1, err2, err3;
5457 INEXACT REAL _i, _j;
5460 adx = (REAL) (pa[0] - pd[0]);
5461 bdx = (REAL) (pb[0] - pd[0]);
5462 cdx = (REAL) (pc[0] - pd[0]);
5463 ady = (REAL) (pa[1] - pd[1]);
5464 bdy = (REAL) (pb[1] - pd[1]);
5465 cdy = (REAL) (pc[1] - pd[1]);
5467 Two_Product(bdx, cdy, bdxcdy1, bdxcdy0);
5468 Two_Product(cdx, bdy, cdxbdy1, cdxbdy0);
5469 Two_Two_Diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0, bc3, bc[2], bc[1], bc[0]);
5471 axbclen = scale_expansion_zeroelim(4, bc, adx, axbc);
5472 axxbclen = scale_expansion_zeroelim(axbclen, axbc, adx, axxbc);
5473 aybclen = scale_expansion_zeroelim(4, bc, ady, aybc);
5474 ayybclen = scale_expansion_zeroelim(aybclen, aybc, ady, ayybc);
5475 alen = fast_expansion_sum_zeroelim(axxbclen, axxbc, ayybclen, ayybc, adet);
5477 Two_Product(cdx, ady, cdxady1, cdxady0);
5478 Two_Product(adx, cdy, adxcdy1, adxcdy0);
5479 Two_Two_Diff(cdxady1, cdxady0, adxcdy1, adxcdy0, ca3, ca[2], ca[1], ca[0]);
5481 bxcalen = scale_expansion_zeroelim(4, ca, bdx, bxca);
5482 bxxcalen = scale_expansion_zeroelim(bxcalen, bxca, bdx, bxxca);
5483 bycalen = scale_expansion_zeroelim(4, ca, bdy, byca);
5484 byycalen = scale_expansion_zeroelim(bycalen, byca, bdy, byyca);
5485 blen = fast_expansion_sum_zeroelim(bxxcalen, bxxca, byycalen, byyca, bdet);
5487 Two_Product(adx, bdy, adxbdy1, adxbdy0);
5488 Two_Product(bdx, ady, bdxady1, bdxady0);
5489 Two_Two_Diff(adxbdy1, adxbdy0, bdxady1, bdxady0, ab3, ab[2], ab[1], ab[0]);
5491 cxablen = scale_expansion_zeroelim(4, ab, cdx, cxab);
5492 cxxablen = scale_expansion_zeroelim(cxablen, cxab, cdx, cxxab);
5493 cyablen = scale_expansion_zeroelim(4, ab, cdy, cyab);
5494 cyyablen = scale_expansion_zeroelim(cyablen, cyab, cdy, cyyab);
5495 clen = fast_expansion_sum_zeroelim(cxxablen, cxxab, cyyablen, cyyab, cdet);
5497 ablen = fast_expansion_sum_zeroelim(alen, adet, blen, bdet, abdet);
5498 finlength = fast_expansion_sum_zeroelim(ablen, abdet, clen, cdet, fin1);
5500 det = estimate(finlength, fin1);
5501 errbound = iccerrboundB * permanent;
5502 if ((det >= errbound) || (-det >= errbound)) {
5506 Two_Diff_Tail(pa[0], pd[0], adx, adxtail);
5507 Two_Diff_Tail(pa[1], pd[1], ady, adytail);
5508 Two_Diff_Tail(pb[0], pd[0], bdx, bdxtail);
5509 Two_Diff_Tail(pb[1], pd[1], bdy, bdytail);
5510 Two_Diff_Tail(pc[0], pd[0], cdx, cdxtail);
5511 Two_Diff_Tail(pc[1], pd[1], cdy, cdytail);
5512 if ((adxtail == 0.0) && (bdxtail == 0.0) && (cdxtail == 0.0)
5513 && (adytail == 0.0) && (bdytail == 0.0) && (cdytail == 0.0)) {
5517 errbound = iccerrboundC * permanent + resulterrbound * Absolute(det);
5518 det += ((adx * adx + ady * ady) * ((bdx * cdytail + cdy * bdxtail)
5519 - (bdy * cdxtail + cdx * bdytail))
5520 + 2.0 * (adx * adxtail + ady * adytail) * (bdx * cdy - bdy * cdx))
5521 + ((bdx * bdx + bdy * bdy) * ((cdx * adytail + ady * cdxtail)
5522 - (cdy * adxtail + adx * cdytail))
5523 + 2.0 * (bdx * bdxtail + bdy * bdytail) * (cdx * ady - cdy * adx))
5524 + ((cdx * cdx + cdy * cdy) * ((adx * bdytail + bdy * adxtail)
5525 - (ady * bdxtail + bdx * adytail))
5526 + 2.0 * (cdx * cdxtail + cdy * cdytail) * (adx * bdy - ady * bdx));
5527 if ((det >= errbound) || (-det >= errbound)) {
5534 if ((bdxtail != 0.0) || (bdytail != 0.0)
5535 || (cdxtail != 0.0) || (cdytail != 0.0)) {
5536 Square(adx, adxadx1, adxadx0);
5537 Square(ady, adyady1, adyady0);
5538 Two_Two_Sum(adxadx1, adxadx0, adyady1, adyady0, aa3, aa[2], aa[1], aa[0]);
5541 if ((cdxtail != 0.0) || (cdytail != 0.0)
5542 || (adxtail != 0.0) || (adytail != 0.0)) {
5543 Square(bdx, bdxbdx1, bdxbdx0);
5544 Square(bdy, bdybdy1, bdybdy0);
5545 Two_Two_Sum(bdxbdx1, bdxbdx0, bdybdy1, bdybdy0, bb3, bb[2], bb[1], bb[0]);
5548 if ((adxtail != 0.0) || (adytail != 0.0)
5549 || (bdxtail != 0.0) || (bdytail != 0.0)) {
5550 Square(cdx, cdxcdx1, cdxcdx0);
5551 Square(cdy, cdycdy1, cdycdy0);
5552 Two_Two_Sum(cdxcdx1, cdxcdx0, cdycdy1, cdycdy0, cc3, cc[2], cc[1], cc[0]);
5556 if (adxtail != 0.0) {
5557 axtbclen = scale_expansion_zeroelim(4, bc, adxtail, axtbc);
5558 temp16alen = scale_expansion_zeroelim(axtbclen, axtbc, 2.0 * adx,
5561 axtcclen = scale_expansion_zeroelim(4, cc, adxtail, axtcc);
5562 temp16blen = scale_expansion_zeroelim(axtcclen, axtcc, bdy, temp16b);
5564 axtbblen = scale_expansion_zeroelim(4, bb, adxtail, axtbb);
5565 temp16clen = scale_expansion_zeroelim(axtbblen, axtbb, -cdy, temp16c);
5567 temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5568 temp16blen, temp16b, temp32a);
5569 temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
5570 temp32alen, temp32a, temp48);
5571 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5573 finswap = finnow; finnow = finother; finother = finswap;
5575 if (adytail != 0.0) {
5576 aytbclen = scale_expansion_zeroelim(4, bc, adytail, aytbc);
5577 temp16alen = scale_expansion_zeroelim(aytbclen, aytbc, 2.0 * ady,
5580 aytbblen = scale_expansion_zeroelim(4, bb, adytail, aytbb);
5581 temp16blen = scale_expansion_zeroelim(aytbblen, aytbb, cdx, temp16b);
5583 aytcclen = scale_expansion_zeroelim(4, cc, adytail, aytcc);
5584 temp16clen = scale_expansion_zeroelim(aytcclen, aytcc, -bdx, temp16c);
5586 temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5587 temp16blen, temp16b, temp32a);
5588 temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
5589 temp32alen, temp32a, temp48);
5590 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5592 finswap = finnow; finnow = finother; finother = finswap;
5594 if (bdxtail != 0.0) {
5595 bxtcalen = scale_expansion_zeroelim(4, ca, bdxtail, bxtca);
5596 temp16alen = scale_expansion_zeroelim(bxtcalen, bxtca, 2.0 * bdx,
5599 bxtaalen = scale_expansion_zeroelim(4, aa, bdxtail, bxtaa);
5600 temp16blen = scale_expansion_zeroelim(bxtaalen, bxtaa, cdy, temp16b);
5602 bxtcclen = scale_expansion_zeroelim(4, cc, bdxtail, bxtcc);
5603 temp16clen = scale_expansion_zeroelim(bxtcclen, bxtcc, -ady, temp16c);
5605 temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5606 temp16blen, temp16b, temp32a);
5607 temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
5608 temp32alen, temp32a, temp48);
5609 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5611 finswap = finnow; finnow = finother; finother = finswap;
5613 if (bdytail != 0.0) {
5614 bytcalen = scale_expansion_zeroelim(4, ca, bdytail, bytca);
5615 temp16alen = scale_expansion_zeroelim(bytcalen, bytca, 2.0 * bdy,
5618 bytcclen = scale_expansion_zeroelim(4, cc, bdytail, bytcc);
5619 temp16blen = scale_expansion_zeroelim(bytcclen, bytcc, adx, temp16b);
5621 bytaalen = scale_expansion_zeroelim(4, aa, bdytail, bytaa);
5622 temp16clen = scale_expansion_zeroelim(bytaalen, bytaa, -cdx, temp16c);
5624 temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5625 temp16blen, temp16b, temp32a);
5626 temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
5627 temp32alen, temp32a, temp48);
5628 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5630 finswap = finnow; finnow = finother; finother = finswap;
5632 if (cdxtail != 0.0) {
5633 cxtablen = scale_expansion_zeroelim(4, ab, cdxtail, cxtab);
5634 temp16alen = scale_expansion_zeroelim(cxtablen, cxtab, 2.0 * cdx,
5637 cxtbblen = scale_expansion_zeroelim(4, bb, cdxtail, cxtbb);
5638 temp16blen = scale_expansion_zeroelim(cxtbblen, cxtbb, ady, temp16b);
5640 cxtaalen = scale_expansion_zeroelim(4, aa, cdxtail, cxtaa);
5641 temp16clen = scale_expansion_zeroelim(cxtaalen, cxtaa, -bdy, temp16c);
5643 temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5644 temp16blen, temp16b, temp32a);
5645 temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
5646 temp32alen, temp32a, temp48);
5647 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5649 finswap = finnow; finnow = finother; finother = finswap;
5651 if (cdytail != 0.0) {
5652 cytablen = scale_expansion_zeroelim(4, ab, cdytail, cytab);
5653 temp16alen = scale_expansion_zeroelim(cytablen, cytab, 2.0 * cdy,
5656 cytaalen = scale_expansion_zeroelim(4, aa, cdytail, cytaa);
5657 temp16blen = scale_expansion_zeroelim(cytaalen, cytaa, bdx, temp16b);
5659 cytbblen = scale_expansion_zeroelim(4, bb, cdytail, cytbb);
5660 temp16clen = scale_expansion_zeroelim(cytbblen, cytbb, -adx, temp16c);
5662 temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5663 temp16blen, temp16b, temp32a);
5664 temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
5665 temp32alen, temp32a, temp48);
5666 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5668 finswap = finnow; finnow = finother; finother = finswap;
5671 if ((adxtail != 0.0) || (adytail != 0.0)) {
5672 if ((bdxtail != 0.0) || (bdytail != 0.0)
5673 || (cdxtail != 0.0) || (cdytail != 0.0)) {
5674 Two_Product(bdxtail, cdy, ti1, ti0);
5675 Two_Product(bdx, cdytail, tj1, tj0);
5676 Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
5679 Two_Product(cdxtail, negate, ti1, ti0);
5681 Two_Product(cdx, negate, tj1, tj0);
5682 Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
5684 bctlen = fast_expansion_sum_zeroelim(4, u, 4, v, bct);
5686 Two_Product(bdxtail, cdytail, ti1, ti0);
5687 Two_Product(cdxtail, bdytail, tj1, tj0);
5688 Two_Two_Diff(ti1, ti0, tj1, tj0, bctt3, bctt[2], bctt[1], bctt[0]);
5698 if (adxtail != 0.0) {
5699 temp16alen = scale_expansion_zeroelim(axtbclen, axtbc, adxtail, temp16a);
5700 axtbctlen = scale_expansion_zeroelim(bctlen, bct, adxtail, axtbct);
5701 temp32alen = scale_expansion_zeroelim(axtbctlen, axtbct, 2.0 * adx,
5703 temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5704 temp32alen, temp32a, temp48);
5705 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5707 finswap = finnow; finnow = finother; finother = finswap;
5708 if (bdytail != 0.0) {
5709 temp8len = scale_expansion_zeroelim(4, cc, adxtail, temp8);
5710 temp16alen = scale_expansion_zeroelim(temp8len, temp8, bdytail,
5712 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
5714 finswap = finnow; finnow = finother; finother = finswap;
5716 if (cdytail != 0.0) {
5717 temp8len = scale_expansion_zeroelim(4, bb, -adxtail, temp8);
5718 temp16alen = scale_expansion_zeroelim(temp8len, temp8, cdytail,
5720 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
5722 finswap = finnow; finnow = finother; finother = finswap;
5725 temp32alen = scale_expansion_zeroelim(axtbctlen, axtbct, adxtail,
5727 axtbcttlen = scale_expansion_zeroelim(bcttlen, bctt, adxtail, axtbctt);
5728 temp16alen = scale_expansion_zeroelim(axtbcttlen, axtbctt, 2.0 * adx,
5730 temp16blen = scale_expansion_zeroelim(axtbcttlen, axtbctt, adxtail,
5732 temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5733 temp16blen, temp16b, temp32b);
5734 temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
5735 temp32blen, temp32b, temp64);
5736 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
5738 finswap = finnow; finnow = finother; finother = finswap;
5740 if (adytail != 0.0) {
5741 temp16alen = scale_expansion_zeroelim(aytbclen, aytbc, adytail, temp16a);
5742 aytbctlen = scale_expansion_zeroelim(bctlen, bct, adytail, aytbct);
5743 temp32alen = scale_expansion_zeroelim(aytbctlen, aytbct, 2.0 * ady,
5745 temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5746 temp32alen, temp32a, temp48);
5747 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5749 finswap = finnow; finnow = finother; finother = finswap;
5752 temp32alen = scale_expansion_zeroelim(aytbctlen, aytbct, adytail,
5754 aytbcttlen = scale_expansion_zeroelim(bcttlen, bctt, adytail, aytbctt);
5755 temp16alen = scale_expansion_zeroelim(aytbcttlen, aytbctt, 2.0 * ady,
5757 temp16blen = scale_expansion_zeroelim(aytbcttlen, aytbctt, adytail,
5759 temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5760 temp16blen, temp16b, temp32b);
5761 temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
5762 temp32blen, temp32b, temp64);
5763 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
5765 finswap = finnow; finnow = finother; finother = finswap;
5768 if ((bdxtail != 0.0) || (bdytail != 0.0)) {
5769 if ((cdxtail != 0.0) || (cdytail != 0.0)
5770 || (adxtail != 0.0) || (adytail != 0.0)) {
5771 Two_Product(cdxtail, ady, ti1, ti0);
5772 Two_Product(cdx, adytail, tj1, tj0);
5773 Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
5776 Two_Product(adxtail, negate, ti1, ti0);
5778 Two_Product(adx, negate, tj1, tj0);
5779 Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
5781 catlen = fast_expansion_sum_zeroelim(4, u, 4, v, cat);
5783 Two_Product(cdxtail, adytail, ti1, ti0);
5784 Two_Product(adxtail, cdytail, tj1, tj0);
5785 Two_Two_Diff(ti1, ti0, tj1, tj0, catt3, catt[2], catt[1], catt[0]);
5795 if (bdxtail != 0.0) {
5796 temp16alen = scale_expansion_zeroelim(bxtcalen, bxtca, bdxtail, temp16a);
5797 bxtcatlen = scale_expansion_zeroelim(catlen, cat, bdxtail, bxtcat);
5798 temp32alen = scale_expansion_zeroelim(bxtcatlen, bxtcat, 2.0 * bdx,
5800 temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5801 temp32alen, temp32a, temp48);
5802 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5804 finswap = finnow; finnow = finother; finother = finswap;
5805 if (cdytail != 0.0) {
5806 temp8len = scale_expansion_zeroelim(4, aa, bdxtail, temp8);
5807 temp16alen = scale_expansion_zeroelim(temp8len, temp8, cdytail,
5809 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
5811 finswap = finnow; finnow = finother; finother = finswap;
5813 if (adytail != 0.0) {
5814 temp8len = scale_expansion_zeroelim(4, cc, -bdxtail, temp8);
5815 temp16alen = scale_expansion_zeroelim(temp8len, temp8, adytail,
5817 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
5819 finswap = finnow; finnow = finother; finother = finswap;
5822 temp32alen = scale_expansion_zeroelim(bxtcatlen, bxtcat, bdxtail,
5824 bxtcattlen = scale_expansion_zeroelim(cattlen, catt, bdxtail, bxtcatt);
5825 temp16alen = scale_expansion_zeroelim(bxtcattlen, bxtcatt, 2.0 * bdx,
5827 temp16blen = scale_expansion_zeroelim(bxtcattlen, bxtcatt, bdxtail,
5829 temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5830 temp16blen, temp16b, temp32b);
5831 temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
5832 temp32blen, temp32b, temp64);
5833 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
5835 finswap = finnow; finnow = finother; finother = finswap;
5837 if (bdytail != 0.0) {
5838 temp16alen = scale_expansion_zeroelim(bytcalen, bytca, bdytail, temp16a);
5839 bytcatlen = scale_expansion_zeroelim(catlen, cat, bdytail, bytcat);
5840 temp32alen = scale_expansion_zeroelim(bytcatlen, bytcat, 2.0 * bdy,
5842 temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5843 temp32alen, temp32a, temp48);
5844 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5846 finswap = finnow; finnow = finother; finother = finswap;
5849 temp32alen = scale_expansion_zeroelim(bytcatlen, bytcat, bdytail,
5851 bytcattlen = scale_expansion_zeroelim(cattlen, catt, bdytail, bytcatt);
5852 temp16alen = scale_expansion_zeroelim(bytcattlen, bytcatt, 2.0 * bdy,
5854 temp16blen = scale_expansion_zeroelim(bytcattlen, bytcatt, bdytail,
5856 temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5857 temp16blen, temp16b, temp32b);
5858 temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
5859 temp32blen, temp32b, temp64);
5860 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
5862 finswap = finnow; finnow = finother; finother = finswap;
5865 if ((cdxtail != 0.0) || (cdytail != 0.0)) {
5866 if ((adxtail != 0.0) || (adytail != 0.0)
5867 || (bdxtail != 0.0) || (bdytail != 0.0)) {
5868 Two_Product(adxtail, bdy, ti1, ti0);
5869 Two_Product(adx, bdytail, tj1, tj0);
5870 Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
5873 Two_Product(bdxtail, negate, ti1, ti0);
5875 Two_Product(bdx, negate, tj1, tj0);
5876 Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
5878 abtlen = fast_expansion_sum_zeroelim(4, u, 4, v, abt);
5880 Two_Product(adxtail, bdytail, ti1, ti0);
5881 Two_Product(bdxtail, adytail, tj1, tj0);
5882 Two_Two_Diff(ti1, ti0, tj1, tj0, abtt3, abtt[2], abtt[1], abtt[0]);
5892 if (cdxtail != 0.0) {
5893 temp16alen = scale_expansion_zeroelim(cxtablen, cxtab, cdxtail, temp16a);
5894 cxtabtlen = scale_expansion_zeroelim(abtlen, abt, cdxtail, cxtabt);
5895 temp32alen = scale_expansion_zeroelim(cxtabtlen, cxtabt, 2.0 * cdx,
5897 temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5898 temp32alen, temp32a, temp48);
5899 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5901 finswap = finnow; finnow = finother; finother = finswap;
5902 if (adytail != 0.0) {
5903 temp8len = scale_expansion_zeroelim(4, bb, cdxtail, temp8);
5904 temp16alen = scale_expansion_zeroelim(temp8len, temp8, adytail,
5906 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
5908 finswap = finnow; finnow = finother; finother = finswap;
5910 if (bdytail != 0.0) {
5911 temp8len = scale_expansion_zeroelim(4, aa, -cdxtail, temp8);
5912 temp16alen = scale_expansion_zeroelim(temp8len, temp8, bdytail,
5914 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
5916 finswap = finnow; finnow = finother; finother = finswap;
5919 temp32alen = scale_expansion_zeroelim(cxtabtlen, cxtabt, cdxtail,
5921 cxtabttlen = scale_expansion_zeroelim(abttlen, abtt, cdxtail, cxtabtt);
5922 temp16alen = scale_expansion_zeroelim(cxtabttlen, cxtabtt, 2.0 * cdx,
5924 temp16blen = scale_expansion_zeroelim(cxtabttlen, cxtabtt, cdxtail,
5926 temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5927 temp16blen, temp16b, temp32b);
5928 temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
5929 temp32blen, temp32b, temp64);
5930 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
5932 finswap = finnow; finnow = finother; finother = finswap;
5934 if (cdytail != 0.0) {
5935 temp16alen = scale_expansion_zeroelim(cytablen, cytab, cdytail, temp16a);
5936 cytabtlen = scale_expansion_zeroelim(abtlen, abt, cdytail, cytabt);
5937 temp32alen = scale_expansion_zeroelim(cytabtlen, cytabt, 2.0 * cdy,
5939 temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5940 temp32alen, temp32a, temp48);
5941 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
5943 finswap = finnow; finnow = finother; finother = finswap;
5946 temp32alen = scale_expansion_zeroelim(cytabtlen, cytabt, cdytail,
5948 cytabttlen = scale_expansion_zeroelim(abttlen, abtt, cdytail, cytabtt);
5949 temp16alen = scale_expansion_zeroelim(cytabttlen, cytabtt, 2.0 * cdy,
5951 temp16blen = scale_expansion_zeroelim(cytabttlen, cytabtt, cdytail,
5953 temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
5954 temp16blen, temp16b, temp32b);
5955 temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
5956 temp32blen, temp32b, temp64);
5957 finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
5959 finswap = finnow; finnow = finother; finother = finswap;
5963 return finnow[finlength - 1];
5966#ifdef ANSI_DECLARATORS
5968 vertex pa, vertex pb, vertex pc, vertex pd)
5970REAL incircle(m, b, pa, pb, pc, pd)
5980 REAL adx, bdx, cdx, ady, bdy, cdy;
5981 REAL bdxcdy, cdxbdy, cdxady, adxcdy, adxbdy, bdxady;
5982 REAL alift, blift, clift;
5984 REAL permanent, errbound;
5988 adx = pa[0] - pd[0];
5989 bdx = pb[0] - pd[0];
5990 cdx = pc[0] - pd[0];
5991 ady = pa[1] - pd[1];
5992 bdy = pb[1] - pd[1];
5993 cdy = pc[1] - pd[1];
5997 alift = adx * adx + ady * ady;
6001 blift = bdx * bdx + bdy * bdy;
6005 clift = cdx * cdx + cdy * cdy;
6007 det = alift * (bdxcdy - cdxbdy)
6008 + blift * (cdxady - adxcdy)
6009 + clift * (adxbdy - bdxady);
6015 permanent = (Absolute(bdxcdy) + Absolute(cdxbdy)) * alift
6016 + (Absolute(cdxady) + Absolute(adxcdy)) * blift
6017 + (Absolute(adxbdy) + Absolute(bdxady)) * clift;
6018 errbound = iccerrboundA * permanent;
6019 if ((det > errbound) || (-det > errbound)) {
6023 return incircleadapt(pa, pb, pc, pd, permanent);
6048#ifdef ANSI_DECLARATORS
6049REAL orient3dadapt(vertex pa, vertex pb, vertex pc, vertex pd,
6050 REAL aheight, REAL bheight, REAL cheight, REAL dheight,
6053REAL orient3dadapt(pa, pb, pc, pd,
6054 aheight, bheight, cheight, dheight, permanent)
6067 INEXACT REAL adx, bdx, cdx, ady, bdy, cdy, adheight, bdheight, cdheight;
6070 INEXACT REAL bdxcdy1, cdxbdy1, cdxady1, adxcdy1, adxbdy1, bdxady1;
6071 REAL bdxcdy0, cdxbdy0, cdxady0, adxcdy0, adxbdy0, bdxady0;
6072 REAL bc[4], ca[4], ab[4];
6073 INEXACT REAL bc3, ca3, ab3;
6074 REAL adet[8], bdet[8], cdet[8];
6075 int alen, blen, clen;
6078 REAL *finnow, *finother, *finswap;
6079 REAL fin1[192], fin2[192];
6082 REAL adxtail, bdxtail, cdxtail;
6083 REAL adytail, bdytail, cdytail;
6084 REAL adheighttail, bdheighttail, cdheighttail;
6085 INEXACT REAL at_blarge, at_clarge;
6086 INEXACT REAL bt_clarge, bt_alarge;
6087 INEXACT REAL ct_alarge, ct_blarge;
6088 REAL at_b[4], at_c[4], bt_c[4], bt_a[4], ct_a[4], ct_b[4];
6089 int at_blen, at_clen, bt_clen, bt_alen, ct_alen, ct_blen;
6090 INEXACT REAL bdxt_cdy1, cdxt_bdy1, cdxt_ady1;
6091 INEXACT REAL adxt_cdy1, adxt_bdy1, bdxt_ady1;
6092 REAL bdxt_cdy0, cdxt_bdy0, cdxt_ady0;
6093 REAL adxt_cdy0, adxt_bdy0, bdxt_ady0;
6094 INEXACT REAL bdyt_cdx1, cdyt_bdx1, cdyt_adx1;
6095 INEXACT REAL adyt_cdx1, adyt_bdx1, bdyt_adx1;
6096 REAL bdyt_cdx0, cdyt_bdx0, cdyt_adx0;
6097 REAL adyt_cdx0, adyt_bdx0, bdyt_adx0;
6098 REAL bct[8], cat[8], abt[8];
6099 int bctlen, catlen, abtlen;
6100 INEXACT REAL bdxt_cdyt1, cdxt_bdyt1, cdxt_adyt1;
6101 INEXACT REAL adxt_cdyt1, adxt_bdyt1, bdxt_adyt1;
6102 REAL bdxt_cdyt0, cdxt_bdyt0, cdxt_adyt0;
6103 REAL adxt_cdyt0, adxt_bdyt0, bdxt_adyt0;
6104 REAL u[4], v[12], w[16];
6106 int vlength, wlength;
6110 REAL avirt, bround, around;
6113 REAL ahi, alo, bhi, blo;
6114 REAL err1, err2, err3;
6115 INEXACT REAL _i, _j, _k;
6118 adx = (REAL) (pa[0] - pd[0]);
6119 bdx = (REAL) (pb[0] - pd[0]);
6120 cdx = (REAL) (pc[0] - pd[0]);
6121 ady = (REAL) (pa[1] - pd[1]);
6122 bdy = (REAL) (pb[1] - pd[1]);
6123 cdy = (REAL) (pc[1] - pd[1]);
6124 adheight = (REAL) (aheight - dheight);
6125 bdheight = (REAL) (bheight - dheight);
6126 cdheight = (REAL) (cheight - dheight);
6128 Two_Product(bdx, cdy, bdxcdy1, bdxcdy0);
6129 Two_Product(cdx, bdy, cdxbdy1, cdxbdy0);
6130 Two_Two_Diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0, bc3, bc[2], bc[1], bc[0]);
6132 alen = scale_expansion_zeroelim(4, bc, adheight, adet);
6134 Two_Product(cdx, ady, cdxady1, cdxady0);
6135 Two_Product(adx, cdy, adxcdy1, adxcdy0);
6136 Two_Two_Diff(cdxady1, cdxady0, adxcdy1, adxcdy0, ca3, ca[2], ca[1], ca[0]);
6138 blen = scale_expansion_zeroelim(4, ca, bdheight, bdet);
6140 Two_Product(adx, bdy, adxbdy1, adxbdy0);
6141 Two_Product(bdx, ady, bdxady1, bdxady0);
6142 Two_Two_Diff(adxbdy1, adxbdy0, bdxady1, bdxady0, ab3, ab[2], ab[1], ab[0]);
6144 clen = scale_expansion_zeroelim(4, ab, cdheight, cdet);
6146 ablen = fast_expansion_sum_zeroelim(alen, adet, blen, bdet, abdet);
6147 finlength = fast_expansion_sum_zeroelim(ablen, abdet, clen, cdet, fin1);
6149 det = estimate(finlength, fin1);
6150 errbound = o3derrboundB * permanent;
6151 if ((det >= errbound) || (-det >= errbound)) {
6155 Two_Diff_Tail(pa[0], pd[0], adx, adxtail);
6156 Two_Diff_Tail(pb[0], pd[0], bdx, bdxtail);
6157 Two_Diff_Tail(pc[0], pd[0], cdx, cdxtail);
6158 Two_Diff_Tail(pa[1], pd[1], ady, adytail);
6159 Two_Diff_Tail(pb[1], pd[1], bdy, bdytail);
6160 Two_Diff_Tail(pc[1], pd[1], cdy, cdytail);
6161 Two_Diff_Tail(aheight, dheight, adheight, adheighttail);
6162 Two_Diff_Tail(bheight, dheight, bdheight, bdheighttail);
6163 Two_Diff_Tail(cheight, dheight, cdheight, cdheighttail);
6165 if ((adxtail == 0.0) && (bdxtail == 0.0) && (cdxtail == 0.0) &&
6166 (adytail == 0.0) && (bdytail == 0.0) && (cdytail == 0.0) &&
6167 (adheighttail == 0.0) &&
6168 (bdheighttail == 0.0) &&
6169 (cdheighttail == 0.0)) {
6173 errbound = o3derrboundC * permanent + resulterrbound * Absolute(det);
6174 det += (adheight * ((bdx * cdytail + cdy * bdxtail) -
6175 (bdy * cdxtail + cdx * bdytail)) +
6176 adheighttail * (bdx * cdy - bdy * cdx)) +
6177 (bdheight * ((cdx * adytail + ady * cdxtail) -
6178 (cdy * adxtail + adx * cdytail)) +
6179 bdheighttail * (cdx * ady - cdy * adx)) +
6180 (cdheight * ((adx * bdytail + bdy * adxtail) -
6181 (ady * bdxtail + bdx * adytail)) +
6182 cdheighttail * (adx * bdy - ady * bdx));
6183 if ((det >= errbound) || (-det >= errbound)) {
6190 if (adxtail == 0.0) {
6191 if (adytail == 0.0) {
6198 Two_Product(negate, bdx, at_blarge, at_b[0]);
6199 at_b[1] = at_blarge;
6201 Two_Product(adytail, cdx, at_clarge, at_c[0]);
6202 at_c[1] = at_clarge;
6206 if (adytail == 0.0) {
6207 Two_Product(adxtail, bdy, at_blarge, at_b[0]);
6208 at_b[1] = at_blarge;
6211 Two_Product(negate, cdy, at_clarge, at_c[0]);
6212 at_c[1] = at_clarge;
6215 Two_Product(adxtail, bdy, adxt_bdy1, adxt_bdy0);
6216 Two_Product(adytail, bdx, adyt_bdx1, adyt_bdx0);
6217 Two_Two_Diff(adxt_bdy1, adxt_bdy0, adyt_bdx1, adyt_bdx0,
6218 at_blarge, at_b[2], at_b[1], at_b[0]);
6219 at_b[3] = at_blarge;
6221 Two_Product(adytail, cdx, adyt_cdx1, adyt_cdx0);
6222 Two_Product(adxtail, cdy, adxt_cdy1, adxt_cdy0);
6223 Two_Two_Diff(adyt_cdx1, adyt_cdx0, adxt_cdy1, adxt_cdy0,
6224 at_clarge, at_c[2], at_c[1], at_c[0]);
6225 at_c[3] = at_clarge;
6229 if (bdxtail == 0.0) {
6230 if (bdytail == 0.0) {
6237 Two_Product(negate, cdx, bt_clarge, bt_c[0]);
6238 bt_c[1] = bt_clarge;
6240 Two_Product(bdytail, adx, bt_alarge, bt_a[0]);
6241 bt_a[1] = bt_alarge;
6245 if (bdytail == 0.0) {
6246 Two_Product(bdxtail, cdy, bt_clarge, bt_c[0]);
6247 bt_c[1] = bt_clarge;
6250 Two_Product(negate, ady, bt_alarge, bt_a[0]);
6251 bt_a[1] = bt_alarge;
6254 Two_Product(bdxtail, cdy, bdxt_cdy1, bdxt_cdy0);
6255 Two_Product(bdytail, cdx, bdyt_cdx1, bdyt_cdx0);
6256 Two_Two_Diff(bdxt_cdy1, bdxt_cdy0, bdyt_cdx1, bdyt_cdx0,
6257 bt_clarge, bt_c[2], bt_c[1], bt_c[0]);
6258 bt_c[3] = bt_clarge;
6260 Two_Product(bdytail, adx, bdyt_adx1, bdyt_adx0);
6261 Two_Product(bdxtail, ady, bdxt_ady1, bdxt_ady0);
6262 Two_Two_Diff(bdyt_adx1, bdyt_adx0, bdxt_ady1, bdxt_ady0,
6263 bt_alarge, bt_a[2], bt_a[1], bt_a[0]);
6264 bt_a[3] = bt_alarge;
6268 if (cdxtail == 0.0) {
6269 if (cdytail == 0.0) {
6276 Two_Product(negate, adx, ct_alarge, ct_a[0]);
6277 ct_a[1] = ct_alarge;
6279 Two_Product(cdytail, bdx, ct_blarge, ct_b[0]);
6280 ct_b[1] = ct_blarge;
6284 if (cdytail == 0.0) {
6285 Two_Product(cdxtail, ady, ct_alarge, ct_a[0]);
6286 ct_a[1] = ct_alarge;
6289 Two_Product(negate, bdy, ct_blarge, ct_b[0]);
6290 ct_b[1] = ct_blarge;
6293 Two_Product(cdxtail, ady, cdxt_ady1, cdxt_ady0);
6294 Two_Product(cdytail, adx, cdyt_adx1, cdyt_adx0);
6295 Two_Two_Diff(cdxt_ady1, cdxt_ady0, cdyt_adx1, cdyt_adx0,
6296 ct_alarge, ct_a[2], ct_a[1], ct_a[0]);
6297 ct_a[3] = ct_alarge;
6299 Two_Product(cdytail, bdx, cdyt_bdx1, cdyt_bdx0);
6300 Two_Product(cdxtail, bdy, cdxt_bdy1, cdxt_bdy0);
6301 Two_Two_Diff(cdyt_bdx1, cdyt_bdx0, cdxt_bdy1, cdxt_bdy0,
6302 ct_blarge, ct_b[2], ct_b[1], ct_b[0]);
6303 ct_b[3] = ct_blarge;
6308 bctlen = fast_expansion_sum_zeroelim(bt_clen, bt_c, ct_blen, ct_b, bct);
6309 wlength = scale_expansion_zeroelim(bctlen, bct, adheight, w);
6310 finlength = fast_expansion_sum_zeroelim(finlength, finnow, wlength, w,
6312 finswap = finnow; finnow = finother; finother = finswap;
6314 catlen = fast_expansion_sum_zeroelim(ct_alen, ct_a, at_clen, at_c, cat);
6315 wlength = scale_expansion_zeroelim(catlen, cat, bdheight, w);
6316 finlength = fast_expansion_sum_zeroelim(finlength, finnow, wlength, w,
6318 finswap = finnow; finnow = finother; finother = finswap;
6320 abtlen = fast_expansion_sum_zeroelim(at_blen, at_b, bt_alen, bt_a, abt);
6321 wlength = scale_expansion_zeroelim(abtlen, abt, cdheight, w);
6322 finlength = fast_expansion_sum_zeroelim(finlength, finnow, wlength, w,
6324 finswap = finnow; finnow = finother; finother = finswap;
6326 if (adheighttail != 0.0) {
6327 vlength = scale_expansion_zeroelim(4, bc, adheighttail, v);
6328 finlength = fast_expansion_sum_zeroelim(finlength, finnow, vlength, v,
6330 finswap = finnow; finnow = finother; finother = finswap;
6332 if (bdheighttail != 0.0) {
6333 vlength = scale_expansion_zeroelim(4, ca, bdheighttail, v);
6334 finlength = fast_expansion_sum_zeroelim(finlength, finnow, vlength, v,
6336 finswap = finnow; finnow = finother; finother = finswap;
6338 if (cdheighttail != 0.0) {
6339 vlength = scale_expansion_zeroelim(4, ab, cdheighttail, v);
6340 finlength = fast_expansion_sum_zeroelim(finlength, finnow, vlength, v,
6342 finswap = finnow; finnow = finother; finother = finswap;
6345 if (adxtail != 0.0) {
6346 if (bdytail != 0.0) {
6347 Two_Product(adxtail, bdytail, adxt_bdyt1, adxt_bdyt0);
6348 Two_One_Product(adxt_bdyt1, adxt_bdyt0, cdheight, u3, u[2], u[1], u[0]);
6350 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6352 finswap = finnow; finnow = finother; finother = finswap;
6353 if (cdheighttail != 0.0) {
6354 Two_One_Product(adxt_bdyt1, adxt_bdyt0, cdheighttail,
6355 u3, u[2], u[1], u[0]);
6357 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6359 finswap = finnow; finnow = finother; finother = finswap;
6362 if (cdytail != 0.0) {
6364 Two_Product(negate, cdytail, adxt_cdyt1, adxt_cdyt0);
6365 Two_One_Product(adxt_cdyt1, adxt_cdyt0, bdheight, u3, u[2], u[1], u[0]);
6367 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6369 finswap = finnow; finnow = finother; finother = finswap;
6370 if (bdheighttail != 0.0) {
6371 Two_One_Product(adxt_cdyt1, adxt_cdyt0, bdheighttail,
6372 u3, u[2], u[1], u[0]);
6374 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6376 finswap = finnow; finnow = finother; finother = finswap;
6380 if (bdxtail != 0.0) {
6381 if (cdytail != 0.0) {
6382 Two_Product(bdxtail, cdytail, bdxt_cdyt1, bdxt_cdyt0);
6383 Two_One_Product(bdxt_cdyt1, bdxt_cdyt0, adheight, u3, u[2], u[1], u[0]);
6385 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6387 finswap = finnow; finnow = finother; finother = finswap;
6388 if (adheighttail != 0.0) {
6389 Two_One_Product(bdxt_cdyt1, bdxt_cdyt0, adheighttail,
6390 u3, u[2], u[1], u[0]);
6392 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6394 finswap = finnow; finnow = finother; finother = finswap;
6397 if (adytail != 0.0) {
6399 Two_Product(negate, adytail, bdxt_adyt1, bdxt_adyt0);
6400 Two_One_Product(bdxt_adyt1, bdxt_adyt0, cdheight, u3, u[2], u[1], u[0]);
6402 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6404 finswap = finnow; finnow = finother; finother = finswap;
6405 if (cdheighttail != 0.0) {
6406 Two_One_Product(bdxt_adyt1, bdxt_adyt0, cdheighttail,
6407 u3, u[2], u[1], u[0]);
6409 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6411 finswap = finnow; finnow = finother; finother = finswap;
6415 if (cdxtail != 0.0) {
6416 if (adytail != 0.0) {
6417 Two_Product(cdxtail, adytail, cdxt_adyt1, cdxt_adyt0);
6418 Two_One_Product(cdxt_adyt1, cdxt_adyt0, bdheight, u3, u[2], u[1], u[0]);
6420 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6422 finswap = finnow; finnow = finother; finother = finswap;
6423 if (bdheighttail != 0.0) {
6424 Two_One_Product(cdxt_adyt1, cdxt_adyt0, bdheighttail,
6425 u3, u[2], u[1], u[0]);
6427 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6429 finswap = finnow; finnow = finother; finother = finswap;
6432 if (bdytail != 0.0) {
6434 Two_Product(negate, bdytail, cdxt_bdyt1, cdxt_bdyt0);
6435 Two_One_Product(cdxt_bdyt1, cdxt_bdyt0, adheight, u3, u[2], u[1], u[0]);
6437 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6439 finswap = finnow; finnow = finother; finother = finswap;
6440 if (adheighttail != 0.0) {
6441 Two_One_Product(cdxt_bdyt1, cdxt_bdyt0, adheighttail,
6442 u3, u[2], u[1], u[0]);
6444 finlength = fast_expansion_sum_zeroelim(finlength, finnow, 4, u,
6446 finswap = finnow; finnow = finother; finother = finswap;
6451 if (adheighttail != 0.0) {
6452 wlength = scale_expansion_zeroelim(bctlen, bct, adheighttail, w);
6453 finlength = fast_expansion_sum_zeroelim(finlength, finnow, wlength, w,
6455 finswap = finnow; finnow = finother; finother = finswap;
6457 if (bdheighttail != 0.0) {
6458 wlength = scale_expansion_zeroelim(catlen, cat, bdheighttail, w);
6459 finlength = fast_expansion_sum_zeroelim(finlength, finnow, wlength, w,
6461 finswap = finnow; finnow = finother; finother = finswap;
6463 if (cdheighttail != 0.0) {
6464 wlength = scale_expansion_zeroelim(abtlen, abt, cdheighttail, w);
6465 finlength = fast_expansion_sum_zeroelim(finlength, finnow, wlength, w,
6467 finswap = finnow; finnow = finother; finother = finswap;
6470 return finnow[finlength - 1];
6473#ifdef ANSI_DECLARATORS
6475 vertex pa, vertex pb, vertex pc, vertex pd,
6476 REAL aheight, REAL bheight, REAL cheight, REAL dheight)
6478REAL orient3d(m, b, pa, pb, pc, pd, aheight, bheight, cheight, dheight)
6492 REAL adx, bdx, cdx, ady, bdy, cdy, adheight, bdheight, cdheight;
6493 REAL bdxcdy, cdxbdy, cdxady, adxcdy, adxbdy, bdxady;
6495 REAL permanent, errbound;
6499 adx = pa[0] - pd[0];
6500 bdx = pb[0] - pd[0];
6501 cdx = pc[0] - pd[0];
6502 ady = pa[1] - pd[1];
6503 bdy = pb[1] - pd[1];
6504 cdy = pc[1] - pd[1];
6505 adheight = aheight - dheight;
6506 bdheight = bheight - dheight;
6507 cdheight = cheight - dheight;
6518 det = adheight * (bdxcdy - cdxbdy)
6519 + bdheight * (cdxady - adxcdy)
6520 + cdheight * (adxbdy - bdxady);
6526 permanent = (Absolute(bdxcdy) + Absolute(cdxbdy)) * Absolute(adheight)
6527 + (Absolute(cdxady) + Absolute(adxcdy)) * Absolute(bdheight)
6528 + (Absolute(adxbdy) + Absolute(bdxady)) * Absolute(cdheight);
6529 errbound = o3derrboundA * permanent;
6530 if ((det > errbound) || (-det > errbound)) {
6534 return orient3dadapt(pa, pb, pc, pd, aheight, bheight, cheight, dheight,
6556#ifdef ANSI_DECLARATORS
6558 vertex pa, vertex pb, vertex pc, vertex pd)
6560REAL nonregular(m, b, pa, pb, pc, pd)
6570 if (b->weighted == 0) {
6571 return incircle(m, b, pa, pb, pc, pd);
6572 }
else if (b->weighted == 1) {
6573 return orient3d(m, b, pa, pb, pc, pd,
6574 pa[0] * pa[0] + pa[1] * pa[1] - pa[2],
6575 pb[0] * pb[0] + pb[1] * pb[1] - pb[2],
6576 pc[0] * pc[0] + pc[1] * pc[1] - pc[2],
6577 pd[0] * pd[0] + pd[1] * pd[1] - pd[2]);
6579 return orient3d(m, b, pa, pb, pc, pd, pa[2], pb[2], pc[2], pd[2]);
6597#ifdef ANSI_DECLARATORS
6598void findcircumcenter(
struct mesh *m,
struct behavior *b,
6599 vertex torg, vertex tdest, vertex tapex,
6600 vertex circumcenter, REAL *xi, REAL *eta,
int offcenter)
6602void findcircumcenter(m, b, torg, tdest, tapex, circumcenter, xi, eta,
6616 REAL xdo, ydo, xao, yao;
6617 REAL dodist, aodist, dadist;
6619 REAL dx, dy, dxoff, dyoff;
6621 m->circumcentercount++;
6624 xdo = tdest[0] - torg[0];
6625 ydo = tdest[1] - torg[1];
6626 xao = tapex[0] - torg[0];
6627 yao = tapex[1] - torg[1];
6628 dodist = xdo * xdo + ydo * ydo;
6629 aodist = xao * xao + yao * yao;
6630 dadist = (tdest[0] - tapex[0]) * (tdest[0] - tapex[0]) +
6631 (tdest[1] - tapex[1]) * (tdest[1] - tapex[1]);
6633 denominator = 0.5 / (xdo * yao - xao * ydo);
6638 denominator = 0.5 / counterclockwise(m, b, tdest, tapex, torg);
6640 m->counterclockcount--;
6642 dx = (yao * dodist - ydo * aodist) * denominator;
6643 dy = (xdo * aodist - xao * dodist) * denominator;
6650 if ((dodist < aodist) && (dodist < dadist)) {
6651 if (offcenter && (b->offconstant > 0.0)) {
6653 dxoff = 0.5 * xdo - b->offconstant * ydo;
6654 dyoff = 0.5 * ydo + b->offconstant * xdo;
6657 if (dxoff * dxoff + dyoff * dyoff < dx * dx + dy * dy) {
6662 }
else if (aodist < dadist) {
6663 if (offcenter && (b->offconstant > 0.0)) {
6664 dxoff = 0.5 * xao + b->offconstant * yao;
6665 dyoff = 0.5 * yao - b->offconstant * xao;
6668 if (dxoff * dxoff + dyoff * dyoff < dx * dx + dy * dy) {
6674 if (offcenter && (b->offconstant > 0.0)) {
6675 dxoff = 0.5 * (tapex[0] - tdest[0]) -
6676 b->offconstant * (tapex[1] - tdest[1]);
6677 dyoff = 0.5 * (tapex[1] - tdest[1]) +
6678 b->offconstant * (tapex[0] - tdest[0]);
6681 if (dxoff * dxoff + dyoff * dyoff <
6682 (dx - xdo) * (dx - xdo) + (dy - ydo) * (dy - ydo)) {
6689 circumcenter[0] = torg[0] + dx;
6690 circumcenter[1] = torg[1] + dy;
6697 *xi = (yao * dx - xao * dy) * (2.0 * denominator);
6698 *eta = (xdo * dy - ydo * dx) * (2.0 * denominator);
6711#ifdef ANSI_DECLARATORS
6712void triangleinit(
struct mesh *m)
6719 poolzero(&m->vertices);
6720 poolzero(&m->triangles);
6721 poolzero(&m->subsegs);
6723 poolzero(&m->badsubsegs);
6724 poolzero(&m->badtriangles);
6725 poolzero(&m->flipstackers);
6726 poolzero(&m->splaynodes);
6728 m->recenttri.tri = (triangle *) NULL;
6731 m->checksegments = 0;
6732 m->checkquality = 0;
6733 m->incirclecount = m->counterclockcount = m->orient3dcount = 0;
6734 m->hyperbolacount = m->circletopcount = m->circumcentercount = 0;
6750#ifdef ANSI_DECLARATORS
6751unsigned long randomnation(
unsigned int choices)
6753unsigned long randomnation(choices)
6754unsigned int choices;
6758 randomseed = (randomseed * 1366l + 150889l) % 714025l;
6759 return randomseed / (714025l / choices + 1);
6774#ifdef ANSI_DECLARATORS
6783 struct otri triangleloop;
6784 struct otri oppotri, oppooppotri;
6785 vertex triorg, tridest, triapex;
6786 vertex oppoorg, oppodest;
6792 saveexact = b->noexact;
6795 printf(
" Checking consistency of mesh...\n");
6799 traversalinit(&m->triangles);
6800 triangleloop.tri = triangletraverse(m);
6801 while (triangleloop.tri != (triangle *) NULL) {
6803 for (triangleloop.orient = 0; triangleloop.orient < 3;
6804 triangleloop.orient++) {
6805 org(triangleloop, triorg);
6806 dest(triangleloop, tridest);
6807 if (triangleloop.orient == 0) {
6809 apex(triangleloop, triapex);
6810 if (counterclockwise(m, b, triorg, tridest, triapex) <= 0.0) {
6811 printf(
" !! !! Inverted ");
6812 printtriangle(m, b, &triangleloop);
6817 sym(triangleloop, oppotri);
6818 if (oppotri.tri != m->dummytri) {
6820 sym(oppotri, oppooppotri);
6821 if ((triangleloop.tri != oppooppotri.tri)
6822 || (triangleloop.orient != oppooppotri.orient)) {
6823 printf(
" !! !! Asymmetric triangle-triangle bond:\n");
6824 if (triangleloop.tri == oppooppotri.tri) {
6825 printf(
" (Right triangle, wrong orientation)\n");
6828 printtriangle(m, b, &triangleloop);
6829 printf(
" Second (nonreciprocating) ");
6830 printtriangle(m, b, &oppotri);
6835 org(oppotri, oppoorg);
6836 dest(oppotri, oppodest);
6837 if ((triorg != oppodest) || (tridest != oppoorg)) {
6838 printf(
" !! !! Mismatched edge coordinates between two triangles:\n"
6840 printf(
" First mismatched ");
6841 printtriangle(m, b, &triangleloop);
6842 printf(
" Second mismatched ");
6843 printtriangle(m, b, &oppotri);
6848 triangleloop.tri = triangletraverse(m);
6852 printf(
" In my studied opinion, the mesh appears to be consistent.\n");
6854 }
else if (horrors == 1) {
6855 printf(
" !! !! !! !! Precisely one festering wound discovered.\n");
6857 printf(
" !! !! !! !! %d abominations witnessed.\n", horrors);
6860 b->noexact = saveexact;
6873#ifdef ANSI_DECLARATORS
6876void checkdelaunay(m, b)
6882 struct otri triangleloop;
6883 struct otri oppotri;
6884 struct osub opposubseg;
6885 vertex triorg, tridest, triapex;
6887 int shouldbedelaunay;
6894 saveexact = b->noexact;
6897 printf(
" Checking Delaunay property of mesh...\n");
6901 traversalinit(&m->triangles);
6902 triangleloop.tri = triangletraverse(m);
6903 while (triangleloop.tri != (triangle *) NULL) {
6905 for (triangleloop.orient = 0; triangleloop.orient < 3;
6906 triangleloop.orient++) {
6907 org(triangleloop, triorg);
6908 dest(triangleloop, tridest);
6909 apex(triangleloop, triapex);
6910 sym(triangleloop, oppotri);
6911 apex(oppotri, oppoapex);
6915 shouldbedelaunay = (oppotri.tri != m->dummytri) &&
6916 !deadtri(oppotri.tri) && (triangleloop.tri < oppotri.tri) &&
6917 (triorg != m->infvertex1) && (triorg != m->infvertex2) &&
6918 (triorg != m->infvertex3) &&
6919 (tridest != m->infvertex1) && (tridest != m->infvertex2) &&
6920 (tridest != m->infvertex3) &&
6921 (triapex != m->infvertex1) && (triapex != m->infvertex2) &&
6922 (triapex != m->infvertex3) &&
6923 (oppoapex != m->infvertex1) && (oppoapex != m->infvertex2) &&
6924 (oppoapex != m->infvertex3);
6925 if (m->checksegments && shouldbedelaunay) {
6928 tspivot(triangleloop, opposubseg);
6929 if (opposubseg.ss != m->dummysub){
6930 shouldbedelaunay = 0;
6933 if (shouldbedelaunay) {
6934 if (nonregular(m, b, triorg, tridest, triapex, oppoapex) > 0.0) {
6936 printf(
" !! !! Non-Delaunay pair of triangles:\n");
6937 printf(
" First non-Delaunay ");
6938 printtriangle(m, b, &triangleloop);
6939 printf(
" Second non-Delaunay ");
6941 printf(
" !! !! Non-regular pair of triangles:\n");
6942 printf(
" First non-regular ");
6943 printtriangle(m, b, &triangleloop);
6944 printf(
" Second non-regular ");
6946 printtriangle(m, b, &oppotri);
6951 triangleloop.tri = triangletraverse(m);
6956 " By virtue of my perceptive intelligence, I declare the mesh Delaunay.\n");
6958 }
else if (horrors == 1) {
6960 " !! !! !! !! Precisely one terrifying transgression identified.\n");
6962 printf(
" !! !! !! !! %d obscenities viewed with horror.\n", horrors);
6965 b->noexact = saveexact;
6983#ifdef ANSI_DECLARATORS
6984void enqueuebadtriang(
struct mesh *m,
struct behavior *b,
6987void enqueuebadtriang(m, b, badtri)
6994 REAL length, multiplier;
6995 int exponent, expincrement;
7000 if (b->verbose > 2) {
7001 printf(
" Queueing bad triangle:\n");
7002 printf(
" (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
7003 badtri->triangorg[0], badtri->triangorg[1],
7004 badtri->triangdest[0], badtri->triangdest[1],
7005 badtri->triangapex[0], badtri->triangapex[1]);
7010 if (badtri->key >= 1.0) {
7011 length = badtri->key;
7016 length = 1.0 / badtri->key;
7022 while (length > 2.0) {
7026 while (length * multiplier * multiplier > 1.0) {
7028 multiplier *= multiplier;
7031 exponent += expincrement;
7032 length *= multiplier;
7035 exponent = 2 * exponent + (length > SQUAREROOTTWO);
7040 queuenumber = 2047 - exponent;
7042 queuenumber = 2048 + exponent;
7046 if (m->queuefront[queuenumber] == (
struct badtriang *) NULL) {
7049 if (queuenumber > m->firstnonemptyq) {
7051 m->nextnonemptyq[queuenumber] = m->firstnonemptyq;
7052 m->firstnonemptyq = queuenumber;
7056 i = queuenumber + 1;
7057 while (m->queuefront[i] == (
struct badtriang *) NULL) {
7061 m->nextnonemptyq[queuenumber] = m->nextnonemptyq[i];
7062 m->nextnonemptyq[i] = queuenumber;
7065 m->queuefront[queuenumber] = badtri;
7068 m->queuetail[queuenumber]->nexttriang = badtri;
7071 m->queuetail[queuenumber] = badtri;
7073 badtri->nexttriang = (
struct badtriang *) NULL;
7089#ifdef ANSI_DECLARATORS
7091 REAL minedge, vertex enqapex, vertex enqorg, vertex enqdest)
7093void enqueuebadtri(m, b, enqtri, minedge, enqapex, enqorg, enqdest)
7107 newbad = (
struct badtriang *) poolalloc(&m->badtriangles);
7108 newbad->poortri = encode(*enqtri);
7109 newbad->key = minedge;
7110 newbad->triangapex = enqapex;
7111 newbad->triangorg = enqorg;
7112 newbad->triangdest = enqdest;
7113 enqueuebadtriang(m, b, newbad);
7126#ifdef ANSI_DECLARATORS
7137 if (m->firstnonemptyq < 0) {
7141 result = m->queuefront[m->firstnonemptyq];
7143 m->queuefront[m->firstnonemptyq] = result->nexttriang;
7146 if (result == m->queuetail[m->firstnonemptyq]) {
7147 m->firstnonemptyq = m->nextnonemptyq[m->firstnonemptyq];
7179#ifdef ANSI_DECLARATORS
7180int checkseg4encroach(
struct mesh *m,
struct behavior *b,
7181 struct osub *testsubseg)
7183int checkseg4encroach(m, b, testsubseg)
7186struct osub *testsubseg;
7190 struct otri neighbortri;
7191 struct osub testsym;
7196 vertex eorg, edest, eapex;
7202 sorg(*testsubseg, eorg);
7203 sdest(*testsubseg, edest);
7205 stpivot(*testsubseg, neighbortri);
7207 if (neighbortri.tri != m->dummytri) {
7210 apex(neighbortri, eapex);
7216 dotproduct = (eorg[0] - eapex[0]) * (edest[0] - eapex[0]) +
7217 (eorg[1] - eapex[1]) * (edest[1] - eapex[1]);
7218 if (dotproduct < 0.0) {
7219 if (b->conformdel ||
7220 (dotproduct * dotproduct >=
7221 (2.0 * b->goodangle - 1.0) * (2.0 * b->goodangle - 1.0) *
7222 ((eorg[0] - eapex[0]) * (eorg[0] - eapex[0]) +
7223 (eorg[1] - eapex[1]) * (eorg[1] - eapex[1])) *
7224 ((edest[0] - eapex[0]) * (edest[0] - eapex[0]) +
7225 (edest[1] - eapex[1]) * (edest[1] - eapex[1])))) {
7231 ssym(*testsubseg, testsym);
7232 stpivot(testsym, neighbortri);
7234 if (neighbortri.tri != m->dummytri) {
7237 apex(neighbortri, eapex);
7240 dotproduct = (eorg[0] - eapex[0]) * (edest[0] - eapex[0]) +
7241 (eorg[1] - eapex[1]) * (edest[1] - eapex[1]);
7242 if (dotproduct < 0.0) {
7243 if (b->conformdel ||
7244 (dotproduct * dotproduct >=
7245 (2.0 * b->goodangle - 1.0) * (2.0 * b->goodangle - 1.0) *
7246 ((eorg[0] - eapex[0]) * (eorg[0] - eapex[0]) +
7247 (eorg[1] - eapex[1]) * (eorg[1] - eapex[1])) *
7248 ((edest[0] - eapex[0]) * (edest[0] - eapex[0]) +
7249 (edest[1] - eapex[1]) * (edest[1] - eapex[1])))) {
7255 if (encroached && (!b->nobisect || ((b->nobisect == 1) && (sides == 2)))) {
7256 if (b->verbose > 2) {
7258 " Queueing encroached subsegment (%.12g, %.12g) (%.12g, %.12g).\n",
7259 eorg[0], eorg[1], edest[0], edest[1]);
7263 encroachedseg = (
struct badsubseg *) poolalloc(&m->badsubsegs);
7264 if (encroached == 1) {
7265 encroachedseg->encsubseg = sencode(*testsubseg);
7266 encroachedseg->subsegorg = eorg;
7267 encroachedseg->subsegdest = edest;
7269 encroachedseg->encsubseg = sencode(testsym);
7270 encroachedseg->subsegorg = edest;
7271 encroachedseg->subsegdest = eorg;
7292#ifdef ANSI_DECLARATORS
7295void testtriangle(m, b, testtri)
7298struct otri *testtri;
7302 struct otri tri1, tri2;
7303 struct osub testsub;
7304 vertex torg, tdest, tapex;
7305 vertex base1, base2;
7306 vertex org1, dest1, org2, dest2;
7308 REAL dxod, dyod, dxda, dyda, dxao, dyao;
7309 REAL dxod2, dyod2, dxda2, dyda2, dxao2, dyao2;
7310 REAL apexlen, orglen, destlen, minedge;
7317 org(*testtri, torg);
7318 dest(*testtri, tdest);
7319 apex(*testtri, tapex);
7320 dxod = torg[0] - tdest[0];
7321 dyod = torg[1] - tdest[1];
7322 dxda = tdest[0] - tapex[0];
7323 dyda = tdest[1] - tapex[1];
7324 dxao = tapex[0] - torg[0];
7325 dyao = tapex[1] - torg[1];
7326 dxod2 = dxod * dxod;
7327 dyod2 = dyod * dyod;
7328 dxda2 = dxda * dxda;
7329 dyda2 = dyda * dyda;
7330 dxao2 = dxao * dxao;
7331 dyao2 = dyao * dyao;
7333 apexlen = dxod2 + dyod2;
7334 orglen = dxda2 + dyda2;
7335 destlen = dxao2 + dyao2;
7337 if ((apexlen < orglen) && (apexlen < destlen)) {
7341 angle = dxda * dxao + dyda * dyao;
7342 angle = angle * angle / (orglen * destlen);
7345 otricopy(*testtri, tri1);
7346 }
else if (orglen < destlen) {
7350 angle = dxod * dxao + dyod * dyao;
7351 angle = angle * angle / (apexlen * destlen);
7354 lnext(*testtri, tri1);
7359 angle = dxod * dxda + dyod * dyda;
7360 angle = angle * angle / (apexlen * orglen);
7363 lprev(*testtri, tri1);
7366 if (b->vararea || b->fixedarea || b->usertest) {
7368 area = 0.5 * (dxod * dyda - dyod * dxda);
7369 if (b->fixedarea && (area > b->maxarea)) {
7371 enqueuebadtri(m, b, testtri, minedge, tapex, torg, tdest);
7376 if ((b->vararea) && (area > areabound(*testtri)) &&
7377 (areabound(*testtri) > 0.0)) {
7379 enqueuebadtri(m, b, testtri, minedge, tapex, torg, tdest);
7385 if (triunsuitable(torg, tdest, tapex, area)) {
7386 enqueuebadtri(m, b, testtri, minedge, tapex, torg, tdest);
7393 if (angle > b->goodangle) {
7403 if ((vertextype(base1) == SEGMENTVERTEX) &&
7404 (vertextype(base2) == SEGMENTVERTEX)) {
7407 tspivot(tri1, testsub);
7408 if (testsub.ss == m->dummysub) {
7410 otricopy(tri1, tri2);
7413 tspivot(tri1, testsub);
7414 }
while (testsub.ss == m->dummysub);
7416 segorg(testsub, org1);
7417 segdest(testsub, dest1);
7421 tspivot(tri2, testsub);
7422 }
while (testsub.ss == m->dummysub);
7424 segorg(testsub, org2);
7425 segdest(testsub, dest2);
7427 joinvertex = (vertex) NULL;
7428 if ((dest1[0] == org2[0]) && (dest1[1] == org2[1])) {
7430 }
else if ((org1[0] == dest2[0]) && (org1[1] == dest2[1])) {
7433 if (joinvertex != (vertex) NULL) {
7436 dist1 = ((base1[0] - joinvertex[0]) * (base1[0] - joinvertex[0]) +
7437 (base1[1] - joinvertex[1]) * (base1[1] - joinvertex[1]));
7438 dist2 = ((base2[0] - joinvertex[0]) * (base2[0] - joinvertex[0]) +
7439 (base2[1] - joinvertex[1]) * (base2[1] - joinvertex[1]));
7441 if ((dist1 < 1.001 * dist2) && (dist1 > 0.999 * dist2)) {
7450 enqueuebadtri(m, b, testtri, minedge, tapex, torg, tdest);
7478#ifdef ANSI_DECLARATORS
7481void makevertexmap(m, b)
7487 struct otri triangleloop;
7491 printf(
" Constructing mapping from vertices to triangles.\n");
7493 traversalinit(&m->triangles);
7494 triangleloop.tri = triangletraverse(m);
7495 while (triangleloop.tri != (triangle *) NULL) {
7497 for (triangleloop.orient = 0; triangleloop.orient < 3;
7498 triangleloop.orient++) {
7499 org(triangleloop, triorg);
7500 setvertex2tri(triorg, encode(triangleloop));
7502 triangleloop.tri = triangletraverse(m);
7573#ifdef ANSI_DECLARATORS
7574enum locateresult preciselocate(
struct mesh *m,
struct behavior *b,
7575 vertex searchpoint,
struct otri *searchtri,
7576 int stopatsubsegment)
7578enum locateresult preciselocate(m, b, searchpoint, searchtri, stopatsubsegment)
7582struct otri *searchtri;
7583int stopatsubsegment;
7587 struct otri backtracktri;
7588 struct osub checkedge;
7589 vertex forg, fdest, fapex;
7590 REAL orgorient, destorient;
7595 if (b->verbose > 2) {
7596 printf(
" Searching for point (%.12g, %.12g).\n",
7597 searchpoint[0], searchpoint[1]);
7600 org(*searchtri, forg);
7601 dest(*searchtri, fdest);
7602 apex(*searchtri, fapex);
7604 if (b->verbose > 2) {
7605 printf(
" At (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
7606 forg[0], forg[1], fdest[0], fdest[1], fapex[0], fapex[1]);
7609 if ((fapex[0] == searchpoint[0]) && (fapex[1] == searchpoint[1])) {
7610 lprevself(*searchtri);
7615 destorient = counterclockwise(m, b, forg, fapex, searchpoint);
7618 orgorient = counterclockwise(m, b, fapex, fdest, searchpoint);
7619 if (destorient > 0.0) {
7620 if (orgorient > 0.0) {
7626 moveleft = (fapex[0] - searchpoint[0]) * (fdest[0] - forg[0]) +
7627 (fapex[1] - searchpoint[1]) * (fdest[1] - forg[1]) > 0.0;
7632 if (orgorient > 0.0) {
7637 if (destorient == 0.0) {
7638 lprevself(*searchtri);
7641 if (orgorient == 0.0) {
7642 lnextself(*searchtri);
7653 lprev(*searchtri, backtracktri);
7656 lnext(*searchtri, backtracktri);
7659 sym(backtracktri, *searchtri);
7661 if (m->checksegments && stopatsubsegment) {
7663 tspivot(backtracktri, checkedge);
7664 if (checkedge.ss != m->dummysub) {
7666 otricopy(backtracktri, *searchtri);
7671 if (searchtri->tri == m->dummytri) {
7673 otricopy(backtracktri, *searchtri);
7677 apex(*searchtri, fapex);
7717#ifdef ANSI_DECLARATORS
7718enum locateresult locate(
struct mesh *m,
struct behavior *b,
7719 vertex searchpoint,
struct otri *searchtri)
7721enum locateresult locate(m, b, searchpoint, searchtri)
7725struct otri *searchtri;
7731 struct otri sampletri;
7733 unsigned long alignptr;
7734 REAL searchdist, dist;
7736 long samplesperblock, totalsamplesleft, samplesleft;
7737 long population, totalpopulation;
7740 if (b->verbose > 2) {
7741 printf(
" Randomly sampling for a triangle near point (%.12g, %.12g).\n",
7742 searchpoint[0], searchpoint[1]);
7746 org(*searchtri, torg);
7747 searchdist = (searchpoint[0] - torg[0]) * (searchpoint[0] - torg[0]) +
7748 (searchpoint[1] - torg[1]) * (searchpoint[1] - torg[1]);
7749 if (b->verbose > 2) {
7750 printf(
" Boundary triangle has origin (%.12g, %.12g).\n",
7756 if (m->recenttri.tri != (triangle *) NULL) {
7757 if (!deadtri(m->recenttri.tri)) {
7758 org(m->recenttri, torg);
7759 if ((torg[0] == searchpoint[0]) && (torg[1] == searchpoint[1])) {
7760 otricopy(m->recenttri, *searchtri);
7763 dist = (searchpoint[0] - torg[0]) * (searchpoint[0] - torg[0]) +
7764 (searchpoint[1] - torg[1]) * (searchpoint[1] - torg[1]);
7765 if (dist < searchdist) {
7766 otricopy(m->recenttri, *searchtri);
7768 if (b->verbose > 2) {
7769 printf(
" Choosing recent triangle with origin (%.12g, %.12g).\n",
7780 while (SAMPLEFACTOR * m->samples * m->samples * m->samples <
7781 m->triangles.items) {
7789 samplesperblock = (m->samples * TRIPERBLOCK - 1) / m->triangles.maxitems + 1;
7792 samplesleft = (m->samples * m->triangles.itemsfirstblock - 1) /
7793 m->triangles.maxitems + 1;
7794 totalsamplesleft = m->samples;
7795 population = m->triangles.itemsfirstblock;
7796 totalpopulation = m->triangles.maxitems;
7797 sampleblock = m->triangles.firstblock;
7798 sampletri.orient = 0;
7799 while (totalsamplesleft > 0) {
7801 if (population > totalpopulation) {
7802 population = totalpopulation;
7805 alignptr = (
unsigned long) (sampleblock + 1);
7806 firsttri = (
char *) (alignptr +
7807 (
unsigned long) m->triangles.alignbytes -
7809 (
unsigned long) m->triangles.alignbytes));
7813 sampletri.tri = (triangle *) (firsttri +
7814 (randomnation((
unsigned int) population) *
7815 m->triangles.itembytes));
7816 if (!deadtri(sampletri.tri)) {
7817 org(sampletri, torg);
7818 dist = (searchpoint[0] - torg[0]) * (searchpoint[0] - torg[0]) +
7819 (searchpoint[1] - torg[1]) * (searchpoint[1] - torg[1]);
7820 if (dist < searchdist) {
7821 otricopy(sampletri, *searchtri);
7823 if (b->verbose > 2) {
7824 printf(
" Choosing triangle with origin (%.12g, %.12g).\n",
7832 }
while ((samplesleft > 0) && (totalsamplesleft > 0));
7834 if (totalsamplesleft > 0) {
7835 sampleblock = (VOID **) *sampleblock;
7836 samplesleft = samplesperblock;
7837 totalpopulation -= population;
7838 population = TRIPERBLOCK;
7843 org(*searchtri, torg);
7844 dest(*searchtri, tdest);
7846 if ((torg[0] == searchpoint[0]) && (torg[1] == searchpoint[1])) {
7849 if ((tdest[0] == searchpoint[0]) && (tdest[1] == searchpoint[1])) {
7850 lnextself(*searchtri);
7854 ahead = counterclockwise(m, b, torg, tdest, searchpoint);
7858 symself(*searchtri);
7859 }
else if (ahead == 0.0) {
7861 if (((torg[0] < searchpoint[0]) == (searchpoint[0] < tdest[0])) &&
7862 ((torg[1] < searchpoint[1]) == (searchpoint[1] < tdest[1]))) {
7866 return preciselocate(m, b, searchpoint, searchtri, 0);
7888#ifdef ANSI_DECLARATORS
7892void insertsubseg(m, b, tri, subsegmark)
7900 struct otri oppotri;
7901 struct osub newsubseg;
7902 vertex triorg, tridest;
7907 dest(*tri, tridest);
7909 if (vertexmark(triorg) == 0) {
7910 setvertexmark(triorg, subsegmark);
7912 if (vertexmark(tridest) == 0) {
7913 setvertexmark(tridest, subsegmark);
7916 tspivot(*tri, newsubseg);
7917 if (newsubseg.ss == m->dummysub) {
7919 makesubseg(m, &newsubseg);
7920 setsorg(newsubseg, tridest);
7921 setsdest(newsubseg, triorg);
7922 setsegorg(newsubseg, tridest);
7923 setsegdest(newsubseg, triorg);
7928 tsbond(*tri, newsubseg);
7930 ssymself(newsubseg);
7931 tsbond(oppotri, newsubseg);
7932 setmark(newsubseg, subsegmark);
7933 if (b->verbose > 2) {
7934 printf(
" Inserting new ");
7935 printsubseg(m, b, &newsubseg);
7938 if (mark(newsubseg) == 0) {
7939 setmark(newsubseg, subsegmark);
7992#ifdef ANSI_DECLARATORS
7995void flip(m, b, flipedge)
7998struct otri *flipedge;
8002 struct otri botleft, botright;
8003 struct otri topleft, topright;
8005 struct otri botlcasing, botrcasing;
8006 struct otri toplcasing, toprcasing;
8007 struct osub botlsubseg, botrsubseg;
8008 struct osub toplsubseg, toprsubseg;
8009 vertex leftvertex, rightvertex, botvertex;
8015 org(*flipedge, rightvertex);
8016 dest(*flipedge, leftvertex);
8017 apex(*flipedge, botvertex);
8018 sym(*flipedge, top);
8020 if (top.tri == m->dummytri) {
8021 printf(
"Internal error in flip(): Attempt to flip on boundary.\n");
8022 lnextself(*flipedge);
8025 if (m->checksegments) {
8026 tspivot(*flipedge, toplsubseg);
8027 if (toplsubseg.ss != m->dummysub) {
8028 printf(
"Internal error in flip(): Attempt to flip a segment.\n");
8029 lnextself(*flipedge);
8034 apex(top, farvertex);
8037 lprev(top, topleft);
8038 sym(topleft, toplcasing);
8039 lnext(top, topright);
8040 sym(topright, toprcasing);
8041 lnext(*flipedge, botleft);
8042 sym(botleft, botlcasing);
8043 lprev(*flipedge, botright);
8044 sym(botright, botrcasing);
8046 bond(topleft, botlcasing);
8047 bond(botleft, botrcasing);
8048 bond(botright, toprcasing);
8049 bond(topright, toplcasing);
8051 if (m->checksegments) {
8053 tspivot(topleft, toplsubseg);
8054 tspivot(botleft, botlsubseg);
8055 tspivot(botright, botrsubseg);
8056 tspivot(topright, toprsubseg);
8057 if (toplsubseg.ss == m->dummysub) {
8058 tsdissolve(topright);
8060 tsbond(topright, toplsubseg);
8062 if (botlsubseg.ss == m->dummysub) {
8063 tsdissolve(topleft);
8065 tsbond(topleft, botlsubseg);
8067 if (botrsubseg.ss == m->dummysub) {
8068 tsdissolve(botleft);
8070 tsbond(botleft, botrsubseg);
8072 if (toprsubseg.ss == m->dummysub) {
8073 tsdissolve(botright);
8075 tsbond(botright, toprsubseg);
8080 setorg(*flipedge, farvertex);
8081 setdest(*flipedge, botvertex);
8082 setapex(*flipedge, rightvertex);
8083 setorg(top, botvertex);
8084 setdest(top, farvertex);
8085 setapex(top, leftvertex);
8086 if (b->verbose > 2) {
8087 printf(
" Edge flip results in left ");
8088 printtriangle(m, b, &top);
8089 printf(
" and right ");
8090 printtriangle(m, b, flipedge);
8127#ifdef ANSI_DECLARATORS
8130void unflip(m, b, flipedge)
8133struct otri *flipedge;
8137 struct otri botleft, botright;
8138 struct otri topleft, topright;
8140 struct otri botlcasing, botrcasing;
8141 struct otri toplcasing, toprcasing;
8142 struct osub botlsubseg, botrsubseg;
8143 struct osub toplsubseg, toprsubseg;
8144 vertex leftvertex, rightvertex, botvertex;
8150 org(*flipedge, rightvertex);
8151 dest(*flipedge, leftvertex);
8152 apex(*flipedge, botvertex);
8153 sym(*flipedge, top);
8155 if (top.tri == m->dummytri) {
8156 printf(
"Internal error in unflip(): Attempt to flip on boundary.\n");
8157 lnextself(*flipedge);
8160 if (m->checksegments) {
8161 tspivot(*flipedge, toplsubseg);
8162 if (toplsubseg.ss != m->dummysub) {
8163 printf(
"Internal error in unflip(): Attempt to flip a subsegment.\n");
8164 lnextself(*flipedge);
8169 apex(top, farvertex);
8172 lprev(top, topleft);
8173 sym(topleft, toplcasing);
8174 lnext(top, topright);
8175 sym(topright, toprcasing);
8176 lnext(*flipedge, botleft);
8177 sym(botleft, botlcasing);
8178 lprev(*flipedge, botright);
8179 sym(botright, botrcasing);
8181 bond(topleft, toprcasing);
8182 bond(botleft, toplcasing);
8183 bond(botright, botlcasing);
8184 bond(topright, botrcasing);
8186 if (m->checksegments) {
8188 tspivot(topleft, toplsubseg);
8189 tspivot(botleft, botlsubseg);
8190 tspivot(botright, botrsubseg);
8191 tspivot(topright, toprsubseg);
8192 if (toplsubseg.ss == m->dummysub) {
8193 tsdissolve(botleft);
8195 tsbond(botleft, toplsubseg);
8197 if (botlsubseg.ss == m->dummysub) {
8198 tsdissolve(botright);
8200 tsbond(botright, botlsubseg);
8202 if (botrsubseg.ss == m->dummysub) {
8203 tsdissolve(topright);
8205 tsbond(topright, botrsubseg);
8207 if (toprsubseg.ss == m->dummysub) {
8208 tsdissolve(topleft);
8210 tsbond(topleft, toprsubseg);
8215 setorg(*flipedge, botvertex);
8216 setdest(*flipedge, farvertex);
8217 setapex(*flipedge, leftvertex);
8218 setorg(top, farvertex);
8219 setdest(top, botvertex);
8220 setapex(top, rightvertex);
8221 if (b->verbose > 2) {
8222 printf(
" Edge unflip results in left ");
8223 printtriangle(m, b, flipedge);
8224 printf(
" and right ");
8225 printtriangle(m, b, &top);
8276#ifdef ANSI_DECLARATORS
8277enum insertvertexresult insertvertex(
struct mesh *m,
struct behavior *b,
8278 vertex newvertex,
struct otri *searchtri,
8279 struct osub *splitseg,
8280 int segmentflaws,
int triflaws)
8282enum insertvertexresult insertvertex(m, b, newvertex, searchtri, splitseg,
8283 segmentflaws, triflaws)
8287struct otri *searchtri;
8288struct osub *splitseg;
8296 struct otri botleft, botright;
8297 struct otri topleft, topright;
8298 struct otri newbotleft, newbotright;
8299 struct otri newtopright;
8300 struct otri botlcasing, botrcasing;
8301 struct otri toplcasing={NULL, 0}, toprcasing={NULL, 0};
8302 struct otri testtri;
8303 struct osub botlsubseg, botrsubseg;
8304 struct osub toplsubseg, toprsubseg;
8305 struct osub brokensubseg;
8306 struct osub checksubseg;
8307 struct osub rightsubseg;
8308 struct osub newsubseg;
8312 vertex leftvertex, rightvertex, botvertex, topvertex, farvertex;
8313 vertex segmentorg, segmentdest;
8316 enum insertvertexresult success;
8317 enum locateresult intersect;
8325 if (b->verbose > 1) {
8326 printf(
" Inserting (%.12g, %.12g).\n", newvertex[0], newvertex[1]);
8329 if (splitseg == (
struct osub *) NULL) {
8332 if (searchtri->tri == m->dummytri) {
8334 horiz.tri = m->dummytri;
8338 intersect = locate(m, b, newvertex, &horiz);
8341 otricopy(*searchtri, horiz);
8342 intersect = preciselocate(m, b, newvertex, &horiz, 1);
8347 otricopy(*searchtri, horiz);
8351 if (intersect == ONVERTEX) {
8354 otricopy(horiz, *searchtri);
8355 otricopy(horiz, m->recenttri);
8356 return DUPLICATEVERTEX;
8358 if ((intersect == ONEDGE) || (intersect == OUTSIDE)) {
8360 if (m->checksegments && (splitseg == (
struct osub *) NULL)) {
8362 tspivot(horiz, brokensubseg);
8363 if (brokensubseg.ss != m->dummysub) {
8366 enq = b->nobisect != 2;
8367 if (enq && (b->nobisect == 1)) {
8370 sym(horiz, testtri);
8371 enq = testtri.tri != m->dummytri;
8375 encroached = (
struct badsubseg *) poolalloc(&m->badsubsegs);
8376 encroached->encsubseg = sencode(brokensubseg);
8377 sorg(brokensubseg, encroached->subsegorg);
8378 sdest(brokensubseg, encroached->subsegdest);
8379 if (b->verbose > 2) {
8381 " Queueing encroached subsegment (%.12g, %.12g) (%.12g, %.12g).\n",
8382 encroached->subsegorg[0], encroached->subsegorg[1],
8383 encroached->subsegdest[0], encroached->subsegdest[1]);
8389 otricopy(horiz, *searchtri);
8390 otricopy(horiz, m->recenttri);
8391 return VIOLATINGVERTEX;
8397 lprev(horiz, botright);
8398 sym(botright, botrcasing);
8399 sym(horiz, topright);
8401 mirrorflag = topright.tri != m->dummytri;
8403 lnextself(topright);
8404 sym(topright, toprcasing);
8405 maketriangle(m, b, &newtopright);
8410 maketriangle(m, b, &newbotright);
8413 org(horiz, rightvertex);
8414 dest(horiz, leftvertex);
8415 apex(horiz, botvertex);
8416 setorg(newbotright, botvertex);
8417 setdest(newbotright, rightvertex);
8418 setapex(newbotright, newvertex);
8419 setorg(horiz, newvertex);
8420 for (i = 0; i < m->eextras; i++) {
8422 setelemattribute(newbotright, i, elemattribute(botright, i));
8426 setareabound(newbotright, areabound(botright));
8429 dest(topright, topvertex);
8430 setorg(newtopright, rightvertex);
8431 setdest(newtopright, topvertex);
8432 setapex(newtopright, newvertex);
8433 setorg(topright, newvertex);
8434 for (i = 0; i < m->eextras; i++) {
8436 setelemattribute(newtopright, i, elemattribute(topright, i));
8440 setareabound(newtopright, areabound(topright));
8446 if (m->checksegments) {
8447 tspivot(botright, botrsubseg);
8448 if (botrsubseg.ss != m->dummysub) {
8449 tsdissolve(botright);
8450 tsbond(newbotright, botrsubseg);
8453 tspivot(topright, toprsubseg);
8454 if (toprsubseg.ss != m->dummysub) {
8455 tsdissolve(topright);
8456 tsbond(newtopright, toprsubseg);
8462 bond(newbotright, botrcasing);
8463 lprevself(newbotright);
8464 bond(newbotright, botright);
8465 lprevself(newbotright);
8467 bond(newtopright, toprcasing);
8468 lnextself(newtopright);
8469 bond(newtopright, topright);
8470 lnextself(newtopright);
8471 bond(newtopright, newbotright);
8474 if (splitseg != (
struct osub *) NULL) {
8476 setsdest(*splitseg, newvertex);
8477 segorg(*splitseg, segmentorg);
8478 segdest(*splitseg, segmentdest);
8479 ssymself(*splitseg);
8480 spivot(*splitseg, rightsubseg);
8481 insertsubseg(m, b, &newbotright, mark(*splitseg));
8482 tspivot(newbotright, newsubseg);
8483 setsegorg(newsubseg, segmentorg);
8484 setsegdest(newsubseg, segmentdest);
8485 sbond(*splitseg, newsubseg);
8486 ssymself(newsubseg);
8487 sbond(newsubseg, rightsubseg);
8488 ssymself(*splitseg);
8491 if (vertexmark(newvertex) == 0) {
8492 setvertexmark(newvertex, mark(*splitseg));
8496 if (m->checkquality) {
8497 poolrestart(&m->flipstackers);
8498 m->lastflip = (
struct flipstacker *) poolalloc(&m->flipstackers);
8499 m->lastflip->flippedtri = encode(horiz);
8500 m->lastflip->prevflip = (
struct flipstacker *) &insertvertex;
8504 if (counterclockwise(m, b, rightvertex, leftvertex, botvertex) < 0.0) {
8505 printf(
"Internal error in insertvertex():\n");
8507 " Clockwise triangle prior to edge vertex insertion (bottom).\n");
8510 if (counterclockwise(m, b, leftvertex, rightvertex, topvertex) < 0.0) {
8511 printf(
"Internal error in insertvertex():\n");
8512 printf(
" Clockwise triangle prior to edge vertex insertion (top).\n");
8514 if (counterclockwise(m, b, rightvertex, topvertex, newvertex) < 0.0) {
8515 printf(
"Internal error in insertvertex():\n");
8517 " Clockwise triangle after edge vertex insertion (top right).\n");
8519 if (counterclockwise(m, b, topvertex, leftvertex, newvertex) < 0.0) {
8520 printf(
"Internal error in insertvertex():\n");
8522 " Clockwise triangle after edge vertex insertion (top left).\n");
8525 if (counterclockwise(m, b, leftvertex, botvertex, newvertex) < 0.0) {
8526 printf(
"Internal error in insertvertex():\n");
8528 " Clockwise triangle after edge vertex insertion (bottom left).\n");
8530 if (counterclockwise(m, b, botvertex, rightvertex, newvertex) < 0.0) {
8531 printf(
"Internal error in insertvertex():\n");
8533 " Clockwise triangle after edge vertex insertion (bottom right).\n");
8536 if (b->verbose > 2) {
8537 printf(
" Updating bottom left ");
8538 printtriangle(m, b, &botright);
8540 printf(
" Updating top left ");
8541 printtriangle(m, b, &topright);
8542 printf(
" Creating top right ");
8543 printtriangle(m, b, &newtopright);
8545 printf(
" Creating bottom right ");
8546 printtriangle(m, b, &newbotright);
8554 lnext(horiz, botleft);
8555 lprev(horiz, botright);
8556 sym(botleft, botlcasing);
8557 sym(botright, botrcasing);
8558 maketriangle(m, b, &newbotleft);
8559 maketriangle(m, b, &newbotright);
8562 org(horiz, rightvertex);
8563 dest(horiz, leftvertex);
8564 apex(horiz, botvertex);
8565 setorg(newbotleft, leftvertex);
8566 setdest(newbotleft, botvertex);
8567 setapex(newbotleft, newvertex);
8568 setorg(newbotright, botvertex);
8569 setdest(newbotright, rightvertex);
8570 setapex(newbotright, newvertex);
8571 setapex(horiz, newvertex);
8572 for (i = 0; i < m->eextras; i++) {
8574 attrib = elemattribute(horiz, i);
8575 setelemattribute(newbotleft, i, attrib);
8576 setelemattribute(newbotright, i, attrib);
8580 area = areabound(horiz);
8581 setareabound(newbotleft, area);
8582 setareabound(newbotright, area);
8587 if (m->checksegments) {
8588 tspivot(botleft, botlsubseg);
8589 if (botlsubseg.ss != m->dummysub) {
8590 tsdissolve(botleft);
8591 tsbond(newbotleft, botlsubseg);
8593 tspivot(botright, botrsubseg);
8594 if (botrsubseg.ss != m->dummysub) {
8595 tsdissolve(botright);
8596 tsbond(newbotright, botrsubseg);
8601 bond(newbotleft, botlcasing);
8602 bond(newbotright, botrcasing);
8603 lnextself(newbotleft);
8604 lprevself(newbotright);
8605 bond(newbotleft, newbotright);
8606 lnextself(newbotleft);
8607 bond(botleft, newbotleft);
8608 lprevself(newbotright);
8609 bond(botright, newbotright);
8611 if (m->checkquality) {
8612 poolrestart(&m->flipstackers);
8613 m->lastflip = (
struct flipstacker *) poolalloc(&m->flipstackers);
8614 m->lastflip->flippedtri = encode(horiz);
8615 m->lastflip->prevflip = (
struct flipstacker *) NULL;
8619 if (counterclockwise(m, b, rightvertex, leftvertex, botvertex) < 0.0) {
8620 printf(
"Internal error in insertvertex():\n");
8621 printf(
" Clockwise triangle prior to vertex insertion.\n");
8623 if (counterclockwise(m, b, rightvertex, leftvertex, newvertex) < 0.0) {
8624 printf(
"Internal error in insertvertex():\n");
8625 printf(
" Clockwise triangle after vertex insertion (top).\n");
8627 if (counterclockwise(m, b, leftvertex, botvertex, newvertex) < 0.0) {
8628 printf(
"Internal error in insertvertex():\n");
8629 printf(
" Clockwise triangle after vertex insertion (left).\n");
8631 if (counterclockwise(m, b, botvertex, rightvertex, newvertex) < 0.0) {
8632 printf(
"Internal error in insertvertex():\n");
8633 printf(
" Clockwise triangle after vertex insertion (right).\n");
8636 if (b->verbose > 2) {
8637 printf(
" Updating top ");
8638 printtriangle(m, b, &horiz);
8639 printf(
" Creating left ");
8640 printtriangle(m, b, &newbotleft);
8641 printf(
" Creating right ");
8642 printtriangle(m, b, &newbotright);
8648 success = SUCCESSFULVERTEX;
8654 rightvertex = first;
8655 dest(horiz, leftvertex);
8661 if (m->checksegments) {
8663 tspivot(horiz, checksubseg);
8664 if (checksubseg.ss != m->dummysub) {
8670 if (checkseg4encroach(m, b, &checksubseg)) {
8671 success = ENCROACHINGVERTEX;
8681 if (top.tri == m->dummytri) {
8686 apex(top, farvertex);
8692 if ((leftvertex == m->infvertex1) || (leftvertex == m->infvertex2) ||
8693 (leftvertex == m->infvertex3)) {
8698 doflip = counterclockwise(m, b, newvertex, rightvertex, farvertex)
8700 }
else if ((rightvertex == m->infvertex1) ||
8701 (rightvertex == m->infvertex2) ||
8702 (rightvertex == m->infvertex3)) {
8707 doflip = counterclockwise(m, b, farvertex, leftvertex, newvertex)
8709 }
else if ((farvertex == m->infvertex1) ||
8710 (farvertex == m->infvertex2) ||
8711 (farvertex == m->infvertex3)) {
8717 doflip = incircle(m, b, leftvertex, newvertex, rightvertex,
8724 lprev(top, topleft);
8725 sym(topleft, toplcasing);
8726 lnext(top, topright);
8727 sym(topright, toprcasing);
8728 lnext(horiz, botleft);
8729 sym(botleft, botlcasing);
8730 lprev(horiz, botright);
8731 sym(botright, botrcasing);
8733 bond(topleft, botlcasing);
8734 bond(botleft, botrcasing);
8735 bond(botright, toprcasing);
8736 bond(topright, toplcasing);
8737 if (m->checksegments) {
8739 tspivot(topleft, toplsubseg);
8740 tspivot(botleft, botlsubseg);
8741 tspivot(botright, botrsubseg);
8742 tspivot(topright, toprsubseg);
8743 if (toplsubseg.ss == m->dummysub) {
8744 tsdissolve(topright);
8746 tsbond(topright, toplsubseg);
8748 if (botlsubseg.ss == m->dummysub) {
8749 tsdissolve(topleft);
8751 tsbond(topleft, botlsubseg);
8753 if (botrsubseg.ss == m->dummysub) {
8754 tsdissolve(botleft);
8756 tsbond(botleft, botrsubseg);
8758 if (toprsubseg.ss == m->dummysub) {
8759 tsdissolve(botright);
8761 tsbond(botright, toprsubseg);
8765 setorg(horiz, farvertex);
8766 setdest(horiz, newvertex);
8767 setapex(horiz, rightvertex);
8768 setorg(top, newvertex);
8769 setdest(top, farvertex);
8770 setapex(top, leftvertex);
8771 for (i = 0; i < m->eextras; i++) {
8773 attrib = 0.5 * (elemattribute(top, i) + elemattribute(horiz, i));
8774 setelemattribute(top, i, attrib);
8775 setelemattribute(horiz, i, attrib);
8778 if ((areabound(top) <= 0.0) || (areabound(horiz) <= 0.0)) {
8784 area = 0.5 * (areabound(top) + areabound(horiz));
8786 setareabound(top, area);
8787 setareabound(horiz, area);
8790 if (m->checkquality) {
8791 newflip = (
struct flipstacker *) poolalloc(&m->flipstackers);
8792 newflip->flippedtri = encode(horiz);
8793 newflip->prevflip = m->lastflip;
8794 m->lastflip = newflip;
8798 if (newvertex != (vertex) NULL) {
8799 if (counterclockwise(m, b, leftvertex, newvertex, rightvertex) <
8801 printf(
"Internal error in insertvertex():\n");
8802 printf(
" Clockwise triangle prior to edge flip (bottom).\n");
8814 if (counterclockwise(m, b, farvertex, leftvertex, newvertex) <
8816 printf(
"Internal error in insertvertex():\n");
8817 printf(
" Clockwise triangle after edge flip (left).\n");
8819 if (counterclockwise(m, b, newvertex, rightvertex, farvertex) <
8821 printf(
"Internal error in insertvertex():\n");
8822 printf(
" Clockwise triangle after edge flip (right).\n");
8826 if (b->verbose > 2) {
8827 printf(
" Edge flip results in left ");
8829 printtriangle(m, b, &topleft);
8830 printf(
" and right ");
8831 printtriangle(m, b, &horiz);
8837 leftvertex = farvertex;
8846 testtriangle(m, b, &horiz);
8851 sym(horiz, testtri);
8855 if ((leftvertex == first) || (testtri.tri == m->dummytri)) {
8857 lnext(horiz, *searchtri);
8858 lnext(horiz, m->recenttri);
8862 lnext(testtri, horiz);
8863 rightvertex = leftvertex;
8864 dest(horiz, leftvertex);
8933#ifdef ANSI_DECLARATORS
8934void triangulatepolygon(
struct mesh *m,
struct behavior *b,
8935 struct otri *firstedge,
struct otri *lastedge,
8936 int edgecount,
int doflip,
int triflaws)
8938void triangulatepolygon(m, b, firstedge, lastedge, edgecount, doflip, triflaws)
8941struct otri *firstedge;
8942struct otri *lastedge;
8949 struct otri testtri;
8950 struct otri besttri;
8951 struct otri tempedge;
8952 vertex leftbasevertex, rightbasevertex;
8960 apex(*lastedge, leftbasevertex);
8961 dest(*firstedge, rightbasevertex);
8962 if (b->verbose > 2) {
8963 printf(
" Triangulating interior polygon at edge\n");
8964 printf(
" (%.12g, %.12g) (%.12g, %.12g)\n", leftbasevertex[0],
8965 leftbasevertex[1], rightbasevertex[0], rightbasevertex[1]);
8968 onext(*firstedge, besttri);
8969 dest(besttri, bestvertex);
8970 otricopy(besttri, testtri);
8972 for (i = 2; i <= edgecount - 2; i++) {
8974 dest(testtri, testvertex);
8976 if (incircle(m, b, leftbasevertex, rightbasevertex, bestvertex,
8977 testvertex) > 0.0) {
8978 otricopy(testtri, besttri);
8979 bestvertex = testvertex;
8983 if (b->verbose > 2) {
8984 printf(
" Connecting edge to (%.12g, %.12g)\n", bestvertex[0],
8987 if (bestnumber > 1) {
8989 oprev(besttri, tempedge);
8990 triangulatepolygon(m, b, firstedge, &tempedge, bestnumber + 1, 1,
8993 if (bestnumber < edgecount - 2) {
8995 sym(besttri, tempedge);
8996 triangulatepolygon(m, b, &besttri, lastedge, edgecount - bestnumber, 1,
8999 sym(tempedge, besttri);
9003 flip(m, b, &besttri);
9007 sym(besttri, testtri);
9008 testtriangle(m, b, &testtri);
9013 otricopy(besttri, *lastedge);
9032#ifdef ANSI_DECLARATORS
9035void deletevertex(m, b, deltri)
9042 struct otri countingtri;
9043 struct otri firstedge, lastedge;
9044 struct otri deltriright;
9045 struct otri lefttri, righttri;
9046 struct otri leftcasing, rightcasing;
9047 struct osub leftsubseg, rightsubseg;
9054 org(*deltri, delvertex);
9055 if (b->verbose > 1) {
9056 printf(
" Deleting (%.12g, %.12g).\n", delvertex[0], delvertex[1]);
9058 vertexdealloc(m, delvertex);
9061 onext(*deltri, countingtri);
9063 while (!otriequal(*deltri, countingtri)) {
9065 if (countingtri.tri == m->dummytri) {
9066 printf(
"Internal error in deletevertex():\n");
9067 printf(
" Attempt to delete boundary vertex.\n");
9072 onextself(countingtri);
9076 if (edgecount < 3) {
9077 printf(
"Internal error in deletevertex():\n Vertex has degree %d.\n",
9082 if (edgecount > 3) {
9086 onext(*deltri, firstedge);
9087 oprev(*deltri, lastedge);
9088 triangulatepolygon(m, b, &firstedge, &lastedge, edgecount, 0,
9092 lprev(*deltri, deltriright);
9093 dnext(*deltri, lefttri);
9094 sym(lefttri, leftcasing);
9095 oprev(deltriright, righttri);
9096 sym(righttri, rightcasing);
9097 bond(*deltri, leftcasing);
9098 bond(deltriright, rightcasing);
9099 tspivot(lefttri, leftsubseg);
9100 if (leftsubseg.ss != m->dummysub) {
9101 tsbond(*deltri, leftsubseg);
9103 tspivot(righttri, rightsubseg);
9104 if (rightsubseg.ss != m->dummysub) {
9105 tsbond(deltriright, rightsubseg);
9109 org(lefttri, neworg);
9110 setorg(*deltri, neworg);
9112 testtriangle(m, b, deltri);
9116 triangledealloc(m, lefttri.tri);
9117 triangledealloc(m, righttri.tri);
9135#ifdef ANSI_DECLARATORS
9138void undovertex(m, b)
9144 struct otri fliptri;
9145 struct otri botleft, botright, topright;
9146 struct otri botlcasing, botrcasing, toprcasing;
9147 struct otri gluetri;
9148 struct osub botlsubseg, botrsubseg, toprsubseg;
9149 vertex botvertex, rightvertex;
9155 while (m->lastflip != (
struct flipstacker *) NULL) {
9157 decode(m->lastflip->flippedtri, fliptri);
9163 if (m->lastflip->prevflip == (
struct flipstacker *) NULL) {
9166 dprev(fliptri, botleft);
9168 onext(fliptri, botright);
9169 lprevself(botright);
9170 sym(botleft, botlcasing);
9171 sym(botright, botrcasing);
9172 dest(botleft, botvertex);
9174 setapex(fliptri, botvertex);
9176 bond(fliptri, botlcasing);
9177 tspivot(botleft, botlsubseg);
9178 tsbond(fliptri, botlsubseg);
9180 bond(fliptri, botrcasing);
9181 tspivot(botright, botrsubseg);
9182 tsbond(fliptri, botrsubseg);
9185 triangledealloc(m, botleft.tri);
9186 triangledealloc(m, botright.tri);
9187 }
else if (m->lastflip->prevflip == (
struct flipstacker *) &insertvertex) {
9190 lprev(fliptri, gluetri);
9191 sym(gluetri, botright);
9192 lnextself(botright);
9193 sym(botright, botrcasing);
9194 dest(botright, rightvertex);
9196 setorg(fliptri, rightvertex);
9197 bond(gluetri, botrcasing);
9198 tspivot(botright, botrsubseg);
9199 tsbond(gluetri, botrsubseg);
9202 triangledealloc(m, botright.tri);
9204 sym(fliptri, gluetri);
9205 if (gluetri.tri != m->dummytri) {
9207 dnext(gluetri, topright);
9208 sym(topright, toprcasing);
9210 setorg(gluetri, rightvertex);
9211 bond(gluetri, toprcasing);
9212 tspivot(topright, toprsubseg);
9213 tsbond(gluetri, toprsubseg);
9216 triangledealloc(m, topright.tri);
9220 m->lastflip->prevflip = (
struct flipstacker *) NULL;
9223 unflip(m, b, &fliptri);
9227 m->lastflip = m->lastflip->prevflip;
9282#ifdef ANSI_DECLARATORS
9283void vertexsort(vertex *sortarray,
int arraysize)
9285void vertexsort(sortarray, arraysize)
9293 REAL pivotx, pivoty;
9296 if (arraysize == 2) {
9298 if ((sortarray[0][0] > sortarray[1][0]) ||
9299 ((sortarray[0][0] == sortarray[1][0]) &&
9300 (sortarray[0][1] > sortarray[1][1]))) {
9301 temp = sortarray[1];
9302 sortarray[1] = sortarray[0];
9303 sortarray[0] = temp;
9308 pivot = (int) randomnation((
unsigned int) arraysize);
9309 pivotx = sortarray[pivot][0];
9310 pivoty = sortarray[pivot][1];
9314 while (left < right) {
9318 }
while ((left <= right) && ((sortarray[left][0] < pivotx) ||
9319 ((sortarray[left][0] == pivotx) &&
9320 (sortarray[left][1] < pivoty))));
9324 }
while ((left <= right) && ((sortarray[right][0] > pivotx) ||
9325 ((sortarray[right][0] == pivotx) &&
9326 (sortarray[right][1] > pivoty))));
9329 temp = sortarray[left];
9330 sortarray[left] = sortarray[right];
9331 sortarray[right] = temp;
9336 vertexsort(sortarray, left);
9338 if (right < arraysize - 2) {
9340 vertexsort(&sortarray[right + 1], arraysize - right - 1);
9356#ifdef ANSI_DECLARATORS
9357void vertexmedian(vertex *sortarray,
int arraysize,
int median,
int axis)
9359void vertexmedian(sortarray, arraysize, median, axis)
9369 REAL pivot1, pivot2;
9372 if (arraysize == 2) {
9374 if ((sortarray[0][axis] > sortarray[1][axis]) ||
9375 ((sortarray[0][axis] == sortarray[1][axis]) &&
9376 (sortarray[0][1 - axis] > sortarray[1][1 - axis]))) {
9377 temp = sortarray[1];
9378 sortarray[1] = sortarray[0];
9379 sortarray[0] = temp;
9384 pivot = (int) randomnation((
unsigned int) arraysize);
9385 pivot1 = sortarray[pivot][axis];
9386 pivot2 = sortarray[pivot][1 - axis];
9390 while (left < right) {
9394 }
while ((left <= right) && ((sortarray[left][axis] < pivot1) ||
9395 ((sortarray[left][axis] == pivot1) &&
9396 (sortarray[left][1 - axis] < pivot2))));
9400 }
while ((left <= right) && ((sortarray[right][axis] > pivot1) ||
9401 ((sortarray[right][axis] == pivot1) &&
9402 (sortarray[right][1 - axis] > pivot2))));
9405 temp = sortarray[left];
9406 sortarray[left] = sortarray[right];
9407 sortarray[right] = temp;
9412 if (left > median) {
9414 vertexmedian(sortarray, left, median, axis);
9416 if (right < median - 1) {
9418 vertexmedian(&sortarray[right + 1], arraysize - right - 1,
9419 median - right - 1, axis);
9434#ifdef ANSI_DECLARATORS
9435void alternateaxes(vertex *sortarray,
int arraysize,
int axis)
9437void alternateaxes(sortarray, arraysize, axis)
9446 divider = arraysize >> 1;
9447 if (arraysize <= 3) {
9453 vertexmedian(sortarray, arraysize, divider, axis);
9455 if (arraysize - divider >= 2) {
9457 alternateaxes(sortarray, divider, 1 - axis);
9459 alternateaxes(&sortarray[divider], arraysize - divider, 1 - axis);
9498#ifdef ANSI_DECLARATORS
9500 struct otri *innerleft,
struct otri *innerright,
9501 struct otri *farright,
int axis)
9503void mergehulls(m, b, farleft, innerleft, innerright, farright, axis)
9506struct otri *farleft;
9507struct otri *innerleft;
9508struct otri *innerright;
9509struct otri *farright;
9514 struct otri leftcand, rightcand;
9515 struct otri baseedge;
9516 struct otri nextedge;
9517 struct otri sidecasing, topcasing, outercasing;
9518 struct otri checkedge;
9519 vertex innerleftdest;
9520 vertex innerrightorg;
9521 vertex innerleftapex, innerrightapex;
9522 vertex farleftpt, farrightpt;
9523 vertex farleftapex, farrightapex;
9524 vertex lowerleft, lowerright;
9525 vertex upperleft, upperright;
9530 int leftfinished, rightfinished;
9533 dest(*innerleft, innerleftdest);
9534 apex(*innerleft, innerleftapex);
9535 org(*innerright, innerrightorg);
9536 apex(*innerright, innerrightapex);
9538 if (b->dwyer && (axis == 1)) {
9539 org(*farleft, farleftpt);
9540 apex(*farleft, farleftapex);
9541 dest(*farright, farrightpt);
9542 apex(*farright, farrightapex);
9546 while (farleftapex[1] < farleftpt[1]) {
9547 lnextself(*farleft);
9549 farleftpt = farleftapex;
9550 apex(*farleft, farleftapex);
9552 sym(*innerleft, checkedge);
9553 apex(checkedge, checkvertex);
9554 while (checkvertex[1] > innerleftdest[1]) {
9555 lnext(checkedge, *innerleft);
9556 innerleftapex = innerleftdest;
9557 innerleftdest = checkvertex;
9558 sym(*innerleft, checkedge);
9559 apex(checkedge, checkvertex);
9561 while (innerrightapex[1] < innerrightorg[1]) {
9562 lnextself(*innerright);
9563 symself(*innerright);
9564 innerrightorg = innerrightapex;
9565 apex(*innerright, innerrightapex);
9567 sym(*farright, checkedge);
9568 apex(checkedge, checkvertex);
9569 while (checkvertex[1] > farrightpt[1]) {
9570 lnext(checkedge, *farright);
9571 farrightapex = farrightpt;
9572 farrightpt = checkvertex;
9573 sym(*farright, checkedge);
9574 apex(checkedge, checkvertex);
9581 if (counterclockwise(m, b, innerleftdest, innerleftapex, innerrightorg) >
9583 lprevself(*innerleft);
9584 symself(*innerleft);
9585 innerleftdest = innerleftapex;
9586 apex(*innerleft, innerleftapex);
9590 if (counterclockwise(m, b, innerrightapex, innerrightorg, innerleftdest) >
9592 lnextself(*innerright);
9593 symself(*innerright);
9594 innerrightorg = innerrightapex;
9595 apex(*innerright, innerrightapex);
9598 }
while (changemade);
9600 sym(*innerleft, leftcand);
9601 sym(*innerright, rightcand);
9603 maketriangle(m, b, &baseedge);
9605 bond(baseedge, *innerleft);
9606 lnextself(baseedge);
9607 bond(baseedge, *innerright);
9608 lnextself(baseedge);
9609 setorg(baseedge, innerrightorg);
9610 setdest(baseedge, innerleftdest);
9612 if (b->verbose > 2) {
9613 printf(
" Creating base bounding ");
9614 printtriangle(m, b, &baseedge);
9617 org(*farleft, farleftpt);
9618 if (innerleftdest == farleftpt) {
9619 lnext(baseedge, *farleft);
9621 dest(*farright, farrightpt);
9622 if (innerrightorg == farrightpt) {
9623 lprev(baseedge, *farright);
9626 lowerleft = innerleftdest;
9627 lowerright = innerrightorg;
9629 apex(leftcand, upperleft);
9630 apex(rightcand, upperright);
9637 leftfinished = counterclockwise(m, b, upperleft, lowerleft, lowerright) <=
9639 rightfinished = counterclockwise(m, b, upperright, lowerleft, lowerright)
9641 if (leftfinished && rightfinished) {
9643 maketriangle(m, b, &nextedge);
9644 setorg(nextedge, lowerleft);
9645 setdest(nextedge, lowerright);
9648 bond(nextedge, baseedge);
9649 lnextself(nextedge);
9650 bond(nextedge, rightcand);
9651 lnextself(nextedge);
9652 bond(nextedge, leftcand);
9653 if (b->verbose > 2) {
9654 printf(
" Creating top bounding ");
9655 printtriangle(m, b, &nextedge);
9658 if (b->dwyer && (axis == 1)) {
9659 org(*farleft, farleftpt);
9660 apex(*farleft, farleftapex);
9661 dest(*farright, farrightpt);
9662 apex(*farright, farrightapex);
9663 sym(*farleft, checkedge);
9664 apex(checkedge, checkvertex);
9668 while (checkvertex[0] < farleftpt[0]) {
9669 lprev(checkedge, *farleft);
9670 farleftapex = farleftpt;
9671 farleftpt = checkvertex;
9672 sym(*farleft, checkedge);
9673 apex(checkedge, checkvertex);
9675 while (farrightapex[0] > farrightpt[0]) {
9676 lprevself(*farright);
9678 farrightpt = farrightapex;
9679 apex(*farright, farrightapex);
9685 if (!leftfinished) {
9687 lprev(leftcand, nextedge);
9689 apex(nextedge, nextapex);
9692 if (nextapex != (vertex) NULL) {
9694 badedge = incircle(m, b, lowerleft, lowerright, upperleft, nextapex) >
9699 lnextself(nextedge);
9700 sym(nextedge, topcasing);
9701 lnextself(nextedge);
9702 sym(nextedge, sidecasing);
9703 bond(nextedge, topcasing);
9704 bond(leftcand, sidecasing);
9705 lnextself(leftcand);
9706 sym(leftcand, outercasing);
9707 lprevself(nextedge);
9708 bond(nextedge, outercasing);
9710 setorg(leftcand, lowerleft);
9711 setdest(leftcand, NULL);
9712 setapex(leftcand, nextapex);
9713 setorg(nextedge, NULL);
9714 setdest(nextedge, upperleft);
9715 setapex(nextedge, nextapex);
9717 upperleft = nextapex;
9719 otricopy(sidecasing, nextedge);
9720 apex(nextedge, nextapex);
9721 if (nextapex != (vertex) NULL) {
9723 badedge = incircle(m, b, lowerleft, lowerright, upperleft,
9733 if (!rightfinished) {
9735 lnext(rightcand, nextedge);
9737 apex(nextedge, nextapex);
9740 if (nextapex != (vertex) NULL) {
9742 badedge = incircle(m, b, lowerleft, lowerright, upperright, nextapex) >
9747 lprevself(nextedge);
9748 sym(nextedge, topcasing);
9749 lprevself(nextedge);
9750 sym(nextedge, sidecasing);
9751 bond(nextedge, topcasing);
9752 bond(rightcand, sidecasing);
9753 lprevself(rightcand);
9754 sym(rightcand, outercasing);
9755 lnextself(nextedge);
9756 bond(nextedge, outercasing);
9758 setorg(rightcand, NULL);
9759 setdest(rightcand, lowerright);
9760 setapex(rightcand, nextapex);
9761 setorg(nextedge, upperright);
9762 setdest(nextedge, NULL);
9763 setapex(nextedge, nextapex);
9765 upperright = nextapex;
9767 otricopy(sidecasing, nextedge);
9768 apex(nextedge, nextapex);
9769 if (nextapex != (vertex) NULL) {
9771 badedge = incircle(m, b, lowerleft, lowerright, upperright,
9780 if (leftfinished || (!rightfinished &&
9781 (incircle(m, b, upperleft, lowerleft, lowerright, upperright) >
9785 bond(baseedge, rightcand);
9786 lprev(rightcand, baseedge);
9787 setdest(baseedge, lowerleft);
9788 lowerright = upperright;
9789 sym(baseedge, rightcand);
9790 apex(rightcand, upperright);
9794 bond(baseedge, leftcand);
9795 lnext(leftcand, baseedge);
9796 setorg(baseedge, lowerright);
9797 lowerleft = upperleft;
9798 sym(baseedge, leftcand);
9799 apex(leftcand, upperleft);
9801 if (b->verbose > 2) {
9802 printf(
" Connecting ");
9803 printtriangle(m, b, &baseedge);
9825#ifdef ANSI_DECLARATORS
9826void divconqrecurse(
struct mesh *m,
struct behavior *b, vertex *sortarray,
9827 int vertices,
int axis,
9828 struct otri *farleft,
struct otri *farright)
9830void divconqrecurse(m, b, sortarray, vertices, axis, farleft, farright)
9836struct otri *farleft;
9837struct otri *farright;
9841 struct otri midtri, tri1, tri2, tri3;
9842 struct otri innerleft, innerright;
9846 if (b->verbose > 2) {
9847 printf(
" Triangulating %d vertices.\n", vertices);
9849 if (vertices == 2) {
9852 maketriangle(m, b, farleft);
9853 setorg(*farleft, sortarray[0]);
9854 setdest(*farleft, sortarray[1]);
9856 maketriangle(m, b, farright);
9857 setorg(*farright, sortarray[1]);
9858 setdest(*farright, sortarray[0]);
9860 bond(*farleft, *farright);
9861 lprevself(*farleft);
9862 lnextself(*farright);
9863 bond(*farleft, *farright);
9864 lprevself(*farleft);
9865 lnextself(*farright);
9866 bond(*farleft, *farright);
9867 if (b->verbose > 2) {
9868 printf(
" Creating ");
9869 printtriangle(m, b, farleft);
9870 printf(
" Creating ");
9871 printtriangle(m, b, farright);
9874 lprev(*farright, *farleft);
9876 }
else if (vertices == 3) {
9880 maketriangle(m, b, &midtri);
9881 maketriangle(m, b, &tri1);
9882 maketriangle(m, b, &tri2);
9883 maketriangle(m, b, &tri3);
9884 area = counterclockwise(m, b, sortarray[0], sortarray[1], sortarray[2]);
9887 setorg(midtri, sortarray[0]);
9888 setdest(midtri, sortarray[1]);
9889 setorg(tri1, sortarray[1]);
9890 setdest(tri1, sortarray[0]);
9891 setorg(tri2, sortarray[2]);
9892 setdest(tri2, sortarray[1]);
9893 setorg(tri3, sortarray[1]);
9894 setdest(tri3, sortarray[2]);
9911 otricopy(tri1, *farleft);
9913 otricopy(tri2, *farright);
9917 setorg(midtri, sortarray[0]);
9918 setdest(tri1, sortarray[0]);
9919 setorg(tri3, sortarray[0]);
9923 setdest(midtri, sortarray[1]);
9924 setorg(tri1, sortarray[1]);
9925 setdest(tri2, sortarray[1]);
9926 setapex(midtri, sortarray[2]);
9927 setorg(tri2, sortarray[2]);
9928 setdest(tri3, sortarray[2]);
9931 setdest(midtri, sortarray[2]);
9932 setorg(tri1, sortarray[2]);
9933 setdest(tri2, sortarray[2]);
9934 setapex(midtri, sortarray[1]);
9935 setorg(tri2, sortarray[1]);
9936 setdest(tri3, sortarray[1]);
9954 otricopy(tri1, *farleft);
9957 otricopy(tri2, *farright);
9959 lnext(*farleft, *farright);
9962 if (b->verbose > 2) {
9963 printf(
" Creating ");
9964 printtriangle(m, b, &midtri);
9965 printf(
" Creating ");
9966 printtriangle(m, b, &tri1);
9967 printf(
" Creating ");
9968 printtriangle(m, b, &tri2);
9969 printf(
" Creating ");
9970 printtriangle(m, b, &tri3);
9975 divider = vertices >> 1;
9977 divconqrecurse(m, b, sortarray, divider, 1 - axis, farleft, &innerleft);
9978 divconqrecurse(m, b, &sortarray[divider], vertices - divider, 1 - axis,
9979 &innerright, farright);
9980 if (b->verbose > 1) {
9981 printf(
" Joining triangulations with %d and %d vertices.\n", divider,
9982 vertices - divider);
9985 mergehulls(m, b, farleft, &innerleft, &innerright, farright, axis);
9989#ifdef ANSI_DECLARATORS
9990long removeghosts(
struct mesh *m,
struct behavior *b,
struct otri *startghost)
9992long removeghosts(m, b, startghost)
9995struct otri *startghost;
9999 struct otri searchedge;
10000 struct otri dissolveedge;
10001 struct otri deadtriangle;
10007 printf(
" Removing ghost triangles.\n");
10010 lprev(*startghost, searchedge);
10011 symself(searchedge);
10012 m->dummytri[0] = encode(searchedge);
10014 otricopy(*startghost, dissolveedge);
10018 lnext(dissolveedge, deadtriangle);
10019 lprevself(dissolveedge);
10020 symself(dissolveedge);
10025 if (dissolveedge.tri != m->dummytri) {
10026 org(dissolveedge, markorg);
10027 if (vertexmark(markorg) == 0) {
10028 setvertexmark(markorg, 1);
10033 dissolve(dissolveedge);
10035 sym(deadtriangle, dissolveedge);
10037 triangledealloc(m, deadtriangle.tri);
10038 }
while (!otriequal(dissolveedge, *startghost));
10052#ifdef ANSI_DECLARATORS
10053long divconqdelaunay(
struct mesh *m,
struct behavior *b)
10055long divconqdelaunay(m, b)
10062 struct otri hullleft, hullright;
10067 printf(
" Sorting vertices.\n");
10071 sortarray = (vertex *) trimalloc(m->invertices * (
int)
sizeof(vertex));
10072 traversalinit(&m->vertices);
10073 for (i = 0; i < m->invertices; i++) {
10074 sortarray[i] = vertextraverse(m);
10077 vertexsort(sortarray, m->invertices);
10080 for (j = 1; j < m->invertices; j++) {
10081 if ((sortarray[i][0] == sortarray[j][0])
10082 && (sortarray[i][1] == sortarray[j][1])) {
10085"Warning: A duplicate vertex at (%.12g, %.12g) appeared and was ignored.\n",
10086 sortarray[j][0], sortarray[j][1]);
10088 setvertextype(sortarray[j], UNDEADVERTEX);
10092 sortarray[i] = sortarray[j];
10099 if (i - divider >= 2) {
10100 if (divider >= 2) {
10101 alternateaxes(sortarray, divider, 1);
10103 alternateaxes(&sortarray[divider], i - divider, 1);
10108 printf(
" Forming triangulation.\n");
10112 divconqrecurse(m, b, sortarray, i, 0, &hullleft, &hullright);
10113 trifree((VOID *) sortarray);
10115 return removeghosts(m, b, &hullleft);
10139#ifdef ANSI_DECLARATORS
10142void boundingbox(m, b)
10148 struct otri inftri;
10152 printf(
" Creating triangular bounding box.\n");
10155 width = m->xmax - m->xmin;
10156 if (m->ymax - m->ymin > width) {
10157 width = m->ymax - m->ymin;
10159 if (width == 0.0) {
10163 m->infvertex1 = (vertex) trimalloc(m->vertices.itembytes);
10164 m->infvertex2 = (vertex) trimalloc(m->vertices.itembytes);
10165 m->infvertex3 = (vertex) trimalloc(m->vertices.itembytes);
10166 m->infvertex1[0] = m->xmin - 50.0 * width;
10167 m->infvertex1[1] = m->ymin - 40.0 * width;
10168 m->infvertex2[0] = m->xmax + 50.0 * width;
10169 m->infvertex2[1] = m->ymin - 40.0 * width;
10170 m->infvertex3[0] = 0.5 * (m->xmin + m->xmax);
10171 m->infvertex3[1] = m->ymax + 60.0 * width;
10174 maketriangle(m, b, &inftri);
10175 setorg(inftri, m->infvertex1);
10176 setdest(inftri, m->infvertex2);
10177 setapex(inftri, m->infvertex3);
10180 m->dummytri[0] = (triangle) inftri.tri;
10181 if (b->verbose > 2) {
10182 printf(
" Creating ");
10183 printtriangle(m, b, &inftri);
10205#ifdef ANSI_DECLARATORS
10208long removebox(m, b)
10214 struct otri deadtriangle;
10215 struct otri searchedge;
10216 struct otri checkedge;
10217 struct otri nextedge, finaledge, dissolveedge;
10223 printf(
" Removing triangular bounding box.\n");
10226 nextedge.tri = m->dummytri;
10227 nextedge.orient = 0;
10230 lprev(nextedge, finaledge);
10231 lnextself(nextedge);
10235 lprev(nextedge, searchedge);
10236 symself(searchedge);
10239 lnext(nextedge, checkedge);
10240 symself(checkedge);
10241 if (checkedge.tri == m->dummytri) {
10245 lprevself(searchedge);
10246 symself(searchedge);
10250 m->dummytri[0] = encode(searchedge);
10252 while (!otriequal(nextedge, finaledge)) {
10254 lprev(nextedge, dissolveedge);
10255 symself(dissolveedge);
10263 if (dissolveedge.tri != m->dummytri) {
10264 org(dissolveedge, markorg);
10265 if (vertexmark(markorg) == 0) {
10266 setvertexmark(markorg, 1);
10271 dissolve(dissolveedge);
10272 lnext(nextedge, deadtriangle);
10273 sym(deadtriangle, nextedge);
10275 triangledealloc(m, deadtriangle.tri);
10277 if (nextedge.tri == m->dummytri) {
10279 otricopy(dissolveedge, nextedge);
10282 triangledealloc(m, finaledge.tri);
10284 trifree((VOID *) m->infvertex1);
10285 trifree((VOID *) m->infvertex2);
10286 trifree((VOID *) m->infvertex3);
10304#ifdef ANSI_DECLARATORS
10305long incrementaldelaunay(
struct mesh *m,
struct behavior *b)
10307long incrementaldelaunay(m, b)
10313 struct otri starttri;
10319 printf(
" Incrementally inserting vertices.\n");
10321 traversalinit(&m->vertices);
10322 vertexloop = vertextraverse(m);
10323 while (vertexloop != (vertex) NULL) {
10324 starttri.tri = m->dummytri;
10325 if (insertvertex(m, b, vertexloop, &starttri, (
struct osub *) NULL, 0, 0)
10326 == DUPLICATEVERTEX) {
10329"Warning: A duplicate vertex at (%.12g, %.12g) appeared and was ignored.\n",
10330 vertexloop[0], vertexloop[1]);
10332 setvertextype(vertexloop, UNDEADVERTEX);
10335 vertexloop = vertextraverse(m);
10338 return removebox(m, b);
10353#ifdef ANSI_DECLARATORS
10354void eventheapinsert(
struct event **heap,
int heapsize,
struct event *newevent)
10356void eventheapinsert(heap, heapsize, newevent)
10357struct event **heap;
10359struct event *newevent;
10363 REAL eventx, eventy;
10368 eventx = newevent->xkey;
10369 eventy = newevent->ykey;
10370 eventnum = heapsize;
10371 notdone = eventnum > 0;
10373 parent = (eventnum - 1) >> 1;
10374 if ((heap[parent]->ykey < eventy) ||
10375 ((heap[parent]->ykey == eventy)
10376 && (heap[parent]->xkey <= eventx))) {
10379 heap[eventnum] = heap[parent];
10380 heap[eventnum]->heapposition = eventnum;
10383 notdone = eventnum > 0;
10386 heap[eventnum] = newevent;
10387 newevent->heapposition = eventnum;
10394#ifdef ANSI_DECLARATORS
10395void eventheapify(
struct event **heap,
int heapsize,
int eventnum)
10397void eventheapify(heap, heapsize, eventnum)
10398struct event **heap;
10404 struct event *thisevent;
10405 REAL eventx, eventy;
10406 int leftchild, rightchild;
10410 thisevent = heap[eventnum];
10411 eventx = thisevent->xkey;
10412 eventy = thisevent->ykey;
10413 leftchild = 2 * eventnum + 1;
10414 notdone = leftchild < heapsize;
10416 if ((heap[leftchild]->ykey < eventy) ||
10417 ((heap[leftchild]->ykey == eventy)
10418 && (heap[leftchild]->xkey < eventx))) {
10419 smallest = leftchild;
10421 smallest = eventnum;
10423 rightchild = leftchild + 1;
10424 if (rightchild < heapsize) {
10425 if ((heap[rightchild]->ykey < heap[smallest]->ykey) ||
10426 ((heap[rightchild]->ykey == heap[smallest]->ykey)
10427 && (heap[rightchild]->xkey < heap[smallest]->xkey))) {
10428 smallest = rightchild;
10431 if (smallest == eventnum) {
10434 heap[eventnum] = heap[smallest];
10435 heap[eventnum]->heapposition = eventnum;
10436 heap[smallest] = thisevent;
10437 thisevent->heapposition = smallest;
10439 eventnum = smallest;
10440 leftchild = 2 * eventnum + 1;
10441 notdone = leftchild < heapsize;
10450#ifdef ANSI_DECLARATORS
10451void eventheapdelete(
struct event **heap,
int heapsize,
int eventnum)
10453void eventheapdelete(heap, heapsize, eventnum)
10454struct event **heap;
10460 struct event *moveevent;
10461 REAL eventx, eventy;
10465 moveevent = heap[heapsize - 1];
10466 if (eventnum > 0) {
10467 eventx = moveevent->xkey;
10468 eventy = moveevent->ykey;
10470 parent = (eventnum - 1) >> 1;
10471 if ((heap[parent]->ykey < eventy) ||
10472 ((heap[parent]->ykey == eventy)
10473 && (heap[parent]->xkey <= eventx))) {
10476 heap[eventnum] = heap[parent];
10477 heap[eventnum]->heapposition = eventnum;
10480 notdone = eventnum > 0;
10484 heap[eventnum] = moveevent;
10485 moveevent->heapposition = eventnum;
10486 eventheapify(heap, heapsize - 1, eventnum);
10493#ifdef ANSI_DECLARATORS
10494void createeventheap(
struct mesh *m,
struct event ***eventheap,
10495 struct event **events,
struct event **freeevents)
10497void createeventheap(m, eventheap, events, freeevents)
10499struct event ***eventheap;
10500struct event **events;
10501struct event **freeevents;
10509 maxevents = (3 * m->invertices) / 2;
10510 *eventheap = (
struct event **) trimalloc(maxevents *
10511 (
int)
sizeof(
struct event *));
10512 *events = (
struct event *) trimalloc(maxevents * (
int)
sizeof(
struct event));
10513 traversalinit(&m->vertices);
10514 for (i = 0; i < m->invertices; i++) {
10515 thisvertex = vertextraverse(m);
10516 (*events)[i].eventptr = (VOID *) thisvertex;
10517 (*events)[i].xkey = thisvertex[0];
10518 (*events)[i].ykey = thisvertex[1];
10519 eventheapinsert(*eventheap, i, *events + i);
10521 *freeevents = (
struct event *) NULL;
10522 for (i = maxevents - 1; i >= m->invertices; i--) {
10523 (*events)[i].eventptr = (VOID *) *freeevents;
10524 *freeevents = *events + i;
10532#ifdef ANSI_DECLARATORS
10533int rightofhyperbola(
struct mesh *m,
struct otri *fronttri, vertex newsite)
10535int rightofhyperbola(m, fronttri, newsite)
10537struct otri *fronttri;
10542 vertex leftvertex, rightvertex;
10543 REAL dxa, dya, dxb, dyb;
10545 m->hyperbolacount++;
10547 dest(*fronttri, leftvertex);
10548 apex(*fronttri, rightvertex);
10549 if ((leftvertex[1] < rightvertex[1]) ||
10550 ((leftvertex[1] == rightvertex[1]) &&
10551 (leftvertex[0] < rightvertex[0]))) {
10552 if (newsite[0] >= rightvertex[0]) {
10556 if (newsite[0] <= leftvertex[0]) {
10560 dxa = leftvertex[0] - newsite[0];
10561 dya = leftvertex[1] - newsite[1];
10562 dxb = rightvertex[0] - newsite[0];
10563 dyb = rightvertex[1] - newsite[1];
10564 return dya * (dxb * dxb + dyb * dyb) > dyb * (dxa * dxa + dya * dya);
10571#ifdef ANSI_DECLARATORS
10572REAL circletop(
struct mesh *m, vertex pa, vertex pb, vertex pc, REAL ccwabc)
10574REAL circletop(m, pa, pb, pc, ccwabc)
10583 REAL xac, yac, xbc, ybc, xab, yab;
10584 REAL aclen2, bclen2, ablen2;
10586 m->circletopcount++;
10588 xac = pa[0] - pc[0];
10589 yac = pa[1] - pc[1];
10590 xbc = pb[0] - pc[0];
10591 ybc = pb[1] - pc[1];
10592 xab = pa[0] - pb[0];
10593 yab = pa[1] - pb[1];
10594 aclen2 = xac * xac + yac * yac;
10595 bclen2 = xbc * xbc + ybc * ybc;
10596 ablen2 = xab * xab + yab * yab;
10597 return pc[1] + (xac * bclen2 - xbc * aclen2 + sqrt(aclen2 * bclen2 * ablen2))
10605#ifdef ANSI_DECLARATORS
10606void check4deadevent(
struct otri *checktri,
struct event **freeevents,
10607 struct event **eventheap,
int *heapsize)
10609void check4deadevent(checktri, freeevents, eventheap, heapsize)
10610struct otri *checktri;
10611struct event **freeevents;
10612struct event **eventheap;
10617 struct event *deadevent;
10618 vertex eventvertex;
10621 org(*checktri, eventvertex);
10622 if (eventvertex != (vertex) NULL) {
10623 deadevent = (
struct event *) eventvertex;
10624 eventnum = deadevent->heapposition;
10625 deadevent->eventptr = (VOID *) *freeevents;
10626 *freeevents = deadevent;
10627 eventheapdelete(eventheap, *heapsize, eventnum);
10629 setorg(*checktri, NULL);
10637#ifdef ANSI_DECLARATORS
10639 vertex searchpoint,
struct otri *searchtri)
10641struct splaynode *splay(m, splaytree, searchpoint, searchtri)
10645struct otri *searchtri;
10650 struct splaynode *lefttree, *righttree;
10652 vertex checkvertex;
10653 int rightofroot, rightofchild;
10655 if (splaytree == (
struct splaynode *) NULL) {
10658 dest(splaytree->keyedge, checkvertex);
10659 if (checkvertex == splaytree->keydest) {
10660 rightofroot = rightofhyperbola(m, &splaytree->keyedge, searchpoint);
10662 otricopy(splaytree->keyedge, *searchtri);
10663 child = splaytree->rchild;
10665 child = splaytree->lchild;
10667 if (child == (
struct splaynode *) NULL) {
10670 dest(child->keyedge, checkvertex);
10671 if (checkvertex != child->keydest) {
10672 child = splay(m, child, searchpoint, searchtri);
10673 if (child == (
struct splaynode *) NULL) {
10675 splaytree->rchild = (
struct splaynode *) NULL;
10677 splaytree->lchild = (
struct splaynode *) NULL;
10682 rightofchild = rightofhyperbola(m, &child->keyedge, searchpoint);
10683 if (rightofchild) {
10684 otricopy(child->keyedge, *searchtri);
10685 grandchild = splay(m, child->rchild, searchpoint, searchtri);
10686 child->rchild = grandchild;
10688 grandchild = splay(m, child->lchild, searchpoint, searchtri);
10689 child->lchild = grandchild;
10691 if (grandchild == (
struct splaynode *) NULL) {
10693 splaytree->rchild = child->lchild;
10694 child->lchild = splaytree;
10696 splaytree->lchild = child->rchild;
10697 child->rchild = splaytree;
10701 if (rightofchild) {
10703 splaytree->rchild = child->lchild;
10704 child->lchild = splaytree;
10706 splaytree->lchild = grandchild->rchild;
10707 grandchild->rchild = splaytree;
10709 child->rchild = grandchild->lchild;
10710 grandchild->lchild = child;
10713 splaytree->rchild = grandchild->lchild;
10714 grandchild->lchild = splaytree;
10716 splaytree->lchild = child->rchild;
10717 child->rchild = splaytree;
10719 child->lchild = grandchild->rchild;
10720 grandchild->rchild = child;
10724 lefttree = splay(m, splaytree->lchild, searchpoint, searchtri);
10725 righttree = splay(m, splaytree->rchild, searchpoint, searchtri);
10727 pooldealloc(&m->splaynodes, (VOID *) splaytree);
10728 if (lefttree == (
struct splaynode *) NULL) {
10730 }
else if (righttree == (
struct splaynode *) NULL) {
10732 }
else if (lefttree->rchild == (
struct splaynode *) NULL) {
10733 lefttree->rchild = righttree->lchild;
10734 righttree->lchild = lefttree;
10736 }
else if (righttree->lchild == (
struct splaynode *) NULL) {
10737 righttree->lchild = lefttree->rchild;
10738 lefttree->rchild = righttree;
10742 leftright = lefttree->rchild;
10743 while (leftright->rchild != (
struct splaynode *) NULL) {
10744 leftright = leftright->rchild;
10746 leftright->rchild = righttree;
10756#ifdef ANSI_DECLARATORS
10758 struct otri *newkey, vertex searchpoint)
10760struct splaynode *splayinsert(m, splayroot, newkey, searchpoint)
10763struct otri *newkey;
10770 newsplaynode = (
struct splaynode *) poolalloc(&m->splaynodes);
10771 otricopy(*newkey, newsplaynode->keyedge);
10772 dest(*newkey, newsplaynode->keydest);
10773 if (splayroot == (
struct splaynode *) NULL) {
10774 newsplaynode->lchild = (
struct splaynode *) NULL;
10775 newsplaynode->rchild = (
struct splaynode *) NULL;
10776 }
else if (rightofhyperbola(m, &splayroot->keyedge, searchpoint)) {
10777 newsplaynode->lchild = splayroot;
10778 newsplaynode->rchild = splayroot->rchild;
10779 splayroot->rchild = (
struct splaynode *) NULL;
10781 newsplaynode->lchild = splayroot->lchild;
10782 newsplaynode->rchild = splayroot;
10783 splayroot->lchild = (
struct splaynode *) NULL;
10785 return newsplaynode;
10792#ifdef ANSI_DECLARATORS
10795 struct otri *newkey,
10796 vertex pa, vertex pb, vertex pc, REAL topy)
10798struct splaynode *circletopinsert(m, b, splayroot, newkey, pa, pb, pc, topy)
10802struct otri *newkey;
10811 REAL xac, yac, xbc, ybc;
10812 REAL aclen2, bclen2;
10813 REAL searchpoint[2];
10814 struct otri dummytri;
10816 ccwabc = counterclockwise(m, b, pa, pb, pc);
10817 xac = pa[0] - pc[0];
10818 yac = pa[1] - pc[1];
10819 xbc = pb[0] - pc[0];
10820 ybc = pb[1] - pc[1];
10821 aclen2 = xac * xac + yac * yac;
10822 bclen2 = xbc * xbc + ybc * ybc;
10823 searchpoint[0] = pc[0] - (yac * bclen2 - ybc * aclen2) / (2.0 * ccwabc);
10824 searchpoint[1] = topy;
10825 return splayinsert(m, splay(m, splayroot, (vertex) searchpoint, &dummytri),
10826 newkey, (vertex) searchpoint);
10833#ifdef ANSI_DECLARATORS
10835 struct otri *bottommost, vertex searchvertex,
10836 struct otri *searchtri,
int *farright)
10838struct splaynode *frontlocate(m, splayroot, bottommost, searchvertex,
10839 searchtri, farright)
10842struct otri *bottommost;
10843vertex searchvertex;
10844struct otri *searchtri;
10852 otricopy(*bottommost, *searchtri);
10853 splayroot = splay(m, splayroot, searchvertex, searchtri);
10856 while (!farrightflag && rightofhyperbola(m, searchtri, searchvertex)) {
10857 onextself(*searchtri);
10858 farrightflag = otriequal(*searchtri, *bottommost);
10860 *farright = farrightflag;
10868#ifdef ANSI_DECLARATORS
10869long sweeplinedelaunay(
struct mesh *m,
struct behavior *b)
10871long sweeplinedelaunay(m, b)
10877 struct event **eventheap;
10878 struct event *events;
10879 struct event *freeevents;
10880 struct event *nextevent;
10881 struct event *newevent;
10883 struct otri bottommost;
10884 struct otri searchtri;
10885 struct otri fliptri;
10886 struct otri lefttri, righttri, farlefttri, farrighttri;
10887 struct otri inserttri;
10888 vertex firstvertex, secondvertex;
10889 vertex nextvertex, lastvertex;
10890 vertex connectvertex;
10891 vertex leftvertex, midvertex, rightvertex;
10892 REAL lefttest, righttest;
10894 int check4events, farrightflag;
10897 poolinit(&m->splaynodes,
sizeof(
struct splaynode), SPLAYNODEPERBLOCK,
10898 SPLAYNODEPERBLOCK, 0);
10902 printf(
" Placing vertices in event heap.\n");
10904 createeventheap(m, &eventheap, &events, &freeevents);
10905 heapsize = m->invertices;
10908 printf(
" Forming triangulation.\n");
10910 maketriangle(m, b, &lefttri);
10911 maketriangle(m, b, &righttri);
10912 bond(lefttri, righttri);
10913 lnextself(lefttri);
10914 lprevself(righttri);
10915 bond(lefttri, righttri);
10916 lnextself(lefttri);
10917 lprevself(righttri);
10918 bond(lefttri, righttri);
10919 firstvertex = (vertex) eventheap[0]->eventptr;
10920 eventheap[0]->eventptr = (VOID *) freeevents;
10921 freeevents = eventheap[0];
10922 eventheapdelete(eventheap, heapsize, 0);
10925 if (heapsize == 0) {
10926 printf(
"Error: Input vertices are all identical.\n");
10929 secondvertex = (vertex) eventheap[0]->eventptr;
10930 eventheap[0]->eventptr = (VOID *) freeevents;
10931 freeevents = eventheap[0];
10932 eventheapdelete(eventheap, heapsize, 0);
10934 if ((firstvertex[0] == secondvertex[0]) &&
10935 (firstvertex[1] == secondvertex[1])) {
10938"Warning: A duplicate vertex at (%.12g, %.12g) appeared and was ignored.\n",
10939 secondvertex[0], secondvertex[1]);
10941 setvertextype(secondvertex, UNDEADVERTEX);
10944 }
while ((firstvertex[0] == secondvertex[0]) &&
10945 (firstvertex[1] == secondvertex[1]));
10946 setorg(lefttri, firstvertex);
10947 setdest(lefttri, secondvertex);
10948 setorg(righttri, secondvertex);
10949 setdest(righttri, firstvertex);
10950 lprev(lefttri, bottommost);
10951 lastvertex = secondvertex;
10952 while (heapsize > 0) {
10953 nextevent = eventheap[0];
10954 eventheapdelete(eventheap, heapsize, 0);
10957 if (nextevent->xkey < m->xmin) {
10958 decode(nextevent->eventptr, fliptri);
10959 oprev(fliptri, farlefttri);
10960 check4deadevent(&farlefttri, &freeevents, eventheap, &heapsize);
10961 onext(fliptri, farrighttri);
10962 check4deadevent(&farrighttri, &freeevents, eventheap, &heapsize);
10964 if (otriequal(farlefttri, bottommost)) {
10965 lprev(fliptri, bottommost);
10967 flip(m, b, &fliptri);
10968 setapex(fliptri, NULL);
10969 lprev(fliptri, lefttri);
10970 lnext(fliptri, righttri);
10971 sym(lefttri, farlefttri);
10973 if (randomnation(SAMPLERATE) == 0) {
10975 dest(fliptri, leftvertex);
10976 apex(fliptri, midvertex);
10977 org(fliptri, rightvertex);
10978 splayroot = circletopinsert(m, b, splayroot, &lefttri, leftvertex,
10979 midvertex, rightvertex, nextevent->ykey);
10982 nextvertex = (vertex) nextevent->eventptr;
10983 if ((nextvertex[0] == lastvertex[0]) &&
10984 (nextvertex[1] == lastvertex[1])) {
10987"Warning: A duplicate vertex at (%.12g, %.12g) appeared and was ignored.\n",
10988 nextvertex[0], nextvertex[1]);
10990 setvertextype(nextvertex, UNDEADVERTEX);
10994 lastvertex = nextvertex;
10996 splayroot = frontlocate(m, splayroot, &bottommost, nextvertex,
10997 &searchtri, &farrightflag);
11007 check4deadevent(&searchtri, &freeevents, eventheap, &heapsize);
11009 otricopy(searchtri, farrighttri);
11010 sym(searchtri, farlefttri);
11011 maketriangle(m, b, &lefttri);
11012 maketriangle(m, b, &righttri);
11013 dest(farrighttri, connectvertex);
11014 setorg(lefttri, connectvertex);
11015 setdest(lefttri, nextvertex);
11016 setorg(righttri, nextvertex);
11017 setdest(righttri, connectvertex);
11018 bond(lefttri, righttri);
11019 lnextself(lefttri);
11020 lprevself(righttri);
11021 bond(lefttri, righttri);
11022 lnextself(lefttri);
11023 lprevself(righttri);
11024 bond(lefttri, farlefttri);
11025 bond(righttri, farrighttri);
11026 if (!farrightflag && otriequal(farrighttri, bottommost)) {
11027 otricopy(lefttri, bottommost);
11030 if (randomnation(SAMPLERATE) == 0) {
11031 splayroot = splayinsert(m, splayroot, &lefttri, nextvertex);
11032 }
else if (randomnation(SAMPLERATE) == 0) {
11033 lnext(righttri, inserttri);
11034 splayroot = splayinsert(m, splayroot, &inserttri, nextvertex);
11038 nextevent->eventptr = (VOID *) freeevents;
11039 freeevents = nextevent;
11041 if (check4events) {
11042 apex(farlefttri, leftvertex);
11043 dest(lefttri, midvertex);
11044 apex(lefttri, rightvertex);
11045 lefttest = counterclockwise(m, b, leftvertex, midvertex, rightvertex);
11046 if (lefttest > 0.0) {
11047 newevent = freeevents;
11048 freeevents = (
struct event *) freeevents->eventptr;
11049 newevent->xkey = m->xminextreme;
11050 newevent->ykey = circletop(m, leftvertex, midvertex, rightvertex,
11052 newevent->eventptr = (VOID *) encode(lefttri);
11053 eventheapinsert(eventheap, heapsize, newevent);
11055 setorg(lefttri, newevent);
11057 apex(righttri, leftvertex);
11058 org(righttri, midvertex);
11059 apex(farrighttri, rightvertex);
11060 righttest = counterclockwise(m, b, leftvertex, midvertex, rightvertex);
11061 if (righttest > 0.0) {
11062 newevent = freeevents;
11063 freeevents = (
struct event *) freeevents->eventptr;
11064 newevent->xkey = m->xminextreme;
11065 newevent->ykey = circletop(m, leftvertex, midvertex, rightvertex,
11067 newevent->eventptr = (VOID *) encode(farrighttri);
11068 eventheapinsert(eventheap, heapsize, newevent);
11070 setorg(farrighttri, newevent);
11075 pooldeinit(&m->splaynodes);
11076 lprevself(bottommost);
11077 return removeghosts(m, b, &bottommost);
11096#ifdef ANSI_DECLARATORS
11108 initializetrisubpools(m, b);
11113 "Constructing Delaunay triangulation by divide-and-conquer method.\n");
11115 hulledges = divconqdelaunay(m, b);
11118 printf(
"Constructing Delaunay triangulation ");
11119 if (b->incremental) {
11120 printf(
"by incremental method.\n");
11121 }
else if (b->sweepline) {
11122 printf(
"by sweepline method.\n");
11124 printf(
"by divide-and-conquer method.\n");
11127 if (b->incremental) {
11128 hulledges = incrementaldelaunay(m, b);
11129 }
else if (b->sweepline) {
11130 hulledges = sweeplinedelaunay(m, b);
11132 hulledges = divconqdelaunay(m, b);
11136 if (m->triangles.items == 0) {
11173#ifdef ANSI_DECLARATORS
11174int reconstruct(
struct mesh *m,
struct behavior *b,
int *trianglelist,
11175 REAL *triangleattriblist, REAL *trianglearealist,
11176 int elements,
int corners,
int attribs,
11177 int *segmentlist,
int *segmentmarkerlist,
int numberofsegments)
11179int reconstruct(m, b, trianglelist, triangleattriblist, trianglearealist,
11180 elements, corners, attribs, segmentlist, segmentmarkerlist,
11185REAL *triangleattriblist;
11186REAL *trianglearealist;
11191int *segmentmarkerlist;
11192int numberofsegments;
11197#ifdef ANSI_DECLARATORS
11198long reconstruct(
struct mesh *m,
struct behavior *b,
char *elefilename,
11199 char *areafilename,
char* polyfilename, FILE* polyfile)
11201long reconstruct(m, b, elefilename, areafilename, polyfilename, polyfile)
11219 char inputline[INPUTLINESIZE];
11223 struct otri triangleloop;
11224 struct otri triangleleft;
11225 struct otri checktri;
11226 struct otri checkleft;
11227 struct otri checkneighbor;
11228 struct osub subsegloop;
11229 triangle *vertexarray;
11230 triangle* prevlink;
11232 vertex tdest, tapex;
11233 vertex checkdest, checkapex;
11236 vertex segmentorg, segmentdest;
11240 int killvertexindex;
11242 int segmentmarkers;
11247 long elementnumber, segmentnumber;
11252 m->inelements = elements;
11253 incorners = corners;
11254 if (incorners < 3) {
11255 printf(
"Error: Triangles must have at least 3 vertices.\n");
11258 m->eextras = attribs;
11262 printf(
"Opening %s.\n", elefilename);
11264 elefile = fopen(elefilename,
"r");
11265 if (elefile == (FILE *) NULL) {
11266 printf(
" Error: Cannot access file %s.\n", elefilename);
11271 stringptr = readline(inputline, elefile, elefilename);
11272 m->inelements = (int) strtol(stringptr, &stringptr, 0);
11273 stringptr = findfield(stringptr);
11274 if (*stringptr ==
'\0') {
11277 incorners = (int) strtol(stringptr, &stringptr, 0);
11278 if (incorners < 3) {
11279 printf(
"Error: Triangles in %s must have at least 3 vertices.\n",
11284 stringptr = findfield(stringptr);
11285 if (*stringptr ==
'\0') {
11288 m->eextras = (int) strtol(stringptr, &stringptr, 0);
11292 initializetrisubpools(m, b);
11295 for (elementnumber = 1; elementnumber <= m->inelements; elementnumber++) {
11296 maketriangle(m, b, &triangleloop);
11298 triangleloop.tri[3] = (triangle) triangleloop.tri;
11301 segmentmarkers = 0;
11304 m->insegments = numberofsegments;
11305 segmentmarkers = segmentmarkerlist != (
int *) NULL;
11309 stringptr = readline(inputline, polyfile, b->inpolyfilename);
11310 m->insegments = (int) strtol(stringptr, &stringptr, 0);
11311 stringptr = findfield(stringptr);
11312 if (*stringptr !=
'\0') {
11313 segmentmarkers = (int) strtol(stringptr, &stringptr, 0);
11318 for (segmentnumber = 1; segmentnumber <= m->insegments; segmentnumber++) {
11319 makesubseg(m, &subsegloop);
11321 subsegloop.ss[2] = (subseg) subsegloop.ss;
11332 printf(
"Opening %s.\n", areafilename);
11334 areafile = fopen(areafilename,
"r");
11335 if (areafile == (FILE *) NULL) {
11336 printf(
" Error: Cannot access file %s.\n", areafilename);
11339 stringptr = readline(inputline, areafile, areafilename);
11340 areaelements = (int) strtol(stringptr, &stringptr, 0);
11341 if (areaelements != m->inelements) {
11342 printf(
"Error: %s and %s disagree on number of triangles.\n",
11343 elefilename, areafilename);
11350 printf(
"Reconstructing mesh.\n");
11355 vertexarray = (triangle *) trimalloc(m->vertices.items *
11356 (
int)
sizeof(triangle));
11358 for (i = 0; i < m->vertices.items; i++) {
11359 vertexarray[i] = (triangle) m->dummytri;
11363 printf(
" Assembling triangles.\n");
11367 traversalinit(&m->triangles);
11368 triangleloop.tri = triangletraverse(m);
11369 elementnumber = b->firstnumber;
11370 while (triangleloop.tri != (triangle *) NULL) {
11373 for (j = 0; j < 3; j++) {
11374 corner[j] = trianglelist[vertexindex++];
11375 if ((corner[j] < b->firstnumber) ||
11376 (corner[j] >= b->firstnumber + m->invertices)) {
11377 printf(
"Error: Triangle %ld has an invalid vertex index.\n",
11384 stringptr = readline(inputline, elefile, elefilename);
11385 for (j = 0; j < 3; j++) {
11386 stringptr = findfield(stringptr);
11387 if (*stringptr ==
'\0') {
11388 printf(
"Error: Triangle %ld is missing vertex %d in %s.\n",
11389 elementnumber, j + 1, elefilename);
11392 corner[j] = (int) strtol(stringptr, &stringptr, 0);
11393 if ((corner[j] < b->firstnumber) ||
11394 (corner[j] >= b->firstnumber + m->invertices)) {
11395 printf(
"Error: Triangle %ld has an invalid vertex index.\n",
11404 for (j = 3; j < incorners; j++) {
11406 killvertexindex = trianglelist[vertexindex++];
11408 stringptr = findfield(stringptr);
11409 if (*stringptr !=
'\0') {
11410 killvertexindex = (int) strtol(stringptr, &stringptr, 0);
11412 if ((killvertexindex >= b->firstnumber) &&
11413 (killvertexindex < b->firstnumber + m->invertices)) {
11415 killvertex = getvertex(m, b, killvertexindex);
11416 if (vertextype(killvertex) != DEADVERTEX) {
11417 vertexdealloc(m, killvertex);
11426 for (j = 0; j < m->eextras; j++) {
11428 setelemattribute(triangleloop, j, triangleattriblist[attribindex++]);
11430 stringptr = findfield(stringptr);
11431 if (*stringptr ==
'\0') {
11432 setelemattribute(triangleloop, j, 0);
11434 setelemattribute(triangleloop, j,
11435 (REAL) strtod(stringptr, &stringptr));
11442 area = trianglearealist[elementnumber - b->firstnumber];
11445 stringptr = readline(inputline, areafile, areafilename);
11446 stringptr = findfield(stringptr);
11447 if (*stringptr ==
'\0') {
11450 area = (REAL) strtod(stringptr, &stringptr);
11453 setareabound(triangleloop, area);
11457 triangleloop.orient = 0;
11458 setorg(triangleloop, getvertex(m, b, corner[0]));
11459 setdest(triangleloop, getvertex(m, b, corner[1]));
11460 setapex(triangleloop, getvertex(m, b, corner[2]));
11462 for (triangleloop.orient = 0; triangleloop.orient < 3;
11463 triangleloop.orient++) {
11465 aroundvertex = corner[triangleloop.orient];
11467 nexttri = vertexarray[aroundvertex - b->firstnumber];
11469 triangleloop.tri[6 + triangleloop.orient] = nexttri;
11471 vertexarray[aroundvertex - b->firstnumber] = encode(triangleloop);
11472 decode(nexttri, checktri);
11473 if (checktri.tri != m->dummytri) {
11474 dest(triangleloop, tdest);
11475 apex(triangleloop, tapex);
11478 dest(checktri, checkdest);
11479 apex(checktri, checkapex);
11480 if (tapex == checkdest) {
11482 lprev(triangleloop, triangleleft);
11483 bond(triangleleft, checktri);
11485 if (tdest == checkapex) {
11487 lprev(checktri, checkleft);
11488 bond(triangleloop, checkleft);
11491 nexttri = checktri.tri[6 + checktri.orient];
11492 decode(nexttri, checktri);
11493 }
while (checktri.tri != m->dummytri);
11496 triangleloop.tri = triangletraverse(m);
11512 printf(
" Marking segments in triangulation.\n");
11517 traversalinit(&m->subsegs);
11518 subsegloop.ss = subsegtraverse(m);
11519 segmentnumber = b->firstnumber;
11520 while (subsegloop.ss != (subseg *) NULL) {
11522 end[0] = segmentlist[vertexindex++];
11523 end[1] = segmentlist[vertexindex++];
11524 if (segmentmarkers) {
11525 boundmarker = segmentmarkerlist[segmentnumber - b->firstnumber];
11529 stringptr = readline(inputline, polyfile, b->inpolyfilename);
11531 stringptr = findfield(stringptr);
11532 if (*stringptr ==
'\0') {
11533 printf(
"Error: Segment %ld has no endpoints in %s.\n", segmentnumber,
11537 end[0] = (int) strtol(stringptr, &stringptr, 0);
11539 stringptr = findfield(stringptr);
11540 if (*stringptr ==
'\0') {
11541 printf(
"Error: Segment %ld is missing its second endpoint in %s.\n",
11542 segmentnumber, polyfilename);
11545 end[1] = (int) strtol(stringptr, &stringptr, 0);
11547 if (segmentmarkers) {
11548 stringptr = findfield(stringptr);
11549 if (*stringptr ==
'\0') {
11552 boundmarker = (int) strtol(stringptr, &stringptr, 0);
11556 for (j = 0; j < 2; j++) {
11557 if ((end[j] < b->firstnumber) ||
11558 (end[j] >= b->firstnumber + m->invertices)) {
11559 printf(
"Error: Segment %ld has an invalid vertex index.\n",
11566 subsegloop.ssorient = 0;
11567 segmentorg = getvertex(m, b, end[0]);
11568 segmentdest = getvertex(m, b, end[1]);
11569 setsorg(subsegloop, segmentorg);
11570 setsdest(subsegloop, segmentdest);
11571 setsegorg(subsegloop, segmentorg);
11572 setsegdest(subsegloop, segmentdest);
11573 setmark(subsegloop, boundmarker);
11575 for (subsegloop.ssorient = 0; subsegloop.ssorient < 2;
11576 subsegloop.ssorient++) {
11578 aroundvertex = end[1 - subsegloop.ssorient];
11580 prevlink = &vertexarray[aroundvertex - b->firstnumber];
11581 nexttri = vertexarray[aroundvertex - b->firstnumber];
11582 decode(nexttri, checktri);
11583 sorg(subsegloop, shorg);
11592 while (notfound && (checktri.tri != m->dummytri)) {
11593 dest(checktri, checkdest);
11594 if (shorg == checkdest) {
11596 * prevlink = checktri.tri[6 + checktri.orient];
11598 tsbond(checktri, subsegloop);
11600 sym(checktri, checkneighbor);
11601 if (checkneighbor.tri == m->dummytri) {
11605 insertsubseg(m, b, &checktri, 1);
11611 prevlink = &checktri.tri[6 + checktri.orient];
11612 nexttri = checktri.tri[6 + checktri.orient];
11613 decode(nexttri, checktri);
11616 subsegloop.ss = subsegtraverse(m);
11623 for (i = 0; i < m->vertices.items; i++) {
11625 nexttri = vertexarray[i];
11626 decode(nexttri, checktri);
11627 while (checktri.tri != m->dummytri) {
11630 nexttri = checktri.tri[6 + checktri.orient];
11632 tsdissolve(checktri);
11633 sym(checktri, checkneighbor);
11634 if (checkneighbor.tri == m->dummytri) {
11635 insertsubseg(m, b, &checktri, 1);
11638 decode(nexttri, checktri);
11642 trifree((VOID *) vertexarray);
11673#ifdef ANSI_DECLARATORS
11674enum finddirectionresult finddirection(
struct mesh *m,
struct behavior *b,
11675 struct otri *searchtri,
11676 vertex searchpoint)
11678enum finddirectionresult finddirection(m, b, searchtri, searchpoint)
11681struct otri *searchtri;
11686 struct otri checktri;
11687 vertex startvertex;
11688 vertex leftvertex, rightvertex;
11689 REAL leftccw, rightccw;
11690 int leftflag, rightflag;
11693 org(*searchtri, startvertex);
11694 dest(*searchtri, rightvertex);
11695 apex(*searchtri, leftvertex);
11697 leftccw = counterclockwise(m, b, searchpoint, startvertex, leftvertex);
11698 leftflag = leftccw > 0.0;
11700 rightccw = counterclockwise(m, b, startvertex, searchpoint, rightvertex);
11701 rightflag = rightccw > 0.0;
11702 if (leftflag && rightflag) {
11705 onext(*searchtri, checktri);
11706 if (checktri.tri == m->dummytri) {
11714 onextself(*searchtri);
11715 if (searchtri->tri == m->dummytri) {
11716 printf(
"Internal error in finddirection(): Unable to find a\n");
11717 printf(
" triangle leading from (%.12g, %.12g) to", startvertex[0],
11719 printf(
" (%.12g, %.12g).\n", searchpoint[0], searchpoint[1]);
11722 apex(*searchtri, leftvertex);
11723 rightccw = leftccw;
11724 leftccw = counterclockwise(m, b, searchpoint, startvertex, leftvertex);
11725 leftflag = leftccw > 0.0;
11727 while (rightflag) {
11729 oprevself(*searchtri);
11730 if (searchtri->tri == m->dummytri) {
11731 printf(
"Internal error in finddirection(): Unable to find a\n");
11732 printf(
" triangle leading from (%.12g, %.12g) to", startvertex[0],
11734 printf(
" (%.12g, %.12g).\n", searchpoint[0], searchpoint[1]);
11737 dest(*searchtri, rightvertex);
11738 leftccw = rightccw;
11739 rightccw = counterclockwise(m, b, startvertex, searchpoint, rightvertex);
11740 rightflag = rightccw > 0.0;
11742 if (leftccw == 0.0) {
11743 return LEFTCOLLINEAR;
11744 }
else if (rightccw == 0.0) {
11745 return RIGHTCOLLINEAR;
11768#ifdef ANSI_DECLARATORS
11769void segmentintersection(
struct mesh *m,
struct behavior *b,
11770 struct otri *splittri,
struct osub *splitsubseg,
11773void segmentintersection(m, b, splittri, splitsubseg, endpoint2)
11776struct otri *splittri;
11777struct osub *splitsubseg;
11782 struct osub opposubseg;
11784 vertex torg, tdest;
11785 vertex leftvertex, rightvertex;
11787 enum insertvertexresult success;
11798 apex(*splittri, endpoint1);
11799 org(*splittri, torg);
11800 dest(*splittri, tdest);
11802 tx = tdest[0] - torg[0];
11803 ty = tdest[1] - torg[1];
11804 ex = endpoint2[0] - endpoint1[0];
11805 ey = endpoint2[1] - endpoint1[1];
11806 etx = torg[0] - endpoint2[0];
11807 ety = torg[1] - endpoint2[1];
11808 denom = ty * ex - tx * ey;
11809 if (denom == 0.0) {
11810 printf(
"Internal error in segmentintersection():");
11811 printf(
" Attempt to find intersection of parallel segments.\n");
11814 split = (ey * etx - ex * ety) / denom;
11816 newvertex = (vertex) poolalloc(&m->vertices);
11818 for (i = 0; i < 2 + m->nextras; i++) {
11819 newvertex[i] = torg[i] + split * (tdest[i] - torg[i]);
11821 setvertexmark(newvertex, mark(*splitsubseg));
11822 setvertextype(newvertex, INPUTVERTEX);
11823 if (b->verbose > 1) {
11825 " Splitting subsegment (%.12g, %.12g) (%.12g, %.12g) at (%.12g, %.12g).\n",
11826 torg[0], torg[1], tdest[0], tdest[1], newvertex[0], newvertex[1]);
11829 success = insertvertex(m, b, newvertex, splittri, splitsubseg, 0, 0);
11830 if (success != SUCCESSFULVERTEX) {
11831 printf(
"Internal error in segmentintersection():\n");
11832 printf(
" Failure to split a segment.\n");
11836 setvertex2tri(newvertex, encode(*splittri));
11837 if (m->steinerleft > 0) {
11842 ssymself(*splitsubseg);
11843 spivot(*splitsubseg, opposubseg);
11844 sdissolve(*splitsubseg);
11845 sdissolve(opposubseg);
11847 setsegorg(*splitsubseg, newvertex);
11848 snextself(*splitsubseg);
11849 }
while (splitsubseg->ss != m->dummysub);
11851 setsegorg(opposubseg, newvertex);
11852 snextself(opposubseg);
11853 }
while (opposubseg.ss != m->dummysub);
11857 finddirection(m, b, splittri, endpoint1);
11858 dest(*splittri, rightvertex);
11859 apex(*splittri, leftvertex);
11860 if ((leftvertex[0] == endpoint1[0]) && (leftvertex[1] == endpoint1[1])) {
11861 onextself(*splittri);
11862 }
else if ((rightvertex[0] != endpoint1[0]) ||
11863 (rightvertex[1] != endpoint1[1])) {
11864 printf(
"Internal error in segmentintersection():\n");
11865 printf(
" Topological inconsistency after splitting a segment.\n");
11895#ifdef ANSI_DECLARATORS
11896int scoutsegment(
struct mesh *m,
struct behavior *b,
struct otri *searchtri,
11897 vertex endpoint2,
int newmark)
11899int scoutsegment(m, b, searchtri, endpoint2, newmark)
11902struct otri *searchtri;
11908 struct otri crosstri;
11909 struct osub crosssubseg;
11910 vertex leftvertex, rightvertex;
11911 enum finddirectionresult collinear;
11914 collinear = finddirection(m, b, searchtri, endpoint2);
11915 dest(*searchtri, rightvertex);
11916 apex(*searchtri, leftvertex);
11917 if (((leftvertex[0] == endpoint2[0]) && (leftvertex[1] == endpoint2[1])) ||
11918 ((rightvertex[0] == endpoint2[0]) && (rightvertex[1] == endpoint2[1]))) {
11920 if ((leftvertex[0] == endpoint2[0]) && (leftvertex[1] == endpoint2[1])) {
11921 lprevself(*searchtri);
11924 insertsubseg(m, b, searchtri, newmark);
11926 }
else if (collinear == LEFTCOLLINEAR) {
11929 lprevself(*searchtri);
11930 insertsubseg(m, b, searchtri, newmark);
11932 return scoutsegment(m, b, searchtri, endpoint2, newmark);
11933 }
else if (collinear == RIGHTCOLLINEAR) {
11935 insertsubseg(m, b, searchtri, newmark);
11937 lnextself(*searchtri);
11939 return scoutsegment(m, b, searchtri, endpoint2, newmark);
11941 lnext(*searchtri, crosstri);
11942 tspivot(crosstri, crosssubseg);
11944 if (crosssubseg.ss == m->dummysub) {
11948 segmentintersection(m, b, &crosstri, &crosssubseg, endpoint2);
11949 otricopy(crosstri, *searchtri);
11950 insertsubseg(m, b, searchtri, newmark);
11952 return scoutsegment(m, b, searchtri, endpoint2, newmark);
11979#ifdef ANSI_DECLARATORS
11980void conformingedge(
struct mesh *m,
struct behavior *b,
11981 vertex endpoint1, vertex endpoint2,
int newmark)
11983void conformingedge(m, b, endpoint1, endpoint2, newmark)
11992 struct otri searchtri1, searchtri2;
11993 struct osub brokensubseg;
11995 vertex midvertex1, midvertex2;
11996 enum insertvertexresult success;
12000 if (b->verbose > 2) {
12001 printf(
"Forcing segment into triangulation by recursive splitting:\n");
12002 printf(
" (%.12g, %.12g) (%.12g, %.12g)\n", endpoint1[0], endpoint1[1],
12003 endpoint2[0], endpoint2[1]);
12006 newvertex = (vertex) poolalloc(&m->vertices);
12008 for (i = 0; i < 2 + m->nextras; i++) {
12009 newvertex[i] = 0.5 * (endpoint1[i] + endpoint2[i]);
12011 setvertexmark(newvertex, newmark);
12012 setvertextype(newvertex, SEGMENTVERTEX);
12014 searchtri1.tri = m->dummytri;
12016 success = insertvertex(m, b, newvertex, &searchtri1, (
struct osub *) NULL,
12018 if (success == DUPLICATEVERTEX) {
12019 if (b->verbose > 2) {
12020 printf(
" Segment intersects existing vertex (%.12g, %.12g).\n",
12021 newvertex[0], newvertex[1]);
12024 vertexdealloc(m, newvertex);
12025 org(searchtri1, newvertex);
12027 if (success == VIOLATINGVERTEX) {
12028 if (b->verbose > 2) {
12029 printf(
" Two segments intersect at (%.12g, %.12g).\n",
12030 newvertex[0], newvertex[1]);
12033 tspivot(searchtri1, brokensubseg);
12034 success = insertvertex(m, b, newvertex, &searchtri1, &brokensubseg,
12036 if (success != SUCCESSFULVERTEX) {
12037 printf(
"Internal error in conformingedge():\n");
12038 printf(
" Failure to split a segment.\n");
12043 if (m->steinerleft > 0) {
12047 otricopy(searchtri1, searchtri2);
12053 finddirection(m, b, &searchtri2, endpoint2);
12054 if (!scoutsegment(m, b, &searchtri1, endpoint1, newmark)) {
12057 org(searchtri1, midvertex1);
12058 conformingedge(m, b, midvertex1, endpoint1, newmark);
12060 if (!scoutsegment(m, b, &searchtri2, endpoint2, newmark)) {
12063 org(searchtri2, midvertex2);
12064 conformingedge(m, b, midvertex2, endpoint2, newmark);
12109#ifdef ANSI_DECLARATORS
12111 struct otri *fixuptri,
int leftside)
12113void delaunayfixup(m, b, fixuptri, leftside)
12116struct otri *fixuptri;
12121 struct otri neartri;
12122 struct otri fartri;
12123 struct osub faredge;
12124 vertex nearvertex, leftvertex, rightvertex, farvertex;
12128 lnext(*fixuptri, neartri);
12129 sym(neartri, fartri);
12131 if (fartri.tri == m->dummytri) {
12134 tspivot(neartri, faredge);
12135 if (faredge.ss != m->dummysub) {
12139 apex(neartri, nearvertex);
12140 org(neartri, leftvertex);
12141 dest(neartri, rightvertex);
12142 apex(fartri, farvertex);
12145 if (counterclockwise(m, b, nearvertex, leftvertex, farvertex) <= 0.0) {
12151 if (counterclockwise(m, b, farvertex, rightvertex, nearvertex) <= 0.0) {
12157 if (counterclockwise(m, b, rightvertex, leftvertex, farvertex) > 0.0) {
12162 if (incircle(m, b, leftvertex, farvertex, rightvertex, nearvertex) <=
12168 flip(m, b, &neartri);
12169 lprevself(*fixuptri);
12171 delaunayfixup(m, b, fixuptri, leftside);
12172 delaunayfixup(m, b, &fartri, leftside);
12229#ifdef ANSI_DECLARATORS
12230void constrainededge(
struct mesh *m,
struct behavior *b,
12231 struct otri *starttri, vertex endpoint2,
int newmark)
12233void constrainededge(m, b, starttri, endpoint2, newmark)
12236struct otri *starttri;
12242 struct otri fixuptri, fixuptri2;
12243 struct osub crosssubseg;
12252 org(*starttri, endpoint1);
12253 lnext(*starttri, fixuptri);
12254 flip(m, b, &fixuptri);
12260 org(fixuptri, farvertex);
12263 if ((farvertex[0] == endpoint2[0]) && (farvertex[1] == endpoint2[1])) {
12264 oprev(fixuptri, fixuptri2);
12266 delaunayfixup(m, b, &fixuptri, 0);
12267 delaunayfixup(m, b, &fixuptri2, 1);
12273 area = counterclockwise(m, b, endpoint1, endpoint2, farvertex);
12277 oprev(fixuptri, fixuptri2);
12279 delaunayfixup(m, b, &fixuptri, 0);
12280 delaunayfixup(m, b, &fixuptri2, 1);
12284 oprev(fixuptri, fixuptri2);
12287 delaunayfixup(m, b, &fixuptri2, 1);
12291 lprevself(fixuptri);
12293 delaunayfixup(m, b, &fixuptri, 0);
12297 oprevself(fixuptri);
12300 tspivot(fixuptri, crosssubseg);
12301 if (crosssubseg.ss == m->dummysub) {
12302 flip(m, b, &fixuptri);
12307 segmentintersection(m, b, &fixuptri, &crosssubseg, endpoint2);
12314 insertsubseg(m, b, &fixuptri, newmark);
12319 if (!scoutsegment(m, b, &fixuptri, endpoint2, newmark)) {
12320 constrainededge(m, b, &fixuptri, endpoint2, newmark);
12331#ifdef ANSI_DECLARATORS
12333 vertex endpoint1, vertex endpoint2,
int newmark)
12335void insertsegment(m, b, endpoint1, endpoint2, newmark)
12344 struct otri searchtri1, searchtri2;
12345 triangle encodedtri;
12346 vertex checkvertex;
12349 if (b->verbose > 1) {
12350 printf(
" Connecting (%.12g, %.12g) to (%.12g, %.12g).\n",
12351 endpoint1[0], endpoint1[1], endpoint2[0], endpoint2[1]);
12355 checkvertex = (vertex) NULL;
12356 encodedtri = vertex2tri(endpoint1);
12357 if (encodedtri != (triangle) NULL) {
12358 decode(encodedtri, searchtri1);
12359 org(searchtri1, checkvertex);
12361 if (checkvertex != endpoint1) {
12363 searchtri1.tri = m->dummytri;
12364 searchtri1.orient = 0;
12365 symself(searchtri1);
12367 if (locate(m, b, endpoint1, &searchtri1) != ONVERTEX) {
12369 "Internal error in insertsegment(): Unable to locate PSLG vertex\n");
12370 printf(
" (%.12g, %.12g) in triangulation.\n",
12371 endpoint1[0], endpoint1[1]);
12376 otricopy(searchtri1, m->recenttri);
12379 if (scoutsegment(m, b, &searchtri1, endpoint2, newmark)) {
12385 org(searchtri1, endpoint1);
12388 checkvertex = (vertex) NULL;
12389 encodedtri = vertex2tri(endpoint2);
12390 if (encodedtri != (triangle) NULL) {
12391 decode(encodedtri, searchtri2);
12392 org(searchtri2, checkvertex);
12394 if (checkvertex != endpoint2) {
12396 searchtri2.tri = m->dummytri;
12397 searchtri2.orient = 0;
12398 symself(searchtri2);
12400 if (locate(m, b, endpoint2, &searchtri2) != ONVERTEX) {
12402 "Internal error in insertsegment(): Unable to locate PSLG vertex\n");
12403 printf(
" (%.12g, %.12g) in triangulation.\n",
12404 endpoint2[0], endpoint2[1]);
12409 otricopy(searchtri2, m->recenttri);
12412 if (scoutsegment(m, b, &searchtri2, endpoint1, newmark)) {
12418 org(searchtri2, endpoint2);
12424 conformingedge(m, b, endpoint1, endpoint2, newmark);
12429 constrainededge(m, b, &searchtri1, endpoint2, newmark);
12443#ifdef ANSI_DECLARATORS
12452 struct otri hulltri;
12453 struct otri nexttri;
12454 struct otri starttri;
12458 hulltri.tri = m->dummytri;
12459 hulltri.orient = 0;
12462 otricopy(hulltri, starttri);
12466 insertsubseg(m, b, &hulltri, 1);
12468 lnextself(hulltri);
12469 oprev(hulltri, nexttri);
12470 while (nexttri.tri != m->dummytri) {
12471 otricopy(nexttri, hulltri);
12472 oprev(hulltri, nexttri);
12474 }
while (!otriequal(hulltri, starttri));
12489#ifdef ANSI_DECLARATORS
12490void formskeleton(
struct mesh *m,
struct behavior *b,
int *segmentlist,
12491 int *segmentmarkerlist,
int numberofsegments)
12493void formskeleton(m, b, segmentlist, segmentmarkerlist, numberofsegments)
12497int *segmentmarkerlist;
12498int numberofsegments;
12503#ifdef ANSI_DECLARATORS
12505 FILE* polyfile,
char* polyfilename)
12507void formskeleton(m, b, polyfile, polyfilename)
12518 char polyfilename[6];
12521 char inputline[INPUTLINESIZE];
12524 vertex endpoint1, endpoint2;
12525 int segmentmarkers;
12532 printf(
"Recovering segments in Delaunay triangulation.\n");
12535 strcpy(polyfilename,
"input");
12536 m->insegments = numberofsegments;
12537 segmentmarkers = segmentmarkerlist != (
int *) NULL;
12542 stringptr = readline(inputline, polyfile, polyfilename);
12543 m->insegments = (int) strtol(stringptr, &stringptr, 0);
12544 stringptr = findfield(stringptr);
12545 if (*stringptr ==
'\0') {
12546 segmentmarkers = 0;
12548 segmentmarkers = (int) strtol(stringptr, &stringptr, 0);
12553 if (m->triangles.items == 0) {
12559 if (m->insegments > 0) {
12560 makevertexmap(m, b);
12562 printf(
" Recovering PSLG segments.\n");
12568 for (i = 0; i < m->insegments; i++) {
12570 end1 = segmentlist[index++];
12571 end2 = segmentlist[index++];
12572 if (segmentmarkers) {
12573 boundmarker = segmentmarkerlist[i];
12576 stringptr = readline(inputline, polyfile, b->inpolyfilename);
12577 stringptr = findfield(stringptr);
12578 if (*stringptr ==
'\0') {
12579 printf(
"Error: Segment %d has no endpoints in %s.\n",
12580 b->firstnumber + i, polyfilename);
12583 end1 = (int) strtol(stringptr, &stringptr, 0);
12585 stringptr = findfield(stringptr);
12586 if (*stringptr ==
'\0') {
12587 printf(
"Error: Segment %d is missing its second endpoint in %s.\n",
12588 b->firstnumber + i, polyfilename);
12591 end2 = (int) strtol(stringptr, &stringptr, 0);
12593 if (segmentmarkers) {
12594 stringptr = findfield(stringptr);
12595 if (*stringptr ==
'\0') {
12598 boundmarker = (int) strtol(stringptr, &stringptr, 0);
12602 if ((end1 < b->firstnumber) ||
12603 (end1 >= b->firstnumber + m->invertices)) {
12605 printf(
"Warning: Invalid first endpoint of segment %d in %s.\n",
12606 b->firstnumber + i, polyfilename);
12608 }
else if ((end2 < b->firstnumber) ||
12609 (end2 >= b->firstnumber + m->invertices)) {
12611 printf(
"Warning: Invalid second endpoint of segment %d in %s.\n",
12612 b->firstnumber + i, polyfilename);
12616 endpoint1 = getvertex(m, b, end1);
12617 endpoint2 = getvertex(m, b, end2);
12618 if ((endpoint1[0] == endpoint2[0]) && (endpoint1[1] == endpoint2[1])) {
12620 printf(
"Warning: Endpoints of segment %d are coincident in %s.\n",
12621 b->firstnumber + i, polyfilename);
12624 insertsegment(m, b, endpoint1, endpoint2, boundmarker);
12631 if (b->convex || !b->poly) {
12634 printf(
" Enclosing convex hull with segments.\n");
12656#ifdef ANSI_DECLARATORS
12659void infecthull(m, b)
12665 struct otri hulltri;
12666 struct otri nexttri;
12667 struct otri starttri;
12668 struct osub hullsubseg;
12669 triangle **deadtriangle;
12670 vertex horg, hdest;
12675 printf(
" Marking concavities (external triangles) for elimination.\n");
12678 hulltri.tri = m->dummytri;
12679 hulltri.orient = 0;
12682 otricopy(hulltri, starttri);
12686 if (!infected(hulltri)) {
12688 tspivot(hulltri, hullsubseg);
12689 if (hullsubseg.ss == m->dummysub) {
12691 if (!infected(hulltri)) {
12693 deadtriangle = (triangle **) poolalloc(&m->viri);
12694 *deadtriangle = hulltri.tri;
12698 if (mark(hullsubseg) == 0) {
12699 setmark(hullsubseg, 1);
12700 org(hulltri, horg);
12701 dest(hulltri, hdest);
12702 if (vertexmark(horg) == 0) {
12703 setvertexmark(horg, 1);
12705 if (vertexmark(hdest) == 0) {
12706 setvertexmark(hdest, 1);
12712 lnextself(hulltri);
12713 oprev(hulltri, nexttri);
12714 while (nexttri.tri != m->dummytri) {
12715 otricopy(nexttri, hulltri);
12716 oprev(hulltri, nexttri);
12718 }
while (!otriequal(hulltri, starttri));
12738#ifdef ANSI_DECLARATORS
12747 struct otri testtri;
12748 struct otri neighbor;
12749 triangle **virusloop;
12750 triangle **deadtriangle;
12751 struct osub neighborsubseg;
12753 vertex norg, ndest;
12754 vertex deadorg, deaddest, deadapex;
12760 printf(
" Marking neighbors of marked triangles.\n");
12764 traversalinit(&m->viri);
12765 virusloop = (triangle **) traverse(&m->viri);
12766 while (virusloop != (triangle **) NULL) {
12767 testtri.tri = *virusloop;
12773 if (b->verbose > 2) {
12776 testtri.orient = 0;
12777 org(testtri, deadorg);
12778 dest(testtri, deaddest);
12779 apex(testtri, deadapex);
12780 printf(
" Checking (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
12781 deadorg[0], deadorg[1], deaddest[0], deaddest[1],
12782 deadapex[0], deadapex[1]);
12785 for (testtri.orient = 0; testtri.orient < 3; testtri.orient++) {
12787 sym(testtri, neighbor);
12789 tspivot(testtri, neighborsubseg);
12791 if ((neighbor.tri == m->dummytri) || infected(neighbor)) {
12792 if (neighborsubseg.ss != m->dummysub) {
12796 subsegdealloc(m, neighborsubseg.ss);
12797 if (neighbor.tri != m->dummytri) {
12800 uninfect(neighbor);
12801 tsdissolve(neighbor);
12806 if (neighborsubseg.ss == m->dummysub) {
12809 if (b->verbose > 2) {
12810 org(neighbor, deadorg);
12811 dest(neighbor, deaddest);
12812 apex(neighbor, deadapex);
12814 " Marking (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
12815 deadorg[0], deadorg[1], deaddest[0], deaddest[1],
12816 deadapex[0], deadapex[1]);
12820 deadtriangle = (triangle **) poolalloc(&m->viri);
12821 *deadtriangle = neighbor.tri;
12824 stdissolve(neighborsubseg);
12826 if (mark(neighborsubseg) == 0) {
12827 setmark(neighborsubseg, 1);
12829 org(neighbor, norg);
12830 dest(neighbor, ndest);
12831 if (vertexmark(norg) == 0) {
12832 setvertexmark(norg, 1);
12834 if (vertexmark(ndest) == 0) {
12835 setvertexmark(ndest, 1);
12843 virusloop = (triangle **) traverse(&m->viri);
12847 printf(
" Deleting marked triangles.\n");
12850 traversalinit(&m->viri);
12851 virusloop = (triangle **) traverse(&m->viri);
12852 while (virusloop != (triangle **) NULL) {
12853 testtri.tri = *virusloop;
12858 for (testtri.orient = 0; testtri.orient < 3; testtri.orient++) {
12859 org(testtri, testvertex);
12861 if (testvertex != (vertex) NULL) {
12864 setorg(testtri, NULL);
12866 onext(testtri, neighbor);
12868 while ((neighbor.tri != m->dummytri) &&
12869 (!otriequal(neighbor, testtri))) {
12870 if (infected(neighbor)) {
12872 setorg(neighbor, NULL);
12878 onextself(neighbor);
12881 if (neighbor.tri == m->dummytri) {
12883 oprev(testtri, neighbor);
12885 while (neighbor.tri != m->dummytri) {
12886 if (infected(neighbor)) {
12888 setorg(neighbor, NULL);
12894 oprevself(neighbor);
12898 if (b->verbose > 1) {
12899 printf(
" Deleting vertex (%.12g, %.12g)\n",
12900 testvertex[0], testvertex[1]);
12902 setvertextype(testvertex, UNDEADVERTEX);
12910 for (testtri.orient = 0; testtri.orient < 3; testtri.orient++) {
12911 sym(testtri, neighbor);
12912 if (neighbor.tri == m->dummytri) {
12919 dissolve(neighbor);
12926 triangledealloc(m, testtri.tri);
12927 virusloop = (triangle **) traverse(&m->viri);
12930 poolrestart(&m->viri);
12948#ifdef ANSI_DECLARATORS
12950 REAL attribute, REAL area)
12952void regionplague(m, b, attribute, area)
12960 struct otri testtri;
12961 struct otri neighbor;
12962 triangle **virusloop;
12963 triangle **regiontri;
12964 struct osub neighborsubseg;
12965 vertex regionorg, regiondest, regionapex;
12969 if (b->verbose > 1) {
12970 printf(
" Marking neighbors of marked triangles.\n");
12975 traversalinit(&m->viri);
12976 virusloop = (triangle **) traverse(&m->viri);
12977 while (virusloop != (triangle **) NULL) {
12978 testtri.tri = *virusloop;
12984 if (b->regionattrib) {
12986 setelemattribute(testtri, m->eextras, attribute);
12990 setareabound(testtri, area);
12992 if (b->verbose > 2) {
12995 testtri.orient = 0;
12996 org(testtri, regionorg);
12997 dest(testtri, regiondest);
12998 apex(testtri, regionapex);
12999 printf(
" Checking (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
13000 regionorg[0], regionorg[1], regiondest[0], regiondest[1],
13001 regionapex[0], regionapex[1]);
13004 for (testtri.orient = 0; testtri.orient < 3; testtri.orient++) {
13006 sym(testtri, neighbor);
13008 tspivot(testtri, neighborsubseg);
13011 if ((neighbor.tri != m->dummytri) && !infected(neighbor)
13012 && (neighborsubseg.ss == m->dummysub)) {
13013 if (b->verbose > 2) {
13014 org(neighbor, regionorg);
13015 dest(neighbor, regiondest);
13016 apex(neighbor, regionapex);
13017 printf(
" Marking (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
13018 regionorg[0], regionorg[1], regiondest[0], regiondest[1],
13019 regionapex[0], regionapex[1]);
13024 regiontri = (triangle **) poolalloc(&m->viri);
13025 *regiontri = neighbor.tri;
13031 virusloop = (triangle **) traverse(&m->viri);
13035 if (b->verbose > 1) {
13036 printf(
" Unmarking marked triangles.\n");
13038 traversalinit(&m->viri);
13039 virusloop = (triangle **) traverse(&m->viri);
13040 while (virusloop != (triangle **) NULL) {
13041 testtri.tri = *virusloop;
13043 virusloop = (triangle **) traverse(&m->viri);
13046 poolrestart(&m->viri);
13061#ifdef ANSI_DECLARATORS
13062void carveholes(
struct mesh *m,
struct behavior *b, REAL *holelist,
int holes,
13063 REAL *regionlist,
int regions)
13065void carveholes(m, b, holelist, holes, regionlist, regions)
13075 struct otri searchtri;
13076 struct otri triangleloop;
13077 struct otri *regiontris;
13078 triangle **holetri;
13079 triangle **regiontri;
13080 vertex searchorg, searchdest;
13081 enum locateresult intersect;
13085 if (!(b->quiet || (b->noholes && b->convex))) {
13086 printf(
"Removing unwanted triangles.\n");
13087 if (b->verbose && (holes > 0)) {
13088 printf(
" Marking holes for elimination.\n");
13094 regiontris = (
struct otri *) trimalloc(regions *
13095 (
int)
sizeof(
struct otri));
13097 regiontris = (
struct otri *) NULL;
13100 if (((holes > 0) && !b->noholes) || !b->convex || (regions > 0)) {
13103 poolinit(&m->viri,
sizeof(triangle *), VIRUSPERBLOCK, VIRUSPERBLOCK, 0);
13112 if ((holes > 0) && !b->noholes) {
13114 for (i = 0; i < 2 * holes; i += 2) {
13116 if ((holelist[i] >= m->xmin) && (holelist[i] <= m->xmax)
13117 && (holelist[i + 1] >= m->ymin) && (holelist[i + 1] <= m->ymax)) {
13119 searchtri.tri = m->dummytri;
13120 searchtri.orient = 0;
13121 symself(searchtri);
13125 org(searchtri, searchorg);
13126 dest(searchtri, searchdest);
13127 if (counterclockwise(m, b, searchorg, searchdest, &holelist[i]) >
13130 intersect = locate(m, b, &holelist[i], &searchtri);
13131 if ((intersect != OUTSIDE) && (!infected(searchtri))) {
13135 holetri = (triangle **) poolalloc(&m->viri);
13136 *holetri = searchtri.tri;
13151 for (i = 0; i < regions; i++) {
13152 regiontris[i].tri = m->dummytri;
13154 if ((regionlist[4 * i] >= m->xmin) && (regionlist[4 * i] <= m->xmax) &&
13155 (regionlist[4 * i + 1] >= m->ymin) &&
13156 (regionlist[4 * i + 1] <= m->ymax)) {
13158 searchtri.tri = m->dummytri;
13159 searchtri.orient = 0;
13160 symself(searchtri);
13164 org(searchtri, searchorg);
13165 dest(searchtri, searchdest);
13166 if (counterclockwise(m, b, searchorg, searchdest, ®ionlist[4 * i]) >
13169 intersect = locate(m, b, ®ionlist[4 * i], &searchtri);
13170 if ((intersect != OUTSIDE) && (!infected(searchtri))) {
13173 otricopy(searchtri, regiontris[i]);
13180 if (m->viri.items > 0) {
13188 if (b->regionattrib) {
13190 printf(
"Spreading regional attributes and area constraints.\n");
13192 printf(
"Spreading regional attributes.\n");
13195 printf(
"Spreading regional area constraints.\n");
13198 if (b->regionattrib && !b->refine) {
13200 traversalinit(&m->triangles);
13201 triangleloop.orient = 0;
13202 triangleloop.tri = triangletraverse(m);
13203 while (triangleloop.tri != (triangle *) NULL) {
13204 setelemattribute(triangleloop, m->eextras, 0.0);
13205 triangleloop.tri = triangletraverse(m);
13208 for (i = 0; i < regions; i++) {
13209 if (regiontris[i].tri != m->dummytri) {
13212 if (!deadtri(regiontris[i].tri)) {
13214 infect(regiontris[i]);
13215 regiontri = (triangle **) poolalloc(&m->viri);
13216 *regiontri = regiontris[i].tri;
13218 regionplague(m, b, regionlist[4 * i + 2], regionlist[4 * i + 3]);
13223 if (b->regionattrib && !b->refine) {
13230 if (((holes > 0) && !b->noholes) || !b->convex || (regions > 0)) {
13231 pooldeinit(&m->viri);
13234 trifree((VOID *) regiontris);
13255#ifdef ANSI_DECLARATORS
13258void tallyencs(m, b)
13264 struct osub subsegloop;
13267 traversalinit(&m->subsegs);
13268 subsegloop.ssorient = 0;
13269 subsegloop.ss = subsegtraverse(m);
13270 while (subsegloop.ss != (subseg *) NULL) {
13273 checkseg4encroach(m, b, &subsegloop);
13274 subsegloop.ss = subsegtraverse(m);
13288void precisionerror()
13290 printf(
"Try increasing the area criterion and/or reducing the minimum\n");
13291 printf(
" allowable angle so that tiny triangles are not created.\n");
13293 printf(
"Alternatively, try recompiling me with double precision\n");
13294 printf(
" arithmetic (by removing \"#define SINGLE\" from the\n");
13295 printf(
" source file or \"-DSINGLE\" from the makefile).\n");
13317#ifdef ANSI_DECLARATORS
13318void splitencsegs(
struct mesh *m,
struct behavior *b,
int triflaws)
13320void splitencsegs(m, b, triflaws)
13327 struct otri enctri;
13328 struct otri testtri;
13329 struct osub testsh;
13330 struct osub currentenc;
13332 vertex eorg, edest, eapex;
13334 enum insertvertexresult success;
13335 REAL segmentlength, nearestpoweroftwo;
13337 REAL multiplier, divisor;
13338 int acuteorg, acuteorg2, acutedest, acutedest2;
13346 while ((m->badsubsegs.items > 0) && (m->steinerleft != 0)) {
13347 traversalinit(&m->badsubsegs);
13348 encloop = badsubsegtraverse(m);
13349 while ((encloop != (
struct badsubseg *) NULL) && (m->steinerleft != 0)) {
13350 sdecode(encloop->encsubseg, currentenc);
13351 sorg(currentenc, eorg);
13352 sdest(currentenc, edest);
13357 if (!deadsubseg(currentenc.ss) &&
13358 (eorg == encloop->subsegorg) && (edest == encloop->subsegdest)) {
13375 stpivot(currentenc, enctri);
13376 lnext(enctri, testtri);
13377 tspivot(testtri, testsh);
13378 acuteorg = testsh.ss != m->dummysub;
13380 lnextself(testtri);
13381 tspivot(testtri, testsh);
13382 acutedest = testsh.ss != m->dummysub;
13387 if (!b->conformdel && !acuteorg && !acutedest) {
13388 apex(enctri, eapex);
13389 while ((vertextype(eapex) == FREEVERTEX) &&
13390 ((eorg[0] - eapex[0]) * (edest[0] - eapex[0]) +
13391 (eorg[1] - eapex[1]) * (edest[1] - eapex[1]) < 0.0)) {
13392 deletevertex(m, b, &testtri);
13393 stpivot(currentenc, enctri);
13394 apex(enctri, eapex);
13395 lprev(enctri, testtri);
13401 sym(enctri, testtri);
13402 if (testtri.tri != m->dummytri) {
13404 lnextself(testtri);
13405 tspivot(testtri, testsh);
13406 acutedest2 = testsh.ss != m->dummysub;
13407 acutedest = acutedest || acutedest2;
13409 lnextself(testtri);
13410 tspivot(testtri, testsh);
13411 acuteorg2 = testsh.ss != m->dummysub;
13412 acuteorg = acuteorg || acuteorg2;
13415 if (!b->conformdel && !acuteorg2 && !acutedest2) {
13416 org(testtri, eapex);
13417 while ((vertextype(eapex) == FREEVERTEX) &&
13418 ((eorg[0] - eapex[0]) * (edest[0] - eapex[0]) +
13419 (eorg[1] - eapex[1]) * (edest[1] - eapex[1]) < 0.0)) {
13420 deletevertex(m, b, &testtri);
13421 sym(enctri, testtri);
13422 apex(testtri, eapex);
13423 lprevself(testtri);
13430 if (acuteorg || acutedest) {
13431 segmentlength = sqrt((edest[0] - eorg[0]) * (edest[0] - eorg[0]) +
13432 (edest[1] - eorg[1]) * (edest[1] - eorg[1]));
13435 nearestpoweroftwo = 1.0;
13436 while (segmentlength > 3.0 * nearestpoweroftwo) {
13437 nearestpoweroftwo *= 2.0;
13439 while (segmentlength < 1.5 * nearestpoweroftwo) {
13440 nearestpoweroftwo *= 0.5;
13443 split = nearestpoweroftwo / segmentlength;
13445 split = 1.0 - split;
13454 newvertex = (vertex) poolalloc(&m->vertices);
13456 for (i = 0; i < 2 + m->nextras; i++) {
13457 newvertex[i] = eorg[i] + split * (edest[i] - eorg[i]);
13464 multiplier = counterclockwise(m, b, eorg, edest, newvertex);
13465 divisor = ((eorg[0] - edest[0]) * (eorg[0] - edest[0]) +
13466 (eorg[1] - edest[1]) * (eorg[1] - edest[1]));
13467 if ((multiplier != 0.0) && (divisor != 0.0)) {
13468 multiplier = multiplier / divisor;
13470 if (multiplier == multiplier) {
13471 newvertex[0] += multiplier * (edest[1] - eorg[1]);
13472 newvertex[1] += multiplier * (eorg[0] - edest[0]);
13477 setvertexmark(newvertex, mark(currentenc));
13478 setvertextype(newvertex, SEGMENTVERTEX);
13479 if (b->verbose > 1) {
13481 " Splitting subsegment (%.12g, %.12g) (%.12g, %.12g) at (%.12g, %.12g).\n",
13482 eorg[0], eorg[1], edest[0], edest[1],
13483 newvertex[0], newvertex[1]);
13486 if (((newvertex[0] == eorg[0]) && (newvertex[1] == eorg[1])) ||
13487 ((newvertex[0] == edest[0]) && (newvertex[1] == edest[1]))) {
13488 printf(
"Error: Ran out of precision at (%.12g, %.12g).\n",
13489 newvertex[0], newvertex[1]);
13490 printf(
"I attempted to split a segment to a smaller size than\n");
13491 printf(
" can be accommodated by the finite precision of\n");
13492 printf(
" floating point arithmetic.\n");
13497 success = insertvertex(m, b, newvertex, &enctri, ¤tenc,
13499 if ((success != SUCCESSFULVERTEX) && (success != ENCROACHINGVERTEX)) {
13500 printf(
"Internal error in splitencsegs():\n");
13501 printf(
" Failure to split a segment.\n");
13504 if (m->steinerleft > 0) {
13509 checkseg4encroach(m, b, ¤tenc);
13510 snextself(currentenc);
13512 checkseg4encroach(m, b, ¤tenc);
13515 badsubsegdealloc(m, encloop);
13516 encloop = badsubsegtraverse(m);
13531#ifdef ANSI_DECLARATORS
13534void tallyfaces(m, b)
13540 struct otri triangleloop;
13543 printf(
" Making a list of bad triangles.\n");
13545 traversalinit(&m->triangles);
13546 triangleloop.orient = 0;
13547 triangleloop.tri = triangletraverse(m);
13548 while (triangleloop.tri != (triangle *) NULL) {
13550 testtriangle(m, b, &triangleloop);
13551 triangleloop.tri = triangletraverse(m);
13567#ifdef ANSI_DECLARATORS
13571void splittriangle(m, b, badtri)
13578 struct otri badotri;
13579 vertex borg, bdest, bapex;
13582 enum insertvertexresult success;
13586 decode(badtri->poortri, badotri);
13587 org(badotri, borg);
13588 dest(badotri, bdest);
13589 apex(badotri, bapex);
13593 if (!deadtri(badotri.tri) && (borg == badtri->triangorg) &&
13594 (bdest == badtri->triangdest) && (bapex == badtri->triangapex)) {
13595 if (b->verbose > 1) {
13596 printf(
" Splitting this triangle at its circumcenter:\n");
13597 printf(
" (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n", borg[0],
13598 borg[1], bdest[0], bdest[1], bapex[0], bapex[1]);
13603 newvertex = (vertex) poolalloc(&m->vertices);
13604 findcircumcenter(m, b, borg, bdest, bapex, newvertex, &xi, &eta, 1);
13607 if (((newvertex[0] == borg[0]) && (newvertex[1] == borg[1])) ||
13608 ((newvertex[0] == bdest[0]) && (newvertex[1] == bdest[1])) ||
13609 ((newvertex[0] == bapex[0]) && (newvertex[1] == bapex[1]))) {
13612 "Warning: New vertex (%.12g, %.12g) falls on existing vertex.\n",
13613 newvertex[0], newvertex[1]);
13616 vertexdealloc(m, newvertex);
13618 for (i = 2; i < 2 + m->nextras; i++) {
13620 newvertex[i] = borg[i] + xi * (bdest[i] - borg[i])
13621 + eta * (bapex[i] - borg[i]);
13625 setvertexmark(newvertex, 0);
13626 setvertextype(newvertex, FREEVERTEX);
13636 lprevself(badotri);
13641 success = insertvertex(m, b, newvertex, &badotri, (
struct osub *) NULL,
13643 if (success == SUCCESSFULVERTEX) {
13644 if (m->steinerleft > 0) {
13647 }
else if (success == ENCROACHINGVERTEX) {
13651 if (b->verbose > 1) {
13652 printf(
" Rejecting (%.12g, %.12g).\n", newvertex[0], newvertex[1]);
13654 vertexdealloc(m, newvertex);
13655 }
else if (success == VIOLATINGVERTEX) {
13658 vertexdealloc(m, newvertex);
13663 "Warning: New vertex (%.12g, %.12g) falls on existing vertex.\n",
13664 newvertex[0], newvertex[1]);
13667 vertexdealloc(m, newvertex);
13672 printf(
" The new vertex is at the circumcenter of triangle\n");
13673 printf(
" (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
13674 borg[0], borg[1], bdest[0], bdest[1], bapex[0], bapex[1]);
13676 printf(
"This probably means that I am trying to refine triangles\n");
13677 printf(
" to a smaller size than can be accommodated by the finite\n");
13678 printf(
" precision of floating point arithmetic. (You can be\n");
13679 printf(
" sure of this if I fail to terminate.)\n");
13696#ifdef ANSI_DECLARATORS
13697void enforcequality(
struct mesh *m,
struct behavior *b)
13699void enforcequality(m, b)
13709 printf(
"Adding Steiner points to enforce quality.\n");
13712 poolinit(&m->badsubsegs,
sizeof(
struct badsubseg), BADSUBSEGPERBLOCK,
13713 BADSUBSEGPERBLOCK, 0);
13715 printf(
" Looking for encroached subsegments.\n");
13719 if (b->verbose && (m->badsubsegs.items > 0)) {
13720 printf(
" Splitting encroached subsegments.\n");
13723 splitencsegs(m, b, 0);
13728 if ((b->minangle > 0.0) || b->vararea || b->fixedarea || b->usertest) {
13730 poolinit(&m->badtriangles,
sizeof(
struct badtriang), BADTRIPERBLOCK,
13731 BADTRIPERBLOCK, 0);
13733 for (i = 0; i < 4096; i++) {
13734 m->queuefront[i] = (
struct badtriang *) NULL;
13736 m->firstnonemptyq = -1;
13740 poolinit(&m->flipstackers,
sizeof(
struct flipstacker), FLIPSTACKERPERBLOCK,
13741 FLIPSTACKERPERBLOCK, 0);
13742 m->checkquality = 1;
13744 printf(
" Splitting bad triangles.\n");
13746 while ((m->badtriangles.items > 0) && (m->steinerleft != 0)) {
13748 badtri = dequeuebadtriang(m);
13749 splittriangle(m, b, badtri);
13750 if (m->badsubsegs.items > 0) {
13752 enqueuebadtriang(m, b, badtri);
13755 splitencsegs(m, b, 1);
13758 pooldealloc(&m->badtriangles, (VOID *) badtri);
13767 if (!b->quiet && b->conformdel && (m->badsubsegs.items > 0) &&
13768 (m->steinerleft == 0)) {
13769 printf(
"\nWarning: I ran out of Steiner points, but the mesh has\n");
13770 if (m->badsubsegs.items == 1) {
13771 printf(
" one encroached subsegment, and therefore might not be truly\n"
13774 printf(
" %ld encroached subsegments, and therefore might not be truly\n"
13775 , m->badsubsegs.items);
13777 printf(
" Delaunay. If the Delaunay property is important to you,\n");
13778 printf(
" try increasing the number of Steiner points (controlled by\n");
13779 printf(
" the -S switch) slightly and try again.\n\n");
13795#ifdef ANSI_DECLARATORS
13798void highorder(m, b)
13804 struct otri triangleloop, trisym;
13805 struct osub checkmark;
13807 vertex torg, tdest;
13813 printf(
"Adding vertices for second-order triangles.\n");
13820 m->vertices.deaditemstack = (VOID *) NULL;
13822 traversalinit(&m->triangles);
13823 triangleloop.tri = triangletraverse(m);
13830 while (triangleloop.tri != (triangle *) NULL) {
13831 for (triangleloop.orient = 0; triangleloop.orient < 3;
13832 triangleloop.orient++) {
13833 sym(triangleloop, trisym);
13834 if ((triangleloop.tri < trisym.tri) || (trisym.tri == m->dummytri)) {
13835 org(triangleloop, torg);
13836 dest(triangleloop, tdest);
13839 newvertex = (vertex) poolalloc(&m->vertices);
13840 for (i = 0; i < 2 + m->nextras; i++) {
13841 newvertex[i] = 0.5 * (torg[i] + tdest[i]);
13845 setvertexmark(newvertex, trisym.tri == m->dummytri);
13846 setvertextype(newvertex,
13847 trisym.tri == m->dummytri ? FREEVERTEX : SEGMENTVERTEX);
13848 if (b->usesegments) {
13849 tspivot(triangleloop, checkmark);
13851 if (checkmark.ss != m->dummysub) {
13852 setvertexmark(newvertex, mark(checkmark));
13853 setvertextype(newvertex, SEGMENTVERTEX);
13856 if (b->verbose > 1) {
13857 printf(
" Creating (%.12g, %.12g).\n", newvertex[0], newvertex[1]);
13860 triangleloop.tri[m->highorderindex + triangleloop.orient] =
13861 (triangle) newvertex;
13862 if (trisym.tri != m->dummytri) {
13863 trisym.tri[m->highorderindex + trisym.orient] = (triangle) newvertex;
13867 triangleloop.tri = triangletraverse(m);
13886#ifdef ANSI_DECLARATORS
13887char *readline(
char *
string, FILE *infile,
char *infilename)
13889char *readline(
string, infile, infilename)
13900 result = fgets(
string, INPUTLINESIZE, infile);
13901 if (result == (
char *) NULL) {
13902 printf(
" Error: Unexpected end of file in %s.\n", infilename);
13907 while ((*result !=
'\0') && (*result !=
'#')
13908 && (*result !=
'.') && (*result !=
'+') && (*result !=
'-')
13909 && ((*result <
'0') || (*result >
'9'))) {
13913 }
while ((*result ==
'#') || (*result ==
'\0'));
13930#ifdef ANSI_DECLARATORS
13931char *findfield(
char *
string)
13933char *findfield(
string)
13942 while ((*result !=
'\0') && (*result !=
'#')
13943 && (*result !=
' ') && (*result !=
'\t')) {
13948 while ((*result !=
'\0') && (*result !=
'#')
13949 && (*result !=
'.') && (*result !=
'+') && (*result !=
'-')
13950 && ((*result <
'0') || (*result >
'9'))) {
13954 if (*result ==
'#') {
13971#ifdef ANSI_DECLARATORS
13972void readnodes(
struct mesh *m,
struct behavior *b,
char *nodefilename,
13973 char* polyfilename, FILE **polyfile)
13975void readnodes(m, b, nodefilename, polyfilename, polyfile)
13986 char inputline[INPUTLINESIZE];
13998 printf(
"Opening %s.\n", polyfilename);
14000 * polyfile = fopen(polyfilename,
"r");
14001 if (*polyfile == (FILE *) NULL) {
14002 printf(
" Error: Cannot access file %s.\n", polyfilename);
14007 stringptr = readline(inputline,* polyfile, polyfilename);
14008 m->invertices = (int) strtol(stringptr, &stringptr, 0);
14009 stringptr = findfield(stringptr);
14010 if (*stringptr ==
'\0') {
14013 m->mesh_dim = (int) strtol(stringptr, &stringptr, 0);
14015 stringptr = findfield(stringptr);
14016 if (*stringptr ==
'\0') {
14019 m->nextras = (int) strtol(stringptr, &stringptr, 0);
14021 stringptr = findfield(stringptr);
14022 if (*stringptr ==
'\0') {
14025 nodemarkers = (int) strtol(stringptr, &stringptr, 0);
14027 if (m->invertices > 0) {
14028 infile =* polyfile;
14029 infilename = polyfilename;
14030 m->readnodefile = 0;
14034 m->readnodefile = 1;
14035 infilename = nodefilename;
14038 m->readnodefile = 1;
14039 infilename = nodefilename;
14040 * polyfile = (FILE *) NULL;
14043 if (m->readnodefile) {
14046 printf(
"Opening %s.\n", nodefilename);
14048 infile = fopen(nodefilename,
"r");
14049 if (infile == (FILE *) NULL) {
14050 printf(
" Error: Cannot access file %s.\n", nodefilename);
14055 stringptr = readline(inputline, infile, nodefilename);
14056 m->invertices = (int) strtol(stringptr, &stringptr, 0);
14057 stringptr = findfield(stringptr);
14058 if (*stringptr ==
'\0') {
14061 m->mesh_dim = (int) strtol(stringptr, &stringptr, 0);
14063 stringptr = findfield(stringptr);
14064 if (*stringptr ==
'\0') {
14067 m->nextras = (int) strtol(stringptr, &stringptr, 0);
14069 stringptr = findfield(stringptr);
14070 if (*stringptr ==
'\0') {
14073 nodemarkers = (int) strtol(stringptr, &stringptr, 0);
14077 if (m->invertices < 3) {
14078 printf(
"Error: Input must have at least three input vertices.\n");
14081 if (m->mesh_dim != 2) {
14082 printf(
"Error: Triangle only works with two-dimensional meshes.\n");
14085 if (m->nextras == 0) {
14089 initializevertexpool(m, b);
14092 for (i = 0; i < m->invertices; i++) {
14093 vertexloop = (vertex) poolalloc(&m->vertices);
14094 stringptr = readline(inputline, infile, infilename);
14096 firstnode = (int) strtol(stringptr, &stringptr, 0);
14097 if ((firstnode == 0) || (firstnode == 1)) {
14098 b->firstnumber = firstnode;
14101 stringptr = findfield(stringptr);
14102 if (*stringptr ==
'\0') {
14103 printf(
"Error: Vertex %d has no x coordinate.\n", b->firstnumber + i);
14106 x = (REAL) strtod(stringptr, &stringptr);
14107 stringptr = findfield(stringptr);
14108 if (*stringptr ==
'\0') {
14109 printf(
"Error: Vertex %d has no y coordinate.\n", b->firstnumber + i);
14112 y = (REAL) strtod(stringptr, &stringptr);
14116 for (j = 2; j < 2 + m->nextras; j++) {
14117 stringptr = findfield(stringptr);
14118 if (*stringptr ==
'\0') {
14119 vertexloop[j] = 0.0;
14121 vertexloop[j] = (REAL) strtod(stringptr, &stringptr);
14126 stringptr = findfield(stringptr);
14127 if (*stringptr ==
'\0') {
14128 setvertexmark(vertexloop, 0);
14130 currentmarker = (int) strtol(stringptr, &stringptr, 0);
14131 setvertexmark(vertexloop, currentmarker);
14135 setvertexmark(vertexloop, 0);
14137 setvertextype(vertexloop, INPUTVERTEX);
14140 m->xmin = m->xmax = x;
14141 m->ymin = m->ymax = y;
14143 m->xmin = (x < m->xmin) ? x : m->xmin;
14144 m->xmax = (x > m->xmax) ? x : m->xmax;
14145 m->ymin = (y < m->ymin) ? y : m->ymin;
14146 m->ymax = (y > m->ymax) ? y : m->ymax;
14149 if (m->readnodefile) {
14155 m->xminextreme = 10 * m->xmin - 9 * m->xmax;
14168#ifdef ANSI_DECLARATORS
14169void transfernodes(
struct mesh *m,
struct behavior *b, REAL* pointlist,
14170 REAL* pointattriblist,
int* pointmarkerlist,
14171 int numberofpoints,
int numberofpointattribs)
14173void transfernodes(m, b, pointlist, pointattriblist, pointmarkerlist,
14174 numberofpoints, numberofpointattribs)
14178REAL* pointattriblist;
14179int* pointmarkerlist;
14181int numberofpointattribs;
14191 m->invertices = numberofpoints;
14193 m->nextras = numberofpointattribs;
14194 m->readnodefile = 0;
14195 if (m->invertices < 3) {
14196 printf(
"Error: Input must have at least three input vertices.\n");
14199 if (m->nextras == 0) {
14203 initializevertexpool(m, b);
14208 for (i = 0; i < m->invertices; i++) {
14209 vertexloop = (vertex) poolalloc(&m->vertices);
14211 x = vertexloop[0] = pointlist[coordindex++];
14212 y = vertexloop[1] = pointlist[coordindex++];
14214 for (j = 0; j < numberofpointattribs; j++) {
14215 vertexloop[2 + j] = pointattriblist[attribindex++];
14217 if (pointmarkerlist != (
int *) NULL) {
14219 setvertexmark(vertexloop, pointmarkerlist[i]);
14222 setvertexmark(vertexloop, 0);
14224 setvertextype(vertexloop, INPUTVERTEX);
14227 m->xmin = m->xmax = x;
14228 m->ymin = m->ymax = y;
14230 m->xmin = (x < m->xmin) ? x : m->xmin;
14231 m->xmax = (x > m->xmax) ? x : m->xmax;
14232 m->ymin = (y < m->ymin) ? y : m->ymin;
14233 m->ymax = (y > m->ymax) ? y : m->ymax;
14239 m->xminextreme = 10 * m->xmin - 9 * m->xmax;
14253#ifdef ANSI_DECLARATORS
14255 FILE* polyfile,
char* polyfilename, REAL **hlist,
int *holes,
14256 REAL **rlist,
int *regions)
14258void readholes(m, b, polyfile, polyfilename, hlist, holes, rlist, regions)
14272 char inputline[INPUTLINESIZE];
14278 stringptr = readline(inputline, polyfile, polyfilename);
14279 *holes = (int) strtol(stringptr, &stringptr, 0);
14281 holelist = (REAL *) trimalloc(2 * *holes * (
int)
sizeof(REAL));
14283 for (i = 0; i < 2 * *holes; i += 2) {
14284 stringptr = readline(inputline, polyfile, polyfilename);
14285 stringptr = findfield(stringptr);
14286 if (*stringptr ==
'\0') {
14287 printf(
"Error: Hole %d has no x coordinate.\n",
14288 b->firstnumber + (i >> 1));
14291 holelist[i] = (REAL) strtod(stringptr, &stringptr);
14293 stringptr = findfield(stringptr);
14294 if (*stringptr ==
'\0') {
14295 printf(
"Error: Hole %d has no y coordinate.\n",
14296 b->firstnumber + (i >> 1));
14299 holelist[i + 1] = (REAL) strtod(stringptr, &stringptr);
14303 *hlist = (REAL *) NULL;
14307 if ((b->regionattrib || b->vararea) && !b->refine) {
14309 stringptr = readline(inputline, polyfile, polyfilename);
14310 *regions = (int) strtol(stringptr, &stringptr, 0);
14311 if (*regions > 0) {
14312 regionlist = (REAL *) trimalloc(4 * *regions * (
int)
sizeof(REAL));
14313 *rlist = regionlist;
14315 for (i = 0; i < *regions; i++) {
14316 stringptr = readline(inputline, polyfile, polyfilename);
14317 stringptr = findfield(stringptr);
14318 if (*stringptr ==
'\0') {
14319 printf(
"Error: Region %d has no x coordinate.\n",
14320 b->firstnumber + i);
14323 regionlist[index++] = (REAL) strtod(stringptr, &stringptr);
14325 stringptr = findfield(stringptr);
14326 if (*stringptr ==
'\0') {
14327 printf(
"Error: Region %d has no y coordinate.\n",
14328 b->firstnumber + i);
14331 regionlist[index++] = (REAL) strtod(stringptr, &stringptr);
14333 stringptr = findfield(stringptr);
14334 if (*stringptr ==
'\0') {
14336 "Error: Region %d has no region attribute or area constraint.\n",
14337 b->firstnumber + i);
14340 regionlist[index++] = (REAL) strtod(stringptr, &stringptr);
14342 stringptr = findfield(stringptr);
14343 if (*stringptr ==
'\0') {
14344 regionlist[index] = regionlist[index - 1];
14346 regionlist[index] = (REAL) strtod(stringptr, &stringptr);
14354 *rlist = (REAL *) NULL;
14372#ifdef ANSI_DECLARATORS
14373void finishfile(FILE *outfile,
int argc,
char **argv)
14375void finishfile(outfile, argc, argv)
14384 fprintf(outfile,
"# Generated by");
14385 for (i = 0; i < argc; i++) {
14386 fprintf(outfile,
" ");
14387 fputs(argv[i], outfile);
14389 fprintf(outfile,
"\n");
14406#ifdef ANSI_DECLARATORS
14407void writenodes(
struct mesh *m,
struct behavior *b, REAL **pointlist,
14408 REAL **pointattriblist,
int **pointmarkerlist)
14410void writenodes(m, b, pointlist, pointattriblist, pointmarkerlist)
14414REAL **pointattriblist;
14415int **pointmarkerlist;
14420#ifdef ANSI_DECLARATORS
14421void writenodes(
struct mesh *m,
struct behavior *b,
char *nodefilename,
14422 int argc,
char **argv)
14424void writenodes(m, b, nodefilename, argc, argv)
14450 outvertices = m->vertices.items - m->undeads;
14452 outvertices = m->vertices.items;
14457 printf(
"Writing vertices.\n");
14460 if (*pointlist == (REAL *) NULL) {
14461 * pointlist = (REAL *) trimalloc((
int) (outvertices * 2 *
sizeof(REAL)));
14464 if ((m->nextras > 0) && (*pointattriblist == (REAL *) NULL)) {
14465 * pointattriblist = (REAL *) trimalloc((
int) (outvertices * m->nextras *
14469 if (!b->nobound && (*pointmarkerlist == (
int *) NULL)) {
14470 * pointmarkerlist = (
int *) trimalloc((
int) (outvertices *
sizeof(int)));
14472 plist =* pointlist;
14473 palist =* pointattriblist;
14474 pmlist =* pointmarkerlist;
14479 printf(
"Writing %s.\n", nodefilename);
14481 outfile = fopen(nodefilename,
"w");
14482 if (outfile == (FILE *) NULL) {
14483 printf(
" Error: Cannot create file %s.\n", nodefilename);
14488 fprintf(outfile,
"%ld %d %d %d\n", outvertices, m->mesh_dim,
14489 m->nextras, 1 - b->nobound);
14492 traversalinit(&m->vertices);
14493 vertexnumber = b->firstnumber;
14494 vertexloop = vertextraverse(m);
14495 while (vertexloop != (vertex) NULL) {
14496 if (!b->jettison || (vertextype(vertexloop) != UNDEADVERTEX)) {
14499 plist[coordindex++] = vertexloop[0];
14500 plist[coordindex++] = vertexloop[1];
14502 for (i = 0; i < m->nextras; i++) {
14503 palist[attribindex++] = vertexloop[2 + i];
14507 pmlist[vertexnumber - b->firstnumber] = vertexmark(vertexloop);
14511 fprintf(outfile,
"%4d %.17g %.17g", vertexnumber, vertexloop[0],
14513 for (i = 0; i < m->nextras; i++) {
14515 fprintf(outfile,
" %.17g", vertexloop[i + 2]);
14518 fprintf(outfile,
"\n");
14521 fprintf(outfile,
" %d\n", vertexmark(vertexloop));
14525 setvertexmark(vertexloop, vertexnumber);
14528 vertexloop = vertextraverse(m);
14532 finishfile(outfile, argc, argv);
14546#ifdef ANSI_DECLARATORS
14549void numbernodes(m, b)
14558 traversalinit(&m->vertices);
14559 vertexnumber = b->firstnumber;
14560 vertexloop = vertextraverse(m);
14561 while (vertexloop != (vertex) NULL) {
14562 setvertexmark(vertexloop, vertexnumber);
14563 if (!b->jettison || (vertextype(vertexloop) != UNDEADVERTEX)) {
14566 vertexloop = vertextraverse(m);
14578#ifdef ANSI_DECLARATORS
14580 int **trianglelist, REAL **triangleattriblist)
14582void writeelements(m, b, trianglelist, triangleattriblist)
14586REAL **triangleattriblist;
14591#ifdef ANSI_DECLARATORS
14592void writeelements(
struct mesh *m,
struct behavior *b,
char *elefilename,
14593 int argc,
char **argv)
14595void writeelements(m, b, elefilename, argc, argv)
14614 struct otri triangleloop;
14616 vertex mid1, mid2, mid3;
14617 long elementnumber;
14622 printf(
"Writing triangles.\n");
14625 if (*trianglelist == (
int *) NULL) {
14626 *trianglelist = (
int *) trimalloc((
int) (m->triangles.items *
14627 ((b->order + 1) * (b->order + 2) /
14628 2) *
sizeof(
int)));
14631 if ((m->eextras > 0) && (*triangleattriblist == (REAL *) NULL)) {
14632 *triangleattriblist = (REAL *) trimalloc((
int) (m->triangles.items *
14636 tlist = *trianglelist;
14637 talist = *triangleattriblist;
14642 printf(
"Writing %s.\n", elefilename);
14644 outfile = fopen(elefilename,
"w");
14645 if (outfile == (FILE *) NULL) {
14646 printf(
" Error: Cannot create file %s.\n", elefilename);
14650 fprintf(outfile,
"%ld %d %d\n", m->triangles.items,
14651 (b->order + 1) * (b->order + 2) / 2, m->eextras);
14654 traversalinit(&m->triangles);
14655 triangleloop.tri = triangletraverse(m);
14656 triangleloop.orient = 0;
14657 elementnumber = b->firstnumber;
14658 while (triangleloop.tri != (triangle *) NULL) {
14659 org(triangleloop, p1);
14660 dest(triangleloop, p2);
14661 apex(triangleloop, p3);
14662 if (b->order == 1) {
14664 tlist[vertexindex++] = vertexmark(p1);
14665 tlist[vertexindex++] = vertexmark(p2);
14666 tlist[vertexindex++] = vertexmark(p3);
14669 fprintf(outfile,
"%4ld %4d %4d %4d", elementnumber,
14670 vertexmark(p1), vertexmark(p2), vertexmark(p3));
14673 mid1 = (vertex) triangleloop.tri[m->highorderindex + 1];
14674 mid2 = (vertex) triangleloop.tri[m->highorderindex + 2];
14675 mid3 = (vertex) triangleloop.tri[m->highorderindex];
14677 tlist[vertexindex++] = vertexmark(p1);
14678 tlist[vertexindex++] = vertexmark(p2);
14679 tlist[vertexindex++] = vertexmark(p3);
14680 tlist[vertexindex++] = vertexmark(mid1);
14681 tlist[vertexindex++] = vertexmark(mid2);
14682 tlist[vertexindex++] = vertexmark(mid3);
14685 fprintf(outfile,
"%4ld %4d %4d %4d %4d %4d %4d", elementnumber,
14686 vertexmark(p1), vertexmark(p2), vertexmark(p3), vertexmark(mid1),
14687 vertexmark(mid2), vertexmark(mid3));
14692 for (i = 0; i < m->eextras; i++) {
14693 talist[attribindex++] = elemattribute(triangleloop, i);
14696 for (i = 0; i < m->eextras; i++) {
14697 fprintf(outfile,
" %.17g", elemattribute(triangleloop, i));
14699 fprintf(outfile,
"\n");
14702 triangleloop.tri = triangletraverse(m);
14707 finishfile(outfile, argc, argv);
14719#ifdef ANSI_DECLARATORS
14721 int **segmentlist,
int **segmentmarkerlist)
14723void writepoly(m, b, segmentlist, segmentmarkerlist)
14727int **segmentmarkerlist;
14732#ifdef ANSI_DECLARATORS
14733void writepoly(
struct mesh *m,
struct behavior *b,
char* polyfilename,
14734 REAL *holelist,
int holes, REAL *regionlist,
int regions,
14735 int argc,
char **argv)
14737void writepoly(m, b, polyfilename, holelist, holes, regionlist, regions,
14759 long holenumber, regionnumber;
14761 struct osub subsegloop;
14762 vertex endpoint1, endpoint2;
14767 printf(
"Writing segments.\n");
14770 if (*segmentlist == (
int *) NULL) {
14771 *segmentlist = (
int *) trimalloc((
int) (m->subsegs.items * 2 *
14775 if (!b->nobound && (*segmentmarkerlist == (
int *) NULL)) {
14776 *segmentmarkerlist = (
int *) trimalloc((
int) (m->subsegs.items *
14779 slist = *segmentlist;
14780 smlist = *segmentmarkerlist;
14784 printf(
"Writing %s.\n", polyfilename);
14786 outfile = fopen(polyfilename,
"w");
14787 if (outfile == (FILE *) NULL) {
14788 printf(
" Error: Cannot create file %s.\n", polyfilename);
14794 fprintf(outfile,
"%d %d %d %d\n", 0, m->mesh_dim, m->nextras,
14797 fprintf(outfile,
"%ld %d\n", m->subsegs.items, 1 - b->nobound);
14800 traversalinit(&m->subsegs);
14801 subsegloop.ss = subsegtraverse(m);
14802 subsegloop.ssorient = 0;
14803 subsegnumber = b->firstnumber;
14804 while (subsegloop.ss != (subseg *) NULL) {
14805 sorg(subsegloop, endpoint1);
14806 sdest(subsegloop, endpoint2);
14809 slist[index++] = vertexmark(endpoint1);
14810 slist[index++] = vertexmark(endpoint2);
14813 smlist[subsegnumber - b->firstnumber] = mark(subsegloop);
14818 fprintf(outfile,
"%4ld %4d %4d\n", subsegnumber,
14819 vertexmark(endpoint1), vertexmark(endpoint2));
14821 fprintf(outfile,
"%4ld %4d %4d %4d\n", subsegnumber,
14822 vertexmark(endpoint1), vertexmark(endpoint2), mark(subsegloop));
14826 subsegloop.ss = subsegtraverse(m);
14832 fprintf(outfile,
"%d\n", holes);
14834 for (holenumber = 0; holenumber < holes; holenumber++) {
14836 fprintf(outfile,
"%4ld %.17g %.17g\n", b->firstnumber + holenumber,
14837 holelist[2 * holenumber], holelist[2 * holenumber + 1]);
14841 fprintf(outfile,
"%d\n", regions);
14842 for (regionnumber = 0; regionnumber < regions; regionnumber++) {
14844 fprintf(outfile,
"%4ld %.17g %.17g %.17g %.17g\n",
14845 b->firstnumber + regionnumber,
14846 regionlist[4 * regionnumber], regionlist[4 * regionnumber + 1],
14847 regionlist[4 * regionnumber + 2],
14848 regionlist[4 * regionnumber + 3]);
14853 finishfile(outfile, argc, argv);
14865#ifdef ANSI_DECLARATORS
14867 int **edgelist,
int **edgemarkerlist)
14869void writeedges(m, b, edgelist, edgemarkerlist)
14873int **edgemarkerlist;
14878#ifdef ANSI_DECLARATORS
14879void writeedges(
struct mesh *m,
struct behavior *b,
char *edgefilename,
14880 int argc,
char **argv)
14882void writeedges(m, b, edgefilename, argc, argv)
14900 struct otri triangleloop, trisym;
14901 struct osub checkmark;
14909 printf(
"Writing edges.\n");
14912 if (*edgelist == (
int *) NULL) {
14913 *edgelist = (
int *) trimalloc((
int) (m->edges * 2 *
sizeof(int)));
14916 if (!b->nobound && (*edgemarkerlist == (
int *) NULL)) {
14917 *edgemarkerlist = (
int *) trimalloc((
int) (m->edges *
sizeof(int)));
14920 emlist = *edgemarkerlist;
14924 printf(
"Writing %s.\n", edgefilename);
14926 outfile = fopen(edgefilename,
"w");
14927 if (outfile == (FILE *) NULL) {
14928 printf(
" Error: Cannot create file %s.\n", edgefilename);
14932 fprintf(outfile,
"%ld %d\n", m->edges, 1 - b->nobound);
14935 traversalinit(&m->triangles);
14936 triangleloop.tri = triangletraverse(m);
14937 edgenumber = b->firstnumber;
14944 while (triangleloop.tri != (triangle *) NULL) {
14945 for (triangleloop.orient = 0; triangleloop.orient < 3;
14946 triangleloop.orient++) {
14947 sym(triangleloop, trisym);
14948 if ((triangleloop.tri < trisym.tri) || (trisym.tri == m->dummytri)) {
14949 org(triangleloop, p1);
14950 dest(triangleloop, p2);
14952 elist[index++] = vertexmark(p1);
14953 elist[index++] = vertexmark(p2);
14958 fprintf(outfile,
"%4ld %d %d\n", edgenumber,
14959 vertexmark(p1), vertexmark(p2));
14964 if (b->usesegments) {
14965 tspivot(triangleloop, checkmark);
14966 if (checkmark.ss == m->dummysub) {
14968 emlist[edgenumber - b->firstnumber] = 0;
14970 fprintf(outfile,
"%4ld %d %d %d\n", edgenumber,
14971 vertexmark(p1), vertexmark(p2), 0);
14975 emlist[edgenumber - b->firstnumber] = mark(checkmark);
14977 fprintf(outfile,
"%4ld %d %d %d\n", edgenumber,
14978 vertexmark(p1), vertexmark(p2), mark(checkmark));
14983 emlist[edgenumber - b->firstnumber] = trisym.tri == m->dummytri;
14985 fprintf(outfile,
"%4ld %d %d %d\n", edgenumber,
14986 vertexmark(p1), vertexmark(p2), trisym.tri == m->dummytri);
14993 triangleloop.tri = triangletraverse(m);
14997 finishfile(outfile, argc, argv);
15019#ifdef ANSI_DECLARATORS
15020void writevoronoi(
struct mesh *m,
struct behavior *b, REAL **vpointlist,
15021 REAL **vpointattriblist,
int **vpointmarkerlist,
15022 int **vedgelist,
int **vedgemarkerlist, REAL **vnormlist)
15024void writevoronoi(m, b, vpointlist, vpointattriblist, vpointmarkerlist,
15025 vedgelist, vedgemarkerlist, vnormlist)
15029REAL **vpointattriblist;
15030int **vpointmarkerlist;
15032int **vedgemarkerlist;
15038#ifdef ANSI_DECLARATORS
15039void writevoronoi(
struct mesh *m,
struct behavior *b,
char *vnodefilename,
15040 char *vedgefilename,
int argc,
char **argv)
15042void writevoronoi(m, b, vnodefilename, vedgefilename, argc, argv)
15045char *vnodefilename;
15046char *vedgefilename;
15064 struct otri triangleloop, trisym;
15065 vertex torg, tdest, tapex;
15066 REAL circumcenter[2];
15068 long vnodenumber, vedgenumber;
15075 printf(
"Writing Voronoi vertices.\n");
15078 if (*vpointlist == (REAL *) NULL) {
15079 *vpointlist = (REAL *) trimalloc((
int) (m->triangles.items * 2 *
15083 if (*vpointattriblist == (REAL *) NULL) {
15084 *vpointattriblist = (REAL *) trimalloc((
int) (m->triangles.items *
15085 m->nextras *
sizeof(REAL)));
15087 *vpointmarkerlist = (
int *) NULL;
15088 plist = *vpointlist;
15089 palist = *vpointattriblist;
15094 printf(
"Writing %s.\n", vnodefilename);
15096 outfile = fopen(vnodefilename,
"w");
15097 if (outfile == (FILE *) NULL) {
15098 printf(
" Error: Cannot create file %s.\n", vnodefilename);
15103 fprintf(outfile,
"%ld %d %d %d\n", m->triangles.items, 2, m->nextras, 0);
15106 traversalinit(&m->triangles);
15107 triangleloop.tri = triangletraverse(m);
15108 triangleloop.orient = 0;
15109 vnodenumber = b->firstnumber;
15110 while (triangleloop.tri != (triangle *) NULL) {
15111 org(triangleloop, torg);
15112 dest(triangleloop, tdest);
15113 apex(triangleloop, tapex);
15114 findcircumcenter(m, b, torg, tdest, tapex, circumcenter, &xi, &eta, 0);
15117 plist[coordindex++] = circumcenter[0];
15118 plist[coordindex++] = circumcenter[1];
15119 for (i = 2; i < 2 + m->nextras; i++) {
15121 palist[attribindex++] = torg[i] + xi * (tdest[i] - torg[i])
15122 + eta * (tapex[i] - torg[i]);
15126 fprintf(outfile,
"%4ld %.17g %.17g", vnodenumber, circumcenter[0],
15128 for (i = 2; i < 2 + m->nextras; i++) {
15130 fprintf(outfile,
" %.17g", torg[i] + xi * (tdest[i] - torg[i])
15131 + eta * (tapex[i] - torg[i]));
15133 fprintf(outfile,
"\n");
15136 * (
int *) (triangleloop.tri + 6) = (int) vnodenumber;
15137 triangleloop.tri = triangletraverse(m);
15142 finishfile(outfile, argc, argv);
15147 printf(
"Writing Voronoi edges.\n");
15150 if (*vedgelist == (
int *) NULL) {
15151 *vedgelist = (
int *) trimalloc((
int) (m->edges * 2 *
sizeof(int)));
15153 *vedgemarkerlist = (
int *) NULL;
15155 if (*vnormlist == (REAL *) NULL) {
15156 *vnormlist = (REAL *) trimalloc((
int) (m->edges * 2 *
sizeof(REAL)));
15158 elist = *vedgelist;
15159 normlist = *vnormlist;
15163 printf(
"Writing %s.\n", vedgefilename);
15165 outfile = fopen(vedgefilename,
"w");
15166 if (outfile == (FILE *) NULL) {
15167 printf(
" Error: Cannot create file %s.\n", vedgefilename);
15171 fprintf(outfile,
"%ld %d\n", m->edges, 0);
15174 traversalinit(&m->triangles);
15175 triangleloop.tri = triangletraverse(m);
15176 vedgenumber = b->firstnumber;
15183 while (triangleloop.tri != (triangle *) NULL) {
15184 for (triangleloop.orient = 0; triangleloop.orient < 3;
15185 triangleloop.orient++) {
15186 sym(triangleloop, trisym);
15187 if ((triangleloop.tri < trisym.tri) || (trisym.tri == m->dummytri)) {
15189 p1 = * (
int *) (triangleloop.tri + 6);
15190 if (trisym.tri == m->dummytri) {
15191 org(triangleloop, torg);
15192 dest(triangleloop, tdest);
15195 elist[coordindex] = p1;
15196 normlist[coordindex++] = tdest[1] - torg[1];
15197 elist[coordindex] = -1;
15198 normlist[coordindex++] = torg[0] - tdest[0];
15203 fprintf(outfile,
"%4ld %d %d %.17g %.17g\n", vedgenumber,
15204 p1, -1, tdest[1] - torg[1], torg[0] - tdest[0]);
15208 p2 = * (
int *) (trisym.tri + 6);
15211 elist[coordindex] = p1;
15212 normlist[coordindex++] = 0.0;
15213 elist[coordindex] = p2;
15214 normlist[coordindex++] = 0.0;
15216 fprintf(outfile,
"%4ld %d %d\n", vedgenumber, p1, p2);
15222 triangleloop.tri = triangletraverse(m);
15226 finishfile(outfile, argc, argv);
15232#ifdef ANSI_DECLARATORS
15233void writeneighbors(
struct mesh *m,
struct behavior *b,
int **neighborlist)
15235void writeneighbors(m, b, neighborlist)
15243#ifdef ANSI_DECLARATORS
15244void writeneighbors(
struct mesh *m,
struct behavior *b,
char *neighborfilename,
15245 int argc,
char **argv)
15247void writeneighbors(m, b, neighborfilename, argc, argv)
15250char *neighborfilename;
15264 struct otri triangleloop, trisym;
15265 long elementnumber;
15266 int neighbor1, neighbor2, neighbor3;
15271 printf(
"Writing neighbors.\n");
15274 if (*neighborlist == (
int *) NULL) {
15275 *neighborlist = (
int *) trimalloc((
int) (m->triangles.items * 3 *
15278 nlist = *neighborlist;
15282 printf(
"Writing %s.\n", neighborfilename);
15284 outfile = fopen(neighborfilename,
"w");
15285 if (outfile == (FILE *) NULL) {
15286 printf(
" Error: Cannot create file %s.\n", neighborfilename);
15290 fprintf(outfile,
"%ld %d\n", m->triangles.items, 3);
15293 traversalinit(&m->triangles);
15294 triangleloop.tri = triangletraverse(m);
15295 triangleloop.orient = 0;
15296 elementnumber = b->firstnumber;
15297 while (triangleloop.tri != (triangle *) NULL) {
15298 * (
int *) (triangleloop.tri + 6) = (int) elementnumber;
15299 triangleloop.tri = triangletraverse(m);
15302 * (
int *) (m->dummytri + 6) = -1;
15304 traversalinit(&m->triangles);
15305 triangleloop.tri = triangletraverse(m);
15306 elementnumber = b->firstnumber;
15307 while (triangleloop.tri != (triangle *) NULL) {
15308 triangleloop.orient = 1;
15309 sym(triangleloop, trisym);
15310 neighbor1 = * (
int *) (trisym.tri + 6);
15311 triangleloop.orient = 2;
15312 sym(triangleloop, trisym);
15313 neighbor2 = * (
int *) (trisym.tri + 6);
15314 triangleloop.orient = 0;
15315 sym(triangleloop, trisym);
15316 neighbor3 = * (
int *) (trisym.tri + 6);
15318 nlist[index++] = neighbor1;
15319 nlist[index++] = neighbor2;
15320 nlist[index++] = neighbor3;
15323 fprintf(outfile,
"%4ld %d %d %d\n", elementnumber,
15324 neighbor1, neighbor2, neighbor3);
15327 triangleloop.tri = triangletraverse(m);
15332 finishfile(outfile, argc, argv);
15347#ifdef ANSI_DECLARATORS
15348void writeoff(
struct mesh *m,
struct behavior *b,
char *offfilename,
15349 int argc,
char **argv)
15351void writeoff(m, b, offfilename, argc, argv)
15361 struct otri triangleloop;
15367 printf(
"Writing %s.\n", offfilename);
15371 outvertices = m->vertices.items - m->undeads;
15373 outvertices = m->vertices.items;
15376 outfile = fopen(offfilename,
"w");
15377 if (outfile == (FILE *) NULL) {
15378 printf(
" Error: Cannot create file %s.\n", offfilename);
15382 fprintf(outfile,
"OFF\n%ld %ld %ld\n", outvertices, m->triangles.items,
15386 traversalinit(&m->vertices);
15387 vertexloop = vertextraverse(m);
15388 while (vertexloop != (vertex) NULL) {
15389 if (!b->jettison || (vertextype(vertexloop) != UNDEADVERTEX)) {
15391 fprintf(outfile,
" %.17g %.17g %.17g\n", vertexloop[0], vertexloop[1],
15394 vertexloop = vertextraverse(m);
15398 traversalinit(&m->triangles);
15399 triangleloop.tri = triangletraverse(m);
15400 triangleloop.orient = 0;
15401 while (triangleloop.tri != (triangle *) NULL) {
15402 org(triangleloop, p1);
15403 dest(triangleloop, p2);
15404 apex(triangleloop, p3);
15406 fprintf(outfile,
" 3 %4d %4d %4d\n", vertexmark(p1) - b->firstnumber,
15407 vertexmark(p2) - b->firstnumber, vertexmark(p3) - b->firstnumber);
15408 triangleloop.tri = triangletraverse(m);
15410 finishfile(outfile, argc, argv);
15425#ifdef ANSI_DECLARATORS
15426void quality_statistics(
struct mesh *m,
struct behavior *b)
15428void quality_statistics(m, b)
15434 struct otri triangleloop;
15436 REAL cossquaretable[8];
15437 REAL ratiotable[16];
15439 REAL edgelength[3];
15443 REAL shortest, longest;
15445 REAL smallestarea, biggestarea;
15446 REAL triminaltitude2;
15450 REAL smallestangle, biggestangle;
15451 REAL radconst, degconst;
15452 int angletable[18];
15453 int aspecttable[16];
15459 printf(
"Mesh quality statistics:\n\n");
15460 radconst = PI / 18.0;
15461 degconst = 180.0 / PI;
15462 for (i = 0; i < 8; i++) {
15463 cossquaretable[i] = cos(radconst * (REAL) (i + 1));
15464 cossquaretable[i] = cossquaretable[i] * cossquaretable[i];
15466 for (i = 0; i < 18; i++) {
15470 ratiotable[0] = 1.5; ratiotable[1] = 2.0;
15471 ratiotable[2] = 2.5; ratiotable[3] = 3.0;
15472 ratiotable[4] = 4.0; ratiotable[5] = 6.0;
15473 ratiotable[6] = 10.0; ratiotable[7] = 15.0;
15474 ratiotable[8] = 25.0; ratiotable[9] = 50.0;
15475 ratiotable[10] = 100.0; ratiotable[11] = 300.0;
15476 ratiotable[12] = 1000.0; ratiotable[13] = 10000.0;
15477 ratiotable[14] = 100000.0; ratiotable[15] = 0.0;
15478 for (i = 0; i < 16; i++) {
15479 aspecttable[i] = 0;
15483 minaltitude = m->xmax - m->xmin + m->ymax - m->ymin;
15484 minaltitude = minaltitude * minaltitude;
15485 shortest = minaltitude;
15487 smallestarea = minaltitude;
15490 smallestangle = 0.0;
15491 biggestangle = 2.0;
15494 traversalinit(&m->triangles);
15495 triangleloop.tri = triangletraverse(m);
15496 triangleloop.orient = 0;
15497 while (triangleloop.tri != (triangle *) NULL) {
15498 org(triangleloop, p[0]);
15499 dest(triangleloop, p[1]);
15500 apex(triangleloop, p[2]);
15503 for (i = 0; i < 3; i++) {
15506 dx[i] = p[j][0] - p[k][0];
15507 dy[i] = p[j][1] - p[k][1];
15508 edgelength[i] = dx[i] * dx[i] + dy[i] * dy[i];
15509 if (edgelength[i] > trilongest2) {
15510 trilongest2 = edgelength[i];
15512 if (edgelength[i] > longest) {
15513 longest = edgelength[i];
15515 if (edgelength[i] < shortest) {
15516 shortest = edgelength[i];
15520 triarea = counterclockwise(m, b, p[0], p[1], p[2]);
15521 if (triarea < smallestarea) {
15522 smallestarea = triarea;
15524 if (triarea > biggestarea) {
15525 biggestarea = triarea;
15527 triminaltitude2 = triarea * triarea / trilongest2;
15528 if (triminaltitude2 < minaltitude) {
15529 minaltitude = triminaltitude2;
15531 triaspect2 = trilongest2 / triminaltitude2;
15532 if (triaspect2 > worstaspect) {
15533 worstaspect = triaspect2;
15536 while ((triaspect2 > ratiotable[aspectindex] * ratiotable[aspectindex])
15537 && (aspectindex < 15)) {
15540 aspecttable[aspectindex]++;
15542 for (i = 0; i < 3; i++) {
15545 dotproduct = dx[j] * dx[k] + dy[j] * dy[k];
15546 cossquare = dotproduct * dotproduct / (edgelength[j] * edgelength[k]);
15548 for (ii = 7; ii >= 0; ii--) {
15549 if (cossquare > cossquaretable[ii]) {
15553 if (dotproduct <= 0.0) {
15554 angletable[tendegree]++;
15555 if (cossquare > smallestangle) {
15556 smallestangle = cossquare;
15558 if (acutebiggest && (cossquare < biggestangle)) {
15559 biggestangle = cossquare;
15562 angletable[17 - tendegree]++;
15563 if (acutebiggest || (cossquare > biggestangle)) {
15564 biggestangle = cossquare;
15569 triangleloop.tri = triangletraverse(m);
15572 shortest = sqrt(shortest);
15573 longest = sqrt(longest);
15574 minaltitude = sqrt(minaltitude);
15575 worstaspect = sqrt(worstaspect);
15576 smallestarea *= 0.5;
15577 biggestarea *= 0.5;
15578 if (smallestangle >= 1.0) {
15579 smallestangle = 0.0;
15581 smallestangle = degconst * acos(sqrt(smallestangle));
15583 if (biggestangle >= 1.0) {
15584 biggestangle = 180.0;
15586 if (acutebiggest) {
15587 biggestangle = degconst * acos(sqrt(biggestangle));
15589 biggestangle = 180.0 - degconst * acos(sqrt(biggestangle));
15593 printf(
" Smallest area: %16.5g | Largest area: %16.5g\n",
15594 smallestarea, biggestarea);
15595 printf(
" Shortest edge: %16.5g | Longest edge: %16.5g\n",
15596 shortest, longest);
15597 printf(
" Shortest altitude: %12.5g | Largest aspect ratio: %8.5g\n\n",
15598 minaltitude, worstaspect);
15600 printf(
" Triangle aspect ratio histogram:\n");
15601 printf(
" 1.1547 - %-6.6g : %8d | %6.6g - %-6.6g : %8d\n",
15602 ratiotable[0], aspecttable[0], ratiotable[7], ratiotable[8],
15604 for (i = 1; i < 7; i++) {
15605 printf(
" %6.6g - %-6.6g : %8d | %6.6g - %-6.6g : %8d\n",
15606 ratiotable[i - 1], ratiotable[i], aspecttable[i],
15607 ratiotable[i + 7], ratiotable[i + 8], aspecttable[i + 8]);
15609 printf(
" %6.6g - %-6.6g : %8d | %6.6g - : %8d\n",
15610 ratiotable[6], ratiotable[7], aspecttable[7], ratiotable[14],
15612 printf(
" (Aspect ratio is longest edge divided by shortest altitude)\n\n");
15614 printf(
" Smallest angle: %15.5g | Largest angle: %15.5g\n\n",
15615 smallestangle, biggestangle);
15617 printf(
" Angle histogram:\n");
15618 for (i = 0; i < 9; i++) {
15619 printf(
" %3d - %3d degrees: %8d | %3d - %3d degrees: %8d\n",
15620 i * 10, i * 10 + 10, angletable[i],
15621 i * 10 + 90, i * 10 + 100, angletable[i + 9]);
15632#ifdef ANSI_DECLARATORS
15635void statistics(m, b)
15641 printf(
"\nStatistics:\n\n");
15642 printf(
" Input vertices: %d\n", m->invertices);
15644 printf(
" Input triangles: %d\n", m->inelements);
15647 printf(
" Input segments: %d\n", m->insegments);
15649 printf(
" Input holes: %d\n", m->holes);
15653 printf(
"\n Mesh vertices: %ld\n", m->vertices.items - m->undeads);
15654 printf(
" Mesh triangles: %ld\n", m->triangles.items);
15655 printf(
" Mesh edges: %ld\n", m->edges);
15656 printf(
" Mesh exterior boundary edges: %ld\n", m->hullsize);
15657 if (b->poly || b->refine) {
15658 printf(
" Mesh interior boundary edges: %ld\n",
15659 m->subsegs.items - m->hullsize);
15660 printf(
" Mesh subsegments (constrained edges): %ld\n",
15666 quality_statistics(m, b);
15667 printf(
"Memory allocation statistics:\n\n");
15668 printf(
" Maximum number of vertices: %ld\n", m->vertices.maxitems);
15669 printf(
" Maximum number of triangles: %ld\n", m->triangles.maxitems);
15670 if (m->subsegs.maxitems > 0) {
15671 printf(
" Maximum number of subsegments: %ld\n", m->subsegs.maxitems);
15673 if (m->viri.maxitems > 0) {
15674 printf(
" Maximum number of viri: %ld\n", m->viri.maxitems);
15676 if (m->badsubsegs.maxitems > 0) {
15677 printf(
" Maximum number of encroached subsegments: %ld\n",
15678 m->badsubsegs.maxitems);
15680 if (m->badtriangles.maxitems > 0) {
15681 printf(
" Maximum number of bad triangles: %ld\n",
15682 m->badtriangles.maxitems);
15684 if (m->flipstackers.maxitems > 0) {
15685 printf(
" Maximum number of stacked triangle flips: %ld\n",
15686 m->flipstackers.maxitems);
15688 if (m->splaynodes.maxitems > 0) {
15689 printf(
" Maximum number of splay tree nodes: %ld\n",
15690 m->splaynodes.maxitems);
15692 printf(
" Approximate heap memory use (bytes): %ld\n\n",
15693 m->vertices.maxitems * m->vertices.itembytes +
15694 m->triangles.maxitems * m->triangles.itembytes +
15695 m->subsegs.maxitems * m->subsegs.itembytes +
15696 m->viri.maxitems * m->viri.itembytes +
15697 m->badsubsegs.maxitems * m->badsubsegs.itembytes +
15698 m->badtriangles.maxitems * m->badtriangles.itembytes +
15699 m->flipstackers.maxitems * m->flipstackers.itembytes +
15700 m->splaynodes.maxitems * m->splaynodes.itembytes);
15702 printf(
"Algorithmic statistics:\n\n");
15703 if (!b->weighted) {
15704 printf(
" Number of incircle tests: %ld\n", m->incirclecount);
15706 printf(
" Number of 3D orientation tests: %ld\n", m->orient3dcount);
15708 printf(
" Number of 2D orientation tests: %ld\n", m->counterclockcount);
15709 if (m->hyperbolacount > 0) {
15710 printf(
" Number of right-of-hyperbola tests: %ld\n",
15711 m->hyperbolacount);
15713 if (m->circletopcount > 0) {
15714 printf(
" Number of circle top computations: %ld\n",
15715 m->circletopcount);
15717 if (m->circumcentercount > 0) {
15718 printf(
" Number of triangle circumcenter computations: %ld\n",
15719 m->circumcentercount);
15752#ifdef ANSI_DECLARATORS
15753void triangulate(
char *triswitches,
struct triangulateio *in,
15754 struct triangulateio *out,
struct triangulateio *vorout)
15756void triangulate(triswitches, in, out, vorout)
15758struct triangulateio *in;
15759struct triangulateio *out;
15760struct triangulateio *vorout;
15765#ifdef ANSI_DECLARATORS
15766int main(
int argc,
char **argv)
15768int main(argc, argv)
15786 struct timeval tv0, tv1, tv2, tv3, tv4, tv5, tv6;
15787 struct timezone tz;
15791 gettimeofday(&tv0, &tz);
15796 parsecommandline(1, &triswitches, &b);
15798 parsecommandline(argc, argv, &b);
15800 m.steinerleft = b.steiner;
15803 transfernodes(&m, &b, in->pointlist, in->pointattributelist,
15804 in->pointmarkerlist, in->numberofpoints,
15805 in->numberofpointattributes);
15807 readnodes(&m, &b, b.innodefilename, b.inpolyfilename, &polyfile);
15812 gettimeofday(&tv1, &tz);
15817 m.hullsize = delaunay(&m, &b);
15822 m.hullsize = reconstruct(&m, &b, in->trianglelist,
15823 in->triangleattributelist, in->trianglearealist,
15824 in->numberoftriangles, in->numberofcorners,
15825 in->numberoftriangleattributes,
15826 in->segmentlist, in->segmentmarkerlist,
15827 in->numberofsegments);
15829 m.hullsize = reconstruct(&m, &b, b.inelefilename, b.areafilename,
15830 b.inpolyfilename, polyfile);
15833 m.hullsize = delaunay(&m, &b);
15839 gettimeofday(&tv2, &tz);
15841 printf(
"Mesh reconstruction");
15843 printf(
"Delaunay");
15845 printf(
" milliseconds: %ld\n", 1000l * (tv2.tv_sec - tv1.tv_sec) +
15846 (tv2.tv_usec - tv1.tv_usec) / 1000l);
15852 m.infvertex1 = (vertex) NULL;
15853 m.infvertex2 = (vertex) NULL;
15854 m.infvertex3 = (vertex) NULL;
15856 if (b.usesegments) {
15857 m.checksegments = 1;
15861 formskeleton(&m, &b, in->segmentlist,
15862 in->segmentmarkerlist, in->numberofsegments);
15864 formskeleton(&m, &b, polyfile, b.inpolyfilename);
15871 gettimeofday(&tv3, &tz);
15872 if (b.usesegments && !b.refine) {
15873 printf(
"Segment milliseconds: %ld\n",
15874 1000l * (tv3.tv_sec - tv2.tv_sec) +
15875 (tv3.tv_usec - tv2.tv_usec) / 1000l);
15880 if (b.poly && (m.triangles.items > 0)) {
15882 holearray = in->holelist;
15883 m.holes = in->numberofholes;
15884 regionarray = in->regionlist;
15885 m.regions = in->numberofregions;
15887 readholes(&m, &b, polyfile, b.inpolyfilename, &holearray, &m.holes,
15888 ®ionarray, &m.regions);
15892 carveholes(&m, &b, holearray, m.holes, regionarray, m.regions);
15904 gettimeofday(&tv4, &tz);
15905 if (b.poly && !b.refine) {
15906 printf(
"Hole milliseconds: %ld\n", 1000l * (tv4.tv_sec - tv3.tv_sec) +
15907 (tv4.tv_usec - tv3.tv_usec) / 1000l);
15913 if (b.quality && (m.triangles.items > 0)) {
15914 enforcequality(&m, &b);
15920 gettimeofday(&tv5, &tz);
15923 printf(
"Quality milliseconds: %ld\n",
15924 1000l * (tv5.tv_sec - tv4.tv_sec) +
15925 (tv5.tv_usec - tv4.tv_usec) / 1000l);
15932 m.edges = (3l * m.triangles.items + m.hullsize) / 2l;
15943 out->numberofpoints = m.vertices.items - m.undeads;
15945 out->numberofpoints = m.vertices.items;
15947 out->numberofpointattributes = m.nextras;
15948 out->numberoftriangles = m.triangles.items;
15949 out->numberofcorners = (b.order + 1) * (b.order + 2) / 2;
15950 out->numberoftriangleattributes = m.eextras;
15951 out->numberofedges = m.edges;
15952 if (b.usesegments) {
15953 out->numberofsegments = m.subsegs.items;
15955 out->numberofsegments = m.hullsize;
15957 if (vorout != (
struct triangulateio *) NULL) {
15958 vorout->numberofpoints = m.triangles.items;
15959 vorout->numberofpointattributes = m.nextras;
15960 vorout->numberofedges = m.edges;
15965 if (b.nonodewritten || (b.noiterationnum && m.readnodefile)) {
15968 printf(
"NOT writing vertices.\n");
15970 printf(
"NOT writing a .node file.\n");
15973 numbernodes(&m, &b);
15977 writenodes(&m, &b, &out->pointlist, &out->pointattributelist,
15978 &out->pointmarkerlist);
15980 writenodes(&m, &b, b.outnodefilename, argc, argv);
15983 if (b.noelewritten) {
15986 printf(
"NOT writing triangles.\n");
15988 printf(
"NOT writing an .ele file.\n");
15993 writeelements(&m, &b, &out->trianglelist, &out->triangleattributelist);
15995 writeelements(&m, &b, b.outelefilename, argc, argv);
16000 if (b.poly || b.convex) {
16002 if (b.nopolywritten || b.noiterationnum) {
16005 printf(
"NOT writing segments.\n");
16007 printf(
"NOT writing a .poly file.\n");
16012 writepoly(&m, &b, &out->segmentlist, &out->segmentmarkerlist);
16013 out->numberofholes = m.holes;
16014 out->numberofregions = m.regions;
16016 out->holelist = in->holelist;
16017 out->regionlist = in->regionlist;
16019 out->holelist = (REAL *) NULL;
16020 out->regionlist = (REAL *) NULL;
16023 writepoly(&m, &b, b.outpolyfilename, holearray, m.holes, regionarray,
16024 m.regions, argc, argv);
16030 if (m.regions > 0) {
16031 trifree((VOID *) regionarray);
16035 trifree((VOID *) holearray);
16038 writeoff(&m, &b, b.offfilename, argc, argv);
16043 writeedges(&m, &b, &out->edgelist, &out->edgemarkerlist);
16045 writeedges(&m, &b, b.edgefilename, argc, argv);
16050 writevoronoi(&m, &b, &vorout->pointlist, &vorout->pointattributelist,
16051 &vorout->pointmarkerlist, &vorout->edgelist,
16052 &vorout->edgemarkerlist, &vorout->normlist);
16054 writevoronoi(&m, &b, b.vnodefilename, b.vedgefilename, argc, argv);
16059 writeneighbors(&m, &b, &out->neighborlist);
16061 writeneighbors(&m, &b, b.neighborfilename, argc, argv);
16067 gettimeofday(&tv6, &tz);
16068 printf(
"\nOutput milliseconds: %ld\n",
16069 1000l * (tv6.tv_sec - tv5.tv_sec) +
16070 (tv6.tv_usec - tv5.tv_usec) / 1000l);
16071 printf(
"Total running milliseconds: %ld\n",
16072 1000l * (tv6.tv_sec - tv0.tv_sec) +
16073 (tv6.tv_usec - tv0.tv_usec) / 1000l);
16076 statistics(&m, &b);
16082 checkdelaunay(&m, &b);
16086 triangledeinit(&m, &b);