///////////////////////////////////////////////////////////////////////////////1// //2// TetGen //3// //4// A Quality Tetrahedral Mesh Generator and A 3D Delaunay Triangulator //5// //6// Version 1.5 //7// May 31, 2014 //8// //9// Copyright (C) 2002--2014 //10// //11// TetGen is freely available through the website: http://www.tetgen.org. //12// It may be copied, modified, and redistributed for non-commercial use. //13// Please consult the file LICENSE for the detailed copyright notices. //14// //15///////////////////////////////////////////////////////////////////////////////161718#ifndef tetgenH19#define tetgenH2021// To compile TetGen as a library instead of an executable program, define22// the TETLIBRARY symbol.2324// #define TETLIBRARY252627// TetGen default uses the double precision (64 bit) for a real number.28// Alternatively, one can use the single precision (32 bit) 'float' if the29// memory is limited.3031#define REAL double // #define REAL float3233// Maximum number of characters in a file name (including the null).3435#define FILENAMESIZE 10243637// Maximum number of chars in a line read from a file (including the null).3839#define INPUTLINESIZE 20484041// TetGen only uses the C standard library.4243#include <QtGlobal>44#include <stdio.h>45#include <stdlib.h>46#include <string.h>47#include <math.h>48#include <time.h>4950// 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, Microsoft53// Visual C++ 2003 -- 2008 (Visual C++ 7.1 - 9) doesn't ship with this header54// file. In such case, we can define them by ourself.55// Update (learned from Stack Overflow): Visual Studio 2010 and Visual C++ 201056// Express both have stdint.h5758// The following piece of code was provided by Steven Johnson (MIT). Define the59// symbol _MSC_VER if you are using Microsoft Visual C++. Moreover, define60// the _WIN64 symbol if you are running TetGen on Win64 systems.6162#ifdef _MSC_VER // Microsoft Visual C++63# ifdef _WIN6464typedef __int64 intptr_t;65typedef unsigned __int64 uintptr_t;66# else // not _WIN6467typedef int intptr_t;68typedef unsigned int uintptr_t;69# endif70#else // not Visual C++71# include <stdint.h>72#endif7374///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////109110class tetgenio {111112public:113114// A "polygon" describes a simple polygon (no holes). It is not necessarily115// convex. Each polygon contains a number of corners (points) and the same116// number of sides (edges). The points of the polygon must be given in117// either counterclockwise or clockwise order and they form a ring, so118// every two consecutive points forms an edge of the polygon.119typedef struct {120int *vertexlist;121int numberofvertices;122} polygon;123124// A "facet" describes a polygonal region possibly with holes, edges, and125// points floating in it. Each facet consists of a list of polygons and126// a list of hole points (which lie strictly inside holes).127typedef struct {128polygon *polygonlist;129int numberofpolygons;130REAL *holelist;131int numberofholes;132} facet;133134// A "voroedge" is an edge of the Voronoi diagram. It corresponds to a135// Delaunay face. Each voroedge is either a line segment connecting136// two Voronoi vertices or a ray starting from a Voronoi vertex to an137// "infinite vertex". 'v1' and 'v2' are two indices pointing to the138// list of Voronoi vertices. 'v1' must be non-negative, while 'v2' may139// be -1 if it is a ray, in this case, the unit normal of this ray is140// given in 'vnormal'.141typedef struct {142int v1, v2;143REAL vnormal[3];144} voroedge;145146// A "vorofacet" is an facet of the Voronoi diagram. It corresponds to a147// Delaunay edge. Each Voronoi facet is a convex polygon formed by a148// list of Voronoi edges, it may not be closed. 'c1' and 'c2' are two149// indices pointing into the list of Voronoi cells, i.e., the two cells150// share this facet. 'elist' is an array of indices pointing into the151// list of Voronoi edges, 'elist[0]' saves the number of Voronoi edges152// (including rays) of this facet.153typedef struct {154int c1, c2;155int *elist;156} vorofacet;157158159// Additional parameters associated with an input (or mesh) vertex.160// These informations are provided by CAD libraries.161typedef struct {162REAL uv[2];163int tag;164int type; // 0, 1, or 2.165} pointparam;166167// Callback functions for meshing PSCs.168typedef REAL (* GetVertexParamOnEdge)(void*, int, int);169typedef void (* GetSteinerOnEdge)(void*, int, REAL, REAL*);170typedef void (* GetVertexParamOnFace)(void*, int, int, REAL*);171typedef void (* GetEdgeSteinerParamOnFace)(void*, int, REAL, int, REAL*);172typedef void (* GetSteinerOnFace)(void*, int, REAL*, REAL*);173174// A callback function for mesh refinement.175typedef bool (* TetSizeFunc)(REAL*, REAL*, REAL*, REAL*, REAL*, REAL);176177// Items are numbered starting from 'firstnumber' (0 or 1), default is 0.178int firstnumber;179180// Dimension of the mesh (2 or 3), default is 3.181int mesh_dim;182183// Does the lines in .node file contain index or not, default is 1.184int useindex;185186// 'pointlist': An array of point coordinates. The first point's x187// coordinate is at index [0] and its y coordinate at index [1], its188// z coordinate is at index [2], followed by the coordinates of the189// remaining points. Each point occupies three REALs.190// 'pointattributelist': An array of point attributes. Each point's191// attributes occupy 'numberofpointattributes' REALs.192// 'pointmtrlist': An array of metric tensors at points. Each point's193// 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.196REAL *pointlist;197REAL *pointattributelist;198REAL *pointmtrlist;199int *pointmarkerlist;200int *point2tetlist;201pointparam *pointparamlist;202int numberofpoints;203int numberofpointattributes;204int numberofpointmtrs;205206// 'tetrahedronlist': An array of tetrahedron corners. The first207// tetrahedron's first corner is at index [0], followed by its other208// corners, followed by six nodes on the edges of the tetrahedron if the209// second order option (-o2) is applied. Each tetrahedron occupies210// 'numberofcorners' ints. The second order nodes are ouput only.211// 'tetrahedronattributelist': An array of tetrahedron attributes. Each212// tetrahedron's attributes occupy 'numberoftetrahedronattributes' REALs.213// 'tetrahedronvolumelist': An array of constraints, i.e. tetrahedron's214// 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.218int *tetrahedronlist;219REAL *tetrahedronattributelist;220REAL *tetrahedronvolumelist;221int *neighborlist;222int *tet2facelist;223int *tet2edgelist;224int numberoftetrahedra;225int numberofcorners;226int numberoftetrahedronattributes;227228// 'facetlist': An array of facets. Each entry is a structure of facet.229// 'facetmarkerlist': An array of facet markers; one int per facet.230facet *facetlist;231int *facetmarkerlist;232int numberoffacets;233234// 'holelist': An array of holes (in volume). Each hole is given by a235// seed (point) which lies strictly inside it. The first seed's x, y and z236// coordinates are at indices [0], [1] and [2], followed by the237// remaining seeds. Three REALs per hole.238REAL *holelist;239int numberofholes;240241// 'regionlist': An array of regions (subdomains). Each region is given by242// a seed (point) which lies strictly inside it. The first seed's x, y and243// z coordinates are at indices [0], [1] and [2], followed by the regional244// 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 the248// 'a' switch (with no number following).249REAL *regionlist;250int numberofregions;251252// 'facetconstraintlist': An array of facet constraints. Each constraint253// specifies a maximum area bound on the subfaces of that facet. The254// first facet constraint is given by a facet marker at index [0] and its255// maximum area bound at index [1], followed by the remaining facet con-256// straints. Two REALs per facet constraint. Note: the facet marker is257// actually an integer.258REAL *facetconstraintlist;259int numberoffacetconstraints;260261// 'segmentconstraintlist': An array of segment constraints. Each constraint262// specifies a maximum length bound on the subsegments of that segment.263// The first constraint is given by the two endpoints of the segment at264// index [0] and [1], and the maximum length bound at index [2], followed265// by the remaining segment constraints. Three REALs per constraint.266// Note the segment endpoints are actually integers.267REAL *segmentconstraintlist;268int numberofsegmentconstraints;269270271// 'trifacelist': An array of face (triangle) corners. The first face's272// three corners are at indices [0], [1] and [2], followed by the remaining273// 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. The277// 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.281int *trifacelist;282int *trifacemarkerlist;283int *o2facelist;284int *face2tetlist;285int *face2edgelist;286int numberoftrifaces;287288// 'edgelist': An array of edge endpoints. The first edge's endpoints289// 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 the293// second order option (-o2) is applied. One int per edge.294// 'edge2tetlist': An array of tetrahedra indices. One int per edge.295int *edgelist;296int *edgemarkerlist;297int *o2edgelist;298int *edge2tetlist;299int numberofedges;300301// '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 of305// indices pointing into 'vfacetlist'. The 0th entry is used to store306// the length of this array.307REAL *vpointlist;308voroedge *vedgelist;309vorofacet *vfacetlist;310int **vcelllist;311int numberofvpoints;312int numberofvedges;313int numberofvfacets;314int numberofvcells;315316317// Variable (and callback functions) for meshing PSCs.318void *geomhandle;319GetVertexParamOnEdge getvertexparamonedge;320GetSteinerOnEdge getsteineronedge;321GetVertexParamOnFace getvertexparamonface;322GetEdgeSteinerParamOnFace getedgesteinerparamonface;323GetSteinerOnFace getsteineronface;324325// A callback function.326TetSizeFunc tetunsuitable;327328// Input & output routines.329virtual bool load_node_call(FILE* infile, int markers, int uvflag, char*);330virtual bool load_node(char*);331virtual bool load_edge(char*);332virtual bool load_face(char*);333virtual bool load_tet(char*);334virtual bool load_vol(char*);335virtual bool load_var(char*);336virtual bool load_mtr(char*);337// virtual bool load_pbc(char*);338virtual bool load_poly(char*);339virtual bool load_off(char*);340virtual bool load_ply(char*);341virtual bool load_stl(char*);342virtual bool load_vtk(char*);343virtual bool load_medit(char*, int);344virtual bool load_plc(char*, int);345virtual bool load_tetmesh(char*, int);346virtual void save_nodes(char*);347virtual void save_elements(char*);348virtual void save_faces(char*);349virtual void save_edges(char*);350virtual void save_neighbors(char*);351virtual void save_poly(char*);352virtual void save_faces2smesh(char*);353354// Read line and parse string functions.355virtual char *readline(char* string, FILE* infile, int *linenumber);356virtual char *findnextfield(char* string);357virtual char *readnumberline(char* string, FILE* infile, char* infilename);358virtual char *findnextnumber(char* string);359360virtual void init(polygon* p) {361p->vertexlist = (int *) NULL;362p->numberofvertices = 0;363}364365virtual void init(facet* f) {366f->polygonlist = (polygon *) NULL;367f->numberofpolygons = 0;368f->holelist = (REAL *) NULL;369f->numberofholes = 0;370}371372// Initialize routine.373virtual void initialize()374{375firstnumber = 0;376mesh_dim = 3;377useindex = 1;378379pointlist = (REAL *) NULL;380pointattributelist = (REAL *) NULL;381pointmtrlist = (REAL *) NULL;382pointmarkerlist = (int *) NULL;383point2tetlist = (int *) NULL;384pointparamlist = (pointparam *) NULL;385numberofpoints = 0;386numberofpointattributes = 0;387numberofpointmtrs = 0;388389tetrahedronlist = (int *) NULL;390tetrahedronattributelist = (REAL *) NULL;391tetrahedronvolumelist = (REAL *) NULL;392neighborlist = (int *) NULL;393tet2facelist = (int *) NULL;394tet2edgelist = (int *) NULL;395numberoftetrahedra = 0;396numberofcorners = 4;397numberoftetrahedronattributes = 0;398399trifacelist = (int *) NULL;400trifacemarkerlist = (int *) NULL;401o2facelist = (int *) NULL;402face2tetlist = (int *) NULL;403face2edgelist = (int *) NULL;404numberoftrifaces = 0;405406edgelist = (int *) NULL;407edgemarkerlist = (int *) NULL;408o2edgelist = (int *) NULL;409edge2tetlist = (int *) NULL;410numberofedges = 0;411412facetlist = (facet *) NULL;413facetmarkerlist = (int *) NULL;414numberoffacets = 0;415416holelist = (REAL *) NULL;417numberofholes = 0;418419regionlist = (REAL *) NULL;420numberofregions = 0;421422facetconstraintlist = (REAL *) NULL;423numberoffacetconstraints = 0;424segmentconstraintlist = (REAL *) NULL;425numberofsegmentconstraints = 0;426427428vpointlist = (REAL *) NULL;429vedgelist = (voroedge *) NULL;430vfacetlist = (vorofacet *) NULL;431vcelllist = (int **) NULL;432numberofvpoints = 0;433numberofvedges = 0;434numberofvfacets = 0;435numberofvcells = 0;436437438tetunsuitable = NULL;439440geomhandle = NULL;441getvertexparamonedge = NULL;442getsteineronedge = NULL;443getvertexparamonface = NULL;444getedgesteinerparamonface = NULL;445getsteineronface = NULL;446}447448// Free the memory allocated in 'tetgenio'. Note that it assumes that the449// memory was allocated by the "new" operator (C++).450virtual void deinitialize()451{452int i, j;453454if (pointlist != (REAL *) NULL) {455delete [] pointlist;456}457if (pointattributelist != (REAL *) NULL) {458delete [] pointattributelist;459}460if (pointmtrlist != (REAL *) NULL) {461delete [] pointmtrlist;462}463if (pointmarkerlist != (int *) NULL) {464delete [] pointmarkerlist;465}466if (point2tetlist != (int *) NULL) {467delete [] point2tetlist;468}469if (pointparamlist != (pointparam *) NULL) {470delete [] pointparamlist;471}472473if (tetrahedronlist != (int *) NULL) {474delete [] tetrahedronlist;475}476if (tetrahedronattributelist != (REAL *) NULL) {477delete [] tetrahedronattributelist;478}479if (tetrahedronvolumelist != (REAL *) NULL) {480delete [] tetrahedronvolumelist;481}482if (neighborlist != (int *) NULL) {483delete [] neighborlist;484}485if (tet2facelist != (int *) NULL) {486delete [] tet2facelist;487}488if (tet2edgelist != (int *) NULL) {489delete [] tet2edgelist;490}491492if (trifacelist != (int *) NULL) {493delete [] trifacelist;494}495if (trifacemarkerlist != (int *) NULL) {496delete [] trifacemarkerlist;497}498if (o2facelist != (int *) NULL) {499delete [] o2facelist;500}501if (face2tetlist != (int *) NULL) {502delete [] face2tetlist;503}504if (face2edgelist != (int *) NULL) {505delete [] face2edgelist;506}507508if (edgelist != (int *) NULL) {509delete [] edgelist;510}511if (edgemarkerlist != (int *) NULL) {512delete [] edgemarkerlist;513}514if (o2edgelist != (int *) NULL) {515delete [] o2edgelist;516}517if (edge2tetlist != (int *) NULL) {518delete [] edge2tetlist;519}520521if (facetlist != (facet *) NULL) {522facet *f;523polygon *p;524for (i = 0; i < numberoffacets; i++) {525f = &facetlist[i];526for (j = 0; j < f->numberofpolygons; j++) {527p = &f->polygonlist[j];528delete [] p->vertexlist;529}530delete [] f->polygonlist;531if (f->holelist != (REAL *) NULL) {532delete [] f->holelist;533}534}535delete [] facetlist;536}537if (facetmarkerlist != (int *) NULL) {538delete [] facetmarkerlist;539}540541if (holelist != (REAL *) NULL) {542delete [] holelist;543}544if (regionlist != (REAL *) NULL) {545delete [] regionlist;546}547if (facetconstraintlist != (REAL *) NULL) {548delete [] facetconstraintlist;549}550if (segmentconstraintlist != (REAL *) NULL) {551delete [] segmentconstraintlist;552}553if (vpointlist != (REAL *) NULL) {554delete [] vpointlist;555}556if (vedgelist != (voroedge *) NULL) {557delete [] vedgelist;558}559if (vfacetlist != (vorofacet *) NULL) {560for (i = 0; i < numberofvfacets; i++) {561delete [] vfacetlist[i].elist;562}563delete [] vfacetlist;564}565if (vcelllist != (int **) NULL) {566for (i = 0; i < numberofvcells; i++) {567delete [] vcelllist[i];568}569delete [] vcelllist;570}571}572573// Constructor & destructor.574tetgenio() {initialize();}575virtual ~tetgenio() {deinitialize();}576577}; // class tetgenio578579///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////595596class tetgenbehavior {597598public:599600// Switches of TetGen.601int plc; // '-p', 0.602int psc; // '-s', 0.603int refine; // '-r', 0.604int quality; // '-q', 0.605int nobisect; // '-Y', 0.606int coarsen; // '-R', 0.607int weighted; // '-w', 0.608int brio_hilbert; // '-b', 1.609int incrflip; // '-l', 0.610int flipinsert; // '-L', 0.611int metric; // '-m', 0.612int varvolume; // '-a', 0.613int fixedvolume; // '-a', 0.614int regionattrib; // '-A', 0.615int cdtrefine; // '-D', 0.616int insertaddpoints; // '-i', 0.617int diagnose; // '-d', 0.618int convex; // '-c', 0.619int nomergefacet; // '-M', 0.620int nomergevertex; // '-M', 0.621int noexact; // '-X', 0.622int nostaticfilter; // '-X', 0.623int zeroindex; // '-z', 0.624int facesout; // '-f', 0.625int edgesout; // '-e', 0.626int neighout; // '-n', 0.627int voroout; // '-v', 0.628int meditview; // '-g', 0.629int vtkview; // '-k', 0.630int nobound; // '-B', 0.631int nonodewritten; // '-N', 0.632int noelewritten; // '-E', 0.633int nofacewritten; // '-F', 0.634int noiterationnum; // '-I', 0.635int nojettison; // '-J', 0.636int docheck; // '-C', 0.637int quiet; // '-Q', 0.638int verbose; // '-V', 0.639640// Parameters of TetGen.641int vertexperblock; // '-x', 4092.642int tetrahedraperblock; // '-x', 8188.643int shellfaceperblock; // '-x', 2044.644int nobisect_nomerge; // '-Y', 1.645int supsteiner_level; // '-Y/', 2.646int addsteiner_algo; // '-Y//', 1.647int coarsen_param; // '-R', 0.648int weighted_param; // '-w', 0.649int fliplinklevel; // -1.650int flipstarsize; // -1.651int fliplinklevelinc; // 1.652int reflevel; // '-D', 3.653int optlevel; // '-O', 2.654int optscheme; // '-O', 7.655int delmaxfliplevel; // 1.656int order; // '-o', 1.657int reversetetori; // '-o/', 0.658int steinerleft; // '-S', 0.659int no_sort; // 0.660int hilbert_order; // '-b///', 52.661int hilbert_limit; // '-b//' 8.662int brio_threshold; // '-b' 64.663REAL brio_ratio; // '-b/' 0.125.664REAL facet_separate_ang_tol; // '-p', 179.9.665REAL facet_overlap_ang_tol; // '-p/', 0.1.666REAL facet_small_ang_tol; // '-p//', 15.0.667REAL maxvolume; // '-a', -1.0.668REAL minratio; // '-q', 0.0.669REAL mindihedral; // '-q', 5.0.670REAL optmaxdihedral; // 165.0.671REAL optminsmtdihed; // 179.0.672REAL optminslidihed; // 179.0.673REAL epsilon; // '-T', 1.0e-8.674REAL coarsen_percent; // -R1/#, 1.0.675676// Strings of command line arguments and input/output file names.677char commandline[1024];678char infilename[1024];679char outfilename[1024];680char addinfilename[1024];681char bgmeshfilename[1024];682683// The input object of TetGen. They are recognized by either the input684// 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 switch694// (-p or -r) implies the object.695enum objecttype {NODES, POLY, OFF, PLY, STL, MEDIT, VTK, MESH} object;696697698virtual void syntax();699virtual void usage();700701// Command line parse routine.702virtual bool parse_commandline(int argc, char **argv);703virtual bool parse_commandline(char *switches) {704return parse_commandline(0, &switches);705}706707// Initialize all variables.708tetgenbehavior()709{710plc = 0;711psc = 0;712refine = 0;713quality = 0;714nobisect = 0;715coarsen = 0;716metric = 0;717weighted = 0;718brio_hilbert = 1;719incrflip = 0;720flipinsert = 0;721varvolume = 0;722fixedvolume = 0;723noexact = 0;724nostaticfilter = 0;725insertaddpoints = 0;726regionattrib = 0;727cdtrefine = 0;728diagnose = 0;729convex = 0;730zeroindex = 0;731facesout = 0;732edgesout = 0;733neighout = 0;734voroout = 0;735meditview = 0;736vtkview = 0;737nobound = 0;738nonodewritten = 0;739noelewritten = 0;740nofacewritten = 0;741noiterationnum = 0;742nomergefacet = 0;743nomergevertex = 0;744nojettison = 0;745docheck = 0;746quiet = 0;747verbose = 0;748749vertexperblock = 4092;750tetrahedraperblock = 8188;751shellfaceperblock = 4092;752nobisect_nomerge = 1;753supsteiner_level = 2;754addsteiner_algo = 1;755coarsen_param = 0;756weighted_param = 0;757fliplinklevel = -1;758flipstarsize = -1;759fliplinklevelinc = 1;760reflevel = 3;761optscheme = 7;762optlevel = 2;763delmaxfliplevel = 1;764order = 1;765reversetetori = 0;766steinerleft = -1;767no_sort = 0;768hilbert_order = 52; //-1;769hilbert_limit = 8;770brio_threshold = 64;771brio_ratio = 0.125;772facet_separate_ang_tol = 179.9;773facet_overlap_ang_tol = 0.1;774facet_small_ang_tol = 15.0;775maxvolume = -1.0;776minratio = 2.0;777mindihedral = 0.0;778optmaxdihedral = 165.00;779optminsmtdihed = 179.00;780optminslidihed = 179.00;781epsilon = 1.0e-8;782coarsen_percent = 1.0;783object = NODES;784785commandline[0] = '\0';786infilename[0] = '\0';787outfilename[0] = '\0';788addinfilename[0] = '\0';789bgmeshfilename[0] = '\0';790791}792793}; // class tetgenbehavior794795///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////821822void exactinit(int, int, int, REAL, REAL, REAL);823REAL orient3d(REAL *pa, REAL *pb, REAL *pc, REAL *pd);824REAL insphere(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe);825REAL orient4d(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe,826REAL ah, REAL bh, REAL ch, REAL dh, REAL eh);827828///////////////////////////////////////////////////////////////////////////////829// //830// tetgenmesh //831// //832// A structure for creating and updating tetrahedral meshes. //833// //834///////////////////////////////////////////////////////////////////////////////835836class tetgenmesh {837838public:839840///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////870871// 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 size880// (the length of the array) is determined at runtime.881882typedef REAL **tetrahedron;883884// 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.892893typedef REAL **shellface;894895// 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 acutal size is908// determined at the runtime.909910typedef REAL *point;911912///////////////////////////////////////////////////////////////////////////////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 orieation | 0 2 4 //954// cw orieation | 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///////////////////////////////////////////////////////////////////////////////960961class triface {962public:963tetrahedron *tet;964int ver; // Range from 0 to 11.965triface() : tet(0), ver(0) {}966triface& operator=(const triface& t) {967tet = t.tet; ver = t.ver;968return *this;969}970};971972class face {973public:974shellface *sh;975int shver; // Range from 0 to 5.976face() : sh(0), shver(0) {}977face& operator=(const face& s) {978sh = s.sh; shver = s.shver;979return *this;980}981};982983///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////10001001class arraypool {10021003public:10041005int objectbytes;1006int objectsperblock;1007int log2objectsperblock;1008int objectsperblockmark;1009int toparraylen;1010char **toparray;1011long objects;1012unsigned long totalmemory;10131014void restart();1015void poolinit(int sizeofobject, int log2objperblk);1016char* getblock(int objectindex);1017void* lookup(int objectindex);1018int newindex(void **newptr);10191020arraypool(int sizeofobject, int log2objperblk);1021~arraypool();1022};10231024// fastlookup() -- A fast, unsafe operation. Return the pointer to the object1025// with a given index. Note: The object's block must have been allocated,1026// i.e., by the function newindex().10271028#define fastlookup(pool, index) \1029(void *) ((pool)->toparray[(index) >> (pool)->log2objectsperblock] + \1030((index) & (pool)->objectsperblockmark) * (pool)->objectbytes)10311032///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////10521053class memorypool {10541055public:10561057void **firstblock, **nowblock;1058void *nextitem;1059void *deaditemstack;1060void **pathblock;1061void *pathitem;1062int alignbytes;1063int itembytes, itemwords;1064int itemsperblock;1065long items, maxitems;1066int unallocateditems;1067int pathitemsleft;10681069memorypool();1070memorypool(int, int, int, int);1071~memorypool();10721073void poolinit(int, int, int, int);1074void restart();1075void *alloc();1076void dealloc(void*);1077void traversalinit();1078void *traverse();1079};10801081///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////10941095class badface {1096public:1097triface tt;1098face ss;1099REAL key, cent[6]; // circumcenter or cos(dihedral angles) at 6 edges.1100point forg, fdest, fapex, foppo, noppo;1101badface *nextitem;1102badface() : key(0), forg(0), fdest(0), fapex(0), foppo(0), noppo(0),1103nextitem(0) {}1104};11051106///////////////////////////////////////////////////////////////////////////////1107// //1108// insertvertexflags //1109// //1110// A collection of flags that pass to the routine insertvertex(). //1111// //1112///////////////////////////////////////////////////////////////////////////////11131114class insertvertexflags {11151116public:11171118int iloc; // input/output.1119int bowywat, lawson;1120int splitbdflag, validflag, respectbdflag;1121int rejflag, chkencflag, cdtflag;1122int assignmeshsize;1123int sloc, sbowywat;11241125// Used by Delaunay refinement.1126int refineflag; // 0, 1, 2, 31127triface refinetet;1128face refinesh;1129int smlenflag; // for useinsertradius.1130REAL smlen; // for useinsertradius.1131point parentpt;11321133insertvertexflags() {1134iloc = bowywat = lawson = 0;1135splitbdflag = validflag = respectbdflag = 0;1136rejflag = chkencflag = cdtflag = 0;1137assignmeshsize = 0;1138sloc = sbowywat = 0;11391140refineflag = 0;1141refinetet.tet = NULL;1142refinesh.sh = NULL;1143smlenflag = 0;1144smlen = 0.0;1145}1146};11471148///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////11561157class flipconstraints {11581159public:11601161// Elementary flip flags.1162int enqflag; // (= flipflag)1163int chkencflag;11641165// Control flags1166int unflip; // Undo the performed flips.1167int collectnewtets; // Collect the new tets created by flips.1168int collectencsegflag;11691170// Optimization flags.1171int remove_ndelaunay_edge; // Remove a non-Delaunay edge.1172REAL bak_tetprism_vol; // The value to be minimized.1173REAL tetprism_vol_sum;1174int remove_large_angle; // Remove a large dihedral angle at edge.1175REAL cosdihed_in; // The input cosine of the dihedral angle (> 0).1176REAL cosdihed_out; // The improved cosine of the dihedral angle.11771178// Boundary recovery flags.1179int checkflipeligibility;1180point seg[2]; // A constraining edge to be recovered.1181point fac[3]; // A constraining face to be recovered.1182point remvert; // A vertex to be removed.118311841185flipconstraints() {1186enqflag = 0;1187chkencflag = 0;11881189unflip = 0;1190collectnewtets = 0;1191collectencsegflag = 0;11921193remove_ndelaunay_edge = 0;1194bak_tetprism_vol = 0.0;1195tetprism_vol_sum = 0.0;1196remove_large_angle = 0;1197cosdihed_in = 0.0;1198cosdihed_out = 0.0;11991200checkflipeligibility = 0;1201seg[0] = NULL;1202fac[0] = NULL;1203remvert = NULL;1204}1205};12061207///////////////////////////////////////////////////////////////////////////////1208// //1209// optparameters //1210// //1211// Optimization options and parameters. //1212// //1213///////////////////////////////////////////////////////////////////////////////12141215class optparameters {12161217public:12181219// The one of goals of optimization.1220int max_min_volume; // Maximize the minimum volume.1221int min_max_aspectratio; // Minimize the maximum aspect ratio.1222int min_max_dihedangle; // Minimize the maximum dihedral angle.12231224// The initial and improved value.1225REAL initval, imprval;12261227int numofsearchdirs;1228REAL searchstep;1229int maxiter; // Maximum smoothing iterations (disabled by -1).1230int smthiter; // Performed iterations.123112321233optparameters() {1234max_min_volume = 0;1235min_max_aspectratio = 0;1236min_max_dihedangle = 0;12371238initval = imprval = 0.0;12391240numofsearchdirs = 10;1241searchstep = 0.01;1242maxiter = -1; // Unlimited smoothing iterations.1243smthiter = 0;12441245}1246};124712481249///////////////////////////////////////////////////////////////////////////////1250// //1251// Labels (enumeration declarations) used by TetGen. //1252// //1253///////////////////////////////////////////////////////////////////////////////12541255// Labels that signify the type of a vertex.1256enum verttype {UNUSEDVERTEX, DUPLICATEDVERTEX, RIDGEVERTEX, ACUTEVERTEX,1257FACETVERTEX, VOLVERTEX, FREESEGVERTEX, FREEFACETVERTEX,1258FREEVOLVERTEX, NREGULARVERTEX, DEADVERTEX};12591260// Labels that signify the result of triangle-triangle intersection test.1261enum interresult {DISJOINT, INTERSECT, SHAREVERT, SHAREEDGE, SHAREFACE,1262TOUCHEDGE, TOUCHFACE, ACROSSVERT, ACROSSEDGE, ACROSSFACE};12631264// Labels that signify the result of point location.1265enum locateresult {TGUNKNOWN, OUTSIDE, INTETRAHEDRON, ONFACE, ONEDGE, ONVERTEX,1266ENCVERTEX, ENCSEGMENT, ENCSUBFACE, NEARVERTEX, NONREGULAR,1267INSTAR, BADELEMENT};12681269///////////////////////////////////////////////////////////////////////////////1270// //1271// Variables of TetGen //1272// //1273///////////////////////////////////////////////////////////////////////////////12741275// Pointer to the input data (a set of nodes, a PLC, or a mesh).1276tetgenio *in, *addin;12771278// Pointer to the switches and parameters.1279tetgenbehavior *b;12801281// Pointer to a background mesh (contains size specification map).1282tetgenmesh *bgm;12831284// Memorypools to store mesh elements (points, tetrahedra, subfaces, and1285// segments) and extra pointers between tetrahedra, subfaces, and segments.1286memorypool *tetrahedrons, *subfaces, *subsegs, *points;1287memorypool *tet2subpool, *tet2segpool;12881289// Memorypools to store bad-quality (or encroached) elements.1290memorypool *badtetrahedrons, *badsubfacs, *badsubsegs;12911292// A memorypool to store faces to be flipped.1293memorypool *flippool;1294arraypool *unflipqueue;1295badface *flipstack;12961297// Arrays used for point insertion (the Bowyer-Watson algorithm).1298arraypool *cavetetlist, *cavebdrylist, *caveoldtetlist;1299arraypool *cavetetshlist, *cavetetseglist, *cavetetvertlist;1300arraypool *caveencshlist, *caveencseglist;1301arraypool *caveshlist, *caveshbdlist, *cavesegshlist;13021303// Stacks used for CDT construction and boundary recovery.1304arraypool *subsegstack, *subfacstack, *subvertstack;13051306// Arrays of encroached segments and subfaces (for mesh refinement).1307arraypool *encseglist, *encshlist;13081309// The map between facets to their vertices (for mesh refinement).1310int *idx2facetlist;1311point *facetverticeslist;13121313// The map between segments to their endpoints (for mesh refinement).1314point *segmentendpointslist;13151316// The infinite vertex.1317point dummypoint;1318// The recently visited tetrahedron, subface.1319triface recenttet;1320face recentsh;13211322// PI is the ratio of a circle's circumference to its diameter.1323static REAL PI;13241325// Array (size = numberoftetrahedra * 6) for storing high-order nodes of1326// tetrahedra (only used when -o2 switch is selected).1327point *highordertable;13281329// Various variables.1330int numpointattrib; // Number of point attributes.1331int numelemattrib; // Number of tetrahedron attributes.1332int sizeoftensor; // Number of REALs per metric tensor.1333int pointmtrindex; // Index to find the metric tensor of a point.1334int pointparamindex; // Index to find the u,v coordinates of a point.1335int point2simindex; // Index to find a simplex adjacent to a point.1336int pointmarkindex; // Index to find boundary marker of a point.1337int pointinsradiusindex; // Index to find the insertion radius of a point.1338int elemattribindex; // Index to find attributes of a tetrahedron.1339int volumeboundindex; // Index to find volume bound of a tetrahedron.1340int elemmarkerindex; // Index to find marker of a tetrahedron.1341int shmarkindex; // Index to find boundary marker of a subface.1342int areaboundindex; // Index to find area bound of a subface.1343int checksubsegflag; // Are there segments in the tetrahedralization yet?1344int checksubfaceflag; // Are there subfaces in the tetrahedralization yet?1345int checkconstraints; // Are there variant (node, seg, facet) constraints?1346int nonconvex; // Is current mesh non-convex?1347int autofliplinklevel; // The increase of link levels, default is 1.1348int useinsertradius; // Save the insertion radius for Steiner points.1349long samples; // Number of random samples for point location.1350unsigned long randomseed; // Current random number seed.1351REAL cosmaxdihed, cosmindihed; // The cosine values of max/min dihedral.1352REAL cossmtdihed; // The cosine value of a bad dihedral to be smoothed.1353REAL cosslidihed; // The cosine value of the max dihedral of a sliver.1354REAL minfaceang, minfacetdihed; // The minimum input (dihedral) angles.1355REAL tetprism_vol_sum; // The total volume of tetrahedral-prisms (in 4D).1356REAL longest; // The longest possible edge length.1357REAL minedgelength; // = longest * b->epsion.1358REAL xmax, xmin, ymax, ymin, zmax, zmin; // Bounding box of points.13591360// Counters.1361long insegments; // Number of input segments.1362long hullsize; // Number of exterior boundary faces.1363long meshedges; // Number of mesh edges.1364long meshhulledges; // Number of boundary mesh edges.1365long steinerleft; // Number of Steiner points not yet used.1366long dupverts; // Are there duplicated vertices?1367long unuverts; // Are there unused vertices?1368long nonregularcount; // Are there non-regular vertices?1369long st_segref_count, st_facref_count, st_volref_count; // Steiner points.1370long fillregioncount, cavitycount, cavityexpcount;1371long flip14count, flip26count, flipn2ncount;1372long flip23count, flip32count, flip44count, flip41count;1373long flip31count, flip22count;1374unsigned long totalworkmemory; // Total memory used by working arrays.137513761377///////////////////////////////////////////////////////////////////////////////1378// //1379// Mesh manipulation primitives //1380// //1381///////////////////////////////////////////////////////////////////////////////13821383// Fast lookup tables for mesh manipulation primitives.1384static int bondtbl[12][12], fsymtbl[12][12];1385static int esymtbl[12], enexttbl[12], eprevtbl[12];1386static int enextesymtbl[12], eprevesymtbl[12];1387static int eorgoppotbl[12], edestoppotbl[12];1388static int facepivot1[12], facepivot2[12][12];1389static int orgpivot[12], destpivot[12], apexpivot[12], oppopivot[12];1390static int tsbondtbl[12][6], stbondtbl[12][6];1391static int tspivottbl[12][6], stpivottbl[12][6];1392static int ver2edge[12], edge2ver[6], epivot[12];1393static int sorgpivot [6], sdestpivot[6], sapexpivot[6];1394static int snextpivot[6];13951396void inittables();13971398// Primitives for tetrahedra.1399inline tetrahedron encode(triface& t);1400inline tetrahedron encode2(tetrahedron* ptr, int ver);1401inline void decode(tetrahedron ptr, triface& t);1402inline void bond(triface& t1, triface& t2);1403inline void dissolve(triface& t);1404inline void esym(triface& t1, triface& t2);1405inline void esymself(triface& t);1406inline void enext(triface& t1, triface& t2);1407inline void enextself(triface& t);1408inline void eprev(triface& t1, triface& t2);1409inline void eprevself(triface& t);1410inline void enextesym(triface& t1, triface& t2);1411inline void enextesymself(triface& t);1412inline void eprevesym(triface& t1, triface& t2);1413inline void eprevesymself(triface& t);1414inline void eorgoppo(triface& t1, triface& t2);1415inline void eorgoppoself(triface& t);1416inline void edestoppo(triface& t1, triface& t2);1417inline void edestoppoself(triface& t);1418inline void fsym(triface& t1, triface& t2);1419inline void fsymself(triface& t);1420inline void fnext(triface& t1, triface& t2);1421inline void fnextself(triface& t);1422inline point org (triface& t);1423inline point dest(triface& t);1424inline point apex(triface& t);1425inline point oppo(triface& t);1426inline void setorg (triface& t, point p);1427inline void setdest(triface& t, point p);1428inline void setapex(triface& t, point p);1429inline void setoppo(triface& t, point p);1430inline REAL elemattribute(tetrahedron* ptr, int attnum);1431inline void setelemattribute(tetrahedron* ptr, int attnum, REAL value);1432inline REAL volumebound(tetrahedron* ptr);1433inline void setvolumebound(tetrahedron* ptr, REAL value);1434inline int elemindex(tetrahedron* ptr);1435inline void setelemindex(tetrahedron* ptr, int value);1436inline int elemmarker(tetrahedron* ptr);1437inline void setelemmarker(tetrahedron* ptr, int value);1438inline void infect(triface& t);1439inline void uninfect(triface& t);1440inline bool infected(triface& t);1441inline void marktest(triface& t);1442inline void unmarktest(triface& t);1443inline bool marktested(triface& t);1444inline void markface(triface& t);1445inline void unmarkface(triface& t);1446inline bool facemarked(triface& t);1447inline void markedge(triface& t);1448inline void unmarkedge(triface& t);1449inline bool edgemarked(triface& t);1450inline void marktest2(triface& t);1451inline void unmarktest2(triface& t);1452inline bool marktest2ed(triface& t);1453inline int elemcounter(triface& t);1454inline void setelemcounter(triface& t, int value);1455inline void increaseelemcounter(triface& t);1456inline void decreaseelemcounter(triface& t);1457inline bool ishulltet(triface& t);1458inline bool isdeadtet(triface& t);14591460// Primitives for subfaces and subsegments.1461inline void sdecode(shellface sptr, face& s);1462inline shellface sencode(face& s);1463inline shellface sencode2(shellface *sh, int shver);1464inline void spivot(face& s1, face& s2);1465inline void spivotself(face& s);1466inline void sbond(face& s1, face& s2);1467inline void sbond1(face& s1, face& s2);1468inline void sdissolve(face& s);1469inline point sorg(face& s);1470inline point sdest(face& s);1471inline point sapex(face& s);1472inline void setsorg(face& s, point pointptr);1473inline void setsdest(face& s, point pointptr);1474inline void setsapex(face& s, point pointptr);1475inline void sesym(face& s1, face& s2);1476inline void sesymself(face& s);1477inline void senext(face& s1, face& s2);1478inline void senextself(face& s);1479inline void senext2(face& s1, face& s2);1480inline void senext2self(face& s);1481inline REAL areabound(face& s);1482inline void setareabound(face& s, REAL value);1483inline int shellmark(face& s);1484inline void setshellmark(face& s, int value);1485inline void sinfect(face& s);1486inline void suninfect(face& s);1487inline bool sinfected(face& s);1488inline void smarktest(face& s);1489inline void sunmarktest(face& s);1490inline bool smarktested(face& s);1491inline void smarktest2(face& s);1492inline void sunmarktest2(face& s);1493inline bool smarktest2ed(face& s);1494inline void smarktest3(face& s);1495inline void sunmarktest3(face& s);1496inline bool smarktest3ed(face& s);1497inline void setfacetindex(face& f, int value);1498inline int getfacetindex(face& f);14991500// Primitives for interacting tetrahedra and subfaces.1501inline void tsbond(triface& t, face& s);1502inline void tsdissolve(triface& t);1503inline void stdissolve(face& s);1504inline void tspivot(triface& t, face& s);1505inline void stpivot(face& s, triface& t);15061507// Primitives for interacting tetrahedra and segments.1508inline void tssbond1(triface& t, face& seg);1509inline void sstbond1(face& s, triface& t);1510inline void tssdissolve1(triface& t);1511inline void sstdissolve1(face& s);1512inline void tsspivot1(triface& t, face& s);1513inline void sstpivot1(face& s, triface& t);15141515// Primitives for interacting subfaces and segments.1516inline void ssbond(face& s, face& edge);1517inline void ssbond1(face& s, face& edge);1518inline void ssdissolve(face& s);1519inline void sspivot(face& s, face& edge);15201521// Primitives for points.1522inline int pointmark(point pt);1523inline void setpointmark(point pt, int value);1524inline enum verttype pointtype(point pt);1525inline void setpointtype(point pt, enum verttype value);1526inline int pointgeomtag(point pt);1527inline void setpointgeomtag(point pt, int value);1528inline REAL pointgeomuv(point pt, int i);1529inline void setpointgeomuv(point pt, int i, REAL value);1530inline void pinfect(point pt);1531inline void puninfect(point pt);1532inline bool pinfected(point pt);1533inline void pmarktest(point pt);1534inline void punmarktest(point pt);1535inline bool pmarktested(point pt);1536inline void pmarktest2(point pt);1537inline void punmarktest2(point pt);1538inline bool pmarktest2ed(point pt);1539inline void pmarktest3(point pt);1540inline void punmarktest3(point pt);1541inline bool pmarktest3ed(point pt);1542inline tetrahedron point2tet(point pt);1543inline void setpoint2tet(point pt, tetrahedron value);1544inline shellface point2sh(point pt);1545inline void setpoint2sh(point pt, shellface value);1546inline point point2ppt(point pt);1547inline void setpoint2ppt(point pt, point value);1548inline tetrahedron point2bgmtet(point pt);1549inline void setpoint2bgmtet(point pt, tetrahedron value);1550inline void setpointinsradius(point pt, REAL value);1551inline REAL getpointinsradius(point pt);1552inline bool issteinerpoint(point pt);15531554// Advanced primitives.1555inline void point2tetorg(point pt, triface& t);1556inline void point2shorg(point pa, face& s);1557inline point farsorg(face& seg);1558inline point farsdest(face& seg);15591560///////////////////////////////////////////////////////////////////////////////1561// //1562// Memory managment //1563// //1564///////////////////////////////////////////////////////////////////////////////15651566void tetrahedrondealloc(tetrahedron*);1567tetrahedron *tetrahedrontraverse();1568tetrahedron *alltetrahedrontraverse();1569void shellfacedealloc(memorypool*, shellface*);1570shellface *shellfacetraverse(memorypool*);1571void pointdealloc(point);1572point pointtraverse();15731574void makeindex2pointmap(point*&);1575void makepoint2submap(memorypool*, int*&, face*&);1576void maketetrahedron(triface*);1577void makeshellface(memorypool*, face*);1578void makepoint(point*, enum verttype);15791580void initializepools();15811582///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////16051606// Symbolic perturbations (robust)1607REAL insphere_s(REAL*, REAL*, REAL*, REAL*, REAL*);1608REAL orient4d_s(REAL*, REAL*, REAL*, REAL*, REAL*,1609REAL, REAL, REAL, REAL, REAL);16101611// Triangle-edge intersection test (robust)1612int tri_edge_2d(point, point, point, point, point, point, int, int*, int*);1613int tri_edge_tail(point, point, point, point, point, point, REAL, REAL, int,1614int*, int*);1615int tri_edge_test(point, point, point, point, point, point, int, int*, int*);16161617// Triangle-triangle intersection test (robust)1618int tri_edge_inter_tail(point, point, point, point, point, REAL, REAL);1619int tri_tri_inter(point, point, point, point, point, point);16201621// Linear algebra functions1622inline REAL dot(REAL* v1, REAL* v2);1623inline void cross(REAL* v1, REAL* v2, REAL* n);1624bool lu_decmp(REAL lu[4][4], int n, int* ps, REAL* d, int N);1625void lu_solve(REAL lu[4][4], int n, int* ps, REAL* b, int N);16261627// An embedded 2-dimensional geometric predicate (non-robust)1628REAL incircle3d(point pa, point pb, point pc, point pd);16291630// Geometric calculations (non-robust)1631REAL orient3dfast(REAL *pa, REAL *pb, REAL *pc, REAL *pd);1632inline REAL norm2(REAL x, REAL y, REAL z);1633inline REAL distance(REAL* p1, REAL* p2);1634void facenormal(point pa, point pb, point pc, REAL *n, int pivot, REAL *lav);1635REAL shortdistance(REAL* p, REAL* e1, REAL* e2);1636REAL triarea(REAL* pa, REAL* pb, REAL* pc);1637REAL interiorangle(REAL* o, REAL* p1, REAL* p2, REAL* n);1638void projpt2edge(REAL* p, REAL* e1, REAL* e2, REAL* prj);1639void projpt2face(REAL* p, REAL* f1, REAL* f2, REAL* f3, REAL* prj);1640bool tetalldihedral(point, point, point, point, REAL*, REAL*, REAL*);1641void tetallnormal(point, point, point, point, REAL N[4][3], REAL* volume);1642REAL tetaspectratio(point, point, point, point);1643bool circumsphere(REAL*, REAL*, REAL*, REAL*, REAL* cent, REAL* radius);1644bool orthosphere(REAL*,REAL*,REAL*,REAL*,REAL,REAL,REAL,REAL,REAL*,REAL*);1645void planelineint(REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*);1646int linelineint(REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*);1647REAL tetprismvol(REAL* pa, REAL* pb, REAL* pc, REAL* pd);1648bool calculateabovepoint(arraypool*, point*, point*, point*);1649void calculateabovepoint4(point, point, point, point);16501651// PLC error reports.1652void report_overlapping_facets(face*, face*, REAL dihedang = 0.0);1653int report_selfint_edge(point, point, face* sedge, triface* searchtet,1654enum interresult);1655int report_selfint_face(point, point, point, face* sface, triface* iedge,1656int intflag, int* types, int* poss);16571658///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////16891690// The elementary flips.1691void flip23(triface*, int, flipconstraints* fc);1692void flip32(triface*, int, flipconstraints* fc);1693void flip41(triface*, int, flipconstraints* fc);16941695// A generalized edge flip.1696int flipnm(triface*, int n, int level, int, flipconstraints* fc);1697int flipnm_post(triface*, int n, int nn, int, flipconstraints* fc);16981699// Point insertion.1700int insertpoint(point, triface*, face*, face*, insertvertexflags*);1701void insertpoint_abort(face*, insertvertexflags*);17021703///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////17591760void transfernodes();17611762// Point sorting.1763int transgc[8][3][8], tsb1mod3[8];1764void hilbert_init(int n);1765int hilbert_split(point* vertexarray, int arraysize, int gc0, int gc1,1766REAL, REAL, REAL, REAL, REAL, REAL);1767void hilbert_sort3(point* vertexarray, int arraysize, int e, int d,1768REAL, REAL, REAL, REAL, REAL, REAL, int depth);1769void brio_multiscale_sort(point*,int,int threshold,REAL ratio,int* depth);17701771// Point location.1772unsigned long randomnation(unsigned int choices);1773void randomsample(point searchpt, triface *searchtet);1774enum locateresult locate(point searchpt, triface *searchtet,1775int chkencflag = 0);17761777// Incremental flips.1778void flippush(badface*&, triface*);1779int incrementalflip(point newpt, int, flipconstraints *fc);17801781// Incremental Delaunay construction.1782void initialdelaunay(point pa, point pb, point pc, point pd);1783void incrementaldelaunay(clock_t&);17841785///////////////////////////////////////////////////////////////////////////////1786// //1787// Surface triangulation //1788// //1789///////////////////////////////////////////////////////////////////////////////17901791void flipshpush(face*);1792void flip22(face*, int, int);1793void flip31(face*, int);1794long lawsonflip();1795int sinsertvertex(point newpt, face*, face*, int iloc, int bowywat, int);1796int sremovevertex(point delpt, face*, face*, int lawson);17971798enum locateresult slocate(point, face*, int, int, int);1799enum interresult sscoutsegment(face*, point, int, int, int);1800void scarveholes(int, REAL*);1801int triangulate(int, arraypool*, arraypool*, int, REAL*);18021803void unifysegments();1804void identifyinputedges(point*);1805void mergefacets();1806void meshsurface();18071808void interecursive(shellface** subfacearray, int arraysize, int axis,1809REAL, REAL, REAL, REAL, REAL, REAL, int* internum);1810void detectinterfaces();18111812///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////18861887void makesegmentendpointsmap();18881889enum interresult finddirection(triface* searchtet, point endpt);1890enum interresult scoutsegment(point, point, face*, triface*, point*,1891arraypool*);1892int getsteinerptonsegment(face* seg, point refpt, point steinpt);1893void delaunizesegments();18941895int scoutsubface(face* searchsh,triface* searchtet,int shflag);1896void formregion(face*, arraypool*, arraypool*, arraypool*);1897int scoutcrossedge(triface& crosstet, arraypool*, arraypool*);1898bool formcavity(triface*, arraypool*, arraypool*, arraypool*, arraypool*,1899arraypool*, arraypool*);1900// Facet recovery by cavity re-triangulation [Si and Gaertner 2011].1901void delaunizecavity(arraypool*, arraypool*, arraypool*, arraypool*,1902arraypool*, arraypool*);1903bool fillcavity(arraypool*, arraypool*, arraypool*, arraypool*,1904arraypool*, arraypool*, triface* crossedge);1905void carvecavity(arraypool*, arraypool*, arraypool*);1906void restorecavity(arraypool*, arraypool*, arraypool*, arraypool*);1907// Facet recovery by flips [Shewchuk 2003].1908void flipcertify(triface *chkface, badface **pqueue, point, point, point);1909void flipinsertfacet(arraypool*, arraypool*, arraypool*, arraypool*);19101911int insertpoint_cdt(point, triface*, face*, face*, insertvertexflags*,1912arraypool*, arraypool*, arraypool*, arraypool*,1913arraypool*, arraypool*);1914void refineregion(face&, arraypool*, arraypool*, arraypool*, arraypool*,1915arraypool*, arraypool*);1916void constrainedfacets();19171918void constraineddelaunay(clock_t&);19191920///////////////////////////////////////////////////////////////////////////////1921// //1922// Constrained tetrahedralizations. //1923// //1924///////////////////////////////////////////////////////////////////////////////19251926int checkflipeligibility(int fliptype, point, point, point, point, point,1927int level, int edgepivot, flipconstraints* fc);19281929int removeedgebyflips(triface*, flipconstraints*);1930int removefacebyflips(triface*, flipconstraints*);19311932int recoveredgebyflips(point, point, face*, triface*, int fullsearch);1933int add_steinerpt_in_schoenhardtpoly(triface*, int, int chkencflag);1934int add_steinerpt_in_segment(face*, int searchlevel);1935int addsteiner4recoversegment(face*, int);1936int recoversegments(arraypool*, int fullsearch, int steinerflag);19371938int recoverfacebyflips(point, point, point, face*, triface*);1939int recoversubfaces(arraypool*, int steinerflag);19401941int getvertexstar(int, point searchpt, arraypool*, arraypool*, arraypool*);1942int getedge(point, point, triface*);1943int reduceedgesatvertex(point startpt, arraypool* endptlist);1944int removevertexbyflips(point steinerpt);19451946int suppressbdrysteinerpoint(point steinerpt);1947int suppresssteinerpoints();19481949void recoverboundary(clock_t&);19501951///////////////////////////////////////////////////////////////////////////////1952// //1953// Mesh reconstruction //1954// //1955///////////////////////////////////////////////////////////////////////////////19561957void carveholes();19581959void reconstructmesh();19601961int scoutpoint(point, triface*, int randflag);1962REAL getpointmeshsize(point, triface*, int iloc);1963void interpolatemeshsize();19641965void insertconstrainedpoints(point *insertarray, int arylen, int rejflag);1966void insertconstrainedpoints(tetgenio *addio);19671968void collectremovepoints(arraypool *remptlist);1969void meshcoarsening();19701971///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////20062007void makefacetverticesmap();2008int segsegadjacent(face *, face *);2009int segfacetadjacent(face *checkseg, face *checksh);2010int facetfacetadjacent(face *, face *);2011void save_segmentpoint_insradius(point segpt, point parentpt, REAL r);2012void save_facetpoint_insradius(point facpt, point parentpt, REAL r);2013void enqueuesubface(memorypool*, face*);2014void enqueuetetrahedron(triface*);20152016int checkseg4encroach(point pa, point pb, point checkpt);2017int checkseg4split(face *chkseg, point&, int&);2018int splitsegment(face *splitseg, point encpt, REAL, point, point, int, int);2019void repairencsegs(int chkencflag);20202021int checkfac4encroach(point, point, point, point checkpt, REAL*, REAL*);2022int checkfac4split(face *chkfac, point& encpt, int& qflag, REAL *ccent);2023int splitsubface(face *splitfac, point, point, int qflag, REAL *ccent, int);2024void repairencfacs(int chkencflag);20252026int checktet4split(triface *chktet, int& qflag, REAL *ccent);2027int splittetrahedron(triface* splittet,int qflag,REAL *ccent, int);2028void repairbadtets(int chkencflag);20292030void delaunayrefinement();20312032///////////////////////////////////////////////////////////////////////////////2033// //2034// Mesh optimization //2035// //2036///////////////////////////////////////////////////////////////////////////////20372038long lawsonflip3d(flipconstraints *fc);2039void recoverdelaunay();20402041int gettetrahedron(point, point, point, point, triface *);2042long improvequalitybyflips();20432044int smoothpoint(point smtpt, arraypool*, int ccw, optparameters *opm);2045long improvequalitybysmoothing(optparameters *opm);20462047int splitsliver(triface *, REAL, int);2048long removeslivers(int);20492050void optimizemesh();20512052///////////////////////////////////////////////////////////////////////////////2053// //2054// Mesh check and statistics //2055// //2056///////////////////////////////////////////////////////////////////////////////20572058// Mesh validations.2059int checkmesh(int topoflag);2060int checkshells();2061int checksegments();2062int checkdelaunay(int perturb = 1);2063int checkregular(int);2064int checkconforming(int);20652066// Mesh statistics.2067void printfcomma(unsigned long n);2068void qualitystatistics();2069void memorystatistics();2070void statistics();20712072///////////////////////////////////////////////////////////////////////////////2073// //2074// Mesh output //2075// //2076///////////////////////////////////////////////////////////////////////////////20772078void jettisonnodes();2079void highorder();2080void indexelements();2081void numberedges();2082void outnodes(tetgenio*);2083void outmetrics(tetgenio*);2084void outelements(tetgenio*);2085void outfaces(tetgenio*);2086void outhullfaces(tetgenio*);2087void outsubfaces(tetgenio*);2088void outedges(tetgenio*);2089void outsubsegments(tetgenio*);2090void outneighbors(tetgenio*);2091void outvoronoi(tetgenio*);2092void outsmesh(char*);2093void outmesh2medit(char*);2094void outmesh2vtk(char*);2095209620972098///////////////////////////////////////////////////////////////////////////////2099// //2100// Constructor & destructor //2101// //2102///////////////////////////////////////////////////////////////////////////////21032104void initializetetgenmesh()2105{2106in = addin = NULL;2107b = NULL;2108bgm = NULL;21092110tetrahedrons = subfaces = subsegs = points = NULL;2111badtetrahedrons = badsubfacs = badsubsegs = NULL;2112tet2segpool = tet2subpool = NULL;2113flippool = NULL;21142115dummypoint = NULL;2116flipstack = NULL;2117unflipqueue = NULL;21182119cavetetlist = cavebdrylist = caveoldtetlist = NULL;2120cavetetshlist = cavetetseglist = cavetetvertlist = NULL;2121caveencshlist = caveencseglist = NULL;2122caveshlist = caveshbdlist = cavesegshlist = NULL;21232124subsegstack = subfacstack = subvertstack = NULL;2125encseglist = encshlist = NULL;2126idx2facetlist = NULL;2127facetverticeslist = NULL;2128segmentendpointslist = NULL;21292130highordertable = NULL;21312132numpointattrib = numelemattrib = 0;2133sizeoftensor = 0;2134pointmtrindex = 0;2135pointparamindex = 0;2136pointmarkindex = 0;2137point2simindex = 0;2138pointinsradiusindex = 0;2139elemattribindex = 0;2140volumeboundindex = 0;2141shmarkindex = 0;2142areaboundindex = 0;2143checksubsegflag = 0;2144checksubfaceflag = 0;2145checkconstraints = 0;2146nonconvex = 0;2147autofliplinklevel = 1;2148useinsertradius = 0;2149samples = 0l;2150randomseed = 1l;2151minfaceang = minfacetdihed = PI;2152tetprism_vol_sum = 0.0;2153longest = minedgelength = 0.0;2154xmax = xmin = ymax = ymin = zmax = zmin = 0.0;21552156insegments = 0l;2157hullsize = 0l;2158meshedges = meshhulledges = 0l;2159steinerleft = -1;2160dupverts = 0l;2161unuverts = 0l;2162nonregularcount = 0l;2163st_segref_count = st_facref_count = st_volref_count = 0l;2164fillregioncount = cavitycount = cavityexpcount = 0l;2165flip14count = flip26count = flipn2ncount = 0l;2166flip23count = flip32count = flip44count = flip41count = 0l;2167flip22count = flip31count = 0l;2168totalworkmemory = 0l;216921702171} // tetgenmesh()21722173void freememory()2174{2175if (bgm != NULL) {2176delete bgm;2177}21782179if (points != (memorypool *) NULL) {2180delete points;2181delete [] dummypoint;2182}2183if (tetrahedrons != (memorypool *) NULL) {2184delete tetrahedrons;2185}2186if (subfaces != (memorypool *) NULL) {2187delete subfaces;2188delete subsegs;2189}2190if (tet2segpool != NULL) {2191delete tet2segpool;2192delete tet2subpool;2193}21942195if (badtetrahedrons) {2196delete badtetrahedrons;2197}2198if (badsubfacs) {2199delete badsubfacs;2200}2201if (badsubsegs) {2202delete badsubsegs;2203}2204if (encseglist) {2205delete encseglist;2206}2207if (encshlist) {2208delete encshlist;2209}22102211if (flippool != NULL) {2212delete flippool;2213delete unflipqueue;2214}22152216if (cavetetlist != NULL) {2217delete cavetetlist;2218delete cavebdrylist;2219delete caveoldtetlist;2220delete cavetetvertlist;2221}22222223if (caveshlist != NULL) {2224delete caveshlist;2225delete caveshbdlist;2226delete cavesegshlist;2227delete cavetetshlist;2228delete cavetetseglist;2229delete caveencshlist;2230delete caveencseglist;2231}22322233if (subsegstack != NULL) {2234delete subsegstack;2235delete subfacstack;2236delete subvertstack;2237}22382239if (idx2facetlist != NULL) {2240delete [] idx2facetlist;2241delete [] facetverticeslist;2242}22432244if (segmentendpointslist != NULL) {2245delete [] segmentendpointslist;2246}22472248if (highordertable != NULL) {2249delete [] highordertable;2250}22512252initializetetgenmesh();2253}22542255tetgenmesh()2256{2257initializetetgenmesh();2258}22592260~tetgenmesh()2261{2262freememory();2263} // ~tetgenmesh()22642265}; // End of class tetgenmesh.22662267///////////////////////////////////////////////////////////////////////////////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///////////////////////////////////////////////////////////////////////////////22812282void tetrahedralize(tetgenbehavior *b, tetgenio *in, tetgenio *out,2283tetgenio *addin = NULL, tetgenio *bgmin = NULL);22842285#ifdef TETLIBRARY2286void tetrahedralize(char *switches, tetgenio *in, tetgenio *out,2287tetgenio *addin = NULL, tetgenio *bgmin = NULL);22882289extern "C"2290#ifdef Q_WS_WIN2291__declspec(dllexport)2292#endif2293void delegate_tetrahedralize(int bs, tetgenbehavior *b, char *switches,2294tetgenio *in, tetgenio *out, tetgenio *addin = NULL, tetgenio *bgmin = NULL);22952296#endif // #ifdef TETLIBRARY22972298///////////////////////////////////////////////////////////////////////////////2299// //2300// terminatetetgen() Terminate TetGen with a given exit code. //2301// //2302///////////////////////////////////////////////////////////////////////////////23032304inline void terminatetetgen(tetgenmesh *m, int x)2305{2306#ifdef TETLIBRARY2307throw x;2308#else2309switch (x) {2310case 1: // Out of memory.2311printf("Error: Out of memory.\n");2312break;2313case 2: // Encounter an internal error.2314printf("Please report this bug to [email protected]. Include\n");2315printf(" the message above, your input data set, and the exact\n");2316printf(" command line you used to run this program, thank you.\n");2317break;2318case 3:2319printf("A self-intersection was detected. Program stopped.\n");2320printf("Hint: use -d option to detect all self-intersections.\n");2321break;2322case 4:2323printf("A very small input feature size was detected. Program stopped.\n");2324if (m) {2325printf("Hint: use -T option to set a smaller tolerance. Current is %g\n",2326m->b->epsilon);2327}2328break;2329case 5:2330printf("Two very close input facets were detected. Program stopped.\n");2331printf("Hint: use -Y option to avoid adding Steiner points in boundary.\n");2332break;2333case 10:2334printf("An input error was detected. Program stopped.\n");2335break;2336} // switch (x)2337exit(x);2338#endif // #ifdef TETLIBRARY2339}23402341///////////////////////////////////////////////////////////////////////////////2342// //2343// Primitives for tetrahedra //2344// //2345///////////////////////////////////////////////////////////////////////////////23462347// encode() compress a handle into a single pointer. It relies on the2348// assumption that all addresses of tetrahedra are aligned to sixteen-2349// byte boundaries, so that the last four significant bits are zero.23502351inline tetgenmesh::tetrahedron tetgenmesh::encode(triface& t) {2352return (tetrahedron) ((uintptr_t) (t).tet | (uintptr_t) (t).ver);2353}23542355inline tetgenmesh::tetrahedron tetgenmesh::encode2(tetrahedron* ptr, int ver) {2356return (tetrahedron) ((uintptr_t) (ptr) | (uintptr_t) (ver));2357}23582359// decode() converts a pointer to a handle. The version is extracted from2360// the four least significant bits of the pointer.23612362inline 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}23662367// bond() connects two tetrahedra together. (t1,v1) and (t2,v2) must2368// refer to the same face and the same edge.23692370inline void tetgenmesh::bond(triface& t1, triface& t2) {2371t1.tet[t1.ver & 3] = encode2(t2.tet, bondtbl[t1.ver][t2.ver]);2372t2.tet[t2.ver & 3] = encode2(t1.tet, bondtbl[t2.ver][t1.ver]);2373}237423752376// dissolve() a bond (from one side).23772378inline void tetgenmesh::dissolve(triface& t) {2379t.tet[t.ver & 3] = NULL;2380}23812382// enext() finds the next edge (counterclockwise) in the same face.23832384inline void tetgenmesh::enext(triface& t1, triface& t2) {2385t2.tet = t1.tet;2386t2.ver = enexttbl[t1.ver];2387}23882389inline void tetgenmesh::enextself(triface& t) {2390t.ver = enexttbl[t.ver];2391}23922393// eprev() finds the next edge (clockwise) in the same face.23942395inline void tetgenmesh::eprev(triface& t1, triface& t2) {2396t2.tet = t1.tet;2397t2.ver = eprevtbl[t1.ver];2398}23992400inline void tetgenmesh::eprevself(triface& t) {2401t.ver = eprevtbl[t.ver];2402}24032404// esym() finds the reversed edge. It is in the other face of the2405// same tetrahedron.24062407inline void tetgenmesh::esym(triface& t1, triface& t2) {2408(t2).tet = (t1).tet;2409(t2).ver = esymtbl[(t1).ver];2410}24112412inline void tetgenmesh::esymself(triface& t) {2413(t).ver = esymtbl[(t).ver];2414}24152416// enextesym() finds the reversed edge of the next edge. It is in the other2417// face of the same tetrahedron. It is the combination esym() * enext().24182419inline void tetgenmesh::enextesym(triface& t1, triface& t2) {2420t2.tet = t1.tet;2421t2.ver = enextesymtbl[t1.ver];2422}24232424inline void tetgenmesh::enextesymself(triface& t) {2425t.ver = enextesymtbl[t.ver];2426}24272428// eprevesym() finds the reversed edge of the previous edge.24292430inline void tetgenmesh::eprevesym(triface& t1, triface& t2) {2431t2.tet = t1.tet;2432t2.ver = eprevesymtbl[t1.ver];2433}24342435inline void tetgenmesh::eprevesymself(triface& t) {2436t.ver = eprevesymtbl[t.ver];2437}24382439// eorgoppo() Finds the opposite face of the origin of the current edge.2440// Return the opposite edge of the current edge.24412442inline void tetgenmesh::eorgoppo(triface& t1, triface& t2) {2443t2.tet = t1.tet;2444t2.ver = eorgoppotbl[t1.ver];2445}24462447inline void tetgenmesh::eorgoppoself(triface& t) {2448t.ver = eorgoppotbl[t.ver];2449}24502451// edestoppo() Finds the opposite face of the destination of the current2452// edge. Return the opposite edge of the current edge.24532454inline void tetgenmesh::edestoppo(triface& t1, triface& t2) {2455t2.tet = t1.tet;2456t2.ver = edestoppotbl[t1.ver];2457}24582459inline void tetgenmesh::edestoppoself(triface& t) {2460t.ver = edestoppotbl[t.ver];2461}24622463// fsym() finds the adjacent tetrahedron at the same face and the same edge.24642465inline void tetgenmesh::fsym(triface& t1, triface& t2) {2466decode((t1).tet[(t1).ver & 3], t2);2467t2.ver = fsymtbl[t1.ver][t2.ver];2468}246924702471#define fsymself(t) \2472t1ver = (t).ver; \2473decode((t).tet[(t).ver & 3], (t));\2474(t).ver = fsymtbl[t1ver][(t).ver]24752476// fnext() finds the next face while rotating about an edge according to2477// a right-hand rule. The face is in the adjacent tetrahedron. It is2478// the combination: fsym() * esym().24792480inline void tetgenmesh::fnext(triface& t1, triface& t2) {2481decode(t1.tet[facepivot1[t1.ver]], t2);2482t2.ver = facepivot2[t1.ver][t2.ver];2483}248424852486#define fnextself(t) \2487t1ver = (t).ver; \2488decode((t).tet[facepivot1[(t).ver]], (t)); \2489(t).ver = facepivot2[t1ver][(t).ver]249024912492// The following primtives get or set the origin, destination, face apex,2493// or face opposite of an ordered tetrahedron.24942495inline tetgenmesh::point tetgenmesh::org(triface& t) {2496return (point) (t).tet[orgpivot[(t).ver]];2497}24982499inline tetgenmesh::point tetgenmesh:: dest(triface& t) {2500return (point) (t).tet[destpivot[(t).ver]];2501}25022503inline tetgenmesh::point tetgenmesh:: apex(triface& t) {2504return (point) (t).tet[apexpivot[(t).ver]];2505}25062507inline tetgenmesh::point tetgenmesh:: oppo(triface& t) {2508return (point) (t).tet[oppopivot[(t).ver]];2509}25102511inline void tetgenmesh:: setorg(triface& t, point p) {2512(t).tet[orgpivot[(t).ver]] = (tetrahedron) (p);2513}25142515inline void tetgenmesh:: setdest(triface& t, point p) {2516(t).tet[destpivot[(t).ver]] = (tetrahedron) (p);2517}25182519inline void tetgenmesh:: setapex(triface& t, point p) {2520(t).tet[apexpivot[(t).ver]] = (tetrahedron) (p);2521}25222523inline void tetgenmesh:: setoppo(triface& t, point p) {2524(t).tet[oppopivot[(t).ver]] = (tetrahedron) (p);2525}25262527#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)25322533// Check or set a tetrahedron's attributes.25342535inline REAL tetgenmesh::elemattribute(tetrahedron* ptr, int attnum) {2536return ((REAL *) (ptr))[elemattribindex + attnum];2537}25382539inline void tetgenmesh::setelemattribute(tetrahedron* ptr, int attnum,2540REAL value) {2541((REAL *) (ptr))[elemattribindex + attnum] = value;2542}25432544// Check or set a tetrahedron's maximum volume bound.25452546inline REAL tetgenmesh::volumebound(tetrahedron* ptr) {2547return ((REAL *) (ptr))[volumeboundindex];2548}25492550inline void tetgenmesh::setvolumebound(tetrahedron* ptr, REAL value) {2551((REAL *) (ptr))[volumeboundindex] = value;2552}25532554// Get or set a tetrahedron's index (only used for output).2555// These two routines use the reserved slot ptr[10].25562557inline int tetgenmesh::elemindex(tetrahedron* ptr) {2558int *iptr = (int *) &(ptr[10]);2559return iptr[0];2560}25612562inline void tetgenmesh::setelemindex(tetrahedron* ptr, int value) {2563int *iptr = (int *) &(ptr[10]);2564iptr[0] = value;2565}25662567// Get or set a tetrahedron's marker.2568// Set 'value = 0' cleans all the face/edge flags.25692570inline int tetgenmesh::elemmarker(tetrahedron* ptr) {2571return ((int *) (ptr))[elemmarkerindex];2572}25732574inline void tetgenmesh::setelemmarker(tetrahedron* ptr, int value) {2575((int *) (ptr))[elemmarkerindex] = value;2576}25772578// infect(), infected(), uninfect() -- primitives to flag or unflag a2579// tetrahedron. The last bit of the element marker is flagged (1)2580// or unflagged (0).25812582inline void tetgenmesh::infect(triface& t) {2583((int *) (t.tet))[elemmarkerindex] |= 1;2584}25852586inline void tetgenmesh::uninfect(triface& t) {2587((int *) (t.tet))[elemmarkerindex] &= ~1;2588}25892590inline bool tetgenmesh::infected(triface& t) {2591return (((int *) (t.tet))[elemmarkerindex] & 1) != 0;2592}25932594// marktest(), marktested(), unmarktest() -- primitives to flag or unflag a2595// tetrahedron. Use the second lowerest bit of the element marker.25962597inline void tetgenmesh::marktest(triface& t) {2598((int *) (t.tet))[elemmarkerindex] |= 2;2599}26002601inline void tetgenmesh::unmarktest(triface& t) {2602((int *) (t.tet))[elemmarkerindex] &= ~2;2603}26042605inline bool tetgenmesh::marktested(triface& t) {2606return (((int *) (t.tet))[elemmarkerindex] & 2) != 0;2607}26082609// markface(), unmarkface(), facemarked() -- primitives to flag or unflag a2610// face of a tetrahedron. From the last 3rd to 6th bits are used for2611// face markers, e.g., the last third bit corresponds to loc = 0.26122613inline void tetgenmesh::markface(triface& t) {2614((int *) (t.tet))[elemmarkerindex] |= (4 << (t.ver & 3));2615}26162617inline void tetgenmesh::unmarkface(triface& t) {2618((int *) (t.tet))[elemmarkerindex] &= ~(4 << (t.ver & 3));2619}26202621inline bool tetgenmesh::facemarked(triface& t) {2622return (((int *) (t.tet))[elemmarkerindex] & (4 << (t.ver & 3))) != 0;2623}26242625// markedge(), unmarkedge(), edgemarked() -- primitives to flag or unflag an2626// edge of a tetrahedron. From the last 7th to 12th bits are used for2627// 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.26292630inline void tetgenmesh::markedge(triface& t) {2631((int *) (t.tet))[elemmarkerindex] |= (int) (64 << ver2edge[(t).ver]);2632}26332634inline void tetgenmesh::unmarkedge(triface& t) {2635((int *) (t.tet))[elemmarkerindex] &= ~(int) (64 << ver2edge[(t).ver]);2636}26372638inline bool tetgenmesh::edgemarked(triface& t) {2639return (((int *) (t.tet))[elemmarkerindex] &2640(int) (64 << ver2edge[(t).ver])) != 0;2641}26422643// marktest2(), unmarktest2(), marktest2ed() -- primitives to flag and unflag2644// a tetrahedron. The 13th bit (2^12 = 4096) is used for this flag.26452646inline void tetgenmesh::marktest2(triface& t) {2647((int *) (t.tet))[elemmarkerindex] |= (int) (4096);2648}26492650inline void tetgenmesh::unmarktest2(triface& t) {2651((int *) (t.tet))[elemmarkerindex] &= ~(int) (4096);2652}26532654inline bool tetgenmesh::marktest2ed(triface& t) {2655return (((int *) (t.tet))[elemmarkerindex] & (int) (4096)) != 0;2656}26572658// elemcounter(), setelemcounter() -- primitives to read or ser a (small)2659// integer counter in this tet. It is saved from the 16th bit. On 32 bit2660// system, the range of the counter is [0, 2^15 = 32768].26612662inline int tetgenmesh::elemcounter(triface& t) {2663return (((int *) (t.tet))[elemmarkerindex]) >> 16;2664}26652666inline void tetgenmesh::setelemcounter(triface& t, int value) {2667int c = ((int *) (t.tet))[elemmarkerindex];2668// Clear the old counter while keep the other flags.2669c &= 65535; // sum_{i=0^15} 2^i2670c |= (value << 16);2671((int *) (t.tet))[elemmarkerindex] = c;2672}26732674inline void tetgenmesh::increaseelemcounter(triface& t) {2675int c = elemcounter(t);2676setelemcounter(t, c + 1);2677}26782679inline void tetgenmesh::decreaseelemcounter(triface& t) {2680int c = elemcounter(t);2681setelemcounter(t, c - 1);2682}26832684// ishulltet() tests if t is a hull tetrahedron.26852686inline bool tetgenmesh::ishulltet(triface& t) {2687return (point) (t).tet[7] == dummypoint;2688}26892690// isdeadtet() tests if t is a tetrahedron is dead.26912692inline bool tetgenmesh::isdeadtet(triface& t) {2693return ((t.tet == NULL) || (t.tet[4] == NULL));2694}26952696///////////////////////////////////////////////////////////////////////////////2697// //2698// Primitives for subfaces and subsegments //2699// //2700///////////////////////////////////////////////////////////////////////////////27012702// Each subface contains three pointers to its neighboring subfaces, with2703// edge versions. To save memory, both information are kept in a single2704// pointer. To make this possible, all subfaces are aligned to eight-byte2705// boundaries, so that the last three bits of each pointer are zeros. An2706// edge version (in the range 0 to 5) is compressed into the last three2707// bits of each pointer by 'sencode()'. 'sdecode()' decodes a pointer,2708// extracting an edge version and a pointer to the beginning of a subface.27092710inline void tetgenmesh::sdecode(shellface sptr, face& s) {2711s.shver = (int) ((uintptr_t) (sptr) & (uintptr_t) 7);2712s.sh = (shellface *) ((uintptr_t) (sptr) ^ (uintptr_t) (s.shver));2713}27142715inline tetgenmesh::shellface tetgenmesh::sencode(face& s) {2716return (shellface) ((uintptr_t) s.sh | (uintptr_t) s.shver);2717}27182719inline tetgenmesh::shellface tetgenmesh::sencode2(shellface *sh, int shver) {2720return (shellface) ((uintptr_t) sh | (uintptr_t) shver);2721}27222723// sbond() bonds two subfaces (s1) and (s2) together. s1 and s2 must refer2724// to the same edge. No requirement is needed on their orientations.27252726inline void tetgenmesh::sbond(face& s1, face& s2)2727{2728s1.sh[s1.shver >> 1] = sencode(s2);2729s2.sh[s2.shver >> 1] = sencode(s1);2730}27312732// 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.27352736inline void tetgenmesh::sbond1(face& s1, face& s2)2737{2738s1.sh[s1.shver >> 1] = sencode(s2);2739}27402741// Dissolve a subface bond (from one side). Note that the other subface2742// will still think it's connected to this subface.27432744inline void tetgenmesh::sdissolve(face& s)2745{2746s.sh[s.shver >> 1] = NULL;2747}27482749// spivot() finds the adjacent subface (s2) for a given subface (s1).2750// s1 and s2 share at the same edge.27512752inline void tetgenmesh::spivot(face& s1, face& s2)2753{2754shellface sptr = s1.sh[s1.shver >> 1];2755sdecode(sptr, s2);2756}27572758inline void tetgenmesh::spivotself(face& s)2759{2760shellface sptr = s.sh[s.shver >> 1];2761sdecode(sptr, s);2762}27632764// These primitives determine or set the origin, destination, or apex2765// of a subface with respect to the edge version.27662767inline tetgenmesh::point tetgenmesh::sorg(face& s)2768{2769return (point) s.sh[sorgpivot[s.shver]];2770}27712772inline tetgenmesh::point tetgenmesh::sdest(face& s)2773{2774return (point) s.sh[sdestpivot[s.shver]];2775}27762777inline tetgenmesh::point tetgenmesh::sapex(face& s)2778{2779return (point) s.sh[sapexpivot[s.shver]];2780}27812782inline void tetgenmesh::setsorg(face& s, point pointptr)2783{2784s.sh[sorgpivot[s.shver]] = (shellface) pointptr;2785}27862787inline void tetgenmesh::setsdest(face& s, point pointptr)2788{2789s.sh[sdestpivot[s.shver]] = (shellface) pointptr;2790}27912792inline void tetgenmesh::setsapex(face& s, point pointptr)2793{2794s.sh[sapexpivot[s.shver]] = (shellface) pointptr;2795}27962797#define setshvertices(s, pa, pb, pc)\2798setsorg(s, pa);\2799setsdest(s, pb);\2800setsapex(s, pc)28012802// sesym() reserves the direction of the lead edge.28032804inline void tetgenmesh::sesym(face& s1, face& s2)2805{2806s2.sh = s1.sh;2807s2.shver = (s1.shver ^ 1); // Inverse the last bit.2808}28092810inline void tetgenmesh::sesymself(face& s)2811{2812s.shver ^= 1;2813}28142815// senext() finds the next edge (counterclockwise) in the same orientation2816// of this face.28172818inline void tetgenmesh::senext(face& s1, face& s2)2819{2820s2.sh = s1.sh;2821s2.shver = snextpivot[s1.shver];2822}28232824inline void tetgenmesh::senextself(face& s)2825{2826s.shver = snextpivot[s.shver];2827}28282829inline void tetgenmesh::senext2(face& s1, face& s2)2830{2831s2.sh = s1.sh;2832s2.shver = snextpivot[snextpivot[s1.shver]];2833}28342835inline void tetgenmesh::senext2self(face& s)2836{2837s.shver = snextpivot[snextpivot[s.shver]];2838}283928402841// Check or set a subface's maximum area bound.28422843inline REAL tetgenmesh::areabound(face& s)2844{2845return ((REAL *) (s.sh))[areaboundindex];2846}28472848inline void tetgenmesh::setareabound(face& s, REAL value)2849{2850((REAL *) (s.sh))[areaboundindex] = value;2851}28522853// These two primitives read or set a shell marker. Shell markers are used2854// to hold user boundary information.28552856inline int tetgenmesh::shellmark(face& s)2857{2858return ((int *) (s.sh))[shmarkindex];2859}28602861inline void tetgenmesh::setshellmark(face& s, int value)2862{2863((int *) (s.sh))[shmarkindex] = value;2864}2865286628672868// sinfect(), sinfected(), suninfect() -- primitives to flag or unflag a2869// subface. The last bit of ((int *) ((s).sh))[shmarkindex+1] is flagged.28702871inline void tetgenmesh::sinfect(face& s)2872{2873((int *) ((s).sh))[shmarkindex+1] =2874(((int *) ((s).sh))[shmarkindex+1] | (int) 1);2875}28762877inline void tetgenmesh::suninfect(face& s)2878{2879((int *) ((s).sh))[shmarkindex+1] =2880(((int *) ((s).sh))[shmarkindex+1] & ~(int) 1);2881}28822883// Test a subface for viral infection.28842885inline bool tetgenmesh::sinfected(face& s)2886{2887return (((int *) ((s).sh))[shmarkindex+1] & (int) 1) != 0;2888}28892890// smarktest(), smarktested(), sunmarktest() -- primitives to flag or unflag2891// a subface. The last 2nd bit of the integer is flagged.28922893inline void tetgenmesh::smarktest(face& s)2894{2895((int *) ((s).sh))[shmarkindex+1] =2896(((int *)((s).sh))[shmarkindex+1] | (int) 2);2897}28982899inline void tetgenmesh::sunmarktest(face& s)2900{2901((int *) ((s).sh))[shmarkindex+1] =2902(((int *)((s).sh))[shmarkindex+1] & ~(int)2);2903}29042905inline bool tetgenmesh::smarktested(face& s)2906{2907return ((((int *) ((s).sh))[shmarkindex+1] & (int) 2) != 0);2908}29092910// smarktest2(), smarktest2ed(), sunmarktest2() -- primitives to flag or2911// unflag a subface. The last 3rd bit of the integer is flagged.29122913inline void tetgenmesh::smarktest2(face& s)2914{2915((int *) ((s).sh))[shmarkindex+1] =2916(((int *)((s).sh))[shmarkindex+1] | (int) 4);2917}29182919inline void tetgenmesh::sunmarktest2(face& s)2920{2921((int *) ((s).sh))[shmarkindex+1] =2922(((int *)((s).sh))[shmarkindex+1] & ~(int)4);2923}29242925inline bool tetgenmesh::smarktest2ed(face& s)2926{2927return ((((int *) ((s).sh))[shmarkindex+1] & (int) 4) != 0);2928}29292930// The last 4th bit of ((int *) ((s).sh))[shmarkindex+1] is flagged.29312932inline void tetgenmesh::smarktest3(face& s)2933{2934((int *) ((s).sh))[shmarkindex+1] =2935(((int *)((s).sh))[shmarkindex+1] | (int) 8);2936}29372938inline void tetgenmesh::sunmarktest3(face& s)2939{2940((int *) ((s).sh))[shmarkindex+1] =2941(((int *)((s).sh))[shmarkindex+1] & ~(int)8);2942}29432944inline bool tetgenmesh::smarktest3ed(face& s)2945{2946return ((((int *) ((s).sh))[shmarkindex+1] & (int) 8) != 0);2947}294829492950// Each facet has a unique index (automatically indexed). Starting from '0'.2951// We save this index in the same field of the shell type.29522953inline void tetgenmesh::setfacetindex(face& s, int value)2954{2955((int *) (s.sh))[shmarkindex + 2] = value;2956}29572958inline int tetgenmesh::getfacetindex(face& s)2959{2960return ((int *) (s.sh))[shmarkindex + 2];2961}29622963///////////////////////////////////////////////////////////////////////////////2964// //2965// Primitives for interacting between tetrahedra and subfaces //2966// //2967///////////////////////////////////////////////////////////////////////////////29682969// 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 bond2973// the edge in s which corresponds to t's 0th edge, and vice versa.29742975inline void tetgenmesh::tsbond(triface& t, face& s)2976{2977if ((t).tet[9] == NULL) {2978// Allocate space for this tet.2979(t).tet[9] = (tetrahedron) tet2subpool->alloc();2980// Initialize.2981for (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] =2987sencode2((s).sh, tsbondtbl[t.ver][s.shver]);2988// Bond s <== t.2989s.sh[9 + ((s).shver & 1)] =2990(shellface) encode2((t).tet, stbondtbl[t.ver][s.shver]);2991}29922993// tspivot() finds a subface (s) abutting on the given tetrahdera (t).2994// Return s.sh = NULL if there is no subface at t. Otherwise, return2995// the subface s, and s and t must be at the same edge wth the same2996// orientation.29972998inline void tetgenmesh::tspivot(triface& t, face& s)2999{3000if ((t).tet[9] == NULL) {3001(s).sh = NULL;3002return;3003}3004// Get the attached subface s.3005sdecode(((shellface *) (t).tet[9])[(t).ver & 3], (s));3006(s).shver = tspivottbl[t.ver][s.shver];3007}30083009// Quickly check if the handle (t, v) is a subface.3010#define issubface(t) \3011((t).tet[9] && ((t).tet[9])[(t).ver & 3])30123013// stpivot() finds a tetrahedron (t) abutting a given subface (s).3014// Return the t (if it exists) with the same edge and the same3015// orientation of s.30163017inline void tetgenmesh::stpivot(face& s, triface& t)3018{3019decode((tetrahedron) s.sh[9 + (s.shver & 1)], t);3020if ((t).tet == NULL) {3021return;3022}3023(t).ver = stpivottbl[t.ver][s.shver];3024}30253026// Quickly check if this subface is attached to a tetrahedron.30273028#define isshtet(s) \3029((s).sh[9 + ((s).shver & 1)])30303031// tsdissolve() dissolve a bond (from the tetrahedron side).30323033inline void tetgenmesh::tsdissolve(triface& t)3034{3035if ((t).tet[9] != NULL) {3036((shellface *) (t).tet[9])[(t).ver & 3] = NULL;3037}3038}30393040// stdissolve() dissolve a bond (from the subface side).30413042inline void tetgenmesh::stdissolve(face& s)3043{3044(s).sh[9] = NULL;3045(s).sh[10] = NULL;3046}30473048///////////////////////////////////////////////////////////////////////////////3049// //3050// Primitives for interacting between subfaces and segments //3051// //3052///////////////////////////////////////////////////////////////////////////////30533054// ssbond() bond a subface to a subsegment.30553056inline void tetgenmesh::ssbond(face& s, face& edge)3057{3058s.sh[6 + (s.shver >> 1)] = sencode(edge);3059edge.sh[0] = sencode(s);3060}30613062inline void tetgenmesh::ssbond1(face& s, face& edge)3063{3064s.sh[6 + (s.shver >> 1)] = sencode(edge);3065//edge.sh[0] = sencode(s);3066}30673068// ssdisolve() dissolve a bond (from the subface side)30693070inline void tetgenmesh::ssdissolve(face& s)3071{3072s.sh[6 + (s.shver >> 1)] = NULL;3073}30743075// sspivot() finds a subsegment abutting a subface.30763077inline void tetgenmesh::sspivot(face& s, face& edge)3078{3079sdecode((shellface) s.sh[6 + (s.shver >> 1)], edge);3080}30813082// Quickly check if the edge is a subsegment.30833084#define isshsubseg(s) \3085((s).sh[6 + ((s).shver >> 1)])30863087///////////////////////////////////////////////////////////////////////////////3088// //3089// Primitives for interacting between tetrahedra and segments //3090// //3091///////////////////////////////////////////////////////////////////////////////30923093inline void tetgenmesh::tssbond1(triface& t, face& s)3094{3095if ((t).tet[8] == NULL) {3096// Allocate space for this tet.3097(t).tet[8] = (tetrahedron) tet2segpool->alloc();3098// Initialization.3099for (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}31053106inline void tetgenmesh::sstbond1(face& s, triface& t)3107{3108((tetrahedron *) (s).sh)[9] = encode(t);3109}31103111inline void tetgenmesh::tssdissolve1(triface& t)3112{3113if ((t).tet[8] != NULL) {3114((shellface *) (t).tet[8])[ver2edge[(t).ver]] = NULL;3115}3116}31173118inline void tetgenmesh::sstdissolve1(face& s)3119{3120((tetrahedron *) (s).sh)[9] = NULL;3121}31223123inline void tetgenmesh::tsspivot1(triface& t, face& s)3124{3125if ((t).tet[8] != NULL) {3126sdecode(((shellface *) (t).tet[8])[ver2edge[(t).ver]], s);3127} else {3128(s).sh = NULL;3129}3130}31313132// Quickly check whether 't' is a segment or not.31333134#define issubseg(t) \3135((t).tet[8] && ((t).tet[8])[ver2edge[(t).ver]])31363137inline void tetgenmesh::sstpivot1(face& s, triface& t)3138{3139decode((tetrahedron) s.sh[9], t);3140}31413142///////////////////////////////////////////////////////////////////////////////3143// //3144// Primitives for points //3145// //3146///////////////////////////////////////////////////////////////////////////////31473148inline int tetgenmesh::pointmark(point pt) {3149return ((int *) (pt))[pointmarkindex];3150}31513152inline void tetgenmesh::setpointmark(point pt, int value) {3153((int *) (pt))[pointmarkindex] = value;3154}315531563157// These two primitives set and read the type of the point.31583159inline enum tetgenmesh::verttype tetgenmesh::pointtype(point pt) {3160return (enum verttype) (((int *) (pt))[pointmarkindex + 1] >> (int) 8);3161}31623163inline void tetgenmesh::setpointtype(point pt, enum verttype value) {3164((int *) (pt))[pointmarkindex + 1] =3165((int) value << 8) + (((int *) (pt))[pointmarkindex + 1] & (int) 255);3166}31673168// Read and set the geometry tag of the point (used by -s option).31693170inline int tetgenmesh::pointgeomtag(point pt) {3171return ((int *) (pt))[pointmarkindex + 2];3172}31733174inline void tetgenmesh::setpointgeomtag(point pt, int value) {3175((int *) (pt))[pointmarkindex + 2] = value;3176}31773178// Read and set the u,v coordinates of the point (used by -s option).31793180inline REAL tetgenmesh::pointgeomuv(point pt, int i) {3181return pt[pointparamindex + i];3182}31833184inline void tetgenmesh::setpointgeomuv(point pt, int i, REAL value) {3185pt[pointparamindex + i] = value;3186}31873188// pinfect(), puninfect(), pinfected() -- primitives to flag or unflag3189// a point. The last bit of the integer '[pointindex+1]' is flagged.31903191inline void tetgenmesh::pinfect(point pt) {3192((int *) (pt))[pointmarkindex + 1] |= (int) 1;3193}31943195inline void tetgenmesh::puninfect(point pt) {3196((int *) (pt))[pointmarkindex + 1] &= ~(int) 1;3197}31983199inline bool tetgenmesh::pinfected(point pt) {3200return (((int *) (pt))[pointmarkindex + 1] & (int) 1) != 0;3201}32023203// pmarktest(), punmarktest(), pmarktested() -- more primitives to3204// flag or unflag a point.32053206inline void tetgenmesh::pmarktest(point pt) {3207((int *) (pt))[pointmarkindex + 1] |= (int) 2;3208}32093210inline void tetgenmesh::punmarktest(point pt) {3211((int *) (pt))[pointmarkindex + 1] &= ~(int) 2;3212}32133214inline bool tetgenmesh::pmarktested(point pt) {3215return (((int *) (pt))[pointmarkindex + 1] & (int) 2) != 0;3216}32173218inline void tetgenmesh::pmarktest2(point pt) {3219((int *) (pt))[pointmarkindex + 1] |= (int) 4;3220}32213222inline void tetgenmesh::punmarktest2(point pt) {3223((int *) (pt))[pointmarkindex + 1] &= ~(int) 4;3224}32253226inline bool tetgenmesh::pmarktest2ed(point pt) {3227return (((int *) (pt))[pointmarkindex + 1] & (int) 4) != 0;3228}32293230inline void tetgenmesh::pmarktest3(point pt) {3231((int *) (pt))[pointmarkindex + 1] |= (int) 8;3232}32333234inline void tetgenmesh::punmarktest3(point pt) {3235((int *) (pt))[pointmarkindex + 1] &= ~(int) 8;3236}32373238inline bool tetgenmesh::pmarktest3ed(point pt) {3239return (((int *) (pt))[pointmarkindex + 1] & (int) 8) != 0;3240}32413242// These following primitives set and read a pointer to a tetrahedron3243// a subface/subsegment, a point, or a tet of background mesh.32443245inline tetgenmesh::tetrahedron tetgenmesh::point2tet(point pt) {3246return ((tetrahedron *) (pt))[point2simindex];3247}32483249inline void tetgenmesh::setpoint2tet(point pt, tetrahedron value) {3250((tetrahedron *) (pt))[point2simindex] = value;3251}32523253inline tetgenmesh::point tetgenmesh::point2ppt(point pt) {3254return (point) ((tetrahedron *) (pt))[point2simindex + 1];3255}32563257inline void tetgenmesh::setpoint2ppt(point pt, point value) {3258((tetrahedron *) (pt))[point2simindex + 1] = (tetrahedron) value;3259}32603261inline tetgenmesh::shellface tetgenmesh::point2sh(point pt) {3262return (shellface) ((tetrahedron *) (pt))[point2simindex + 2];3263}32643265inline void tetgenmesh::setpoint2sh(point pt, shellface value) {3266((tetrahedron *) (pt))[point2simindex + 2] = (tetrahedron) value;3267}326832693270inline tetgenmesh::tetrahedron tetgenmesh::point2bgmtet(point pt) {3271return ((tetrahedron *) (pt))[point2simindex + 3];3272}32733274inline void tetgenmesh::setpoint2bgmtet(point pt, tetrahedron value) {3275((tetrahedron *) (pt))[point2simindex + 3] = value;3276}327732783279// The primitives for saving and getting the insertion radius.3280inline void tetgenmesh::setpointinsradius(point pt, REAL value)3281{3282pt[pointinsradiusindex] = value;3283}32843285inline REAL tetgenmesh::getpointinsradius(point pt)3286{3287return pt[pointinsradiusindex];3288}32893290inline bool tetgenmesh::issteinerpoint(point pt) {3291return (pointtype(pt) == FREESEGVERTEX) || (pointtype(pt) == FREEFACETVERTEX)3292|| (pointtype(pt) == FREEVOLVERTEX);3293}32943295// point2tetorg() Get the tetrahedron whose origin is the point.32963297inline void tetgenmesh::point2tetorg(point pa, triface& searchtet)3298{3299decode(point2tet(pa), searchtet);3300if ((point) searchtet.tet[4] == pa) {3301searchtet.ver = 11;3302} else if ((point) searchtet.tet[5] == pa) {3303searchtet.ver = 3;3304} else if ((point) searchtet.tet[6] == pa) {3305searchtet.ver = 7;3306} else {3307searchtet.ver = 0;3308}3309}33103311// point2shorg() Get the subface/segment whose origin is the point.33123313inline void tetgenmesh::point2shorg(point pa, face& searchsh)3314{3315sdecode(point2sh(pa), searchsh);3316if ((point) searchsh.sh[3] == pa) {3317searchsh.shver = 0;3318} else if ((point) searchsh.sh[4] == pa) {3319searchsh.shver = (searchsh.sh[5] != NULL ? 2 : 1);3320} else {3321searchsh.shver = 4;3322}3323}33243325// farsorg() Return the origin of the subsegment.3326// farsdest() Return the destination of the subsegment.33273328inline tetgenmesh::point tetgenmesh::farsorg(face& s)3329{3330face travesh, neighsh;33313332travesh = s;3333while (1) {3334senext2(travesh, neighsh);3335spivotself(neighsh);3336if (neighsh.sh == NULL) break;3337if (sorg(neighsh) != sorg(travesh)) sesymself(neighsh);3338senext2(neighsh, travesh);3339}3340return sorg(travesh);3341}33423343inline tetgenmesh::point tetgenmesh::farsdest(face& s)3344{3345face travesh, neighsh;33463347travesh = s;3348while (1) {3349senext(travesh, neighsh);3350spivotself(neighsh);3351if (neighsh.sh == NULL) break;3352if (sdest(neighsh) != sdest(travesh)) sesymself(neighsh);3353senext(neighsh, travesh);3354}3355return sdest(travesh);3356}33573358///////////////////////////////////////////////////////////////////////////////3359// //3360// Linear algebra operators. //3361// //3362///////////////////////////////////////////////////////////////////////////////33633364// dot() returns the dot product: v1 dot v2.3365inline REAL tetgenmesh::dot(REAL* v1, REAL* v2)3366{3367return v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2];3368}33693370// cross() computes the cross product: n = v1 cross v2.3371inline void tetgenmesh::cross(REAL* v1, REAL* v2, REAL* n)3372{3373n[0] = v1[1] * v2[2] - v2[1] * v1[2];3374n[1] = -(v1[0] * v2[2] - v2[0] * v1[2]);3375n[2] = v1[0] * v2[1] - v2[0] * v1[1];3376}33773378// distance() computes the Euclidean distance between two points.3379inline REAL tetgenmesh::distance(REAL* p1, REAL* p2)3380{3381return 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}33853386inline REAL tetgenmesh::norm2(REAL x, REAL y, REAL z)3387{3388return (x) * (x) + (y) * (y) + (z) * (z);3389}339033913392#endif // #ifndef tetgenH3393339433953396