~ [ source navigation ] ~ [ diff markup ] ~ [ identifier search ] ~ [ freetext search ] ~ [ file search ] ~

Linux Cross Reference
Tina5/tina-tools/tinatool/tlvision/tlvisTrn_triangle.c

Version: ~
Architecture: ~ [ i386 ] ~ [ alpha ] ~ [ m68k ] ~ [ mips ] ~ [ ppc ] ~ [ sparc ] ~ [ sparc64 ] ~

  1 /*****************************************************************************/
  2 /*                                                                           */
  3 /*      888888888        ,o,                          / 888                  */
  4 /*         888    88o88o  "    o8888o  88o8888o o88888o 888  o88888o         */
  5 /*         888    888    888       88b 888  888 888 888 888 d888  88b        */
  6 /*         888    888    888  o88^o888 888  888 "88888" 888 8888oo888        */
  7 /*         888    888    888 C888  888 888  888  /      888 q888             */
  8 /*         888    888    888  "88o^888 888  888 Cb      888  "88oooo"        */
  9 /*                                              "8oo8D                       */
 10 /*                                                                           */
 11 /*  A Two-Dimensional Quality Mesh Generator and Delaunay Triangulator.      */
 12 /*  (triangle.c)                                                             */
 13 /*                                                                           */
 14 /*  Version 1.3                                                              */
 15 /*  July 19, 1996                                                            */
 16 /*                                                                           */
 17 /*  Copyright 1996                                                           */
 18 /*  Jonathan Richard Shewchuk                                                */
 19 /*  School of Computer Science                                               */
 20 /*  Carnegie Mellon University                                               */
 21 /*  5000 Forbes Avenue                                                       */
 22 /*  Pittsburgh, Pennsylvania  15213-3891                                     */
 23 /*  jrs@cs.cmu.edu                                                           */
 24 /*                                                                           */
 25 /*  This program may be freely redistributed under the condition that the    */
 26 /*    copyright notices (including this entire header and the copyright      */
 27 /*    notice printed when the `-h' switch is selected) are not removed, and  */
 28 /*    no compensation is received.  Private, research, and institutional     */
 29 /*    use is free.  You may distribute modified versions of this code UNDER  */
 30 /*    THE CONDITION THAT THIS CODE AND ANY MODIFICATIONS MADE TO IT IN THE   */
 31 /*    SAME FILE REMAIN UNDER COPYRIGHT OF THE ORIGINAL AUTHOR, BOTH SOURCE   */
 32 /*    AND OBJECT CODE ARE MADE FREELY AVAILABLE WITHOUT CHARGE, AND CLEAR    */
 33 /*    NOTICE IS GIVEN OF THE MODIFICATIONS.  Distribution of this code as    */
 34 /*    part of a commercial system is permissible ONLY BY DIRECT ARRANGEMENT  */
 35 /*    WITH THE AUTHOR.  (If you are not directly supplying this code to a    */
 36 /*    customer, and you are instead telling them how they can obtain it for  */
 37 /*    free, then you are not required to make any arrangement with me.)      */
 38 /*                                                                           */
 39 /*  Hypertext instructions for Triangle are available on the Web at          */
 40 /*                                                                           */
 41 /*      http://www.cs.cmu.edu/~quake/triangle.html                           */
 42 /*                                                                           */
 43 /*  Some of the references listed below are marked [*].  These are available */
 44 /*    for downloading from the Web page                                      */
 45 /*                                                                           */
 46 /*      http://www.cs.cmu.edu/~quake/triangle.research.html                  */
 47 /*                                                                           */
 48 /*  A paper discussing some aspects of Triangle is available.  See Jonathan  */
 49 /*    Richard Shewchuk, "Triangle:  Engineering a 2D Quality Mesh Generator  */
 50 /*    and Delaunay Triangulator," First Workshop on Applied Computational    */
 51 /*    Geometry, ACM, May 1996.  [*]                                          */
 52 /*                                                                           */
 53 /*  Triangle was created as part of the Archimedes project in the School of  */
 54 /*    Computer Science at Carnegie Mellon University.  Archimedes is a       */
 55 /*    system for compiling parallel finite element solvers.  For further     */
 56 /*    information, see Anja Feldmann, Omar Ghattas, John R. Gilbert, Gary L. */
 57 /*    Miller, David R. O'Hallaron, Eric J. Schwabe, Jonathan R. Shewchuk,    */
 58 /*    and Shang-Hua Teng, "Automated Parallel Solution of Unstructured PDE   */
 59 /*    Problems."  To appear in Communications of the ACM, we hope.           */
 60 /*                                                                           */
 61 /*  The quality mesh generation algorithm is due to Jim Ruppert, "A          */
 62 /*    Delaunay Refinement Algorithm for Quality 2-Dimensional Mesh           */
 63 /*    Generation," Journal of Algorithms 18(3):548-585, May 1995.  [*]       */
 64 /*                                                                           */
 65 /*  My implementation of the divide-and-conquer and incremental Delaunay     */
 66 /*    triangulation algorithms follows closely the presentation of Guibas    */
 67 /*    and Stolfi, even though I use a triangle-based data structure instead  */
 68 /*    of their quad-edge data structure.  (In fact, I originally implemented */
 69 /*    Triangle using the quad-edge data structure, but switching to a        */
 70 /*    triangle-based data structure sped Triangle by a factor of two.)  The  */
 71 /*    mesh manipulation primitives and the two aforementioned Delaunay       */
 72 /*    triangulation algorithms are described by Leonidas J. Guibas and Jorge */
 73 /*    Stolfi, "Primitives for the Manipulation of General Subdivisions and   */
 74 /*    the Computation of Voronoi Diagrams," ACM Transactions on Graphics     */
 75 /*    4(2):74-123, April 1985.                                               */
 76 /*                                                                           */
 77 /*  Their O(n log n) divide-and-conquer algorithm is adapted from Der-Tsai   */
 78 /*    Lee and Bruce J. Schachter, "Two Algorithms for Constructing the       */
 79 /*    Delaunay Triangulation," International Journal of Computer and         */
 80 /*    Information Science 9(3):219-242, 1980.  The idea to improve the       */
 81 /*    divide-and-conquer algorithm by alternating between vertical and       */
 82 /*    horizontal cuts was introduced by Rex A. Dwyer, "A Faster Divide-and-  */
 83 /*    Conquer Algorithm for Constructing Delaunay Triangulations,"           */
 84 /*    Algorithmica 2(2):137-151, 1987.                                       */
 85 /*                                                                           */
 86 /*  The incremental insertion algorithm was first proposed by C. L. Lawson,  */
 87 /*    "Software for C1 Surface Interpolation," in Mathematical Software III, */
 88 /*    John R. Rice, editor, Academic Press, New York, pp. 161-194, 1977.     */
 89 /*    For point location, I use the algorithm of Ernst P. Mucke, Isaac       */
 90 /*    Saias, and Binhai Zhu, "Fast Randomized Point Location Without         */
 91 /*    Preprocessing in Two- and Three-dimensional Delaunay Triangulations,"  */
 92 /*    Proceedings of the Twelfth Annual Symposium on Computational Geometry, */
 93 /*    ACM, May 1996.  [*]  If I were to randomize the order of point         */
 94 /*    insertion (I currently don't bother), their result combined with the   */
 95 /*    result of Leonidas J. Guibas, Donald E. Knuth, and Micha Sharir,       */
 96 /*    "Randomized Incremental Construction of Delaunay and Voronoi           */
 97 /*    Diagrams," Algorithmica 7(4):381-413, 1992, would yield an expected    */
 98 /*    O(n^{4/3}) bound on running time.                                      */
 99 /*                                                                           */
