Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/ElmerGUI/Application/plugins/tetgen.h
3203 views
1
///////////////////////////////////////////////////////////////////////////////
2
// //
3
// TetGen //
4
// //
5
// A Quality Tetrahedral Mesh Generator and A 3D Delaunay Triangulator //
6
// //
7
// Version 1.5 //
8
// May 31, 2014 //
9
// //
10
// Copyright (C) 2002--2014 //
11
// //
12
// TetGen is freely available through the website: http://www.tetgen.org. //
13
// It may be copied, modified, and redistributed for non-commercial use. //
14
// Please consult the file LICENSE for the detailed copyright notices. //
15
// //
16
///////////////////////////////////////////////////////////////////////////////
17
18
19
#ifndef tetgenH
20
#define tetgenH
21
22
// To compile TetGen as a library instead of an executable program, define
23
// the TETLIBRARY symbol.
24
25
// #define TETLIBRARY
26
27
28
// TetGen default uses the double precision (64 bit) for a real number.
29
// Alternatively, one can use the single precision (32 bit) 'float' if the
30
// memory is limited.
31
32
#define REAL double // #define REAL float
33
34
// Maximum number of characters in a file name (including the null).
35
36
#define FILENAMESIZE 1024
37
38
// Maximum number of chars in a line read from a file (including the null).
39
40
#define INPUTLINESIZE 2048
41
42
// TetGen only uses the C standard library.
43
44
#include <stdio.h>
45
#include <stdlib.h>
46
#include <string.h>
47
#include <math.h>
48
#include <time.h>
49
50
// The types 'intptr_t' and 'uintptr_t' are signed and unsigned integer types,
51
// respectively. They are guaranteed to be the same width as a pointer.
52
// They are defined in <stdint.h> by the C99 Standard. However, Microsoft
53
// Visual C++ 2003 -- 2008 (Visual C++ 7.1 - 9) doesn't ship with this header
54
// file. In such case, we can define them by ourself.
55
// Update (learned from Stack Overflow): Visual Studio 2010 and Visual C++ 2010
56
// Express both have stdint.h
57
58
// The following piece of code was provided by Steven Johnson (MIT). Define the
59
// symbol _MSC_VER if you are using Microsoft Visual C++. Moreover, define
60
// the _WIN64 symbol if you are running TetGen on Win64 systems.
61
62
#ifdef _MSC_VER // Microsoft Visual C++
63
# ifdef _WIN64
64
typedef __int64 intptr_t;
65
typedef unsigned __int64 uintptr_t;
66
# else // not _WIN64
67
typedef int intptr_t;
68
typedef unsigned int uintptr_t;
69
# endif
70
#else // not Visual C++
71
# include <stdint.h>
72
#endif
73
74
///////////////////////////////////////////////////////////////////////////////
75
// //
76
// tetgenio //
77
// //
78
// A structure for transferring data into and out of TetGen's mesh structure,//
79
// 'tetgenmesh' (declared below). //
80
// //
81
// The input of TetGen is either a 3D point set, or a 3D piecewise linear //
82
// complex (PLC), or a tetrahedral mesh. Depending on the input object and //
83
// the specified options, the output of TetGen is either a Delaunay (or wei- //
84
// ghted Delaunay) tetrahedralization, or a constrained (Delaunay) tetrahed- //
85
// ralization, or a quality tetrahedral mesh. //
86
// //
87
// A piecewise linear complex (PLC) represents a 3D polyhedral domain with //
88
// possibly internal boundaries(subdomains). It is introduced in [Miller et //
89
// al, 1996]. Basically it is a set of "cells", i.e., vertices, edges, poly- //
90
// gons, and polyhedra, and the intersection of any two of its cells is the //
91
// union of other cells of it. //
92
// //
93
// TetGen uses a set of files to describe the inputs and outputs. Each file //
94
// is identified from its file extension (.node, .ele, .face, .edge, etc). //
95
// //
96
// The 'tetgenio' structure is a collection of arrays of data, i.e., points, //
97
// facets, tetrahedra, and so forth. It contains functions to read and write //
98
// (input and output) files of TetGen as well as other supported mesh files. //
99
// //
100
// Once an object of tetgenio is declared, no array is created. One has to //
101
// allocate enough memory for them. On deletion of this object, the memory //
102
// occupied by these arrays needs to be freed. The routine deinitialize() //
103
// will be automatically called. It frees the memory for an array if it is //
104
// not a NULL. Note that it assumes that the memory is allocated by the C++ //
105
// "new" operator. Otherwise, the user is responsible to free them and all //
106
// pointers must be NULL before the call of the destructor. //
107
// //
108
///////////////////////////////////////////////////////////////////////////////
109
110
class tetgenio {
111
112
public:
113
114
// A "polygon" describes a simple polygon (no holes). It is not necessarily
115
// convex. Each polygon contains a number of corners (points) and the same
116
// number of sides (edges). The points of the polygon must be given in
117
// either counterclockwise or clockwise order and they form a ring, so
118
// every two consecutive points forms an edge of the polygon.
119
typedef struct {
120
int *vertexlist;
121
int numberofvertices;
122
} polygon;
123
124
// A "facet" describes a polygonal region possibly with holes, edges, and
125
// points floating in it. Each facet consists of a list of polygons and
126
// a list of hole points (which lie strictly inside holes).
127
typedef struct {
128
polygon *polygonlist;
129
int numberofpolygons;
130
REAL *holelist;
131
int numberofholes;
132
} facet;
133
134
// A "voroedge" is an edge of the Voronoi diagram. It corresponds to a
135
// Delaunay face. Each voroedge is either a line segment connecting
136
// two Voronoi vertices or a ray starting from a Voronoi vertex to an
137
// "infinite vertex". 'v1' and 'v2' are two indices pointing to the
138
// list of Voronoi vertices. 'v1' must be non-negative, while 'v2' may
139
// be -1 if it is a ray, in this case, the unit normal of this ray is
140
// given in 'vnormal'.
141
typedef struct {
142
int v1, v2;
143
REAL vnormal[3];
144
} voroedge;
145
146
// A "vorofacet" is an facet of the Voronoi diagram. It corresponds to a
147
// Delaunay edge. Each Voronoi facet is a convex polygon formed by a
148
// list of Voronoi edges, it may not be closed. 'c1' and 'c2' are two
149
// indices pointing into the list of Voronoi cells, i.e., the two cells
150
// share this facet. 'elist' is an array of indices pointing into the
151
// list of Voronoi edges, 'elist[0]' saves the number of Voronoi edges
152
// (including rays) of this facet.
153
typedef struct {
154
int c1, c2;
155
int *elist;
156
} vorofacet;
157
158
159
// Additional parameters associated with an input (or mesh) vertex.
160
// This information is provided by CAD libraries.
161
typedef struct {
162
REAL uv[2];
163
int tag;
164
int type; // 0, 1, or 2.
165
} pointparam;
166
167
// Callback functions for meshing PSCs.
168
typedef REAL (* GetVertexParamOnEdge)(void*, int, int);
169
typedef void (* GetSteinerOnEdge)(void*, int, REAL, REAL*);
170
typedef void (* GetVertexParamOnFace)(void*, int, int, REAL*);
171
typedef void (* GetEdgeSteinerParamOnFace)(void*, int, REAL, int, REAL*);
172
typedef void (* GetSteinerOnFace)(void*, int, REAL*, REAL*);
173
174
// A callback function for mesh refinement.
175
typedef bool (* TetSizeFunc)(REAL*, REAL*, REAL*, REAL*, REAL*, REAL);
176
177
// Items are numbered starting from 'firstnumber' (0 or 1), default is 0.
178
int firstnumber;
179
180
// Dimension of the mesh (2 or 3), default is 3.
181
int mesh_dim;
182
183
// Does the lines in .node file contain index or not, default is 1.
184
int useindex;
185
186
// 'pointlist': An array of point coordinates. The first point's x
187
// coordinate is at index [0] and its y coordinate at index [1], its
188
// z coordinate is at index [2], followed by the coordinates of the
189
// remaining points. Each point occupies three REALs.
190
// 'pointattributelist': An array of point attributes. Each point's
191
// attributes occupy 'numberofpointattributes' REALs.
192
// 'pointmtrlist': An array of metric tensors at points. Each point's
193
// tensor occupies 'numberofpointmtr' REALs.
194
// 'pointmarkerlist': An array of point markers; one integer per point.
195
// 'point2tetlist': An array of tetrahedra indices; one integer per point.
196
REAL *pointlist;
197
REAL *pointattributelist;
198
REAL *pointmtrlist;
199
int *pointmarkerlist;
200
int *point2tetlist;
201
pointparam *pointparamlist;
202
int numberofpoints;
203
int numberofpointattributes;
204
int numberofpointmtrs;
205
206
// 'tetrahedronlist': An array of tetrahedron corners. The first
207
// tetrahedron's first corner is at index [0], followed by its other
208
// corners, followed by six nodes on the edges of the tetrahedron if the
209
// second order option (-o2) is applied. Each tetrahedron occupies
210
// 'numberofcorners' ints. The second order nodes are output only.
211
// 'tetrahedronattributelist': An array of tetrahedron attributes. Each
212
// tetrahedron's attributes occupy 'numberoftetrahedronattributes' REALs.
213
// 'tetrahedronvolumelist': An array of constraints, i.e. tetrahedron's
214
// volume; one REAL per element. Input only.
215
// 'neighborlist': An array of tetrahedron neighbors; 4 ints per element.
216
// 'tet2facelist': An array of tetrahedron face indices; 4 ints per element.
217
// 'tet2edgelist': An array of tetrahedron edge indices; 6 ints per element.
218
int *tetrahedronlist;
219
REAL *tetrahedronattributelist;
220
REAL *tetrahedronvolumelist;
221
int *neighborlist;
222
int *tet2facelist;
223
int *tet2edgelist;
224
int numberoftetrahedra;
225
int numberofcorners;
226
int numberoftetrahedronattributes;
227
228
// 'facetlist': An array of facets. Each entry is a structure of facet.
229
// 'facetmarkerlist': An array of facet markers; one int per facet.
230
facet *facetlist;
231
int *facetmarkerlist;
232
int numberoffacets;
233
234
// 'holelist': An array of holes (in volume). Each hole is given by a
235
// seed (point) which lies strictly inside it. The first seed's x, y and z
236
// coordinates are at indices [0], [1] and [2], followed by the
237
// remaining seeds. Three REALs per hole.
238
REAL *holelist;
239
int numberofholes;
240
241
// 'regionlist': An array of regions (subdomains). Each region is given by
242
// a seed (point) which lies strictly inside it. The first seed's x, y and
243
// z coordinates are at indices [0], [1] and [2], followed by the regional
244
// attribute at index [3], followed by the maximum volume at index [4].
245
// Five REALs per region.
246
// Note that each regional attribute is used only if you select the 'A'
247
// switch, and each volume constraint is used only if you select the
248
// 'a' switch (with no number following).
249
REAL *regionlist;
250
int numberofregions;
251
252
// 'facetconstraintlist': An array of facet constraints. Each constraint
253
// specifies a maximum area bound on the subfaces of that facet. The
254
// first facet constraint is given by a facet marker at index [0] and its
255
// maximum area bound at index [1], followed by the remaining facet con-
256
// straints. Two REALs per facet constraint. Note: the facet marker is
257
// actually an integer.
258
REAL *facetconstraintlist;
259
int numberoffacetconstraints;
260
261
// 'segmentconstraintlist': An array of segment constraints. Each constraint
262
// specifies a maximum length bound on the subsegments of that segment.
263
// The first constraint is given by the two endpoints of the segment at
264
// index [0] and [1], and the maximum length bound at index [2], followed
265
// by the remaining segment constraints. Three REALs per constraint.
266
// Note the segment endpoints are actually integers.
267
REAL *segmentconstraintlist;
268
int numberofsegmentconstraints;
269
270
271
// 'trifacelist': An array of face (triangle) corners. The first face's
272
// three corners are at indices [0], [1] and [2], followed by the remaining
273
// faces. Three ints per face.
274
// 'trifacemarkerlist': An array of face markers; one int per face.
275
// 'o2facelist': An array of second order nodes (on the edges) of the face.
276
// It is output only if the second order option (-o2) is applied. The
277
// first face's three second order nodes are at [0], [1], and [2],
278
// followed by the remaining faces. Three ints per face.
279
// 'face2tetlist': An array of tetrahedra indices; 2 ints per face.
280
// 'face2edgelist': An array of edge indices; 3 ints per face.
281
int *trifacelist;
282
int *trifacemarkerlist;
283
int *o2facelist;
284
int *face2tetlist;
285
int *face2edgelist;
286
int numberoftrifaces;
287
288
// 'edgelist': An array of edge endpoints. The first edge's endpoints
289
// are at indices [0] and [1], followed by the remaining edges.
290
// Two ints per edge.
291
// 'edgemarkerlist': An array of edge markers; one int per edge.
292
// 'o2edgelist': An array of midpoints of edges. It is output only if the
293
// second order option (-o2) is applied. One int per edge.
294
// 'edge2tetlist': An array of tetrahedra indices. One int per edge.
295
int *edgelist;
296
int *edgemarkerlist;
297
int *o2edgelist;
298
int *edge2tetlist;
299
int numberofedges;
300
301
// 'vpointlist': An array of Voronoi vertex coordinates (like pointlist).
302
// 'vedgelist': An array of Voronoi edges. Each entry is a 'voroedge'.
303
// 'vfacetlist': An array of Voronoi facets. Each entry is a 'vorofacet'.
304
// 'vcelllist': An array of Voronoi cells. Each entry is an array of
305
// indices pointing into 'vfacetlist'. The 0th entry is used to store
306
// the length of this array.
307
REAL *vpointlist;
308
voroedge *vedgelist;
309
vorofacet *vfacetlist;
310
int **vcelllist;
311
int numberofvpoints;
312
int numberofvedges;
313
int numberofvfacets;
314
int numberofvcells;
315
316
317
// Variable (and callback functions) for meshing PSCs.
318
void *geomhandle;
319
GetVertexParamOnEdge getvertexparamonedge;
320
GetSteinerOnEdge getsteineronedge;
321
GetVertexParamOnFace getvertexparamonface;
322
GetEdgeSteinerParamOnFace getedgesteinerparamonface;
323
GetSteinerOnFace getsteineronface;
324
325
// A callback function.
326
TetSizeFunc tetunsuitable;
327
328
// Input & output routines.
329
virtual bool load_node_call(FILE* infile, int markers, int uvflag, char*);
330
virtual bool load_node(char*);
331
virtual bool load_edge(char*);
332
virtual bool load_face(char*);
333
virtual bool load_tet(char*);
334
virtual bool load_vol(char*);
335
virtual bool load_var(char*);
336
virtual bool load_mtr(char*);
337
// virtual bool load_pbc(char*);
338
virtual bool load_poly(char*);
339
virtual bool load_off(char*);
340
virtual bool load_ply(char*);
341
virtual bool load_stl(char*);
342
virtual bool load_vtk(char*);
343
virtual bool load_medit(char*, int);
344
virtual bool load_plc(char*, int);
345
virtual bool load_tetmesh(char*, int);
346
virtual void save_nodes(char*);
347
virtual void save_elements(char*);
348
virtual void save_faces(char*);
349
virtual void save_edges(char*);
350
virtual void save_neighbors(char*);
351
virtual void save_poly(char*);
352
virtual void save_faces2smesh(char*);
353
354
// Read line and parse string functions.
355
virtual char *readline(char* string, FILE* infile, int *linenumber);
356
virtual char *findnextfield(char* string);
357
virtual char *readnumberline(char* string, FILE* infile, char* infilename);
358
virtual char *findnextnumber(char* string);
359
360
virtual void init(polygon* p) {
361
p->vertexlist = (int *) NULL;
362
p->numberofvertices = 0;
363
}
364
365
virtual void init(facet* f) {
366
f->polygonlist = (polygon *) NULL;
367
f->numberofpolygons = 0;
368
f->holelist = (REAL *) NULL;
369
f->numberofholes = 0;
370
}
371
372
// Initialize routine.
373
virtual void initialize()
374
{
375
firstnumber = 0;
376
mesh_dim = 3;
377
useindex = 1;
378
379
pointlist = (REAL *) NULL;
380
pointattributelist = (REAL *) NULL;
381
pointmtrlist = (REAL *) NULL;
382
pointmarkerlist = (int *) NULL;
383
point2tetlist = (int *) NULL;
384
pointparamlist = (pointparam *) NULL;
385
numberofpoints = 0;
386
numberofpointattributes = 0;
387
numberofpointmtrs = 0;
388
389
tetrahedronlist = (int *) NULL;
390
tetrahedronattributelist = (REAL *) NULL;
391
tetrahedronvolumelist = (REAL *) NULL;
392
neighborlist = (int *) NULL;
393
tet2facelist = (int *) NULL;
394
tet2edgelist = (int *) NULL;
395
numberoftetrahedra = 0;
396
numberofcorners = 4;
397
numberoftetrahedronattributes = 0;
398
399
trifacelist = (int *) NULL;
400
trifacemarkerlist = (int *) NULL;
401
o2facelist = (int *) NULL;
402
face2tetlist = (int *) NULL;
403
face2edgelist = (int *) NULL;
404
numberoftrifaces = 0;
405
406
edgelist = (int *) NULL;
407
edgemarkerlist = (int *) NULL;
408
o2edgelist = (int *) NULL;
409
edge2tetlist = (int *) NULL;
410
numberofedges = 0;
411
412
facetlist = (facet *) NULL;
413
facetmarkerlist = (int *) NULL;
414
numberoffacets = 0;
415
416
holelist = (REAL *) NULL;
417
numberofholes = 0;
418
419
regionlist = (REAL *) NULL;
420
numberofregions = 0;
421
422
facetconstraintlist = (REAL *) NULL;
423
numberoffacetconstraints = 0;
424
segmentconstraintlist = (REAL *) NULL;
425
numberofsegmentconstraints = 0;
426
427
428
vpointlist = (REAL *) NULL;
429
vedgelist = (voroedge *) NULL;
430
vfacetlist = (vorofacet *) NULL;
431
vcelllist = (int **) NULL;
432
numberofvpoints = 0;
433
numberofvedges = 0;
434
numberofvfacets = 0;
435
numberofvcells = 0;
436
437
438
tetunsuitable = NULL;
439
440
geomhandle = NULL;
441
getvertexparamonedge = NULL;
442
getsteineronedge = NULL;
443
getvertexparamonface = NULL;
444
getedgesteinerparamonface = NULL;
445
getsteineronface = NULL;
446
}
447
448
// Free the memory allocated in 'tetgenio'. Note that it assumes that the
449
// memory was allocated by the "new" operator (C++).
450
virtual void deinitialize()
451
{
452
int i, j;
453
454
if (pointlist != (REAL *) NULL) {
455
delete [] pointlist;
456
}
457
if (pointattributelist != (REAL *) NULL) {
458
delete [] pointattributelist;
459
}
460
if (pointmtrlist != (REAL *) NULL) {
461
delete [] pointmtrlist;
462
}
463
if (pointmarkerlist != (int *) NULL) {
464
delete [] pointmarkerlist;
465
}
466
if (point2tetlist != (int *) NULL) {
467
delete [] point2tetlist;
468
}
469
if (pointparamlist != (pointparam *) NULL) {
470
delete [] pointparamlist;
471
}
472
473
if (tetrahedronlist != (int *) NULL) {
474
delete [] tetrahedronlist;
475
}
476
if (tetrahedronattributelist != (REAL *) NULL) {
477
delete [] tetrahedronattributelist;
478
}
479
if (tetrahedronvolumelist != (REAL *) NULL) {
480
delete [] tetrahedronvolumelist;
481
}
482
if (neighborlist != (int *) NULL) {
483
delete [] neighborlist;
484
}
485
if (tet2facelist != (int *) NULL) {
486
delete [] tet2facelist;
487
}
488
if (tet2edgelist != (int *) NULL) {
489
delete [] tet2edgelist;
490
}
491
492
if (trifacelist != (int *) NULL) {
493
delete [] trifacelist;
494
}
495
if (trifacemarkerlist != (int *) NULL) {
496
delete [] trifacemarkerlist;
497
}
498
if (o2facelist != (int *) NULL) {
499
delete [] o2facelist;
500
}
501
if (face2tetlist != (int *) NULL) {
502
delete [] face2tetlist;
503
}
504
if (face2edgelist != (int *) NULL) {
505
delete [] face2edgelist;
506
}
507
508
if (edgelist != (int *) NULL) {
509
delete [] edgelist;
510
}
511
if (edgemarkerlist != (int *) NULL) {
512
delete [] edgemarkerlist;
513
}
514
if (o2edgelist != (int *) NULL) {
515
delete [] o2edgelist;
516
}
517
if (edge2tetlist != (int *) NULL) {
518
delete [] edge2tetlist;
519
}
520
521
if (facetlist != (facet *) NULL) {
522
facet *f;
523
polygon *p;
524
for (i = 0; i < numberoffacets; i++) {
525
f = &facetlist[i];
526
for (j = 0; j < f->numberofpolygons; j++) {
527
p = &f->polygonlist[j];
528
delete [] p->vertexlist;
529
}
530
delete [] f->polygonlist;
531
if (f->holelist != (REAL *) NULL) {
532
delete [] f->holelist;
533
}
534
}
535
delete [] facetlist;
536
}
537
if (facetmarkerlist != (int *) NULL) {
538
delete [] facetmarkerlist;
539
}
540
541
if (holelist != (REAL *) NULL) {
542
delete [] holelist;
543
}
544
if (regionlist != (REAL *) NULL) {
545
delete [] regionlist;
546
}
547
if (facetconstraintlist != (REAL *) NULL) {
548
delete [] facetconstraintlist;
549
}
550
if (segmentconstraintlist != (REAL *) NULL) {
551
delete [] segmentconstraintlist;
552
}
553
if (vpointlist != (REAL *) NULL) {
554
delete [] vpointlist;
555
}
556
if (vedgelist != (voroedge *) NULL) {
557
delete [] vedgelist;
558
}
559
if (vfacetlist != (vorofacet *) NULL) {
560
for (i = 0; i < numberofvfacets; i++) {
561
delete [] vfacetlist[i].elist;
562
}
563
delete [] vfacetlist;
564
}
565
if (vcelllist != (int **) NULL) {
566
for (i = 0; i < numberofvcells; i++) {
567
delete [] vcelllist[i];
568
}
569
delete [] vcelllist;
570
}
571
}
572
573
// Constructor & destructor.
574
tetgenio() {initialize();}
575
virtual ~tetgenio() {deinitialize();}
576
577
}; // class tetgenio
578
579
///////////////////////////////////////////////////////////////////////////////
580
// //
581
// tetgenbehavior //
582
// //
583
// A structure for maintaining the switches and parameters used by TetGen's //
584
// mesh data structure and algorithms. //
585
// //
586
// All switches and parameters are initialized with default values. They can //
587
// be set by the command line arguments (a list of strings) of TetGen. //
588
// //
589
// NOTE: Some of the switches are incompatible. While some may depend on //
590
// other switches. The routine parse_commandline() sets the switches from //
591
// the command line (a list of strings) and checks the consistency of the //
592
// applied switches. //
593
// //
594
///////////////////////////////////////////////////////////////////////////////
595
596
class tetgenbehavior {
597
598
public:
599
600
// Switches of TetGen.
601
int plc; // '-p', 0.
602
int psc; // '-s', 0.
603
int refine; // '-r', 0.
604
int quality; // '-q', 0.
605
int nobisect; // '-Y', 0.
606
int coarsen; // '-R', 0.
607
int weighted; // '-w', 0.
608
int brio_hilbert; // '-b', 1.
609
int incrflip; // '-l', 0.
610
int flipinsert; // '-L', 0.
611
int metric; // '-m', 0.
612
int varvolume; // '-a', 0.
613
int fixedvolume; // '-a', 0.
614
int regionattrib; // '-A', 0.
615
int cdtrefine; // '-D', 0.
616
int insertaddpoints; // '-i', 0.
617
int diagnose; // '-d', 0.
618
int convex; // '-c', 0.
619
int nomergefacet; // '-M', 0.
620
int nomergevertex; // '-M', 0.
621
int noexact; // '-X', 0.
622
int nostaticfilter; // '-X', 0.
623
int zeroindex; // '-z', 0.
624
int facesout; // '-f', 0.
625
int edgesout; // '-e', 0.
626
int neighout; // '-n', 0.
627
int voroout; // '-v', 0.
628
int meditview; // '-g', 0.
629
int vtkview; // '-k', 0.
630
int nobound; // '-B', 0.
631
int nonodewritten; // '-N', 0.
632
int noelewritten; // '-E', 0.
633
int nofacewritten; // '-F', 0.
634
int noiterationnum; // '-I', 0.
635
int nojettison; // '-J', 0.
636
int docheck; // '-C', 0.
637
int quiet; // '-Q', 0.
638
int verbose; // '-V', 0.
639
640
// Parameters of TetGen.
641
int vertexperblock; // '-x', 4092.
642
int tetrahedraperblock; // '-x', 8188.
643
int shellfaceperblock; // '-x', 2044.
644
int nobisect_nomerge; // '-Y', 1.
645
int supsteiner_level; // '-Y/', 2.
646
int addsteiner_algo; // '-Y//', 1.
647
int coarsen_param; // '-R', 0.
648
int weighted_param; // '-w', 0.
649
int fliplinklevel; // -1.
650
int flipstarsize; // -1.
651
int fliplinklevelinc; // 1.
652
int reflevel; // '-D', 3.
653
int optlevel; // '-O', 2.
654
int optscheme; // '-O', 7.
655
int delmaxfliplevel; // 1.
656
int order; // '-o', 1.
657
int reversetetori; // '-o/', 0.
658
int steinerleft; // '-S', 0.
659
int no_sort; // 0.
660
int hilbert_order; // '-b///', 52.
661
int hilbert_limit; // '-b//' 8.
662
int brio_threshold; // '-b' 64.
663
REAL brio_ratio; // '-b/' 0.125.
664
REAL facet_separate_ang_tol; // '-p', 179.9.
665
REAL facet_overlap_ang_tol; // '-p/', 0.1.
666
REAL facet_small_ang_tol; // '-p//', 15.0.
667
REAL maxvolume; // '-a', -1.0.
668
REAL minratio; // '-q', 0.0.
669
REAL mindihedral; // '-q', 5.0.
670
REAL optmaxdihedral; // 165.0.
671
REAL optminsmtdihed; // 179.0.
672
REAL optminslidihed; // 179.0.
673
REAL epsilon; // '-T', 1.0e-8.
674
REAL coarsen_percent; // -R1/#, 1.0.
675
676
// Strings of command line arguments and input/output file names.
677
char commandline[1024];
678
char infilename[1024];
679
char outfilename[1024];
680
char addinfilename[1024];
681
char bgmeshfilename[1024];
682
683
// The input object of TetGen. They are recognized by either the input
684
// file extensions or by the specified options.
685
// Currently the following objects are supported:
686
// - NODES, a list of nodes (.node);
687
// - POLY, a piecewise linear complex (.poly or .smesh);
688
// - OFF, a polyhedron (.off, Geomview's file format);
689
// - PLY, a polyhedron (.ply, file format from gatech, only ASCII);
690
// - STL, a surface mesh (.stl, stereolithography format);
691
// - MEDIT, a surface mesh (.mesh, Medit's file format);
692
// - MESH, a tetrahedral mesh (.ele).
693
// If no extension is available, the imposed command line switch
694
// (-p or -r) implies the object.
695
enum objecttype {NODES, POLY, OFF, PLY, STL, MEDIT, VTK, MESH} object;
696
697
698
virtual void syntax();
699
virtual void usage();
700
701
// Command line parse routine.
702
virtual bool parse_commandline(int argc, char **argv);
703
virtual bool parse_commandline(char *switches) {
704
return parse_commandline(0, &switches);
705
}
706
707
// Initialize all variables.
708
tetgenbehavior()
709
{
710
plc = 0;
711
psc = 0;
712
refine = 0;
713
quality = 0;
714
nobisect = 0;
715
coarsen = 0;
716
metric = 0;
717
weighted = 0;
718
brio_hilbert = 1;
719
incrflip = 0;
720
flipinsert = 0;
721
varvolume = 0;
722
fixedvolume = 0;
723
noexact = 0;
724
nostaticfilter = 0;
725
insertaddpoints = 0;
726
regionattrib = 0;
727
cdtrefine = 0;
728
diagnose = 0;
729
convex = 0;
730
zeroindex = 0;
731
facesout = 0;
732
edgesout = 0;
733
neighout = 0;
734
voroout = 0;
735
meditview = 0;
736
vtkview = 0;
737
nobound = 0;
738
nonodewritten = 0;
739
noelewritten = 0;
740
nofacewritten = 0;
741
noiterationnum = 0;
742
nomergefacet = 0;
743
nomergevertex = 0;
744
nojettison = 0;
745
docheck = 0;
746
quiet = 0;
747
verbose = 0;
748
749
vertexperblock = 4092;
750
tetrahedraperblock = 8188;
751
shellfaceperblock = 4092;
752
nobisect_nomerge = 1;
753
supsteiner_level = 2;
754
addsteiner_algo = 1;
755
coarsen_param = 0;
756
weighted_param = 0;
757
fliplinklevel = -1;
758
flipstarsize = -1;
759
fliplinklevelinc = 1;
760
reflevel = 3;
761
optscheme = 7;
762
optlevel = 2;
763
delmaxfliplevel = 1;
764
order = 1;
765
reversetetori = 0;
766
steinerleft = -1;
767
no_sort = 0;
768
hilbert_order = 52; //-1;
769
hilbert_limit = 8;
770
brio_threshold = 64;
771
brio_ratio = 0.125;
772
facet_separate_ang_tol = 179.9;
773
facet_overlap_ang_tol = 0.1;
774
facet_small_ang_tol = 15.0;
775
maxvolume = -1.0;
776
minratio = 2.0;
777
mindihedral = 0.0;
778
optmaxdihedral = 165.00;
779
optminsmtdihed = 179.00;
780
optminslidihed = 179.00;
781
epsilon = 1.0e-8;
782
coarsen_percent = 1.0;
783
object = NODES;
784
785
commandline[0] = '\0';
786
infilename[0] = '\0';
787
outfilename[0] = '\0';
788
addinfilename[0] = '\0';
789
bgmeshfilename[0] = '\0';
790
791
}
792
793
}; // class tetgenbehavior
794
795
///////////////////////////////////////////////////////////////////////////////
796
// //
797
// Robust Geometric predicates //
798
// //
799
// Geometric predicates are simple tests of spatial relations of a set of d- //
800
// dimensional points, such as the orientation test and the point-in-sphere //
801
// test. Each of these tests is performed by evaluating the sign of a deter- //
802
// minant of a matrix whose entries are the coordinates of these points. If //
803
// the computation is performed by using the floating-point numbers, e.g., //
804
// the single or double precision numbers in C/C++, roundoff error may cause //
805
// an incorrect result. This may either lead to a wrong result or eventually //
806
// lead to a failure of the program. Computing the predicates exactly will //
807
// avoid the error and make the program robust. //
808
// //
809
// The following routines are the robust geometric predicates for 3D orient- //
810
// ation test and point-in-sphere test. They were implemented by Shewchuk. //
811
// The source code are generously provided by him in the public domain, //
812
// http://www.cs.cmu.edu/~quake/robust.html. predicates.cxx is a C++ version //
813
// of the original C code. //
814
// //
815
// The original predicates of Shewchuk only use "dynamic filters", i.e., it //
816
// computes the error at run time step by step. TetGen first adds a "static //
817
// filter" in each predicate. It estimates the maximal possible error in all //
818
// cases. So it can safely and quickly answer many easy cases. //
819
// //
820
///////////////////////////////////////////////////////////////////////////////
821
822
void exactinit(int, int, int, REAL, REAL, REAL);
823
REAL orient3d(REAL *pa, REAL *pb, REAL *pc, REAL *pd);
824
REAL insphere(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe);
825
REAL orient4d(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe,
826
REAL ah, REAL bh, REAL ch, REAL dh, REAL eh);
827
828
///////////////////////////////////////////////////////////////////////////////
829
// //
830
// tetgenmesh //
831
// //
832
// A structure for creating and updating tetrahedral meshes. //
833
// //
834
///////////////////////////////////////////////////////////////////////////////
835
836
class tetgenmesh {
837
838
public:
839
840
///////////////////////////////////////////////////////////////////////////////
841
// //
842
// Mesh data structure //
843
// //
844
// A tetrahedral mesh T of a 3D piecewise linear complex (PLC) X is a 3D //
845
// simplicial complex whose underlying space is equal to the space of X. T //
846
// contains a 2D subcomplex S which is a triangular mesh of the boundary of //
847
// X. S contains a 1D subcomplex L which is a linear mesh of the boundary of //
848
// S. Faces and edges in S and L are respectively called subfaces and segme- //
849
// nts to distinguish them from others in T. //
850
// //
851
// TetGen stores the tetrahedra and vertices of T. The basic structure of a //
852
// tetrahedron contains pointers to its vertices and adjacent tetrahedra. A //
853
// vertex stores its x-, y-, and z-coordinates, and a pointer to a tetrahed- //
854
// ron containing it. Both tetrahedra and vertices may contain user data. //
855
// //
856
// Each face of T belongs to either two tetrahedra or one tetrahedron. In //
857
// the latter case, the face is an exterior boundary face of T. TetGen adds //
858
// fictitious tetrahedra (one-to-one) at such faces, and connects them to an //
859
// "infinite vertex" (which has no geometric coordinates). One can imagine //
860
// such a vertex lies in 4D space and is visible by all exterior boundary //
861
// faces. The extended set of tetrahedra (including the infinite vertex) is //
862
// a tetrahedralization of a 3-pseudomanifold without boundary. It has the //
863
// property that every face is shared by exactly two tetrahedra. //
864
// //
865
// The current version of TetGen stores explicitly the subfaces and segments //
866
// (which are in surface mesh S and the linear mesh L), respectively. Extra //
867
// pointers are allocated in tetrahedra and subfaces to point each others. //
868
// //
869
///////////////////////////////////////////////////////////////////////////////
870
871
// The tetrahedron data structure. It includes the following fields:
872
// - a list of four adjoining tetrahedra;
873
// - a list of four vertices;
874
// - a pointer to a list of four subfaces (optional, for -p switch);
875
// - a pointer to a list of six segments (optional, for -p switch);
876
// - a list of user-defined floating-point attributes (optional);
877
// - a volume constraint (optional, for -a switch);
878
// - an integer of element marker (and flags);
879
// The structure of a tetrahedron is an array of pointers. Its actual size
880
// (the length of the array) is determined at runtime.
881
882
typedef REAL **tetrahedron;
883
884
// The subface data structure. It includes the following fields:
885
// - a list of three adjoining subfaces;
886
// - a list of three vertices;
887
// - a list of three adjoining segments;
888
// - two adjoining tetrahedra;
889
// - an area constraint (optional, for -q switch);
890
// - an integer for boundary marker;
891
// - an integer for type, flags, etc.
892
893
typedef REAL **shellface;
894
895
// The point data structure. It includes the following fields:
896
// - x, y and z coordinates;
897
// - a list of user-defined point attributes (optional);
898
// - u, v coordinates (optional, for -s switch);
899
// - a metric tensor (optional, for -q or -m switch);
900
// - a pointer to an adjacent tetrahedron;
901
// - a pointer to a parent (or a duplicate) point;
902
// - a pointer to an adjacent subface or segment (optional, -p switch);
903
// - a pointer to a tet in background mesh (optional, for -m switch);
904
// - an integer for boundary marker (point index);
905
// - an integer for point type (and flags).
906
// - an integer for geometry tag (optional, for -s switch).
907
// The structure of a point is an array of REALs. Its actual size is
908
// determined at the runtime.
909
910
typedef REAL *point;
911
912
///////////////////////////////////////////////////////////////////////////////
913
// //
914
// Handles //
915
// //
916
// Navigation and manipulation in a tetrahedralization are accomplished by //
917
// operating on structures referred as ``handles". A handle is a pair (t,v), //
918
// where t is a pointer to a tetrahedron, and v is a 4-bit integer, in the //
919
// range from 0 to 11. v is called the ``version'' of a tetrahedron, it rep- //
920
// resents a directed edge of a specific face of the tetrahedron. //
921
// //
922
// There are 12 even permutations of the four vertices, each of them corres- //
923
// ponds to a directed edge (a version) of the tetrahedron. The 12 versions //
924
// can be grouped into 4 distinct ``edge rings'' in 4 ``oriented faces'' of //
925
// this tetrahedron. One can encode each version (a directed edge) into a //
926
// 4-bit integer such that the two upper bits encode the index (from 0 to 2) //
927
// of this edge in the edge ring, and the two lower bits encode the index ( //
928
// from 0 to 3) of the oriented face which contains this edge. //
929
// //
930
// The four vertices of a tetrahedron are indexed from 0 to 3 (according to //
931
// their storage in the data structure). Give each face the same index as //
932
// the node opposite it in the tetrahedron. Denote the edge connecting face //
933
// i to face j as i/j. We number the twelve versions as follows: //
934
// //
935
// | edge 0 edge 1 edge 2 //
936
// --------|-------------------------------- //
937
// face 0 | 0 (0/1) 4 (0/3) 8 (0/2) //
938
// face 1 | 1 (1/2) 5 (1/3) 9 (1/0) //
939
// face 2 | 2 (2/3) 6 (2/1) 10 (2/0) //
940
// face 3 | 3 (3/0) 7 (3/1) 11 (3/2) //
941
// //
942
// Similarly, navigation and manipulation in a (boundary) triangulation are //
943
// done by using handles of triangles. Each handle is a pair (s, v), where s //
944
// is a pointer to a triangle, and v is a version in the range from 0 to 5. //
945
// Each version corresponds to a directed edge of this triangle. //
946
// //
947
// Number the three vertices of a triangle from 0 to 2 (according to their //
948
// storage in the data structure). Give each edge the same index as the node //
949
// opposite it in the triangle. The six versions of a triangle are: //
950
// //
951
// | edge 0 edge 1 edge 2 //
952
// ------------------|-------------------------- //
953
// ccw orientation | 0 2 4 //
954
// cw orientation | 1 3 5 //
955
// //
956
// In the following, a 'triface' is a handle of tetrahedron, and a 'face' is //
957
// a handle of a triangle. //
958
// //
959
///////////////////////////////////////////////////////////////////////////////
960
961
class triface {
962
public:
963
tetrahedron *tet;
964
int ver; // Range from 0 to 11.
965
triface() : tet(0), ver(0) {}
966
triface& operator=(const triface& t) {
967
tet = t.tet; ver = t.ver;
968
return *this;
969
}
970
};
971
972
class face {
973
public:
974
shellface *sh;
975
int shver; // Range from 0 to 5.
976
face() : sh(0), shver(0) {}
977
face& operator=(const face& s) {
978
sh = s.sh; shver = s.shver;
979
return *this;
980
}
981
};
982
983
///////////////////////////////////////////////////////////////////////////////
984
// //
985
// Arraypool //
986
// //
987
// A dynamic linear array. (It is written by J. Shewchuk) //
988
// //
989
// Each arraypool contains an array of pointers to a number of blocks. Each //
990
// block contains the same fixed number of objects. Each index of the array //
991
// addresses a particular object in the pool. The most significant bits add- //
992
// ress the index of the block containing the object. The less significant //
993
// bits address this object within the block. //
994
// //
995
// 'objectbytes' is the size of one object in blocks; 'log2objectsperblock' //
996
// is the base-2 logarithm of 'objectsperblock'; 'objects' counts the number //
997
// of allocated objects; 'totalmemory' is the total memory in bytes. //
998
// //
999
///////////////////////////////////////////////////////////////////////////////
1000
1001
class arraypool {
1002
1003
public:
1004
1005
int objectbytes;
1006
int objectsperblock;
1007
int log2objectsperblock;
1008
int objectsperblockmark;
1009
int toparraylen;
1010
char **toparray;
1011
long objects;
1012
unsigned long totalmemory;
1013
1014
void restart();
1015
void poolinit(int sizeofobject, int log2objperblk);
1016
char* getblock(int objectindex);
1017
void* lookup(int objectindex);
1018
int newindex(void **newptr);
1019
1020
arraypool(int sizeofobject, int log2objperblk);
1021
~arraypool();
1022
};
1023
1024
// fastlookup() -- A fast, unsafe operation. Return the pointer to the object
1025
// with a given index. Note: The object's block must have been allocated,
1026
// i.e., by the function newindex().
1027
1028
#define fastlookup(pool, index) \
1029
(void *) ((pool)->toparray[(index) >> (pool)->log2objectsperblock] + \
1030
((index) & (pool)->objectsperblockmark) * (pool)->objectbytes)
1031
1032
///////////////////////////////////////////////////////////////////////////////
1033
// //
1034
// Memorypool //
1035
// //
1036
// A structure for memory allocation. (It is written by J. Shewchuk) //
1037
// //
1038
// firstblock is the first block of items. nowblock is the block from which //
1039
// items are currently being allocated. nextitem points to the next slab //
1040
// of free memory for an item. deaditemstack is the head of a linked list //
1041
// (stack) of deallocated items that can be recycled. unallocateditems is //
1042
// the number of items that remain to be allocated from nowblock. //
1043
// //
1044
// Traversal is the process of walking through the entire list of items, and //
1045
// is separate from allocation. Note that a traversal will visit items on //
1046
// the "deaditemstack" stack as well as live items. pathblock points to //
1047
// the block currently being traversed. pathitem points to the next item //
1048
// to be traversed. pathitemsleft is the number of items that remain to //
1049
// be traversed in pathblock. //
1050
// //
1051
///////////////////////////////////////////////////////////////////////////////
1052
1053
class memorypool {
1054
1055
public:
1056
1057
void **firstblock, **nowblock;
1058
void *nextitem;
1059
void *deaditemstack;
1060
void **pathblock;
1061
void *pathitem;
1062
int alignbytes;
1063
int itembytes, itemwords;
1064
int itemsperblock;
1065
long items, maxitems;
1066
int unallocateditems;
1067
int pathitemsleft;
1068
1069
memorypool();
1070
memorypool(int, int, int, int);
1071
~memorypool();
1072
1073
void poolinit(int, int, int, int);
1074
void restart();
1075
void *alloc();
1076
void dealloc(void*);
1077
void traversalinit();
1078
void *traverse();
1079
};
1080
1081
///////////////////////////////////////////////////////////////////////////////
1082
// //
1083
// badface //
1084
// //
1085
// Despite of its name, a 'badface' can be used to represent one of the //
1086
// following objects: //
1087
// - a face of a tetrahedron which is (possibly) non-Delaunay; //
1088
// - an encroached subsegment or subface; //
1089
// - a bad-quality tetrahedron, i.e, has too large radius-edge ratio; //
1090
// - a sliver, i.e., has good radius-edge ratio but nearly zero volume; //
1091
// - a recently flipped face (saved for undoing the flip later). //
1092
// //
1093
///////////////////////////////////////////////////////////////////////////////
1094
1095
class badface {
1096
public:
1097
triface tt;
1098
face ss;
1099
REAL key, cent[6]; // circumcenter or cos(dihedral angles) at 6 edges.
1100
point forg, fdest, fapex, foppo, noppo;
1101
badface *nextitem;
1102
badface() : key(0), forg(0), fdest(0), fapex(0), foppo(0), noppo(0),
1103
nextitem(0) {}
1104
};
1105
1106
///////////////////////////////////////////////////////////////////////////////
1107
// //
1108
// insertvertexflags //
1109
// //
1110
// A collection of flags that pass to the routine insertvertex(). //
1111
// //
1112
///////////////////////////////////////////////////////////////////////////////
1113
1114
class insertvertexflags {
1115
1116
public:
1117
1118
int iloc; // input/output.
1119
int bowywat, lawson;
1120
int splitbdflag, validflag, respectbdflag;
1121
int rejflag, chkencflag, cdtflag;
1122
int assignmeshsize;
1123
int sloc, sbowywat;
1124
1125
// Used by Delaunay refinement.
1126
int refineflag; // 0, 1, 2, 3
1127
triface refinetet;
1128
face refinesh;
1129
int smlenflag; // for useinsertradius.
1130
REAL smlen; // for useinsertradius.
1131
point parentpt;
1132
1133
insertvertexflags() {
1134
iloc = bowywat = lawson = 0;
1135
splitbdflag = validflag = respectbdflag = 0;
1136
rejflag = chkencflag = cdtflag = 0;
1137
assignmeshsize = 0;
1138
sloc = sbowywat = 0;
1139
1140
refineflag = 0;
1141
refinetet.tet = NULL;
1142
refinesh.sh = NULL;
1143
smlenflag = 0;
1144
smlen = 0.0;
1145
}
1146
};
1147
1148
///////////////////////////////////////////////////////////////////////////////
1149
// //
1150
// flipconstraints //
1151
// //
1152
// A structure of a collection of data (options and parameters) which pass //
1153
// to the edge flip function flipnm(). //
1154
// //
1155
///////////////////////////////////////////////////////////////////////////////
1156
1157
class flipconstraints {
1158
1159
public:
1160
1161
// Elementary flip flags.
1162
int enqflag; // (= flipflag)
1163
int chkencflag;
1164
1165
// Control flags
1166
int unflip; // Undo the performed flips.
1167
int collectnewtets; // Collect the new tets created by flips.
1168
int collectencsegflag;
1169
1170
// Optimization flags.
1171
int remove_ndelaunay_edge; // Remove a non-Delaunay edge.
1172
REAL bak_tetprism_vol; // The value to be minimized.
1173
REAL tetprism_vol_sum;
1174
int remove_large_angle; // Remove a large dihedral angle at edge.
1175
REAL cosdihed_in; // The input cosine of the dihedral angle (> 0).
1176
REAL cosdihed_out; // The improved cosine of the dihedral angle.
1177
1178
// Boundary recovery flags.
1179
int checkflipeligibility;
1180
point seg[2]; // A constraining edge to be recovered.
1181
point fac[3]; // A constraining face to be recovered.
1182
point remvert; // A vertex to be removed.
1183
1184
1185
flipconstraints() {
1186
enqflag = 0;
1187
chkencflag = 0;
1188
1189
unflip = 0;
1190
collectnewtets = 0;
1191
collectencsegflag = 0;
1192
1193
remove_ndelaunay_edge = 0;
1194
bak_tetprism_vol = 0.0;
1195
tetprism_vol_sum = 0.0;
1196
remove_large_angle = 0;
1197
cosdihed_in = 0.0;
1198
cosdihed_out = 0.0;
1199
1200
checkflipeligibility = 0;
1201
seg[0] = NULL;
1202
fac[0] = NULL;
1203
remvert = NULL;
1204
}
1205
};
1206
1207
///////////////////////////////////////////////////////////////////////////////
1208
// //
1209
// optparameters //
1210
// //
1211
// Optimization options and parameters. //
1212
// //
1213
///////////////////////////////////////////////////////////////////////////////
1214
1215
class optparameters {
1216
1217
public:
1218
1219
// The one of goals of optimization.
1220
int max_min_volume; // Maximize the minimum volume.
1221
int min_max_aspectratio; // Minimize the maximum aspect ratio.
1222
int min_max_dihedangle; // Minimize the maximum dihedral angle.
1223
1224
// The initial and improved value.
1225
REAL initval, imprval;
1226
1227
int numofsearchdirs;
1228
REAL searchstep;
1229
int maxiter; // Maximum smoothing iterations (disabled by -1).
1230
int smthiter; // Performed iterations.
1231
1232
1233
optparameters() {
1234
max_min_volume = 0;
1235
min_max_aspectratio = 0;
1236
min_max_dihedangle = 0;
1237
1238
initval = imprval = 0.0;
1239
1240
numofsearchdirs = 10;
1241
searchstep = 0.01;
1242
maxiter = -1; // Unlimited smoothing iterations.
1243
smthiter = 0;
1244
1245
}
1246
};
1247
1248
1249
///////////////////////////////////////////////////////////////////////////////
1250
// //
1251
// Labels (enumeration declarations) used by TetGen. //
1252
// //
1253
///////////////////////////////////////////////////////////////////////////////
1254
1255
// Labels that signify the type of a vertex.
1256
enum verttype {UNUSEDVERTEX, DUPLICATEDVERTEX, RIDGEVERTEX, ACUTEVERTEX,
1257
FACETVERTEX, VOLVERTEX, FREESEGVERTEX, FREEFACETVERTEX,
1258
FREEVOLVERTEX, NREGULARVERTEX, DEADVERTEX};
1259
1260
// Labels that signify the result of triangle-triangle intersection test.
1261
enum interresult {DISJOINT, INTERSECT, SHAREVERT, SHAREEDGE, SHAREFACE,
1262
TOUCHEDGE, TOUCHFACE, ACROSSVERT, ACROSSEDGE, ACROSSFACE};
1263
1264
// Labels that signify the result of point location.
1265
enum locateresult {TGUNKNOWN, OUTSIDE, INTETRAHEDRON, ONFACE, ONEDGE, ONVERTEX,
1266
ENCVERTEX, ENCSEGMENT, ENCSUBFACE, NEARVERTEX, NONREGULAR,
1267
INSTAR, BADELEMENT};
1268
1269
///////////////////////////////////////////////////////////////////////////////
1270
// //
1271
// Variables of TetGen //
1272
// //
1273
///////////////////////////////////////////////////////////////////////////////
1274
1275
// Pointer to the input data (a set of nodes, a PLC, or a mesh).
1276
tetgenio *in, *addin;
1277
1278
// Pointer to the switches and parameters.
1279
tetgenbehavior *b;
1280
1281
// Pointer to a background mesh (contains size specification map).
1282
tetgenmesh *bgm;
1283
1284
// Memorypools to store mesh elements (points, tetrahedra, subfaces, and
1285
// segments) and extra pointers between tetrahedra, subfaces, and segments.
1286
memorypool *tetrahedrons, *subfaces, *subsegs, *points;
1287
memorypool *tet2subpool, *tet2segpool;
1288
1289
// Memorypools to store bad-quality (or encroached) elements.
1290
memorypool *badtetrahedrons, *badsubfacs, *badsubsegs;
1291
1292
// A memorypool to store faces to be flipped.
1293
memorypool *flippool;
1294
arraypool *unflipqueue;
1295
badface *flipstack;
1296
1297
// Arrays used for point insertion (the Bowyer-Watson algorithm).
1298
arraypool *cavetetlist, *cavebdrylist, *caveoldtetlist;
1299
arraypool *cavetetshlist, *cavetetseglist, *cavetetvertlist;
1300
arraypool *caveencshlist, *caveencseglist;
1301
arraypool *caveshlist, *caveshbdlist, *cavesegshlist;
1302
1303
// Stacks used for CDT construction and boundary recovery.
1304
arraypool *subsegstack, *subfacstack, *subvertstack;
1305
1306
// Arrays of encroached segments and subfaces (for mesh refinement).
1307
arraypool *encseglist, *encshlist;
1308
1309
// The map between facets to their vertices (for mesh refinement).
1310
int *idx2facetlist;
1311
point *facetverticeslist;
1312
1313
// The map between segments to their endpoints (for mesh refinement).
1314
point *segmentendpointslist;
1315
1316
// The infinite vertex.
1317
point dummypoint;
1318
// The recently visited tetrahedron, subface.
1319
triface recenttet;
1320
face recentsh;
1321
1322
// PI is the ratio of a circle's circumference to its diameter.
1323
static REAL PI;
1324
1325
// Array (size = numberoftetrahedra * 6) for storing high-order nodes of
1326
// tetrahedra (only used when -o2 switch is selected).
1327
point *highordertable;
1328
1329
// Various variables.
1330
int numpointattrib; // Number of point attributes.
1331
int numelemattrib; // Number of tetrahedron attributes.
1332
int sizeoftensor; // Number of REALs per metric tensor.
1333
int pointmtrindex; // Index to find the metric tensor of a point.
1334
int pointparamindex; // Index to find the u,v coordinates of a point.
1335
int point2simindex; // Index to find a simplex adjacent to a point.
1336
int pointmarkindex; // Index to find boundary marker of a point.
1337
int pointinsradiusindex; // Index to find the insertion radius of a point.
1338
int elemattribindex; // Index to find attributes of a tetrahedron.
1339
int volumeboundindex; // Index to find volume bound of a tetrahedron.
1340
int elemmarkerindex; // Index to find marker of a tetrahedron.
1341
int shmarkindex; // Index to find boundary marker of a subface.
1342
int areaboundindex; // Index to find area bound of a subface.
1343
int checksubsegflag; // Are there segments in the tetrahedralization yet?
1344
int checksubfaceflag; // Are there subfaces in the tetrahedralization yet?
1345
int checkconstraints; // Are there variant (node, seg, facet) constraints?
1346
int nonconvex; // Is current mesh non-convex?
1347
int autofliplinklevel; // The increase of link levels, default is 1.
1348
int useinsertradius; // Save the insertion radius for Steiner points.
1349
long samples; // Number of random samples for point location.
1350
unsigned long randomseed; // Current random number seed.
1351
REAL cosmaxdihed, cosmindihed; // The cosine values of max/min dihedral.
1352
REAL cossmtdihed; // The cosine value of a bad dihedral to be smoothed.
1353
REAL cosslidihed; // The cosine value of the max dihedral of a sliver.
1354
REAL minfaceang, minfacetdihed; // The minimum input (dihedral) angles.
1355
REAL tetprism_vol_sum; // The total volume of tetrahedral-prisms (in 4D).
1356
REAL longest; // The longest possible edge length.
1357
REAL minedgelength; // = longest * b->epsion.
1358
REAL xmax, xmin, ymax, ymin, zmax, zmin; // Bounding box of points.
1359
1360
// Counters.
1361
long insegments; // Number of input segments.
1362
long hullsize; // Number of exterior boundary faces.
1363
long meshedges; // Number of mesh edges.
1364
long meshhulledges; // Number of boundary mesh edges.
1365
long steinerleft; // Number of Steiner points not yet used.
1366
long dupverts; // Are there duplicated vertices?
1367
long unuverts; // Are there unused vertices?
1368
long nonregularcount; // Are there non-regular vertices?
1369
long st_segref_count, st_facref_count, st_volref_count; // Steiner points.
1370
long fillregioncount, cavitycount, cavityexpcount;
1371
long flip14count, flip26count, flipn2ncount;
1372
long flip23count, flip32count, flip44count, flip41count;
1373
long flip31count, flip22count;
1374
unsigned long totalworkmemory; // Total memory used by working arrays.
1375
1376
1377
///////////////////////////////////////////////////////////////////////////////
1378
// //
1379
// Mesh manipulation primitives //
1380
// //
1381
///////////////////////////////////////////////////////////////////////////////
1382
1383
// Fast lookup tables for mesh manipulation primitives.
1384
static int bondtbl[12][12], fsymtbl[12][12];
1385
static int esymtbl[12], enexttbl[12], eprevtbl[12];
1386
static int enextesymtbl[12], eprevesymtbl[12];
1387
static int eorgoppotbl[12], edestoppotbl[12];
1388
static int facepivot1[12], facepivot2[12][12];
1389
static int orgpivot[12], destpivot[12], apexpivot[12], oppopivot[12];
1390
static int tsbondtbl[12][6], stbondtbl[12][6];
1391
static int tspivottbl[12][6], stpivottbl[12][6];
1392
static int ver2edge[12], edge2ver[6], epivot[12];
1393
static int sorgpivot [6], sdestpivot[6], sapexpivot[6];
1394
static int snextpivot[6];
1395
1396
void inittables();
1397
1398
// Primitives for tetrahedra.
1399
inline tetrahedron encode(triface& t);
1400
inline tetrahedron encode2(tetrahedron* ptr, int ver);
1401
inline void decode(tetrahedron ptr, triface& t);
1402
inline void bond(triface& t1, triface& t2);
1403
inline void dissolve(triface& t);
1404
inline void esym(triface& t1, triface& t2);
1405
inline void esymself(triface& t);
1406
inline void enext(triface& t1, triface& t2);
1407
inline void enextself(triface& t);
1408
inline void eprev(triface& t1, triface& t2);
1409
inline void eprevself(triface& t);
1410
inline void enextesym(triface& t1, triface& t2);
1411
inline void enextesymself(triface& t);
1412
inline void eprevesym(triface& t1, triface& t2);
1413
inline void eprevesymself(triface& t);
1414
inline void eorgoppo(triface& t1, triface& t2);
1415
inline void eorgoppoself(triface& t);
1416
inline void edestoppo(triface& t1, triface& t2);
1417
inline void edestoppoself(triface& t);
1418
inline void fsym(triface& t1, triface& t2);
1419
inline void fsymself(triface& t);
1420
inline void fnext(triface& t1, triface& t2);
1421
inline void fnextself(triface& t);
1422
inline point org (triface& t);
1423
inline point dest(triface& t);
1424
inline point apex(triface& t);
1425
inline point oppo(triface& t);
1426
inline void setorg (triface& t, point p);
1427
inline void setdest(triface& t, point p);
1428
inline void setapex(triface& t, point p);
1429
inline void setoppo(triface& t, point p);
1430
inline REAL elemattribute(tetrahedron* ptr, int attnum);
1431
inline void setelemattribute(tetrahedron* ptr, int attnum, REAL value);
1432
inline REAL volumebound(tetrahedron* ptr);
1433
inline void setvolumebound(tetrahedron* ptr, REAL value);
1434
inline int elemindex(tetrahedron* ptr);
1435
inline void setelemindex(tetrahedron* ptr, int value);
1436
inline int elemmarker(tetrahedron* ptr);
1437
inline void setelemmarker(tetrahedron* ptr, int value);
1438
inline void infect(triface& t);
1439
inline void uninfect(triface& t);
1440
inline bool infected(triface& t);
1441
inline void marktest(triface& t);
1442
inline void unmarktest(triface& t);
1443
inline bool marktested(triface& t);
1444
inline void markface(triface& t);
1445
inline void unmarkface(triface& t);
1446
inline bool facemarked(triface& t);
1447
inline void markedge(triface& t);
1448
inline void unmarkedge(triface& t);
1449
inline bool edgemarked(triface& t);
1450
inline void marktest2(triface& t);
1451
inline void unmarktest2(triface& t);
1452
inline bool marktest2ed(triface& t);
1453
inline int elemcounter(triface& t);
1454
inline void setelemcounter(triface& t, int value);
1455
inline void increaseelemcounter(triface& t);
1456
inline void decreaseelemcounter(triface& t);
1457
inline bool ishulltet(triface& t);
1458
inline bool isdeadtet(triface& t);
1459
1460
// Primitives for subfaces and subsegments.
1461
inline void sdecode(shellface sptr, face& s);
1462
inline shellface sencode(face& s);
1463
inline shellface sencode2(shellface *sh, int shver);
1464
inline void spivot(face& s1, face& s2);
1465
inline void spivotself(face& s);
1466
inline void sbond(face& s1, face& s2);
1467
inline void sbond1(face& s1, face& s2);
1468
inline void sdissolve(face& s);
1469
inline point sorg(face& s);
1470
inline point sdest(face& s);
1471
inline point sapex(face& s);
1472
inline void setsorg(face& s, point pointptr);
1473
inline void setsdest(face& s, point pointptr);
1474
inline void setsapex(face& s, point pointptr);
1475
inline void sesym(face& s1, face& s2);
1476
inline void sesymself(face& s);
1477
inline void senext(face& s1, face& s2);
1478
inline void senextself(face& s);
1479
inline void senext2(face& s1, face& s2);
1480
inline void senext2self(face& s);
1481
inline REAL areabound(face& s);
1482
inline void setareabound(face& s, REAL value);
1483
inline int shellmark(face& s);
1484
inline void setshellmark(face& s, int value);
1485
inline void sinfect(face& s);
1486
inline void suninfect(face& s);
1487
inline bool sinfected(face& s);
1488
inline void smarktest(face& s);
1489
inline void sunmarktest(face& s);
1490
inline bool smarktested(face& s);
1491
inline void smarktest2(face& s);
1492
inline void sunmarktest2(face& s);
1493
inline bool smarktest2ed(face& s);
1494
inline void smarktest3(face& s);
1495
inline void sunmarktest3(face& s);
1496
inline bool smarktest3ed(face& s);
1497
inline void setfacetindex(face& f, int value);
1498
inline int getfacetindex(face& f);
1499
1500
// Primitives for interacting tetrahedra and subfaces.
1501
inline void tsbond(triface& t, face& s);
1502
inline void tsdissolve(triface& t);
1503
inline void stdissolve(face& s);
1504
inline void tspivot(triface& t, face& s);
1505
inline void stpivot(face& s, triface& t);
1506
1507
// Primitives for interacting tetrahedra and segments.
1508
inline void tssbond1(triface& t, face& seg);
1509
inline void sstbond1(face& s, triface& t);
1510
inline void tssdissolve1(triface& t);
1511
inline void sstdissolve1(face& s);
1512
inline void tsspivot1(triface& t, face& s);
1513
inline void sstpivot1(face& s, triface& t);
1514
1515
// Primitives for interacting subfaces and segments.
1516
inline void ssbond(face& s, face& edge);
1517
inline void ssbond1(face& s, face& edge);
1518
inline void ssdissolve(face& s);
1519
inline void sspivot(face& s, face& edge);
1520
1521
// Primitives for points.
1522
inline int pointmark(point pt);
1523
inline void setpointmark(point pt, int value);
1524
inline enum verttype pointtype(point pt);
1525
inline void setpointtype(point pt, enum verttype value);
1526
inline int pointgeomtag(point pt);
1527
inline void setpointgeomtag(point pt, int value);
1528
inline REAL pointgeomuv(point pt, int i);
1529
inline void setpointgeomuv(point pt, int i, REAL value);
1530
inline void pinfect(point pt);
1531
inline void puninfect(point pt);
1532
inline bool pinfected(point pt);
1533
inline void pmarktest(point pt);
1534
inline void punmarktest(point pt);
1535
inline bool pmarktested(point pt);
1536
inline void pmarktest2(point pt);
1537
inline void punmarktest2(point pt);
1538
inline bool pmarktest2ed(point pt);
1539
inline void pmarktest3(point pt);
1540
inline void punmarktest3(point pt);
1541
inline bool pmarktest3ed(point pt);
1542
inline tetrahedron point2tet(point pt);
1543
inline void setpoint2tet(point pt, tetrahedron value);
1544
inline shellface point2sh(point pt);
1545
inline void setpoint2sh(point pt, shellface value);
1546
inline point point2ppt(point pt);
1547
inline void setpoint2ppt(point pt, point value);
1548
inline tetrahedron point2bgmtet(point pt);
1549
inline void setpoint2bgmtet(point pt, tetrahedron value);
1550
inline void setpointinsradius(point pt, REAL value);
1551
inline REAL getpointinsradius(point pt);
1552
inline bool issteinerpoint(point pt);
1553
1554
// Advanced primitives.
1555
inline void point2tetorg(point pt, triface& t);
1556
inline void point2shorg(point pa, face& s);
1557
inline point farsorg(face& seg);
1558
inline point farsdest(face& seg);
1559
1560
///////////////////////////////////////////////////////////////////////////////
1561
// //
1562
// Memory management //
1563
// //
1564
///////////////////////////////////////////////////////////////////////////////
1565
1566
void tetrahedrondealloc(tetrahedron*);
1567
tetrahedron *tetrahedrontraverse();
1568
tetrahedron *alltetrahedrontraverse();
1569
void shellfacedealloc(memorypool*, shellface*);
1570
shellface *shellfacetraverse(memorypool*);
1571
void pointdealloc(point);
1572
point pointtraverse();
1573
1574
void makeindex2pointmap(point*&);
1575
void makepoint2submap(memorypool*, int*&, face*&);
1576
void maketetrahedron(triface*);
1577
void makeshellface(memorypool*, face*);
1578
void makepoint(point*, enum verttype);
1579
1580
void initializepools();
1581
1582
///////////////////////////////////////////////////////////////////////////////
1583
// //
1584
// Advanced geometric predicates and calculations //
1585
// //
1586
// TetGen uses a simplified symbolic perturbation scheme from Edelsbrunner, //
1587
// et al [*]. Hence the point-in-sphere test never returns a zero. The idea //
1588
// is to perturb the weights of vertices in the fourth dimension. TetGen //
1589
// uses the indices of the vertices decide the amount of perturbation. It is //
1590
// implemented in the routine insphere_s().
1591
// //
1592
// The routine tri_edge_test() determines whether or not a triangle and an //
1593
// edge intersect in 3D. If they intersect, their intersection type is also //
1594
// reported. This test is a combination of n 3D orientation tests (n is bet- //
1595
// ween 3 and 9). It uses the robust orient3d() test to make the branch dec- //
1596
// isions. The routine tri_tri_test() determines whether or not two triang- //
1597
// les intersect in 3D. It also uses the robust orient3d() test. //
1598
// //
1599
// There are a number of routines to calculate geometrical quantities, e.g., //
1600
// circumcenters, angles, dihedral angles, face normals, face areas, etc. //
1601
// They are so far done by the default floating-point arithmetics which are //
1602
// non-robust. They should be improved in the future. //
1603
// //
1604
///////////////////////////////////////////////////////////////////////////////
1605
1606
// Symbolic perturbations (robust)
1607
REAL insphere_s(REAL*, REAL*, REAL*, REAL*, REAL*);
1608
REAL orient4d_s(REAL*, REAL*, REAL*, REAL*, REAL*,
1609
REAL, REAL, REAL, REAL, REAL);
1610
1611
// Triangle-edge intersection test (robust)
1612
int tri_edge_2d(point, point, point, point, point, point, int, int*, int*);
1613
int tri_edge_tail(point, point, point, point, point, point, REAL, REAL, int,
1614
int*, int*);
1615
int tri_edge_test(point, point, point, point, point, point, int, int*, int*);
1616
1617
// Triangle-triangle intersection test (robust)
1618
int tri_edge_inter_tail(point, point, point, point, point, REAL, REAL);
1619
int tri_tri_inter(point, point, point, point, point, point);
1620
1621
// Linear algebra functions
1622
inline REAL dot(REAL* v1, REAL* v2);
1623
inline void cross(REAL* v1, REAL* v2, REAL* n);
1624
bool lu_decmp(REAL lu[4][4], int n, int* ps, REAL* d, int N);
1625
void lu_solve(REAL lu[4][4], int n, int* ps, REAL* b, int N);
1626
1627
// An embedded 2-dimensional geometric predicate (non-robust)
1628
REAL incircle3d(point pa, point pb, point pc, point pd);
1629
1630
// Geometric calculations (non-robust)
1631
REAL orient3dfast(REAL *pa, REAL *pb, REAL *pc, REAL *pd);
1632
inline REAL norm2(REAL x, REAL y, REAL z);
1633
inline REAL distance(REAL* p1, REAL* p2);
1634
void facenormal(point pa, point pb, point pc, REAL *n, int pivot, REAL *lav);
1635
REAL shortdistance(REAL* p, REAL* e1, REAL* e2);
1636
REAL triarea(REAL* pa, REAL* pb, REAL* pc);
1637
REAL interiorangle(REAL* o, REAL* p1, REAL* p2, REAL* n);
1638
void projpt2edge(REAL* p, REAL* e1, REAL* e2, REAL* prj);
1639
void projpt2face(REAL* p, REAL* f1, REAL* f2, REAL* f3, REAL* prj);
1640
bool tetalldihedral(point, point, point, point, REAL*, REAL*, REAL*);
1641
void tetallnormal(point, point, point, point, REAL N[4][3], REAL* volume);
1642
REAL tetaspectratio(point, point, point, point);
1643
bool circumsphere(REAL*, REAL*, REAL*, REAL*, REAL* cent, REAL* radius);
1644
bool orthosphere(REAL*,REAL*,REAL*,REAL*,REAL,REAL,REAL,REAL,REAL*,REAL*);
1645
void planelineint(REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*);
1646
int linelineint(REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*);
1647
REAL tetprismvol(REAL* pa, REAL* pb, REAL* pc, REAL* pd);
1648
bool calculateabovepoint(arraypool*, point*, point*, point*);
1649
void calculateabovepoint4(point, point, point, point);
1650
1651
// PLC error reports.
1652
void report_overlapping_facets(face*, face*, REAL dihedang = 0.0);
1653
int report_selfint_edge(point, point, face* sedge, triface* searchtet,
1654
enum interresult);
1655
int report_selfint_face(point, point, point, face* sface, triface* iedge,
1656
int intflag, int* types, int* poss);
1657
1658
///////////////////////////////////////////////////////////////////////////////
1659
// //
1660
// Local mesh transformations //
1661
// //
1662
// A local transformation replaces a small set of tetrahedra with another //
1663
// set of tetrahedra which fills the same space and the same boundaries. //
1664
// In 3D, the most simplest local transformations are the elementary flips //
1665
// performed within the convex hull of five vertices: 2-to-3, 3-to-2, 1-to-4,//
1666
// and 4-to-1 flips, where the numbers indicate the number of tetrahedra //
1667
// before and after each flip. The 1-to-4 and 4-to-1 flip involve inserting //
1668
// or deleting a vertex, respectively. //
1669
// There are complex local transformations which can be decomposed as a //
1670
// combination of elementary flips. For example,a 4-to-4 flip which replaces //
1671
// two coplanar edges can be regarded by a 2-to-3 flip and a 3-to-2 flip. //
1672
// Note that the first 2-to-3 flip will temporarily create a degenerate tet- //
1673
// rahedron which is removed immediately by the followed 3-to-2 flip. More //
1674
// generally, a n-to-m flip, where n > 3, m = (n - 2) * 2, which removes an //
1675
// edge can be done by first performing a sequence of (n - 3) 2-to-3 flips //
1676
// followed by a 3-to-2 flip. //
1677
// //
1678
// The routines flip23(), flip32(), and flip41() perform the three element- //
1679
// ray flips. The flip14() is available inside the routine insertpoint(). //
1680
// //
1681
// The routines flipnm() and flipnm_post() implement a generalized edge flip //
1682
// algorithm which uses a combination of elementary flips. //
1683
// //
1684
// The routine insertpoint() implements a variant of Bowyer-Watson's cavity //
1685
// algorithm to insert a vertex. It works for arbitrary tetrahedralization, //
1686
// either Delaunay, or constrained Delaunay, or non-Delaunay. //
1687
// //
1688
///////////////////////////////////////////////////////////////////////////////
1689
1690
// The elementary flips.
1691
void flip23(triface*, int, flipconstraints* fc);
1692
void flip32(triface*, int, flipconstraints* fc);
1693
void flip41(triface*, int, flipconstraints* fc);
1694
1695
// A generalized edge flip.
1696
int flipnm(triface*, int n, int level, int, flipconstraints* fc);
1697
int flipnm_post(triface*, int n, int nn, int, flipconstraints* fc);
1698
1699
// Point insertion.
1700
int insertpoint(point, triface*, face*, face*, insertvertexflags*);
1701
void insertpoint_abort(face*, insertvertexflags*);
1702
1703
///////////////////////////////////////////////////////////////////////////////
1704
// //
1705
// Delaunay tetrahedralization //
1706
// //
1707
// The routine incrementaldelaunay() implemented two incremental algorithms //
1708
// for constructing Delaunay tetrahedralizations (DTs): the Bowyer-Watson //
1709
// (B-W) algorithm and the incremental flip algorithm of Edelsbrunner and //
1710
// Shah, "Incremental topological flipping works for regular triangulation," //
1711
// Algorithmica, 15:233-241, 1996. //
1712
// //
1713
// The routine incrementalflip() implements the flip algorithm of [Edelsbru- //
1714
// nner and Shah, 1996]. It flips a queue of locally non-Delaunay faces (in //
1715
// an arbitrary order). The success is guaranteed when the Delaunay tetrah- //
1716
// edralization is constructed incrementally by adding one vertex at a time. //
1717
// //
1718
// The routine locate() finds a tetrahedron contains a new point in current //
1719
// DT. It uses a simple stochastic walk algorithm: starting from an arbitr- //
1720
// ary tetrahedron in DT, it finds the destination by visit one tetrahedron //
1721
// at a time, randomly chooses a tetrahedron if there are more than one //
1722
// choices. This algorithm terminates due to Edelsbrunner's acyclic theorem. //
1723
// Choose a good starting tetrahedron is crucial to the speed of the walk. //
1724
// TetGen originally uses the "jump-and-walk" algorithm of Muecke, E.P., //
1725
// Saias, I., and Zhu, B. "Fast Randomized Point Location Without Preproces- //
1726
// sing." In Proceedings of the 12th ACM Symposium on Computational Geometry,//
1727
// 274-283, 1996. It first randomly samples several tetrahedra in the DT //
1728
// and then choosing the closet one to start walking. //
1729
// The above algorithm slows download dramatically as the number of points //
1730
// grows -- reported in Amenta, N., Choi, S. and Rote, G., "Incremental //
1731
// construction con {BRIO}," In Proceedings of 19th ACM Symposium on //
1732
// Computational Geometry, 211-219, 2003. On the other hand, Liu and //
1733
// Snoeyink showed that the point location can be made in constant time if //
1734
// the points are pre-sorted so that the nearby points in space have nearby //
1735
// indices, then adding the points in this order. They sorted the points //
1736
// along the 3D Hilbert curve. //
1737
// //
1738
// The routine hilbert_sort3() sorts a set of 3D points along the 3D Hilbert //
1739
// curve. It recursively splits a point set according to the Hilbert indices //
1740
// mapped to the subboxes of the bounding box of the point set. //
1741
// The Hilbert indices is calculated by Butz's algorithm in 1971. A nice //
1742
// exposition of this algorithm can be found in the paper of Hamilton, C., //
1743
// "Compact Hilbert Indices", Technical Report CS-2006-07, Computer Science, //
1744
// Dalhousie University, 2006 (the Section 2). My implementation also refer- //
1745
// enced Steven Witham's implementation of "Hilbert walk" (hopefully, it is //
1746
// still available at: http://www.tiac.net/~sw/2008/10/Hilbert/). //
1747
// //
1748
// TetGen sorts the points using the method in the paper of Boissonnat,J.-D.,//
1749
// Devillers, O. and Hornus, S. "Incremental Construction of the Delaunay //
1750
// Triangulation and the Delaunay Graph in Medium Dimension," In Proceedings //
1751
// of the 25th ACM Symposium on Computational Geometry, 2009. //
1752
// It first randomly sorts the points into subgroups using the Biased Rand-//
1753
// omized Insertion Ordering (BRIO) of Amenta et al 2003, then sorts the //
1754
// points in each subgroup along the 3D Hilbert curve. Inserting points in //
1755
// this order ensures a randomized "sprinkling" of the points over the //
1756
// domain, while sorting of each subset ensures locality. //
1757
// //
1758
///////////////////////////////////////////////////////////////////////////////
1759
1760
void transfernodes();
1761
1762
// Point sorting.
1763
int transgc[8][3][8], tsb1mod3[8];
1764
void hilbert_init(int n);
1765
int hilbert_split(point* vertexarray, int arraysize, int gc0, int gc1,
1766
REAL, REAL, REAL, REAL, REAL, REAL);
1767
void hilbert_sort3(point* vertexarray, int arraysize, int e, int d,
1768
REAL, REAL, REAL, REAL, REAL, REAL, int depth);
1769
void brio_multiscale_sort(point*,int,int threshold,REAL ratio,int* depth);
1770
1771
// Point location.
1772
unsigned long randomnation(unsigned int choices);
1773
void randomsample(point searchpt, triface *searchtet);
1774
enum locateresult locate(point searchpt, triface *searchtet,
1775
int chkencflag = 0);
1776
1777
// Incremental flips.
1778
void flippush(badface*&, triface*);
1779
int incrementalflip(point newpt, int, flipconstraints *fc);
1780
1781
// Incremental Delaunay construction.
1782
void initialdelaunay(point pa, point pb, point pc, point pd);
1783
void incrementaldelaunay(clock_t&);
1784
1785
///////////////////////////////////////////////////////////////////////////////
1786
// //
1787
// Surface triangulation //
1788
// //
1789
///////////////////////////////////////////////////////////////////////////////
1790
1791
void flipshpush(face*);
1792
void flip22(face*, int, int);
1793
void flip31(face*, int);
1794
long lawsonflip();
1795
int sinsertvertex(point newpt, face*, face*, int iloc, int bowywat, int);
1796
int sremovevertex(point delpt, face*, face*, int lawson);
1797
1798
enum locateresult slocate(point, face*, int, int, int);
1799
enum interresult sscoutsegment(face*, point, int, int, int);
1800
void scarveholes(int, REAL*);
1801
int triangulate(int, arraypool*, arraypool*, int, REAL*);
1802
1803
void unifysegments();
1804
void identifyinputedges(point*);
1805
void mergefacets();
1806
void meshsurface();
1807
1808
void interecursive(shellface** subfacearray, int arraysize, int axis,
1809
REAL, REAL, REAL, REAL, REAL, REAL, int* internum);
1810
void detectinterfaces();
1811
1812
///////////////////////////////////////////////////////////////////////////////
1813
// //
1814
// Constrained Delaunay tetrahedralization //
1815
// //
1816
// A constrained Delaunay tetrahedralization (CDT) is a variation of a Dela- //
1817
// unay tetrahedralization (DT) that is constrained to respect the boundary //
1818
// of a 3D PLC (domain). In a CDT of a 3D PLC, every vertex or edge of the //
1819
// PLC is also a vertex or an edge of the CDT, every polygon of the PLC is a //
1820
// union of triangles of the CDT. A crucial difference between a CDT and a //
1821
// DT is that triangles in the PLC's polygons are not required to be locally //
1822
// Delaunay, which frees the CDT to better respect the PLC's polygons. CDTs //
1823
// have optimal properties similar to those of DTs. //
1824
// //
1825
// Steiner Points and Steiner CDTs. It is known that even a simple 3D polyh- //
1826
// edron may not have a tetrahedralization which only uses its own vertices. //
1827
// Some extra points, so-called "Steiner points" are needed in order to form //
1828
// a tetrahedralization of such polyhedron. It is true for tetrahedralizing //
1829
// a 3D PLC as well. A Steiner CDT of a 3D PLC is a CDT containing Steiner //
1830
// points. The CDT algorithms of TetGen in general create Steiner CDTs. //
1831
// Almost all of the Steiner points are added in the edges of the PLC. They //
1832
// guarantee the existence of a CDT of the modified PLC. //
1833
// //
1834
// The routine constraineddelaunay() starts from a DT of the vertices of a //
1835
// PLC and creates a (Steiner) CDT of the PLC (including Steiner points). It //
1836
// is constructed by two steps, (1) segment recovery and (2) facet (polygon) //
1837
// recovery. Each step is accomplished by its own algorithm. //
1838
// //
1839
// The routine delaunizesegments() implements the segment recovery algorithm //
1840
// of Si, H. and Gaertner, K. "Meshing Piecewise Linear Complexes by Constr- //
1841
// ained Delaunay Tetrahedralizations," In Proceedings of the 14th Internat- //
1842
// ional Meshing Roundtable, 147--163, 2005. It adds Steiner points into //
1843
// non-Delaunay segments until all subsegments appear together in a DT. The //
1844
// running time of this algorithm is proportional to the number of added //
1845
// Steiner points. //
1846
// //
1847
// There are two incremental facet recovery algorithms: the cavity re-trian- //
1848
// gulation algorithm of Si, H. and Gaertner, K. "3D Boundary Recovery by //
1849
// Constrained Delaunay Tetrahedralization," International Journal for Numer-//
1850
// ical Methods in Engineering, 85:1341-1364, 2011, and the flip algorithm //
1851
// of Shewchuk, J. "Updating and Constructing Constrained Delaunay and //
1852
// Constrained Regular Triangulations by Flips." In Proceedings of the 19th //
1853
// ACM Symposium on Computational Geometry, 86-95, 2003. //
1854
// //
1855
// It is guaranteed in theory, no Steiner point is needed in both algorithms //
1856
// However, a facet with non-coplanar vertices might cause the additions of //
1857
// Steiner points. It is discussed in the paper of Si, H., and Shewchuk, J.,//
1858
// "Incrementally Constructing and Updating Constrained Delaunay //
1859
// Tetrahedralizations with Finite Precision Coordinates." In Proceedings of //
1860
// the 21th International Meshing Roundtable, 2012. //
1861
// //
1862
// Our implementation of the facet recovery algorithms recover a "missing //
1863
// region" at a time. Each missing region is a subset of connected interiors //
1864
// of a polygon. The routine formcavity() creates the cavity of crossing //
1865
// tetrahedra of the missing region. //
1866
// //
1867
// The cavity re-triangulation algorithm is implemented by three subroutines,//
1868
// delaunizecavity(), fillcavity(), and carvecavity(). Since it may fail due //
1869
// to non-coplanar vertices, the subroutine restorecavity() is used to rest- //
1870
// ore the original cavity. //
1871
// //
1872
// The routine flipinsertfacet() implements the flip algorithm. The subrout- //
1873
// ine flipcertify() is used to maintain the priority queue of flips. //
1874
// //
1875
// The routine refineregion() is called when the facet recovery algorithm //
1876
// fail to recover a missing region. It inserts Steiner points to refine the //
1877
// missing region. In order to avoid inserting Steiner points very close to //
1878
// existing segments. The classical encroachment rules of the Delaunay //
1879
// refinement algorithm are used to choose the Steiner points. //
1880
// //
1881
// The routine constrainedfacets() does the facet recovery by using either //
1882
// the cavity re-triangulation algorithm (default) or the flip algorithm. It //
1883
// results a CDT of the (modified) PLC (including Steiner points). //
1884
// //
1885
///////////////////////////////////////////////////////////////////////////////
1886
1887
void makesegmentendpointsmap();
1888
1889
enum interresult finddirection(triface* searchtet, point endpt);
1890
enum interresult scoutsegment(point, point, face*, triface*, point*,
1891
arraypool*);
1892
int getsteinerptonsegment(face* seg, point refpt, point steinpt);
1893
void delaunizesegments();
1894
1895
int scoutsubface(face* searchsh,triface* searchtet,int shflag);
1896
void formregion(face*, arraypool*, arraypool*, arraypool*);
1897
int scoutcrossedge(triface& crosstet, arraypool*, arraypool*);
1898
bool formcavity(triface*, arraypool*, arraypool*, arraypool*, arraypool*,
1899
arraypool*, arraypool*);
1900
// Facet recovery by cavity re-triangulation [Si and Gaertner 2011].
1901
void delaunizecavity(arraypool*, arraypool*, arraypool*, arraypool*,
1902
arraypool*, arraypool*);
1903
bool fillcavity(arraypool*, arraypool*, arraypool*, arraypool*,
1904
arraypool*, arraypool*, triface* crossedge);
1905
void carvecavity(arraypool*, arraypool*, arraypool*);
1906
void restorecavity(arraypool*, arraypool*, arraypool*, arraypool*);
1907
// Facet recovery by flips [Shewchuk 2003].
1908
void flipcertify(triface *chkface, badface **pqueue, point, point, point);
1909
void flipinsertfacet(arraypool*, arraypool*, arraypool*, arraypool*);
1910
1911
int insertpoint_cdt(point, triface*, face*, face*, insertvertexflags*,
1912
arraypool*, arraypool*, arraypool*, arraypool*,
1913
arraypool*, arraypool*);
1914
void refineregion(face&, arraypool*, arraypool*, arraypool*, arraypool*,
1915
arraypool*, arraypool*);
1916
void constrainedfacets();
1917
1918
void constraineddelaunay(clock_t&);
1919
1920
///////////////////////////////////////////////////////////////////////////////
1921
// //
1922
// Constrained tetrahedralizations. //
1923
// //
1924
///////////////////////////////////////////////////////////////////////////////
1925
1926
int checkflipeligibility(int fliptype, point, point, point, point, point,
1927
int level, int edgepivot, flipconstraints* fc);
1928
1929
int removeedgebyflips(triface*, flipconstraints*);
1930
int removefacebyflips(triface*, flipconstraints*);
1931
1932
int recoveredgebyflips(point, point, face*, triface*, int fullsearch);
1933
int add_steinerpt_in_schoenhardtpoly(triface*, int, int chkencflag);
1934
int add_steinerpt_in_segment(face*, int searchlevel);
1935
int addsteiner4recoversegment(face*, int);
1936
int recoversegments(arraypool*, int fullsearch, int steinerflag);
1937
1938
int recoverfacebyflips(point, point, point, face*, triface*);
1939
int recoversubfaces(arraypool*, int steinerflag);
1940
1941
int getvertexstar(int, point searchpt, arraypool*, arraypool*, arraypool*);
1942
int getedge(point, point, triface*);
1943
int reduceedgesatvertex(point startpt, arraypool* endptlist);
1944
int removevertexbyflips(point steinerpt);
1945
1946
int suppressbdrysteinerpoint(point steinerpt);
1947
int suppresssteinerpoints();
1948
1949
void recoverboundary(clock_t&);
1950
1951
///////////////////////////////////////////////////////////////////////////////
1952
// //
1953
// Mesh reconstruction //
1954
// //
1955
///////////////////////////////////////////////////////////////////////////////
1956
1957
void carveholes();
1958
1959
void reconstructmesh();
1960
1961
int scoutpoint(point, triface*, int randflag);
1962
REAL getpointmeshsize(point, triface*, int iloc);
1963
void interpolatemeshsize();
1964
1965
void insertconstrainedpoints(point *insertarray, int arylen, int rejflag);
1966
void insertconstrainedpoints(tetgenio *addio);
1967
1968
void collectremovepoints(arraypool *remptlist);
1969
void meshcoarsening();
1970
1971
///////////////////////////////////////////////////////////////////////////////
1972
// //
1973
// Mesh refinement //
1974
// //
1975
// The purpose of mesh refinement is to obtain a tetrahedral mesh with well- //
1976
// -shaped tetrahedra and appropriate mesh size. It is necessary to insert //
1977
// new Steiner points to achieve this property. The questions are (1) how to //
1978
// choose the Steiner points? and (2) how to insert them? //
1979
// //
1980
// Delaunay refinement is a technique first developed by Chew [1989] and //
1981
// Ruppert [1993, 1995] to generate quality triangular meshes in the plane. //
1982
// It provides guarantee on the smallest angle of the triangles. Rupper's //
1983
// algorithm guarantees that the mesh is size-optimal (to within a constant //
1984
// factor) among all meshes with the same quality. //
1985
// Shewchuk generalized Ruppert's algorithm into 3D in his PhD thesis //
1986
// [Shewchuk 1997]. A short version of his algorithm appears in "Tetrahedral //
1987
// Mesh Generation by Delaunay Refinement," In Proceedings of the 14th ACM //
1988
// Symposium on Computational Geometry, 86-95, 1998. It guarantees that all //
1989
// tetrahedra of the output mesh have a "radius-edge ratio" (equivalent to //
1990
// the minimal face angle) bounded. However, it does not remove slivers, a //
1991
// type of very flat tetrahedra which can have no small face angles but have //
1992
// very small (and large) dihedral angles. Moreover, it may not terminate if //
1993
// the input PLC contains "sharp features", e.g., two edges (or two facets) //
1994
// meet at an acute angle (or dihedral angle). //
1995
// //
1996
// TetGen uses the basic Delaunay refinement scheme to insert Steiner points.//
1997
// While it always maintains a constrained Delaunay mesh. The algorithm is //
1998
// described in Si, H., "Adaptive Constrained Delaunay Mesh Generation," //
1999
// International Journal for Numerical Methods in Engineering, 75:856-880. //
2000
// This algorithm always terminates and sharp features are easily preserved. //
2001
// The mesh has good quality (same as Shewchuk's Delaunay refinement algori- //
2002
// thm) in the bulk of the mesh domain. Moreover, it supports the generation //
2003
// of adaptive mesh according to a (isotropic) mesh sizing function. //
2004
// //
2005
///////////////////////////////////////////////////////////////////////////////
2006
2007
void makefacetverticesmap();
2008
int segsegadjacent(face *, face *);
2009
int segfacetadjacent(face *checkseg, face *checksh);
2010
int facetfacetadjacent(face *, face *);
2011
void save_segmentpoint_insradius(point segpt, point parentpt, REAL r);
2012
void save_facetpoint_insradius(point facpt, point parentpt, REAL r);
2013
void enqueuesubface(memorypool*, face*);
2014
void enqueuetetrahedron(triface*);
2015
2016
int checkseg4encroach(point pa, point pb, point checkpt);
2017
int checkseg4split(face *chkseg, point&, int&);
2018
int splitsegment(face *splitseg, point encpt, REAL, point, point, int, int);
2019
void repairencsegs(int chkencflag);
2020
2021
int checkfac4encroach(point, point, point, point checkpt, REAL*, REAL*);
2022
int checkfac4split(face *chkfac, point& encpt, int& qflag, REAL *ccent);
2023
int splitsubface(face *splitfac, point, point, int qflag, REAL *ccent, int);
2024
void repairencfacs(int chkencflag);
2025
2026
int checktet4split(triface *chktet, int& qflag, REAL *ccent);
2027
int splittetrahedron(triface* splittet,int qflag,REAL *ccent, int);
2028
void repairbadtets(int chkencflag);
2029
2030
void delaunayrefinement();
2031
2032
///////////////////////////////////////////////////////////////////////////////
2033
// //
2034
// Mesh optimization //
2035
// //
2036
///////////////////////////////////////////////////////////////////////////////
2037
2038
long lawsonflip3d(flipconstraints *fc);
2039
void recoverdelaunay();
2040
2041
int gettetrahedron(point, point, point, point, triface *);
2042
long improvequalitybyflips();
2043
2044
int smoothpoint(point smtpt, arraypool*, int ccw, optparameters *opm);
2045
long improvequalitybysmoothing(optparameters *opm);
2046
2047
int splitsliver(triface *, REAL, int);
2048
long removeslivers(int);
2049
2050
void optimizemesh();
2051
2052
///////////////////////////////////////////////////////////////////////////////
2053
// //
2054
// Mesh check and statistics //
2055
// //
2056
///////////////////////////////////////////////////////////////////////////////
2057
2058
// Mesh validations.
2059
int checkmesh(int topoflag);
2060
int checkshells();
2061
int checksegments();
2062
int checkdelaunay(int perturb = 1);
2063
int checkregular(int);
2064
int checkconforming(int);
2065
2066
// Mesh statistics.
2067
void printfcomma(unsigned long n);
2068
void qualitystatistics();
2069
void memorystatistics();
2070
void statistics();
2071
2072
///////////////////////////////////////////////////////////////////////////////
2073
// //
2074
// Mesh output //
2075
// //
2076
///////////////////////////////////////////////////////////////////////////////
2077
2078
void jettisonnodes();
2079
void highorder();
2080
void indexelements();
2081
void numberedges();
2082
void outnodes(tetgenio*);
2083
void outmetrics(tetgenio*);
2084
void outelements(tetgenio*);
2085
void outfaces(tetgenio*);
2086
void outhullfaces(tetgenio*);
2087
void outsubfaces(tetgenio*);
2088
void outedges(tetgenio*);
2089
void outsubsegments(tetgenio*);
2090
void outneighbors(tetgenio*);
2091
void outvoronoi(tetgenio*);
2092
void outsmesh(char*);
2093
void outmesh2medit(char*);
2094
void outmesh2vtk(char*);
2095
2096
2097
2098
///////////////////////////////////////////////////////////////////////////////
2099
// //
2100
// Constructor & destructor //
2101
// //
2102
///////////////////////////////////////////////////////////////////////////////
2103
2104
void initializetetgenmesh()
2105
{
2106
in = addin = NULL;
2107
b = NULL;
2108
bgm = NULL;
2109
2110
tetrahedrons = subfaces = subsegs = points = NULL;
2111
badtetrahedrons = badsubfacs = badsubsegs = NULL;
2112
tet2segpool = tet2subpool = NULL;
2113
flippool = NULL;
2114
2115
dummypoint = NULL;
2116
flipstack = NULL;
2117
unflipqueue = NULL;
2118
2119
cavetetlist = cavebdrylist = caveoldtetlist = NULL;
2120
cavetetshlist = cavetetseglist = cavetetvertlist = NULL;
2121
caveencshlist = caveencseglist = NULL;
2122
caveshlist = caveshbdlist = cavesegshlist = NULL;
2123
2124
subsegstack = subfacstack = subvertstack = NULL;
2125
encseglist = encshlist = NULL;
2126
idx2facetlist = NULL;
2127
facetverticeslist = NULL;
2128
segmentendpointslist = NULL;
2129
2130
highordertable = NULL;
2131
2132
numpointattrib = numelemattrib = 0;
2133
sizeoftensor = 0;
2134
pointmtrindex = 0;
2135
pointparamindex = 0;
2136
pointmarkindex = 0;
2137
point2simindex = 0;
2138
pointinsradiusindex = 0;
2139
elemattribindex = 0;
2140
volumeboundindex = 0;
2141
shmarkindex = 0;
2142
areaboundindex = 0;
2143
checksubsegflag = 0;
2144
checksubfaceflag = 0;
2145
checkconstraints = 0;
2146
nonconvex = 0;
2147
autofliplinklevel = 1;
2148
useinsertradius = 0;
2149
samples = 0l;
2150
randomseed = 1l;
2151
minfaceang = minfacetdihed = PI;
2152
tetprism_vol_sum = 0.0;
2153
longest = minedgelength = 0.0;
2154
xmax = xmin = ymax = ymin = zmax = zmin = 0.0;
2155
2156
insegments = 0l;
2157
hullsize = 0l;
2158
meshedges = meshhulledges = 0l;
2159
steinerleft = -1;
2160
dupverts = 0l;
2161
unuverts = 0l;
2162
nonregularcount = 0l;
2163
st_segref_count = st_facref_count = st_volref_count = 0l;
2164
fillregioncount = cavitycount = cavityexpcount = 0l;
2165
flip14count = flip26count = flipn2ncount = 0l;
2166
flip23count = flip32count = flip44count = flip41count = 0l;
2167
flip22count = flip31count = 0l;
2168
totalworkmemory = 0l;
2169
2170
2171
} // tetgenmesh()
2172
2173
void freememory()
2174
{
2175
if (bgm != NULL) {
2176
delete bgm;
2177
}
2178
2179
if (points != (memorypool *) NULL) {
2180
delete points;
2181
delete [] dummypoint;
2182
}
2183
if (tetrahedrons != (memorypool *) NULL) {
2184
delete tetrahedrons;
2185
}
2186
if (subfaces != (memorypool *) NULL) {
2187
delete subfaces;
2188
delete subsegs;
2189
}
2190
if (tet2segpool != NULL) {
2191
delete tet2segpool;
2192
delete tet2subpool;
2193
}
2194
2195
if (badtetrahedrons) {
2196
delete badtetrahedrons;
2197
}
2198
if (badsubfacs) {
2199
delete badsubfacs;
2200
}
2201
if (badsubsegs) {
2202
delete badsubsegs;
2203
}
2204
if (encseglist) {
2205
delete encseglist;
2206
}
2207
if (encshlist) {
2208
delete encshlist;
2209
}
2210
2211
if (flippool != NULL) {
2212
delete flippool;
2213
delete unflipqueue;
2214
}
2215
2216
if (cavetetlist != NULL) {
2217
delete cavetetlist;
2218
delete cavebdrylist;
2219
delete caveoldtetlist;
2220
delete cavetetvertlist;
2221
}
2222
2223
if (caveshlist != NULL) {
2224
delete caveshlist;
2225
delete caveshbdlist;
2226
delete cavesegshlist;
2227
delete cavetetshlist;
2228
delete cavetetseglist;
2229
delete caveencshlist;
2230
delete caveencseglist;
2231
}
2232
2233
if (subsegstack != NULL) {
2234
delete subsegstack;
2235
delete subfacstack;
2236
delete subvertstack;
2237
}
2238
2239
if (idx2facetlist != NULL) {
2240
delete [] idx2facetlist;
2241
delete [] facetverticeslist;
2242
}
2243
2244
if (segmentendpointslist != NULL) {
2245
delete [] segmentendpointslist;
2246
}
2247
2248
if (highordertable != NULL) {
2249
delete [] highordertable;
2250
}
2251
2252
initializetetgenmesh();
2253
}
2254
2255
tetgenmesh()
2256
{
2257
initializetetgenmesh();
2258
}
2259
2260
~tetgenmesh()
2261
{
2262
freememory();
2263
} // ~tetgenmesh()
2264
2265
}; // End of class tetgenmesh.
2266
2267
///////////////////////////////////////////////////////////////////////////////
2268
// //
2269
// tetrahedralize() Interface for using TetGen's library to generate //
2270
// Delaunay tetrahedralizations, constrained Delaunay //
2271
// tetrahedralizations, quality tetrahedral meshes. //
2272
// //
2273
// 'in' is an object of 'tetgenio' which contains a PLC you want to tetrahed-//
2274
// ralize or a previously generated tetrahedral mesh you want to refine. It //
2275
// must not be a NULL. 'out' is another object of 'tetgenio' for storing the //
2276
// generated tetrahedral mesh. It can be a NULL. If so, the output will be //
2277
// saved to file(s). If 'bgmin' != NULL, it contains a background mesh which //
2278
// defines a mesh size function. //
2279
// //
2280
///////////////////////////////////////////////////////////////////////////////
2281
2282
void tetrahedralize(tetgenbehavior *b, tetgenio *in, tetgenio *out,
2283
tetgenio *addin = NULL, tetgenio *bgmin = NULL);
2284
2285
#ifdef TETLIBRARY
2286
void tetrahedralize(char *switches, tetgenio *in, tetgenio *out,
2287
tetgenio *addin = NULL, tetgenio *bgmin = NULL);
2288
2289
extern "C"
2290
#ifdef WIN32
2291
__declspec(dllexport)
2292
#endif
2293
void delegate_tetrahedralize(int bs, tetgenbehavior *b, char *switches,
2294
tetgenio *in, tetgenio *out, tetgenio *addin = NULL, tetgenio *bgmin = NULL);
2295
2296
#endif // #ifdef TETLIBRARY
2297
2298
///////////////////////////////////////////////////////////////////////////////
2299
// //
2300
// terminatetetgen() Terminate TetGen with a given exit code. //
2301
// //
2302
///////////////////////////////////////////////////////////////////////////////
2303
2304
inline void terminatetetgen(tetgenmesh *m, int x)
2305
{
2306
#ifdef TETLIBRARY
2307
throw x;
2308
#else
2309
switch (x) {
2310
case 1: // Out of memory.
2311
printf("Error: Out of memory.\n");
2312
break;
2313
case 2: // Encounter an internal error.
2314
printf("Please report this bug to [email protected]. Include\n");
2315
printf(" the message above, your input data set, and the exact\n");
2316
printf(" command line you used to run this program, thank you.\n");
2317
break;
2318
case 3:
2319
printf("A self-intersection was detected. Program stopped.\n");
2320
printf("Hint: use -d option to detect all self-intersections.\n");
2321
break;
2322
case 4:
2323
printf("A very small input feature size was detected. Program stopped.\n");
2324
if (m) {
2325
printf("Hint: use -T option to set a smaller tolerance. Current is %g\n",
2326
m->b->epsilon);
2327
}
2328
break;
2329
case 5:
2330
printf("Two very close input facets were detected. Program stopped.\n");
2331
printf("Hint: use -Y option to avoid adding Steiner points in boundary.\n");
2332
break;
2333
case 10:
2334
printf("An input error was detected. Program stopped.\n");
2335
break;
2336
} // switch (x)
2337
exit(x);
2338
#endif // #ifdef TETLIBRARY
2339
}
2340
2341
///////////////////////////////////////////////////////////////////////////////
2342
// //
2343
// Primitives for tetrahedra //
2344
// //
2345
///////////////////////////////////////////////////////////////////////////////
2346
2347
// encode() compress a handle into a single pointer. It relies on the
2348
// assumption that all addresses of tetrahedra are aligned to sixteen-
2349
// byte boundaries, so that the last four significant bits are zero.
2350
2351
inline tetgenmesh::tetrahedron tetgenmesh::encode(triface& t) {
2352
return (tetrahedron) ((uintptr_t) (t).tet | (uintptr_t) (t).ver);
2353
}
2354
2355
inline tetgenmesh::tetrahedron tetgenmesh::encode2(tetrahedron* ptr, int ver) {
2356
return (tetrahedron) ((uintptr_t) (ptr) | (uintptr_t) (ver));
2357
}
2358
2359
// decode() converts a pointer to a handle. The version is extracted from
2360
// the four least significant bits of the pointer.
2361
2362
inline void tetgenmesh::decode(tetrahedron ptr, triface& t) {
2363
(t).ver = (int) ((uintptr_t) (ptr) & (uintptr_t) 15);
2364
(t).tet = (tetrahedron *) ((uintptr_t) (ptr) ^ (uintptr_t) (t).ver);
2365
}
2366
2367
// bond() connects two tetrahedra together. (t1,v1) and (t2,v2) must
2368
// refer to the same face and the same edge.
2369
2370
inline void tetgenmesh::bond(triface& t1, triface& t2) {
2371
t1.tet[t1.ver & 3] = encode2(t2.tet, bondtbl[t1.ver][t2.ver]);
2372
t2.tet[t2.ver & 3] = encode2(t1.tet, bondtbl[t2.ver][t1.ver]);
2373
}
2374
2375
2376
// dissolve() a bond (from one side).
2377
2378
inline void tetgenmesh::dissolve(triface& t) {
2379
t.tet[t.ver & 3] = NULL;
2380
}
2381
2382
// enext() finds the next edge (counterclockwise) in the same face.
2383
2384
inline void tetgenmesh::enext(triface& t1, triface& t2) {
2385
t2.tet = t1.tet;
2386
t2.ver = enexttbl[t1.ver];
2387
}
2388
2389
inline void tetgenmesh::enextself(triface& t) {
2390
t.ver = enexttbl[t.ver];
2391
}
2392
2393
// eprev() finds the next edge (clockwise) in the same face.
2394
2395
inline void tetgenmesh::eprev(triface& t1, triface& t2) {
2396
t2.tet = t1.tet;
2397
t2.ver = eprevtbl[t1.ver];
2398
}
2399
2400
inline void tetgenmesh::eprevself(triface& t) {
2401
t.ver = eprevtbl[t.ver];
2402
}
2403
2404
// esym() finds the reversed edge. It is in the other face of the
2405
// same tetrahedron.
2406
2407
inline void tetgenmesh::esym(triface& t1, triface& t2) {
2408
(t2).tet = (t1).tet;
2409
(t2).ver = esymtbl[(t1).ver];
2410
}
2411
2412
inline void tetgenmesh::esymself(triface& t) {
2413
(t).ver = esymtbl[(t).ver];
2414
}
2415
2416
// enextesym() finds the reversed edge of the next edge. It is in the other
2417
// face of the same tetrahedron. It is the combination esym() * enext().
2418
2419
inline void tetgenmesh::enextesym(triface& t1, triface& t2) {
2420
t2.tet = t1.tet;
2421
t2.ver = enextesymtbl[t1.ver];
2422
}
2423
2424
inline void tetgenmesh::enextesymself(triface& t) {
2425
t.ver = enextesymtbl[t.ver];
2426
}
2427
2428
// eprevesym() finds the reversed edge of the previous edge.
2429
2430
inline void tetgenmesh::eprevesym(triface& t1, triface& t2) {
2431
t2.tet = t1.tet;
2432
t2.ver = eprevesymtbl[t1.ver];
2433
}
2434
2435
inline void tetgenmesh::eprevesymself(triface& t) {
2436
t.ver = eprevesymtbl[t.ver];
2437
}
2438
2439
// eorgoppo() Finds the opposite face of the origin of the current edge.
2440
// Return the opposite edge of the current edge.
2441
2442
inline void tetgenmesh::eorgoppo(triface& t1, triface& t2) {
2443
t2.tet = t1.tet;
2444
t2.ver = eorgoppotbl[t1.ver];
2445
}
2446
2447
inline void tetgenmesh::eorgoppoself(triface& t) {
2448
t.ver = eorgoppotbl[t.ver];
2449
}
2450
2451
// edestoppo() Finds the opposite face of the destination of the current
2452
// edge. Return the opposite edge of the current edge.
2453
2454
inline void tetgenmesh::edestoppo(triface& t1, triface& t2) {
2455
t2.tet = t1.tet;
2456
t2.ver = edestoppotbl[t1.ver];
2457
}
2458
2459
inline void tetgenmesh::edestoppoself(triface& t) {
2460
t.ver = edestoppotbl[t.ver];
2461
}
2462
2463
// fsym() finds the adjacent tetrahedron at the same face and the same edge.
2464
2465
inline void tetgenmesh::fsym(triface& t1, triface& t2) {
2466
decode((t1).tet[(t1).ver & 3], t2);
2467
t2.ver = fsymtbl[t1.ver][t2.ver];
2468
}
2469
2470
2471
#define fsymself(t) \
2472
t1ver = (t).ver; \
2473
decode((t).tet[(t).ver & 3], (t));\
2474
(t).ver = fsymtbl[t1ver][(t).ver]
2475
2476
// fnext() finds the next face while rotating about an edge according to
2477
// a right-hand rule. The face is in the adjacent tetrahedron. It is
2478
// the combination: fsym() * esym().
2479
2480
inline void tetgenmesh::fnext(triface& t1, triface& t2) {
2481
decode(t1.tet[facepivot1[t1.ver]], t2);
2482
t2.ver = facepivot2[t1.ver][t2.ver];
2483
}
2484
2485
2486
#define fnextself(t) \
2487
t1ver = (t).ver; \
2488
decode((t).tet[facepivot1[(t).ver]], (t)); \
2489
(t).ver = facepivot2[t1ver][(t).ver]
2490
2491
2492
// The following primtives get or set the origin, destination, face apex,
2493
// or face opposite of an ordered tetrahedron.
2494
2495
inline tetgenmesh::point tetgenmesh::org(triface& t) {
2496
return (point) (t).tet[orgpivot[(t).ver]];
2497
}
2498
2499
inline tetgenmesh::point tetgenmesh:: dest(triface& t) {
2500
return (point) (t).tet[destpivot[(t).ver]];
2501
}
2502
2503
inline tetgenmesh::point tetgenmesh:: apex(triface& t) {
2504
return (point) (t).tet[apexpivot[(t).ver]];
2505
}
2506
2507
inline tetgenmesh::point tetgenmesh:: oppo(triface& t) {
2508
return (point) (t).tet[oppopivot[(t).ver]];
2509
}
2510
2511
inline void tetgenmesh:: setorg(triface& t, point p) {
2512
(t).tet[orgpivot[(t).ver]] = (tetrahedron) (p);
2513
}
2514
2515
inline void tetgenmesh:: setdest(triface& t, point p) {
2516
(t).tet[destpivot[(t).ver]] = (tetrahedron) (p);
2517
}
2518
2519
inline void tetgenmesh:: setapex(triface& t, point p) {
2520
(t).tet[apexpivot[(t).ver]] = (tetrahedron) (p);
2521
}
2522
2523
inline void tetgenmesh:: setoppo(triface& t, point p) {
2524
(t).tet[oppopivot[(t).ver]] = (tetrahedron) (p);
2525
}
2526
2527
#define setvertices(t, torg, tdest, tapex, toppo) \
2528
(t).tet[orgpivot[(t).ver]] = (tetrahedron) (torg);\
2529
(t).tet[destpivot[(t).ver]] = (tetrahedron) (tdest); \
2530
(t).tet[apexpivot[(t).ver]] = (tetrahedron) (tapex); \
2531
(t).tet[oppopivot[(t).ver]] = (tetrahedron) (toppo)
2532
2533
// Check or set a tetrahedron's attributes.
2534
2535
inline REAL tetgenmesh::elemattribute(tetrahedron* ptr, int attnum) {
2536
return ((REAL *) (ptr))[elemattribindex + attnum];
2537
}
2538
2539
inline void tetgenmesh::setelemattribute(tetrahedron* ptr, int attnum,
2540
REAL value) {
2541
((REAL *) (ptr))[elemattribindex + attnum] = value;
2542
}
2543
2544
// Check or set a tetrahedron's maximum volume bound.
2545
2546
inline REAL tetgenmesh::volumebound(tetrahedron* ptr) {
2547
return ((REAL *) (ptr))[volumeboundindex];
2548
}
2549
2550
inline void tetgenmesh::setvolumebound(tetrahedron* ptr, REAL value) {
2551
((REAL *) (ptr))[volumeboundindex] = value;
2552
}
2553
2554
// Get or set a tetrahedron's index (only used for output).
2555
// These two routines use the reserved slot ptr[10].
2556
2557
inline int tetgenmesh::elemindex(tetrahedron* ptr) {
2558
int *iptr = (int *) &(ptr[10]);
2559
return iptr[0];
2560
}
2561
2562
inline void tetgenmesh::setelemindex(tetrahedron* ptr, int value) {
2563
int *iptr = (int *) &(ptr[10]);
2564
iptr[0] = value;
2565
}
2566
2567
// Get or set a tetrahedron's marker.
2568
// Set 'value = 0' cleans all the face/edge flags.
2569
2570
inline int tetgenmesh::elemmarker(tetrahedron* ptr) {
2571
return ((int *) (ptr))[elemmarkerindex];
2572
}
2573
2574
inline void tetgenmesh::setelemmarker(tetrahedron* ptr, int value) {
2575
((int *) (ptr))[elemmarkerindex] = value;
2576
}
2577
2578
// infect(), infected(), uninfect() -- primitives to flag or unflag a
2579
// tetrahedron. The last bit of the element marker is flagged (1)
2580
// or unflagged (0).
2581
2582
inline void tetgenmesh::infect(triface& t) {
2583
((int *) (t.tet))[elemmarkerindex] |= 1;
2584
}
2585
2586
inline void tetgenmesh::uninfect(triface& t) {
2587
((int *) (t.tet))[elemmarkerindex] &= ~1;
2588
}
2589
2590
inline bool tetgenmesh::infected(triface& t) {
2591
return (((int *) (t.tet))[elemmarkerindex] & 1) != 0;
2592
}
2593
2594
// marktest(), marktested(), unmarktest() -- primitives to flag or unflag a
2595
// tetrahedron. Use the second lowerest bit of the element marker.
2596
2597
inline void tetgenmesh::marktest(triface& t) {
2598
((int *) (t.tet))[elemmarkerindex] |= 2;
2599
}
2600
2601
inline void tetgenmesh::unmarktest(triface& t) {
2602
((int *) (t.tet))[elemmarkerindex] &= ~2;
2603
}
2604
2605
inline bool tetgenmesh::marktested(triface& t) {
2606
return (((int *) (t.tet))[elemmarkerindex] & 2) != 0;
2607
}
2608
2609
// markface(), unmarkface(), facemarked() -- primitives to flag or unflag a
2610
// face of a tetrahedron. From the last 3rd to 6th bits are used for
2611
// face markers, e.g., the last third bit corresponds to loc = 0.
2612
2613
inline void tetgenmesh::markface(triface& t) {
2614
((int *) (t.tet))[elemmarkerindex] |= (4 << (t.ver & 3));
2615
}
2616
2617
inline void tetgenmesh::unmarkface(triface& t) {
2618
((int *) (t.tet))[elemmarkerindex] &= ~(4 << (t.ver & 3));
2619
}
2620
2621
inline bool tetgenmesh::facemarked(triface& t) {
2622
return (((int *) (t.tet))[elemmarkerindex] & (4 << (t.ver & 3))) != 0;
2623
}
2624
2625
// markedge(), unmarkedge(), edgemarked() -- primitives to flag or unflag an
2626
// edge of a tetrahedron. From the last 7th to 12th bits are used for
2627
// edge markers, e.g., the last 7th bit corresponds to the 0th edge, etc.
2628
// Remark: The last 7th bit is marked by 2^6 = 64.
2629
2630
inline void tetgenmesh::markedge(triface& t) {
2631
((int *) (t.tet))[elemmarkerindex] |= (int) (64 << ver2edge[(t).ver]);
2632
}
2633
2634
inline void tetgenmesh::unmarkedge(triface& t) {
2635
((int *) (t.tet))[elemmarkerindex] &= ~(int) (64 << ver2edge[(t).ver]);
2636
}
2637
2638
inline bool tetgenmesh::edgemarked(triface& t) {
2639
return (((int *) (t.tet))[elemmarkerindex] &
2640
(int) (64 << ver2edge[(t).ver])) != 0;
2641
}
2642
2643
// marktest2(), unmarktest2(), marktest2ed() -- primitives to flag and unflag
2644
// a tetrahedron. The 13th bit (2^12 = 4096) is used for this flag.
2645
2646
inline void tetgenmesh::marktest2(triface& t) {
2647
((int *) (t.tet))[elemmarkerindex] |= (int) (4096);
2648
}
2649
2650
inline void tetgenmesh::unmarktest2(triface& t) {
2651
((int *) (t.tet))[elemmarkerindex] &= ~(int) (4096);
2652
}
2653
2654
inline bool tetgenmesh::marktest2ed(triface& t) {
2655
return (((int *) (t.tet))[elemmarkerindex] & (int) (4096)) != 0;
2656
}
2657
2658
// elemcounter(), setelemcounter() -- primitives to read or set a (small)
2659
// integer counter in this tet. It is saved from the 16th bit. On 32 bit
2660
// system, the range of the counter is [0, 2^15 = 32768].
2661
2662
inline int tetgenmesh::elemcounter(triface& t) {
2663
return (((int *) (t.tet))[elemmarkerindex]) >> 16;
2664
}
2665
2666
inline void tetgenmesh::setelemcounter(triface& t, int value) {
2667
int c = ((int *) (t.tet))[elemmarkerindex];
2668
// Clear the old counter while keep the other flags.
2669
c &= 65535; // sum_{i=0^15} 2^i
2670
c |= (value << 16);
2671
((int *) (t.tet))[elemmarkerindex] = c;
2672
}
2673
2674
inline void tetgenmesh::increaseelemcounter(triface& t) {
2675
int c = elemcounter(t);
2676
setelemcounter(t, c + 1);
2677
}
2678
2679
inline void tetgenmesh::decreaseelemcounter(triface& t) {
2680
int c = elemcounter(t);
2681
setelemcounter(t, c - 1);
2682
}
2683
2684
// ishulltet() tests if t is a hull tetrahedron.
2685
2686
inline bool tetgenmesh::ishulltet(triface& t) {
2687
return (point) (t).tet[7] == dummypoint;
2688
}
2689
2690
// isdeadtet() tests if t is a tetrahedron is dead.
2691
2692
inline bool tetgenmesh::isdeadtet(triface& t) {
2693
return ((t.tet == NULL) || (t.tet[4] == NULL));
2694
}
2695
2696
///////////////////////////////////////////////////////////////////////////////
2697
// //
2698
// Primitives for subfaces and subsegments //
2699
// //
2700
///////////////////////////////////////////////////////////////////////////////
2701
2702
// Each subface contains three pointers to its neighboring subfaces, with
2703
// edge versions. To save memory, both information are kept in a single
2704
// pointer. To make this possible, all subfaces are aligned to eight-byte
2705
// boundaries, so that the last three bits of each pointer are zeros. An
2706
// edge version (in the range 0 to 5) is compressed into the last three
2707
// bits of each pointer by 'sencode()'. 'sdecode()' decodes a pointer,
2708
// extracting an edge version and a pointer to the beginning of a subface.
2709
2710
inline void tetgenmesh::sdecode(shellface sptr, face& s) {
2711
s.shver = (int) ((uintptr_t) (sptr) & (uintptr_t) 7);
2712
s.sh = (shellface *) ((uintptr_t) (sptr) ^ (uintptr_t) (s.shver));
2713
}
2714
2715
inline tetgenmesh::shellface tetgenmesh::sencode(face& s) {
2716
return (shellface) ((uintptr_t) s.sh | (uintptr_t) s.shver);
2717
}
2718
2719
inline tetgenmesh::shellface tetgenmesh::sencode2(shellface *sh, int shver) {
2720
return (shellface) ((uintptr_t) sh | (uintptr_t) shver);
2721
}
2722
2723
// sbond() bonds two subfaces (s1) and (s2) together. s1 and s2 must refer
2724
// to the same edge. No requirement is needed on their orientations.
2725
2726
inline void tetgenmesh::sbond(face& s1, face& s2)
2727
{
2728
s1.sh[s1.shver >> 1] = sencode(s2);
2729
s2.sh[s2.shver >> 1] = sencode(s1);
2730
}
2731
2732
// sbond1() bonds s1 <== s2, i.e., after bonding, s1 is pointing to s2,
2733
// but s2 is not pointing to s1. s1 and s2 must refer to the same edge.
2734
// No requirement is needed on their orientations.
2735
2736
inline void tetgenmesh::sbond1(face& s1, face& s2)
2737
{
2738
s1.sh[s1.shver >> 1] = sencode(s2);
2739
}
2740
2741
// Dissolve a subface bond (from one side). Note that the other subface
2742
// will still think it's connected to this subface.
2743
2744
inline void tetgenmesh::sdissolve(face& s)
2745
{
2746
s.sh[s.shver >> 1] = NULL;
2747
}
2748
2749
// spivot() finds the adjacent subface (s2) for a given subface (s1).
2750
// s1 and s2 share at the same edge.
2751
2752
inline void tetgenmesh::spivot(face& s1, face& s2)
2753
{
2754
shellface sptr = s1.sh[s1.shver >> 1];
2755
sdecode(sptr, s2);
2756
}
2757
2758
inline void tetgenmesh::spivotself(face& s)
2759
{
2760
shellface sptr = s.sh[s.shver >> 1];
2761
sdecode(sptr, s);
2762
}
2763
2764
// These primitives determine or set the origin, destination, or apex
2765
// of a subface with respect to the edge version.
2766
2767
inline tetgenmesh::point tetgenmesh::sorg(face& s)
2768
{
2769
return (point) s.sh[sorgpivot[s.shver]];
2770
}
2771
2772
inline tetgenmesh::point tetgenmesh::sdest(face& s)
2773
{
2774
return (point) s.sh[sdestpivot[s.shver]];
2775
}
2776
2777
inline tetgenmesh::point tetgenmesh::sapex(face& s)
2778
{
2779
return (point) s.sh[sapexpivot[s.shver]];
2780
}
2781
2782
inline void tetgenmesh::setsorg(face& s, point pointptr)
2783
{
2784
s.sh[sorgpivot[s.shver]] = (shellface) pointptr;
2785
}
2786
2787
inline void tetgenmesh::setsdest(face& s, point pointptr)
2788
{
2789
s.sh[sdestpivot[s.shver]] = (shellface) pointptr;
2790
}
2791
2792
inline void tetgenmesh::setsapex(face& s, point pointptr)
2793
{
2794
s.sh[sapexpivot[s.shver]] = (shellface) pointptr;
2795
}
2796
2797
#define setshvertices(s, pa, pb, pc)\
2798
setsorg(s, pa);\
2799
setsdest(s, pb);\
2800
setsapex(s, pc)
2801
2802
// sesym() reserves the direction of the lead edge.
2803
2804
inline void tetgenmesh::sesym(face& s1, face& s2)
2805
{
2806
s2.sh = s1.sh;
2807
s2.shver = (s1.shver ^ 1); // Inverse the last bit.
2808
}
2809
2810
inline void tetgenmesh::sesymself(face& s)
2811
{
2812
s.shver ^= 1;
2813
}
2814
2815
// senext() finds the next edge (counterclockwise) in the same orientation
2816
// of this face.
2817
2818
inline void tetgenmesh::senext(face& s1, face& s2)
2819
{
2820
s2.sh = s1.sh;
2821
s2.shver = snextpivot[s1.shver];
2822
}
2823
2824
inline void tetgenmesh::senextself(face& s)
2825
{
2826
s.shver = snextpivot[s.shver];
2827
}
2828
2829
inline void tetgenmesh::senext2(face& s1, face& s2)
2830
{
2831
s2.sh = s1.sh;
2832
s2.shver = snextpivot[snextpivot[s1.shver]];
2833
}
2834
2835
inline void tetgenmesh::senext2self(face& s)
2836
{
2837
s.shver = snextpivot[snextpivot[s.shver]];
2838
}
2839
2840
2841
// Check or set a subface's maximum area bound.
2842
2843
inline REAL tetgenmesh::areabound(face& s)
2844
{
2845
return ((REAL *) (s.sh))[areaboundindex];
2846
}
2847
2848
inline void tetgenmesh::setareabound(face& s, REAL value)
2849
{
2850
((REAL *) (s.sh))[areaboundindex] = value;
2851
}
2852
2853
// These two primitives read or set a shell marker. Shell markers are used
2854
// to hold user boundary information.
2855
2856
inline int tetgenmesh::shellmark(face& s)
2857
{
2858
return ((int *) (s.sh))[shmarkindex];
2859
}
2860
2861
inline void tetgenmesh::setshellmark(face& s, int value)
2862
{
2863
((int *) (s.sh))[shmarkindex] = value;
2864
}
2865
2866
2867
2868
// sinfect(), sinfected(), suninfect() -- primitives to flag or unflag a
2869
// subface. The last bit of ((int *) ((s).sh))[shmarkindex+1] is flagged.
2870
2871
inline void tetgenmesh::sinfect(face& s)
2872
{
2873
((int *) ((s).sh))[shmarkindex+1] =
2874
(((int *) ((s).sh))[shmarkindex+1] | (int) 1);
2875
}
2876
2877
inline void tetgenmesh::suninfect(face& s)
2878
{
2879
((int *) ((s).sh))[shmarkindex+1] =
2880
(((int *) ((s).sh))[shmarkindex+1] & ~(int) 1);
2881
}
2882
2883
// Test a subface for viral infection.
2884
2885
inline bool tetgenmesh::sinfected(face& s)
2886
{
2887
return (((int *) ((s).sh))[shmarkindex+1] & (int) 1) != 0;
2888
}
2889
2890
// smarktest(), smarktested(), sunmarktest() -- primitives to flag or unflag
2891
// a subface. The last 2nd bit of the integer is flagged.
2892
2893
inline void tetgenmesh::smarktest(face& s)
2894
{
2895
((int *) ((s).sh))[shmarkindex+1] =
2896
(((int *)((s).sh))[shmarkindex+1] | (int) 2);
2897
}
2898
2899
inline void tetgenmesh::sunmarktest(face& s)
2900
{
2901
((int *) ((s).sh))[shmarkindex+1] =
2902
(((int *)((s).sh))[shmarkindex+1] & ~(int)2);
2903
}
2904
2905
inline bool tetgenmesh::smarktested(face& s)
2906
{
2907
return ((((int *) ((s).sh))[shmarkindex+1] & (int) 2) != 0);
2908
}
2909
2910
// smarktest2(), smarktest2ed(), sunmarktest2() -- primitives to flag or
2911
// unflag a subface. The last 3rd bit of the integer is flagged.
2912
2913
inline void tetgenmesh::smarktest2(face& s)
2914
{
2915
((int *) ((s).sh))[shmarkindex+1] =
2916
(((int *)((s).sh))[shmarkindex+1] | (int) 4);
2917
}
2918
2919
inline void tetgenmesh::sunmarktest2(face& s)
2920
{
2921
((int *) ((s).sh))[shmarkindex+1] =
2922
(((int *)((s).sh))[shmarkindex+1] & ~(int)4);
2923
}
2924
2925
inline bool tetgenmesh::smarktest2ed(face& s)
2926
{
2927
return ((((int *) ((s).sh))[shmarkindex+1] & (int) 4) != 0);
2928
}
2929
2930
// The last 4th bit of ((int *) ((s).sh))[shmarkindex+1] is flagged.
2931
2932
inline void tetgenmesh::smarktest3(face& s)
2933
{
2934
((int *) ((s).sh))[shmarkindex+1] =
2935
(((int *)((s).sh))[shmarkindex+1] | (int) 8);
2936
}
2937
2938
inline void tetgenmesh::sunmarktest3(face& s)
2939
{
2940
((int *) ((s).sh))[shmarkindex+1] =
2941
(((int *)((s).sh))[shmarkindex+1] & ~(int)8);
2942
}
2943
2944
inline bool tetgenmesh::smarktest3ed(face& s)
2945
{
2946
return ((((int *) ((s).sh))[shmarkindex+1] & (int) 8) != 0);
2947
}
2948
2949
2950
// Each facet has a unique index (automatically indexed). Starting from '0'.
2951
// We save this index in the same field of the shell type.
2952
2953
inline void tetgenmesh::setfacetindex(face& s, int value)
2954
{
2955
((int *) (s.sh))[shmarkindex + 2] = value;
2956
}
2957
2958
inline int tetgenmesh::getfacetindex(face& s)
2959
{
2960
return ((int *) (s.sh))[shmarkindex + 2];
2961
}
2962
2963
///////////////////////////////////////////////////////////////////////////////
2964
// //
2965
// Primitives for interacting between tetrahedra and subfaces //
2966
// //
2967
///////////////////////////////////////////////////////////////////////////////
2968
2969
// tsbond() bond a tetrahedron (t) and a subface (s) together.
2970
// Note that t and s must be the same face and the same edge. Moreover,
2971
// t and s have the same orientation.
2972
// Since the edge number in t and in s can be any number in {0,1,2}. We bond
2973
// the edge in s which corresponds to t's 0th edge, and vice versa.
2974
2975
inline void tetgenmesh::tsbond(triface& t, face& s)
2976
{
2977
if ((t).tet[9] == NULL) {
2978
// Allocate space for this tet.
2979
(t).tet[9] = (tetrahedron) tet2subpool->alloc();
2980
// Initialize.
2981
for (int i = 0; i < 4; i++) {
2982
((shellface *) (t).tet[9])[i] = NULL;
2983
}
2984
}
2985
// Bond t <== s.
2986
((shellface *) (t).tet[9])[(t).ver & 3] =
2987
sencode2((s).sh, tsbondtbl[t.ver][s.shver]);
2988
// Bond s <== t.
2989
s.sh[9 + ((s).shver & 1)] =
2990
(shellface) encode2((t).tet, stbondtbl[t.ver][s.shver]);
2991
}
2992
2993
// tspivot() finds a subface (s) abutting on the given tetrahdera (t).
2994
// Return s.sh = NULL if there is no subface at t. Otherwise, return
2995
// the subface s, and s and t must be at the same edge with the same
2996
// orientation.
2997
2998
inline void tetgenmesh::tspivot(triface& t, face& s)
2999
{
3000
if ((t).tet[9] == NULL) {
3001
(s).sh = NULL;
3002
return;
3003
}
3004
// Get the attached subface s.
3005
sdecode(((shellface *) (t).tet[9])[(t).ver & 3], (s));
3006
(s).shver = tspivottbl[t.ver][s.shver];
3007
}
3008
3009
// Quickly check if the handle (t, v) is a subface.
3010
#define issubface(t) \
3011
((t).tet[9] && ((t).tet[9])[(t).ver & 3])
3012
3013
// stpivot() finds a tetrahedron (t) abutting a given subface (s).
3014
// Return the t (if it exists) with the same edge and the same
3015
// orientation of s.
3016
3017
inline void tetgenmesh::stpivot(face& s, triface& t)
3018
{
3019
decode((tetrahedron) s.sh[9 + (s.shver & 1)], t);
3020
if ((t).tet == NULL) {
3021
return;
3022
}
3023
(t).ver = stpivottbl[t.ver][s.shver];
3024
}
3025
3026
// Quickly check if this subface is attached to a tetrahedron.
3027
3028
#define isshtet(s) \
3029
((s).sh[9 + ((s).shver & 1)])
3030
3031
// tsdissolve() dissolve a bond (from the tetrahedron side).
3032
3033
inline void tetgenmesh::tsdissolve(triface& t)
3034
{
3035
if ((t).tet[9] != NULL) {
3036
((shellface *) (t).tet[9])[(t).ver & 3] = NULL;
3037
}
3038
}
3039
3040
// stdissolve() dissolve a bond (from the subface side).
3041
3042
inline void tetgenmesh::stdissolve(face& s)
3043
{
3044
(s).sh[9] = NULL;
3045
(s).sh[10] = NULL;
3046
}
3047
3048
///////////////////////////////////////////////////////////////////////////////
3049
// //
3050
// Primitives for interacting between subfaces and segments //
3051
// //
3052
///////////////////////////////////////////////////////////////////////////////
3053
3054
// ssbond() bond a subface to a subsegment.
3055
3056
inline void tetgenmesh::ssbond(face& s, face& edge)
3057
{
3058
s.sh[6 + (s.shver >> 1)] = sencode(edge);
3059
edge.sh[0] = sencode(s);
3060
}
3061
3062
inline void tetgenmesh::ssbond1(face& s, face& edge)
3063
{
3064
s.sh[6 + (s.shver >> 1)] = sencode(edge);
3065
//edge.sh[0] = sencode(s);
3066
}
3067
3068
// ssdisolve() dissolve a bond (from the subface side)
3069
3070
inline void tetgenmesh::ssdissolve(face& s)
3071
{
3072
s.sh[6 + (s.shver >> 1)] = NULL;
3073
}
3074
3075
// sspivot() finds a subsegment abutting a subface.
3076
3077
inline void tetgenmesh::sspivot(face& s, face& edge)
3078
{
3079
sdecode((shellface) s.sh[6 + (s.shver >> 1)], edge);
3080
}
3081
3082
// Quickly check if the edge is a subsegment.
3083
3084
#define isshsubseg(s) \
3085
((s).sh[6 + ((s).shver >> 1)])
3086
3087
///////////////////////////////////////////////////////////////////////////////
3088
// //
3089
// Primitives for interacting between tetrahedra and segments //
3090
// //
3091
///////////////////////////////////////////////////////////////////////////////
3092
3093
inline void tetgenmesh::tssbond1(triface& t, face& s)
3094
{
3095
if ((t).tet[8] == NULL) {
3096
// Allocate space for this tet.
3097
(t).tet[8] = (tetrahedron) tet2segpool->alloc();
3098
// Initialization.
3099
for (int i = 0; i < 6; i++) {
3100
((shellface *) (t).tet[8])[i] = NULL;
3101
}
3102
}
3103
((shellface *) (t).tet[8])[ver2edge[(t).ver]] = sencode((s));
3104
}
3105
3106
inline void tetgenmesh::sstbond1(face& s, triface& t)
3107
{
3108
((tetrahedron *) (s).sh)[9] = encode(t);
3109
}
3110
3111
inline void tetgenmesh::tssdissolve1(triface& t)
3112
{
3113
if ((t).tet[8] != NULL) {
3114
((shellface *) (t).tet[8])[ver2edge[(t).ver]] = NULL;
3115
}
3116
}
3117
3118
inline void tetgenmesh::sstdissolve1(face& s)
3119
{
3120
((tetrahedron *) (s).sh)[9] = NULL;
3121
}
3122
3123
inline void tetgenmesh::tsspivot1(triface& t, face& s)
3124
{
3125
if ((t).tet[8] != NULL) {
3126
sdecode(((shellface *) (t).tet[8])[ver2edge[(t).ver]], s);
3127
} else {
3128
(s).sh = NULL;
3129
}
3130
}
3131
3132
// Quickly check whether 't' is a segment or not.
3133
3134
#define issubseg(t) \
3135
((t).tet[8] && ((t).tet[8])[ver2edge[(t).ver]])
3136
3137
inline void tetgenmesh::sstpivot1(face& s, triface& t)
3138
{
3139
decode((tetrahedron) s.sh[9], t);
3140
}
3141
3142
///////////////////////////////////////////////////////////////////////////////
3143
// //
3144
// Primitives for points //
3145
// //
3146
///////////////////////////////////////////////////////////////////////////////
3147
3148
inline int tetgenmesh::pointmark(point pt) {
3149
return ((int *) (pt))[pointmarkindex];
3150
}
3151
3152
inline void tetgenmesh::setpointmark(point pt, int value) {
3153
((int *) (pt))[pointmarkindex] = value;
3154
}
3155
3156
3157
// These two primitives set and read the type of the point.
3158
3159
inline enum tetgenmesh::verttype tetgenmesh::pointtype(point pt) {
3160
return (enum verttype) (((int *) (pt))[pointmarkindex + 1] >> (int) 8);
3161
}
3162
3163
inline void tetgenmesh::setpointtype(point pt, enum verttype value) {
3164
((int *) (pt))[pointmarkindex + 1] =
3165
((int) value << 8) + (((int *) (pt))[pointmarkindex + 1] & (int) 255);
3166
}
3167
3168
// Read and set the geometry tag of the point (used by -s option).
3169
3170
inline int tetgenmesh::pointgeomtag(point pt) {
3171
return ((int *) (pt))[pointmarkindex + 2];
3172
}
3173
3174
inline void tetgenmesh::setpointgeomtag(point pt, int value) {
3175
((int *) (pt))[pointmarkindex + 2] = value;
3176
}
3177
3178
// Read and set the u,v coordinates of the point (used by -s option).
3179
3180
inline REAL tetgenmesh::pointgeomuv(point pt, int i) {
3181
return pt[pointparamindex + i];
3182
}
3183
3184
inline void tetgenmesh::setpointgeomuv(point pt, int i, REAL value) {
3185
pt[pointparamindex + i] = value;
3186
}
3187
3188
// pinfect(), puninfect(), pinfected() -- primitives to flag or unflag
3189
// a point. The last bit of the integer '[pointindex+1]' is flagged.
3190
3191
inline void tetgenmesh::pinfect(point pt) {
3192
((int *) (pt))[pointmarkindex + 1] |= (int) 1;
3193
}
3194
3195
inline void tetgenmesh::puninfect(point pt) {
3196
((int *) (pt))[pointmarkindex + 1] &= ~(int) 1;
3197
}
3198
3199
inline bool tetgenmesh::pinfected(point pt) {
3200
return (((int *) (pt))[pointmarkindex + 1] & (int) 1) != 0;
3201
}
3202
3203
// pmarktest(), punmarktest(), pmarktested() -- more primitives to
3204
// flag or unflag a point.
3205
3206
inline void tetgenmesh::pmarktest(point pt) {
3207
((int *) (pt))[pointmarkindex + 1] |= (int) 2;
3208
}
3209
3210
inline void tetgenmesh::punmarktest(point pt) {
3211
((int *) (pt))[pointmarkindex + 1] &= ~(int) 2;
3212
}
3213
3214
inline bool tetgenmesh::pmarktested(point pt) {
3215
return (((int *) (pt))[pointmarkindex + 1] & (int) 2) != 0;
3216
}
3217
3218
inline void tetgenmesh::pmarktest2(point pt) {
3219
((int *) (pt))[pointmarkindex + 1] |= (int) 4;
3220
}
3221
3222
inline void tetgenmesh::punmarktest2(point pt) {
3223
((int *) (pt))[pointmarkindex + 1] &= ~(int) 4;
3224
}
3225
3226
inline bool tetgenmesh::pmarktest2ed(point pt) {
3227
return (((int *) (pt))[pointmarkindex + 1] & (int) 4) != 0;
3228
}
3229
3230
inline void tetgenmesh::pmarktest3(point pt) {
3231
((int *) (pt))[pointmarkindex + 1] |= (int) 8;
3232
}
3233
3234
inline void tetgenmesh::punmarktest3(point pt) {
3235
((int *) (pt))[pointmarkindex + 1] &= ~(int) 8;
3236
}
3237
3238
inline bool tetgenmesh::pmarktest3ed(point pt) {
3239
return (((int *) (pt))[pointmarkindex + 1] & (int) 8) != 0;
3240
}
3241
3242
// These following primitives set and read a pointer to a tetrahedron
3243
// a subface/subsegment, a point, or a tet of background mesh.
3244
3245
inline tetgenmesh::tetrahedron tetgenmesh::point2tet(point pt) {
3246
return ((tetrahedron *) (pt))[point2simindex];
3247
}
3248
3249
inline void tetgenmesh::setpoint2tet(point pt, tetrahedron value) {
3250
((tetrahedron *) (pt))[point2simindex] = value;
3251
}
3252
3253
inline tetgenmesh::point tetgenmesh::point2ppt(point pt) {
3254
return (point) ((tetrahedron *) (pt))[point2simindex + 1];
3255
}
3256
3257
inline void tetgenmesh::setpoint2ppt(point pt, point value) {
3258
((tetrahedron *) (pt))[point2simindex + 1] = (tetrahedron) value;
3259
}
3260
3261
inline tetgenmesh::shellface tetgenmesh::point2sh(point pt) {
3262
return (shellface) ((tetrahedron *) (pt))[point2simindex + 2];
3263
}
3264
3265
inline void tetgenmesh::setpoint2sh(point pt, shellface value) {
3266
((tetrahedron *) (pt))[point2simindex + 2] = (tetrahedron) value;
3267
}
3268
3269
3270
inline tetgenmesh::tetrahedron tetgenmesh::point2bgmtet(point pt) {
3271
return ((tetrahedron *) (pt))[point2simindex + 3];
3272
}
3273
3274
inline void tetgenmesh::setpoint2bgmtet(point pt, tetrahedron value) {
3275
((tetrahedron *) (pt))[point2simindex + 3] = value;
3276
}
3277
3278
3279
// The primitives for saving and getting the insertion radius.
3280
inline void tetgenmesh::setpointinsradius(point pt, REAL value)
3281
{
3282
pt[pointinsradiusindex] = value;
3283
}
3284
3285
inline REAL tetgenmesh::getpointinsradius(point pt)
3286
{
3287
return pt[pointinsradiusindex];
3288
}
3289
3290
inline bool tetgenmesh::issteinerpoint(point pt) {
3291
return (pointtype(pt) == FREESEGVERTEX) || (pointtype(pt) == FREEFACETVERTEX)
3292
|| (pointtype(pt) == FREEVOLVERTEX);
3293
}
3294
3295
// point2tetorg() Get the tetrahedron whose origin is the point.
3296
3297
inline void tetgenmesh::point2tetorg(point pa, triface& searchtet)
3298
{
3299
decode(point2tet(pa), searchtet);
3300
if ((point) searchtet.tet[4] == pa) {
3301
searchtet.ver = 11;
3302
} else if ((point) searchtet.tet[5] == pa) {
3303
searchtet.ver = 3;
3304
} else if ((point) searchtet.tet[6] == pa) {
3305
searchtet.ver = 7;
3306
} else {
3307
searchtet.ver = 0;
3308
}
3309
}
3310
3311
// point2shorg() Get the subface/segment whose origin is the point.
3312
3313
inline void tetgenmesh::point2shorg(point pa, face& searchsh)
3314
{
3315
sdecode(point2sh(pa), searchsh);
3316
if ((point) searchsh.sh[3] == pa) {
3317
searchsh.shver = 0;
3318
} else if ((point) searchsh.sh[4] == pa) {
3319
searchsh.shver = (searchsh.sh[5] != NULL ? 2 : 1);
3320
} else {
3321
searchsh.shver = 4;
3322
}
3323
}
3324
3325
// farsorg() Return the origin of the subsegment.
3326
// farsdest() Return the destination of the subsegment.
3327
3328
inline tetgenmesh::point tetgenmesh::farsorg(face& s)
3329
{
3330
face travesh, neighsh;
3331
3332
travesh = s;
3333
while (1) {
3334
senext2(travesh, neighsh);
3335
spivotself(neighsh);
3336
if (neighsh.sh == NULL) break;
3337
if (sorg(neighsh) != sorg(travesh)) sesymself(neighsh);
3338
senext2(neighsh, travesh);
3339
}
3340
return sorg(travesh);
3341
}
3342
3343
inline tetgenmesh::point tetgenmesh::farsdest(face& s)
3344
{
3345
face travesh, neighsh;
3346
3347
travesh = s;
3348
while (1) {
3349
senext(travesh, neighsh);
3350
spivotself(neighsh);
3351
if (neighsh.sh == NULL) break;
3352
if (sdest(neighsh) != sdest(travesh)) sesymself(neighsh);
3353
senext(neighsh, travesh);
3354
}
3355
return sdest(travesh);
3356
}
3357
3358
///////////////////////////////////////////////////////////////////////////////
3359
// //
3360
// Linear algebra operators. //
3361
// //
3362
///////////////////////////////////////////////////////////////////////////////
3363
3364
// dot() returns the dot product: v1 dot v2.
3365
inline REAL tetgenmesh::dot(REAL* v1, REAL* v2)
3366
{
3367
return v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2];
3368
}
3369
3370
// cross() computes the cross product: n = v1 cross v2.
3371
inline void tetgenmesh::cross(REAL* v1, REAL* v2, REAL* n)
3372
{
3373
n[0] = v1[1] * v2[2] - v2[1] * v1[2];
3374
n[1] = -(v1[0] * v2[2] - v2[0] * v1[2]);
3375
n[2] = v1[0] * v2[1] - v2[0] * v1[1];
3376
}
3377
3378
// distance() computes the Euclidean distance between two points.
3379
inline REAL tetgenmesh::distance(REAL* p1, REAL* p2)
3380
{
3381
return sqrt((p2[0] - p1[0]) * (p2[0] - p1[0]) +
3382
(p2[1] - p1[1]) * (p2[1] - p1[1]) +
3383
(p2[2] - p1[2]) * (p2[2] - p1[2]));
3384
}
3385
3386
inline REAL tetgenmesh::norm2(REAL x, REAL y, REAL z)
3387
{
3388
return (x) * (x) + (y) * (y) + (z) * (z);
3389
}
3390
3391
3392
#endif // #ifndef tetgenH
3393
3394
3395