forked from libigl/tetgen
-
Notifications
You must be signed in to change notification settings - Fork 0
/
tetgen.h
3613 lines (3061 loc) · 143 KB
/
tetgen.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
//============================================================================//
// //
// TetGen //
// //
// A Quality Tetrahedral Mesh Generator and A 3D Delaunay Triangulator //
// //
// Version 1.6.0 //
// August 31, 2020 //
// //
// Copyright (C) 2002--2020 //
// //
// Hang Si //
// Research Group: Numerical Mathematics and Scientific Computing //
// Weierstrass Institute for Applied Analysis and Stochastics (WIAS) //
// Mohrenstr. 39, 10117 Berlin, Germany //
// si@wias-berlin.de //
// //
// TetGen is a tetrahedral mesh generator. It creates 3d triangulations of //
// polyhedral domains. It generates meshes with well-shaped elements whose //
// sizes are adapted to the geometric features or user-provided sizing //
// functions. It has applications in various applications in scientific //
// computing, such as computer graphics (CG), computer-aided design (CAD), //
// geometry processing (parametrizations and computer animation), and //
// physical simulations (finite element analysis). //
// //
// TetGen computes (weighted) Delaunay triangulations for three-dimensional //
// (weighted) point sets, and constrained Delaunay triangulations for //
// three-dimensional polyhedral domains. In the latter case, input edges //
// and triangles can be completely preserved in the output meshes. TetGen //
// can refine or coarsen an existing mesh to result in good quality and //
// size-adapted mesh according to the geometric features and user-defined //
// mesh sizing functions. //
// //
// TetGen implements theoretically proven algorithms for computing the //
// Delaunay and constrained Delaunay tetrahedralizations. TetGen achieves //
// robustness and efficiency by using advanced techniques in computational //
// geometry. A technical paper describes the algorithms and methods //
// implemented in TetGen is available in ACM-TOMS, Hang Si ``TetGen, a //
// Delaunay-Based Quality Tetrahedral Mesh Generator", ACM Transactions on //
// Mathematical Software, February 2015, https://doi.org/10.1145/2629697. //
// //
// TetGen is freely available through the website: http://www.tetgen.org. //
// It may be copied, modified, and redistributed for non-commercial use. //
// Please consult the file LICENSE for the detailed copyright notices. //
// //
//============================================================================//
#ifndef tetgenH
#define tetgenH
// To compile TetGen as a library instead of an executable program, define
// the TETLIBRARY symbol.
// #define TETLIBRARY
// TetGen default uses the double-precision (64 bit) for a real number.
// Alternatively, one can use the single-precision (32 bit) 'float' if the
// memory is limited.
#define REAL double // #define REAL float
// The maximum number of characters in a file name (including the null).
#define FILENAMESIZE 1024
// The maximum number of chars in a line read from a file (including the null).
#define INPUTLINESIZE 2048
// C standard libraries to perform Input/output operations, general utililities,
// manipulate strings and arrays, compute common mathematical operations,
// get date and time information.
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include <time.h>
// The types 'intptr_t' and 'uintptr_t' are signed and unsigned integer types,
// respectively. They are guaranteed to be the same width as a pointer.
// They are defined in <stdint.h> by the C99 Standard.
#include <stdint.h>
//============================================================================//
// //
// tetgenio //
// //
// A structure for transferring input/output data between the user and //
// TetGen's internal data structure (class tetgenmesh). //
// //
// This data structure contains a collection of arrays, i.e., points, facets, //
// tetrahedra. It contains functions to read input data from files (.node, //
// .poly, .face, .edge, .ele) as well as write output data into files. //
// //
// Once an object of tetgenio is declared, no array is created. One has to //
// allocate enough memory for them. On the deletion of this object, the //
// memory occupied by these arrays needs to be freed. The routine //
// deinitialize() will be automatically called. It frees the memory for //
// an array if it is not a NULL. Note that it assumes that the memory is //
// allocated by the C++ "new" operator. Otherwise, the user is responsible //
// to free them and all pointers must be NULL. //
// //
//============================================================================//
class tetgenio {
public:
// A "polygon" describes a simple polygon (no holes). It is not necessarily
// convex. Each polygon contains a number of corners (points) and the same
// number of sides (edges). The points of the polygon must be given in
// either counterclockwise or clockwise order and they form a ring, so
// every two consecutive points forms an edge of the polygon.
typedef struct {
int *vertexlist;
int numberofvertices;
} polygon;
// A "facet" describes a polygonal region possibly with holes, edges, and
// points floating in it. Each facet consists of a list of polygons and
// a list of hole points (which lie strictly inside holes).
typedef struct {
polygon *polygonlist;
int numberofpolygons;
REAL *holelist;
int numberofholes;
} facet;
// A "voroedge" is an edge of the Voronoi diagram. It corresponds to a
// Delaunay face. Each voroedge is either a line segment connecting
// two Voronoi vertices or a ray starting from a Voronoi vertex to an
// "infinite vertex". 'v1' and 'v2' are two indices pointing to the
// list of Voronoi vertices. 'v1' must be non-negative, while 'v2' may
// be -1 if it is a ray, in this case, the unit normal of this ray is
// given in 'vnormal'.
typedef struct {
int v1, v2;
REAL vnormal[3];
} voroedge;
// A "vorofacet" is an facet of the Voronoi diagram. It corresponds to a
// Delaunay edge. Each Voronoi facet is a convex polygon formed by a
// list of Voronoi edges, it may not be closed. 'c1' and 'c2' are two
// indices pointing into the list of Voronoi cells, i.e., the two cells
// share this facet. 'elist' is an array of indices pointing into the
// list of Voronoi edges, 'elist[0]' saves the number of Voronoi edges
// (including rays) of this facet.
typedef struct {
int c1, c2;
int *elist;
} vorofacet;
// Additional parameters associated with an input (or mesh) vertex.
// These informations are provided by CAD libraries.
typedef struct {
REAL uv[2];
int tag;
int type; // 0, 1, or 2.
} pointparam;
// Callback functions for meshing PSCs.
typedef REAL (* GetVertexParamOnEdge)(void*, int, int);
typedef void (* GetSteinerOnEdge)(void*, int, REAL, REAL*);
typedef void (* GetVertexParamOnFace)(void*, int, int, REAL*);
typedef void (* GetEdgeSteinerParamOnFace)(void*, int, REAL, int, REAL*);
typedef void (* GetSteinerOnFace)(void*, int, REAL*, REAL*);
// A callback function for mesh refinement.
typedef bool (* TetSizeFunc)(REAL*, REAL*, REAL*, REAL*, REAL*, REAL);
// Items are numbered starting from 'firstnumber' (0 or 1), default is 0.
int firstnumber;
// Dimension of the mesh (2 or 3), default is 3.
int mesh_dim;
// Does the lines in .node file contain index or not, default is 1.
int useindex;
// 'pointlist': An array of point coordinates. The first point's x
// coordinate is at index [0] and its y coordinate at index [1], its
// z coordinate is at index [2], followed by the coordinates of the
// remaining points. Each point occupies three REALs.
// 'pointattributelist': An array of point attributes. Each point's
// attributes occupy 'numberofpointattributes' REALs.
// 'pointmtrlist': An array of metric tensors at points. Each point's
// tensor occupies 'numberofpointmtr' REALs.
// 'pointmarkerlist': An array of point markers; one integer per point.
// 'point2tetlist': An array of tetrahedra indices; one integer per point.
REAL *pointlist;
REAL *pointattributelist;
REAL *pointmtrlist;
int *pointmarkerlist;
int *point2tetlist;
pointparam *pointparamlist;
int numberofpoints;
int numberofpointattributes;
int numberofpointmtrs;
// 'tetrahedronlist': An array of tetrahedron corners. The first
// tetrahedron's first corner is at index [0], followed by its other
// corners, followed by six nodes on the edges of the tetrahedron if the
// second order option (-o2) is applied. Each tetrahedron occupies
// 'numberofcorners' ints. The second order nodes are ouput only.
// 'tetrahedronattributelist': An array of tetrahedron attributes. Each
// tetrahedron's attributes occupy 'numberoftetrahedronattributes' REALs.
// 'tetrahedronvolumelist': An array of constraints, i.e. tetrahedron's
// volume; one REAL per element. Input only.
// 'neighborlist': An array of tetrahedron neighbors; 4 ints per element.
// 'tet2facelist': An array of tetrahedron face indices; 4 ints per element.
// 'tet2edgelist': An array of tetrahedron edge indices; 6 ints per element.
int *tetrahedronlist;
REAL *tetrahedronattributelist;
REAL *tetrahedronvolumelist;
int *neighborlist;
int *tet2facelist;
int *tet2edgelist;
int numberoftetrahedra;
int numberofcorners;
int numberoftetrahedronattributes;
// 'facetlist': An array of facets. Each entry is a structure of facet.
// 'facetmarkerlist': An array of facet markers; one int per facet.
facet *facetlist;
int *facetmarkerlist;
int numberoffacets;
// 'holelist': An array of holes (in volume). Each hole is given by a
// seed (point) which lies strictly inside it. The first seed's x, y and z
// coordinates are at indices [0], [1] and [2], followed by the
// remaining seeds. Three REALs per hole.
REAL *holelist;
int numberofholes;
// 'regionlist': An array of regions (subdomains). Each region is given by
// a seed (point) which lies strictly inside it. The first seed's x, y and
// z coordinates are at indices [0], [1] and [2], followed by the regional
// attribute at index [3], followed by the maximum volume at index [4].
// Five REALs per region.
// Note that each regional attribute is used only if you select the 'A'
// switch, and each volume constraint is used only if you select the
// 'a' switch (with no number following).
REAL *regionlist;
int numberofregions;
// 'refine_elem_list': An array of tetrahedra to be refined. The first
// tetrahedron's first corner is at index [0], followed by its other
// corners. Four integers per element.
// 'refine_elem_vol_list': An array of constraints, i.e. tetrahedron's
// volume; one REAL per element.
int *refine_elem_list;
REAL *refine_elem_vol_list;
int numberofrefineelems;
// 'facetconstraintlist': An array of facet constraints. Each constraint
// specifies a maximum area bound on the subfaces of that facet. The
// first facet constraint is given by a facet marker at index [0] and its
// maximum area bound at index [1], followed by the remaining facet con-
// straints. Two REALs per facet constraint. Note: the facet marker is
// actually an integer.
REAL *facetconstraintlist;
int numberoffacetconstraints;
// 'segmentconstraintlist': An array of segment constraints. Each constraint
// specifies a maximum length bound on the subsegments of that segment.
// The first constraint is given by the two endpoints of the segment at
// index [0] and [1], and the maximum length bound at index [2], followed
// by the remaining segment constraints. Three REALs per constraint.
// Note the segment endpoints are actually integers.
REAL *segmentconstraintlist;
int numberofsegmentconstraints;
// 'trifacelist': An array of face (triangle) corners. The first face's
// three corners are at indices [0], [1] and [2], followed by the remaining
// faces. Three ints per face.
// 'trifacemarkerlist': An array of face markers; one int per face.
// 'o2facelist': An array of second order nodes (on the edges) of the face.
// It is output only if the second order option (-o2) is applied. The
// first face's three second order nodes are at [0], [1], and [2],
// followed by the remaining faces. Three ints per face.
// 'face2tetlist': An array of tetrahedra indices; 2 ints per face.
// 'face2edgelist': An array of edge indices; 3 ints per face.
int *trifacelist;
int *trifacemarkerlist;
int *o2facelist;
int *face2tetlist;
int *face2edgelist;
int numberoftrifaces;
// 'edgelist': An array of edge endpoints. The first edge's endpoints
// are at indices [0] and [1], followed by the remaining edges.
// Two ints per edge.
// 'edgemarkerlist': An array of edge markers; one int per edge.
// 'o2edgelist': An array of midpoints of edges. It is output only if the
// second order option (-o2) is applied. One int per edge.
// 'edge2tetlist': An array of tetrahedra indices. One int per edge.
int *edgelist;
int *edgemarkerlist;
int *o2edgelist;
int *edge2tetlist;
int numberofedges;
// 'vpointlist': An array of Voronoi vertex coordinates (like pointlist).
// 'vedgelist': An array of Voronoi edges. Each entry is a 'voroedge'.
// 'vfacetlist': An array of Voronoi facets. Each entry is a 'vorofacet'.
// 'vcelllist': An array of Voronoi cells. Each entry is an array of
// indices pointing into 'vfacetlist'. The 0th entry is used to store
// the length of this array.
REAL *vpointlist;
voroedge *vedgelist;
vorofacet *vfacetlist;
int **vcelllist;
int numberofvpoints;
int numberofvedges;
int numberofvfacets;
int numberofvcells;
// Variable (and callback functions) for meshing PSCs.
void *geomhandle;
GetVertexParamOnEdge getvertexparamonedge;
GetSteinerOnEdge getsteineronedge;
GetVertexParamOnFace getvertexparamonface;
GetEdgeSteinerParamOnFace getedgesteinerparamonface;
GetSteinerOnFace getsteineronface;
// A callback function.
TetSizeFunc tetunsuitable;
// Input & output routines.
bool load_node_call(FILE* infile, int markers, int uvflag, char*);
bool load_node(char*);
bool load_edge(char*);
bool load_face(char*);
bool load_tet(char*);
bool load_vol(char*);
bool load_var(char*);
bool load_mtr(char*);
bool load_elem(char*);
bool load_poly(char*);
bool load_off(char*);
bool load_ply(char*);
bool load_stl(char*);
bool load_vtk(char*);
bool load_medit(char*, int);
bool load_neumesh(char*, int);
bool load_plc(char*, int);
bool load_tetmesh(char*, int);
void save_nodes(const char*);
void save_elements(const char*);
void save_faces(const char*);
void save_edges(char*);
void save_neighbors(char*);
void save_poly(const char*);
void save_faces2smesh(char*);
// Read line and parse string functions.
char *readline(char* string, FILE* infile, int *linenumber);
char *findnextfield(char* string);
char *readnumberline(char* string, FILE* infile, char* infilename);
char *findnextnumber(char* string);
static void init(polygon* p) {
p->vertexlist = (int *) NULL;
p->numberofvertices = 0;
}
static void init(facet* f) {
f->polygonlist = (polygon *) NULL;
f->numberofpolygons = 0;
f->holelist = (REAL *) NULL;
f->numberofholes = 0;
}
// Initialize routine.
void initialize()
{
firstnumber = 0;
mesh_dim = 3;
useindex = 1;
pointlist = (REAL *) NULL;
pointattributelist = (REAL *) NULL;
pointmtrlist = (REAL *) NULL;
pointmarkerlist = (int *) NULL;
point2tetlist = (int *) NULL;
pointparamlist = (pointparam *) NULL;
numberofpoints = 0;
numberofpointattributes = 0;
numberofpointmtrs = 0;
tetrahedronlist = (int *) NULL;
tetrahedronattributelist = (REAL *) NULL;
tetrahedronvolumelist = (REAL *) NULL;
neighborlist = (int *) NULL;
tet2facelist = (int *) NULL;
tet2edgelist = (int *) NULL;
numberoftetrahedra = 0;
numberofcorners = 4;
numberoftetrahedronattributes = 0;
trifacelist = (int *) NULL;
trifacemarkerlist = (int *) NULL;
o2facelist = (int *) NULL;
face2tetlist = (int *) NULL;
face2edgelist = (int *) NULL;
numberoftrifaces = 0;
edgelist = (int *) NULL;
edgemarkerlist = (int *) NULL;
o2edgelist = (int *) NULL;
edge2tetlist = (int *) NULL;
numberofedges = 0;
facetlist = (facet *) NULL;
facetmarkerlist = (int *) NULL;
numberoffacets = 0;
holelist = (REAL *) NULL;
numberofholes = 0;
regionlist = (REAL *) NULL;
numberofregions = 0;
refine_elem_list = (int *) NULL;
refine_elem_vol_list = (REAL *) NULL;
numberofrefineelems = 0;
facetconstraintlist = (REAL *) NULL;
numberoffacetconstraints = 0;
segmentconstraintlist = (REAL *) NULL;
numberofsegmentconstraints = 0;
vpointlist = (REAL *) NULL;
vedgelist = (voroedge *) NULL;
vfacetlist = (vorofacet *) NULL;
vcelllist = (int **) NULL;
numberofvpoints = 0;
numberofvedges = 0;
numberofvfacets = 0;
numberofvcells = 0;
tetunsuitable = NULL;
geomhandle = NULL;
getvertexparamonedge = NULL;
getsteineronedge = NULL;
getvertexparamonface = NULL;
getedgesteinerparamonface = NULL;
getsteineronface = NULL;
}
// Free the memory allocated in 'tetgenio'. Note that it assumes that the
// memory was allocated by the "new" operator (C++).
void clean_memory()
{
int i, j;
if (pointlist != (REAL *) NULL) {
delete [] pointlist;
}
if (pointattributelist != (REAL *) NULL) {
delete [] pointattributelist;
}
if (pointmtrlist != (REAL *) NULL) {
delete [] pointmtrlist;
}
if (pointmarkerlist != (int *) NULL) {
delete [] pointmarkerlist;
}
if (point2tetlist != (int *) NULL) {
delete [] point2tetlist;
}
if (pointparamlist != (pointparam *) NULL) {
delete [] pointparamlist;
}
if (tetrahedronlist != (int *) NULL) {
delete [] tetrahedronlist;
}
if (tetrahedronattributelist != (REAL *) NULL) {
delete [] tetrahedronattributelist;
}
if (tetrahedronvolumelist != (REAL *) NULL) {
delete [] tetrahedronvolumelist;
}
if (neighborlist != (int *) NULL) {
delete [] neighborlist;
}
if (tet2facelist != (int *) NULL) {
delete [] tet2facelist;
}
if (tet2edgelist != (int *) NULL) {
delete [] tet2edgelist;
}
if (trifacelist != (int *) NULL) {
delete [] trifacelist;
}
if (trifacemarkerlist != (int *) NULL) {
delete [] trifacemarkerlist;
}
if (o2facelist != (int *) NULL) {
delete [] o2facelist;
}
if (face2tetlist != (int *) NULL) {
delete [] face2tetlist;
}
if (face2edgelist != (int *) NULL) {
delete [] face2edgelist;
}
if (edgelist != (int *) NULL) {
delete [] edgelist;
}
if (edgemarkerlist != (int *) NULL) {
delete [] edgemarkerlist;
}
if (o2edgelist != (int *) NULL) {
delete [] o2edgelist;
}
if (edge2tetlist != (int *) NULL) {
delete [] edge2tetlist;
}
if (facetlist != (facet *) NULL) {
facet *f;
polygon *p;
for (i = 0; i < numberoffacets; i++) {
f = &facetlist[i];
for (j = 0; j < f->numberofpolygons; j++) {
p = &f->polygonlist[j];
delete [] p->vertexlist;
}
delete [] f->polygonlist;
if (f->holelist != (REAL *) NULL) {
delete [] f->holelist;
}
}
delete [] facetlist;
}
if (facetmarkerlist != (int *) NULL) {
delete [] facetmarkerlist;
}
if (holelist != (REAL *) NULL) {
delete [] holelist;
}
if (regionlist != (REAL *) NULL) {
delete [] regionlist;
}
if (refine_elem_list != (int *) NULL) {
delete [] refine_elem_list;
if (refine_elem_vol_list != (REAL *) NULL) {
delete [] refine_elem_vol_list;
}
}
if (facetconstraintlist != (REAL *) NULL) {
delete [] facetconstraintlist;
}
if (segmentconstraintlist != (REAL *) NULL) {
delete [] segmentconstraintlist;
}
if (vpointlist != (REAL *) NULL) {
delete [] vpointlist;
}
if (vedgelist != (voroedge *) NULL) {
delete [] vedgelist;
}
if (vfacetlist != (vorofacet *) NULL) {
for (i = 0; i < numberofvfacets; i++) {
delete [] vfacetlist[i].elist;
}
delete [] vfacetlist;
}
if (vcelllist != (int **) NULL) {
for (i = 0; i < numberofvcells; i++) {
delete [] vcelllist[i];
}
delete [] vcelllist;
}
}
// Constructor & destructor.
tetgenio() {initialize();}
~tetgenio() {clean_memory();}
}; // class tetgenio
//============================================================================//
// //
// tetgenbehavior //
// //
// A structure for maintaining the switches and parameters used by TetGen's //
// internal data structure and algorithms. //
// //
// All switches and parameters are initialized with default values. They are //
// set by the command line arguments (argc, argv). //
// //
// NOTE: Some switches are incompatible with others. While some may depend //
// on other switches. The routine parse_commandline() sets the switches from //
// the command line (a list of strings) and checks the consistency of the //
// applied switches. //
// //
//============================================================================//
class tetgenbehavior {
public:
// Switches of TetGen.
int plc; // '-p', 0.
int psc; // '-s', 0.
int refine; // '-r', 0.
int quality; // '-q', 0.
int nobisect; // '-Y', 0.
int cdt; // '-D', 0.
int cdtrefine; // '-D#', 7.
int coarsen; // '-R', 0.
int weighted; // '-w', 0.
int brio_hilbert; // '-b', 1.
int flipinsert; // '-L', 0.
int metric; // '-m', 0.
int varvolume; // '-a', 0.
int fixedvolume; // '-a', 0.
int regionattrib; // '-A', 0.
int insertaddpoints; // '-i', 0.
int diagnose; // '-d', 0.
int convex; // '-c', 0.
int nomergefacet; // '-M', 0.
int nomergevertex; // '-M', 0.
int noexact; // '-X', 0.
int nostaticfilter; // '-X', 0.
int zeroindex; // '-z', 0.
int facesout; // '-f', 0.
int edgesout; // '-e', 0.
int neighout; // '-n', 0.
int voroout; // '-v', 0.
int meditview; // '-g', 0.
int vtkview; // '-k', 0.
int vtksurfview; // '-k', 0.
int nobound; // '-B', 0.
int nonodewritten; // '-N', 0.
int noelewritten; // '-E', 0.
int nofacewritten; // '-F', 0.
int noiterationnum; // '-I', 0.
int nojettison; // '-J', 0.
int docheck; // '-C', 0.
int quiet; // '-Q', 0.
int nowarning; // '-W', 0.
int verbose; // '-V', 0.
// Parameters of TetGen.
int vertexperblock; // '-x', 4092.
int tetrahedraperblock; // '-x', 8188.
int shellfaceperblock; // '-x', 2044.
int supsteiner_level; // '-Y/', 2.
int addsteiner_algo; // '-Y//', 1.
int coarsen_param; // '-R', 0.
int weighted_param; // '-w', 0.
int fliplinklevel; // -1.
int flipstarsize; // -1.
int fliplinklevelinc; // 1.
int opt_max_flip_level; // '-O', 3.
int opt_scheme; // '-O/#', 7.
int opt_iterations; // -O//#, 3.
int smooth_cirterion; // -s, 1.
int smooth_maxiter; // -s, 7.
int delmaxfliplevel; // 1.
int order; // '-o', 1.
int reversetetori; // '-o/', 0.
int steinerleft; // '-S', 0.
int unflip_queue_limit; // '-U#', 1000.
int no_sort; // 0.
int hilbert_order; // '-b///', 52.
int hilbert_limit; // '-b//' 8.
int brio_threshold; // '-b' 64.
REAL brio_ratio; // '-b/' 0.125.
REAL epsilon; // '-T', 1.0e-8.
REAL facet_separate_ang_tol; // '-p', 179.9.
REAL collinear_ang_tol; // '-p/', 179.9.
REAL facet_small_ang_tol; // '-p//', 15.0.
REAL maxvolume; // '-a', -1.0.
REAL maxvolume_length; // '-a', -1.0.
REAL minratio; // '-q', 0.0.
REAL opt_max_asp_ratio; // 1000.0.
REAL opt_max_edge_ratio; // 100.0.
REAL mindihedral; // '-q', 5.0.
REAL optmaxdihedral; // -o/# 177.0.
REAL metric_scale; // -m#, 1.0.
REAL smooth_alpha; // '-s', 0.3.
REAL coarsen_percent; // -R1/#, 1.0.
REAL elem_growth_ratio; // Growth ratio of # elements, -r#, 0.0.
REAL refine_progress_ratio; // -r/#, 0.333.
// Strings of command line arguments and input/output file names.
char commandline[1024];
char infilename[1024];
char outfilename[1024];
char addinfilename[1024];
char bgmeshfilename[1024];
// Read an additional tetrahedral mesh and treat it as holes [2018-07-30].
int hole_mesh; // '-H', 0.
char hole_mesh_filename[1024];
// The input object of TetGen. They are recognized by either the input
// file extensions or by the specified options.
// Currently the following objects are supported:
// - NODES, a list of nodes (.node);
// - POLY, a piecewise linear complex (.poly or .smesh);
// - OFF, a polyhedron (.off, Geomview's file format);
// - PLY, a polyhedron (.ply, file format from gatech, only ASCII);
// - STL, a surface mesh (.stl, stereolithography format);
// - MEDIT, a surface mesh (.mesh, Medit's file format);
// - MESH, a tetrahedral mesh (.ele).
// If no extension is available, the imposed command line switch
// (-p or -r) implies the object.
enum objecttype {NODES, POLY, OFF, PLY, STL, MEDIT, VTK, MESH, NEU_MESH} object;
void syntax();
void usage();
// Command line parse routine.
bool parse_commandline(int argc, char **argv);
bool parse_commandline(char *switches) {
return parse_commandline(0, &switches);
}
// Initialize all variables.
tetgenbehavior()
{
plc = 0;
psc = 0;
refine = 0;
quality = 0;
nobisect = 0;
cdt = 0; // set by -D (without a number following it)
cdtrefine = 7; // default, set by -D#
coarsen = 0;
metric = 0;
weighted = 0;
brio_hilbert = 1;
flipinsert = 0;
varvolume = 0;
fixedvolume = 0;
noexact = 0;
nostaticfilter = 0;
insertaddpoints = 0;
regionattrib = 0;
diagnose = 0;
convex = 0;
zeroindex = 0;
facesout = 0;
edgesout = 0;
neighout = 0;
voroout = 0;
meditview = 0;
vtkview = 0;
vtksurfview = 0;
nobound = 0;
nonodewritten = 0;
noelewritten = 0;
nofacewritten = 0;
noiterationnum = 0;
nomergefacet = 0;
nomergevertex = 0;
nojettison = 0;
docheck = 0;
quiet = 0;
nowarning = 0;
verbose = 0;
vertexperblock = 4092;
tetrahedraperblock = 8188;
shellfaceperblock = 4092;
supsteiner_level = 2;
addsteiner_algo = 1;
coarsen_param = 0;
weighted_param = 0;
fliplinklevel = -1;
flipstarsize = -1;
fliplinklevelinc = 1;
opt_scheme = 7;
opt_max_flip_level = 3;
opt_iterations = 3;
delmaxfliplevel = 1;
order = 1;
reversetetori = 0;
steinerleft = -1;
unflip_queue_limit = 1000;
no_sort = 0;
hilbert_order = 52; //-1;
hilbert_limit = 8;
brio_threshold = 64;
brio_ratio = 0.125;
facet_separate_ang_tol = 179.9;
collinear_ang_tol = 179.9;
facet_small_ang_tol = 15.0;
maxvolume = -1.0;
maxvolume_length = -1.0;
minratio = 2.0;
opt_max_asp_ratio = 1000.;
opt_max_edge_ratio = 100.;
mindihedral = 3.5;
optmaxdihedral = 177.00;
epsilon = 1.0e-8;
coarsen_percent = 1.0;
metric_scale = 1.0; // -m#
elem_growth_ratio = 0.0; // -r#
refine_progress_ratio = 0.333; // -r/#
object = NODES;
smooth_cirterion = 3; // -s# default smooth surface and volume vertices.
smooth_maxiter = 7; // set by -s#/7
smooth_alpha = 0.3; // relax parameter, set by -s#/#/0.3
commandline[0] = '\0';
infilename[0] = '\0';
outfilename[0] = '\0';
addinfilename[0] = '\0';
bgmeshfilename[0] = '\0';
hole_mesh = 0;
hole_mesh_filename[0] = '\0';
}
}; // class tetgenbehavior
//============================================================================//
// //
// Robust Geometric predicates //
// //
// The following routines are the robust geometric predicates for orientation //
// test and point-in-sphere test implemented by Jonathan Shewchuk. //
// He generously provided the source code in the public domain, //
// http://www.cs.cmu.edu/~quake/robust.html. //
// predicates.cxx is a C++ version of the original C code. //
// //
// The original predicates of Shewchuk only use "dynamic filters", i.e., it //
// computes the error at runtime step by step. TetGen first uses a "static //
// filter" in each predicate. It estimates the maximal possible error in all //
// cases. It safely and quickly "filters" many easy cases. //
// //
//============================================================================//
void exactinit(int, int, int, REAL, REAL, REAL);
REAL orient3d(REAL *pa, REAL *pb, REAL *pc, REAL *pd);
REAL insphere(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe);
REAL orient4d(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe,
REAL ah, REAL bh, REAL ch, REAL dh, REAL eh);
REAL orient2dexact(REAL *pa, REAL *pb, REAL *pc);
REAL orient3dexact(REAL *pa, REAL *pb, REAL *pc, REAL *pd);
REAL orient4dexact(REAL* pa, REAL* pb, REAL* pc, REAL* pd, REAL* pe,
REAL ah, REAL bh, REAL ch, REAL dh, REAL eh);
//============================================================================//
// //
// tetgenmesh TetGen's internal mesh data structure. //
// //
// It uses a tetrahedron-based mesh data structure. It implements elementary //
// flip operations to locally modify the mesh. It implements basic meshing //
// algorithms to create Delaunay tetrahedraliations, to perform boundary //
// recovery, to place Steiner points in the mesh domain, and to optimize the //
// quality of the mesh. //
// //
//============================================================================//
class tetgenmesh {
public:
//============================================================================//
// //
// Mesh data structure //
// //
// A tetrahedral mesh T of a 3D piecewise linear complex (PLC) X is a 3D //
// simplicial complex whose underlying space is equal to the space of X. T //
// contains a 2D subcomplex S which is a triangular mesh of the boundary of //
// X. S contains a 1D subcomplex L which is a linear mesh of the boundary of //
// S. Faces and edges in S and L are respectively called subfaces and segme- //
// nts to distinguish them from others in T. //
// //
// TetGen uses a tetrahedron-based data structure. It stores tetrahedra and //
// vertices. This data structure is pointer-based. Each tetrahedron contains //
// pointers to its vertices and adjacent tetrahedra. Each vertex holds its x-,//
// y-, z-coordinates, and a pointer to one of the tetrahedra having it. Both //
// tetrahedra and vertices may contain user data. //
// //
// Let T be a tetrahedralization. Each triangular face of T belongs to either //
// two or one tetrahedron. In the latter case, it is an exterior boundary //
// face of T. TetGen attaches tetrahedra (one-to-one) to such faces. All such //
// tetrahedra contain an "infinite vertex" (which has no geometric coordinates//
// ). One can imagine such a vertex lies in 4D space and is visible by all //
// exterior boundary faces simultaneously. This extended set of tetrahedra //
// (including the infinite vertex) becomes a tetrahedralization of a 3-sphere //
// that has no boundary in 3d. It has the nice property that every triangular //
// face is shared by exactly two tetrahedra. //
// //
// The current version of TetGen stores explicitly the subfaces and segments //
// (which are in surface mesh S and the linear mesh L), respectively. Extra //
// pointers are allocated in tetrahedra and subfaces to point each other. //
// //
//============================================================================//
// The tetrahedron data structure. It includes the following fields:
// - a list of four adjoining tetrahedra;
// - a list of four vertices;
// - a pointer to a list of four subfaces (optional, for -p switch);
// - a pointer to a list of six segments (optional, for -p switch);
// - a list of user-defined floating-point attributes (optional);
// - a volume constraint (optional, for -a switch);
// - an integer of element marker (and flags);
// The structure of a tetrahedron is an array of pointers. Its actual size
// (the length of the array) is determined at runtime.
typedef REAL **tetrahedron;
// The subface data structure. It includes the following fields:
// - a list of three adjoining subfaces;
// - a list of three vertices;
// - a list of three adjoining segments;
// - two adjoining tetrahedra;
// - an area constraint (optional, for -q switch);
// - an integer for boundary marker;
// - an integer for type, flags, etc.
typedef REAL **shellface;
// The point data structure. It includes the following fields:
// - x, y and z coordinates;
// - a list of user-defined point attributes (optional);
// - u, v coordinates (optional, for -s switch);
// - a metric tensor (optional, for -q or -m switch);
// - a pointer to an adjacent tetrahedron;
// - a pointer to a parent (or a duplicate) point;
// - a pointer to an adjacent subface or segment (optional, -p switch);
// - a pointer to a tet in background mesh (optional, for -m switch);
// - an integer for boundary marker (point index);
// - an integer for point type (and flags).
// - an integer for geometry tag (optional, for -s switch).
// The structure of a point is an array of REALs. Its acutal size is
// determined at the runtime.
typedef REAL *point;
//============================================================================//
// //
// Handles //
// //
// Navigation and manipulation in a tetrahedralization are accomplished by //
// operating on structures referred as ``handles". A handle is a pair (t,v), //
// where t is a pointer to a tetrahedron, and v is a 4-bit integer, in the //
// range from 0 to 11. v is called the ``version'' of a tetrahedron, it rep- //
// resents a directed edge of a specific face of the tetrahedron. //
// //
// There are 12 even permutations of the four vertices, each of them corres- //
// ponds to a directed edge (a version) of the tetrahedron. The 12 versions //
// can be grouped into 4 distinct ``edge rings'' in 4 ``oriented faces'' of //
// this tetrahedron. One can encode each version (a directed edge) into a //
// 4-bit integer such that the two upper bits encode the index (from 0 to 2) //
// of this edge in the edge ring, and the two lower bits encode the index ( //
// from 0 to 3) of the oriented face which contains this edge. //
// //
// The four vertices of a tetrahedron are indexed from 0 to 3 (according to //
// their storage in the data structure). Give each face the same index as //
// the node opposite it in the tetrahedron. Denote the edge connecting face //
// i to face j as i/j. We number the twelve versions as follows: //
// //
// | edge 0 edge 1 edge 2 //
// --------|-------------------------------- //
// face 0 | 0 (0/1) 4 (0/3) 8 (0/2) //
// face 1 | 1 (1/2) 5 (1/3) 9 (1/0) //
// face 2 | 2 (2/3) 6 (2/1) 10 (2/0) //
// face 3 | 3 (3/0) 7 (3/1) 11 (3/2) //
// //
// Similarly, navigation and manipulation in a (boundary) triangulation are //
// done by using handles of triangles. Each handle is a pair (s, v), where s //
// is a pointer to a triangle, and v is a version in the range from 0 to 5. //
// Each version corresponds to a directed edge of this triangle. //
// //
// Number the three vertices of a triangle from 0 to 2 (according to their //
// storage in the data structure). Give each edge the same index as the node //
// opposite it in the triangle. The six versions of a triangle are: //
// //
// | edge 0 edge 1 edge 2 //
// ---------------|-------------------------- //