100 /*  The O(n log n) sweepline Delaunay triangulation algorithm is taken from  */
101 /*    Steven Fortune, "A Sweepline Algorithm for Voronoi Diagrams",          */
102 /*    Algorithmica 2(2):153-174, 1987.  A random sample of edges on the      */
103 /*    boundary of the triangulation are maintained in a splay tree for the   */
104 /*    purpose of point location.  Splay trees are described by Daniel        */
105 /*    Dominic Sleator and Robert Endre Tarjan, "Self-Adjusting Binary Search */
106 /*    Trees," Journal of the ACM 32(3):652-686, July 1985.                   */
107 /*                                                                           */
108 /*  The algorithms for exact computation of the signs of determinants are    */
109 /*    described in Jonathan Richard Shewchuk, "Adaptive Precision Floating-  */
110 /*    Point Arithmetic and Fast Robust Geometric Predicates," Technical      */
111 /*    Report CMU-CS-96-140, School of Computer Science, Carnegie Mellon      */
112 /*    University, Pittsburgh, Pennsylvania, May 1996.  [*]  (Submitted to    */
113 /*    Discrete & Computational Geometry.)  An abbreviated version appears as */
114 /*    Jonathan Richard Shewchuk, "Robust Adaptive Floating-Point Geometric   */
115 /*    Predicates," Proceedings of the Twelfth Annual Symposium on Computa-   */
116 /*    tional Geometry, ACM, May 1996.  [*]  Many of the ideas for my exact   */
117 /*    arithmetic routines originate with Douglas M. Priest, "Algorithms for  */
118 /*    Arbitrary Precision Floating Point Arithmetic," Tenth Symposium on     */
119 /*    Computer Arithmetic, 132-143, IEEE Computer Society Press, 1991.  [*]  */
120 /*    Many of the ideas for the correct evaluation of the signs of           */
121 /*    determinants are taken from Steven Fortune and Christopher J. Van Wyk, */
122 /*    "Efficient Exact Arithmetic for Computational Geometry," Proceedings   */
123 /*    of the Ninth Annual Symposium on Computational Geometry, ACM,          */
124 /*    pp. 163-172, May 1993, and from Steven Fortune, "Numerical Stability   */
125 /*    of Algorithms for 2D Delaunay Triangulations," International Journal   */
126 /*    of Computational Geometry & Applications 5(1-2):193-213, March-June    */
127 /*    1995.                                                                  */
128 /*                                                                           */
129 /*  For definitions of and results involving Delaunay triangulations,        */
130 /*    constrained and conforming versions thereof, and other aspects of      */
131 /*    triangular mesh generation, see the excellent survey by Marshall Bern  */
132 /*    and David Eppstein, "Mesh Generation and Optimal Triangulation," in    */
133 /*    Computing and Euclidean Geometry, Ding-Zhu Du and Frank Hwang,         */
134 /*    editors, World Scientific, Singapore, pp. 23-90, 1992.                 */
135 /*                                                                           */
136 /*  The time for incrementally adding PSLG (planar straight line graph)      */
137 /*    segments to create a constrained Delaunay triangulation is probably    */
138 /*    O(n^2) per segment in the worst case and O(n) per edge in the common   */
139 /*    case, where n is the number of triangles that intersect the segment    */
140 /*    before it is inserted.  This doesn't count point location, which can   */
141 /*    be much more expensive.  (This note does not apply to conforming       */
142 /*    Delaunay triangulations, for which a different method is used to       */
143 /*    insert segments.)                                                      */
144 /*                                                                           */
145 /*  The time for adding segments to a conforming Delaunay triangulation is   */
146 /*    not clear, but does not depend upon n alone.  In some cases, very      */
147 /*    small features (like a point lying next to a segment) can cause a      */
148 /*    single segment to be split an arbitrary number of times.  Of course,   */
149 /*    floating-point precision is a practical barrier to how much this can   */
150 /*    happen.                                                                */
151 /*                                                                           */
152 /*  The time for deleting a point from a Delaunay triangulation is O(n^2) in */
153 /*    the worst case and O(n) in the common case, where n is the degree of   */
154 /*    the point being deleted.  I could improve this to expected O(n) time   */
155 /*    by "inserting" the neighboring vertices in random order, but n is      */
156 /*    usually quite small, so it's not worth the bother.  (The O(n) time     */
157 /*    for random insertion follows from L. Paul Chew, "Building Voronoi      */
158 /*    Diagrams for Convex Polygons in Linear Expected Time," Technical       */
159 /*    Report PCS-TR90-147, Department of Mathematics and Computer Science,   */
160 /*    Dartmouth College, 1990.                                               */
161 /*                                                                           */
162 /*  Ruppert's Delaunay refinement algorithm typically generates triangles    */
163 /*    at a linear rate (constant time per triangle) after the initial        */
164 /*    triangulation is formed.  There may be pathological cases where more   */
165 /*    time is required, but these never arise in practice.                   */
166 /*                                                                           */
167 /*  The segment intersection formulae are straightforward.  If you want to   */
168 /*    see them derived, see Franklin Antonio.  "Faster Line Segment          */
169 /*    Intersection."  In Graphics Gems III (David Kirk, editor), pp. 199-    */
170 /*    202.  Academic Press, Boston, 1992.                                    */
171 /*                                                                           */
172 /*  If you make any improvements to this code, please please please let me   */
173 /*    know, so that I may obtain the improvements.  Even if you don't change */
174 /*    the code, I'd still love to hear what it's being used for.             */
175 /*                                                                           */
176 /*  Disclaimer:  Neither I nor Carnegie Mellon warrant this code in any way  */
177 /*    whatsoever.  This code is provided "as-is".  Use at your own risk.     */
178 /*                                                                           */
179 /*****************************************************************************/
180 
181 /* For single precision (which will save some memory and reduce paging),     */
182 /*   define the symbol SINGLE by using the -DSINGLE compiler switch or by    */
183 /*   writing "#define SINGLE" below.                                         */
184 /*                                                                           */
185 /* For double precision (which will allow you to refine meshes to a smaller  */
186 /*   edge length), leave SINGLE undefined.                                   */
187 /*                                                                           */
188 /* Double precision uses more memory, but improves the resolution of the     */
189 /*   meshes you can generate with Triangle.  It also reduces the likelihood  */
190 /*   of a floating exception due to overflow.  Finally, it is much faster    */
191 /*   than single precision on 64-bit architectures like the DEC Alpha.  I    */
192 /*   recommend double precision unless you want to generate a mesh for which */
193 /*   you do not have enough memory.                                          */
194 
195 /* #define SINGLE */
196 
197 #ifdef SINGLE
198 #define REAL float
199 #else /* not SINGLE */
200 #define REAL double
201 #endif /* not SINGLE */
202 
203 /* If yours is not a Unix system, define the NO_TIMER compiler switch to     */
204 /*   remove the Unix-specific timing code.                                   */
205 
206 #define NO_TIMER 1
207 
208 /* To insert lots of self-checks for internal errors, define the SELF_CHECK  */
209 /*   symbol.  This will slow down the program significantly.  It is best to  */
210 /*   define the symbol using the -DSELF_CHECK compiler switch, but you could */
211 /*   write "#define SELF_CHECK" below.  If you are modifying this code, I    */
212 /*   recommend you turn self-checks on.                                      */
213 
214 /* #define SELF_CHECK */
215 
216 /* To compile Triangle as a callable object library (triangle.o), define the */
217 /*   TRILIBRARY symbol.  Read the file triangle.h for details on how to call */
218 /*   the procedure triangulate() that results.                               */
219 
220 #define TRILIBRARY
221 
222 /* It is possible to generate a smaller version of Triangle using one or     */
223 /*   both of the following symbols.  Define the REDUCED symbol to eliminate  */
224 /*   all features that are primarily of research interest; specifically, the */
225 /*   -i, -F, -s, and -C switches.  Define the CDT_ONLY symbol to eliminate   */
226 /*   all meshing algorithms above and beyond constrained Delaunay            */
227 /*   triangulation; specifically, the -r, -q, -a, -S, and -s switches.       */
228 /*   These reductions are most likely to be useful when generating an object */
229 /*   library (triangle.o) by defining the TRILIBRARY symbol.                 */
230 
231 #define REDUCED
232 #define CDT_ONLY
233 
234 /* On some machines, the exact arithmetic routines might be defeated by the  */
235 /*   use of internal extended precision floating-point registers.  Sometimes */
236 /*   this problem can be fixed by defining certain values to be volatile,    */
237 /*   thus forcing them to be stored to memory and rounded off.  This isn't   */
238 /*   a great solution, though, as it slows Triangle down.                    */
239 /*                                                                           */
240 /* To try this out, write "#define INEXACT volatile" below.  Normally,       */
241 /*   however, INEXACT should be defined to be nothing.  ("#define INEXACT".) */
242 
243 #define INEXACT /* Nothing */
244 /* #define INEXACT volatile */
245 
246 /* Maximum number of characters in a file name (including the null).         */
247 
248 #define FILENAMESIZE 512
249 
250 /* Maximum number of characters in a line read from a file (including the    */
251 /*   null).                                                                  */
252 
253 #define INPUTLINESIZE 512
254 
255 /* For efficiency, a variety of data structures are allocated in bulk.  The  */
256 /*   following constants determine how many of each structure is allocated   */
257 /*   at once.                                                                */
258 
259 #define TRIPERBLOCK 4092           /* Number of triangles allocated at once. */
260 #define SHELLEPERBLOCK 508       /* Number of shell edges allocated at once. */
261 #define POINTPERBLOCK 4092            /* Number of points allocated at once. */
262 #define VIRUSPERBLOCK 1020   /* Number of virus triangles allocated at once. */
263 /* Number of encroached segments allocated at once. */
264 #define BADSEGMENTPERBLOCK 252
265 /* Number of skinny triangles allocated at once. */
266 #define BADTRIPERBLOCK 4092
267 /* Number of splay tree nodes allocated at once. */
268 #define SPLAYNODEPERBLOCK 508
269 
270 /* The point marker DEADPOINT is an arbitrary number chosen large enough to  */
271 /*   (hopefully) not conflict with user boundary markers.  Make sure that it */
272 /*   is small enough to fit into your machine's integer size.                */
273 
274 #define DEADPOINT -1073741824
275 
276 /* The next line is used to outsmart some very stupid compilers.  If your    */
277 /*   compiler is smarter, feel free to replace the "int" with "void".        */
278 /*   Not that it matters.                                                    */
279 
280 #define VOID int
281 
282 /* Two constants for algorithms based on random sampling.  Both constants    */
283 /*   have been chosen empirically to optimize their respective algorithms.   */
284 
285 /* Used for the point location scheme of Mucke, Saias, and Zhu, to decide    */
286 /*   how large a random sample of triangles to inspect.                      */
287 #define SAMPLEFACTOR 11
288 /* Used in Fortune's sweepline Delaunay algorithm to determine what fraction */
289 /*   of boundary edges should be maintained in the splay tree for point      */
290 /*   location on the front.                                                  */
291 #define SAMPLERATE 10
292 
293 /* A number that speaks for itself, every kissable digit.                    */
294 
295 #define PI 3.141592653589793238462643383279502884197169399375105820974944592308
296 
297 /* Another fave.                                                             */
298 
299 #define SQUAREROOTTWO 1.4142135623730950488016887242096980785696718753769480732
300 
301 /* And here's one for those of you who are intimidated by math.              */
302 
303 #define ONETHIRD 0.333333333333333333333333333333333333333333333333333333333333
304 
305 #include <stdio.h>
306 #include <string.h>
307 #include <math.h>
308 #ifndef NO_TIMER
309 #include <sys/time.h>
310 #endif /* NO_TIMER */
311 #ifdef TRILIBRARY
312 #include <tinatool/tlvision/tlvisTrn_triangle.h>
313 #endif /* TRILIBRARY */
314 
315 /* The following obscenity seems to be necessary to ensure that this program */
316 /* will port to Dec Alphas running OSF/1, because their stdio.h file commits */
317 /* the unpardonable sin of including stdlib.h.  Hence, malloc(), free(), and */
318 /* exit() may or may not already be defined at this point.  I declare these  */
319 /* functions explicitly because some non-ANSI C compilers lack stdlib.h.     */
320 
321 #ifndef _STDLIB_H_
322 extern void *malloc();
323 extern void free();
324 extern void exit();
325 extern double strtod();
326 extern long strtol();
327 #endif /* _STDLIB_H_ */
328 
329 /* A few forward declarations.                                               */
330 
331 void poolrestart();
332 #ifndef TRILIBRARY
333 char *freadline();
334 char *findfield();
335 #endif /* not TRILIBRARY */
336 
337 /* Labels that signify whether a record consists primarily of pointers or of */
338 /*   floating-point words.  Used to make decisions about data alignment.     */
339 
340 enum wordtype {POINTER, FLOATINGPOINT};
341 
342 /* Labels that signify the result of point location.  The result of a        */
343 /*   search indicates that the point falls in the interior of a triangle, on */
344 /*   an edge, on a vertex, or outside the mesh.                              */
345 
346 enum locateresult {INTRIANGLE, ONEDGE, ONVERTEX, OUTSIDE};
347 
348 /* Labels that signify the result of site insertion.  The result indicates   */
349 /*   that the point was inserted with complete success, was inserted but     */
350 /*   encroaches on a segment, was not inserted because it lies on a segment, */
351 /*   or was not inserted because another point occupies the same location.   */
352 
353 enum insertsiteresult {SUCCESSFULPOINT, ENCROACHINGPOINT, VIOLATINGPOINT,
354                        DUPLICATEPOINT};
355 
356 /* Labels that signify the result of direction finding.  The result          */
357 /*   indicates that a segment connecting the two query points falls within   */
358 /*   the direction triangle, along the left edge of the direction triangle,  */
359 /*   or along the right edge of the direction triangle.                      */
360 
361 enum finddirectionresult {WITHIN, LEFTCOLLINEAR, RIGHTCOLLINEAR};
362 
363 /* Labels that signify the result of the circumcenter computation routine.   */
364 /*   The return value indicates which edge of the triangle is shortest.      */
365 
366 enum circumcenterresult {OPPOSITEORG, OPPOSITEDEST, OPPOSITEAPEX};
367 
368 /*****************************************************************************/
369 /*                                                                           */
370 /*  The basic mesh data structures                                           */
371 /*                                                                           */
372 /*  There are three:  points, triangles, and shell edges (abbreviated        */
373 /*  `shelle').  These three data structures, linked by pointers, comprise    */
374 /*  the mesh.  A point simply represents a point in space and its properties.*/
375 /*  A triangle is a triangle.  A shell edge is a special data structure used */
376 /*  to represent impenetrable segments in the mesh (including the outer      */
377 /*  boundary, boundaries of holes, and internal boundaries separating two    */
378 /*  triangulated regions).  Shell edges represent boundaries defined by the  */
379 /*  user that triangles may not lie across.                                  */
380 /*                                                                           */
381 /*  A triangle consists of a list of three vertices, a list of three         */
382 /*  adjoining triangles, a list of three adjoining shell edges (when shell   */
383 /*  edges are used), an arbitrary number of optional user-defined floating-  */
384 /*  point attributes, and an optional area constraint.  The latter is an     */
385 /*  upper bound on the permissible area of each triangle in a region, used   */
386 /*  for mesh refinement.                                                     */
387 /*                                                                           */
388 /*  For a triangle on a boundary of the mesh, some or all of the neighboring */
389 /*  triangles may not be present.  For a triangle in the interior of the     */
390 /*  mesh, often no neighboring shell edges are present.  Such absent         */
391 /*  triangles and shell edges are never represented by NULL pointers; they   */
392 /*  are represented by two special records:  `dummytri', the triangle that   */
393 /*  fills "outer space", and `dummysh', the omnipresent shell edge.          */
394 /*  `dummytri' and `dummysh' are used for several reasons; for instance,     */
395 /*  they can be dereferenced and their contents examined without causing the */
396 /*  memory protection exception that would occur if NULL were dereferenced.  */
397 /*                                                                           */
398 /*  However, it is important to understand that a triangle includes other    */
399 /*  information as well.  The pointers to adjoining vertices, triangles, and */
400 /*  shell edges are ordered in a way that indicates their geometric relation */
401 /*  to each other.  Furthermore, each of these pointers contains orientation */
402 /*  information.  Each pointer to an adjoining triangle indicates which face */
403 /*  of that triangle is contacted.  Similarly, each pointer to an adjoining  */
404 /*  shell edge indicates which side of that shell edge is contacted, and how */
405 /*  the shell edge is oriented relative to the triangle.                     */
406 /*                                                                           */
407 /*  Shell edges are found abutting edges of triangles; either sandwiched     */
408 /*  between two triangles, or resting against one triangle on an exterior    */
409 /*  boundary or hole boundary.                                               */
410 /*                                                                           */
411 /*  A shell edge consists of a list of two vertices, a list of two           */
412 /*  adjoining shell edges, and a list of two adjoining triangles.  One of    */
413 /*  the two adjoining triangles may not be present (though there should      */
414 /*  always be one), and neighboring shell edges might not be present.        */
415 /*  Shell edges also store a user-defined integer "boundary marker".         */
416 /*  Typically, this integer is used to indicate what sort of boundary        */
417 /*  conditions are to be applied at that location in a finite element        */
418 /*  simulation.                                                              */
419 /*                                                                           */
420 /*  Like triangles, shell edges maintain information about the relative      */
421 /*  orientation of neighboring objects.                                      */
422 /*                                                                           */
423 /*  Points are relatively simple.  A point is a list of floating point       */
424 /*  numbers, starting with the x, and y coordinates, followed by an          */
425 /*  arbitrary number of optional user-defined floating-point attributes,     */
426 /*  followed by an integer boundary marker.  During the segment insertion    */
427 /*  phase, there is also a pointer from each point to a triangle that may    */
428 /*  contain it.  Each pointer is not always correct, but when one is, it     */
429 /*  speeds up segment insertion.  These pointers are assigned values once    */
430 /*  at the beginning of the segment insertion phase, and are not used or     */
431 /*  updated at any other time.  Edge swapping during segment insertion will  */
432 /*  render some of them incorrect.  Hence, don't rely upon them for          */
433 /*  anything.  For the most part, points do not have any information about   */
434 /*  what triangles or shell edges they are linked to.                        */
435 /*                                                                           */
436 /*****************************************************************************/
437 
438 /*****************************************************************************/
439 /*                                                                           */
440 /*  Handles                                                                  */
441 /*                                                                           */
442 /*  The oriented triangle (`triedge') and oriented shell edge (`edge') data  */
443 /*  structures defined below do not themselves store any part of the mesh.   */
444 /*  The mesh itself is made of `triangle's, `shelle's, and `point's.         */
445 /*                                                                           */
446 /*  Oriented triangles and oriented shell edges will usually be referred to  */
447 /*  as "handles".  A handle is essentially a pointer into the mesh; it       */
448 /*  allows you to "hold" one particular part of the mesh.  Handles are used  */
449 /*  to specify the regions in which one is traversing and modifying the mesh.*/
450 /*  A single `triangle' may be held by many handles, or none at all.  (The   */
451 /*  latter case is not a memory leak, because the triangle is still          */
452 /*  connected to other triangles in the mesh.)                               */
453 /*                                                                           */
454 /*  A `triedge' is a handle that holds a triangle.  It holds a specific side */
455 /*  of the triangle.  An `edge' is a handle that holds a shell edge.  It     */
456 /*  holds either the left or right side of the edge.                         */
457 /*                                                                           */
458 /*  Navigation about the mesh is accomplished through a set of mesh          */
459 /*  manipulation primitives, further below.  Many of these primitives take   */
460 /*  a handle and produce a new handle that holds the mesh near the first     */
461 /*  handle.  Other primitives take two handles and glue the corresponding    */
462 /*  parts of the mesh together.  The exact position of the handles is        */
463 /*  important.  For instance, when two triangles are glued together by the   */
464 /*  bond() primitive, they are glued by the sides on which the handles lie.  */
465 /*                                                                           */
466 /*  Because points have no information about which triangles they are        */
467 /*  attached to, I commonly represent a point by use of a handle whose       */
468 /*  origin is the point.  A single handle can simultaneously represent a     */
469 /*  triangle, an edge, and a point.                                          */
470 /*                                                                           */
471 /*****************************************************************************/
472 
473 /* The triangle data structure.  Each triangle contains three pointers to    */
474 /*   adjoining triangles, plus three pointers to vertex points, plus three   */
475 /*   pointers to shell edges (defined below; these pointers are usually      */
476 /*   `dummysh').  It may or may not also contain user-defined attributes     */
477 /*   and/or a floating-point "area constraint".  It may also contain extra   */
478 /*   pointers for nodes, when the user asks for high-order elements.         */
479 /*   Because the size and structure of a `triangle' is not decided until     */
480 /*   runtime, I haven't simply defined the type `triangle' to be a struct.   */
481 
482 typedef REAL **triangle;            /* Really:  typedef triangle *triangle   */
483 
484 /* An oriented triangle:  includes a pointer to a triangle and orientation.  */
485 /*   The orientation denotes an edge of the triangle.  Hence, there are      */
486 /*   three possible orientations.  By convention, each edge is always        */
487 /*   directed to point counterclockwise about the corresponding triangle.    */
488 
489 struct triedge {
490   triangle *tri;
491   int orient;                                         /* Ranges from 0 to 2. */
492 };
493 
494 /* The shell data structure.  Each shell edge contains two pointers to       */
495 /*   adjoining shell edges, plus two pointers to vertex points, plus two     */
496 /*   pointers to adjoining triangles, plus one shell marker.                 */
497 
498 typedef REAL **shelle;                  /* Really:  typedef shelle *shelle   */
499 
500 /* An oriented shell edge:  includes a pointer to a shell edge and an        */
501 /*   orientation.  The orientation denotes a side of the edge.  Hence, there */
502 /*   are two possible orientations.  By convention, the edge is always       */
503 /*   directed so that the "side" denoted is the right side of the edge.      */
504 
505 struct edge {
506   shelle *sh;
507   int shorient;                                       /* Ranges from 0 to 1. */
508 };
509 
510 /* The point data structure.  Each point is actually an array of REALs.      */
511 /*   The number of REALs is unknown until runtime.  An integer boundary      */
512 /*   marker, and sometimes a pointer to a triangle, is appended after the    */
513 /*   REALs.                                                                  */
514 
515 typedef REAL *point;
516 
517 /* A queue used to store encroached segments.  Each segment's vertices are   */
518 /*   stored so that one can check whether a segment is still the same.       */
519 
520 struct badsegment {
521   struct edge encsegment;                          /* An encroached segment. */
522   point segorg, segdest;                                /* The two vertices. */
523   struct badsegment *nextsegment;     /* Pointer to next encroached segment. */
524 };
525 
526 /* A queue used to store bad triangles.  The key is the square of the cosine */
527 /*   of the smallest angle of the triangle.  Each triangle's vertices are    */
528 /*   stored so that one can check whether a triangle is still the same.      */
529 
530 struct badface {
531   struct triedge badfacetri;                              /* A bad triangle. */
532   REAL key;                             /* cos^2 of smallest (apical) angle. */
533   point faceorg, facedest, faceapex;                  /* The three vertices. */
534   struct badface *nextface;                 /* Pointer to next bad triangle. */
535 };
536 
537 /* A node in a heap used to store events for the sweepline Delaunay          */
538 /*   algorithm.  Nodes do not point directly to their parents or children in */
539 /*   the heap.  Instead, each node knows its position in the heap, and can   */
540 /*   look up its parent and children in a separate array.  The `eventptr'    */
541 /*   points either to a `point' or to a triangle (in encoded format, so that */
542 /*   an orientation is included).  In the latter case, the origin of the     */
543 /*   oriented triangle is the apex of a "circle event" of the sweepline      */
544 /*   algorithm.  To distinguish site events from circle events, all circle   */
545 /*   events are given an invalid (smaller than `xmin') x-coordinate `xkey'.  */
546 
547 struct event {
548   REAL xkey, ykey;                              /* Coordinates of the event. */
549   VOID *eventptr;       /* Can be a point or the location of a circle event. */
550   int heapposition;              /* Marks this event's position in the heap. */
551 };
552 
553 /* A node in the splay tree.  Each node holds an oriented ghost triangle     */
554 /*   that represents a boundary edge of the growing triangulation.  When a   */
555 /*   circle event covers two boundary edges with a triangle, so that they    */
556 /*   are no longer boundary edges, those edges are not immediately deleted   */
557 /*   from the tree; rather, they are lazily deleted when they are next       */
558 /*   encountered.  (Since only a random sample of boundary edges are kept    */
559 /*   in the tree, lazy deletion is faster.)  `keydest' is used to verify     */
560 /*   that a triangle is still the same as when it entered the splay tree; if */
561 /*   it has been rotated (due to a circle event), it no longer represents a  */
562 /*   boundary edge and should be deleted.                                    */
563 
564 struct splaynode {
565   struct triedge keyedge;                  /* Lprev of an edge on the front. */
566   point keydest;            /* Used to verify that splay node is still live. */
567   struct splaynode *lchild, *rchild;              /* Children in splay tree. */
568 };
569 
570 /* A type used to allocate memory.  firstblock is the first block of items.  */
571 /*   nowblock is the block from which items are currently being allocated.   */
572 /*   nextitem points to the next slab of free memory for an item.            */
573 /*   deaditemstack is the head of a linked list (stack) of deallocated items */
574 /*   that can be recycled.  unallocateditems is the number of items that     */
575 /*   remain to be allocated from nowblock.                                   */
576 /*                                                                           */
577 /* Traversal is the process of walking through the entire list of items, and */
578 /*   is separate from allocation.  Note that a traversal will visit items on */
579 /*   the "deaditemstack" stack as well as live items.  pathblock points to   */
580 /*   the block currently being traversed.  pathitem points to the next item  */
581 /*   to be traversed.  pathitemsleft is the number of items that remain to   */
582 /*   be traversed in pathblock.                                              */
583 /*                                                                           */
584 /* itemwordtype is set to POINTER or FLOATINGPOINT, and is used to suggest   */
585 /*   what sort of word the record is primarily made up of.  alignbytes       */
586 /*   determines how new records should be aligned in memory.  itembytes and  */
587 /*   itemwords are the length of a record in bytes (after rounding up) and   */
588 /*   words.  itemsperblock is the number of items allocated at once in a     */
589 /*   single block.  items is the number of currently allocated items.        */
590 /*   maxitems is the maximum number of items that have been allocated at     */
591 /*   once; it is the current number of items plus the number of records kept */
592 /*   on deaditemstack.                                                       */
593 
594 struct memorypool {
595   VOID **firstblock, **nowblock;
596   VOID *nextitem;
597   VOID *deaditemstack;
598   VOID **pathblock;
599   VOID *pathitem;
600   enum wordtype itemwordtype;
601   int alignbytes;
602   int itembytes, itemwords;
603   int itemsperblock;
604   long items, maxitems;
605   int unallocateditems;
606   int pathitemsleft;
607 };
608 
609 /* Variables used to allocate memory for triangles, shell edges, points,     */
610 /*   viri (triangles being eaten), bad (encroached) segments, bad (skinny    */
611 /*   or too large) triangles, and splay tree nodes.                          */
612 
613 struct memorypool triangles;
614 struct memorypool shelles;
615 struct memorypool points;
616 struct memorypool viri;
617 struct memorypool badsegments;
618 struct memorypool badtriangles;
619 struct memorypool splaynodes;
620 
621 /* Variables that maintain the bad triangle queues.  The tails are pointers  */
622 /*   to the pointers that have to be filled in to enqueue an item.           */
623 
624 struct badface *queuefront[64];
625 struct badface **queuetail[64];
626 
627 REAL xmin, xmax, ymin, ymax;                              /* x and y bounds. */
628 REAL xminextreme;        /* Nonexistent x value used as a flag in sweepline. */
629 int inpoints;                                     /* Number of input points. */
630 int inelements;                                /* Number of input triangles. */
631 int insegments;                                 /* Number of input segments. */
632 int holes;                                         /* Number of input holes. */
633 int regions;                                     /* Number of input regions. */
634 long edges;                                       /* Number of output edges. */
635 int mesh_dim;                                  /* Dimension (ought to be 2). */
636 int nextras;                              /* Number of attributes per point. */
637 int eextras;                           /* Number of attributes per triangle. */
638 long hullsize;                            /* Number of edges of convex hull. */
639 int triwords;                                   /* Total words per triangle. */
640 int shwords;                                  /* Total words per shell edge. */
641 int pointmarkindex;             /* Index to find boundary marker of a point. */
642 int point2triindex;         /* Index to find a triangle adjacent to a point. */
643 int highorderindex;    /* Index to find extra nodes for high-order elements. */
644 int elemattribindex;              /* Index to find attributes of a triangle. */
645 int areaboundindex;               /* Index to find area bound of a triangle. */
646 int checksegments;           /* Are there segments in the triangulation yet? */
647 int readnodefile;                             /* Has a .node file been read? */
648 long samples;                /* Number of random samples for point location. */
649 unsigned long randomseed;                     /* Current random number seed. */
650 
651 REAL splitter;       /* Used to split REAL factors for exact multiplication. */
652 REAL epsilon;                             /* Floating-point machine epsilon. */
653 REAL resulterrbound;
654 REAL ccwerrboundA, ccwerrboundB, ccwerrboundC;
655 REAL iccerrboundA, iccerrboundB, iccerrboundC;
656 
657 long incirclecount;                   /* Number of incircle tests performed. */
658 long counterclockcount;       /* Number of counterclockwise tests performed. */
659 long hyperbolacount;        /* Number of right-of-hyperbola tests performed. */
660 long circumcentercount;    /* Number of circumcenter calculations performed. */
661 long circletopcount;         /* Number of circle top calculations performed. */
662 
663 /* Switches for the triangulator.                                            */
664 /*   poly: -p switch.  refine: -r switch.                                    */
665 /*   quality: -q switch.                                                     */
666 /*     minangle: minimum angle bound, specified after -q switch.             */
667 /*     goodangle: cosine squared of minangle.                                */
668 /*   vararea: -a switch without number.                                      */
669 /*   fixedarea: -a switch with number.                                       */
670 /*     maxarea: maximum area bound, specified after -a switch.               */
671 /*   regionattrib: -A switch.  convex: -c switch.                            */
672 /*   firstnumber: inverse of -z switch.  All items are numbered starting     */
673 /*     from firstnumber.                                                     */
674 /*   edgesout: -e switch.  voronoi: -v switch.                               */
675 /*   neighbors: -n switch.  geomview: -g switch.                             */
676 /*   nobound: -B switch.  nopolywritten: -P switch.                          */
677 /*   nonodewritten: -N switch.  noelewritten: -E switch.                     */
678 /*   noiterationnum: -I switch.  noholes: -O switch.                         */
679 /*   noexact: -X switch.                                                     */
680 /*   order: element order, specified after -o switch.                        */
681 /*   nobisect: count of how often -Y switch is selected.                     */
682 /*   steiner: maximum number of Steiner points, specified after -S switch.   */
683 /*     steinerleft: number of Steiner points not yet used.                   */
684 /*   incremental: -i switch.  sweepline: -F switch.                          */
685 /*   dwyer: inverse of -l switch.                                            */
686 /*   splitseg: -s switch.                                                    */
687 /*   docheck: -C switch.                                                     */
688 /*   quiet: -Q switch.  verbose: count of how often -V switch is selected.   */
689 /*   useshelles: -p, -r, -q, or -c switch; determines whether shell edges    */
690 /*     are used at all.                                                      */
691 /*                                                                           */
692 /* Read the instructions to find out the meaning of these switches.          */
693 
694 int poly, refine, quality, vararea, fixedarea, regionattrib, convex;
695 int firstnumber;
696 int edgesout, voronoi, neighbors, geomview;
697 int nobound, nopolywritten, nonodewritten, noelewritten, noiterationnum;
698 int noholes, noexact;
699 int incremental, sweepline, dwyer;
700 int splitseg;
701 int docheck;
702 int quiet, verbose;
703 int useshelles;
704 int order;
705 int nobisect;
706 int steiner, steinerleft;
707 REAL minangle, goodangle;
708 REAL maxarea;
709 
710 /* Variables for file names.                                                 */
711 
712 #ifndef TRILIBRARY
713 char innodefilename[FILENAMESIZE];
714 char inelefilename[FILENAMESIZE];
715 char inpolyfilename[FILENAMESIZE];
716 char areafilename[FILENAMESIZE];
717 char outnodefilename[FILENAMESIZE];
718 char outelefilename[FILENAMESIZE];
719 char outpolyfilename[FILENAMESIZE];
720 char edgefilename[FILENAMESIZE];
721 char vnodefilename[FILENAMESIZE];
722 char vedgefilename[FILENAMESIZE];
723 char neighborfilename[FILENAMESIZE];
724 char offfilename[FILENAMESIZE];
725 #endif /* not TRILIBRARY */
726 
727 /* Triangular bounding box points.                                           */
728 
729 point infpoint1, infpoint2, infpoint3;
730 
731 /* Pointer to the `triangle' that occupies all of "outer space".             */
732 
733 triangle *dummytri;
734 triangle *dummytribase;      /* Keep base address so we can free() it later. */
735 
736 /* Pointer to the omnipresent shell edge.  Referenced by any triangle or     */
737 /*   shell edge that isn't really connected to a shell edge at that          */
738 /*   location.                                                               */
739 
740 shelle *dummysh;
741 shelle *dummyshbase;         /* Keep base address so we can free() it later. */
742 
743 /* Pointer to a recently visited triangle.  Improves point location if       */
744 /*   proximate points are inserted sequentially.                             */
745 
746 struct triedge recenttri;
747 
748 /*****************************************************************************/
749 /*                                                                           */
750 /*  Mesh manipulation primitives.  Each triangle contains three pointers to  */
751 /*  other triangles, with orientations.  Each pointer points not to the      */
752 /*  first byte of a triangle, but to one of the first three bytes of a       */
753 /*  triangle.  It is necessary to extract both the triangle itself and the   */
754 /*  orientation.  To save memory, I keep both pieces of information in one   */
755 /*  pointer.  To make this possible, I assume that all triangles are aligned */
756 /*  to four-byte boundaries.  The `decode' routine below decodes a pointer,  */
757 /*  extracting an orientation (in the range 0 to 2) and a pointer to the     */
758 /*  beginning of a triangle.  The `encode' routine compresses a pointer to a */
759 /*  triangle and an orientation into a single pointer.  My assumptions that  */
760 /*  triangles are four-byte-aligned and that the `unsigned long' type is     */
761 /*  long enough to hold a pointer are two of the few kludges in this program.*/
762 /*                                                                           */
763 /*  Shell edges are manipulated similarly.  A pointer to a shell edge        */
764 /*  carries both an address and an orientation in the range 0 to 1.          */
765 /*                                                                           */
766 /*  The other primitives take an oriented triangle or oriented shell edge,   */
767 /*  and return an oriented triangle or oriented shell edge or point; or they */
768 /*  change the connections in the data structure.                            */
769 /*                                                                           */
770 /*****************************************************************************/
771 
772 /********* Mesh manipulation primitives begin here                   *********/
773 /**                                                                         **/
774 /**                                                                         **/
775 
776 /* Fast lookup arrays to speed some of the mesh manipulation primitives.     */
777 
778 int plus1mod3[3] = {1, 2, 0};
779 int minus1mod3[3] = {2, 0, 1};
780 
781 /********* Primitives for triangles                                  *********/
782 /*                                                                           */
783 /*                                                                           */
784 
785 /* decode() converts a pointer to an oriented triangle.  The orientation is  */
786 /*   extracted from the two least significant bits of the pointer.           */
787 
788 #define decode(ptr, triedge)                                                  \
789   (triedge).orient = (int) ((unsigned long) (ptr) & (unsigned long) 3l);      \
790   (triedge).tri = (triangle *)                                                \
791                   ((unsigned long) (ptr) ^ (unsigned long) (triedge).orient)
792 
793 /* encode() compresses an oriented triangle into a single pointer.  It       */
794 /*   relies on the assumption that all triangles are aligned to four-byte    */
795 /*   boundaries, so the two least significant bits of (triedge).tri are zero.*/
796 
797 #define encode(triedge)                                                       \
798   (triangle) ((unsigned long) (triedge).tri | (unsigned long) (triedge).orient)
799 
800 /* The following edge manipulation primitives are all described by Guibas    */
801 /*   and Stolfi.  However, they use an edge-based data structure, whereas I  */
802 /*   am using a triangle-based data structure.                               */
803 
804 /* sym() finds the abutting triangle, on the same edge.  Note that the       */
805 /*   edge direction is necessarily reversed, because triangle/edge handles   */
806 /*   are always directed counterclockwise around the triangle.               */
807 
808 #define sym(triedge1, triedge2)                                               \
809   ptr = (triedge1).tri[(triedge1).orient];                                    \
810   decode(ptr, triedge2);
811 
812 #define symself(triedge)                                                      \
813   ptr = (triedge).tri[(triedge).orient];                                      \
814   decode(ptr, triedge);
815 
816 /* lnext() finds the next edge (counterclockwise) of a triangle.             */
817 
818 #define lnext(triedge1, triedge2)                                             \
819   (triedge2).tri = (triedge1).tri;                                            \
820   (triedge2).orient = plus1mod3[(triedge1).orient]
821 
822 #define lnextself(triedge)                                                    \
823   (triedge).orient = plus1mod3[(triedge).orient]
824 
825 /* lprev() finds the previous edge (clockwise) of a triangle.                */
826 
827 #define lprev(triedge1, triedge2)                                             \
828   (triedge2).tri = (triedge1).tri;                                            \
829   (triedge2).orient = minus1mod3[(triedge1).orient]
830 
831 #define lprevself(triedge)                                                    \
832   (triedge).orient = minus1mod3[(triedge).orient]
833 
834 /* onext() spins counterclockwise around a point; that is, it finds the next */
835 /*   edge with the same origin in the counterclockwise direction.  This edge */
836 /*   will be part of a different triangle.                                   */
837 
838 #define onext(triedge1, triedge2)                                             \
839   lprev(triedge1, triedge2);                                                  \
840   symself(triedge2);
841 
842 #define onextself(triedge)                                                    \
843   lprevself(triedge);                                                         \
844   symself(triedge);
845 
846 /* oprev() spins clockwise around a point; that is, it finds the next edge   */
847 /*   with the same origin in the clockwise direction.  This edge will be     */
848 /*   part of a different triangle.                                           */
849 
850 #define oprev(triedge1, triedge2)                                             \
851   sym(triedge1, triedge2);                                                    \
852   lnextself(triedge2);
853 
854 #define oprevself(triedge)                                                    \
855   symself(triedge);                                                           \
856   lnextself(triedge);
857 
858 /* dnext() spins counterclockwise around a point; that is, it finds the next */
859 /*   edge with the same destination in the counterclockwise direction.  This */
860 /*   edge will be part of a different triangle.                              */
861 
862 #define dnext(triedge1, triedge2)                                             \
863   sym(triedge1, triedge2);                                                    \
864   lprevself(triedge2);
865 
866 #define dnextself(triedge)                                                    \
867   symself(triedge);                                                           \
868   lprevself(triedge);
869 
870 /* dprev() spins clockwise around a point; that is, it finds the next edge   */
871 /*   with the same destination in the clockwise direction.  This edge will   */
872 /*   be part of a different triangle.                                        */
873 
874 #define dprev(triedge1, triedge2)                                             \
875   lnext(triedge1, triedge2);                                                  \
876   symself(triedge2);
877 
878 #define dprevself(triedge)                                                    \
879   lnextself(triedge);                                                         \
880   symself(triedge);
881 
882 /* rnext() moves one edge counterclockwise about the adjacent triangle.      */
883 /*   (It's best understood by reading Guibas and Stolfi.  It involves        */
884 /*   changing triangles twice.)                                              */
885 
886 #define rnext(triedge1, triedge2)                                             \
887   sym(triedge1, triedge2);                                                    \
888   lnextself(triedge2);                                                        \
889   symself(triedge2);
890 
891 #define rnextself(triedge)                                                    \
892   symself(triedge);                                                           \
893   lnextself(triedge);                                                         \
894   symself(triedge);
895 
896 /* rnext() moves one edge clockwise about the adjacent triangle.             */
897 /*   (It's best understood by reading Guibas and Stolfi.  It involves        */
898 /*   changing triangles twice.)                                              */
899 
900 #define rprev(triedge1, triedge2)                                             \
901   sym(triedge1, triedge2);                                                    \
902   lprevself(triedge2);                                                        \
903   symself(triedge2);
904 
905 #define rprevself(triedge)                                                    \
906   symself(triedge);                                                           \
907   lprevself(triedge);                                                         \
908   symself(triedge);
909 
910 /* These primitives determine or set the origin, destination, or apex of a   */
911 /* triangle.                                                                 */
912 
913 #define org(triedge, pointptr)                                                \
914   pointptr = (point) (triedge).tri[plus1mod3[(triedge).orient] + 3]
915 
916 #define dest(triedge, pointptr)                                               \
917   pointptr = (point) (triedge).tri[minus1mod3[(triedge).orient] + 3]
918 
919 #define apex(triedge, pointptr)                                               \
920   pointptr = (point) (triedge).tri[(triedge).orient + 3]
921 
922 #define setorg(triedge, pointptr)                                             \
923   (triedge).tri[plus1mod3[(triedge).orient] + 3] = (triangle) pointptr
924 
925 #define setdest(triedge, pointptr)                                            \
926   (triedge).tri[minus1mod3[(triedge).orient] + 3] = (triangle) pointptr
927 
928 #define setapex(triedge, pointptr)                                            \
929   (triedge).tri[(triedge).orient + 3] = (triangle) pointptr
930 
931 #define setvertices2null(triedge)                                             \
932   (triedge).tri[3] = (triangle) NULL;                                         \
933   (triedge).tri[4] = (triangle) NULL;                                         \
934   (triedge).tri[5] = (triangle) NULL;
935 
936 /* Bond two triangles together.                                              */
937 
938 #define bond(triedge1, triedge2)                                              \
939   (triedge1).tri[(triedge1).orient] = encode(triedge2);                       \
940   (triedge2).tri[(triedge2).orient] = encode(triedge1)
941 
942 /* Dissolve a bond (from one side).  Note that the other triangle will still */
943 /*   think it's connected to this triangle.  Usually, however, the other     */
944 /*   triangle is being deleted entirely, or bonded to another triangle, so   */
945 /*   it doesn't matter.                                                      */
946 
947 #define dissolve(triedge)                                                     \
948   (triedge).tri[(triedge).orient] = (triangle) dummytri
949 
950 /* Copy a triangle/edge handle.                                              */
951 
952 #define triedgecopy(triedge1, triedge2)                                       \
953   (triedge2).tri = (triedge1).tri;                                            \
954   (triedge2).orient = (triedge1).orient
955 
956 /* Test for equality of triangle/edge handles.                               */
957 
958 #define triedgeequal(triedge1, triedge2)                                      \
959   (((triedge1).tri == (triedge2).tri) &&                                      \
960    ((triedge1).orient == (triedge2).orient))
961 
962 /* Primitives to infect or cure a triangle with the virus.  These rely on    */
963 /*   the assumption that all shell edges are aligned to four-byte boundaries.*/
964 
965 #define infect(triedge)                                                       \
966   (triedge).tri[6] = (triangle)                                               \
967                      ((unsigned long) (triedge).tri[6] | (unsigned long) 2l)
968 
969 #define uninfect(triedge)                                                     \
970   (triedge).tri[6] = (triangle)                                               \
971                      ((unsigned long) (triedge).tri[6] & ~ (unsigned long) 2l)
972 
973 /* Test a triangle for viral infection.                                      */
974 
975 #define infected(triedge)                                                     \
976   (((unsigned long) (triedge).tri[6] & (unsigned long) 2l) != 0)
977 
978 /* Check or set a triangle's attributes.                                     */
979 
980 #define elemattribute(triedge, attnum)                                        \
981   ((REAL *) (triedge).tri)[elemattribindex + (attnum)]
982 
983 #define setelemattribute(triedge, attnum, value)                              \
984   ((REAL *) (triedge).tri)[elemattribindex + (attnum)] = value
985 
986 /* Check or set a triangle's maximum area bound.                             */
987 
988 #define areabound(triedge)  ((REAL *) (triedge).tri)[areaboundindex]
989 
990 #define setareabound(triedge, value)                                          \
991   ((REAL *) (triedge).tri)[areaboundindex] = value
992 
993 /********* Primitives for shell edges                                *********/
994 /*                                                                           */
995 /*                                                                           */
996 
997 /* sdecode() converts a pointer to an oriented shell edge.  The orientation  */
998 /*   is extracted from the least significant bit of the pointer.  The two    */
999 /*   least significant bits (one for orientation, one for viral infection)   */
1000 /*   are masked out to produce the real pointer.                             */
1001 
1002 #define sdecode(sptr, edge)                                                   \
1003   (edge).shorient = (int) ((unsigned long) (sptr) & (unsigned long) 1l);      \
1004   (edge).sh = (shelle *)                                                      \
1005               ((unsigned long) (sptr) & ~ (unsigned long) 3l)
1006 
1007 /* sencode() compresses an oriented shell edge into a single pointer.  It    */
1008 /*   relies on the assumption that all shell edges are aligned to two-byte   */
1009 /*   boundaries, so the least significant bit of (edge).sh is zero.          */
1010 
1011 #define sencode(edge)                                                         \
1012   (shelle) ((unsigned long) (edge).sh | (unsigned long) (edge).shorient)
1013 
1014 /* ssym() toggles the orientation of a shell edge.                           */
1015 
1016 #define ssym(edge1, edge2)                                                    \
1017   (edge2).sh = (edge1).sh;                                                    \
1018   (edge2).shorient = 1 - (edge1).shorient
1019 
1020 #define ssymself(edge)                                                        \
1021   (edge).shorient = 1 - (edge).shorient
1022 
1023 /* spivot() finds the other shell edge (from the same segment) that shares   */
1024 /*   the same origin.                                                        */
1025 
1026 #define spivot(edge1, edge2)                                                  \
1027   sptr = (edge1).sh[(edge1).shorient];                                        \
1028   sdecode(sptr, edge2)
1029 
1030 #define spivotself(edge)                                                      \
1031   sptr = (edge).sh[(edge).shorient];                                          \
1032   sdecode(sptr, edge)
1033 
1034 /* snext() finds the next shell edge (from the same segment) in sequence;    */
1035 /*   one whose origin is the input shell edge's destination.                 */
1036 
1037 #define snext(edge1, edge2)                                                   \
1038   sptr = (edge1).sh[1 - (edge1).shorient];                                    \
1039   sdecode(sptr, edge2)
1040 
1041 #define snextself(edge)                                                       \
1042   sptr = (edge).sh[1 - (edge).shorient];                                      \
1043   sdecode(sptr, edge)
1044 
1045 /* These primitives determine or set the origin or destination of a shell    */
1046 /*   edge.                                                                   */
1047 
1048 #define sorg(edge, pointptr)                                                  \
1049   pointptr = (point) (edge).sh[2 + (edge).shorient]
1050 
1051 #define sdest(edge, pointptr)                                                 \
1052   pointptr = (point) (edge).sh[3 - (edge).shorient]
1053 
1054 #define setsorg(edge, pointptr)                                               \
1055   (edge).sh[2 + (edge).shorient] = (shelle) pointptr
1056 
1057 #define setsdest(edge, pointptr)                                              \
1058   (edge).sh[3 - (edge).shorient] = (shelle) pointptr
1059 
1060 /* These primitives read or set a shell marker.  Shell markers are used to   */
1061 /*   hold user boundary information.                                         */
1062 
1063 #define mark(edge)  (* (int *) ((edge).sh + 6))
1064 
1065 #define setmark(edge, value)                                                  \
1066   * (int *) ((edge).sh + 6) = value
1067 
1068 /* Bond two shell edges together.                                            */
1069 
1070 #define sbond(edge1, edge2)                                                   \
1071   (edge1).sh[(edge1).shorient] = sencode(edge2);                              \
1072   (edge2).sh[(edge2).shorient] = sencode(edge1)
1073 
1074 /* Dissolve a shell edge bond (from one side).  Note that the other shell    */
1075 /*   edge will still think it's connected to this shell edge.                */
1076 
1077 #define sdissolve(edge)                                                       \
1078   (edge).sh[(edge).shorient] = (shelle) dummysh
1079 
1080 /* Copy a shell edge.                                                        */
1081 
1082 #define shellecopy(edge1, edge2)                                              \
1083   (edge2).sh = (edge1).sh;                                                    \
1084   (edge2).shorient = (edge1).shorient
1085 
1086 /* Test for equality of shell edges.                                         */
1087 
1088 #define shelleequal(edge1, edge2)                                             \
1089   (((edge1).sh == (edge2).sh) &&                                              \
1090    ((edge1).shorient == (edge2).shorient))
1091 
1092 /********* Primitives for interacting triangles and shell edges      *********/
1093 /*                                                                           */
1094 /*                                                                           */
1095 
1096 /* tspivot() finds a shell edge abutting a triangle.                         */
1097 
1098 #define tspivot(triedge, edge)                                                \
1099   sptr = (shelle) (triedge).tri[6 + (triedge).orient];                        \
1100   sdecode(sptr, edge)
1101 
1102 /* stpivot() finds a triangle abutting a shell edge.  It requires that the   */
1103 /*   variable `ptr' of type `triangle' be defined.                           */
1104 
1105 #define stpivot(edge, triedge)                                                \
1106   ptr = (triangle) (edge).sh[4 + (edge).shorient];                            \
1107   decode(ptr, triedge)
1108 
1109 /* Bond a triangle to a shell edge.                                          */
1110 
1111 #define tsbond(triedge, edge)                                                 \
1112   (triedge).tri[6 + (triedge).orient] = (triangle) sencode(edge);             \
1113   (edge).sh[4 + (edge).shorient] = (shelle) encode(triedge)
1114 
1115 /* Dissolve a bond (from the triangle side).                                 */
1116 
1117 #define tsdissolve(triedge)                                                   \
1118   (triedge).tri[6 + (triedge).orient] = (triangle) dummysh
1119 
1120 /* Dissolve a bond (from the shell edge side).                               */
1121 
1122 #define stdissolve(edge)                                                      \
1123   (edge).sh[4 + (edge).shorient] = (shelle) dummytri
1124 
1125 /********* Primitives for points                                     *********/
1126 /*                                                                           */
1127 /*                                                                           */
1128 
1129 #define pointmark(pt)  ((int *) (pt))[pointmarkindex]
1130 
1131 #define setpointmark(pt, value)                                               \
1132   ((int *) (pt))[pointmarkindex] = value
1133 
1134 #define point2tri(pt)  ((triangle *) (pt))[point2triindex]
1135 
1136 #define setpoint2tri(pt, value)                                               \
1137   ((triangle *) (pt))[point2triindex] = value
1138 
1139 /**                                                                         **/
1140 /**                                                                         **/
1141 /********* Mesh manipulation primitives end here                     *********/
1142 
1143 /********* User interaction routines begin here                      *********/
1144 /**                                                                         **/
1145 /**                                                                         **/
1146 
1147 /*****************************************************************************/
1148 /*                                                                           */
1149 /*  syntax()   Print list of command line switches.                          */
1150 /*                                                                           */
1151 /*****************************************************************************/
1152 
1153 #ifndef TRILIBRARY
1154 
1155 void syntax()
1156 {
1157 #ifdef CDT_ONLY
1158 #ifdef REDUCED
1159   printf("triangle [-pAcevngBPNEIOXzo_lQVh] input_file\n");
1160 #else /* not REDUCED */
1161   printf("triangle [-pAcevngBPNEIOXzo_iFlCQVh] input_file\n");
1162 #endif /* not REDUCED */
1163 #else /* not CDT_ONLY */
1164 #ifdef REDUCED
1165   printf("triangle [-prq__a__AcevngBPNEIOXzo_YS__lQVh] input_file\n");
1166 #else /* not REDUCED */
1167   printf("triangle [-prq__a__AcevngBPNEIOXzo_YS__iFlsCQVh] input_file\n");
1168 #endif /* not REDUCED */
1169 #endif /* not CDT_ONLY */
1170 
1171   printf("    -p  Triangulates a Planar Straight Line Graph (.poly file).\n");
1172 #ifndef CDT_ONLY
1173   printf("    -r  Refines a previously generated mesh.\n");
1174   printf(
1175     "    -q  Quality mesh generation.  A minimum angle may be specified.\n");
1176   printf("    -a  Applies a maximum triangle area constraint.\n");
1177 #endif /* not CDT_ONLY */
1178   printf(
1179     "    -A  Applies attributes to identify elements in certain regions.\n");
1180   printf("    -c  Encloses the convex hull with segments.\n");
1181   printf("    -e  Generates an edge list.\n");
1182   printf("    -v  Generates a Voronoi diagram.\n");
1183   printf("    -n  Generates a list of triangle neighbors.\n");
1184   printf("    -g  Generates an .off file for Geomview.\n");
1185   printf("    -B  Suppresses output of boundary information.\n");
1186   printf("    -P  Suppresses output of .poly file.\n");
1187   printf("    -N  Suppresses output of .node file.\n");
1188   printf("    -E  Suppresses output of .ele file.\n");
1189   printf("    -I  Suppresses mesh iteration numbers.\n");
1190   printf("    -O  Ignores holes in .poly file.\n");
1191   printf("    -X  Suppresses use of exact arithmetic.\n");
1192   printf("    -z  Numbers all items starting from zero (rather than one).\n");
1193   printf("    -o2 Generates second-order subparametric elements.\n");
1194 #ifndef CDT_ONLY
1195   printf("    -Y  Suppresses boundary segment splitting.\n");
1196   printf("    -S  Specifies maximum number of added Steiner points.\n");
1197 #endif /* not CDT_ONLY */
1198 #ifndef REDUCED
1199   printf("    -i  Uses incremental method, rather than divide-and-conquer.\n");
1200   printf("    -F  Uses Fortune's sweepline algorithm, rather than d-and-c.\n");
1201 #endif /* not REDUCED */
1202   printf("    -l  Uses vertical cuts only, rather than alternating cuts.\n");
1203 #ifndef REDUCED
1204 #ifndef CDT_ONLY
1205   printf(
1206     "    -s  Force segments into mesh by splitting (instead of using CDT).\n");
1207 #endif /* not CDT_ONLY */
1208   printf("    -C  Check consistency of final mesh.\n");
1209 #endif /* not REDUCED */
1210   printf("    -Q  Quiet:  No terminal output except errors.\n");
1211   printf("    -V  Verbose:  Detailed information on what I'm doing.\n");
1212   printf("    -h  Help:  Detailed instructions for Triangle.\n");
1213   exit(0);
1214 }
1215 
1216 #endif /* not TRILIBRARY */
1217 
1218 /*****************************************************************************/
1219 /*                                                                           */
1220 /*  info()   Print out complete instructions.                                */
1221 /*                                                                           */
1222 /*****************************************************************************/
1223 
1224 #ifndef TRILIBRARY
1225 
1226 void info()
1227 {
1228   printf("Triangle\n");
1229   printf(
1230 "A Two-Dimensional Quality Mesh Generator and Delaunay Triangulator.\n");
1231   printf("Version 1.3\n\n");
1232   printf(
1233 "Copyright 1996 Jonathan Richard Shewchuk  (bugs/comments to jrs@cs.cmu.edu)\n"
1234 );
1235   printf("School of Computer Science / Carnegie Mellon University\n");
1236   printf("5000 Forbes Avenue / Pittsburgh, Pennsylvania  15213-3891\n");
1237   printf(
1238 "Created as part of the Archimedes project (tools for parallel FEM).\n");
1239   printf(
1240 "Supported in part by NSF Grant CMS-9318163 and an NSERC 1967 Scholarship.\n");
1241   printf("There is no warranty whatsoever.  Use at your own risk.\n");
1242 #ifdef SINGLE
1243   printf("This executable is compiled for single precision arithmetic.\n\n\n");
1244 #else /* not SINGLE */
1245   printf("This executable is compiled for double precision arithmetic.\n\n\n");
1246 #endif /* not SINGLE */
1247   printf(
1248 "Triangle generates exact Delaunay triangulations, constrained Delaunay\n");
1249   printf(
1250 "triangulations, and quality conforming Delaunay triangulations.  The latter\n"
1251 );
1252   printf(
1253 "can be generated with no small angles, and are thus suitable for finite\n");
1254   printf(
1255 "element analysis.  If no command line switches are specified, your .node\n");
1256   printf(
1257 "input file will be read, and the Delaunay triangulation will be returned in\n"
1258 );
1259   printf(".node and .ele output files.  The command syntax is:\n\n");
1260 #ifdef CDT_ONLY
1261 #ifdef REDUCED
1262   printf("triangle [-pAcevngBPNEIOXzo_lQVh] input_file\n\n");
1263 #else /* not REDUCED */
1264   printf("triangle [-pAcevngBPNEIOXzo_iFlCQVh] input_file\n\n");
1265 #endif /* not REDUCED */
1266 #else /* not CDT_ONLY */
1267 #ifdef REDUCED
1268   printf("triangle [-prq__a__AcevngBPNEIOXzo_YS__lQVh] input_file\n\n");
1269 #else /* not REDUCED */
1270   printf("triangle [-prq__a__AcevngBPNEIOXzo_YS__iFlsCQVh] input_file\n\n");
1271 #endif /* not REDUCED */
1272 #endif /* not CDT_ONLY */
1273   printf(
1274 "Underscores indicate that numbers may optionally follow certain switches;\n");
1275   printf(
1276 "do not leave any space between a switch and its numeric parameter.\n");
1277   printf(
1278 "input_file must be a file with extension .node, or extension .poly if the\n");
1279   printf(
1280 "-p switch is used.  If -r is used, you must supply .node and .ele files,\n");
1281   printf(
1282 "and possibly a .poly file and .area file as well.  The formats of these\n");
1283   printf("files are described below.\n\n");
1284   printf("Command Line Switches:\n\n");
1285   printf(
1286 "    -p  Reads a Planar Straight Line Graph (.poly file), which can specify\n"
1287 );
1288   printf(
1289 "        points, segments, holes, and regional attributes and area\n");
1290   printf(
1291 "        constraints.  Will generate a constrained Delaunay triangulation\n");
1292   printf(
1293 "        fitting the input; or, if -s, -q, or -a is used, a conforming\n");
1294   printf(
1295 "        Delaunay triangulation.  If -p is not used, Triangle reads a .node\n"
1296 );
1297   printf("        file by default.\n");
1298   printf(
1299 "    -r  Refines a previously generated mesh.  The mesh is read from a .node\n"
1300 );
1301   printf(
1302 "        file and an .ele file.  If -p is also used, a .poly file is read\n");
1303   printf(
1304 "        and used to constrain edges in the mesh.  Further details on\n");
1305   printf("        refinement are given below.\n");
1306   printf(
1307 "    -q  Quality mesh generation by Jim Ruppert's Delaunay refinement\n");
1308   printf(
1309 "        algorithm.  Adds points to the mesh to ensure that no angles\n");
1310   printf(
1311 "        smaller than 20 degrees occur.  An alternative minimum angle may be\n"
1312 );
1313   printf(
1314 "        specified after the `q'.  If the minimum angle is 20.7 degrees or\n");
1315   printf(
1316 "        smaller, the triangulation algorithm is theoretically guaranteed to\n"
1317 );
1318   printf(
1319 "        terminate (assuming infinite precision arithmetic - Triangle may\n");
1320   printf(
1321 "        fail to terminate if you run out of precision).  In practice, the\n");
1322   printf(
1323 "        algorithm often succeeds for minimum angles up to 33.8 degrees.\n");
1324   printf(
1325 "        For highly refined meshes, however, it may be necessary to reduce\n");
1326   printf(
1327 "        the minimum angle to well below 20 to avoid problems associated\n");
1328   printf(
1329 "        with insufficient floating-point precision.  The specified angle\n");
1330   printf("        may include a decimal point.\n");
1331   printf(
1332 "    -a  Imposes a maximum triangle area.  If a number follows the `a', no\n");
1333   printf(
1334 "        triangle will be generated whose area is larger than that number.\n");
1335   printf(
1336 "        If no number is specified, an .area file (if -r is used) or .poly\n");
1337   printf(
1338 "        file (if -r is not used) specifies a number of maximum area\n");
1339   printf(
1340 "        constraints.  An .area file contains a separate area constraint for\n"
1341 );
1342   printf(
1343 "        each triangle, and is useful for refining a finite element mesh\n");
1344   printf(
1345 "        based on a posteriori error estimates.  A .poly file can optionally\n"
1346 );
1347   printf(
1348 "        contain an area constraint for each segment-bounded region, thereby\n"
1349 );
1350   printf(
1351 "        enforcing triangle densities in a first triangulation.  You can\n");
1352   printf(
1353 "        impose both a fixed area constraint and a varying area constraint\n");
1354   printf(
1355 "        by invoking the -a switch twice, once with and once without a\n");
1356   printf(
1357 "        number following.  Each area specified may include a decimal point.\n"
1358 );
1359   printf(
1360 "    -A  Assigns an additional attribute to each triangle that identifies\n");
1361   printf(
1362 "        what segment-bounded region each triangle belongs to.  Attributes\n");
1363   printf(
1364 "        are assigned to regions by the .poly file.  If a region is not\n");
1365   printf(
1366 "        explicitly marked by the .poly file, triangles in that region are\n");
1367   printf(
1368 "        assigned an attribute of zero.  The -A switch has an effect only\n");
1369   printf("        when the -p switch is used and the -r switch is not.\n");
1370   printf(
1371 "    -c  Creates segments on the convex hull of the triangulation.  If you\n");
1372   printf(
1373 "        are triangulating a point set, this switch causes a .poly file to\n");
1374   printf(
1375 "        be written, containing all edges in the convex hull.  (By default,\n"
1376 );
1377   printf(
1378 "        a .poly file is written only if a .poly file is read.)  If you are\n"
1379 );
1380   printf(
1381 "        triangulating a PSLG, this switch specifies that the interior of\n");
1382   printf(
1383 "        the convex hull of the PSLG should be triangulated.  If you do not\n"
1384 );
1385   printf(
1386 "        use this switch when triangulating a PSLG, it is assumed that you\n");
1387   printf(
1388 "        have identified the region to be triangulated by surrounding it\n");
1389   printf(
1390 "        with segments of the input PSLG.  Beware:  if you are not careful,\n"
1391 );
1392   printf(
1393 "        this switch can cause the introduction of an extremely thin angle\n");
1394   printf(
1395 "        between a PSLG segment and a convex hull segment, which can cause\n");
1396   printf(
1397 "        overrefinement or failure if Triangle runs out of precision.  If\n");
1398   printf(
1399 "        you are refining a mesh, the -c switch works differently; it\n");
1400   printf(
1401 "        generates the set of boundary edges of the mesh, rather than the\n");
1402   printf("        convex hull.\n");
1403   printf(
1404 "    -e  Outputs (to an .edge file) a list of edges of the triangulation.\n");
1405   printf(
1406 "    -v  Outputs the Voronoi diagram associated with the triangulation.\n");
1407   printf("        Does not attempt to detect degeneracies.\n");
1408   printf(
1409 "    -n  Outputs (to a .neigh file) a list of triangles neighboring each\n");
1410   printf("        triangle.\n");
1411   printf(
1412 "    -g  Outputs the mesh to an Object File Format (.off) file, suitable for\n"
1413 );
1414   printf("        viewing with the Geometry Center's Geomview package.\n");
1415   printf(
1416 "    -B  No boundary markers in the output .node, .poly, and .edge output\n");
1417   printf(
1418 "        files.  See the detailed discussion of boundary markers below.\n");
1419   printf(
1420 "    -P  No output .poly file.  Saves disk space, but you lose the ability\n");
1421   printf(
1422 "        to impose segment constraints on later refinements of the mesh.\n");
1423   printf("    -N  No output .node file.\n");
1424   printf("    -E  No output .ele file.\n");
1425   printf(
1426 "    -I  No iteration numbers.  Suppresses the output of .node and .poly\n");
1427   printf(
1428 "        files, so your input files won't be overwritten.  (If your input is\n"
1429 );
1430   printf(
1431 "        a .poly file only, a .node file will be written.)  Cannot be used\n");
1432   printf(
1433 "        with the -r switch, because that would overwrite your input .ele\n");
1434   printf(
1435 "        file.  Shouldn't be used with the -s, -q, or -a switch if you are\n");
1436   printf(
1437 "        using a .node file for input, because no .node file will be\n");
1438   printf("        written, so there will be no record of any added points.\n");
1439   printf("    -O  No holes.  Ignores the holes in the .poly file.\n");
1440   printf(
1441 "    -X  No exact arithmetic.  Normally, Triangle uses exact floating-point\n"
1442 );
1443   printf(
1444 "        arithmetic for certain tests if it thinks the inexact tests are not\n"
1445 );
1446   printf(
1447 "        accurate enough.  Exact arithmetic ensures the robustness of the\n");
1448   printf(
1449 "        triangulation algorithms, despite floating-point roundoff error.\n");
1450   printf(
1451 "        Disabling exact arithmetic with the -X switch will cause a small\n");
1452   printf(
1453 "        improvement in speed and create the possibility (albeit small) that\n"
1454 );
1455   printf(
1456 "        Triangle will fail to produce a valid mesh.  Not recommended.\n");
1457   printf(
1458 "    -z  Numbers all items starting from zero (rather than one).  Note that\n"
1459 );
1460   printf(
1461 "        this switch is normally overrided by the value used to number the\n");
1462   printf(
1463 "        first point of the input .node or .poly file.  However, this switch\n"
1464 );
1465   printf("        is useful when calling Triangle from another program.\n");
1466   printf(
1467 "    -o2 Generates second-order subparametric elements with six nodes each.\n"
1468 );
1469   printf(
1470 "    -Y  No new points on the boundary.  This switch is useful when the mesh\n"
1471 );
1472   printf(
1473 "        boundary must be preserved so that it conforms to some adjacent\n");
1474   printf(
1475 "        mesh.  Be forewarned that you will probably sacrifice some of the\n");
1476   printf(
1477 "        quality of the mesh; Triangle will try, but the resulting mesh may\n"
1478 );
1479   printf(
1480 "        contain triangles of poor aspect ratio.  Works well if all the\n");
1481   printf(
1482 "        boundary points are closely spaced.  Specify this switch twice\n");
1483   printf(
1484 "        (`-YY') to prevent all segment splitting, including internal\n");
1485   printf("        boundaries.\n");
1486   printf(
1487 "    -S  Specifies the maximum number of Steiner points (points that are not\n"
1488 );
1489   printf(
1490 "        in the input, but are added to meet the constraints of minimum\n");
1491   printf(
1492 "        angle and maximum area).  The default is to allow an unlimited\n");
1493   printf(
1494 "        number.  If you specify this switch with no number after it,\n");
1495   printf(
1496 "        the limit is set to zero.  Triangle always adds points at segment\n");
1497   printf(
1498 "        intersections, even if it needs to use more points than the limit\n");
1499   printf(
1500 "        you set.  When Triangle inserts segments by splitting (-s), it\n");
1501   printf(
1502 "        always adds enough points to ensure that all the segments appear in\n"
1503 );
1504   printf(
1505 "        the triangulation, again ignoring the limit.  Be forewarned that\n");
1506   printf(
1507 "        the -S switch may result in a conforming triangulation that is not\n"
1508 );
1509   printf(
1510 "        truly Delaunay, because Triangle may be forced to stop adding\n");
1511   printf(
1512 "        points when the mesh is in a state where a segment is non-Delaunay\n"
1513 );
1514   printf(
1515 "        and needs to be split.  If so, Triangle will print a warning.\n");
1516   printf(
1517 "    -i  Uses an incremental rather than divide-and-conquer algorithm to\n");
1518   printf(
1519 "        form a Delaunay triangulation.  Try it if the divide-and-conquer\n");
1520   printf("        algorithm fails.\n");
1521   printf(
1522 "    -F  Uses Steven Fortune's sweepline algorithm to form a Delaunay\n");
1523   printf(
1524 "        triangulation.  Warning:  does not use exact arithmetic for all\n");
1525   printf("        calculations.  An exact result is not guaranteed.\n");
1526   printf(
1527 "    -l  Uses only vertical cuts in the divide-and-conquer algorithm.  By\n");
1528   printf(
1529 "        default, Triangle uses alternating vertical and horizontal cuts,\n");
1530   printf(
1531 "        which usually improve the speed except with point sets that are\n");
1532   printf(
1533 "        small or short and wide.  This switch is primarily of theoretical\n");
1534   printf("        interest.\n");
1535   printf(
1536 "    -s  Specifies that segments should be forced into the triangulation by\n"
1537 );
1538   printf(
1539 "        recursively splitting them at their midpoints, rather than by\n");
1540   printf(
1541 "        generating a constrained Delaunay triangulation.  Segment splitting\n"
1542 );
1543   printf(
1544 "        is true to Ruppert's original algorithm, but can create needlessly\n"
1545 );
1546   printf("        small triangles near external small features.\n");
1547   printf(
1548 "    -C  Check the consistency of the final mesh.  Uses exact arithmetic for\n"
1549 );
1550   printf(
1551 "        checking, even if the -X switch is used.  Useful if you suspect\n");
1552   printf("        Triangle is buggy.\n");
1553   printf(
1554 "    -Q  Quiet: Suppresses all explanation of what Triangle is doing, unless\n"
1555 );
1556   printf("        an error occurs.\n");
1557   printf(
1558 "    -V  Verbose: Gives detailed information about what Triangle is doing.\n");
1559   printf(
1560 "        Add more `V's for increasing amount of detail.  `-V' gives\n");
1561   printf(
1562 "        information on algorithmic progress and more detailed statistics.\n");
1563   printf(
1564 "        `-VV' gives point-by-point details, and will print so much that\n");
1565   printf(
1566 "        Triangle will run much more slowly.  `-VVV' gives information only\n"
1567 );
1568   printf("        a debugger could love.\n");
1569   printf("    -h  Help:  Displays these instructions.\n");
1570   printf("\n");
1571   printf("Definitions:\n");
1572   printf("\n");
1573   printf(
1574 "  A Delaunay triangulation of a point set is a triangulation whose vertices\n"
1575 );
1576   printf(
1577 "  are the point set, having the property that no point in the point set\n");
1578   printf(
1579 "  falls in the interior of the circumcircle (circle that passes through all\n"
1580 );
1581   printf("  three vertices) of any triangle in the triangulation.\n\n");
1582   printf(
1583 "  A Voronoi diagram of a point set is a subdivision of the plane into\n");
1584   printf(
1585 "  polygonal regions (some of which may be infinite), where each region is\n");
1586   printf(
1587 "  the set of points in the plane that are closer to some input point than\n");
1588   printf(
1589 "  to any other input point.  (The Voronoi diagram is the geometric dual of\n"
1590 );
1591   printf("  the Delaunay triangulation.)\n\n");
1592   printf(
1593 "  A Planar Straight Line Graph (PSLG) is a collection of points and\n");
1594   printf(
1595 "  segments.  Segments are simply edges, whose endpoints are points in the\n");
1596   printf(
1597 "  PSLG.  The file format for PSLGs (.poly files) is described below.\n");
1598   printf("\n");
1599   printf(
1600 "  A constrained Delaunay triangulation of a PSLG is similar to a Delaunay\n");
1601   printf(
1602 "  triangulation, but each PSLG segment is present as a single edge in the\n");
1603   printf(
1604 "  triangulation.  (A constrained Delaunay triangulation is not truly a\n");
1605   printf("  Delaunay triangulation.)\n\n");
1606   printf(
1607 "  A conforming Delaunay triangulation of a PSLG is a true Delaunay\n");
1608   printf(
1609 "  triangulation in which each PSLG segment may have been subdivided into\n");
1610   printf(
1611 "  several edges by the insertion of additional points.  These inserted\n");
1612   printf(
1613 "  points are necessary to allow the segments to exist in the mesh while\n");
1614   printf("  maintaining the Delaunay property.\n\n");
1615   printf("File Formats:\n\n");
1616   printf(
1617 "  All files may contain comments prefixed by the character '#'.  Points,\n");
1618   printf(
1619 "  triangles, edges, holes, and maximum area constraints must be numbered\n");
1620   printf(
1621 "  consecutively, starting from either 1 or 0.  Whichever you choose, all\n");
1622   printf(
1623 "  input files must be consistent; if the nodes are numbered from 1, so must\n"
1624 );
1625   printf(
1626 "  be all other objects.  Triangle automatically detects your choice while\n");
1627   printf(
1628 "  reading the .node (or .poly) file.  (When calling Triangle from another\n");
1629   printf(
1630 "  program, use the -z switch if you wish to number objects from zero.)\n");
1631   printf("  Examples of these file formats are given below.\n\n");
1632   printf("  .node files:\n");
1633   printf(
1634 "    First line:  <# of points> <dimension (must be 2)> <# of attributes>\n");
1635   printf(
1636 "                                           <# of boundary markers (0 or 1)>\n"
1637 );
1638   printf(
1639 "    Remaining lines:  <point #> <x> <y> [attributes] [boundary marker]\n");
1640   printf("\n");
1641   printf(
1642 "    The attributes, which are typically floating-point values of physical\n");
1643   printf(
1644 "    quantities (such as mass or conductivity) associated with the nodes of\n"
1645 );
1646   printf(
1647 "    a finite element mesh, are copied unchanged to the output mesh.  If -s,\n"
1648 );
1649   printf(
1650 "    -q, or -a is selected, each new Steiner point added to the mesh will\n");
1651   printf("    have attributes assigned to it by linear interpolation.\n\n");
1652   printf(
1653 "    If the fourth entry of the first line is `1', the last column of the\n");
1654   printf(
1655 "    remainder of the file is assumed to contain boundary markers.  Boundary\n"
1656 );
1657   printf(
1658 "    markers are used to identify boundary points and points resting on PSLG\n"
1659 );
1660   printf(
1661 "    segments; a complete description appears in a section below.  The .node\n"
1662 );
1663   printf(
1664 "    file produced by Triangle will contain boundary markers in the last\n");
1665   printf("    column unless they are suppressed by the -B switch.\n\n");
1666   printf("  .ele files:\n");
1667   printf(
1668 "    First line:  <# of triangles> <points per triangle> <# of attributes>\n");
1669   printf(
1670 "    Remaining lines:  <triangle #> <point> <point> <point> ... [attributes]\n"
1671 );
1672   printf("\n");
1673   printf(
1674 "    Points are indices into the corresponding .node file.  The first three\n"
1675 );
1676   printf(
1677 "    points are the corners, and are listed in counterclockwise order around\n"
1678 );
1679   printf(
1680 "    each triangle.  (The remaining points, if any, depend on the type of\n");
1681   printf(
1682 "    finite element used.)  The attributes are just like those of .node\n");
1683   printf(
1684 "    files.  Because there is no simple mapping from input to output\n");
1685   printf(
1686 "    triangles, an attempt is made to interpolate attributes, which may\n");
1687   printf(
1688 "    result in a good deal of diffusion of attributes among nearby triangles\n"
1689 );
1690   printf(
1691 "    as the triangulation is refined.  Diffusion does not occur across\n");
1692   printf(
1693 "    segments, so attributes used to identify segment-bounded regions remain\n"
1694 );
1695   printf(
1696 "    intact.  In output .ele files, all triangles have three points each\n");
1697   printf(
1698 "    unless the -o2 switch is used, in which case they have six, and the\n");
1699   printf(
1700 "    fourth, fifth, and sixth points lie on the midpoints of the edges\n");
1701   printf("    opposite the first, second, and third corners.\n\n");
1702   printf("  .poly files:\n");
1703   printf(
1704 "    First line:  <# of points> <dimension (must be 2)> <# of attributes>\n");
1705   printf(
1706 "                                           <# of boundary markers (0 or 1)>\n"
1707 );
1708   printf(
1709 "    Following lines:  <point #> <x> <y> [attributes] [boundary marker]\n");
1710   printf("    One line:  <# of segments> <# of boundary markers (0 or 1)>\n");
1711   printf(
1712 "    Following lines:  <segment #> <endpoint> <endpoint> [boundary marker]\n");
1713   printf("    One line:  <# of holes>\n");
1714   printf("    Following lines:  <hole #> <x> <y>\n");
1715   printf(
1716 "    Optional line:  <# of regional attributes and/or area constraints>\n");
1717   printf(
1718 "    Optional following lines:  <constraint #> <x> <y> <attrib> <max area>\n");
1719   printf("\n");
1720   printf(
1721 "    A .poly file represents a PSLG, as well as some additional information.\n"
1722 );
1723   printf(
1724 "    The first section lists all the points, and is identical to the format\n"
1725 );
1726   printf(
1727 "    of .node files.  <# of points> may be set to zero to indicate that the\n"
1728 );
1729   printf(
1730 "    points are listed in a separate .node file; .poly files produced by\n");
1731   printf(
1732 "    Triangle always have this format.  This has the advantage that a point\n"
1733 );
1734   printf(
1735 "    set may easily be triangulated with or without segments.  (The same\n");
1736   printf(
1737 "    effect can be achieved, albeit using more disk space, by making a copy\n"
1738 );
1739   printf(
1740 "    of the .poly file with the extension .node; all sections of the file\n");
1741   printf("    but the first are ignored.)\n\n");
1742   printf(
1743 "    The second section lists the segments.  Segments are edges whose\n");
1744   printf(
1745 "    presence in the triangulation is enforced.  Each segment is specified\n");
1746   printf(
1747 "    by listing the indices of its two endpoints.  This means that you must\n"
1748 );
1749   printf(
1750 "    include its endpoints in the point list.  If -s, -q, and -a are not\n");
1751   printf(
1752 "    selected, Triangle will produce a constrained Delaunay triangulation,\n");
1753   printf(
1754 "    in which each segment appears as a single edge in the triangulation.\n");
1755   printf(
1756 "    If -q or -a is selected, Triangle will produce a conforming Delaunay\n");
1757   printf(
1758 "    triangulation, in which segments may be subdivided into smaller edges.\n"
1759 );
1760   printf("    Each segment, like each point, may have a boundary marker.\n\n");
1761   printf(
1762 "    The third section lists holes (and concavities, if -c is selected) in\n");
1763   printf(
1764 "    the triangulation.  Holes are specified by identifying a point inside\n");
1765   printf(
1766 "    each hole.  After the triangulation is formed, Triangle creates holes\n");
1767   printf(
1768 "    by eating triangles, spreading out from each hole point until its\n");
1769   printf(
1770 "    progress is blocked by PSLG segments; you must be careful to enclose\n");
1771   printf(
1772 "    each hole in segments, or your whole triangulation may be eaten away.\n");
1773   printf(
1774 "    If the two triangles abutting a segment are eaten, the segment itself\n");
1775   printf(
1776 "    is also eaten.  Do not place a hole directly on a segment; if you do,\n");
1777   printf("    Triangle will choose one side of the segment arbitrarily.\n\n");
1778   printf(
1779 "    The optional fourth section lists regional attributes (to be assigned\n");
1780   printf(
1781 "    to all triangles in a region) and regional constraints on the maximum\n");
1782   printf(
1783 "    triangle area.  Triangle will read this section only if the -A switch\n");
1784   printf(
1785 "    is used or the -a switch is used without a number following it, and the\n"
1786 );
1787   printf(
1788 "    -r switch is not used.  Regional attributes and area constraints are\n");
1789   printf(
1790 "    propagated in the same manner as holes; you specify a point for each\n");
1791   printf(
1792 "    attribute and/or constraint, and the attribute and/or constraint will\n");
1793   printf(
1794 "    affect the whole region (bounded by segments) containing the point.  If\n"
1795 );
1796   printf(
1797 "    two values are written on a line after the x and y coordinate, the\n");
1798   printf(
1799 "    former is assumed to be a regional attribute (but will only be applied\n"
1800 );
1801   printf(
1802 "    if the -A switch is selected), and the latter is assumed to be a\n");
1803   printf(
1804 "    regional area constraint (but will only be applied if the -a switch is\n"
1805 );
1806   printf(
1807 "    selected).  You may also specify just one value after the coordinates,\n"
1808 );
1809   printf(
1810 "    which can serve as both an attribute and an area constraint, depending\n"
1811 );
1812   printf(
1813 "    on the choice of switches.  If you are using the -A and -a switches\n");
1814   printf(
1815 "    simultaneously and wish to assign an attribute to some region without\n");
1816   printf("    imposing an area constraint, use a negative maximum area.\n\n");
1817   printf(
1818 "    When a triangulation is created from a .poly file, you must either\n");
1819   printf(
1820 "    enclose the entire region to be triangulated in PSLG segments, or\n");
1821   printf(
1822 "    use the -c switch, which encloses the convex hull of the input point\n");
1823   printf(
1824 "    set.  If you do not use the -c switch, Triangle will eat all triangles\n"
1825 );
1826   printf(
1827 "    on the outer boundary that are not protected by segments; if you are\n");
1828   printf(
1829 "    not careful, your whole triangulation may be eaten away.  If you do\n");
1830   printf(
1831 "    use the -c switch, you can still produce concavities by appropriate\n");
1832   printf("    placement of holes just inside the convex hull.\n\n");
1833   printf(
1834 "    An ideal PSLG has no intersecting segments, nor any points that lie\n");
1835   printf(
1836 "    upon segments (except, of course, the endpoints of each segment.)  You\n"
1837 );
1838   printf(
1839 "    aren't required to make your .poly files ideal, but you should be aware\n"
1840 );
1841   printf(
1842 "    of what can go wrong.  Segment intersections are relatively safe -\n");
1843   printf(
1844 "    Triangle will calculate the intersection points for you and add them to\n"
1845 );
1846   printf(
1847 "    the triangulation - as long as your machine's floating-point precision\n"
1848 );
1849   printf(
1850 "    doesn't become a problem.  You are tempting the fates if you have three\n"
1851 );
1852   printf(
1853 "    segments that cross at the same location, and expect Triangle to figure\n"
1854 );
1855   printf(
1856 "    out where the intersection point is.  Thanks to floating-point roundoff\n"
1857 );
1858   printf(
1859 "    error, Triangle will probably decide that the three segments intersect\n"
1860 );
1861   printf(
1862 "    at three different points, and you will find a minuscule triangle in\n");
1863   printf(
1864 "    your output - unless Triangle tries to refine the tiny triangle, uses\n");
1865   printf(
1866 "    up the last bit of machine precision, and fails to terminate at all.\n");
1867   printf(
1868 "    You're better off putting the intersection point in the input files,\n");
1869   printf(
1870 "    and manually breaking up each segment into two.  Similarly, if you\n");
1871   printf(
1872 "    place a point at the middle of a segment, and hope that Triangle will\n");
1873   printf(
1874 "    break up the segment at that point, you might get lucky.  On the other\n"
1875 );
1876   printf(
1877 "    hand, Triangle might decide that the point doesn't lie precisely on the\n"
1878 );
1879   printf(
1880 "    line, and you'll have a needle-sharp triangle in your output - or a lot\n"
1881 );
1882   printf("    of tiny triangles if you're generating a quality mesh.\n\n");
1883   printf(
1884 "    When Triangle reads a .poly file, it also writes a .poly file, which\n");
1885   printf(
1886 "    includes all edges that are part of input segments.  If the -c switch\n");
1887   printf(
1888 "    is used, the output .poly file will also include all of the edges on\n");
1889   printf(
1890 "    the convex hull.  Hence, the output .poly file is useful for finding\n");
1891   printf(
1892 "    edges associated with input segments and setting boundary conditions in\n"
1893 );
1894   printf(
1895 "    finite element simulations.  More importantly, you will need it if you\n"
1896 );
1897   printf(
1898 "    plan to refine the output mesh, and don't want segments to be missing\n");
1899   printf("    in later triangulations.\n\n");
1900   printf("  .area files:\n");
1901   printf("    First line:  <# of triangles>\n");
1902   printf("    Following lines:  <triangle #> <maximum area>\n\n");
1903   printf(
1904 "    An .area file associates with each triangle a maximum area that is used\n"
1905 );
1906   printf(
1907 "    for mesh refinement.  As with other file formats, every triangle must\n");
1908   printf(
1909 "    be represented, and they must be numbered consecutively.  A triangle\n");
1910   printf(
1911 "    may be left unconstrained by assigning it a negative maximum area.\n");
1912   printf("\n");
1913   printf("  .edge files:\n");
1914   printf("    First line:  <# of edges> <# of boundary markers (0 or 1)>\n");
1915   printf(
1916 "    Following lines:  <edge #> <endpoint> <endpoint> [boundary marker]\n");
1917   printf("\n");
1918   printf(
1919 "    Endpoints are indices into the corresponding .node file.  Triangle can\n"
1920 );
1921   printf(
1922 "    produce .edge files (use the -e switch), but cannot read them.  The\n");
1923   printf(
1924 "    optional column of boundary markers is suppressed by the -B switch.\n");
1925   printf("\n");
1926   printf(
1927 "    In Voronoi diagrams, one also finds a special kind of edge that is an\n");
1928   printf(
1929 "    infinite ray with only one endpoint.  For these edges, a different\n");
1930   printf("    format is used:\n\n");
1931   printf("        <edge #> <endpoint> -1 <direction x> <direction y>\n\n");
1932   printf(
1933 "    The `direction' is a floating-point vector that indicates the direction\n"
1934 );
1935   printf("    of the infinite ray.\n\n");
1936   printf("  .neigh files:\n");
1937   printf(
1938 "    First line:  <# of triangles> <# of neighbors per triangle (always 3)>\n"
1939 );
1940   printf(
1941 "    Following lines:  <triangle #> <neighbor> <neighbor> <neighbor>\n");
1942   printf("\n");
1943   printf(
1944 "    Neighbors are indices into the corresponding .ele file.  An index of -1\n"
1945 );
1946   printf(
1947 "    indicates a mesh boundary, and therefore no neighbor.  Triangle can\n");
1948   printf(
1949 "    produce .neigh files (use the -n switch), but cannot read them.\n");
1950   printf("\n");
1951   printf(
1952 "    The first neighbor of triangle i is opposite the first corner of\n");
1953   printf("    triangle i, and so on.\n\n");
1954   printf("Boundary Markers:\n\n");
1955   printf(
1956 "  Boundary markers are tags used mainly to identify which output points and\n"
1957 );
1958   printf(
1959 "  edges are associated with which PSLG segment, and to identify which\n");
1960   printf(
1961 "  points and edges occur on a boundary of the triangulation.  A common use\n"
1962 );
1963   printf(
1964 "  is to determine where boundary conditions should be applied to a finite\n");
1965   printf(
1966 "  element mesh.  You can prevent boundary markers from being written into\n");
1967   printf("  files produced by Triangle by using the -B switch.\n\n");
1968   printf(
1969 "  The boundary marker associated with each segment in an output .poly file\n"
1970 );
1971   printf("  or edge in an output .edge file is chosen as follows:\n");
1972   printf(
1973 "    - If an output edge is part or all of a PSLG segment with a nonzero\n");
1974   printf(
1975 "      boundary marker, then the edge is assigned the same marker.\n");
1976   printf(
1977 "    - Otherwise, if the edge occurs on a boundary of the triangulation\n");
1978   printf(
1979 "      (including boundaries of holes), then the edge is assigned the marker\n"
1980 );
1981   printf("      one (1).\n");
1982   printf("    - Otherwise, the edge is assigned the marker zero (0).\n");
1983   printf(
1984 "  The boundary marker associated with each point in an output .node file is\n"
1985 );
1986   printf("  chosen as follows:\n");
1987   printf(
1988 "    - If a point is assigned a nonzero boundary marker in the input file,\n");
1989   printf(
1990 "      then it is assigned the same marker in the output .node file.\n");
1991   printf(
1992 "    - Otherwise, if the point lies on a PSLG segment (including the\n");
1993   printf(
1994 "      segment's endpoints) with a nonzero boundary marker, then the point\n");
1995   printf(
1996 "      is assigned the same marker.  If the point lies on several such\n");
1997   printf("      segments, one of the markers is chosen arbitrarily.\n");
1998   printf(
1999 "    - Otherwise, if the point occurs on a boundary of the triangulation,\n");
2000   printf("      then the point is assigned the marker one (1).\n");
2001   printf("    - Otherwise, the point is assigned the marker zero (0).\n");
2002   printf("\n");
2003   printf(
2004 "  If you want Triangle to determine for you which points and edges are on\n");
2005   printf(
2006 "  the boundary, assign them the boundary marker zero (or use no markers at\n"
2007 );
2008   printf(
2009 "  all) in your input files.  Alternatively, you can mark some of them and\n");
2010   printf("  leave others marked zero, allowing Triangle to label them.\n\n");
2011   printf("Triangulation Iteration Numbers:\n\n");
2012   printf(
2013 "  Because Triangle can read and refine its own triangulations, input\n");
2014   printf(
2015 "  and output files have iteration numbers.  For instance, Triangle might\n");
2016   printf(
2017 "  read the files mesh.3.node, mesh.3.ele, and mesh.3.poly, refine the\n");
2018   printf(
2019 "  triangulation, and output the files mesh.4.node, mesh.4.ele, and\n");
2020   printf("  mesh.4.poly.  Files with no iteration number are treated as if\n");
2021   printf(
2022 "  their iteration number is zero; hence, Triangle might read the file\n");
2023   printf(
2024 "  points.node, triangulate it, and produce the files points.1.node and\n");
2025   printf("  points.1.ele.\n\n");
2026   printf(
2027 "  Iteration numbers allow you to create a sequence of successively finer\n");
2028   printf(
2029 "  meshes suitable for multigrid methods.  They also allow you to produce a\n"
2030 );
2031   printf(
2032 "  sequence of meshes using error estimate-driven mesh refinement.\n");
2033   printf("\n");
2034   printf(
2035 "  If you're not using refinement or quality meshing, and you don't like\n");
2036   printf(
2037 "  iteration numbers, use the -I switch to disable them.  This switch will\n");
2038   printf(
2039 "  also disable output of .node and .poly files to prevent your input files\n"
2040 );
2041   printf(
2042 "  from being overwritten.  (If the input is a .poly file that contains its\n"
2043 );
2044   printf("  own points, a .node file will be written.)\n\n");
2045   printf("Examples of How to Use Triangle:\n\n");
2046   printf(
2047 "  `triangle dots' will read points from dots.node, and write their Delaunay\n"
2048 );
2049   printf(
2050 "  triangulation to dots.1.node and dots.1.ele.  (dots.1.node will be\n");
2051   printf(
2052 "  identical to dots.node.)  `triangle -I dots' writes the triangulation to\n"
2053 );
2054   printf(
2055 "  dots.ele instead.  (No additional .node file is needed, so none is\n");
2056   printf("  written.)\n\n");
2057   printf(
2058 "  `triangle -pe object.1' will read a PSLG from object.1.poly (and possibly\n"
2059 );
2060   printf(
2061 "  object.1.node, if the points are omitted from object.1.poly) and write\n");
2062   printf("  their constrained Delaunay triangulation to object.2.node and\n");
2063   printf(
2064 "  object.2.ele.  The segments will be copied to object.2.poly, and all\n");
2065   printf("  edges will be written to object.2.edge.\n\n");
2066   printf(
2067 "  `triangle -pq31.5a.1 object' will read a PSLG from object.poly (and\n");
2068   printf(
2069 "  possibly object.node), generate a mesh whose angles are all greater than\n"
2070 );
2071   printf(
2072 "  31.5 degrees and whose triangles all have area smaller than 0.1, and\n");
2073   printf(
2074 "  write the mesh to object.1.node and object.1.ele.  Each segment may have\n"
2075 );
2076   printf(
2077 "  been broken up into multiple edges; the resulting constrained edges are\n");
2078   printf("  written to object.1.poly.\n\n");
2079   printf(
2080 "  Here is a sample file `box.poly' describing a square with a square hole:\n"
2081 );
2082   printf("\n");
2083   printf(
2084 "    # A box with eight points in 2D, no attributes, one boundary marker.\n");
2085   printf("    8 2 0 1\n");
2086   printf("    # Outer box has these vertices:\n");
2087   printf("     1   0 0   0\n");
2088   printf("     2   0 3   0\n");
2089   printf("     3   3 0   0\n");
2090   printf("     4   3 3   33     # A special marker for this point.\n");
2091   printf("    # Inner square has these vertices:\n");
2092   printf("     5   1 1   0\n");
2093   printf("     6   1 2   0\n");
2094   printf("     7   2 1   0\n");
2095   printf("     8   2 2   0\n");
2096   printf("    # Five segments with boundary markers.\n");
2097   printf("    5 1\n");
2098   printf("     1   1 2   5      # Left side of outer box.\n");
2099   printf("     2   5 7   0      # Segments 2 through 5 enclose the hole.\n");
2100   printf("     3   7 8   0\n");
2101   printf("     4   8 6   10\n");
2102   printf("     5   6 5   0\n");
2103   printf("    # One hole in the middle of the inner square.\n");
2104   printf("    1\n");
2105   printf("     1   1.5 1.5\n\n");
2106   printf(
2107 "  Note that some segments are missing from the outer square, so one must\n");
2108   printf(
2109 "  use the `-c' switch.  After `triangle -pqc box.poly', here is the output\n"
2110 );
2111   printf(
2112 "  file `box.1.node', with twelve points.  The last four points were added\n");
2113   printf(
2114 "  to meet the angle constraint.  Points 1, 2, and 9 have markers from\n");
2115   printf(
2116 "  segment 1.  Points 6 and 8 have markers from segment 4.  All the other\n");
2117   printf(
2118 "  points but 4 have been marked to indicate that they lie on a boundary.\n");
2119   printf("\n");
2120   printf("    12  2  0  1\n");
2121   printf("       1    0   0      5\n");
2122   printf("       2    0   3      5\n");
2123   printf("       3    3   0      1\n");
2124   printf("       4    3   3     33\n");
2125   printf("       5    1   1      1\n");
2126   printf("       6    1   2     10\n");
2127   printf("       7    2   1      1\n");
2128   printf("       8    2   2     10\n");
2129   printf("       9    0   1.5    5\n");
2130   printf("      10    1.5   0    1\n");
2131   printf("      11    3   1.5    1\n");
2132   printf("      12    1.5   3    1\n");
2133   printf("    # Generated by triangle -pqc box.poly\n\n");
2134   printf("  Here is the output file `box.1.ele', with twelve triangles.\n\n");
2135   printf("    12  3  0\n");
2136   printf("       1     5   6   9\n");
2137   printf("       2    10   3   7\n");
2138   printf("       3     6   8  12\n");
2139   printf("       4     9   1   5\n");
2140   printf("       5     6   2   9\n");
2141   printf("       6     7   3  11\n");
2142   printf("       7    11   4   8\n");
2143   printf("       8     7   5  10\n");
2144   printf("       9    12   2   6\n");
2145   printf("      10     8   7  11\n");
2146   printf("      11     5   1  10\n");
2147   printf("      12     8   4  12\n");
2148   printf("    # Generated by triangle -pqc box.poly\n\n");
2149   printf(
2150 "  Here is the output file `box.1.poly'.  Note that segments have been added\n"
2151 );
2152   printf(
2153 "  to represent the convex hull, and some segments have been split by newly\n"
2154 );
2155   printf(
2156 "  added points.  Note also that <# of points> is set to zero to indicate\n");
2157   printf("  that the points should be read from the .node file.\n\n");
2158   printf("    0  2  0  1\n");
2159   printf("    12  1\n");
2160   printf("       1     1   9     5\n");
2161   printf("       2     5   7     1\n");
2162   printf("       3     8   7     1\n");
2163   printf("       4     6   8    10\n");
2164   printf("       5     5   6     1\n");
2165   printf("       6     3  10     1\n");
2166   printf("       7     4  11     1\n");
2167   printf("       8     2  12     1\n");
2168   printf("       9     9   2     5\n");
2169   printf("      10    10   1     1\n");
2170   printf("      11    11   3     1\n");
2171   printf("      12    12   4     1\n");
2172   printf("    1\n");
2173   printf("       1   1.5 1.5\n");
2174   printf("    # Generated by triangle -pqc box.poly\n\n");
2175   printf("Refinement and Area Constraints:\n\n");
2176   printf(
2177 "  The -r switch causes a mesh (.node and .ele files) to be read and\n");
2178   printf(
2179 "  refined.  If the -p switch is also used, a .poly file is read and used to\n"
2180 );
2181   printf(
2182 "  specify edges that are constrained and cannot be eliminated (although\n");
2183   printf(
2184 "  they can be divided into smaller edges) by the refinement process.\n");
2185   printf("\n");
2186   printf(
2187 "  When you refine a mesh, you generally want to impose tighter quality\n");
2188   printf(
2189 "  constraints.  One way to accomplish this is to use -q with a larger\n");
2190   printf(
2191 "  angle, or -a followed by a smaller area than you used to generate the\n");
2192   printf(
2193 "  mesh you are refining.  Another way to do this is to create an .area\n");
2194   printf(
2195 "  file, which specifies a maximum area for each triangle, and use the -a\n");
2196   printf(
2197 "  switch (without a number following).  Each triangle's area constraint is\n"
2198 );
2199   printf(
2200 "  applied to that triangle.  Area constraints tend to diffuse as the mesh\n");
2201   printf(
2202 "  is refined, so if there are large variations in area constraint between\n");
2203   printf("  adjacent triangles, you may not get the results you want.\n\n");
2204   printf(
2205 "  If you are refining a mesh composed of linear (three-node) elements, the\n"
2206 );
2207   printf(
2208 "  output mesh will contain all the nodes present in the input mesh, in the\n"
2209 );
2210   printf(
2211 "  same order, with new nodes added at the end of the .node file.  However,\n"
2212 );
2213   printf(
2214 "  there is no guarantee that each output element is contained in a single\n");
2215   printf(
2216 "  input element.  Often, output elements will overlap two input elements,\n");
2217   printf(
2218 "  and input edges are not present in the output mesh.  Hence, a sequence of\n"
2219 );
2220   printf(
2221 "  refined meshes will form a hierarchy of nodes, but not a hierarchy of\n");
2222   printf(
2223 "  elements.  If you a refining a mesh of higher-order elements, the\n");
2224   printf(
2225 "  hierarchical property applies only to the nodes at the corners of an\n");
2226   printf("  element; other nodes may not be present in the refined mesh.\n\n");
2227   printf(
2228 "  It is important to understand that maximum area constraints in .poly\n");
2229   printf(
2230 "  files are handled differently from those in .area files.  A maximum area\n"
2231 );
2232   printf(
2233 "  in a .poly file applies to the whole (segment-bounded) region in which a\n"
2234 );
2235   printf(
2236 "  point falls, whereas a maximum area in an .area file applies to only one\n"
2237 );
2238   printf(
2239 "  triangle.  Area constraints in .poly files are used only when a mesh is\n");
2240   printf(
2241 "  first generated, whereas area constraints in .area files are used only to\n"
2242 );
2243   printf(
2244 "  refine an existing mesh, and are typically based on a posteriori error\n");
2245   printf(
2246 "  estimates resulting from a finite element simulation on that mesh.\n");
2247   printf("\n");
2248   printf(
2249 "  `triangle -rq25 object.1' will read object.1.node and object.1.ele, then\n"
2250 );
2251   printf(
2252 "  refine the triangulation to enforce a 25 degree minimum angle, and then\n");
2253   printf(
2254 "  write the refined triangulation to object.2.node and object.2.ele.\n");
2255   printf("\n");
2256   printf(
2257 "  `triangle -rpaa6.2 z.3' will read z.3.node, z.3.ele, z.3.poly, and\n");
2258   printf(
2259 "  z.3.area.  After reconstructing the mesh and its segments, Triangle will\n"
2260 );
2261   printf(
2262 "  refine the mesh so that no triangle has area greater than 6.2, and\n");
2263   printf(
2264 "  furthermore the triangles satisfy the maximum area constraints in\n");
2265   printf(
2266 "  z.3.area.  The output is written to z.4.node, z.4.ele, and z.4.poly.\n");
2267   printf("\n");
2268   printf(
2269 "  The sequence `triangle -qa1 x', `triangle -rqa.3 x.1', `triangle -rqa.1\n");
2270   printf(
2271 "  x.2' creates a sequence of successively finer meshes x.1, x.2, and x.3,\n");
2272   printf("  suitable for multigrid.\n\n");
2273   printf("Convex Hulls and Mesh Boundaries:\n\n");
2274   printf(
2275 "  If the input is a point set (rather than a PSLG), Triangle produces its\n");
2276   printf(
2277 "  convex hull as a by-product in the output .poly file if you use the -c\n");
2278   printf(
2279 "  switch.  There are faster algorithms for finding a two-dimensional convex\n"
2280 );
2281   printf(
2282 "  hull than triangulation, of course, but this one comes for free.  If the\n"
2283 );
2284   printf(
2285 "  input is an unconstrained mesh (you are using the -r switch but not the\n");
2286   printf(
2287 "  -p switch), Triangle produces a list of its boundary edges (including\n");
2288   printf("  hole boundaries) as a by-product if you use the -c switch.\n\n");
2289   printf("Voronoi Diagrams:\n\n");
2290   printf(
2291 "  The -v switch produces a Voronoi diagram, in files suffixed .v.node and\n");
2292   printf(
2293 "  .v.edge.  For example, `triangle -v points' will read points.node,\n");
2294   printf(
2295 "  produce its Delaunay triangulation in points.1.node and points.1.ele,\n");
2296   printf(
2297 "  and produce its Voronoi diagram in points.1.v.node and points.1.v.edge.\n");
2298   printf(
2299 "  The .v.node file contains a list of all Voronoi vertices, and the .v.edge\n"
2300 );
2301   printf(
2302 "  file contains a list of all Voronoi edges, some of which may be infinite\n"
2303 );
2304   printf(
2305 "  rays.  (The choice of filenames makes it easy to run the set of Voronoi\n");
2306   printf("  vertices through Triangle, if so desired.)\n\n");
2307   printf(
2308 "  This implementation does not use exact arithmetic to compute the Voronoi\n"
2309 );
2310   printf(
2311 "  vertices, and does not check whether neighboring vertices are identical.\n"
2312 );
2313   printf(
2314 "  Be forewarned that if the Delaunay triangulation is degenerate or\n");
2315   printf(
2316 "  near-degenerate, the Voronoi diagram may have duplicate points, crossing\n"
2317 );
2318   printf(
2319 "  edges, or infinite rays whose direction vector is zero.  Also, if you\n");
2320   printf(
2321 "  generate a constrained (as opposed to conforming) Delaunay triangulation,\n"
2322 );
2323   printf(
2324 "  or if the triangulation has holes, the corresponding Voronoi diagram is\n");
2325   printf("  likely to have crossing edges and unlikely to make sense.\n\n");
2326   printf("Mesh Topology:\n\n");
2327   printf(
2328 "  You may wish to know which triangles are adjacent to a certain Delaunay\n");
2329   printf(
2330 "  edge in an .edge file, which Voronoi regions are adjacent to a certain\n");
2331   printf(
2332 "  Voronoi edge in a .v.edge file, or which Voronoi regions are adjacent to\n"
2333 );
2334   printf(
2335 "  each other.  All of this information can be found by cross-referencing\n");
2336   printf(
2337 "  output files with the recollection that the Delaunay triangulation and\n");
2338   printf("  the Voronoi diagrams are planar duals.\n\n");
2339   printf(
2340 "  Specifically, edge i of an .edge file is the dual of Voronoi edge i of\n");
2341   printf(
2342 "  the corresponding .v.edge file, and is rotated 90 degrees counterclock-\n");
2343   printf(
2344 "  wise from the Voronoi edge.  Triangle j of an .ele file is the dual of\n");
2345   printf(
2346 "  vertex j of the corresponding .v.node file; and Voronoi region k is the\n");
2347   printf("  dual of point k of the corresponding .node file.\n\n");
2348   printf(
2349 "  Hence, to find the triangles adjacent to a Delaunay edge, look at the\n");
2350   printf(
2351 "  vertices of the corresponding Voronoi edge; their dual triangles are on\n");
2352   printf(
2353 "  the left and right of the Delaunay edge, respectively.  To find the\n");
2354   printf(
2355 "  Voronoi regions adjacent to a Voronoi edge, look at the endpoints of the\n"
2356 );
2357   printf(
2358 "  corresponding Delaunay edge; their dual regions are on the right and left\n"
2359 );
2360   printf(
2361 "  of the Voronoi edge, respectively.  To find which Voronoi regions are\n");
2362   printf("  adjacent to each other, just read the list of Delaunay edges.\n");
2363   printf("\n");
2364   printf("Statistics:\n");
2365   printf("\n");
2366   printf(
2367 "  After generating a mesh, Triangle prints a count of the number of points,\n"
2368 );
2369   printf(
2370 "  triangles, edges, boundary edges, and segments in the output mesh.  If\n");
2371   printf(
2372 "  you've forgotten the statistics for an existing mesh, the -rNEP switches\n"
2373 );
2374   printf(
2375 "  (or -rpNEP if you've got a .poly file for the existing mesh) will\n");
2376   printf("  regenerate these statistics without writing any output.\n\n");
2377   printf(
2378 "  The -V switch produces extended statistics, including a rough estimate\n");
2379   printf(
2380 "  of memory use and a histogram of triangle aspect ratios and angles in the\n"
2381 );
2382   printf("  mesh.\n\n");
2383   printf("Exact Arithmetic:\n\n");
2384   printf(
2385 "  Triangle uses adaptive exact arithmetic to perform what computational\n");
2386   printf(
2387 "  geometers call the `orientation' and `incircle' tests.  If the floating-\n"
2388 );
2389   printf(
2390 "  point arithmetic of your machine conforms to the IEEE 754 standard (as\n");
2391   printf(
2392 "  most workstations do), and does not use extended precision internal\n");
2393   printf(
2394 "  registers, then your output is guaranteed to be an absolutely true\n");
2395   printf("  Delaunay or conforming Delaunay triangulation, roundoff error\n");
2396   printf(
2397 "  notwithstanding.  The word `adaptive' implies that these arithmetic\n");
2398   printf(
2399 "  routines compute the result only to the precision necessary to guarantee\n"
2400 );
2401   printf(
2402 "  correctness, so they are usually nearly as fast as their approximate\n");
2403   printf(
2404 "  counterparts.  The exact tests can be disabled with the -X switch.  On\n");
2405   printf(
2406 "  most inputs, this switch will reduce the computation time by about eight\n"
2407 );
2408   printf(
2409 "  percent - it's not worth the risk.  There are rare difficult inputs\n");
2410   printf(
2411 "  (having many collinear and cocircular points), however, for which the\n");
2412   printf(
2413 "  difference could be a factor of two.  These are precisely the inputs most\n"
2414 );
2415   printf("  likely to cause errors if you use the -X switch.\n\n");
2416   printf(
2417 "  Unfortunately, these routines don't solve every numerical problem.  Exact\n"
2418 );
2419   printf(
2420 "  arithmetic is not used to compute the positions of points, because the\n");
2421   printf(
2422 "  bit complexity of point coordinates would grow without bound.  Hence,\n");
2423   printf(
2424 "  segment intersections aren't computed exactly; in very unusual cases,\n");
2425   printf(
2426 "  roundoff error in computing an intersection point might actually lead to\n"
2427 );
2428   printf(
2429 "  an inverted triangle and an invalid triangulation.  (This is one reason\n");
2430   printf(
2431 "  to compute your own intersection points in your .poly files.)  Similarly,\n"
2432 );
2433   printf(
2434 "  exact arithmetic is not used to compute the vertices of the Voronoi\n");
2435   printf("  diagram.\n\n");
2436   printf(
2437 "  Underflow and overflow can also cause difficulties; the exact arithmetic\n"
2438 );
2439   printf(
2440 "  routines do not ameliorate out-of-bounds exponents, which can arise\n");
2441   printf(
2442 "  during the orientation and incircle tests.  As a rule of thumb, you\n");
2443   printf(
2444 "  should ensure that your input values are within a range such that their\n");
2445   printf(
2446 "  third powers can be taken without underflow or overflow.  Underflow can\n");
2447   printf(
2448 "  silently prevent the tests from being performed exactly, while overflow\n");
2449   printf("  will typically cause a floating exception.\n\n");
2450   printf("Calling Triangle from Another Program:\n\n");
2451   printf("  Read the file triangle.h for details.\n\n");
2452   printf("Troubleshooting:\n\n");
2453   printf("  Please read this section before mailing me bugs.\n\n");
2454   printf("  `My output mesh has no triangles!'\n\n");
2455   printf(
2456 "    If you're using a PSLG, you've probably failed to specify a proper set\n"
2457 );
2458   printf(
2459 "    of bounding segments, or forgotten to use the -c switch.  Or you may\n");
2460   printf(
2461 "    have placed a hole badly.  To test these possibilities, try again with\n"
2462 );
2463   printf(
2464 "    the -c and -O switches.  Alternatively, all your input points may be\n");
2465   printf(
2466 "    collinear, in which case you can hardly expect to triangulate them.\n");
2467   printf("\n");
2468   printf("  `Triangle doesn't terminate, or just crashes.'\n");
2469   printf("\n");
2470   printf(
2471 "    Bad things can happen when triangles get so small that the distance\n");
2472   printf(
2473 "    between their vertices isn't much larger than the precision of your\n");
2474   printf(
2475 "    machine's arithmetic.  If you've compiled Triangle for single-precision\n"
2476 );
2477   printf(
2478 "    arithmetic, you might do better by recompiling it for double-precision.\n"
2479 );
2480   printf(
2481 "    Then again, you might just have to settle for more lenient constraints\n"
2482 );
2483   printf(
2484 "    on the minimum angle and the maximum area than you had planned.\n");
2485   printf("\n");
2486   printf(
2487 "    You can minimize precision problems by ensuring that the origin lies\n");
2488   printf(
2489 "    inside your point set, or even inside the densest part of your\n");
2490   printf(
2491 "    mesh.  On the other hand, if you're triangulating an object whose x\n");
2492   printf(
2493 "    coordinates all fall between 6247133 and 6247134, you're not leaving\n");
2494   printf("    much floating-point precision for Triangle to work with.\n\n");
2495   printf(
2496 "    Precision problems can occur covertly if the input PSLG contains two\n");
2497   printf(
2498 "    segments that meet (or intersect) at a very small angle, or if such an\n"
2499 );
2500   printf(
2501 "    angle is introduced by the -c switch, which may occur if a point lies\n");
2502   printf(
2503 "    ever-so-slightly inside the convex hull, and is connected by a PSLG\n");
2504   printf(
2505 "    segment to a point on the convex hull.  If you don't realize that a\n");
2506   printf(
2507 "    small angle is being formed, you might never discover why Triangle is\n");
2508   printf(
2509 "    crashing.  To check for this possibility, use the -S switch (with an\n");
2510   printf(
2511 "    appropriate limit on the number of Steiner points, found by trial-and-\n"
2512 );
2513   printf(
2514 "    error) to stop Triangle early, and view the output .poly file with\n");
2515   printf(
2516 "    Show Me (described below).  Look carefully for small angles between\n");
2517   printf(
2518 "    segments; zoom in closely, as such segments might look like a single\n");
2519   printf("    segment from a distance.\n\n");
2520   printf(
2521 "    If some of the input values are too large, Triangle may suffer a\n");
2522   printf(
2523 "    floating exception due to overflow when attempting to perform an\n");
2524   printf(
2525 "    orientation or incircle test.  (Read the section on exact arithmetic\n");
2526   printf(
2527 "    above.)  Again, I recommend compiling Triangle for double (rather\n");
2528   printf("    than single) precision arithmetic.\n\n");
2529   printf(
2530 "  `The numbering of the output points doesn't match the input points.'\n");
2531   printf("\n");
2532   printf(
2533 "    You may have eaten some of your input points with a hole, or by placing\n"
2534 );
2535   printf("    them outside the area enclosed by segments.\n\n");
2536   printf(
2537 "  `Triangle executes without incident, but when I look at the resulting\n");
2538   printf(
2539 "  mesh, it has overlapping triangles or other geometric inconsistencies.'\n");
2540   printf("\n");
2541   printf(
2542 "    If you select the -X switch, Triangle's divide-and-conquer Delaunay\n");
2543   printf(
2544 "    triangulation algorithm occasionally makes mistakes due to floating-\n");
2545   printf(
2546 "    point roundoff error.  Although these errors are rare, don't use the -X\n"
2547 );
2548   printf("    switch.  If you still have problems, please report the bug.\n");
2549   printf("\n");
2550   printf(
2551 "  Strange things can happen if you've taken liberties with your PSLG.  Do\n");
2552   printf(
2553 "  you have a point lying in the middle of a segment?  Triangle sometimes\n");
2554   printf(
2555 "  copes poorly with that sort of thing.  Do you want to lay out a collinear\n"
2556 );
2557   printf(
2558 "  row of evenly spaced, segment-connected points?  Have you simply defined\n"
2559 );
2560   printf(
2561 "  one long segment connecting the leftmost point to the rightmost point,\n");
2562   printf(
2563 "  and a bunch of points lying along it?  This method occasionally works,\n");
2564   printf(
2565 "  especially with horizontal and vertical lines, but often it doesn't, and\n"
2566 );
2567   printf(
2568 "  you'll have to connect each adjacent pair of points with a separate\n");
2569   printf("  segment.  If you don't like it, tough.\n\n");
2570   printf(
2571 "  Furthermore, if you have segments that intersect other than at their\n");
2572   printf(
2573 "  endpoints, try not to let the intersections fall extremely close to PSLG\n"
2574 );
2575   printf("  points or each other.\n\n");
2576   printf(
2577 "  If you have problems refining a triangulation not produced by Triangle:\n");
2578   printf(
2579 "  Are you sure the triangulation is geometrically valid?  Is it formatted\n");
2580   printf(
2581 "  correctly for Triangle?  Are the triangles all listed so the first three\n"
2582 );
2583   printf("  points are their corners in counterclockwise order?\n\n");
2584   printf("Show Me:\n\n");
2585   printf(
2586 "  Triangle comes with a separate program named `Show Me', whose primary\n");
2587   printf(
2588 "  purpose is to draw meshes on your screen or in PostScript.  Its secondary\n"
2589 );
2590   printf(
2591 "  purpose is to check the validity of your input files, and do so more\n");
2592   printf(
2593 "  thoroughly than Triangle does.  Show Me requires that you have the X\n");
2594   printf(
2595 "  Windows system.  If you didn't receive Show Me with Triangle, complain to\n"
2596 );
2597   printf("  whomever you obtained Triangle from, then send me mail.\n\n");
2598   printf("Triangle on the Web:\n\n");
2599   printf(
2600 "  To see an illustrated, updated version of these instructions, check out\n");
2601   printf("\n");
2602   printf("    http://www.cs.cmu.edu/~quake/triangle.html\n");
2603   printf("\n");
2604   printf("A Brief Plea:\n");
2605   printf("\n");
2606   printf(
2607 "  If you use Triangle, and especially if you use it to accomplish real\n");
2608   printf(
2609 "  work, I would like very much to hear from you.  A short letter or email\n");
2610   printf(
2611 "  (to jrs@cs.cmu.edu) describing how you use Triangle will mean a lot to\n");
2612   printf(
2613 "  me.  The more people I know are using this program, the more easily I can\n"
2614 );
2615   printf(
2616 "  justify spending time on improvements and on the three-dimensional\n");
2617   printf(
2618 "  successor to Triangle, which in turn will benefit you.  Also, I can put\n");
2619   printf(
2620 "  you on a list to receive email whenever a new version of Triangle is\n");
2621   printf("  available.\n\n");
2622   printf(
2623 "  If you use a mesh generated by Triangle in a publication, please include\n"
2624 );
2625   printf("  an acknowledgment as well.\n\n");
2626   printf("Research credit:\n\n");
2627   printf(
2628 "  Of course, I can take credit for only a fraction of the ideas that made\n");
2629   printf(
2630 "  this mesh generator possible.  Triangle owes its existence to the efforts\n"
2631 );
2632   printf(
2633 "  of many fine computational geometers and other researchers, including\n");
2634   printf(
2635 "  Marshall Bern, L. Paul Chew, Boris Delaunay, Rex A. Dwyer, David\n");
2636   printf(
2637 "  Eppstein, Steven Fortune, Leonidas J. Guibas, Donald E. Knuth, C. L.\n");
2638   printf(
2639 "  Lawson, Der-Tsai Lee, Ernst P. Mucke, Douglas M. Priest, Jim Ruppert,\n");
2640   printf(
2641 "  Isaac Saias, Bruce J. Schachter, Micha Sharir, Jorge Stolfi, Christopher\n"
2642 );
2643   printf(
2644 "  J. Van Wyk, David F. Watson, and Binhai Zhu.  See the comments at the\n");
2645   printf("  beginning of the source code for references.\n\n");
2646   exit(0);
2647 }
2648 
2649 #endif /* not TRILIBRARY */
2650 
2651 /*****************************************************************************/
2652 /*                                                                           */
2653 /*  internalerror()   Ask the user to send me the defective product.  Exit.  */
2654 /*                                                                           */
2655 /*****************************************************************************/
2656 
2657 void internalerror()
2658 {
2659   printf("  Please report this bug to jrs@cs.cmu.edu\n");
2660   printf("  Include the message above, your input data set, and the exact\n");
2661   printf("    command line you used to run Triangle.\n");
2662   exit(1);
2663 }
2664 
2665 /*****************************************************************************/
2666 /*                                                                           */
2667 /*  parsecommandline()   Read the command line, identify switches, and set   */
2668 /*                       up options and file names.                          */
2669 /*                                                                           */
2670 /*  The effects of this routine are felt entirely through global variables.  */
2671 /*                                                                           */
2672 /*****************************************************************************/
2673 
2674 void parsecommandline(argc, argv)
2675 int argc;
2676 char **argv;
2677 {
2678 #ifdef TRILIBRARY
2679 #define STARTINDEX 0
2680 #else /* not TRILIBRARY */
2681 #define STARTINDEX 1
2682   int increment;
2683   int meshnumber;
2684 #endif /* not TRILIBRARY */
2685   int i, j, k;
2686   char workstring[FILENAMESIZE];
2687 
2688   poly = refine = quality = vararea = fixedarea = regionattrib = convex = 0;
2689   firstnumber = 1;
2690   edgesout = voronoi = neighbors = geomview = 0;
2691   nobound = nopolywritten = nonodewritten = noelewritten = noiterationnum = 0;
2692   noholes = noexact = 0;
2693   incremental = sweepline = 0;
2694   dwyer = 1;
2695   splitseg = 0;
2696   docheck = 0;
2697   nobisect = 0;
2698   steiner = -1;
2699   order = 1;
2700   minangle = 0.0;
2701   maxarea = -1.0;
2702   quiet = verbose = 0;
2703 #ifndef TRILIBRARY
2704   innodefilename[0] = '\0';
2705 #endif /* not TRILIBRARY */
2706 
2707   for (i = STARTINDEX; i < argc; i++) {
2708 #ifndef TRILIBRARY
2709     if (argv[i][0] == '-') {
2710 #endif /* not TRILIBRARY */
2711       for (j = STARTINDEX; argv[i][j] != '\0'; j++) {
2712         if (argv[i][j] == 'p') {
2713           poly = 1;
2714         }
2715 #ifndef CDT_ONLY
2716         if (argv[i][j] == 'r') {
2717           refine = 1;
2718         }
2719         if (argv[i][j] == 'q') {
2720           quality = 1;
2721           if (((argv[i][j + 1] >= '') && (argv[i][j + 1] <= '9')) ||
2722               (argv[i][j + 1] == '.')) {
2723             k = 0;
2724             while (((argv[i][j + 1] >= '') && (argv[i][j + 1] <= '9')) ||
2725                    (argv[i][j + 1] == '.')) {
2726               j++;
2727               workstring[k] = argv[i][j];
2728               k++;
2729             }
2730             workstring[k] = '\0';
2731             minangle = (REAL) strtod(workstring, (char **) NULL);
2732           } else {
2733             minangle = 20.0;
2734           }
2735         }
2736         if (argv[i][j] == 'a') {
2737           quality = 1;
2738           if (((argv[i][j + 1] >= '') && (argv[i][j + 1] <= '9')) ||
2739               (argv[i][j + 1] == '.')) {
2740             fixedarea = 1;
2741             k = 0;
2742             while (((argv[i][j + 1] >= '') && (argv[i][j + 1] <= '9')) ||
2743                    (argv[i][j + 1] == '.')) {
2744               j++;
2745               workstring[k] = argv[i][j];
2746               k++;
2747             }
2748             workstring[k] = '\0';
2749             maxarea = (REAL) strtod(workstring, (char **) NULL);
2750             if (maxarea <= 0.0) {
2751               printf("Error:  Maximum area must be greater than zero.\n");
2752               exit(1);
2753             }
2754           } else {
2755             vararea = 1;
2756           }
2757         }
2758 #endif /* not CDT_ONLY */
2759         if (argv[i][j] == 'A') {
2760           regionattrib = 1;
2761         }
2762         if (argv[i][j] == 'c') {
2763           convex = 1;
2764         }
2765         if (argv[i][j] == 'z') {
2766           firstnumber = 0;
2767         }
2768         if (argv[i][j] == 'e') {
2769           edgesout = 1;
2770         }
2771         if (argv[i][j] == 'v') {
2772           voronoi = 1;
2773         }
2774         if (argv[i][j] == 'n') {
2775           neighbors = 1;
2776         }
2777         if (argv[i][j] == 'g') {
2778           geomview = 1;
2779         }
2780         if (argv[i][j] == 'B') {
2781           nobound = 1;
2782         }
2783         if (argv[i][j] == 'P') {
2784           nopolywritten = 1;
2785         }
2786         if (argv[i][j] == 'N') {
2787           nonodewritten = 1;
2788         }
2789         if (argv[i][j] == 'E') {
2790           noelewritten = 1;
2791         }
2792 #ifndef TRILIBRARY
2793         if (argv[i][j] == 'I') {
2794           noiterationnum = 1;
2795         }
2796 #endif /* not TRILIBRARY */
2797         if (argv[i][j] == 'O') {
2798           noholes = 1;
2799         }
2800         if (argv[i][j] == 'X') {
2801           noexact = 1;
2802         }
2803         if (argv[i][j] == 'o') {
2804           if (argv[i][j + 1] == '2') {
2805             j++;
2806             order = 2;
2807           }
2808         }
2809 #ifndef CDT_ONLY
2810         if (argv[i][j] == 'Y') {
2811           nobisect++;
2812         }
2813         if (argv[i][j] == 'S') {
2814           steiner = 0;
2815           while ((argv[i][j + 1] >= '') && (argv[i][j + 1] <= '9')) {
2816             j++;
2817             steiner = steiner * 10 + (int) (argv[i][j] - '');
2818           }
2819         }
2820 #endif /* not CDT_ONLY */
2821 #ifndef REDUCED
2822         if (argv[i][j] == 'i') {
2823           incremental = 1;
2824         }
2825         if (argv[i][j] == 'F') {
2826           sweepline = 1;
2827         }
2828 #endif /* not REDUCED */
2829         if (argv[i][j] == 'l') {
2830           dwyer = 0;
2831         }
2832 #ifndef REDUCED
2833 #ifndef CDT_ONLY
2834         if (argv[i][j] == 's') {
2835           splitseg = 1;
2836         }
2837 #endif /* not CDT_ONLY */
2838         if (argv[i][j] == 'C') {
2839           docheck = 1;
2840         }
2841 #endif /* not REDUCED */
2842         if (argv[i][j] == 'Q') {
2843           quiet = 1;
2844         }
2845         if (argv[i][j] == 'V') {
2846           verbose++;
2847         }
2848 #ifndef TRILIBRARY
2849         if ((argv[i][j] == 'h') || (argv[i][j] == 'H') ||
2850             (argv[i][j] == '?')) {
2851           info();
2852         }
2853 #endif /* not TRILIBRARY */
2854       }
2855 #ifndef TRILIBRARY
2856     } else {
2857       strncpy(innodefilename, argv[i], FILENAMESIZE - 1);
2858       innodefilename[FILENAMESIZE - 1] = '\0';
2859     }
2860 #endif /* not TRILIBRARY */
2861   }
2862 #ifndef TRILIBRARY
2863   if (innodefilename[0] == '\0') {
2864     syntax();
2865   }
2866   if (!strcmp(&innodefilename[strlen(innodefilename) - 5], ".node")) {
2867     innodefilename[strlen(innodefilename) - 5] = '\0';
2868   }
2869   if (!strcmp(&innodefilename[strlen(innodefilename) - 5], ".poly")) {
2870     innodefilename[strlen(innodefilename) - 5] = '\0';
2871     poly = 1;
2872   }
2873 #ifndef CDT_ONLY
2874   if (!strcmp(&innodefilename[strlen(innodefilename) - 4], ".ele")) {
2875     innodefilename[strlen(innodefilename) - 4] = '\0';
2876     refine = 1;
2877   }
2878   if (!strcmp(&innodefilename[strlen(innodefilename) - 5], ".area")) {
2879     innodefilename[strlen(innodefilename) - 5] = '\0';
2880     refine = 1;
2881     quality = 1;
2882     vararea = 1;
2883   }
2884 #endif /* not CDT_ONLY */
2885 #endif /* not TRILIBRARY */
2886   steinerleft = steiner;
2887   useshelles = poly || refine || quality || convex;
2888   goodangle = cos(minangle * PI / 180.0);
2889   goodangle *= goodangle;
2890   if (refine && noiterationnum) {
2891     printf(
2892       "Error:  You cannot use the -I switch when refining a triangulation.\n");
2893     exit(1);
2894   }
2895   /* Be careful not to allocate space for element area constraints that */
2896   /*   will never be assigned any value (other than the default -1.0).  */
2897   if (!refine && !poly) {
2898     vararea = 0;
2899   }
2900   /* Be careful not to add an extra attribute to each element unless the */
2901   /*   input supports it (PSLG in, but not refining a preexisting mesh). */
2902   if (refine || !poly) {
2903     regionattrib = 0;
2904   }
2905 
2906 #ifndef TRILIBRARY
2907   strcpy(inpolyfilename, innodefilename);
2908   strcpy(inelefilename, innodefilename);
2909   strcpy(areafilename, innodefilename);
2910   increment = 0;
2911   strcpy(workstring, innodefilename);
2912   j = 1;
2913   while (workstring[j] != '\0') {
2914     if ((workstring[j] == '.') && (workstring[j + 1] != '\0')) {
2915       increment = j + 1;
2916     }
2917     j++;
2918   }
2919   meshnumber = 0;
2920   if (increment > 0) {
2921     j = increment;
2922     do {
2923       if ((workstring[j] >= '') && (workstring[j] <= '9')) {
2924         meshnumber = meshnumber * 10 + (int) (workstring[j] - '');
2925       } else {
2926         increment = 0;
2927       }
2928       j++;
2929     } while (workstring[j] != '\0');
2930   }
2931   if (noiterationnum) {
2932     strcpy(outnodefilename, innodefilename);
2933     strcpy(outelefilename, innodefilename);
2934     strcpy(edgefilename, innodefilename);
2935     strcpy(vnodefilename, innodefilename);
2936     strcpy(vedgefilename, innodefilename);
2937     strcpy(neighborfilename, innodefilename);
2938     strcpy(offfilename, innodefilename);
2939     strcat(outnodefilename, ".node");
2940     strcat(outelefilename, ".ele");
2941     strcat(edgefilename, ".edge");
2942     strcat(vnodefilename, ".v.node");
2943     strcat(vedgefilename, ".v.edge");
2944     strcat(neighborfilename, ".neigh");
2945     strcat(offfilename, ".off");
2946   } else if (increment == 0) {
2947     strcpy(outnodefilename, innodefilename);
2948     strcpy(outpolyfilename, innodefilename);
2949     strcpy(outelefilename, innodefilename);
2950     strcpy(edgefilename, innodefilename);
2951     strcpy(vnodefilename, innodefilename);
2952     strcpy(vedgefilename, innodefilename);
2953     strcpy(neighborfilename, innodefilename);
2954     strcpy(offfilename, innodefilename);
2955     strcat(outnodefilename, ".1.node");
2956     strcat(outpolyfilename, ".1.poly");
2957     strcat(outelefilename, ".1.ele");
2958     strcat(edgefilename, ".1.edge");
2959     strcat(vnodefilename, ".1.v.node");
2960     strcat(vedgefilename, ".1.v.edge");
2961     strcat(neighborfilename, ".1.neigh");
2962     strcat(offfilename, ".1.off");
2963   } else {
2964     workstring[increment] = '%';
2965     workstring[increment + 1] = 'd';
2966     workstring[increment + 2] = '\0';
2967     sprintf(outnodefilename, workstring, meshnumber + 1);
2968     strcpy(outpolyfilename, outnodefilename);
2969     strcpy(outelefilename, outnodefilename);
2970     strcpy(edgefilename, outnodefilename);
2971     strcpy(vnodefilename, outnodefilename);
2972     strcpy(vedgefilename, outnodefilename);
2973     strcpy(neighborfilename, outnodefilename);
2974     strcpy(offfilename, outnodefilename);
2975     strcat(outnodefilename, ".node");
2976     strcat(outpolyfilename, ".poly");
2977     strcat(outelefilename, ".ele");
2978     strcat(edgefilename, ".edge");
2979     strcat(vnodefilename, ".v.node");
2980     strcat(vedgefilename, ".v.edge");
2981     strcat(neighborfilename, ".neigh");
2982     strcat(offfilename, ".off");
2983   }
2984   strcat(innodefilename, ".node");
2985   strcat(inpolyfilename, ".poly");
2986   strcat(inelefilename, ".ele");
2987   strcat(areafilename, ".area");
2988 #endif /* not TRILIBRARY */
2989 }
2990 
2991 /**                                                                         **/
2992 /**                                                                         **/
2993 /********* User interaction routines begin here                      *********/
2994 
2995 /********* Debugging routines begin here                             *********/
2996 /**                                                                         **/
2997 /**                                                                         **/
2998 
2999 /*****************************************************************************/
3000 /*                                                                           */
3001 /*  printtriangle()   Print out the details of a triangle/edge handle.       */
3002 /*                                                                           */
3003 /*  I originally wrote this procedure to simplify debugging; it can be       */
3004 /*  called directly from the debugger, and presents information about a      */
3005 /*  triangle/edge handle in digestible form.  It's also used when the        */
3006 /*  highest level of verbosity (`-VVV') is specified.                        */
3007 /*                                                                           */
3008 /*****************************************************************************/
3009 
3010 void printtriangle(t)
3011 struct triedge *t;
3012 {
3013   struct triedge printtri;
3014   struct edge printsh;
3015   point printpoint;
3016 
3017   printf("triangle x%lx with orientation %d:\n", (unsigned long) t->tri,
3018          t->orient);
3019   decode(t->tri[0], printtri);
3020   if (printtri.tri == dummytri) {
3021     printf("    [0] = Outer space\n");
3022   } else {
3023     printf("    [0] = x%lx  %d\n", (unsigned long) printtri.tri,
3024            printtri.orient);
3025   }
3026   decode(t->tri[1], printtri);
3027   if (printtri.tri == dummytri) {
3028     printf("    [1] = Outer space\n");
3029   } else {
3030     printf("    [1] = x%lx  %d\n", (unsigned long) printtri.tri,
3031            printtri.orient);
3032   }
3033   decode(t->tri[2], printtri);
3034   if (printtri.tri == dummytri) {
3035     printf("    [2] = Outer space\n");
3036   } else {
3037     printf("    [2] = x%lx  %d\n", (unsigned long) printtri.tri,
3038            printtri.orient);
3039   }
3040   org(*t, printpoint);
3041   if (printpoint == (point) NULL)
3042     printf("    Origin[%d] = NULL\n", (t->orient + 1) % 3 + 3);
3043   else
3044     printf("    Origin[%d] = x%lx  (%.12g, %.12g)\n",
3045            (t->orient + 1) % 3 + 3, (unsigned long) printpoint,
3046            printpoint[0], printpoint[1]);
3047   dest(*t, printpoint);
3048   if (printpoint == (point) NULL)
3049     printf("    Dest  [%d] = NULL\n", (t->orient + 2) % 3 + 3);
3050   else
3051     printf("    Dest  [%d] = x%lx  (%.12g, %.12g)\n",
3052            (t->orient + 2) % 3 + 3, (unsigned long) printpoint,
3053            printpoint[0], printpoint[1]);
3054   apex(*t, printpoint);
3055   if (printpoint == (point) NULL)
3056     printf("    Apex  [%d] = NULL\n", t->orient + 3);
3057   else
3058     printf("    Apex  [%d] = x%lx  (%.12g, %.12g)\n",
3059            t->orient + 3, (unsigned long) printpoint,
3060            printpoint[0], printpoint[1]);
3061   if (useshelles) {
3062     sdecode(t->tri[6], printsh);
3063     if (printsh.sh != dummysh) {
3064       printf("    [6] = x%lx  %d\n", (unsigned long) printsh.sh,
3065              printsh.shorient);
3066     }
3067     sdecode(t->tri[7], printsh);
3068     if (printsh.sh != dummysh) {
3069       printf("    [7] = x%lx  %d\n", (unsigned long) printsh.sh,
3070              printsh.shorient);
3071     }
3072     sdecode(t->tri[8], printsh);
3073     if (printsh.sh != dummysh) {
3074       printf("    [8] = x%lx  %d\n", (unsigned long) printsh.sh,
3075              printsh.shorient);
3076     }
3077   }
3078   if (vararea) {
3079     printf("    Area constraint:  %.4g\n", areabound(*t));
3080   }
3081 }
3082 
3083 /*****************************************************************************/
3084 /*                                                                           */
3085 /*  printshelle()   Print out the details of a shell edge handle.            */
3086 /*                                                                           */
3087 /*  I originally wrote this procedure to simplify debugging; it can be       */
3088 /*  called directly from the debugger, and presents information about a      */
3089 /*  shell edge handle in digestible form.  It's also used when the highest   */
3090 /*  level of verbosity (`-VVV') is specified.                                */
3091 /*                                                                           */
3092 /*****************************************************************************/
3093 
3094 void printshelle(s)
3095 struct edge *s;
3096 {
3097   struct edge printsh;
3098   struct triedge printtri;
3099   point printpoint;
3100 
3101   printf("shell edge x%lx with orientation %d and mark %d:\n",
3102          (unsigned long) s->sh, s->shorient, mark(*s));
3103   sdecode(s->sh[0], printsh);
3104   if (printsh.sh == dummysh) {
3105     printf("    [0] = No shell\n");
3106   } else {
3107     printf("    [0] = x%lx  %d\n", (unsigned long) printsh.sh,
3108            printsh.shorient);
3109   }
3110   sdecode(s->sh[1], printsh);
3111   if (printsh.sh == dummysh) {
3112     printf("    [1] = No shell\n");
3113   } else {
3114     printf("    [1] = x%lx  %d\n", (unsigned long) printsh.sh,
3115            printsh.shorient);
3116   }
3117   sorg(*s, printpoint);
3118   if (printpoint == (point) NULL)
3119     printf("    Origin[%d] = NULL\n", 2 + s->shorient);
3120   else
3121     printf("    Origin[%d] = x%lx  (%.12g, %.12g)\n",
3122            2 + s->shorient, (unsigned long) printpoint,
3123            printpoint[0], printpoint[1]);
3124   sdest(*s, printpoint);
3125   if (printpoint == (point) NULL)
3126     printf("    Dest  [%d] = NULL\n", 3 - s->shorient);
3127   else
3128     printf("    Dest  [%d] = x%lx  (%.12g, %.12g)\n",
3129            3 - s->shorient, (unsigned long) printpoint,
3130            printpoint[0], printpoint[1]);
3131   decode(s->sh[4], printtri);
3132   if (printtri.tri == dummytri) {
3133     printf("    [4] = Outer space\n");
3134   } else {
3135     printf("    [4] = x%lx  %d\n", (unsigned long) printtri.tri,
3136            printtri.orient);
3137   }
3138   decode(s->sh[5], printtri);
3139   if (printtri.tri == dummytri) {
3140     printf("    [5] = Outer space\n");
3141   } else {
3142     printf("    [5] = x%lx  %d\n", (unsigned long) printtri.tri,
3143            printtri.orient);
3144   }
3145 }
3146 
3147 /**                                                                         **/
3148 /**                                                                         **/
3149 /********* Debugging routines end here                               *********/
3150 
3151 /********* Memory management routines begin here                     *********/
3152 /**                                                                         **/
3153 /**                                                                         **/
3154 
3155 /*****************************************************************************/
3156 /*                                                                           */
3157 /*  poolinit()   Initialize a pool of memory for allocation of items.        */
3158 /*                                                                           */
3159 /*  This routine initializes the machinery for allocating items.  A `pool'   */
3160 /*  is created whose records have size at least `bytecount'.  Items will be  */
3161 /*  allocated in `itemcount'-item blocks.  Each item is assumed to be a      */
3162 /*  collection of words, and either pointers or floating-point values are    */
3163 /*  assumed to be the "primary" word type.  (The "primary" word type is used */
3164 /*  to determine alignment of items.)  If `alignment' isn't zero, all items  */
3165 /*  will be `alignment'-byte aligned in memory.  `alignment' must be either  */
3166 /*  a multiple or a factor of the primary word size; powers of two are safe. */
3167 /*  `alignment' is normally used to create a few unused bits at the bottom   */
3168 /*  of each item's pointer, in which information may be stored.              */
3169 /*                                                                           */
3170 /*  Don't change this routine unless you understand it.                      */
3171 /*                                                                           */
3172 /*****************************************************************************/
3173 
3174 void poolinit(pool, bytecount, itemcount, wtype, alignment)
3175 struct memorypool *pool;
3176 int bytecount;
3177 int itemcount;
3178 enum wordtype wtype;
3179 int alignment;
3180 {
3181   int wordsize;
3182 
3183   /* Initialize values in the pool. */
3184   pool->itemwordtype = wtype;
3185   wordsize = (pool->itemwordtype == POINTER) ? sizeof(VOID *) : sizeof(REAL);
3186   /* Find the proper alignment, which must be at least as large as:   */
3187   /*   - The parameter `alignment'.                                   */
3188   /*   - The primary word type, to avoid unaligned accesses.          */
3189   /*   - sizeof(VOID *), so the stack of dead items can be maintained */
3190   /*       without unaligned accesses.                                */
3191   if (alignment > wordsize) {
3192     pool->alignbytes = alignment;
3193   } else {
3194     pool->alignbytes = wordsize;
3195   }
3196   if (sizeof(VOID *) > pool->alignbytes) {
3197     pool->alignbytes = sizeof(VOID *);
3198   }
3199   pool->itemwords = ((bytecount + pool->alignbytes - 1) / pool->alignbytes)
3200                   * (pool->alignbytes / wordsize);
3201   pool->itembytes = pool->itemwords * wordsize;
3202   pool->itemsperblock = itemcount;
3203 
3204   /* Allocate a block of items.  Space for `itemsperblock' items and one    */
3205   /*   pointer (to point to the next block) are allocated, as well as space */
3206   /*   to ensure alignment of the items.                                    */
3207   pool->firstblock = (VOID **) malloc(pool->itemsperblock * pool->itembytes
3208                                       + sizeof(VOID *) + pool->alignbytes);
3209   if (pool->firstblock == (VOID **) NULL) {
3210     printf("Error:  Out of memory.\n");
3211     exit(1);
3212   }
3213   /* Set the next block pointer to NULL. */
3214   *(pool->firstblock) = (VOID *) NULL;
3215   poolrestart(pool);
3216 }
3217 
3218 /*****************************************************************************/
3219 /*                                                                           */
3220 /*  poolrestart()   Deallocate all items in a pool.                          */
3221 /*                                                                           */
3222 /*  The pool is returned to its starting state, except that no memory is     */
3223 /*  freed to the operating system.  Rather, the previously allocated blocks  */
3224 /*  are ready to be reused.                                                  */
3225 /*                                                                           */
3226 /*****************************************************************************/
3227 
3228 void poolrestart(pool)
3229 struct memorypool *pool;
3230 {
3231   unsigned long alignptr;
3232 
3233   pool->items = 0;
3234   pool->maxitems = 0;
3235 
3236   /* Set the currently active block. */
3237   pool->nowblock = pool->firstblock;
3238   /* Find the first item in the pool.  Increment by the size of (VOID *). */
3239   alignptr = (unsigned long) (pool->nowblock + 1);
3240   /* Align the item on an `alignbytes'-byte boundary. */
3241   pool->nextitem = (VOID *)
3242     (alignptr + (unsigned long) pool->alignbytes
3243      - (alignptr % (unsigned long) pool->alignbytes));
3244   /* There are lots of unallocated items left in this block. */
3245   pool->unallocateditems = pool->itemsperblock;
3246   /* The stack of deallocated items is empty. */
3247   pool->deaditemstack = (VOID *) NULL;
3248 }
3249 
3250 /*****************************************************************************/
3251 /*                                                                           */
3252 /*  pooldeinit()   Free to the operating system all memory taken by a pool.  */
3253 /*                                                                           */
3254 /*****************************************************************************/
3255 
3256 void pooldeinit(pool)
3257 struct memorypool *pool;
3258 {
3259   while (pool->firstblock != (VOID **) NULL) {
3260     pool->nowblock = (VOID **) *(pool->firstblock);
3261     free(pool->firstblock);
3262     pool->firstblock = pool->nowblock;
3263   }
3264 }
3265 
3266 /*****************************************************************************/
3267 /*                                                                           */
3268 /*  poolalloc()   Allocate space for an item.                                */
3269 /*                                                                           */
3270 /*****************************************************************************/
3271 
3272 VOID *poolalloc(pool)
3273 struct memorypool *pool;
3274 {
3275   VOID *newitem;
3276   VOID **newblock;
3277   unsigned long alignptr;
3278 
3279   /* First check the linked list of dead items.  If the list is not   */
3280   /*   empty, allocate an item from the list rather than a fresh one. */
3281   if (pool->deaditemstack != (VOID *) NULL) {
3282     newitem = pool->deaditemstack;               /* Take first item in list. */
3283     pool->deaditemstack = * (VOID **) pool->deaditemstack;
3284   } else {
3285     /* Check if there are any free items left in the current block. */
3286     if (pool->unallocateditems == 0) {
3287       /* Check if another block must be allocated. */
3288       if (*(pool->nowblock) == (VOID *) NULL) {
3289         /* Allocate a new block of items, pointed to by the previous block. */
3290         newblock = (VOID **) malloc(pool->itemsperblock * pool->itembytes
3291                                     + sizeof(VOID *) + pool->alignbytes);
3292         if (newblock == (VOID **) NULL) {
3293           printf("Error:  Out of memory.\n");
3294           exit(1);
3295         }
3296         *(pool->nowblock) = (VOID *) newblock;
3297         /* The next block pointer is NULL. */
3298         *newblock = (VOID *) NULL;
3299       }
3300       /* Move to the new block. */
3301       pool->nowblock = (VOID **) *(pool->nowblock);
3302       /* Find the first item in the block.    */
3303       /*   Increment by the size of (VOID *). */
3304       alignptr = (unsigned long) (pool->nowblock + 1);
3305       /* Align the item on an `alignbytes'-byte boundary. */
3306       pool->nextitem = (VOID *)
3307         (alignptr + (unsigned long) pool->alignbytes
3308          - (alignptr % (unsigned long) pool->alignbytes));
3309       /* There are lots of unallocated items left in this block. */
3310       pool->unallocateditems = pool->itemsperblock;
3311     }
3312     /* Allocate a new item. */
3313     newitem = pool->nextitem;
3314     /* Advance `nextitem' pointer to next free item in block. */
3315     if (pool->itemwordtype == POINTER) {
3316       pool->nextitem = (VOID *) ((VOID **) pool->nextitem + pool->itemwords);
3317     } else {
3318       pool->nextitem = (VOID *) ((REAL *) pool->nextitem + pool->itemwords);
3319     }
3320     pool->unallocateditems--;
3321     pool->maxitems++;
3322   }
3323   pool->items++;
3324   return newitem;
3325 }
3326 
3327 /*****************************************************************************/
3328 /*                                                                           */
3329 /*  pooldealloc()   Deallocate space for an item.                            */
3330 /*                                                                           */
3331 /*  The deallocated space is stored in a queue for later reuse.              */
3332 /*                                                                           */
3333 /*****************************************************************************/
3334 
3335 void pooldealloc(pool, dyingitem)
3336 struct memorypool *pool;
3337 VOID *dyingitem;
3338 {
3339   /* Push freshly killed item onto stack. */
3340   *((VOID **) dyingitem) = pool->deaditemstack;
3341   pool->deaditemstack = dyingitem;
3342   pool->items--;
3343 }
3344 
3345 /*****************************************************************************/
3346 /*                                                                           */
3347 /*  traversalinit()   Prepare to traverse the entire list of items.          */
3348 /*                                                                           */
3349 /*  This routine is used in conjunction with traverse().                     */
3350 /*                                                                           */
3351 /*****************************************************************************/
3352 
3353 void traversalinit(pool)
3354 struct memorypool *pool;
3355 {
3356   unsigned long alignptr;
3357 
3358   /* Begin the traversal in the first block. */
3359   pool->pathblock = pool->firstblock;
3360   /* Find the first item in the block.  Increment by the size of (VOID *). */
3361   alignptr = (unsigned long) (pool->pathblock + 1);
3362   /* Align with item on an `alignbytes'-byte boundary. */
3363   pool->pathitem = (VOID *)
3364     (alignptr + (unsigned long) pool->alignbytes
3365      - (alignptr % (unsigned long) pool->alignbytes));
3366   /* Set the number of items left in the current block. */
3367   pool->pathitemsleft = pool->itemsperblock;
3368 }
3369 
3370 /*****************************************************************************/
3371 /*                                                                           */
3372 /*  traverse()   Find the next item in the list.                             */
3373 /*                                                                           */
3374 /*  This routine is used in conjunction with traversalinit().  Be forewarned */
3375 /*  that this routine successively returns all items in the list, including  */
3376 /*  deallocated ones on the deaditemqueue.  It's up to you to figure out     */
3377 /*  which ones are actually dead.  Why?  I don't want to allocate extra      */
3378 /*  space just to demarcate dead items.  It can usually be done more         */
3379 /*  space-efficiently by a routine that knows something about the structure  */
3380 /*  of the item.                                                             */
3381 /*                                                                           */
3382 /*****************************************************************************/
3383 
3384 VOID *traverse(pool)
3385 struct memorypool *pool;
3386 {
3387   VOID *newitem;
3388   unsigned long alignptr;
3389 
3390   /* Stop upon exhausting the list of items. */
3391   if (pool->pathitem == pool->nextitem) {
3392     return (VOID *) NULL;
3393   }
3394   /* Check whether any untraversed items remain in the current block. */
3395   if (pool->pathitemsleft == 0) {
3396     /* Find the next block. */
3397     pool->pathblock = (VOID **) *(pool->pathblock);
3398     /* Find the first item in the block.  Increment by the size of (VOID *). */
3399     alignptr = (unsigned long) (pool->pathblock + 1);
3400     /* Align with item on an `alignbytes'-byte boundary. */
3401     pool->pathitem = (VOID *)
3402       (alignptr + (unsigned long) pool->alignbytes
3403        - (alignptr % (unsigned long) pool->alignbytes));
3404     /* Set the number of items left in the current block. */
3405     pool->pathitemsleft = pool->itemsperblock;
3406   }
3407   newitem = pool->pathitem;
3408   /* Find the next item in the block. */
3409   if (pool->itemwordtype == POINTER) {
3410     pool->pathitem = (VOID *) ((VOID **) pool->pathitem + pool->itemwords);
3411   } else {
3412     pool->pathitem = (VOID *) ((REAL *) pool->pathitem + pool->itemwords);
3413   }
3414   pool->pathitemsleft--;
3415   return newitem;
3416 }
3417 
3418 /*****************************************************************************/
3419 /*                                                                           */
3420 /*  dummyinit()   Initialize the triangle that fills "outer space" and the   */
3421 /*                omnipresent shell edge.                                    */
3422 /*                                                                           */
3423 /*  The triangle that fills "outer space", called `dummytri', is pointed to  */
3424 /*  by every triangle and shell edge on a boundary (be it outer or inner) of */
3425 /*  the triangulation.  Also, `dummytri' points to one of the triangles on   */
3426 /*  the convex hull (until the holes and concavities are carved), making it  */
3427 /*  possible to find a starting triangle for point location.                 */
3428 /*                                                                           */
3429 /*  The omnipresent shell edge, `dummysh', is pointed to by every triangle   */
3430 /*  or shell edge that doesn't have a full complement of real shell edges    */
3431 /*  to point to.                                                             */
3432 /*                                                                           */
3433 /*****************************************************************************/
3434 
3435 void dummyinit(trianglewords, shellewords)
3436 int trianglewords;
3437 int shellewords;
3438 {
3439   unsigned long alignptr;
3440 
3441   /* `triwords' and `shwords' are used by the mesh manipulation primitives */
3442   /*   to extract orientations of triangles and shell edges from pointers. */
3443   triwords = trianglewords;       /* Initialize `triwords' once and for all. */
3444   shwords = shellewords;           /* Initialize `shwords' once and for all. */
3445 
3446   /* Set up `dummytri', the `triangle' that occupies "outer space". */
3447   dummytribase = (triangle *) malloc(triwords * sizeof(triangle)
3448                                      + triangles.alignbytes);
3449   if (dummytribase == (triangle *) NULL) {
3450     printf("Error:  Out of memory.\n");
3451     exit(1);
3452   }
3453   /* Align `dummytri' on a `triangles.alignbytes'-byte boundary. */
3454   alignptr = (unsigned long) dummytribase;
3455   dummytri = (triangle *)
3456     (alignptr + (unsigned long) triangles.alignbytes
3457      - (alignptr % (unsigned long) triangles.alignbytes));
3458   /* Initialize the three adjoining triangles to be "outer space".  These  */
3459   /*   will eventually be changed by various bonding operations, but their */
3460   /*   values don't really matter, as long as they can legally be          */
3461   /*   dereferenced.                                                       */
3462   dummytri[0] = (triangle) dummytri;
3463   dummytri[1] = (triangle) dummytri;
3464   dummytri[2] = (triangle) dummytri;
3465   /* Three NULL vertex points. */
3466   dummytri[3] = (triangle) NULL;
3467   dummytri[4] = (triangle) NULL;
3468   dummytri[5] = (triangle) NULL;
3469 
3470   if (useshelles) {
3471     /* Set up `dummysh', the omnipresent "shell edge" pointed to by any      */
3472     /*   triangle side or shell edge end that isn't attached to a real shell */
3473     /*   edge.                                                               */
3474     dummyshbase = (shelle *) malloc(shwords * sizeof(shelle)
3475                                     + shelles.alignbytes);
3476     if (dummyshbase == (shelle *) NULL) {
3477       printf("Error:  Out of memory.\n");
3478       exit(1);
3479     }
3480     /* Align `dummysh' on a `shelles.alignbytes'-byte boundary. */
3481     alignptr = (unsigned long) dummyshbase;
3482     dummysh = (shelle *)
3483       (alignptr + (unsigned long) shelles.alignbytes
3484        - (alignptr % (unsigned long) shelles.alignbytes));
3485     /* Initialize the two adjoining shell edges to be the omnipresent shell */
3486     /*   edge.  These will eventually be changed by various bonding         */
3487     /*   operations, but their values don't really matter, as long as they  */
3488     /*   can legally be dereferenced.                                       */
3489     dummysh[0] = (shelle) dummysh;
3490     dummysh[1] = (shelle) dummysh;
3491     /* Two NULL vertex points. */
3492     dummysh[2] = (shelle) NULL;
3493     dummysh[3] = (shelle) NULL;
3494     /* Initialize the two adjoining triangles to be "outer space". */
3495     dummysh[4] = (shelle) dummytri;
3496     dummysh[5] = (shelle) dummytri;
3497     /* Set the boundary marker to zero. */
3498     * (int *) (dummysh + 6) = 0;
3499 
3500     /* Initialize the three adjoining shell edges of `dummytri' to be */
3501     /*   the omnipresent shell edge.                                  */
3502     dummytri[6] = (triangle) dummysh;
3503     dummytri[7] = (triangle) dummysh;
3504     dummytri[8] = (triangle) dummysh;
3505   }
3506 }
3507 
3508 /*****************************************************************************/
3509 /*                                                                           */
3510 /*  initializepointpool()   Calculate the size of the point data structure   */
3511 /*                          and initialize its memory pool.                  */
3512 /*                                                                           */
3513 /*  This routine also computes the `pointmarkindex' and `point2triindex'     */
3514 /*  indices used to find values within each point.                           */
3515 /*                                                                           */
3516 /*****************************************************************************/
3517 
3518 void initializepointpool()
3519 {
3520   int pointsize;
3521 
3522   /* The index within each point at which the boundary marker is found.  */
3523   /*   Ensure the point marker is aligned to a sizeof(int)-byte address. */
3524   pointmarkindex = ((mesh_dim + nextras) * sizeof(REAL) + sizeof(int) - 1)
3525                  / sizeof(int);
3526   pointsize = (pointmarkindex + 1) * sizeof(int);
3527   if (poly) {
3528     /* The index within each point at which a triangle pointer is found.   */
3529     /*   Ensure the pointer is aligned to a sizeof(triangle)-byte address. */
3530     point2triindex = (pointsize + sizeof(triangle) - 1) / sizeof(triangle);
3531     pointsize = (point2triindex + 1) * sizeof(triangle);
3532   }
3533   /* Initialize the pool of points. */
3534   poolinit(&points, pointsize, POINTPERBLOCK,
3535            (sizeof(REAL) >= sizeof(triangle)) ? FLOATINGPOINT : POINTER, 0);
3536 }
3537 
3538 /*****************************************************************************/
3539 /*                                                                           */
3540 /*  initializetrisegpools()   Calculate the sizes of the triangle and shell  */
3541 /*                            edge data structures and initialize their      */
3542 /*                            memory pools.                                  */
3543 /*                                                                           */
3544 /*  This routine also computes the `highorderindex', `elemattribindex', and  */
3545 /*  `areaboundindex' indices used to find values within each triangle.       */
3546 /*                                                                           */
3547 /*****************************************************************************/
3548 
3549 void initializetrisegpools()
3550 {
3551   int trisize;
3552 
3553   /* The index within each triangle at which the extra nodes (above three)  */
3554   /*   associated with high order elements are found.  There are three      */
3555   /*   pointers to other triangles, three pointers to corners, and possibly */
3556   /*   three pointers to shell edges before the extra nodes.                */
3557   highorderindex = 6 + (useshelles * 3);
3558   /* The number of bytes occupied by a triangle. */
3559   trisize = ((order + 1) * (order + 2) / 2 + (highorderindex - 3)) *
3560             sizeof(triangle);
3561   /* The index within each triangle at which its attributes are found, */
3562   /*   where the index is measured in REALs.                           */
3563   elemattribindex = (trisize + sizeof(REAL) - 1) / sizeof(REAL);
3564   /* The index within each triangle at which the maximum area constraint  */
3565   /*   is found, where the index is measured in REALs.  Note that if the  */
3566   /*   `regionattrib' flag is set, an additional attribute will be added. */
3567   areaboundindex = elemattribindex + eextras + regionattrib;
3568   /* If triangle attributes or an area bound are needed, increase the number */
3569   /*   of bytes occupied by a triangle.                                      */
3570   if (vararea) {
3571     trisize = (areaboundindex + 1) * sizeof(REAL);
3572   } else if (eextras + regionattrib > 0) {
3573     trisize = areaboundindex * sizeof(REAL);
3574   }
3575   /* If a Voronoi diagram or triangle neighbor graph is requested, make    */
3576   /*   sure there's room to store an integer index in each triangle.  This */
3577   /*   integer index can occupy the same space as the shell edges or       */
3578   /*   attributes or area constraint or extra nodes.                       */
3579   if ((voronoi || neighbors) &&
3580       (trisize < 6 * sizeof(triangle) + sizeof(int))) {
3581     trisize = 6 * sizeof(triangle) + sizeof(int);
3582   }
3583   /* Having determined the memory size of a triangle, initialize the pool. */
3584   poolinit(&triangles, trisize, TRIPERBLOCK, POINTER, 4);
3585 
3586   if (useshelles) {
3587     /* Initialize the pool of shell edges. */
3588     poolinit(&shelles, 6 * sizeof(triangle) + sizeof(int), SHELLEPERBLOCK,
3589              POINTER, 4);
3590 
3591     /* Initialize the "outer space" triangle and omnipresent shell edge. */
3592     dummyinit(triangles.itemwords, shelles.itemwords);
3593   } else {
3594     /* Initialize the "outer space" triangle. */
3595     dummyinit(triangles.itemwords, 0);
3596   }
3597 }
3598 
3599 /*****************************************************************************/
3600 /*                                                                           */
3601 /*  triangledealloc()   Deallocate space for a triangle, marking it dead.    */
3602 /*                                                                           */
3603 /*****************************************************************************/
3604 
3605 void triangledealloc(dyingtriangle)
3606 triangle *dyingtriangle;
3607 {
3608   /* Set triangle's vertices to NULL.  This makes it possible to        */
3609   /*   detect dead triangles when traversing the list of all triangles. */
3610   dyingtriangle[3] = (triangle) NULL;
3611   dyingtriangle[4] = (triangle) NULL;
3612   dyingtriangle[5] = (triangle) NULL;
3613   pooldealloc(&triangles, (VOID *) dyingtriangle);
3614 }
3615 
3616 /*****************************************************************************/
3617 /*                                                                           */
3618 /*  triangletraverse()   Traverse the triangles, skipping dead ones.         */
3619 /*                                                                           */
3620 /*****************************************************************************/
3621 
3622 triangle *triangletraverse()
3623 {
3624   triangle *newtriangle;
3625 
3626   do {
3627     newtriangle = (triangle *) traverse(&triangles);
3628     if (newtriangle == (triangle *) NULL) {
3629       return (triangle *) NULL;
3630     }
3631   } while (newtriangle[3] == (triangle) NULL);            /* Skip dead ones. */
3632   return newtriangle;
3633 }
3634 
3635 /*****************************************************************************/
3636 /*                                                                           */
3637 /*  shelledealloc()   Deallocate space for a shell edge, marking it dead.    */
3638 /*                                                                           */
3639 /*****************************************************************************/
3640 
3641 void shelledealloc(dyingshelle)
3642 shelle *dyingshelle;
3643 {
3644   /* Set shell edge's vertices to NULL.  This makes it possible to */
3645   /*   detect dead shells when traversing the list of all shells.  */
3646   dyingshelle[2] = (shelle) NULL;
3647   dyingshelle[3] = (shelle) NULL;
3648   pooldealloc(&shelles, (VOID *) dyingshelle);
3649 }
3650 
3651 /*****************************************************************************/
3652 /*                                                                           */
3653 /*  shelletraverse()   Traverse the shell edges, skipping dead ones.         */
3654 /*                                                                           */
3655 /*****************************************************************************/
3656 
3657 shelle *shelletraverse()
3658 {
3659   shelle *newshelle;
3660 
3661   do {
3662     newshelle = (shelle *) traverse(&shelles);
3663     if (newshelle == (shelle *) NULL) {
3664       return (shelle *) NULL;
3665     }
3666   } while (newshelle[2] == (shelle) NULL);                /* Skip dead ones. */
3667   return newshelle;
3668 }
3669 
3670 /*****************************************************************************/
3671 /*                                                                           */
3672 /*  pointdealloc()   Deallocate space for a point, marking it dead.          */
3673 /*                                                                           */
3674 /*****************************************************************************/
3675 
3676 void pointdealloc(dyingpoint)
3677 point dyingpoint;
3678 {
3679   /* Mark the point as dead.  This makes it possible to detect dead points */
3680   /*   when traversing the list of all points.                             */
3681   setpointmark(dyingpoint, DEADPOINT);
3682   pooldealloc(&points, (VOID *) dyingpoint);
3683 }
3684 
3685 /*****************************************************************************/
3686 /*                                                                           */
3687 /*  pointtraverse()   Traverse the points, skipping dead ones.               */
3688 /*                                                                           */
3689 /*****************************************************************************/
3690 
3691 point pointtraverse()
3692 {
3693   point newpoint;
3694 
3695   do {
3696     newpoint = (point) traverse(&points);
3697     if (newpoint == (point) NULL) {
3698       return (point) NULL;
3699     }
3700   } while (pointmark(newpoint) == DEADPOINT);             /* Skip dead ones. */
3701   return newpoint;
3702 }
3703 
3704 /*****************************************************************************/
3705 /*                                                                           */
3706 /*  badsegmentdealloc()   Deallocate space for a bad segment, marking it     */
3707 /*                        dead.                                              */
3708 /*                                                                           */
3709 /*****************************************************************************/
3710 
3711 #ifndef CDT_ONLY
3712 
3713 void badsegmentdealloc(dyingseg)
3714 struct edge *dyingseg;
3715 {
3716   /* Set segment's orientation to -1.  This makes it possible to      */
3717   /*   detect dead segments when traversing the list of all segments. */
3718   dyingseg->shorient = -1;
3719   pooldealloc(&badsegments, (VOID *) dyingseg);
3720 }
3721 
3722 #endif /* not CDT_ONLY */
3723 
3724 /*****************************************************************************/
3725 /*                                                                           */
3726 /*  badsegmenttraverse()   Traverse the bad segments, skipping dead ones.    */
3727 /*                                                                           */
3728 /*****************************************************************************/
3729 
3730 #ifndef CDT_ONLY
3731 
3732 struct edge *badsegmenttraverse()
3733 {
3734   struct edge *newseg;
3735 
3736   do {
3737     newseg = (struct edge *) traverse(&badsegments);
3738     if (newseg == (struct edge *) NULL) {
3739       return (struct edge *) NULL;
3740     }
3741   } while (newseg->shorient == -1);                       /* Skip dead ones. */
3742   return newseg;
3743 }
3744 
3745 #endif /* not CDT_ONLY */
3746 
3747 /*****************************************************************************/
3748 /*                                                                           */
3749 /*  getpoint()   Get a specific point, by number, from the list.             */
3750 /*                                                                           */
3751 /*  The first point is number 'firstnumber'.                                 */
3752 /*                                                                           */
3753 /*  Note that this takes O(n) time (with a small constant, if POINTPERBLOCK  */
3754 /*  is large).  I don't care to take the trouble to make it work in constant */
3755 /*  time.                                                                    */
3756 /*                                                                           */
3757 /*****************************************************************************/
3758 
3759 point getpoint(number)
3760 int number;
3761 {
3762   VOID **getblock;
3763   point foundpoint;
3764   unsigned long alignptr;
3765   int current;
3766 
3767   getblock = points.firstblock;
3768   current = firstnumber;
3769   /* Find the right block. */
3770   while (current + points.itemsperblock <= number) {
3771     getblock = (VOID **) *getblock;
3772     current += points.itemsperblock;
3773   }
3774   /* Now find the right point. */
3775   alignptr = (unsigned long) (getblock + 1);
3776   foundpoint = (point) (alignptr + (unsigned long) points.alignbytes
3777                         - (alignptr % (unsigned long) points.alignbytes));
3778   while (current < number) {
3779     foundpoint += points.itemwords;
3780     current++;
3781   }
3782   return foundpoint;
3783 }
3784 
3785 /*****************************************************************************/
3786 /*                                                                           */
3787 /*  triangledeinit()   Free all remaining allocated memory.                  */
3788 /*                                                                           */
3789 /*****************************************************************************/
3790 
3791 void triangledeinit()
3792 {
3793   pooldeinit(&triangles);
3794   free(dummytribase);
3795   if (useshelles) {
3796     pooldeinit(&shelles);
3797     free(dummyshbase);
3798   }
3799   pooldeinit(&points);
3800 #ifndef CDT_ONLY
3801   if (quality) {
3802     pooldeinit(&badsegments);
3803     if ((minangle > 0.0) || vararea || fixedarea) {
3804       pooldeinit(&badtriangles);
3805     }
3806   }
3807 #endif /* not CDT_ONLY */
3808 }
3809 
3810 /**                                                                         **/
3811 /**                                                                         **/
3812 /********* Memory management routines end here                       *********/
3813 
3814 /********* Constructors begin here                                   *********/
3815 /**                                                                         **/
3816 /**                                                                         **/
3817 
3818 /*****************************************************************************/
3819 /*                                                                           */
3820 /*  maketriangle()   Create a new triangle with orientation zero.            */
3821 /*                                                                           */
3822 /*****************************************************************************/
3823 
3824 void maketriangle(newtriedge)
3825 struct triedge *newtriedge;
3826 {
3827   int i;
3828 
3829   newtriedge->tri = (triangle *) poolalloc(&triangles);
3830   /* Initialize the three adjoining triangles to be "outer space". */
3831   newtriedge->tri[0] = (triangle) dummytri;
3832   newtriedge->tri[1] = (triangle) dummytri;
3833   newtriedge->tri[2] = (triangle) dummytri;
3834   /* Three NULL vertex points. */
3835   newtriedge->tri[3] = (triangle) NULL;
3836   newtriedge->tri[4] = (triangle) NULL;
3837   newtriedge->tri[5] = (triangle) NULL;
3838   /* Initialize the three adjoining shell edges to be the omnipresent */
3839   /*   shell edge.                                                    */
3840   if (useshelles) {
3841     newtriedge->tri[6] = (triangle) dummysh;
3842     newtriedge->tri[7] = (triangle) dummysh;
3843     newtriedge->tri[8] = (triangle) dummysh;
3844   }
3845   for (i = 0; i < eextras; i++) {
3846     setelemattribute(*newtriedge, i, 0.0);
3847   }
3848   if (vararea) {
3849     setareabound(*newtriedge, -1.0);
3850   }
3851 
3852   newtriedge->orient = 0;
3853 }
3854 
3855 /*****************************************************************************/
3856 /*                                                                           */
3857 /*  makeshelle()   Create a new shell edge with orientation zero.            */
3858 /*                                                                           */
3859 /*****************************************************************************/
3860 
3861 void makeshelle(newedge)
3862 struct edge *newedge;
3863 {
3864   newedge->sh = (shelle *) poolalloc(&shelles);
3865   /* Initialize the two adjoining shell edges to be the omnipresent */
3866   /*   shell edge.                                                  */
3867   newedge->sh[0] = (shelle) dummysh;
3868   newedge->sh[1] = (shelle) dummysh;
3869   /* Two NULL vertex points. */
3870   newedge->sh[2] = (shelle) NULL;
3871   newedge->sh[3] = (shelle) NULL;
3872   /* Initialize the two adjoining triangles to be "outer space". */
3873   newedge->sh[4] = (shelle) dummytri;
3874   newedge->sh[5] = (shelle) dummytri;
3875   /* Set the boundary marker to zero. */
3876   setmark(*newedge, 0);
3877 
3878   newedge->shorient = 0;
3879 }
3880 
3881 /**                                                                         **/
3882 /**                                                                         **/
3883 /********* Constructors end here                                     *********/
3884 
3885 /********* Determinant evaluation routines begin here                *********/
3886 /**                                                                         **/
3887 /**                                                                         **/
3888 
3889 /* The adaptive exact arithmetic geometric predicates implemented herein are */
3890 /*   described in detail in my Technical Report CMU-CS-96-140.  The complete */
3891 /*   reference is given in the header.                                       */
3892 
3893 /* Which of the following two methods of finding the absolute values is      */
3894 /*   fastest is compiler-dependent.  A few compilers can inline and optimize */
3895 /*   the fabs() call; but most will incur the overhead of a function call,   */
3896 /*   which is disastrously slow.  A faster way on IEEE machines might be to  */
3897 /*   mask the appropriate bit, but that's difficult to do in C.              */
3898 
3899 #define Absolute(a)  ((a) >= 0.0 ? (a) : -(a))
3900 /* #define Absolute(a)  fabs(a) */
3901 
3902 /* Many of the operations are broken up into two pieces, a main part that    */
3903 /*   performs an approximate operation, and a "tail" that computes the       */
3904 /*   roundoff error of that operation.                                       */
3905 /*                                                                           */
3906 /* The operations Fast_Two_Sum(), Fast_Two_Diff(), Two_Sum(), Two_Diff(),    */
3907 /*   Split(), and Two_Product() are all implemented as described in the      */
3908 /*   reference.  Each of these macros requires certain variables to be       */
3909 /*   defined in the calling routine.  The variables `bvirt', `c', `abig',    */
3910 /*   `_i', `_j', `_k', `_l', `_m', and `_n' are declared `INEXACT' because   */
3911 /*   they store the result of an operation that may incur roundoff error.    */
3912 /*   The input parameter `x' (or the highest numbered `x_' parameter) must   */
3913 /*   also be declared `INEXACT'.                                             */
3914 
3915 #define Fast_Two_Sum_Tail(a, b, x, y) \
3916   bvirt = x - a; \
3917   y = b - bvirt
3918 
3919 #define Fast_Two_Sum(a, b, x, y) \
3920   x = (REAL) (a + b); \
3921   Fast_Two_Sum_Tail(a, b, x, y)
3922 
3923 #define Two_Sum_Tail(a, b, x, y) \
3924   bvirt = (REAL) (x - a); \
3925   avirt = x - bvirt; \
3926   bround = b - bvirt; \
3927   around = a - avirt; \
3928   y = around + bround
3929 
3930 #define Two_Sum(a, b, x, y) \
3931   x = (REAL) (a + b); \
3932   Two_Sum_Tail(a, b, x, y)
3933 
3934 #define Two_Diff_Tail(a, b, x, y) \
3935   bvirt = (REAL) (a - x); \
3936   avirt = x + bvirt; \
3937   bround = bvirt - b; \
3938   around = a - avirt; \
3939   y = around + bround
3940 
3941 #define Two_Diff(a, b, x, y) \
3942   x = (REAL) (a - b); \
3943   Two_Diff_Tail(a, b, x, y)
3944 
3945 #define Split(a, ahi, alo) \
3946   c = (REAL) (splitter * a); \
3947   abig = (REAL) (c - a); \
3948   ahi = c - abig; \
3949   alo = a - ahi
3950 
3951 #define Two_Product_Tail(a, b, x, y) \
3952   Split(a, ahi, alo); \
3953   Split(b, bhi, blo); \
3954   err1 = x - (ahi * bhi); \
3955   err2 = err1 - (alo * bhi); \
3956   err3 = err2 - (ahi * blo); \
3957   y = (alo * blo) - err3
3958 
3959 #define Two_Product(a, b, x, y) \
3960   x = (REAL) (a * b); \
3961   Two_Product_Tail(a, b, x, y)
3962 
3963 /* Two_Product_Presplit() is Two_Product() where one of the inputs has       */
3964 /*   already been split.  Avoids redundant splitting.                        */
3965 
3966 #define Two_Product_Presplit(a, b, bhi, blo, x, y) \
3967   x = (REAL) (a * b); \
3968   Split(a, ahi, alo); \
3969   err1 = x - (ahi * bhi); \
3970   err2 = err1 - (alo * bhi); \
3971   err3 = err2 - (ahi * blo); \
3972   y = (alo * blo) - err3
3973 
3974 /* Square() can be done more quickly than Two_Product().                     */
3975 
3976 #define Square_Tail(a, x, y) \
3977   Split(a, ahi, alo); \
3978   err1 = x - (ahi * ahi); \
3979   err3 = err1 - ((ahi + ahi) * alo); \
3980   y = (alo * alo) - err3
3981 
3982 #define Square(a, x, y) \
3983   x = (REAL) (a * a); \
3984   Square_Tail(a, x, y)
3985 
3986 /* Macros for summing expansions of various fixed lengths.  These are all    */
3987 /*   unrolled versions of Expansion_Sum().                                   */
3988 
3989 #define Two_One_Sum(a1, a0, b, x2, x1, x0) \
3990   Two_Sum(a0, b , _i, x0); \
3991   Two_Sum(a1, _i, x2, x1)
3992 
3993 #define Two_One_Diff(a1, a0, b, x2, x1, x0) \
3994   Two_Diff(a0, b , _i, x0); \
3995   Two_Sum( a1, _i, x2, x1)
3996 
3997 #define Two_Two_Sum(a1, a0, b1, b0, x3, x2, x1, x0) \
3998   Two_One_Sum(a1, a0, b0, _j, _0, x0); \
3999   Two_One_Sum(_j, _0, b1, x3, x2, x1)
4000 
4001 #define Two_Two_Diff(a1, a0, b1, b0, x3, x2, x1, x0) \
4002   Two_One_Diff(a1, a0, b0, _j, _0, x0); \
4003   Two_One_Diff(_j, _0, b1, x3, x2, x1)
4004 
4005 /*****************************************************************************/
4006 /*                                                                           */
4007 /*  exactinit()   Initialize the variables used for exact arithmetic.        */
4008 /*                                                                           */
4009 /*  `epsilon' is the largest power of two such that 1.0 + epsilon = 1.0 in   */
4010 /*  floating-point arithmetic.  `epsilon' bounds the relative roundoff       */
4011 /*  error.  It is used for floating-point error analysis.                    */
4012 /*                                                                           */
4013 /*  `splitter' is used to split floating-point numbers into two half-        */
4014 /*  length significands for exact multiplication.                            */
4015 /*                                                                           */
4016 /*  I imagine that a highly optimizing compiler might be too smart for its   */
4017 /*  own good, and somehow cause this routine to fail, if it pretends that    */
4018 /*  floating-point arithmetic is too much like real arithmetic.              */
4019 /*                                                                           */
4020 /*  Don't change this routine unless you fully understand it.                */
4021 /*                                                                           */
4022 /*****************************************************************************/
4023 
4024 void exactinit()
4025 {
4026   REAL half;
4027   REAL check, lastcheck;
4028   int every_other;
4029 
4030   every_other = 1;
4031   half = 0.5;
4032   epsilon = 1.0;
4033   splitter = 1.0;
4034   check = 1.0;
4035   /* Repeatedly divide `epsilon' by two until it is too small to add to      */
4036   /*   one without causing roundoff.  (Also check if the sum is equal to     */
4037   /*   the previous sum, for machines that round up instead of using exact   */
4038   /*   rounding.  Not that these routines will work on such machines anyway. */
4039   do {
4040     lastcheck = check;
4041     epsilon *= half;
4042     if (every_other) {
4043       splitter *= 2.0;
4044     }
4045     every_other = !every_other;
4046     check = 1.0 + epsilon;
4047   } while ((check != 1.0) && (check != lastcheck));
4048   splitter += 1.0;
4049   if (verbose > 1) {
4050     printf("Floating point roundoff is of magnitude %.17g\n", epsilon);
4051     printf("Floating point splitter is %.17g\n", splitter);
4052   }
4053   /* Error bounds for orientation and incircle tests. */
4054   resulterrbound = (3.0 + 8.0 * epsilon) * epsilon;
4055   ccwerrboundA = (3.0 + 16.0 * epsilon) * epsilon;
4056   ccwerrboundB = (2.0 + 12.0 * epsilon) * epsilon;
4057   ccwerrboundC = (9.0 + 64.0 * epsilon) * epsilon * epsilon;
4058   iccerrboundA = (10.0 + 96.0 * epsilon) * epsilon;
4059   iccerrboundB = (4.0 + 48.0 * epsilon) * epsilon;
4060   iccerrboundC = (44.0 + 576.0 * epsilon) * epsilon * epsilon;
4061 }
4062 
4063 /*****************************************************************************/
4064 /*                                                                           */
4065 /*  fast_expansion_sum_zeroelim()   Sum two expansions, eliminating zero     */
4066 /*                                  components from the output expansion.    */
4067 /*                                                                           */
4068 /*  Sets h = e + f.  See my Robust Predicates paper for details.             */
4069 /*                                                                           */
4070 /*  If round-to-even is used (as with IEEE 754), maintains the strongly      */
4071 /*  nonoverlapping property.  (That is, if e is strongly nonoverlapping, h   */
4072 /*  will be also.)  Does NOT maintain the nonoverlapping or nonadjacent      */
4073 /*  properties.                                                              */
4074 /*                                                                           */
4075 /*****************************************************************************/
4076 
4077 int fast_expansion_sum_zeroelim(elen, e, flen, f, h)  /* h cannot be e or f. */
4078 int elen;
4079 REAL *e;
4080 int flen;
4081 REAL *f;
4082 REAL *h;
4083 {
4084   REAL Q;
4085   INEXACT REAL Qnew;
4086   INEXACT REAL hh;
4087   INEXACT REAL bvirt;
4088   REAL avirt, bround, around;
4089   int eindex, findex, hindex;
4090   REAL enow, fnow;
4091 
4092   enow = e[0];
4093   fnow = f[0];
4094   eindex = findex = 0;
4095   if ((fnow > enow) == (fnow > -enow)) {
4096     Q = enow;
4097     enow = e[++eindex];
4098   } else {
4099     Q = fnow;
4100     fnow = f[++findex];
4101   }
4102   hindex = 0;
4103   if ((eindex < elen) && (findex < flen)) {
4104     if ((fnow > enow) == (fnow > -enow)) {
4105       Fast_Two_Sum(enow, Q, Qnew, hh);
4106       enow = e[++eindex];
4107     } else {
4108       Fast_Two_Sum(fnow, Q, Qnew, hh);
4109       fnow = f[++findex];
4110     }
4111     Q = Qnew;
4112     if (hh != 0.0) {
4113       h[hindex++] = hh;
4114     }
4115     while ((eindex < elen) && (findex < flen)) {
4116       if ((fnow > enow) == (fnow > -enow)) {
4117         Two_Sum(Q, enow, Qnew, hh);
4118         enow = e[++eindex];
4119       } else {
4120         Two_Sum(Q, fnow, Qnew, hh);
4121         fnow = f[++findex];
4122       }
4123       Q = Qnew;
4124       if (hh != 0.0) {
4125         h[hindex++] = hh;
4126       }
4127     }
4128   }
4129   while (eindex < elen) {
4130     Two_Sum(Q, enow, Qnew, hh);
4131     enow = e[++eindex];
4132     Q = Qnew;
4133     if (hh != 0.0) {
4134       h[hindex++] = hh;
4135     }
4136   }
4137   while (findex < flen) {
4138     Two_Sum(Q, fnow, Qnew, hh);
4139     fnow = f[++findex];
4140     Q = Qnew;
4141     if (hh != 0.0) {
4142       h[hindex++] = hh;
4143     }
4144   }
4145   if ((Q != 0.0) || (hindex == 0)) {
4146     h[hindex++] = Q;
4147   }
4148   return hindex;
4149 }
4150 
4151 /*****************************************************************************/
4152 /*                                                                           */
4153 /*  scale_expansion_zeroelim()   Multiply an expansion by a scalar,          */
4154 /*                               eliminating zero components from the        */
4155 /*                               output expansion.                           */
4156 /*                                                                           */
4157 /*  Sets h = be.  See my Robust Predicates paper for details.                */
4158 /*                                                                           */
4159 /*  Maintains the nonoverlapping property.  If round-to-even is used (as     */
4160 /*  with IEEE 754), maintains the strongly nonoverlapping and nonadjacent    */
4161 /*  properties as well.  (That is, if e has one of these properties, so      */
4162 /*  will h.)                                                                 */
4163 /*                                                                           */
4164 /*****************************************************************************/
4165 
4166 int scale_expansion_zeroelim(elen, e, b, h)   /* e and h cannot be the same. */
4167 int elen;
4168 REAL *e;
4169 REAL b;
4170 REAL *h;
4171 {
4172   INEXACT REAL Q, sum;
4173   REAL hh;
4174   INEXACT REAL product1;
4175   REAL product0;
4176   int eindex, hindex;
4177   REAL enow;
4178   INEXACT REAL bvirt;
4179   REAL avirt, bround, around;
4180   INEXACT REAL c;
4181   INEXACT REAL abig;
4182   REAL ahi, alo, bhi, blo;
4183   REAL err1, err2, err3;
4184 
4185   Split(b, bhi, blo);
4186   Two_Product_Presplit(e[0], b, bhi, blo, Q, hh);
4187   hindex = 0;
4188   if (hh != 0) {
4189     h[hindex++] = hh;
4190   }
4191   for (eindex = 1; eindex < elen; eindex++) {
4192     enow = e[eindex];
4193     Two_Product_Presplit(enow, b, bhi, blo, product1, product0);
4194     Two_Sum(Q, product0, sum, hh);
4195     if (hh != 0) {
4196       h[hindex++] = hh;
4197     }
4198     Fast_Two_Sum(product1, sum, Q, hh);
4199     if (hh != 0) {
4200       h[hindex++] = hh;
4201     }
4202   }
4203   if ((Q != 0.0) || (hindex == 0)) {
4204     h[hindex++] = Q;
4205   }
4206   return hindex;
4207 }
4208 
4209 /*****************************************************************************/
4210 /*                                                                           */
4211 /*  estimate()   Produce a one-word estimate of an expansion's value.        */
4212 /*                                                                           */
4213 /*  See my Robust Predicates paper for details.                              */
4214 /*                                                                           */
4215 /*****************************************************************************/
4216 
4217 REAL estimate(elen, e)
4218 int elen;
4219 REAL *e;
4220 {
4221   REAL Q;
4222   int eindex;
4223 
4224   Q = e[0];
4225   for (eindex = 1; eindex < elen; eindex++) {
4226     Q += e[eindex];
4227   }
4228   return Q;
4229 }
4230 
4231 /*****************************************************************************/
4232 /*                                                                           */
4233 /*  counterclockwise()   Return a positive value if the points pa, pb, and   */
4234 /*                       pc occur in counterclockwise order; a negative      */
4235 /*                       value if they occur in clockwise order; and zero    */
4236 /*                       if they are collinear.  The result is also a rough  */
4237 /*                       approximation of twice the signed area of the       */
4238 /*                       triangle defined by the three points.               */
4239 /*                                                                           */
4240 /*  Uses exact arithmetic if necessary to ensure a correct answer.  The      */
4241 /*  result returned is the determinant of a matrix.  This determinant is     */
4242 /*  computed adaptively, in the sense that exact arithmetic is used only to  */
4243 /*  the degree it is needed to ensure that the returned value has the        */
4244 /*  correct sign.  Hence, this function is usually quite fast, but will run  */
4245 /*  more slowly when the input points are collinear or nearly so.            */
4246 /*                                                                           */
4247 /*  See my Robust Predicates paper for details.                              */
4248 /*                                                                           */
4249 /*****************************************************************************/
4250 
4251 REAL counterclockwiseadapt(pa, pb, pc, detsum)
4252 point pa;
4253 point pb;
4254 point pc;
4255 REAL detsum;
4256 {
4257   INEXACT REAL acx, acy, bcx, bcy;
4258   REAL acxtail, acytail, bcxtail, bcytail;
4259   INEXACT REAL detleft, detright;
4260   REAL detlefttail, detrighttail;
4261   REAL det, errbound;
4262   REAL B[4], C1[8], C2[12], D[16];
4263   INEXACT REAL B3;
4264   int C1length, C2length, Dlength;
4265   REAL u[4];
4266   INEXACT REAL u3;
4267   INEXACT REAL s1, t1;
4268   REAL s0, t0;
4269 
4270   INEXACT REAL bvirt;
4271   REAL avirt, bround, around;
4272   INEXACT REAL c;
4273   INEXACT REAL abig;
4274   REAL ahi, alo, bhi, blo;
4275   REAL err1, err2, err3;
4276   INEXACT REAL _i, _j;
4277   REAL _0;
4278 
4279   acx = (REAL) (pa[0] - pc[0]);
4280   bcx = (REAL) (pb[0] - pc[0]);
4281   acy = (REAL) (pa[1] - pc[1]);
4282   bcy = (REAL) (pb[1] - pc[1]);
4283 
4284   Two_Product(acx, bcy, detleft, detlefttail);
4285   Two_Product(acy, bcx, detright, detrighttail);
4286 
4287   Two_Two_Diff(detleft, detlefttail, detright, detrighttail,
4288                B3, B[2], B[1], B[0]);
4289   B[3] = B3;
4290 
4291   det = estimate(4, B);
4292   errbound = ccwerrboundB * detsum;
4293   if ((det >= errbound) || (-det >= errbound)) {
4294     return det;
4295   }
4296 
4297   Two_Diff_Tail(pa[0], pc[0], acx, acxtail);
4298   Two_Diff_Tail(pb[0], pc[0], bcx, bcxtail);
4299   Two_Diff_Tail(pa[1], pc[1], acy, acytail);
4300   Two_Diff_Tail(pb[1], pc[1], bcy, bcytail);
4301 
4302   if ((acxtail == 0.0) && (acytail == 0.0)
4303       && (bcxtail == 0.0) && (bcytail == 0.0)) {
4304     return det;
4305   }
4306 
4307   errbound = ccwerrboundC * detsum + resulterrbound * Absolute(det);
4308   det += (acx * bcytail + bcy * acxtail)
4309        - (acy * bcxtail + bcx * acytail);
4310   if ((det >= errbound) || (-det >= errbound)) {
4311     return det;
4312   }
4313 
4314   Two_Product(acxtail, bcy, s1, s0);
4315   Two_Product(acytail, bcx, t1, t0);
4316   Two_Two_Diff(s1, s0, t1, t0, u3, u[2], u[1], u[0]);
4317   u[3] = u3;
4318   C1length = fast_expansion_sum_zeroelim(4, B, 4, u, C1);
4319 
4320   Two_Product(acx, bcytail, s1, s0);
4321   Two_Product(acy, bcxtail, t1, t0);
4322   Two_Two_Diff(s1, s0, t1, t0, u3, u[2], u[1], u[0]);
4323   u[3] = u3;
4324   C2length = fast_expansion_sum_zeroelim(C1length, C1, 4, u, C2);
4325 
4326   Two_Product(acxtail, bcytail, s1, s0);
4327   Two_Product(acytail, bcxtail, t1, t0);
4328   Two_Two_Diff(s1, s0, t1, t0, u3, u[2], u[1], u[0]);
4329   u[3] = u3;
4330   Dlength = fast_expansion_sum_zeroelim(C2length, C2, 4, u, D);
4331 
4332   return(D[Dlength - 1]);
4333 }
4334 
4335 REAL counterclockwise(pa, pb, pc)
4336 point pa;
4337 point pb;
4338 point pc;
4339 {
4340   REAL detleft, detright, det;
4341   REAL detsum, errbound;
4342 
4343   counterclockcount++;
4344 
4345   detleft = (pa[0] - pc[0]) * (pb[1] - pc[1]);
4346   detright = (pa[1] - pc[1]) * (pb[0] - pc[0]);
4347   det = detleft - detright;
4348 
4349   if (noexact) {
4350     return det;
4351   }
4352 
4353   if (detleft > 0.0) {
4354     if (detright <= 0.0) {
4355       return det;
4356     } else {
4357       detsum = detleft + detright;
4358     }
4359   } else if (detleft < 0.0) {
4360     if (detright >= 0.0) {
4361       return det;
4362     } else {
4363       detsum = -detleft - detright;
4364     }
4365   } else {
4366     return det;
4367   }
4368 
4369   errbound = ccwerrboundA * detsum;
4370   if ((det >= errbound) || (-det >= errbound)) {
4371     return det;
4372   }
4373 
4374   return counterclockwiseadapt(pa, pb, pc, detsum);
4375 }
4376 
4377 /*****************************************************************************/
4378 /*                                                                           */
4379 /*  incircle()   Return a positive value if the point pd lies inside the     */
4380 /*               circle passing through pa, pb, and pc; a negative value if  */
4381 /*               it lies outside; and zero if the four points are cocircular.*/
4382 /*               The points pa, pb, and pc must be in counterclockwise       */
4383 /*               order, or the sign of the result will be reversed.          */
4384 /*                                                                           */
4385 /*  Uses exact arithmetic if necessary to ensure a correct answer.  The      */
4386 /*  result returned is the determinant of a matrix.  This determinant is     */
4387 /*  computed adaptively, in the sense that exact arithmetic is used only to  */
4388 /*  the degree it is needed to ensure that the returned value has the        */
4389 /*  correct sign.  Hence, this function is usually quite fast, but will run  */
4390 /*  more slowly when the input points are cocircular or nearly so.           */
4391 /*                                                                           */
4392 /*  See my Robust Predicates paper for details.                              */
4393 /*                                                                           */
4394 /*****************************************************************************/
4395 
4396 REAL incircleadapt(pa, pb, pc, pd, permanent)
4397 point pa;
4398 point pb;
4399 point pc;
4400 point pd;
4401 REAL permanent;
4402 {
4403   INEXACT REAL adx, bdx, cdx, ady, bdy, cdy;
4404   REAL det, errbound;
4405 
4406   INEXACT REAL bdxcdy1, cdxbdy1, cdxady1, adxcdy1, adxbdy1, bdxady1;
4407   REAL bdxcdy0, cdxbdy0, cdxady0, adxcdy0, adxbdy0, bdxady0;
4408   REAL bc[4], ca[4], ab[4];
4409   INEXACT REAL bc3, ca3, ab3;
4410   REAL axbc[8], axxbc[16], aybc[8], ayybc[16], adet[32];
4411   int axbclen, axxbclen, aybclen, ayybclen, alen;
4412   REAL bxca[8], bxxca[16], byca[8], byyca[16], bdet[32];
4413   int bxcalen, bxxcalen, bycalen, byycalen, blen;
4414   REAL cxab[8], cxxab[16], cyab[8], cyyab[16], cdet[32];
4415   int cxablen, cxxablen, cyablen, cyyablen, clen;
4416   REAL abdet[64];
4417   int ablen;
4418   REAL fin1[1152], fin2[1152];
4419   REAL *finnow, *finother, *finswap;
4420   int finlength;
4421 
4422   REAL adxtail, bdxtail, cdxtail, adytail, bdytail, cdytail;
4423   INEXACT REAL adxadx1, adyady1, bdxbdx1, bdybdy1, cdxcdx1, cdycdy1;
4424   REAL adxadx0, adyady0, bdxbdx0, bdybdy0, cdxcdx0, cdycdy0;
4425   REAL aa[4], bb[4], cc[4];
4426   INEXACT REAL aa3, bb3, cc3;
4427   INEXACT REAL ti1, tj1;
4428   REAL ti0, tj0;
4429   REAL u[4], v[4];
4430   INEXACT REAL u3, v3;
4431   REAL temp8[8], temp16a[16], temp16b[16], temp16c[16];
4432   REAL temp32a[32], temp32b[32], temp48[48], temp64[64];
4433   int temp8len, temp16alen, temp16blen, temp16clen;
4434   int temp32alen, temp32blen, temp48len, temp64len;
4435   REAL axtbb[8], axtcc[8], aytbb[8], aytcc[8];
4436   int axtbblen, axtcclen, aytbblen, aytcclen;
4437   REAL bxtaa[8], bxtcc[8], bytaa[8], bytcc[8];
4438   int bxtaalen, bxtcclen, bytaalen, bytcclen;
4439   REAL cxtaa[8], cxtbb[8], cytaa[8], cytbb[8];
4440   int cxtaalen, cxtbblen, cytaalen, cytbblen;
4441   REAL axtbc[8], aytbc[8], bxtca[8], bytca[8], cxtab[8], cytab[8];
4442   int axtbclen, aytbclen, bxtcalen, bytcalen, cxtablen, cytablen;
4443   REAL axtbct[16], aytbct[16], bxtcat[16], bytcat[16], cxtabt[16], cytabt[16];
4444   int axtbctlen, aytbctlen, bxtcatlen, bytcatlen, cxtabtlen, cytabtlen;
4445   REAL axtbctt[8], aytbctt[8], bxtcatt[8];
4446   REAL bytcatt[8], cxtabtt[8], cytabtt[8];
4447   int axtbcttlen, aytbcttlen, bxtcattlen, bytcattlen, cxtabttlen, cytabttlen;
4448   REAL abt[8], bct[8], cat[8];
4449   int abtlen, bctlen, catlen;
4450   REAL abtt[4], bctt[4], catt[4];
4451   int abttlen, bcttlen, cattlen;
4452   INEXACT REAL abtt3, bctt3, catt3;
4453   REAL negate;
4454 
4455   INEXACT REAL bvirt;
4456   REAL avirt, bround, around;
4457   INEXACT REAL c;
4458   INEXACT REAL abig;
4459   REAL ahi, alo, bhi, blo;
4460   REAL err1, err2, err3;
4461   INEXACT REAL _i, _j;
4462   REAL _0;
4463 
4464   adx = (REAL) (pa[0] - pd[0]);
4465   bdx = (REAL) (pb[0] - pd[0]);
4466   cdx = (REAL) (pc[0] - pd[0]);
4467   ady = (REAL) (pa[1] - pd[1]);
4468   bdy = (REAL) (pb[1] - pd[1]);
4469   cdy = (REAL) (pc[1] - pd[1]);
4470 
4471   Two_Product(bdx, cdy, bdxcdy1, bdxcdy0);
4472   Two_Product(cdx, bdy, cdxbdy1, cdxbdy0);
4473   Two_Two_Diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0, bc3, bc[2], bc[1], bc[0]);
4474   bc[3] = bc3;
4475   axbclen = scale_expansion_zeroelim(4, bc, adx, axbc);
4476   axxbclen = scale_expansion_zeroelim(axbclen, axbc, adx, axxbc);
4477   aybclen = scale_expansion_zeroelim(4, bc, ady, aybc);
4478   ayybclen = scale_expansion_zeroelim(aybclen, aybc, ady, ayybc);
4479   alen = fast_expansion_sum_zeroelim(axxbclen, axxbc, ayybclen, ayybc, adet);
4480 
4481   Two_Product(cdx, ady, cdxady1, cdxady0);
4482   Two_Product(adx, cdy, adxcdy1, adxcdy0);
4483   Two_Two_Diff(cdxady1, cdxady0, adxcdy1, adxcdy0, ca3, ca[2], ca[1], ca[0]);
4484   ca[3] = ca3;
4485   bxcalen = scale_expansion_zeroelim(4, ca, bdx, bxca);
4486   bxxcalen = scale_expansion_zeroelim(bxcalen, bxca, bdx, bxxca);
4487   bycalen = scale_expansion_zeroelim(4, ca, bdy, byca);
4488   byycalen = scale_expansion_zeroelim(bycalen, byca, bdy, byyca);
4489   blen = fast_expansion_sum_zeroelim(bxxcalen, bxxca, byycalen, byyca, bdet);
4490 
4491   Two_Product(adx, bdy, adxbdy1, adxbdy0);
4492   Two_Product(bdx, ady, bdxady1, bdxady0);
4493   Two_Two_Diff(adxbdy1, adxbdy0, bdxady1, bdxady0, ab3, ab[2], ab[1], ab[0]);
4494   ab[3] = ab3;
4495   cxablen = scale_expansion_zeroelim(4, ab, cdx, cxab);
4496   cxxablen = scale_expansion_zeroelim(cxablen, cxab, cdx, cxxab);
4497   cyablen = scale_expansion_zeroelim(4, ab, cdy, cyab);
4498   cyyablen = scale_expansion_zeroelim(cyablen, cyab, cdy, cyyab);
4499   clen = fast_expansion_sum_zeroelim(cxxablen, cxxab, cyyablen, cyyab, cdet);
4500 
4501   ablen = fast_expansion_sum_zeroelim(alen, adet, blen, bdet, abdet);
4502   finlength = fast_expansion_sum_zeroelim(ablen, abdet, clen, cdet, fin1);
4503 
4504   det = estimate(finlength, fin1);
4505   errbound = iccerrboundB * permanent;
4506   if ((det >= errbound) || (-det >= errbound)) {
4507     return det;
4508   }
4509 
4510   Two_Diff_Tail(pa[0], pd[0], adx, adxtail);
4511   Two_Diff_Tail(pa[1], pd[1], ady, adytail);
4512   Two_Diff_Tail(pb[0], pd[0], bdx, bdxtail);
4513   Two_Diff_Tail(pb[1], pd[1], bdy, bdytail);
4514   Two_Diff_Tail(pc[0], pd[0], cdx, cdxtail);
4515   Two_Diff_Tail(pc[1], pd[1], cdy, cdytail);
4516   if ((adxtail == 0.0) && (bdxtail == 0.0) && (cdxtail == 0.0)
4517       && (adytail == 0.0) && (bdytail == 0.0) && (cdytail == 0.0)) {
4518     return det;
4519   }
4520 
4521   errbound = iccerrboundC * permanent + resulterrbound * Absolute(det);
4522   det += ((adx * adx + ady * ady) * ((bdx * cdytail + cdy * bdxtail)
4523                                      - (bdy * cdxtail + cdx * bdytail))
4524           + 2.0 * (adx * adxtail + ady * adytail) * (bdx * cdy - bdy * cdx))
4525        + ((bdx * bdx + bdy * bdy) * ((cdx * adytail + ady * cdxtail)
4526                                      - (cdy * adxtail + adx * cdytail))
4527           + 2.0 * (bdx * bdxtail + bdy * bdytail) * (cdx * ady - cdy * adx))
4528        + ((cdx * cdx + cdy * cdy) * ((adx * bdytail + bdy * adxtail)
4529                                      - (ady * bdxtail + bdx * adytail))
4530           + 2.0 * (cdx * cdxtail + cdy * cdytail) * (adx * bdy - ady * bdx));
4531   if ((det >= errbound) || (-det >= errbound)) {
4532     return det;
4533   }
4534 
4535   finnow = fin1;
4536   finother = fin2;
4537 
4538   if ((bdxtail != 0.0) || (bdytail != 0.0)
4539       || (cdxtail != 0.0) || (cdytail != 0.0)) {
4540     Square(adx, adxadx1, adxadx0);
4541     Square(ady, adyady1, adyady0);
4542     Two_Two_Sum(adxadx1, adxadx0, adyady1, adyady0, aa3, aa[2], aa[1], aa[0]);
4543     aa[3] = aa3;
4544   }
4545   if ((cdxtail != 0.0) || (cdytail != 0.0)
4546       || (adxtail != 0.0) || (adytail != 0.0)) {
4547     Square(bdx, bdxbdx1, bdxbdx0);
4548     Square(bdy, bdybdy1, bdybdy0);
4549     Two_Two_Sum(bdxbdx1, bdxbdx0, bdybdy1, bdybdy0, bb3, bb[2], bb[1], bb[0]);
4550     bb[3] = bb3;
4551   }
4552   if ((adxtail != 0.0) || (adytail != 0.0)
4553       || (bdxtail != 0.0) || (bdytail != 0.0)) {
4554     Square(cdx, cdxcdx1, cdxcdx0);
4555     Square(cdy, cdycdy1, cdycdy0);
4556     Two_Two_Sum(cdxcdx1, cdxcdx0, cdycdy1, cdycdy0, cc3, cc[2], cc[1], cc[0]);
4557     cc[3] = cc3;
4558   }
4559 
4560   if (adxtail != 0.0) {
4561     axtbclen = scale_expansion_zeroelim(4, bc, adxtail, axtbc);
4562     temp16alen = scale_expansion_zeroelim(axtbclen, axtbc, 2.0 * adx,
4563                                           temp16a);
4564 
4565     axtcclen = scale_expansion_zeroelim(4, cc, adxtail, axtcc);
4566     temp16blen = scale_expansion_zeroelim(axtcclen, axtcc, bdy, temp16b);
4567 
4568     axtbblen = scale_expansion_zeroelim(4, bb, adxtail, axtbb);
4569     temp16clen = scale_expansion_zeroelim(axtbblen, axtbb, -cdy, temp16c);
4570 
4571     temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4572                                             temp16blen, temp16b, temp32a);
4573     temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
4574                                             temp32alen, temp32a, temp48);
4575     finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4576                                             temp48, finother);
4577     finswap = finnow; finnow = finother; finother = finswap;
4578   }
4579   if (adytail != 0.0) {
4580     aytbclen = scale_expansion_zeroelim(4, bc, adytail, aytbc);
4581     temp16alen = scale_expansion_zeroelim(aytbclen, aytbc, 2.0 * ady,
4582                                           temp16a);
4583 
4584     aytbblen = scale_expansion_zeroelim(4, bb, adytail, aytbb);
4585     temp16blen = scale_expansion_zeroelim(aytbblen, aytbb, cdx, temp16b);
4586 
4587     aytcclen = scale_expansion_zeroelim(4, cc, adytail, aytcc);
4588     temp16clen = scale_expansion_zeroelim(aytcclen, aytcc, -bdx, temp16c);
4589 
4590     temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4591                                             temp16blen, temp16b, temp32a);
4592     temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
4593                                             temp32alen, temp32a, temp48);
4594     finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4595                                             temp48, finother);
4596     finswap = finnow; finnow = finother; finother = finswap;
4597   }
4598   if (bdxtail != 0.0) {
4599     bxtcalen = scale_expansion_zeroelim(4, ca, bdxtail, bxtca);
4600     temp16alen = scale_expansion_zeroelim(bxtcalen, bxtca, 2.0 * bdx,
4601                                           temp16a);
4602 
4603     bxtaalen = scale_expansion_zeroelim(4, aa, bdxtail, bxtaa);
4604     temp16blen = scale_expansion_zeroelim(bxtaalen, bxtaa, cdy, temp16b);
4605 
4606     bxtcclen = scale_expansion_zeroelim(4, cc, bdxtail, bxtcc);
4607     temp16clen = scale_expansion_zeroelim(bxtcclen, bxtcc, -ady, temp16c);
4608 
4609     temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4610                                             temp16blen, temp16b, temp32a);
4611     temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
4612                                             temp32alen, temp32a, temp48);
4613     finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4614                                             temp48, finother);
4615     finswap = finnow; finnow = finother; finother = finswap;
4616   }
4617   if (bdytail != 0.0) {
4618     bytcalen = scale_expansion_zeroelim(4, ca, bdytail, bytca);
4619     temp16alen = scale_expansion_zeroelim(bytcalen, bytca, 2.0 * bdy,
4620                                           temp16a);
4621 
4622     bytcclen = scale_expansion_zeroelim(4, cc, bdytail, bytcc);
4623     temp16blen = scale_expansion_zeroelim(bytcclen, bytcc, adx, temp16b);
4624 
4625     bytaalen = scale_expansion_zeroelim(4, aa, bdytail, bytaa);
4626     temp16clen = scale_expansion_zeroelim(bytaalen, bytaa, -cdx, temp16c);
4627 
4628     temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4629                                             temp16blen, temp16b, temp32a);
4630     temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
4631                                             temp32alen, temp32a, temp48);
4632     finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4633                                             temp48, finother);
4634     finswap = finnow; finnow = finother; finother = finswap;
4635   }
4636   if (cdxtail != 0.0) {
4637     cxtablen = scale_expansion_zeroelim(4, ab, cdxtail, cxtab);
4638     temp16alen = scale_expansion_zeroelim(cxtablen, cxtab, 2.0 * cdx,
4639                                           temp16a);
4640 
4641     cxtbblen = scale_expansion_zeroelim(4, bb, cdxtail, cxtbb);
4642     temp16blen = scale_expansion_zeroelim(cxtbblen, cxtbb, ady, temp16b);
4643 
4644     cxtaalen = scale_expansion_zeroelim(4, aa, cdxtail, cxtaa);
4645     temp16clen = scale_expansion_zeroelim(cxtaalen, cxtaa, -bdy, temp16c);
4646 
4647     temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4648                                             temp16blen, temp16b, temp32a);
4649     temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
4650                                             temp32alen, temp32a, temp48);
4651     finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4652                                             temp48, finother);
4653     finswap = finnow; finnow = finother; finother = finswap;
4654   }
4655   if (cdytail != 0.0) {
4656     cytablen = scale_expansion_zeroelim(4, ab, cdytail, cytab);
4657     temp16alen = scale_expansion_zeroelim(cytablen, cytab, 2.0 * cdy,
4658                                           temp16a);
4659 
4660     cytaalen = scale_expansion_zeroelim(4, aa, cdytail, cytaa);
4661     temp16blen = scale_expansion_zeroelim(cytaalen, cytaa, bdx, temp16b);
4662 
4663     cytbblen = scale_expansion_zeroelim(4, bb, cdytail, cytbb);
4664     temp16clen = scale_expansion_zeroelim(cytbblen, cytbb, -adx, temp16c);
4665 
4666     temp32alen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4667                                             temp16blen, temp16b, temp32a);
4668     temp48len = fast_expansion_sum_zeroelim(temp16clen, temp16c,
4669                                             temp32alen, temp32a, temp48);
4670     finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4671                                             temp48, finother);
4672     finswap = finnow; finnow = finother; finother = finswap;
4673   }
4674 
4675   if ((adxtail != 0.0) || (adytail != 0.0)) {
4676     if ((bdxtail != 0.0) || (bdytail != 0.0)
4677         || (cdxtail != 0.0) || (cdytail != 0.0)) {
4678       Two_Product(bdxtail, cdy, ti1, ti0);
4679       Two_Product(bdx, cdytail, tj1, tj0);
4680       Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
4681       u[3] = u3;
4682       negate = -bdy;
4683       Two_Product(cdxtail, negate, ti1, ti0);
4684       negate = -bdytail;
4685       Two_Product(cdx, negate, tj1, tj0);
4686       Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
4687       v[3] = v3;
4688       bctlen = fast_expansion_sum_zeroelim(4, u, 4, v, bct);
4689 
4690       Two_Product(bdxtail, cdytail, ti1, ti0);
4691       Two_Product(cdxtail, bdytail, tj1, tj0);
4692       Two_Two_Diff(ti1, ti0, tj1, tj0, bctt3, bctt[2], bctt[1], bctt[0]);
4693       bctt[3] = bctt3;
4694       bcttlen = 4;
4695     } else {
4696       bct[0] = 0.0;
4697       bctlen = 1;
4698       bctt[0] = 0.0;
4699       bcttlen = 1;
4700     }
4701 
4702     if (adxtail != 0.0) {
4703       temp16alen = scale_expansion_zeroelim(axtbclen, axtbc, adxtail, temp16a);
4704       axtbctlen = scale_expansion_zeroelim(bctlen, bct, adxtail, axtbct);
4705       temp32alen = scale_expansion_zeroelim(axtbctlen, axtbct, 2.0 * adx,
4706                                             temp32a);
4707       temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4708                                               temp32alen, temp32a, temp48);
4709       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4710                                               temp48, finother);
4711       finswap = finnow; finnow = finother; finother = finswap;
4712       if (bdytail != 0.0) {
4713         temp8len = scale_expansion_zeroelim(4, cc, adxtail, temp8);
4714         temp16alen = scale_expansion_zeroelim(temp8len, temp8, bdytail,
4715                                               temp16a);
4716         finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
4717                                                 temp16a, finother);
4718         finswap = finnow; finnow = finother; finother = finswap;
4719       }
4720       if (cdytail != 0.0) {
4721         temp8len = scale_expansion_zeroelim(4, bb, -adxtail, temp8);
4722         temp16alen = scale_expansion_zeroelim(temp8len, temp8, cdytail,
4723                                               temp16a);
4724         finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
4725                                                 temp16a, finother);
4726         finswap = finnow; finnow = finother; finother = finswap;
4727       }
4728 
4729       temp32alen = scale_expansion_zeroelim(axtbctlen, axtbct, adxtail,
4730                                             temp32a);
4731       axtbcttlen = scale_expansion_zeroelim(bcttlen, bctt, adxtail, axtbctt);
4732       temp16alen = scale_expansion_zeroelim(axtbcttlen, axtbctt, 2.0 * adx,
4733                                             temp16a);
4734       temp16blen = scale_expansion_zeroelim(axtbcttlen, axtbctt, adxtail,
4735                                             temp16b);
4736       temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4737                                               temp16blen, temp16b, temp32b);
4738       temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
4739                                               temp32blen, temp32b, temp64);
4740       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
4741                                               temp64, finother);
4742       finswap = finnow; finnow = finother; finother = finswap;
4743     }
4744     if (adytail != 0.0) {
4745       temp16alen = scale_expansion_zeroelim(aytbclen, aytbc, adytail, temp16a);
4746       aytbctlen = scale_expansion_zeroelim(bctlen, bct, adytail, aytbct);
4747       temp32alen = scale_expansion_zeroelim(aytbctlen, aytbct, 2.0 * ady,
4748                                             temp32a);
4749       temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4750                                               temp32alen, temp32a, temp48);
4751       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4752                                               temp48, finother);
4753       finswap = finnow; finnow = finother; finother = finswap;
4754 
4755 
4756       temp32alen = scale_expansion_zeroelim(aytbctlen, aytbct, adytail,
4757                                             temp32a);
4758       aytbcttlen = scale_expansion_zeroelim(bcttlen, bctt, adytail, aytbctt);
4759       temp16alen = scale_expansion_zeroelim(aytbcttlen, aytbctt, 2.0 * ady,
4760                                             temp16a);
4761       temp16blen = scale_expansion_zeroelim(aytbcttlen, aytbctt, adytail,
4762                                             temp16b);
4763       temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4764                                               temp16blen, temp16b, temp32b);
4765       temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
4766                                               temp32blen, temp32b, temp64);
4767       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
4768                                               temp64, finother);
4769       finswap = finnow; finnow = finother; finother = finswap;
4770     }
4771   }
4772   if ((bdxtail != 0.0) || (bdytail != 0.0)) {
4773     if ((cdxtail != 0.0) || (cdytail != 0.0)
4774         || (adxtail != 0.0) || (adytail != 0.0)) {
4775       Two_Product(cdxtail, ady, ti1, ti0);
4776       Two_Product(cdx, adytail, tj1, tj0);
4777       Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
4778       u[3] = u3;
4779       negate = -cdy;
4780       Two_Product(adxtail, negate, ti1, ti0);
4781       negate = -cdytail;
4782       Two_Product(adx, negate, tj1, tj0);
4783       Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
4784       v[3] = v3;
4785       catlen = fast_expansion_sum_zeroelim(4, u, 4, v, cat);
4786 
4787       Two_Product(cdxtail, adytail, ti1, ti0);
4788       Two_Product(adxtail, cdytail, tj1, tj0);
4789       Two_Two_Diff(ti1, ti0, tj1, tj0, catt3, catt[2], catt[1], catt[0]);
4790       catt[3] = catt3;
4791       cattlen = 4;
4792     } else {
4793       cat[0] = 0.0;
4794       catlen = 1;
4795       catt[0] = 0.0;
4796       cattlen = 1;
4797     }
4798 
4799     if (bdxtail != 0.0) {
4800       temp16alen = scale_expansion_zeroelim(bxtcalen, bxtca, bdxtail, temp16a);
4801       bxtcatlen = scale_expansion_zeroelim(catlen, cat, bdxtail, bxtcat);
4802       temp32alen = scale_expansion_zeroelim(bxtcatlen, bxtcat, 2.0 * bdx,
4803                                             temp32a);
4804       temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4805                                               temp32alen, temp32a, temp48);
4806       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4807                                               temp48, finother);
4808       finswap = finnow; finnow = finother; finother = finswap;
4809       if (cdytail != 0.0) {
4810         temp8len = scale_expansion_zeroelim(4, aa, bdxtail, temp8);
4811         temp16alen = scale_expansion_zeroelim(temp8len, temp8, cdytail,
4812                                               temp16a);
4813         finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
4814                                                 temp16a, finother);
4815         finswap = finnow; finnow = finother; finother = finswap;
4816       }
4817       if (adytail != 0.0) {
4818         temp8len = scale_expansion_zeroelim(4, cc, -bdxtail, temp8);
4819         temp16alen = scale_expansion_zeroelim(temp8len, temp8, adytail,
4820                                               temp16a);
4821         finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
4822                                                 temp16a, finother);
4823         finswap = finnow; finnow = finother; finother = finswap;
4824       }
4825 
4826       temp32alen = scale_expansion_zeroelim(bxtcatlen, bxtcat, bdxtail,
4827                                             temp32a);
4828       bxtcattlen = scale_expansion_zeroelim(cattlen, catt, bdxtail, bxtcatt);
4829       temp16alen = scale_expansion_zeroelim(bxtcattlen, bxtcatt, 2.0 * bdx,
4830                                             temp16a);
4831       temp16blen = scale_expansion_zeroelim(bxtcattlen, bxtcatt, bdxtail,
4832                                             temp16b);
4833       temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4834                                               temp16blen, temp16b, temp32b);
4835       temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
4836                                               temp32blen, temp32b, temp64);
4837       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
4838                                               temp64, finother);
4839       finswap = finnow; finnow = finother; finother = finswap;
4840     }
4841     if (bdytail != 0.0) {
4842       temp16alen = scale_expansion_zeroelim(bytcalen, bytca, bdytail, temp16a);
4843       bytcatlen = scale_expansion_zeroelim(catlen, cat, bdytail, bytcat);
4844       temp32alen = scale_expansion_zeroelim(bytcatlen, bytcat, 2.0 * bdy,
4845                                             temp32a);
4846       temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4847                                               temp32alen, temp32a, temp48);
4848       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4849                                               temp48, finother);
4850       finswap = finnow; finnow = finother; finother = finswap;
4851 
4852 
4853       temp32alen = scale_expansion_zeroelim(bytcatlen, bytcat, bdytail,
4854                                             temp32a);
4855       bytcattlen = scale_expansion_zeroelim(cattlen, catt, bdytail, bytcatt);
4856       temp16alen = scale_expansion_zeroelim(bytcattlen, bytcatt, 2.0 * bdy,
4857                                             temp16a);
4858       temp16blen = scale_expansion_zeroelim(bytcattlen, bytcatt, bdytail,
4859                                             temp16b);
4860       temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4861                                               temp16blen, temp16b, temp32b);
4862       temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
4863                                               temp32blen, temp32b, temp64);
4864       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
4865                                               temp64, finother);
4866       finswap = finnow; finnow = finother; finother = finswap;
4867     }
4868   }
4869   if ((cdxtail != 0.0) || (cdytail != 0.0)) {
4870     if ((adxtail != 0.0) || (adytail != 0.0)
4871         || (bdxtail != 0.0) || (bdytail != 0.0)) {
4872       Two_Product(adxtail, bdy, ti1, ti0);
4873       Two_Product(adx, bdytail, tj1, tj0);
4874       Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
4875       u[3] = u3;
4876       negate = -ady;
4877       Two_Product(bdxtail, negate, ti1, ti0);
4878       negate = -adytail;
4879       Two_Product(bdx, negate, tj1, tj0);
4880       Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
4881       v[3] = v3;
4882       abtlen = fast_expansion_sum_zeroelim(4, u, 4, v, abt);
4883 
4884       Two_Product(adxtail, bdytail, ti1, ti0);
4885       Two_Product(bdxtail, adytail, tj1, tj0);
4886       Two_Two_Diff(ti1, ti0, tj1, tj0, abtt3, abtt[2], abtt[1], abtt[0]);
4887       abtt[3] = abtt3;
4888       abttlen = 4;
4889     } else {
4890       abt[0] = 0.0;
4891       abtlen = 1;
4892       abtt[0] = 0.0;
4893       abttlen = 1;
4894     }
4895 
4896     if (cdxtail != 0.0) {
4897       temp16alen = scale_expansion_zeroelim(cxtablen, cxtab, cdxtail, temp16a);
4898       cxtabtlen = scale_expansion_zeroelim(abtlen, abt, cdxtail, cxtabt);
4899       temp32alen = scale_expansion_zeroelim(cxtabtlen, cxtabt, 2.0 * cdx,
4900                                             temp32a);
4901       temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4902                                               temp32alen, temp32a, temp48);
4903       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4904                                               temp48, finother);
4905       finswap = finnow; finnow = finother; finother = finswap;
4906       if (adytail != 0.0) {
4907         temp8len = scale_expansion_zeroelim(4, bb, cdxtail, temp8);
4908         temp16alen = scale_expansion_zeroelim(temp8len, temp8, adytail,
4909                                               temp16a);
4910         finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
4911                                                 temp16a, finother);
4912         finswap = finnow; finnow = finother; finother = finswap;
4913       }
4914       if (bdytail != 0.0) {
4915         temp8len = scale_expansion_zeroelim(4, aa, -cdxtail, temp8);
4916         temp16alen = scale_expansion_zeroelim(temp8len, temp8, bdytail,
4917                                               temp16a);
4918         finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp16alen,
4919                                                 temp16a, finother);
4920         finswap = finnow; finnow = finother; finother = finswap;
4921       }
4922 
4923       temp32alen = scale_expansion_zeroelim(cxtabtlen, cxtabt, cdxtail,
4924                                             temp32a);
4925       cxtabttlen = scale_expansion_zeroelim(abttlen, abtt, cdxtail, cxtabtt);
4926       temp16alen = scale_expansion_zeroelim(cxtabttlen, cxtabtt, 2.0 * cdx,
4927                                             temp16a);
4928       temp16blen = scale_expansion_zeroelim(cxtabttlen, cxtabtt, cdxtail,
4929                                             temp16b);
4930       temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4931                                               temp16blen, temp16b, temp32b);
4932       temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
4933                                               temp32blen, temp32b, temp64);
4934       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
4935                                               temp64, finother);
4936       finswap = finnow; finnow = finother; finother = finswap;
4937     }
4938     if (cdytail != 0.0) {
4939       temp16alen = scale_expansion_zeroelim(cytablen, cytab, cdytail, temp16a);
4940       cytabtlen = scale_expansion_zeroelim(abtlen, abt, cdytail, cytabt);
4941       temp32alen = scale_expansion_zeroelim(cytabtlen, cytabt, 2.0 * cdy,
4942                                             temp32a);
4943       temp48len = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4944                                               temp32alen, temp32a, temp48);
4945       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp48len,
4946                                               temp48, finother);
4947       finswap = finnow; finnow = finother; finother = finswap;
4948 
4949 
4950       temp32alen = scale_expansion_zeroelim(cytabtlen, cytabt, cdytail,
4951                                             temp32a);
4952       cytabttlen = scale_expansion_zeroelim(abttlen, abtt, cdytail, cytabtt);
4953       temp16alen = scale_expansion_zeroelim(cytabttlen, cytabtt, 2.0 * cdy,
4954                                             temp16a);
4955       temp16blen = scale_expansion_zeroelim(cytabttlen, cytabtt, cdytail,
4956                                             temp16b);
4957       temp32blen = fast_expansion_sum_zeroelim(temp16alen, temp16a,
4958                                               temp16blen, temp16b, temp32b);
4959       temp64len = fast_expansion_sum_zeroelim(temp32alen, temp32a,
4960                                               temp32blen, temp32b, temp64);
4961       finlength = fast_expansion_sum_zeroelim(finlength, finnow, temp64len,
4962                                               temp64, finother);
4963       finswap = finnow; finnow = finother; finother = finswap;
4964     }
4965   }
4966 
4967   return finnow[finlength - 1];
4968 }
4969 
4970 REAL incircle(pa, pb, pc, pd)
4971 point pa;
4972 point pb;
4973 point pc;
4974 point pd;
4975 {
4976   REAL adx, bdx, cdx, ady, bdy, cdy;
4977   REAL bdxcdy, cdxbdy, cdxady, adxcdy, adxbdy, bdxady;
4978   REAL alift, blift, clift;
4979   REAL det;
4980   REAL permanent, errbound;
4981 
4982   incirclecount++;
4983 
4984   adx = pa[0] - pd[0];
4985   bdx = pb[0] - pd[0];
4986   cdx = pc[0] - pd[0];
4987   ady = pa[1] - pd[1];
4988   bdy = pb[1] - pd[1];
4989   cdy = pc[1] - pd[1];
4990 
4991   bdxcdy = bdx * cdy;
4992   cdxbdy = cdx * bdy;
4993   alift = adx * adx + ady * ady;
4994 
4995   cdxady = cdx * ady;
4996   adxcdy = adx * cdy;
4997   blift = bdx * bdx + bdy * bdy;
4998 
4999   adxbdy = adx * bdy;
5000   bdxady = bdx * ady;
5001   clift = cdx * cdx + cdy * cdy;
5002 
5003   det = alift * (bdxcdy - cdxbdy)
5004       + blift * (cdxady - adxcdy)
5005       + clift * (adxbdy - bdxady);
5006 
5007   if (noexact) {
5008     return det;
5009   }
5010 
5011   permanent = (Absolute(bdxcdy) + Absolute(cdxbdy)) * alift
5012             + (Absolute(cdxady) + Absolute(adxcdy)) * blift
5013             + (Absolute(adxbdy) + Absolute(bdxady)) * clift;
5014   errbound = iccerrboundA * permanent;
5015   if ((det > errbound) || (-det > errbound)) {
5016     return det;
5017   }
5018 
5019   return incircleadapt(pa, pb, pc, pd, permanent);
5020 }
5021 
5022 /**                                                                         **/
5023 /**                                                                         **/
5024 /********* Determinant evaluation routines end here                  *********/
5025 
5026 /*****************************************************************************/
5027 /*                                                                           */
5028 /*  triangleinit()   Initialize some variables.                              */
5029 /*                                                                           */
5030 /*****************************************************************************/
5031 
5032 void triangleinit()
5033 {
5034   points.maxitems = triangles.maxitems = shelles.maxitems = viri.maxitems =
5035     badsegments.maxitems = badtriangles.maxitems = splaynodes.maxitems = 0l;
5036   points.itembytes = triangles.itembytes = shelles.itembytes = viri.itembytes =
5037     badsegments.itembytes = badtriangles.itembytes = splaynodes.itembytes = 0;
5038   recenttri.tri = (triangle *) NULL;    /* No triangle has been visited yet. */
5039   samples = 1;            /* Point location should take at least one sample. */
5040   checksegments = 0;      /* There are no segments in the triangulation yet. */
5041   incirclecount = counterclockcount = hyperbolacount = 0;
5042   circumcentercount = circletopcount = 0;
5043   randomseed = 1;
5044 
5045   exactinit();                     /* Initialize exact arithmetic constants. */
5046 }
5047 
5048 /*****************************************************************************/
5049 /*                                                                           */
5050 /*  randomnation()   Generate a random number between 0 and `choices' - 1.   */
5051 /*                                                                           */
5052 /*  This is a simple linear congruential random number generator.  Hence, it */
5053 /*  is a bad random number generator, but good enough for most randomized    */
5054 /*  geometric algorithms.                                                    */
5055 /*                                                                           */
5056 /*****************************************************************************/
5057 
5058 unsigned long randomnation(choices)
5059 unsigned int choices;
5060 {
5061   randomseed = (randomseed * 1366l + 150889l) % 714025l;
5062   return randomseed / (714025l / choices + 1);
5063 }
5064 
5065 /********* Mesh quality testing routines begin here                  *********/
5066 /**                                                                         **/
5067 /**                                                                         **/
5068 
5069 /*****************************************************************************/
5070 /*                                                                           */
5071 /*  checkmesh()   Test the mesh for topological consistency.                 */
5072 /*                                                                           */
5073 /*****************************************************************************/
5074 
5075 #ifndef REDUCED
5076 
5077 void checkmesh()
5078 {
5079   struct triedge triangleloop;
5080   struct triedge oppotri, oppooppotri;
5081   point triorg, tridest, triapex;
5082   point oppoorg, oppodest;
5083   int horrors;
5084   int saveexact;
5085   triangle ptr;                         /* Temporary variable used by sym(). */
5086 
5087   /* Temporarily turn on exact arithmetic if it's off. */
5088   saveexact = noexact;
5089   noexact = 0;
5090   if (!quiet) {
5091     printf("  Checking consistency of mesh...\n");
5092   }
5093   horrors = 0;
5094   /* Run through the list of triangles, checking each one. */
5095   traversalinit(&triangles);
5096   triangleloop.tri = triangletraverse();
5097   while (triangleloop.tri != (triangle *) NULL) {
5098     /* Check all three edges of the triangle. */
5099     for (triangleloop.orient = 0; triangleloop.orient < 3;
5100          triangleloop.orient++) {
5101       org(triangleloop, triorg);
5102       dest(triangleloop, tridest);
5103       if (triangleloop.orient == 0) {       /* Only test for inversion once. */
5104         /* Test if the triangle is flat or inverted. */
5105         apex(triangleloop, triapex);
5106         if (counterclockwise(triorg, tridest, triapex) <= 0.0) {
5107           printf("  !! !! Inverted ");
5108           printtriangle(&triangleloop);
5109           horrors++;
5110         }
5111       }
5112       /* Find the neighboring triangle on this edge. */
5113       sym(triangleloop, oppotri);
5114       if (oppotri.tri != dummytri) {
5115         /* Check that the triangle's neighbor knows it's a neighbor. */
5116         sym(oppotri, oppooppotri);
5117         if ((triangleloop.tri != oppooppotri.tri)
5118             || (triangleloop.orient != oppooppotri.orient)) {
5119           printf("  !! !! Asymmetric triangle-triangle bond:\n");
5120           if (triangleloop.tri == oppooppotri.tri) {
5121             printf("   (Right triangle, wrong orientation)\n");
5122           }
5123           printf("    First ");
5124           printtriangle(&triangleloop);
5125           printf("    Second (nonreciprocating) ");
5126           printtriangle(&oppotri);
5127           horrors++;
5128         }
5129         /* Check that both triangles agree on the identities */
5130         /*   of their shared vertices.                       */
5131         org(oppotri, oppoorg);
5132         dest(oppotri, oppodest);
5133         if ((triorg != oppodest) || (tridest != oppoorg)) {
5134           printf("  !! !! Mismatched edge coordinates between two triangles:\n"
5135                  );
5136           printf("    First mismatched ");
5137           printtriangle(&triangleloop);
5138           printf("    Second mismatched ");
5139           printtriangle(&oppotri);
5140           horrors++;
5141         }
5142       }
5143     }
5144     triangleloop.tri = triangletraverse();
5145   }
5146   if (horrors == 0) {
5147     if (!quiet) {
5148       printf("  In my studied opinion, the mesh appears to be consistent.\n");
5149     }
5150   } else if (horrors == 1) {
5151     printf("  !! !! !! !! Precisely one festering wound discovered.\n");
5152   } else {
5153     printf("  !! !! !! !! %d abominations witnessed.\n", horrors);
5154   }
5155   /* Restore the status of exact arithmetic. */
5156   noexact = saveexact;
5157 }
5158 
5159 #endif /* not REDUCED */
5160 
5161 /*****************************************************************************/
5162 /*                                                                           */
5163 /*  checkdelaunay()   Ensure that the mesh is (constrained) Delaunay.        */
5164 /*                                                                           */
5165 /*****************************************************************************/
5166 
5167 #ifndef REDUCED
5168 
5169 void checkdelaunay()
5170 {
5171   struct triedge triangleloop;
5172   struct triedge oppotri;
5173   struct edge opposhelle;
5174   point triorg, tridest, triapex;
5175   point oppoapex;
5176   int shouldbedelaunay;
5177   int horrors;
5178   int saveexact;
5179   triangle ptr;                         /* Temporary variable used by sym(). */
5180   shelle sptr;                      /* Temporary variable used by tspivot(). */
5181 
5182   /* Temporarily turn on exact arithmetic if it's off. */
5183   saveexact = noexact;
5184   noexact = 0;
5185   if (!quiet) {
5186     printf("  Checking Delaunay property of mesh...\n");
5187   }
5188   horrors = 0;
5189   /* Run through the list of triangles, checking each one. */
5190   traversalinit(&triangles);
5191   triangleloop.tri = triangletraverse();
5192   while (triangleloop.tri != (triangle *) NULL) {
5193     /* Check all three edges of the triangle. */
5194     for (triangleloop.orient = 0; triangleloop.orient < 3;
5195          triangleloop.orient++) {
5196       org(triangleloop, triorg);
5197       dest(triangleloop, tridest);
5198       apex(triangleloop, triapex);
5199       sym(triangleloop, oppotri);
5200       apex(oppotri, oppoapex);
5201       /* Only test that the edge is locally Delaunay if there is an   */
5202       /*   adjoining triangle whose pointer is larger (to ensure that */
5203       /*   each pair isn't tested twice).                             */
5204       shouldbedelaunay = (oppotri.tri != dummytri)
5205             && (triapex != (point) NULL) && (oppoapex != (point) NULL)
5206             && (triangleloop.tri < oppotri.tri);
5207       if (checksegments && shouldbedelaunay) {
5208         /* If a shell edge separates the triangles, then the edge is */
5209         /*   constrained, so no local Delaunay test should be done.  */
5210         tspivot(triangleloop, opposhelle);
5211         if (opposhelle.sh != dummysh){
5212           shouldbedelaunay = 0;
5213         }
5214       }
5215       if (shouldbedelaunay) {
5216         if (incircle(triorg, tridest, triapex, oppoapex) > 0.0) {
5217           printf("  !! !! Non-Delaunay pair of triangles:\n");
5218           printf("    First non-Delaunay ");
5219           printtriangle(&triangleloop);
5220           printf("    Second non-Delaunay ");
5221           printtriangle(&oppotri);
5222           horrors++;
5223         }
5224       }
5225     }
5226     triangleloop.tri = triangletraverse();
5227   }
5228   if (horrors == 0) {
5229     if (!quiet) {
5230       printf(
5231   "  By virtue of my perceptive intelligence, I declare the mesh Delaunay.\n");
5232     }
5233   } else if (horrors == 1) {
5234     printf(
5235          "  !! !! !! !! Precisely one terrifying transgression identified.\n");
5236   } else {
5237     printf("  !! !! !! !! %d obscenities viewed with horror.\n", horrors);
5238   }
5239   /* Restore the status of exact arithmetic. */
5240   noexact = saveexact;
5241 }
5242 
5243 #endif /* not REDUCED */
5244 
5245 /*****************************************************************************/
5246 /*                                                                           */
5247 /*  enqueuebadtri()   Add a bad triangle to the end of a queue.              */
5248 /*                                                                           */
5249 /*  The queue is actually a set of 64 queues.  I use multiple queues to give */
5250 /*  priority to smaller angles.  I originally implemented a heap, but the    */
5251 /*  queues are (to my surprise) much faster.                                 */
5252 /*                                                                           */
5253 /*****************************************************************************/
5254 
5255 #ifndef CDT_ONLY
5256 
5257 void enqueuebadtri(instri, angle, insapex, insorg, insdest)
5258 struct triedge *instri;
5259 REAL angle;
5260 point insapex;
5261 point insorg;
5262 point insdest;
5263 {
5264   struct badface *newface;
5265   int queuenumber;
5266 
5267   if (verbose > 2) {
5268     printf("  Queueing bad triangle:\n");
5269     printf("    (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n", insorg[0],
5270            insorg[1], insdest[0], insdest[1], insapex[0], insapex[1]);
5271   }
5272   /* Allocate space for the bad triangle. */
5273   newface = (struct badface *) poolalloc(&badtriangles);
5274   triedgecopy(*instri, newface->badfacetri);
5275   newface->key = angle;
5276   newface->faceapex = insapex;
5277   newface->faceorg = insorg;
5278   newface->facedest = insdest;
5279   newface->nextface = (struct badface *) NULL;
5280   /* Determine the appropriate queue to put the bad triangle into. */
5281   if (angle > 0.6) {
5282     queuenumber = (int) (160.0 * (angle - 0.6));
5283     if (queuenumber > 63) {
5284       queuenumber = 63;
5285     }
5286   } else {
5287     /* It's not a bad angle; put the triangle in the lowest-priority queue. */
5288     queuenumber = 0;
5289   }
5290   /* Add the triangle to the end of a queue. */
5291   *queuetail[queuenumber] = newface;
5292   /* Maintain a pointer to the NULL pointer at the end of the queue. */
5293   queuetail[queuenumber] = &newface->nextface;
5294 }
5295 
5296 #endif /* not CDT_ONLY */
5297 
5298 /*****************************************************************************/
5299 /*                                                                           */
5300 /*  dequeuebadtri()   Remove a triangle from the front of the queue.         */
5301 /*                                                                           */
5302 /*****************************************************************************/
5303 
5304 #ifndef CDT_ONLY
5305 
5306 struct badface *dequeuebadtri()
5307 {
5308   struct badface *result;
5309   int queuenumber;
5310 
5311   /* Look for a nonempty queue. */
5312   for (queuenumber = 63; queuenumber >= 0; queuenumber--) {
5313     result = queuefront[queuenumber];
5314     if (result != (struct badface *) NULL) {
5315       /* Remove the triangle from the queue. */
5316       queuefront[queuenumber] = result->nextface;
5317       /* Maintain a pointer to the NULL pointer at the end of the queue. */
5318       if (queuefront[queuenumber] == (struct badface *) NULL) {
5319         queuetail[queuenumber] = &queuefront[queuenumber];
5320       }
5321       return result;
5322     }
5323   }
5324   return (struct badface *) NULL;
5325 }
5326 
5327 #endif /* not CDT_ONLY */
5328 
5329 /*****************************************************************************/
5330 /*                                                                           */
5331 /*  checkedge4encroach()   Check a segment to see if it is encroached; add   */
5332 /*                         it to the list if it is.                          */
5333 /*                                                                           */
5334 /*  An encroached segment is an unflippable edge that has a point in its     */
5335 /*  diametral circle (that is, it faces an angle greater than 90 degrees).   */
5336 /*  This definition is due to Ruppert.                                       */
5337 /*                                                                           */
5338 /*  Returns a nonzero value if the edge is encroached.                       */
5339 /*                                                                           */
5340 /*****************************************************************************/
5341 
5342 #ifndef CDT_ONLY
5343 
5344 int checkedge4encroach(testedge)
5345 struct edge *testedge;
5346 {
5347   struct triedge neighbortri;
5348   struct edge testsym;
5349   struct edge *badedge;
5350   int addtolist;
5351   int sides;
5352   point eorg, edest, eapex;
5353   triangle ptr;                     /* Temporary variable used by stpivot(). */
5354 
5355   addtolist = 0;
5356   sides = 0;
5357 
5358   sorg(*testedge, eorg);
5359   sdest(*testedge, edest);
5360   /* Check one neighbor of the shell edge. */
5361   stpivot(*testedge, neighbortri);
5362   /* Does the neighbor exist, or is this a boundary edge? */
5363   if (neighbortri.tri != dummytri) {
5364     sides++;
5365     /* Find a vertex opposite this edge. */
5366     apex(neighbortri, eapex);
5367     /* Check whether the vertex is inside the diametral circle of the  */
5368     /*   shell edge.  Pythagoras' Theorem is used to check whether the */
5369     /*   angle at the vertex is greater than 90 degrees.               */
5370     if (eapex[0] * (eorg[0] + edest[0]) + eapex[1] * (eorg[1] + edest[1]) >
5371         eapex[0] * eapex[0] + eorg[0] * edest[0] +
5372         eapex[1] * eapex[1] + eorg[1] * edest[1]) {
5373       addtolist = 1;
5374     }
5375   }
5376   /* Check the other neighbor of the shell edge. */
5377   ssym(*testedge, testsym);
5378   stpivot(testsym, neighbortri);
5379   /* Does the neighbor exist, or is this a boundary edge? */
5380   if (neighbortri.tri != dummytri) {
5381     sides++;
5382     /* Find the other vertex opposite this edge. */
5383     apex(neighbortri, eapex);
5384     /* Check whether the vertex is inside the diametral circle of the  */
5385     /*   shell edge.  Pythagoras' Theorem is used to check whether the */
5386     /*   angle at the vertex is greater than 90 degrees.               */
5387     if (eapex[0] * (eorg[0] + edest[0]) +
5388         eapex[1] * (eorg[1] + edest[1]) >
5389         eapex[0] * eapex[0] + eorg[0] * edest[0] +
5390         eapex[1] * eapex[1] + eorg[1] * edest[1]) {
5391       addtolist += 2;
5392     }
5393   }
5394 
5395   if (addtolist && (!nobisect || ((nobisect == 1) && (sides == 2)))) {
5396     if (verbose > 2) {
5397       printf("  Queueing encroached segment (%.12g, %.12g) (%.12g, %.12g).\n",
5398              eorg[0], eorg[1], edest[0], edest[1]);
5399     }
5400     /* Add the shell edge to the list of encroached segments. */
5401     /*   Be sure to get the orientation right.                */
5402     badedge = (struct edge *) poolalloc(&badsegments);
5403     if (addtolist == 1) {
5404       shellecopy(*testedge, *badedge);
5405     } else {
5406       shellecopy(testsym, *badedge);
5407     }
5408   }
5409   return addtolist;
5410 }
5411 
5412 #endif /* not CDT_ONLY */
5413 
5414 /*****************************************************************************/
5415 /*                                                                           */
5416 /*  testtriangle()   Test a face for quality measures.                       */
5417 /*                                                                           */
5418 /*  Tests a triangle to see if it satisfies the minimum angle condition and  */
5419 /*  the maximum area condition.  Triangles that aren't up to spec are added  */
5420 /*  to the bad triangle queue.                                               */
5421 /*                                                                           */
5422 /*****************************************************************************/
5423 
5424 #ifndef CDT_ONLY
5425 
5426 void testtriangle(testtri)
5427 struct triedge *testtri;
5428 {
5429   struct triedge sametesttri;
5430   struct edge edge1, edge2;
5431   point torg, tdest, tapex;
5432   point anglevertex;
5433   REAL dxod, dyod, dxda, dyda, dxao, dyao;
5434   REAL dxod2, dyod2, dxda2, dyda2, dxao2, dyao2;
5435   REAL apexlen, orglen, destlen;
5436   REAL angle;
5437   REAL area;
5438   shelle sptr;                      /* Temporary variable used by tspivot(). */
5439 
5440   org(*testtri, torg);
5441   dest(*testtri, tdest);
5442   apex(*testtri, tapex);
5443   dxod = torg[0] - tdest[0];
5444   dyod = torg[1] - tdest[1];
5445   dxda = tdest[0] - tapex[0];
5446   dyda = tdest[1] - tapex[1];
5447   dxao = tapex[0] - torg[0];
5448   dyao = tapex[1] - torg[1];
5449   dxod2 = dxod * dxod;
5450   dyod2 = dyod * dyod;
5451   dxda2 = dxda * dxda;
5452   dyda2 = dyda * dyda;
5453   dxao2 = dxao * dxao;
5454   dyao2 = dyao * dyao;
5455   /* Find the lengths of the triangle's three edges. */
5456   apexlen = dxod2 + dyod2;
5457   orglen = dxda2 + dyda2;
5458   destlen = dxao2 + dyao2;
5459   if ((apexlen < orglen) && (apexlen < destlen)) {
5460     /* The edge opposite the apex is shortest. */
5461     /* Find the square of the cosine of the angle at the apex. */
5462     angle = dxda * dxao + dyda * dyao;
5463     angle = angle * angle / (orglen * destlen);
5464     anglevertex = tapex;
5465     lnext(*testtri, sametesttri);
5466     tspivot(sametesttri, edge1);
5467     lnextself(sametesttri);
5468     tspivot(sametesttri, edge2);
5469   } else if (orglen < destlen) {
5470     /* The edge opposite the origin is shortest. */
5471     /* Find the square of the cosine of the angle at the origin. */
5472     angle = dxod * dxao + dyod * dyao;
5473     angle = angle * angle / (apexlen * destlen);
5474     anglevertex = torg;
5475     tspivot(*testtri, edge1);
5476     lprev(*testtri, sametesttri);
5477     tspivot(sametesttri, edge2);
5478   } else {
5479     /* The edge opposite the destination is shortest. */
5480     /* Find the square of the cosine of the angle at the destination. */
5481     angle = dxod * dxda + dyod * dyda;
5482     angle = angle * angle / (apexlen * orglen);
5483     anglevertex = tdest;
5484     tspivot(*testtri, edge1);
5485     lnext(*testtri, sametesttri);
5486     tspivot(sametesttri, edge2);
5487   }
5488   /* Check if both edges that form the angle are segments. */
5489   if ((edge1.sh != dummysh) && (edge2.sh != dummysh)) {
5490     /* The angle is a segment intersection. */
5491     if ((angle > 0.9924) && !quiet) {                  /* Roughly 5 degrees. */
5492       if (angle > 1.0) {
5493         /* Beware of a floating exception in acos(). */
5494         angle = 1.0;
5495       }
5496       /* Find the actual angle in degrees, for printing. */
5497       angle = acos(sqrt(angle)) * (180.0 / PI);
5498       printf(
5499       "Warning:  Small angle (%.4g degrees) between segments at point\n",
5500              angle);
5501       printf("  (%.12g, %.12g)\n", anglevertex[0], anglevertex[1]);
5502     }
5503     /* Don't add this bad triangle to the list; there's nothing that */
5504     /*   can be done about a small angle between two segments.       */
5505     angle = 0.0;
5506   }
5507   /* Check whether the angle is smaller than permitted. */
5508   if (angle > goodangle) {
5509     /* Add this triangle to the list of bad triangles. */
5510     enqueuebadtri(testtri, angle, tapex, torg, tdest);
5511     return;
5512   }
5513   if (vararea || fixedarea) {
5514     /* Check whether the area is larger than permitted. */
5515     area = 0.5 * (dxod * dyda - dyod * dxda);
5516     if (fixedarea && (area > maxarea)) {
5517       /* Add this triangle to the list of bad triangles. */
5518       enqueuebadtri(testtri, angle, tapex, torg, tdest);
5519     } else if (vararea) {
5520       /* Nonpositive area constraints are treated as unconstrained. */
5521       if ((area > areabound(*testtri)) && (areabound(*testtri) > 0.0)) {
5522         /* Add this triangle to the list of bad triangles. */
5523         enqueuebadtri(testtri, angle, tapex, torg, tdest);
5524       }
5525     }
5526   }
5527 }
5528 
5529 #endif /* not CDT_ONLY */
5530 
5531 /**                                                                         **/
5532 /**                                                                         **/
5533 /********* Mesh quality testing routines end here                    *********/
5534 
5535 /********* Point location routines begin here                        *********/
5536 /**                                                                         **/
5537 /**                                                                         **/
5538 
5539 /*****************************************************************************/
5540 /*                                                                           */
5541 /*  makepointmap()   Construct a mapping from points to triangles to improve  */
5542 /*                  the speed of point location for segment insertion.       */
5543 /*                                                                           */
5544 /*  Traverses all the triangles, and provides each corner of each triangle   */
5545 /*  with a pointer to that triangle.  Of course, pointers will be            */
5546 /*  overwritten by other pointers because (almost) each point is a corner    */
5547 /*  of several triangles, but in the end every point will point to some      */
5548 /*  triangle that contains it.                                               */
5549 /*                                                                           */
5550 /*****************************************************************************/
5551 
5552 void makepointmap()
5553 {
5554   struct triedge triangleloop;
5555   point triorg;
5556 
5557   if (verbose) {
5558     printf("    Constructing mapping from points to triangles.\n");
5559   }
5560   traversalinit(&triangles);
5561   triangleloop.tri = triangletraverse();
5562   while (triangleloop.tri != (triangle *) NULL) {
5563     /* Check all three points of the triangle. */
5564     for (triangleloop.orient = 0; triangleloop.orient < 3;
5565          triangleloop.orient++) {
5566       org(triangleloop, triorg);
5567       setpoint2tri(triorg, encode(triangleloop));
5568     }
5569     triangleloop.tri = triangletraverse();
5570   }
5571 }
5572 
5573 /*****************************************************************************/
5574 /*                                                                           */
5575 /*  preciselocate()   Find a triangle or edge containing a given point.      */
5576 /*                                                                           */
5577 /*  Begins its search from `searchtri'.  It is important that `searchtri'    */
5578 /*  be a handle with the property that `searchpoint' is strictly to the left */
5579 /*  of the edge denoted by `searchtri', or is collinear with that edge and   */
5580 /*  does not intersect that edge.  (In particular, `searchpoint' should not  */
5581 /*  be the origin or destination of that edge.)                              */
5582 /*                                                                           */
5583 /*  These conditions are imposed because preciselocate() is normally used in */
5584 /*  one of two situations:                                                   */
5585 /*                                                                           */
5586 /*  (1)  To try to find the location to insert a new point.  Normally, we    */
5587 /*       know an edge that the point is strictly to the left of.  In the     */
5588 /*       incremental Delaunay algorithm, that edge is a bounding box edge.   */
5589 /*       In Ruppert's Delaunay refinement algorithm for quality meshing,     */
5590 /*       that edge is the shortest edge of the triangle whose circumcenter   */
5591 /*       is being inserted.                                                  */
5592 /*                                                                           */
5593 /*  (2)  To try to find an existing point.  In this case, any edge on the    */
5594 /*       convex hull is a good starting edge.  The possibility that the      */
5595 /*       vertex one seeks is an endpoint of the starting edge must be        */
5596 /*       screened out before preciselocate() is called.                      */
5597 /*                                                                           */
5598 /*  On completion, `searchtri' is a triangle that contains `searchpoint'.    */
5599 /*                                                                           */
5600 /*  This implementation differs from that given by Guibas and Stolfi.  It    */
5601 /*  walks from triangle to triangle, crossing an edge only if `searchpoint'  */
5602 /*  is on the other side of the line containing that edge.  After entering   */
5603 /*  a triangle, there are two edges by which one can leave that triangle.    */
5604 /*  If both edges are valid (`searchpoint' is on the other side of both      */
5605 /*  edges), one of the two is chosen by drawing a line perpendicular to      */
5606 /*  the entry edge (whose endpoints are `forg' and `fdest') passing through  */
5607 /*  `fapex'.  Depending on which side of this perpendicular `searchpoint'    */
5608 /*  falls on, an exit edge is chosen.                                        */
5609 /*                                                                           */
5610 /*  This implementation is empirically faster than the Guibas and Stolfi     */
5611 /*  point location routine (which I originally used), which tends to spiral  */
5612 /*  in toward its target.                                                    */
5613 /*                                                                           */
5614 /*  Returns ONVERTEX if the point lies on an existing vertex.  `searchtri'   */
5615 /*  is a handle whose origin is the existing vertex.                         */
5616 /*                                                                           */
5617 /*  Returns ONEDGE if the point lies on a mesh edge.  `searchtri' is a       */
5618 /*  handle whose primary edge is the edge on which the point lies.           */
5619 /*                                                                           */
5620 /*  Returns INTRIANGLE if the point lies strictly within a triangle.         */
5621 /*  `searchtri' is a handle on the triangle that contains the point.         */
5622 /*                                                                           */
5623 /*  Returns OUTSIDE if the point lies outside the mesh.  `searchtri' is a    */
5624 /*  handle whose primary edge the point is to the right of.  This might      */
5625 /*  occur when the circumcenter of a triangle falls just slightly outside    */
5626 /*  the mesh due to floating-point roundoff error.  It also occurs when      */
5627 /*  seeking a hole or region point that a foolish user has placed outside    */
5628 /*  the mesh.                                                                */
5629 /*                                                                           */
5630 /*  WARNING:  This routine is designed for convex triangulations, and will   */
5631 /*  not generally work after the holes and concavities have been carved.     */
5632 /*  However, it can still be used to find the circumcenter of a triangle, as */
5633 /*  long as the search is begun from the triangle in question.               */
5634 /*                                                                           */
5635 /*****************************************************************************/
5636 
5637 enum locateresult preciselocate(searchpoint, searchtri)
5638 point searchpoint;
5639 struct triedge *searchtri;
5640 {
5641   struct triedge backtracktri;
5642   point forg, fdest, fapex;
5643   point swappoint;
5644   REAL orgorient, destorient;
5645   int moveleft;
5646   triangle ptr;                         /* Temporary variable used by sym(). */
5647 
5648   if (verbose > 2) {
5649     printf("  Searching for point (%.12g, %.12g).\n",
5650            searchpoint[0], searchpoint[1]);
5651   }
5652   /* Where are we? */
5653   org(*searchtri, forg);
5654   dest(*searchtri, fdest);
5655   apex(*searchtri, fapex);
5656   while (1) {
5657     if (verbose > 2) {
5658       printf("    At (%.12g, %.12g) (%.12g, %.12g) (%.12g, %.12g)\n",
5659              forg[0], forg[1], fdest[0], fdest[1], fapex[0], fapex[1]);
5660     }
5661     /* Check whether the apex is the point we seek. */
5662     if ((fapex[0] == searchpoint[0]) && (fapex[1] == searchpoint[1])) {
5663       lprevself(*searchtri);
5664       return ONVERTEX;
5665     }
5666     /* Does the point lie on the other side of the line defined by the */
5667     /*   triangle edge opposite the triangle's destination?            */
5668     destorient = counterclockwise(forg, fapex, searchpoint);
5669     /* Does the point lie on the other side of the line defined by the */
5670     /*   triangle edge opposite the triangle's origin?                 */
5671     orgorient = counterclockwise(fapex, fdest, searchpoint);
5672     if (destorient > 0.0) {
5673       if (orgorient > 0.0) {
5674         /* Move left if the inner product of (fapex - searchpoint) and  */
5675         /*   (fdest - forg) is positive.  This is equivalent to drawing */
5676         /*   a line perpendicular to the line (forg, fdest) passing     */
5677         /*   through `fapex', and determining which side of this line   */
5678         /*   `searchpoint' falls on.                                    */
5679         moveleft = (fapex[0] - searchpoint[0]) * (fdest[0] - forg[0]) +
5680                    (fapex[1] - searchpoint[1]) * (fdest[1] - forg[1]) > 0.0;
5681       } else {
5682         moveleft = 1;
5683       }
5684     } else {
5685       if (orgorient > 0.0) {
5686         moveleft = 0;
5687       } else {
5688         /* The point we seek must be on the boundary of or inside this */
5689         /*   triangle.                                                 */
5690         if (destorient == 0.0) {
5691           lprevself(*searchtri);
5692           return ONEDGE;
5693         }
5694         if (orgorient == 0.0) {
5695           lnextself(*searchtri);
5696           return ONEDGE;
5697         }
5698         return INTRIANGLE;
5699       }
5700     }
5701 
5702     /* Move to another triangle.  Leave a trace `backtracktri' in case */
5703     /*   floating-point roundoff or some such bogey causes us to walk  */
5704     /*   off a boundary of the triangulation.  We can just bounce off  */
5705     /*   the boundary as if it were an elastic band.                   */
5706     if (moveleft) {
5707       lprev(*searchtri, backtracktri);
5708       fdest = fapex;
5709     } else {
5710       lnext(*searchtri, backtracktri);
5711       forg = fapex;
5712     }
5713     sym(backtracktri, *searchtri);
5714 
5715     /* Check for walking off the edge. */
5716     if (searchtri->tri == dummytri) {
5717       /* Turn around. */
5718       triedgecopy(backtracktri, *searchtri);
5719       swappoint = forg;
5720       forg = fdest;
5721       fdest = swappoint;
5722       apex(*searchtri, fapex);
5723       /* Check if the point really is beyond the triangulation boundary. */
5724       destorient = counterclockwise(forg, fapex, searchpoint);
5725       orgorient = counterclockwise(fapex, fdest, searchpoint);
5726       if ((orgorient < 0.0) && (destorient < 0.0)) {
5727         return OUTSIDE;
5728       }
5729     } else {
5730       apex(*searchtri, fapex);
5731     }
5732   }
5733 }
5734 
5735 /*****************************************************************************/
5736 /*                                                                           */
5737 /*  locate()   Find a triangle or edge containing a given point.             */
5738 /*                                                                           */
5739 /*  Searching begins from one of:  the input `searchtri', a recently         */
5740 /*  encountered triangle `recenttri', or from a triangle chosen from a       */
5741 /*  random sample.  The choice is made by determining which triangle's       */
5742 /*  origin is closest to the point we are searcing for.  Normally,           */
5743 /*  `searchtri' should be a handle on the convex hull of the triangulation.  */
5744 /*                                                                           */
5745 /*  Details on the random sampling method can be found in the Mucke, Saias,  */
5746 /*  and Zhu paper cited in the header of this code.                          */
5747 /*                                                                           */
5748 /*  On completion, `searchtri' is a triangle that contains `searchpoint'.    */
5749 /*                                                                           */
5750 /*  Returns ONVERTEX if the point lies on an existing vertex.  `searchtri'   */
5751 /*  is a handle whose origin is the existing vertex.                         */
5752 /*                                                                           */
5753 /*  Returns ONEDGE if the point lies on a mesh edge.  `searchtri' is a       */
5754 /*  handle whose primary edge is the edge on which the point lies.           */
5755 /*                                                                           */
5756 /*  Returns INTRIANGLE if the point lies strictly within a triangle.         */
5757 /*  `searchtri' is a handle on the triangle that contains the point.         */
5758 /*                                                                           */
5759 /*  Returns OUTSIDE if the point lies outside the mesh.  `searchtri' is a    */
5760 /*  handle whose primary edge the point is to the right of.  This might      */
5761 /*  occur when the circumcenter of a triangle falls just slightly outside    */
5762 /*  the mesh due to floating-point roundoff error.  It also occurs when      */
5763 /*  seeking a hole or region point that a foolish user has placed outside    */
5764 /*  the mesh.                                                                */
5765 /*                                                                           */
5766 /*  WARNING:  This routine is designed for convex triangulations, and will   */
5767 /*  not generally work after the holes and concavities have been carved.     */
5768 /*                                                                           */
5769 /*****************************************************************************/
5770 
5771 enum locateresult locate(searchpoint, searchtri)
5772 point searchpoint;
5773 struct triedge *searchtri;
5774 {
5775   VOID **sampleblock;
5776   triangle *firsttri;
5777   struct triedge sampletri;
5778   point torg, tdest;
5779   unsigned long alignptr;
5780   REAL searchdist, dist;
5781   REAL ahead;
5782   long sampleblocks, samplesperblock, samplenum;
5783   long triblocks;
5784   long i, j;
5785   triangle ptr;                         /* Temporary variable used by sym(). */
5786 
5787   if (verbose > 2) {
5788     printf("  Randomly sampling for a triangle near point (%.12g, %.12g).\n",
5789            searchpoint[0], searchpoint[1]);
5790   }
5791   /* Record the distance from the suggested starting triangle to the */
5792   /*   point we seek.                                                */
5793   org(*searchtri, torg);
5794   searchdist = (searchpoint[0] - torg[0]) * (searchpoint[0] - torg[0])
5795              + (searchpoint[1] - torg[1]) * (searchpoint[1] - torg[1]);
5796   if (verbose > 2) {
5797     printf("    Boundary triangle has origin (%.12g, %.12g).\n",
5798            torg[0], torg[1]);
5799   }
5800 
5801   /* If a recently encountered triangle has been recorded and has not been */
5802   /*   deallocated, test it as a good starting point.                      */
5803   if (recenttri.tri != (triangle *) NULL) {
5804     if (recenttri.tri[3] != (triangle) NULL) {
5805       org(recenttri, torg);
5806       if ((torg[0] == searchpoint[0]) && (torg[1] == searchpoint[1])) {
5807         triedgecopy(recenttri, *searchtri);
5808         return ONVERTEX;
5809       }
5810       dist = (searchpoint[0] - torg[0]) * (searchpoint[0] - torg[0])
5811            + (searchpoint[1] - torg[1]) * (searchpoint[1] - torg[1]);
5812       if (dist < searchdist) {
5813         triedgecopy(recenttri, *searchtri);
5814         searchdist = dist;
5815         if (verbose > 2) {
5816           printf("    Choosing recent triangle with origin (%.12g, %.12g).\n",
5817                  torg[0], torg[1]);
5818         }
5819       }
5820     }
5821   }
5822 
5823   /* The number of random samples taken is proportional to the cube root of */
5824   /*   the number of triangles in the mesh.  The next bit of code assumes   */
5825   /*   that the number of triangles increases monotonically.                */
5826   while (SAMPLEFACTOR * samples * samples * samples < triangles.items) {
5827     samples++;
5828   }
5829   triblocks = (triangles.maxitems + TRIPERBLOCK - 1) / TRIPERBLOCK;
5830   samplesperblock = 1 + (samples / triblocks);
5831   sampleblocks = samples / samplesperblock;
5832   sampleblock = triangles.firstblock;
5833   sampletri.orient = 0;
5834   for (i = 0; i < sampleblocks; i++) {
5835     alignptr = (unsigned long) (sampleblock + 1);
5836     firsttri = (triangle *) (alignptr + (unsigned long) triangles.alignbytes
5837                           - (alignptr % (unsigned long) triangles.alignbytes));
5838     for (j = 0; j < samplesperblock; j++) {
5839       if (i == triblocks - 1) {
5840         samplenum = randomnation((int)
5841                                  (triangles.maxitems - (i * TRIPERBLOCK)));
5842       } else {
5843         samplenum = randomnation(TRIPERBLOCK);
5844       }
5845       sampletri.tri = (triangle *)
5846                       (firsttri + (samplenum * triangles.itemwords));
5847       if (sampletri.tri[3] != (triangle) NULL) {
5848         org(sampletri, torg);
5849         dist = (searchpoint[0] - torg[0]) * (searchpoint[0] - torg[0])
5850              + (searchpoint[1] - torg[1]) * (searchpoint[1] - torg[1]);
5851         if (dist < searchdist) {
5852           triedgecopy(sampletri, *searchtri);
5853           searchdist = dist;
5854           if (verbose > 2) {
5855             printf("    Choosing triangle with origin (%.12g, %.12g).\n",
5856                    torg[0], torg[1]);
5857           }
5858         }
5859       }
5860     }
5861     sampleblock = (VOID **) *sampleblock;
5862   }
5863   /* Where are we? */
5864   org(*searchtri, torg);
5865   dest(*searchtri, tdest);
5866   /* Check the starting triangle's vertices. */
5867   if ((torg[0] == searchpoint[0]) && (torg[1] == searchpoint[1])) {
5868     return ONVERTEX;
5869   }
5870   if ((tdest[0] == searchpoint[0]) && (tdest[1] == searchpoint[1])) {
5871     lnextself(*searchtri);
5872     return ONVERTEX;
5873   }
5874   /* Orient `searchtri' to fit the preconditions of calling preciselocate(). */
5875   ahead = counterclockwise(torg, tdest, searchpoint);
5876   if (ahead < 0.0) {
5877     /* Turn around so that `searchpoint' is to the left of the */
5878     /*   edge specified by `searchtri'.                        */
5879     symself(*searchtri);
5880   } else if (ahead == 0.0) {
5881     /* Check if `searchpoint' is between `torg' and `tdest'. */
5882     if (((torg[0] < searchpoint[0]) == (searchpoint[0] < tdest[0]))
5883         && ((torg[1] < searchpoint[1]) == (searchpoint[1] < tdest[1]))) {
5884       return ONEDGE;
5885     }
5886   }
5887   return preciselocate(searchpoint, searchtri);
5888 }
5889 
5890 /**                                                                         **/
5891 /**                                                                         **/
5892 /********* Point location routines end here                          *********/
5893 
5894 /********* Mesh transformation routines begin here                   *********/
5895 /**                                                                         **/
5896 /**                                                                         **/
5897 
5898 /*****************************************************************************/
5899 /*                                                                           */
5900 /*  insertshelle()   Create a new shell edge and insert it between two       */
5901 /*                   triangles.                                              */
5902 /*                                                                           */
5903 /*  The new shell edge is inserted at the edge described by the handle       */
5904 /*  `tri'.  Its vertices are properly initialized.  The marker `shellemark'  */
5905 /*  is applied to the shell edge and, if appropriate, its vertices.          */
5906 /*                                                                           */
5907 /*****************************************************************************/
5908 
5909 void insertshelle(tri, shellemark)
5910 struct triedge *tri;          /* Edge at which to insert the new shell edge. */
5911 int shellemark;                            /* Marker for the new shell edge. */
5912 {
5913   struct triedge oppotri;
5914   struct edge newshelle;
5915   point triorg, tridest;
5916   triangle ptr;                         /* Temporary variable used by sym(). */
5917   shelle sptr;                      /* Temporary variable used by tspivot(). */
5918 
5919   /* Mark points if possible. */
5920   org(*tri, triorg);
5921   dest(*tri, tridest);
5922   if (pointmark(triorg) == 0) {
5923     setpointmark(triorg, shellemark);
5924   }
5925   if (pointmark(tridest) == 0) {
5926     setpointmark(tridest, shellemark);
5927   }
5928   /* Check if there's already a shell edge here. */
5929   tspivot(*tri, newshelle);
5930   if (newshelle.sh == dummysh) {
5931     /* Make new shell edge and initialize its vertices. */
5932     makeshelle(&newshelle);
5933     setsorg(newshelle, tridest);
5934     setsdest(newshelle, triorg);
5935     /* Bond new shell edge to the two triangles it is sandwiched between. */
5936     /*   Note that the facing triangle `oppotri' might be equal to        */
5937     /*   `dummytri' (outer space), but the new shell edge is bonded to it */
5938     /*   all the same.                                                    */
5939     tsbond(*tri, newshelle);
5940     sym(*tri, oppotri);
5941     ssymself(newshelle);
5942     tsbond(oppotri, newshelle);
5943     setmark(newshelle, shellemark);
5944     if (verbose > 2) {
5945       printf("  Inserting new ");
5946       printshelle(&newshelle);
5947     }
5948   } else {
5949     if (mark(newshelle) == 0) {
5950       setmark(newshelle, shellemark);
5951     }
5952   }
5953 }
5954 
5955 /*****************************************************************************/
5956 /*                                                                           */
5957 /*  Terminology                                                              */
5958 /*                                                                           */
5959 /*  A "local transformation" replaces a small set of triangles with another  */
5960 /*  set of triangles.  This may or may not involve inserting or deleting a   */
5961 /*  point.                                                                   */
5962 /*                                                                           */
5963 /*  The term "casing" is used to describe the set of triangles that are      */
5964 /*  attached to the triangles being transformed, but are not transformed     */
5965 /*  themselves.  Think of the casing as a fixed hollow structure inside      */
5966 /*  which all the action happens.  A "casing" is only defined relative to    */
5967 /*  a single transformation; each occurrence of a transformation will        */
5968 /*  involve a different casing.                                              */
5969 /*                                                                           */
5970 /*  A "shell" is similar to a "casing".  The term "shell" describes the set  */
5971 /*  of shell edges (if any) that are attached to the triangles being         */
5972 /*  transformed.  However, I sometimes use "shell" to refer to a single      */
5973 /*  shell edge, so don't get confused.                                       */
5974 /*                                                                           */
5975 /*****************************************************************************/
5976 
5977 /*****************************************************************************/
5978 /*                                                                           */
5979 /*  flip()   Transform two triangles to two different triangles by flipping  */
5980 /*           an edge within a quadrilateral.                                 */
5981 /*                                                                           */
5982 /*  Imagine the original triangles, abc and bad, oriented so that the        */
5983 /*  shared edge ab lies in a horizontal plane, with the point b on the left  */
5984 /*  and the point a on the right.  The point c lies below the edge, and the  */
5985 /*  point d lies above the edge.  The `flipedge' handle holds the edge ab    */
5986 /*  of triangle abc, and is directed left, from vertex a to vertex b.        */
5987 /*                                                                           */
5988 /*  The triangles abc and bad are deleted and replaced by the triangles cdb  */
5989 /*  and dca.  The triangles that represent abc and bad are NOT deallocated;  */
5990 /*  they are reused for dca and cdb, respectively.  Hence, any handles that  */
5991 /*  may have held the original triangles are still valid, although not       */
5992 /*  directed as they were before.                                            */
5993 /*                                                                           */
5994 /*  Upon completion of this routine, the `flipedge' handle holds the edge    */
5995 /*  dc of triangle dca, and is directed down, from vertex d to vertex c.     */
5996 /*  (Hence, the two triangles have rotated counterclockwise.)                */
5997 /*                                                                           */
5998 /*  WARNING:  This transformation is geometrically valid only if the         */
5999 /*  quadrilateral adbc is convex.  Furthermore, this transformation is       */
6000 /*  valid only if there is not a shell edge between the triangles abc and    */
6001 /*  bad.  This routine does not check either of these preconditions, and     */
6002 /*  it is the responsibility of the calling routine to ensure that they are  */
6003 /*  met.  If they are not, the streets shall be filled with wailing and      */
6004 /*  gnashing of teeth.                                                       */
6005 /*                                                                           */
6006 /*****************************************************************************/
6007 
6008 void flip(flipedge)
6009 struct triedge *flipedge;                    /* Handle for the triangle abc. */
6010 {
6011   struct triedge botleft, botright;
6012   struct triedge topleft, topright;
6013   struct triedge top;
6014   struct triedge botlcasing, botrcasing;
6015   struct triedge toplcasing, toprcasing;
6016   struct edge botlshelle, botrshelle;
6017   struct edge toplshelle, toprshelle;
6018   point leftpoint, rightpoint, botpoint;
6019   point farpoint;
6020   triangle ptr;                         /* Temporary variable used by sym(). */
6021   shelle sptr;                      /* Temporary variable used by tspivot(). */
6022 
6023   /* Identify the vertices of the quadrilateral. */
6024   org(*flipedge, rightpoint);
6025   dest(*flipedge, leftpoint);
6026   apex(*flipedge, botpoint);
6027   sym(*flipedge, top);
6028 #ifdef SELF_CHECK
6029   if (top.tri == dummytri) {
6030     printf("Internal error in flip():  Attempt to flip on boundary.\n");
6031     lnextself(*flipedge);
6032     return;
6033   }
6034   if (checksegments) {
6035     tspivot(*flipedge, toplshelle);
6036     if (toplshelle.sh != dummysh) {
6037       printf("Internal error in flip():  Attempt to flip a segment.\n");
6038       lnextself(*flipedge);
6039       return;
6040     }
6041   }
6042 #endif /* SELF_CHECK */
6043   apex(top, farpoint);
6044 
6045   /* Identify the casing of the quadrilateral. */
6046   lprev(top, topleft);
6047   sym(topleft, toplcasing);
6048   lnext(top, topright);
6049   sym(topright, toprcasing);
6050   lnext(*flipedge, botleft);
6051   sym(botleft, botlcasing);
6052   lprev(*flipedge, botright);
6053   sym(botright, botrcasing);
6054   /* Rotate the quadrilateral one-quarter turn counterclockwise. */
6055   bond(topleft, botlcasing);
6056   bond(botleft, botrcasing);
6057   bond(botright, toprcasing);
6058   bond(topright, toplcasing);
6059 
6060   if (checksegments) {
6061     /* Check for shell edges and rebond them to the quadrilateral. */
6062     tspivot(topleft, toplshelle);
6063     tspivot(botleft, botlshelle);
6064     tspivot(botright, botrshelle);
6065     tspivot(topright, toprshelle);
6066     if (toplshelle.sh == dummysh) {
6067       tsdissolve(topright);
6068     } else {
6069       tsbond(topright, toplshelle);
6070     }
6071     if (botlshelle.sh == dummysh) {
6072       tsdissolve(topleft);
6073     } else {
6074       tsbond(topleft, botlshelle);
6075     }
6076     if (botrshelle.sh == dummysh) {
6077       tsdissolve(botleft);
6078     } else {
6079       tsbond(botleft, botrshelle);
6080     }
6081     if (toprshelle.sh == dummysh) {
6082       tsdissolve(botright);
6083     } else {
6084       tsbond(botright, toprshelle);
6085     }
6086   }
6087 
6088   /* New point assignments for the rotated quadrilateral. */
6089   setorg(*flipedge, farpoint);
6090   setdest(*flipedge, botpoint);
6091   setapex(*flipedge, rightpoint);
6092   setorg(top, botpoint);
6093   setdest(top, farpoint);
6094   setapex(top, leftpoint);
6095   if (verbose > 2) {
6096     printf("  Edge flip results in left ");
6097     lnextself(topleft);
6098     printtriangle(&topleft);
6099     printf("  and right ");
6100     printtriangle(flipedge);
6101   }
6102 }
6103 
6104 /*****************************************************************************/
6105 /*                                                                           */
6106 /*  insertsite()   Insert a vertex into a Delaunay triangulation,            */
6107 /*                 performing flips as necessary to maintain the Delaunay    */
6108 /*                 property.                                                 */
6109 /*                                                                           */
6110 /*  The point `insertpoint' is located.  If `searchtri.tri' is not NULL,     */
6111 /*  the search for the containing triangle begins from `searchtri'.  If      */
6112 /*  `searchtri.tri' is NULL, a full point location procedure is called.      */
6113 /*  If `insertpoint' is found inside a triangle, the triangle is split into  */
6114 /*  three; if `insertpoint' lies on an edge, the edge is split in two,       */
6115 /*  thereby splitting the two adjacent triangles into four.  Edge flips are  */
6116 /*  used to restore the Delaunay property.  If `insertpoint' lies on an      */
6117 /*  existing vertex, no action is taken, and the value DUPLICATEPOINT is     */
6118 /*  returned.  On return, `searchtri' is set to a handle whose origin is the */
6119 /*  existing vertex.                                                         */
6120 /*                                                                           */
6121 /*  Normally, the parameter `splitedge' is set to NULL, implying that no     */
6122 /*  segment should be split.  In this case, if `insertpoint' is found to     */
6123 /*  lie on a segment, no action is taken, and the value VIOLATINGPOINT is    */
6124 /*  returned.  On return, `searchtri' is set to a handle whose primary edge  */
6125 /*  is the violated segment.                                                 */
6126 /*                                                                           */
6127 /*  If the calling routine wishes to split a segment by inserting a point in */
6128 /*  it, the parameter `splitedge' should be that segment.  In this case,     */
6129 /*  `searchtri' MUST be the triangle handle reached by pivoting from that    */
6130 /*  segment; no point location is done.                                      */
6131 /*                                                                           */
6132 /*  `segmentflaws' and `triflaws' are flags that indicate whether or not     */
6133 /*  there should be checks for the creation of encroached segments or bad    */
6134 /*  quality faces.  If a newly inserted point encroaches upon segments,      */
6135 /*  these segments are added to the list of segments to be split if          */
6136 /*  `segmentflaws' is set.  If bad triangles are created, these are added    */
6137 /*  to the queue if `triflaws' is set.                                       */
6138 /*                                                                           */
6139 /*  If a duplicate point or violated segment does not prevent the point      */
6140 /*  from being inserted, the return value will be ENCROACHINGPOINT if the    */
6141 /*  point encroaches upon a segment (and checking is enabled), or            */
6142 /*  SUCCESSFULPOINT otherwise.  In either case, `searchtri' is set to a      */
6143 /*  handle whose origin is the newly inserted vertex.                        */
6144 /*                                                                           */
6145 /*  insertsite() does not use flip() for reasons of speed; some              */
6146 /*  information can be reused from edge flip to edge flip, like the          */
6147 /*  locations of shell edges.                                                */
6148 /*                                                                           */
6149 /*****************************************************************************/
6150 
6151 enum insertsiteresult insertsite(insertpoint, searchtri, splitedge,
6152                                  segmentflaws, triflaws)
6153 point insertpoint;
6154 struct triedge *searchtri;
6155 struct edge *splitedge;
6156 int segmentflaws;
6157 int triflaws;
6158 {
6159   struct triedge horiz;
6160   struct triedge top;
6161   struct triedge botleft, botright;
6162   struct triedge topleft, topright;
6163   struct triedge newbotleft, newbotright;
6164   struct triedge newtopright;
6165   struct triedge botlcasing, botrcasing;
6166   struct triedge toplcasing, toprcasing;
6167   struct triedge testtri;
6168   struct edge botlshelle, botrshelle;
6169   struct edge toplshelle, toprshelle;
6170   struct edge brokenshelle;
6171   struct edge checkshelle;
6172   struct edge rightedge;
6173   struct edge newedge;
6174   struct edge *encroached;
6175   point first;
6176   point leftpoint, rightpoint, botpoint, toppoint, farpoint;
6177   REAL attrib;
6178   REAL area;
6179   enum insertsiteresult success;
6180   enum locateresult intersect;
6181   int doflip;
6182   int mirrorflag;
6183   int i;
6184   triangle ptr;                         /* Temporary variable used by sym(). */
6185   shelle sptr;         /* Temporary variable used by spivot() and tspivot(). */
6186 
6187   if (verbose > 1) {
6188     printf("  Inserting (%.12g, %.12g).\n", insertpoint[0], insertpoint[1]);
6189   }
6190   if (splitedge == (struct edge *) NULL) {
6191     /* Find the location of the point to be inserted.  Check if a good */
6192     /*   starting triangle has already been provided by the caller.    */
6193     if (searchtri->tri == (triangle *) NULL) {
6194       /* Find a boundary triangle. */
6195       horiz.tri = dummytri;
6196       horiz.orient = 0;
6197       symself(horiz);
6198       /* Search for a triangle containing `insertpoint'. */
6199       intersect = locate(insertpoint, &horiz);
6200     } else {
6201       /* Start searching from the triangle provided by the caller. */
6202       triedgecopy(*searchtri, horiz);
6203       intersect = preciselocate(insertpoint, &horiz);
6204     }
6205   } else {
6206     /* The calling routine provides the edge in which the point is inserted. */
6207     triedgecopy(*searchtri, horiz);
6208     intersect = ONEDGE;
6209   }
6210   if (intersect == ONVERTEX) {
6211     /* There's already a vertex there.  Return in `searchtri' a triangle */
6212     /*   whose origin is the existing vertex.                            */
6213     triedgecopy(horiz, *searchtri);
6214     triedgecopy(horiz, recenttri);
6215     return DUPLICATEPOINT;
6216   }
6217   if ((intersect == ONEDGE) || (intersect == OUTSIDE)) {
6218     /* The vertex falls on an edge or boundary. */
6219     if (checksegments && (splitedge == (struct edge *) NULL)) {
6220       /* Check whether the vertex falls on a shell edge. */
6221       tspivot(horiz, brokenshelle);
6222       if (brokenshelle.sh != dummysh) {
6223         /* The vertex falls on a shell edge. */
6224         if (segmentflaws) {
6225           if (nobisect == 0) {
6226             /* Add the shell edge to the list of encroached segments. */
6227             encroached = (struct edge *) poolalloc(&badsegments);
6228             shellecopy(brokenshelle, *encroached);
6229           } else if ((nobisect == 1) && (intersect == ONEDGE)) {
6230             /* This segment may be split only if it is an internal boundary. */
6231             sym(horiz, testtri);
6232             if (testtri.tri != dummytri) {
6233               /* Add the shell edge to the list of encroached segments. */
6234               encroached = (struct edge *) poolalloc(&badsegments);
6235               shellecopy(brokenshelle, *encroached);
6236             }
6237           }
6238         }
6239         /* Return a handle whose primary edge contains the point, */
6240         /*   which has not been inserted.                         */
6241         triedgecopy(horiz, *searchtri);
6242         triedgecopy(horiz, recenttri);
6243         return VIOLATINGPOINT;
6244       }
6245     }
6246     /* Insert the point on an edge, dividing one triangle into two (if */
6247     /*   the edge lies on a boundary) or two triangles into four.      */
6248     lprev(horiz, botright);
6249     sym(botright, botrcasing);
6250     sym(horiz, topright);
6251     /* Is there a second triangle?  (Or does this edge lie on a boundary?) */
6252     mirrorflag = topright.tri != dummytri;
6253     if (mirrorflag) {
6254       lnextself(topright);
6255       sym(topright, toprcasing);
6256       maketriangle(&newtopright);
6257     } else {
6258       /* Splitting the boundary edge increases the number of boundary edges. */
6259       hullsize++;
6260     }
6261     maketriangle(&newbotright);
6262 
6263     /* Set the vertices of changed and new triangles. */
6264     org(horiz, rightpoint);
6265     dest(horiz, leftpoint);
6266     apex(horiz, botpoint);
6267     setorg(newbotright, botpoint);
6268     setdest(newbotright, rightpoint);
6269     setapex(newbotright, insertpoint);
6270     setorg(horiz, insertpoint);
6271     for (i = 0; i < eextras; i++) {
6272       /* Set the element attributes of a new triangle. */
6273       setelemattribute(newbotright, i, elemattribute(botright, i));
6274     }
6275     if (vararea) {
6276       /* Set the area constraint of a new triangle. */
6277       setareabound(newbotright, areabound(botright));
6278     }
6279     if (mirrorflag) {
6280       dest(topright, toppoint);
6281       setorg(newtopright, rightpoint);
6282       setdest(newtopright, toppoint);
6283       setapex(newtopright, insertpoint);
6284       setorg(topright, insertpoint);
6285       for (i = 0; i < eextras; i++) {
6286         /* Set the element attributes of another new triangle. */
6287         setelemattribute(newtopright, i, elemattribute(topright, i));
6288       }
6289       if (vararea) {
6290         /* Set the area constraint of another new triangle. */
6291         setareabound(newtopright, areabound(topright));
6292       }
6293     }
6294 
6295     /* There may be shell edges that need to be bonded */
6296     /*   to the new triangle(s).                       */
6297     if (checksegments) {
6298       tspivot(botright, botrshelle);
6299       if (botrshelle.sh != dummysh) {
6300         tsdissolve(botright);
6301         tsbond(newbotright, botrshelle);
6302       }
6303       if (mirrorflag) {
6304         tspivot(topright, toprshelle);
6305         if (toprshelle.sh != dummysh) {
6306           tsdissolve(topright);
6307           tsbond(newtopright, toprshelle);
6308         }
6309       }
6310     }
6311 
6312     /* Bond the new triangle(s) to the surrounding triangles. */
6313     bond(newbotright, botrcasing);
6314     lprevself(newbotright);
6315     bond(newbotright, botright);
6316     lprevself(newbotright);
6317     if (mirrorflag) {
6318       bond(newtopright, toprcasing);
6319       lnextself(newtopright);
6320       bond(newtopright, topright);
6321       lnextself(newtopright);
6322       bond(newtopright, newbotright);
6323     }
6324 
6325     if (splitedge != (struct edge *) NULL) {
6326       /* Split the shell edge into two. */
6327       setsdest(*splitedge, insertpoint);
6328       ssymself(*splitedge);
6329       spivot(*splitedge, rightedge);
6330       insertshelle(&newbotright, mark(*splitedge));
6331       tspivot(newbotright, newedge);
6332       sbond(*splitedge, newedge);
6333       ssymself(newedge);
6334       sbond(newedge, rightedge);
6335       ssymself(*splitedge);
6336     }
6337 
6338 #ifdef SELF_CHECK
6339     if (counterclockwise(rightpoint, leftpoint, botpoint) < 0.0) {
6340       printf("Internal error in insertsite():\n");
6341       printf("  Clockwise triangle prior to edge point insertion (bottom).\n");
6342     }
6343     if (mirrorflag) {
6344       if (counterclockwise(leftpoint, rightpoint, toppoint) < 0.0) {
6345         printf("Internal error in insertsite():\n");
6346         printf("  Clockwise triangle prior to edge point insertion (top).\n");
6347       }
6348       if (counterclockwise(rightpoint, toppoint, insertpoint) < 0.0) {
6349         printf("Internal error in insertsite():\n");
6350         printf("  Clockwise triangle after edge point insertion (top right).\n"
6351                );
6352       }
6353       if (counterclockwise(toppoint, leftpoint, insertpoint) < 0.0) {
6354         printf("Internal error in insertsite():\n");
6355         printf("  Clockwise triangle after edge point insertion (top left).\n"
6356                );
6357       }
6358     }
6359     if (counterclockwise(leftpoint, botpoint, insertpoint) < 0.0) {
6360       printf("Internal error in insertsite():\n");
6361       printf("  Clockwise triangle after edge point insertion (bottom left).\n"
6362              );
6363     }
6364     if (counterclockwise(botpoint, rightpoint, insertpoint) < 0.0) {
6365       printf("Internal error in insertsite():\n");
6366       printf(
6367         "  Clockwise triangle after edge point insertion (bottom right).\n");
6368     }
6369 #endif /* SELF_CHECK */
6370     if (verbose > 2) {
6371       printf("  Updating bottom left ");
6372       printtriangle(&botright);
6373       if (mirrorflag) {
6374         printf("  Updating top left ");
6375         printtriangle(&topright);
6376         printf("  Creating top right ");
6377         printtriangle(&newtopright);
6378       }
6379       printf("  Creating bottom right ");
6380       printtriangle(&newbotright);
6381     }
6382 
6383     /* Position `horiz' on the first edge to check for */
6384     /*   the Delaunay property.                        */
6385     lnextself(horiz);
6386   } else {
6387     /* Insert the point in a triangle, splitting it into three. */
6388     lnext(horiz, botleft);
6389     lprev(horiz, botright);
6390     sym(botleft, botlcasing);
6391     sym(botright, botrcasing);
6392     maketriangle(&newbotleft);
6393     maketriangle(&newbotright);
6394 
6395     /* Set the vertices of changed and new triangles. */
6396     org(horiz, rightpoint);
6397     dest(horiz, leftpoint);
6398     apex(horiz, botpoint);
6399     setorg(newbotleft, leftpoint);
6400     setdest(newbotleft, botpoint);
6401     setapex(newbotleft, insertpoint);
6402     setorg(newbotright, botpoint);
6403     setdest(newbotright, rightpoint);
6404     setapex(newbotright, insertpoint);
6405     setapex(horiz, insertpoint);
6406     for (i = 0; i < eextras; i++) {
6407       /* Set the element attributes of the new triangles. */
6408       attrib = elemattribute(horiz, i);
6409       setelemattribute(newbotleft, i, attrib);
6410       setelemattribute(newbotright, i, attrib);
6411     }
6412     if (vararea) {
6413       /* Set the area constraint of the new triangles. */
6414       area = areabound(horiz);
6415       setareabound(newbotleft, area);
6416       setareabound(newbotright, area);
6417     }
6418 
6419     /* There may be shell edges that need to be bonded */
6420     /*   to the new triangles.                         */
6421     if (checksegments) {
6422       tspivot(botleft, botlshelle);
6423       if (botlshelle.sh != dummysh) {
6424         tsdissolve(botleft);
6425         tsbond(newbotleft, botlshelle);
6426       }
6427       tspivot(botright, botrshelle);
6428       if (botrshelle.sh != dummysh) {
6429         tsdissolve(botright);
6430         tsbond(newbotright, botrshelle);
6431       }
6432     }
6433 
6434     /* Bond the new triangles to the surrounding triangles. */
6435     bond(newbotleft, botlcasing);
6436     bond(newbotright, botrcasing);
6437     lnextself(newbotleft);
6438     lprevself(newbotright);
6439     bond(newbotleft, newbotright);
6440     lnextself(newbotleft);
6441     bond(botleft, newbotleft);
6442     lprevself(newbotright);
6443     bond(botright, newbotright);
6444 
6445 #ifdef SELF_CHECK
6446     if (counterclockwise(rightpoint, leftpoint, botpoint) < 0.0) {
6447       printf("Internal error in insertsite():\n");
6448       printf("  Clockwise triangle prior to point insertion.\n");
6449     }
6450     if (counterclockwise(rightpoint, leftpoint, insertpoint) < 0.0) {
6451       printf("Internal error in insertsite():\n");
6452       printf("  Clockwise triangle after point insertion (top).\n");
6453     }
6454     if (counterclockwise(leftpoint, botpoint, insertpoint) < 0.0) {
6455       printf("Internal error in insertsite():\n");
6456       printf("  Clockwise triangle after point insertion (left).\n");
6457     }
6458     if (counterclockwise(botpoint, rightpoint, insertpoint) < 0.0) {
6459       printf("Internal error in insertsite():\n");
6460       printf("  Clockwise triangle after point insertion (right).\n");
6461     }
6462 #endif /* SELF_CHECK */
6463     if (verbose > 2) {
6464       printf("  Updating top ");
6465       printtriangle(&horiz);
6466       printf("  Creating left ");
6467       printtriangle(&newbotleft);
6468       printf("  Creating right ");
6469       printtriangle(&newbotright);
6470     }
6471   }
6472 
6473   /* The insertion is successful by default, unless an encroached */
6474   /*   edge is found.                                             */
6475   success = SUCCESSFULPOINT;
6476   /* Circle around the newly inserted vertex, checking each edge opposite */
6477   /*   it for the Delaunay property.  Non-Delaunay edges are flipped.     */
6478   /*   `horiz' is always the edge being checked.  `first' marks where to  */
6479   /*   stop circling.                                                     */
6480   org(horiz, first);
6481   rightpoint = first;
6482   dest(horiz, leftpoint);
6483   /* Circle until finished. */
6484   while (1) {
6485     /* By default, the edge will be flipped. */
6486     doflip = 1;
6487     if (checksegments) {
6488       /* Check for a segment, which cannot be flipped. */
6489       tspivot(horiz, checkshelle);
6490       if (checkshelle.sh != dummysh) {
6491         /* The edge is a segment and cannot be flipped. */
6492         doflip = 0;
6493 #ifndef CDT_ONLY
6494         if (segmentflaws) {
6495           /* Does the new point encroach upon this segment? */
6496           if (checkedge4encroach(&checkshelle)) {
6497             success = ENCROACHINGPOINT;
6498           }
6499         }
6500 #endif /* not CDT_ONLY */
6501       }
6502     }
6503     if (doflip) {
6504       /* Check if the edge is a boundary edge. */
6505       sym(horiz, top);
6506       if (top.tri == dummytri) {
6507         /* The edge is a boundary edge and cannot be flipped. */
6508         doflip = 0;
6509       } else {
6510         /* Find the point on the other side of the edge. */
6511         apex(top, farpoint);
6512         /* In the incremental Delaunay triangulation algorithm, any of    */
6513         /*   `leftpoint', `rightpoint', and `farpoint' could be vertices  */
6514         /*   of the triangular bounding box.  These vertices must be      */
6515         /*   treated as if they are infinitely distant, even though their */
6516         /*   "coordinates" are not.                                       */
6517         if ((leftpoint == infpoint1) || (leftpoint == infpoint2)
6518                    || (leftpoint == infpoint3)) {
6519           /* `leftpoint' is infinitely distant.  Check the convexity of */
6520           /*   the boundary of the triangulation.  'farpoint' might be  */
6521           /*   infinite as well, but trust me, this same condition      */
6522           /*   should be applied.                                       */
6523           doflip = counterclockwise(insertpoint, rightpoint, farpoint) > 0.0;
6524         } else if ((rightpoint == infpoint1) || (rightpoint == infpoint2)
6525                    || (rightpoint == infpoint3)) {
6526           /* `rightpoint' is infinitely distant.  Check the convexity of */
6527           /*   the boundary of the triangulation.  'farpoint' might be  */
6528           /*   infinite as well, but trust me, this same condition      */
6529           /*   should be applied.                                       */
6530           doflip = counterclockwise(farpoint, leftpoint, insertpoint) > 0.0;
6531         } else if ((farpoint == infpoint1) || (farpoint == infpoint2)
6532             || (farpoint == infpoint3)) {
6533           /* `farpoint' is infinitely distant and cannot be inside */
6534           /*   the circumcircle of the triangle `horiz'.           */
6535           doflip = 0;
6536         } else {
6537           /* Test whether the edge is locally Delaunay. */
6538           doflip = incircle(leftpoint, insertpoint, rightpoint, farpoint)
6539                    > 0.0;
6540         }
6541         if (doflip) {
6542           /* We made it!  Flip the edge `horiz' by rotating its containing */
6543           /*   quadrilateral (the two triangles adjacent to `horiz').      */
6544           /* Identify the casing of the quadrilateral. */
6545           lprev(top, topleft);
6546           sym(topleft, toplcasing);
6547           lnext(top, topright);
6548           sym(topright, toprcasing);
6549           lnext(horiz, botleft);
6550           sym(botleft, botlcasing);
6551           lprev(horiz, botright);
6552           sym(botright, botrcasing);
6553           /* Rotate the quadrilateral one-quarter turn counterclockwise. */
6554           bond(topleft, botlcasing);
6555           bond(botleft, botrcasing);
6556           bond(botright, toprcasing);
6557           bond(topright, toplcasing);
6558           if (checksegments) {
6559             /* Check for shell edges and rebond them to the quadrilateral. */
6560             tspivot(topleft, toplshelle);
6561             tspivot(botleft, botlshelle);
6562             tspivot(botright, botrshelle);
6563             tspivot(topright, toprshelle);
6564             if (toplshelle.sh == dummysh) {
6565               tsdissolve(topright);
6566             } else {
6567               tsbond(topright, toplshelle);
6568             }
6569             if (botlshelle.sh == dummysh) {
6570               tsdissolve(topleft);
6571             } else {
6572               tsbond(topleft, botlshelle);
6573             }
6574             if (botrshelle.sh == dummysh) {
6575               tsdissolve(botleft);
6576             } else {
6577               tsbond(botleft, botrshelle);
6578             }
6579             if (toprshelle.sh == dummysh) {
6580               tsdissolve(botright);
6581             } else {
6582               tsbond(botright, toprshelle);
6583             }
6584           }
6585           /* New point assignments for the rotated quadrilateral. */
6586           setorg(horiz, farpoint);
6587           setdest(horiz, insertpoint);
6588           setapex(horiz, rightpoint);
6589           setorg(top, insertpoint);
6590           setdest(top, farpoint);
6591           setapex(top, leftpoint);
6592           for (i = 0; i < eextras; i++) {
6593             /* Take the average of the two triangles' attributes. */
6594             attrib = 0.5 * (elemattribute(top, i) + elemattribute(horiz, i));
6595             setelemattribute(top, i, attrib);
6596             setelemattribute(horiz, i, attrib);
6597           }
6598           if (vararea) {
6599             if ((areabound(top) <= 0.0) || (areabound(horiz) <= 0.0)) {
6600               area = -1.0;
6601             } else {
6602               /* Take the average of the two triangles' area constraints.    */
6603               /*   This prevents small area constraints from migrating a     */
6604               /*   long, long way from their original location due to flips. */
6605               area = 0.5 * (areabound(top) + areabound(horiz));
6606             }
6607             setareabound(top, area);
6608             setareabound(horiz, area);
6609           }
6610 #ifdef SELF_CHECK
6611           if (insertpoint != (point) NULL) {
6612             if (counterclockwise(leftpoint, insertpoint, rightpoint) < 0.0) {
6613               printf("Internal error in insertsite():\n");
6614               printf("  Clockwise triangle prior to edge flip (bottom).\n");
6615             }
6616             /* The following test has been removed because constrainededge() */
6617             /*   sometimes generates inverted triangles that insertsite()    */
6618             /*   removes.                                                    */
6619 /*
6620             if (counterclockwise(rightpoint, farpoint, leftpoint) < 0.0) {
6621               printf("Internal error in insertsite():\n");
6622               printf("  Clockwise triangle prior to edge flip (top).\n");
6623             }
6624 */
6625             if (counterclockwise(farpoint, leftpoint, insertpoint) < 0.0) {
6626               printf("Internal error in insertsite():\n");
6627               printf("  Clockwise triangle after edge flip (left).\n");
6628             }
6629             if (counterclockwise(insertpoint, rightpoint, farpoint) < 0.0) {
6630               printf("Internal error in insertsite():\n");
6631               printf("  Clockwise triangle after edge flip (right).\n");
6632             }
6633           }
6634 #endif /* SELF_CHECK */
6635           if (verbose > 2) {
6636             printf("  Edge flip results in left ");
6637             lnextself(topleft);
6638             printtriangle(&topleft);
6639             printf("  and right ");
6640             printtriangle(&horiz);
6641           }
6642           /* On the next iterations, consider the two edges that were  */
6643           /*   exposed (this is, are now visible to the newly inserted */
6644           /*   point) by the edge flip.                                */
6645           lprevself(horiz);
6646           leftpoint = farpoint;
6647         }
6648       }
6649     }
6650     if (!doflip) {
6651       /* The handle `horiz' is accepted as locally Delaunay. */
6652 #ifndef CDT_ONLY
6653       if (triflaws) {
6654         /* Check the triangle `horiz' for quality. */
6655         testtriangle(&horiz);
6656       }
6657 #endif /* not CDT_ONLY */
6658       /* Look for the next edge around the newly inserted point. */
6659       lnextself(horiz);
6660       sym(horiz, testtri);
6661       /* Check for finishing a complete revolution about the new point, or */
6662       /*   falling off the edge of the triangulation.  The latter will     */
6663       /*   happen when a point is inserted at a boundary.                  */
6664       if ((leftpoint == first) || (testtri.tri == dummytri)) {
6665         /* We're done.  Return a triangle whose origin is the new point. */
6666         lnext(horiz, *searchtri);
6667         lnext(horiz, recenttri);
6668         return success;
6669       }
6670       /* Finish finding the next edge around the newly inserted point. */
6671       lnext(testtri, horiz);
6672       rightpoint = leftpoint;
6673       dest(horiz, leftpoint);
6674     }
6675   }
6676 }
6677 
6678 /*****************************************************************************/
6679 /*                                                                           */
6680 /*  triangulatepolygon()   Find the Delaunay triangulation of a polygon that */
6681 /*                         has a certain "nice" shape.  This includes the    */
6682 /*                         polygons that result from deletion of a point or  */
6683 /*                         insertion of a segment.                           */
6684 /*                                                                           */
6685 /*  This is a conceptually difficult routine.  The starting assumption is    */
6686 /*  that we have a polygon with n sides.  n - 1 of these sides are currently */
6687 /*  represented as edges in the mesh.  One side, called the "base", need not */
6688 /*  be.                                                                      */
6689 /*                                                                           */
6690 /*  Inside the polygon is a structure I call a "fan", consisting of n - 1    */
6691 /*  triangles that share a common origin.  For each of these triangles, the  */
6692 /*  edge opposite the origin is one of the sides of the polygon.  The        */
6693 /*  primary edge of each triangle is the edge directed from the origin to    */
6694 /*  the destination; note that this is not the same edge that is a side of   */
6695 /*  the polygon.  `firstedge' is the primary edge of the first triangle.     */
6696 /*  From there, the triangles follow in counterclockwise order about the     */
6697 /*  polygon, until `lastedge', the primary edge of the last triangle.        */
6698 /*  `firstedge' and `lastedge' are probably connected to other triangles     */
6699 /*  beyond the extremes of the fan, but their identity is not important, as  */
6700 /*  long as the fan remains connected to them.                               */
6701 /*                                                                           */
6702 /*  Imagine the polygon oriented so that its base is at the bottom.  This    */
6703 /*  puts `firstedge' on the far right, and `lastedge' on the far left.       */
6704 /*  The right vertex of the base is the destination of `firstedge', and the  */
6705 /*  left vertex of the base is the apex of `lastedge'.                       */
6706 /*                                                                           */
6707 /*  The challenge now is to find the right sequence of edge flips to         */
6708 /*  transform the fan into a Delaunay triangulation of the polygon.  Each    */
6709 /*  edge flip effectively removes one triangle from the fan, committing it   */
6710 /*  to the polygon.  The resulting polygon has one fewer edge.  If `doflip'  */
6711 /*  is set, the final flip will be performed, resulting in a fan of one      */
6712 /*  (useless?) triangle.  If `doflip' is not set, the final flip is not      */
6713 /*  performed, resulting in a fan of two triangles, and an unfinished        */
6714 /*  triangular polygon that is not yet filled out with a single triangle.    */
6715 /*  On completion of the routine, `lastedge' is the last remaining triangle, */
6716 /*  or the leftmost of the last two.                                         */
6717 /*                                                                           */
6718 /*  Although the flips are performed in the order described above, the       */
6719 /*  decisions about what flips to perform are made in precisely the reverse  */
6720 /*  order.  The recursive triangulatepolygon() procedure makes a decision,   */
6721 /*  uses up to two recursive calls to triangulate the "subproblems"          */
6722 /*  (polygons with fewer edges), and then performs an edge flip.             */
6723 /*                                                                           */
6724 /*  The "decision" it makes is which vertex of the polygon should be         */
6725 /*  connected to the base.  This decision is made by testing every possible  */
6726 /*  vertex.  Once the best vertex is found, the two edges that connect this  */
6727 /*  vertex to the base become the bases for two smaller polygons.  These     */
6728 /*  are triangulated recursively.  Unfortunately, this approach can take     */
6729 /*  O(n^2) time not only in the worst case, but in many common cases.  It's  */
6730 /*  rarely a big deal for point deletion, where n is rarely larger than ten, */
6731 /*  but it could be a big deal for segment insertion, especially if there's  */
6732 /*  a lot of long segments that each cut many triangles.  I ought to code    */
6733 /*  a faster algorithm some time.                                            */
6734 /*                                                                           */
6735 /*  The `edgecount' parameter is the number of sides of the polygon,         */
6736 /*  including its base.  `triflaws' is a flag that determines whether the    */
6737 /*  new triangles should be tested for quality, and enqueued if they are     */
6738 /*  bad.                                                                     */
6739 /*                                                                           */
6740 /*****************************************************************************/
6741 
6742 void triangulatepolygon(firstedge, lastedge, edge