Path: blob/devel/ElmerGUI/Application/plugins/tetgen.h
3203 views
///////////////////////////////////////////////////////////////////////////////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 <stdio.h>44#include <stdlib.h>45#include <string.h>46#include <math.h>47#include <time.h>4849// The types 'intptr_t' and 'uintptr_t' are signed and unsigned integer types,50// respectively. They are guaranteed to be the same width as a pointer.51// They are defined in <stdint.h> by the C99 Standard. However, Microsoft52// Visual C++ 2003 -- 2008 (Visual C++ 7.1 - 9) doesn't ship with this header53// file. In such case, we can define them by ourself.54// Update (learned from Stack Overflow): Visual Studio 2010 and Visual C++ 201055// Express both have stdint.h5657// The following piece of code was provided by Steven Johnson (MIT). Define the58// symbol _MSC_VER if you are using Microsoft Visual C++. Moreover, define59// the _WIN64 symbol if you are running TetGen on Win64 systems.6061#ifdef _MSC_VER // Microsoft Visual C++62# ifdef _WIN6463typedef __int64 intptr_t;64typedef unsigned __int64 uintptr_t;65# else // not _WIN6466typedef int intptr_t;67typedef unsigned int uintptr_t;68# endif69#else // not Visual C++70# include <stdint.h>71#endif7273///////////////////////////////////////////////////////////////////////////////74// //75// tetgenio //76// //77// A structure for transferring data into and out of TetGen's mesh structure,//78// 'tetgenmesh' (declared below). //79// //80// The input of TetGen is either a 3D point set, or a 3D piecewise linear //81// complex (PLC), or a tetrahedral mesh. Depending on the input object and //82// the specified options, the output of TetGen is either a Delaunay (or wei- //83// ghted Delaunay) tetrahedralization, or a constrained (Delaunay) tetrahed- //84// ralization, or a quality tetrahedral mesh. //85// //86// A piecewise linear complex (PLC) represents a 3D polyhedral domain with //87// possibly internal boundaries(subdomains). It is introduced in [Miller et //88// al, 1996]. Basically it is a set of "cells", i.e., vertices, edges, poly- //89// gons, and polyhedra, and the intersection of any two of its cells is the //90// union of other cells of it. //91// //92// TetGen uses a set of files to describe the inputs and outputs. Each file //93// is identified from its file extension (.node, .ele, .face, .edge, etc). //94// //95// The 'tetgenio' structure is a collection of arrays of data, i.e., points, //96// facets, tetrahedra, and so forth. It contains functions to read and write //97// (input and output) files of TetGen as well as other supported mesh files. //98// //99// Once an object of tetgenio is declared, no array is created. One has to //100// allocate enough memory for them. On deletion of this object, the memory //101// occupied by these arrays needs to be freed. The routine deinitialize() //102// will be automatically called. It frees the memory for an array if it is //103// not a NULL. Note that it assumes that the memory is allocated by the C++ //104// "new" operator. Otherwise, the user is responsible to free them and all //105// pointers must be NULL before the call of the destructor. //106// //107///////////////////////////////////////////////////////////////////////////////108109class tetgenio {110111public:112113// A "polygon" describes a simple polygon (no holes). It is not necessarily114// convex. Each polygon contains a number of corners (points) and the same115// number of sides (edges). The points of the polygon must be given in116// either counterclockwise or clockwise order and they form a ring, so117// every two consecutive points forms an edge of the polygon.118typedef struct {119int *vertexlist;120int numberofvertices;121} polygon;122123// A "facet" describes a polygonal region possibly with holes, edges, and124// points floating in it. Each facet consists of a list of polygons and125// a list of hole points (which lie strictly inside holes).126typedef struct {127polygon *polygonlist;128int numberofpolygons;129REAL *holelist;130int numberofholes;131} facet;132133// A "voroedge" is an edge of the Voronoi diagram. It corresponds to a134// Delaunay face. Each voroedge is either a line segment connecting135// two Voronoi vertices or a ray starting from a Voronoi vertex to an136// "infinite vertex". 'v1' and 'v2' are two indices pointing to the137// list of Voronoi vertices. 'v1' must be non-negative, while 'v2' may138// be -1 if it is a ray, in this case, the unit normal of this ray is139// given in 'vnormal'.140typedef struct {141int v1, v2;142REAL vnormal[3];143} voroedge;144145// A "vorofacet" is an facet of the Voronoi diagram. It corresponds to a146// Delaunay edge. Each Voronoi facet is a convex polygon formed by a147// list of Voronoi edges, it may not be closed. 'c1' and 'c2' are two148// indices pointing into the list of Voronoi cells, i.e., the two cells149// share this facet. 'elist' is an array of indices pointing into the150// list of Voronoi edges, 'elist[0]' saves the number of Voronoi edges151// (including rays) of this facet.152typedef struct {153int c1, c2;154int *elist;155} vorofacet;156157158// Additional parameters associated with an input (or mesh) vertex.159// This information is provided by CAD libraries.160typedef struct {161REAL uv[2];162int tag;163int type; // 0, 1, or 2.164} pointparam;165166// Callback functions for meshing PSCs.167typedef REAL (* GetVertexParamOnEdge)(void*, int, int);168typedef void (* GetSteinerOnEdge)(void*, int, REAL, REAL*);169typedef void (* GetVertexParamOnFace)(void*, int, int, REAL*);170typedef void (* GetEdgeSteinerParamOnFace)(void*, int, REAL, int, REAL*);171typedef void (* GetSteinerOnFace)(void*, int, REAL*, REAL*);172173// A callback function for mesh refinement.174typedef bool (* TetSizeFunc)(REAL*, REAL*, REAL*, REAL*, REAL*, REAL);175176// Items are numbered starting from 'firstnumber' (0 or 1), default is 0.177int firstnumber;178179// Dimension of the mesh (2 or 3), default is 3.180int mesh_dim;181182// Does the lines in .node file contain index or not, default is 1.183int useindex;184185// 'pointlist': An array of point coordinates. The first point's x186// coordinate is at index [0] and its y coordinate at index [1], its187// z coordinate is at index [2], followed by the coordinates of the188// remaining points. Each point occupies three REALs.189// 'pointattributelist': An array of point attributes. Each point's190// attributes occupy 'numberofpointattributes' REALs.191// 'pointmtrlist': An array of metric tensors at points. Each point's192// tensor occupies 'numberofpointmtr' REALs.193// 'pointmarkerlist': An array of point markers; one integer per point.194// 'point2tetlist': An array of tetrahedra indices; one integer per point.195REAL *pointlist;196REAL *pointattributelist;197REAL *pointmtrlist;198int *pointmarkerlist;199int *point2tetlist;200pointparam *pointparamlist;201int numberofpoints;202int numberofpointattributes;203int numberofpointmtrs;204205// 'tetrahedronlist': An array of tetrahedron corners. The first206// tetrahedron's first corner is at index [0], followed by its other207// corners, followed by six nodes on the edges of the tetrahedron if the208// second order option (-o2) is applied. Each tetrahedron occupies209// 'numberofcorners' ints. The second order nodes are output only.210// 'tetrahedronattributelist': An array of tetrahedron attributes. Each211// tetrahedron's attributes occupy 'numberoftetrahedronattributes' REALs.212// 'tetrahedronvolumelist': An array of constraints, i.e. tetrahedron's213// volume; one REAL per element. Input only.214// 'neighborlist': An array of tetrahedron neighbors; 4 ints per element.215// 'tet2facelist': An array of tetrahedron face indices; 4 ints per element.216// 'tet2edgelist': An array of tetrahedron edge indices; 6 ints per element.217int *tetrahedronlist;218REAL *tetrahedronattributelist;219REAL *tetrahedronvolumelist;220int *neighborlist;221int *tet2facelist;222int *tet2edgelist;223int numberoftetrahedra;224int numberofcorners;225int numberoftetrahedronattributes;226227// 'facetlist': An array of facets. Each entry is a structure of facet.228// 'facetmarkerlist': An array of facet markers; one int per facet.229facet *facetlist;230int *facetmarkerlist;231int numberoffacets;232233// 'holelist': An array of holes (in volume). Each hole is given by a234// seed (point) which lies strictly inside it. The first seed's x, y and z235// coordinates are at indices [0], [1] and [2], followed by the236// remaining seeds. Three REALs per hole.237REAL *holelist;238int numberofholes;239240// 'regionlist': An array of regions (subdomains). Each region is given by241// a seed (point) which lies strictly inside it. The first seed's x, y and242// z coordinates are at indices [0], [1] and [2], followed by the regional243// attribute at index [3], followed by the maximum volume at index [4].244// Five REALs per region.245// Note that each regional attribute is used only if you select the 'A'246// switch, and each volume constraint is used only if you select the247// 'a' switch (with no number following).248REAL *regionlist;249int numberofregions;250251// 'facetconstraintlist': An array of facet constraints. Each constraint252// specifies a maximum area bound on the subfaces of that facet. The253// first facet constraint is given by a facet marker at index [0] and its254// maximum area bound at index [1], followed by the remaining facet con-255// straints. Two REALs per facet constraint. Note: the facet marker is256// actually an integer.257REAL *facetconstraintlist;258int numberoffacetconstraints;259260// 'segmentconstraintlist': An array of segment constraints. Each constraint261// specifies a maximum length bound on the subsegments of that segment.262// The first constraint is given by the two endpoints of the segment at263// index [0] and [1], and the maximum length bound at index [2], followed264// by the remaining segment constraints. Three REALs per constraint.265// Note the segment endpoints are actually integers.266REAL *segmentconstraintlist;267int numberofsegmentconstraints;268269270// 'trifacelist': An array of face (triangle) corners. The first face's271// three corners are at indices [0], [1] and [2], followed by the remaining272// faces. Three ints per face.273// 'trifacemarkerlist': An array of face markers; one int per face.274// 'o2facelist': An array of second order nodes (on the edges) of the face.275// It is output only if the second order option (-o2) is applied. The276// first face's three second order nodes are at [0], [1], and [2],277// followed by the remaining faces. Three ints per face.278// 'face2tetlist': An array of tetrahedra indices; 2 ints per face.279// 'face2edgelist': An array of edge indices; 3 ints per face.280int *trifacelist;281int *trifacemarkerlist;282int *o2facelist;283int *face2tetlist;284int *face2edgelist;285int numberoftrifaces;286287// 'edgelist': An array of edge endpoints. The first edge's endpoints288// are at indices [0] and [1], followed by the remaining edges.289// Two ints per edge.290// 'edgemarkerlist': An array of edge markers; one int per edge.291// 'o2edgelist': An array of midpoints of edges. It is output only if the292// second order option (-o2) is applied. One int per edge.293// 'edge2tetlist': An array of tetrahedra indices. One int per edge.294int *edgelist;295int *edgemarkerlist;296int *o2edgelist;297int *edge2tetlist;298int numberofedges;299300// 'vpointlist': An array of Voronoi vertex coordinates (like pointlist).301// 'vedgelist': An array of Voronoi edges. Each entry is a 'voroedge'.302// 'vfacetlist': An array of Voronoi facets. Each entry is a 'vorofacet'.303// 'vcelllist': An array of Voronoi cells. Each entry is an array of304// indices pointing into 'vfacetlist'. The 0th entry is used to store305// the length of this array.306REAL *vpointlist;307voroedge *vedgelist;308vorofacet *vfacetlist;309int **vcelllist;310int numberofvpoints;311int numberofvedges;312int numberofvfacets;313int numberofvcells;314315316// Variable (and callback functions) for meshing PSCs.317void *geomhandle;318GetVertexParamOnEdge getvertexparamonedge;319GetSteinerOnEdge getsteineronedge;320GetVertexParamOnFace getvertexparamonface;321GetEdgeSteinerParamOnFace getedgesteinerparamonface;322GetSteinerOnFace getsteineronface;323324// A callback function.325TetSizeFunc tetunsuitable;326327// Input & output routines.328virtual bool load_node_call(FILE* infile, int markers, int uvflag, char*);329virtual bool load_node(char*);330virtual bool load_edge(char*);331virtual bool load_face(char*);332virtual bool load_tet(char*);333virtual bool load_vol(char*);334virtual bool load_var(char*);335virtual bool load_mtr(char*);336// virtual bool load_pbc(char*);337virtual bool load_poly(char*);338virtual bool load_off(char*);339virtual bool load_ply(char*);340virtual bool load_stl(char*);341virtual bool load_vtk(char*);342virtual bool load_medit(char*, int);343virtual bool load_plc(char*, int);344virtual bool load_tetmesh(char*, int);345virtual void save_nodes(char*);346virtual void save_elements(char*);347virtual void save_faces(char*);348virtual void save_edges(char*);349virtual void save_neighbors(char*);350virtual void save_poly(char*);351virtual void save_faces2smesh(char*);352353// Read line and parse string functions.354virtual char *readline(char* string, FILE* infile, int *linenumber);355virtual char *findnextfield(char* string);356virtual char *readnumberline(char* string, FILE* infile, char* infilename);357virtual char *findnextnumber(char* string);358359virtual void init(polygon* p) {360p->vertexlist = (int *) NULL;361p->numberofvertices = 0;362}363364virtual void init(facet* f) {365f->polygonlist = (polygon *) NULL;366f->numberofpolygons = 0;367f->holelist = (REAL *) NULL;368f->numberofholes = 0;369}370371// Initialize routine.372virtual void initialize()373{374firstnumber = 0;375mesh_dim = 3;376useindex = 1;377378pointlist = (REAL *) NULL;379pointattributelist = (REAL *) NULL;380pointmtrlist = (REAL *) NULL;381pointmarkerlist = (int *) NULL;382point2tetlist = (int *) NULL;383pointparamlist = (pointparam *) NULL;384numberofpoints = 0;385numberofpointattributes = 0;386numberofpointmtrs = 0;387388tetrahedronlist = (int *) NULL;389tetrahedronattributelist = (REAL *) NULL;390tetrahedronvolumelist = (REAL *) NULL;391neighborlist = (int *) NULL;392tet2facelist = (int *) NULL;393tet2edgelist = (int *) NULL;394numberoftetrahedra = 0;395numberofcorners = 4;396numberoftetrahedronattributes = 0;397398trifacelist = (int *) NULL;399trifacemarkerlist = (int *) NULL;400o2facelist = (int *) NULL;401face2tetlist = (int *) NULL;402face2edgelist = (int *) NULL;403numberoftrifaces = 0;404405edgelist = (int *) NULL;406edgemarkerlist = (int *) NULL;407o2edgelist = (int *) NULL;408edge2tetlist = (int *) NULL;409numberofedges = 0;410411facetlist = (facet *) NULL;412facetmarkerlist = (int *) NULL;413numberoffacets = 0;414415holelist = (REAL *) NULL;416numberofholes = 0;417418regionlist = (REAL *) NULL;419numberofregions = 0;420421facetconstraintlist = (REAL *) NULL;422numberoffacetconstraints = 0;423segmentconstraintlist = (REAL *) NULL;424numberofsegmentconstraints = 0;425426427vpointlist = (REAL *) NULL;428vedgelist = (voroedge *) NULL;429vfacetlist = (vorofacet *) NULL;430vcelllist = (int **) NULL;431numberofvpoints = 0;432numberofvedges = 0;433numberofvfacets = 0;434numberofvcells = 0;435436437tetunsuitable = NULL;438439geomhandle = NULL;440getvertexparamonedge = NULL;441getsteineronedge = NULL;442getvertexparamonface = NULL;443getedgesteinerparamonface = NULL;444getsteineronface = NULL;445}446447// Free the memory allocated in 'tetgenio'. Note that it assumes that the448// memory was allocated by the "new" operator (C++).449virtual void deinitialize()450{451int i, j;452453if (pointlist != (REAL *) NULL) {454delete [] pointlist;455}456if (pointattributelist != (REAL *) NULL) {457delete [] pointattributelist;458}459if (pointmtrlist != (REAL *) NULL) {460delete [] pointmtrlist;461}462if (pointmarkerlist != (int *) NULL) {463delete [] pointmarkerlist;464}465if (point2tetlist != (int *) NULL) {466delete [] point2tetlist;467}468if (pointparamlist != (pointparam *) NULL) {469delete [] pointparamlist;470}471472if (tetrahedronlist != (int *) NULL) {473delete [] tetrahedronlist;474}475if (tetrahedronattributelist != (REAL *) NULL) {476delete [] tetrahedronattributelist;477}478if (tetrahedronvolumelist != (REAL *) NULL) {479delete [] tetrahedronvolumelist;480}481if (neighborlist != (int *) NULL) {482delete [] neighborlist;483}484if (tet2facelist != (int *) NULL) {485delete [] tet2facelist;486}487if (tet2edgelist != (int *) NULL) {488delete [] tet2edgelist;489}490491if (trifacelist != (int *) NULL) {492delete [] trifacelist;493}494if (trifacemarkerlist != (int *) NULL) {495delete [] trifacemarkerlist;496}497if (o2facelist != (int *) NULL) {498delete [] o2facelist;499}500if (face2tetlist != (int *) NULL) {501delete [] face2tetlist;502}503if (face2edgelist != (int *) NULL) {504delete [] face2edgelist;505}506507if (edgelist != (int *) NULL) {508delete [] edgelist;509}510if (edgemarkerlist != (int *) NULL) {511delete [] edgemarkerlist;512}513if (o2edgelist != (int *) NULL) {514delete [] o2edgelist;515}516if (edge2tetlist != (int *) NULL) {517delete [] edge2tetlist;518}519520if (facetlist != (facet *) NULL) {521facet *f;522polygon *p;523for (i = 0; i < numberoffacets; i++) {524f = &facetlist[i];525for (j = 0; j < f->numberofpolygons; j++) {526p = &f->polygonlist[j];527delete [] p->vertexlist;528}529delete [] f->polygonlist;530if (f->holelist != (REAL *) NULL) {531delete [] f->holelist;532}533}534delete [] facetlist;535}536if (facetmarkerlist != (int *) NULL) {537delete [] facetmarkerlist;538}539540if (holelist != (REAL *) NULL) {541delete [] holelist;542}543if (regionlist != (REAL *) NULL) {544delete [] regionlist;545}546if (facetconstraintlist != (REAL *) NULL) {547delete [] facetconstraintlist;548}549if (segmentconstraintlist != (REAL *) NULL) {550delete [] segmentconstraintlist;551}552if (vpointlist != (REAL *) NULL) {553delete [] vpointlist;554}555if (vedgelist != (voroedge *) NULL) {556delete [] vedgelist;557}558if (vfacetlist != (vorofacet *) NULL) {559for (i = 0; i < numberofvfacets; i++) {560delete [] vfacetlist[i].elist;561}562delete [] vfacetlist;563}564if (vcelllist != (int **) NULL) {565for (i = 0; i < numberofvcells; i++) {566delete [] vcelllist[i];567}568delete [] vcelllist;569}570}571572// Constructor & destructor.573tetgenio() {initialize();}574virtual ~tetgenio() {deinitialize();}575576}; // class tetgenio577578///////////////////////////////////////////////////////////////////////////////579// //580// tetgenbehavior //581// //582// A structure for maintaining the switches and parameters used by TetGen's //583// mesh data structure and algorithms. //584// //585// All switches and parameters are initialized with default values. They can //586// be set by the command line arguments (a list of strings) of TetGen. //587// //588// NOTE: Some of the switches are incompatible. While some may depend on //589// other switches. The routine parse_commandline() sets the switches from //590// the command line (a list of strings) and checks the consistency of the //591// applied switches. //592// //593///////////////////////////////////////////////////////////////////////////////594595class tetgenbehavior {596597public:598599// Switches of TetGen.600int plc; // '-p', 0.601int psc; // '-s', 0.602int refine; // '-r', 0.603int quality; // '-q', 0.604int nobisect; // '-Y', 0.605int coarsen; // '-R', 0.606int weighted; // '-w', 0.607int brio_hilbert; // '-b', 1.608int incrflip; // '-l', 0.609int flipinsert; // '-L', 0.610int metric; // '-m', 0.611int varvolume; // '-a', 0.612int fixedvolume; // '-a', 0.613int regionattrib; // '-A', 0.614int cdtrefine; // '-D', 0.615int insertaddpoints; // '-i', 0.616int diagnose; // '-d', 0.617int convex; // '-c', 0.618int nomergefacet; // '-M', 0.619int nomergevertex; // '-M', 0.620int noexact; // '-X', 0.621int nostaticfilter; // '-X', 0.622int zeroindex; // '-z', 0.623int facesout; // '-f', 0.624int edgesout; // '-e', 0.625int neighout; // '-n', 0.626int voroout; // '-v', 0.627int meditview; // '-g', 0.628int vtkview; // '-k', 0.629int nobound; // '-B', 0.630int nonodewritten; // '-N', 0.631int noelewritten; // '-E', 0.632int nofacewritten; // '-F', 0.633int noiterationnum; // '-I', 0.634int nojettison; // '-J', 0.635int docheck; // '-C', 0.636int quiet; // '-Q', 0.637int verbose; // '-V', 0.638639// Parameters of TetGen.640int vertexperblock; // '-x', 4092.641int tetrahedraperblock; // '-x', 8188.642int shellfaceperblock; // '-x', 2044.643int nobisect_nomerge; // '-Y', 1.644int supsteiner_level; // '-Y/', 2.645int addsteiner_algo; // '-Y//', 1.646int coarsen_param; // '-R', 0.647int weighted_param; // '-w', 0.648int fliplinklevel; // -1.649int flipstarsize; // -1.650int fliplinklevelinc; // 1.651int reflevel; // '-D', 3.652int optlevel; // '-O', 2.653int optscheme; // '-O', 7.654int delmaxfliplevel; // 1.655int order; // '-o', 1.656int reversetetori; // '-o/', 0.657int steinerleft; // '-S', 0.658int no_sort; // 0.659int hilbert_order; // '-b///', 52.660int hilbert_limit; // '-b//' 8.661int brio_threshold; // '-b' 64.662REAL brio_ratio; // '-b/' 0.125.663REAL facet_separate_ang_tol; // '-p', 179.9.664REAL facet_overlap_ang_tol; // '-p/', 0.1.665REAL facet_small_ang_tol; // '-p//', 15.0.666REAL maxvolume; // '-a', -1.0.667REAL minratio; // '-q', 0.0.668REAL mindihedral; // '-q', 5.0.669REAL optmaxdihedral; // 165.0.670REAL optminsmtdihed; // 179.0.671REAL optminslidihed; // 179.0.672REAL epsilon; // '-T', 1.0e-8.673REAL coarsen_percent; // -R1/#, 1.0.674675// Strings of command line arguments and input/output file names.676char commandline[1024];677char infilename[1024];678char outfilename[1024];679char addinfilename[1024];680char bgmeshfilename[1024];681682// The input object of TetGen. They are recognized by either the input683// file extensions or by the specified options.684// Currently the following objects are supported:685// - NODES, a list of nodes (.node);686// - POLY, a piecewise linear complex (.poly or .smesh);687// - OFF, a polyhedron (.off, Geomview's file format);688// - PLY, a polyhedron (.ply, file format from gatech, only ASCII);689// - STL, a surface mesh (.stl, stereolithography format);690// - MEDIT, a surface mesh (.mesh, Medit's file format);691// - MESH, a tetrahedral mesh (.ele).692// If no extension is available, the imposed command line switch693// (-p or -r) implies the object.694enum objecttype {NODES, POLY, OFF, PLY, STL, MEDIT, VTK, MESH} object;695696697virtual void syntax();698virtual void usage();699700// Command line parse routine.701virtual bool parse_commandline(int argc, char **argv);702virtual bool parse_commandline(char *switches) {703return parse_commandline(0, &switches);704}705706// Initialize all variables.707tetgenbehavior()708{709plc = 0;710psc = 0;711refine = 0;712quality = 0;713nobisect = 0;714coarsen = 0;715metric = 0;716weighted = 0;717brio_hilbert = 1;718incrflip = 0;719flipinsert = 0;720varvolume = 0;721fixedvolume = 0;722noexact = 0;723nostaticfilter = 0;724insertaddpoints = 0;725regionattrib = 0;726cdtrefine = 0;727diagnose = 0;728convex = 0;729zeroindex = 0;730facesout = 0;731edgesout = 0;732neighout = 0;733voroout = 0;734meditview = 0;735vtkview = 0;736nobound = 0;737nonodewritten = 0;738noelewritten = 0;739nofacewritten = 0;740noiterationnum = 0;741nomergefacet = 0;742nomergevertex = 0;743nojettison = 0;744docheck = 0;745quiet = 0;746verbose = 0;747748vertexperblock = 4092;749tetrahedraperblock = 8188;750shellfaceperblock = 4092;751nobisect_nomerge = 1;752supsteiner_level = 2;753addsteiner_algo = 1;754coarsen_param = 0;755weighted_param = 0;756fliplinklevel = -1;757flipstarsize = -1;758fliplinklevelinc = 1;759reflevel = 3;760optscheme = 7;761optlevel = 2;762delmaxfliplevel = 1;763order = 1;764reversetetori = 0;765steinerleft = -1;766no_sort = 0;767hilbert_order = 52; //-1;768hilbert_limit = 8;769brio_threshold = 64;770brio_ratio = 0.125;771facet_separate_ang_tol = 179.9;772facet_overlap_ang_tol = 0.1;773facet_small_ang_tol = 15.0;774maxvolume = -1.0;775minratio = 2.0;776mindihedral = 0.0;777optmaxdihedral = 165.00;778optminsmtdihed = 179.00;779optminslidihed = 179.00;780epsilon = 1.0e-8;781coarsen_percent = 1.0;782object = NODES;783784commandline[0] = '\0';785infilename[0] = '\0';786outfilename[0] = '\0';787addinfilename[0] = '\0';788bgmeshfilename[0] = '\0';789790}791792}; // class tetgenbehavior793794///////////////////////////////////////////////////////////////////////////////795// //796// Robust Geometric predicates //797// //798// Geometric predicates are simple tests of spatial relations of a set of d- //799// dimensional points, such as the orientation test and the point-in-sphere //800// test. Each of these tests is performed by evaluating the sign of a deter- //801// minant of a matrix whose entries are the coordinates of these points. If //802// the computation is performed by using the floating-point numbers, e.g., //803// the single or double precision numbers in C/C++, roundoff error may cause //804// an incorrect result. This may either lead to a wrong result or eventually //805// lead to a failure of the program. Computing the predicates exactly will //806// avoid the error and make the program robust. //807// //808// The following routines are the robust geometric predicates for 3D orient- //809// ation test and point-in-sphere test. They were implemented by Shewchuk. //810// The source code are generously provided by him in the public domain, //811// http://www.cs.cmu.edu/~quake/robust.html. predicates.cxx is a C++ version //812// of the original C code. //813// //814// The original predicates of Shewchuk only use "dynamic filters", i.e., it //815// computes the error at run time step by step. TetGen first adds a "static //816// filter" in each predicate. It estimates the maximal possible error in all //817// cases. So it can safely and quickly answer many easy cases. //818// //819///////////////////////////////////////////////////////////////////////////////820821void exactinit(int, int, int, REAL, REAL, REAL);822REAL orient3d(REAL *pa, REAL *pb, REAL *pc, REAL *pd);823REAL insphere(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe);824REAL orient4d(REAL *pa, REAL *pb, REAL *pc, REAL *pd, REAL *pe,825REAL ah, REAL bh, REAL ch, REAL dh, REAL eh);826827///////////////////////////////////////////////////////////////////////////////828// //829// tetgenmesh //830// //831// A structure for creating and updating tetrahedral meshes. //832// //833///////////////////////////////////////////////////////////////////////////////834835class tetgenmesh {836837public:838839///////////////////////////////////////////////////////////////////////////////840// //841// Mesh data structure //842// //843// A tetrahedral mesh T of a 3D piecewise linear complex (PLC) X is a 3D //844// simplicial complex whose underlying space is equal to the space of X. T //845// contains a 2D subcomplex S which is a triangular mesh of the boundary of //846// X. S contains a 1D subcomplex L which is a linear mesh of the boundary of //847// S. Faces and edges in S and L are respectively called subfaces and segme- //848// nts to distinguish them from others in T. //849// //850// TetGen stores the tetrahedra and vertices of T. The basic structure of a //851// tetrahedron contains pointers to its vertices and adjacent tetrahedra. A //852// vertex stores its x-, y-, and z-coordinates, and a pointer to a tetrahed- //853// ron containing it. Both tetrahedra and vertices may contain user data. //854// //855// Each face of T belongs to either two tetrahedra or one tetrahedron. In //856// the latter case, the face is an exterior boundary face of T. TetGen adds //857// fictitious tetrahedra (one-to-one) at such faces, and connects them to an //858// "infinite vertex" (which has no geometric coordinates). One can imagine //859// such a vertex lies in 4D space and is visible by all exterior boundary //860// faces. The extended set of tetrahedra (including the infinite vertex) is //861// a tetrahedralization of a 3-pseudomanifold without boundary. It has the //862// property that every face is shared by exactly two tetrahedra. //863// //864// The current version of TetGen stores explicitly the subfaces and segments //865// (which are in surface mesh S and the linear mesh L), respectively. Extra //866// pointers are allocated in tetrahedra and subfaces to point each others. //867// //868///////////////////////////////////////////////////////////////////////////////869870// The tetrahedron data structure. It includes the following fields:871// - a list of four adjoining tetrahedra;872// - a list of four vertices;873// - a pointer to a list of four subfaces (optional, for -p switch);874// - a pointer to a list of six segments (optional, for -p switch);875// - a list of user-defined floating-point attributes (optional);876// - a volume constraint (optional, for -a switch);877// - an integer of element marker (and flags);878// The structure of a tetrahedron is an array of pointers. Its actual size879// (the length of the array) is determined at runtime.880881typedef REAL **tetrahedron;882883// The subface data structure. It includes the following fields:884// - a list of three adjoining subfaces;885// - a list of three vertices;886// - a list of three adjoining segments;887// - two adjoining tetrahedra;888// - an area constraint (optional, for -q switch);889// - an integer for boundary marker;890// - an integer for type, flags, etc.891892typedef REAL **shellface;893894// The point data structure. It includes the following fields:895// - x, y and z coordinates;896// - a list of user-defined point attributes (optional);897// - u, v coordinates (optional, for -s switch);898// - a metric tensor (optional, for -q or -m switch);899// - a pointer to an adjacent tetrahedron;900// - a pointer to a parent (or a duplicate) point;901// - a pointer to an adjacent subface or segment (optional, -p switch);902// - a pointer to a tet in background mesh (optional, for -m switch);903// - an integer for boundary marker (point index);904// - an integer for point type (and flags).905// - an integer for geometry tag (optional, for -s switch).906// The structure of a point is an array of REALs. Its actual size is907// determined at the runtime.908909typedef REAL *point;910911///////////////////////////////////////////////////////////////////////////////912// //913// Handles //914// //915// Navigation and manipulation in a tetrahedralization are accomplished by //916// operating on structures referred as ``handles". A handle is a pair (t,v), //917// where t is a pointer to a tetrahedron, and v is a 4-bit integer, in the //918// range from 0 to 11. v is called the ``version'' of a tetrahedron, it rep- //919// resents a directed edge of a specific face of the tetrahedron. //920// //921// There are 12 even permutations of the four vertices, each of them corres- //922// ponds to a directed edge (a version) of the tetrahedron. The 12 versions //923// can be grouped into 4 distinct ``edge rings'' in 4 ``oriented faces'' of //924// this tetrahedron. One can encode each version (a directed edge) into a //925// 4-bit integer such that the two upper bits encode the index (from 0 to 2) //926// of this edge in the edge ring, and the two lower bits encode the index ( //927// from 0 to 3) of the oriented face which contains this edge. //928// //929// The four vertices of a tetrahedron are indexed from 0 to 3 (according to //930// their storage in the data structure). Give each face the same index as //931// the node opposite it in the tetrahedron. Denote the edge connecting face //932// i to face j as i/j. We number the twelve versions as follows: //933// //934// | edge 0 edge 1 edge 2 //935// --------|-------------------------------- //936// face 0 | 0 (0/1) 4 (0/3) 8 (0/2) //937// face 1 | 1 (1/2) 5 (1/3) 9 (1/0) //938// face 2 | 2 (2/3) 6 (2/1) 10 (2/0) //939// face 3 | 3 (3/0) 7 (3/1) 11 (3/2) //940// //941// Similarly, navigation and manipulation in a (boundary) triangulation are //942// done by using handles of triangles. Each handle is a pair (s, v), where s //943// is a pointer to a triangle, and v is a version in the range from 0 to 5. //944// Each version corresponds to a directed edge of this triangle. //945// //946// Number the three vertices of a triangle from 0 to 2 (according to their //947// storage in the data structure). Give each edge the same index as the node //948// opposite it in the triangle. The six versions of a triangle are: //949// //950// | edge 0 edge 1 edge 2 //951// ------------------|-------------------------- //952// ccw orientation | 0 2 4 //953// cw orientation | 1 3 5 //954// //955// In the following, a 'triface' is a handle of tetrahedron, and a 'face' is //956// a handle of a triangle. //957// //958///////////////////////////////////////////////////////////////////////////////959960class triface {961public:962tetrahedron *tet;963int ver; // Range from 0 to 11.964triface() : tet(0), ver(0) {}965triface& operator=(const triface& t) {966tet = t.tet; ver = t.ver;967return *this;968}969};970971class face {972public:973shellface *sh;974int shver; // Range from 0 to 5.975face() : sh(0), shver(0) {}976face& operator=(const face& s) {977sh = s.sh; shver = s.shver;978return *this;979}980};981982///////////////////////////////////////////////////////////////////////////////983// //984// Arraypool //985// //986// A dynamic linear array. (It is written by J. Shewchuk) //987// //988// Each arraypool contains an array of pointers to a number of blocks. Each //989// block contains the same fixed number of objects. Each index of the array //990// addresses a particular object in the pool. The most significant bits add- //991// ress the index of the block containing the object. The less significant //992// bits address this object within the block. //993// //994// 'objectbytes' is the size of one object in blocks; 'log2objectsperblock' //995// is the base-2 logarithm of 'objectsperblock'; 'objects' counts the number //996// of allocated objects; 'totalmemory' is the total memory in bytes. //997// //998///////////////////////////////////////////////////////////////////////////////9991000class arraypool {10011002public:10031004int objectbytes;1005int objectsperblock;1006int log2objectsperblock;1007int objectsperblockmark;1008int toparraylen;1009char **toparray;1010long objects;1011unsigned long totalmemory;10121013void restart();1014void poolinit(int sizeofobject, int log2objperblk);1015char* getblock(int objectindex);1016void* lookup(int objectindex);1017int newindex(void **newptr);10181019arraypool(int sizeofobject, int log2objperblk);1020~arraypool();1021};10221023// fastlookup() -- A fast, unsafe operation. Return the pointer to the object1024// with a given index. Note: The object's block must have been allocated,1025// i.e., by the function newindex().10261027#define fastlookup(pool, index) \1028(void *) ((pool)->toparray[(index) >> (pool)->log2objectsperblock] + \1029((index) & (pool)->objectsperblockmark) * (pool)->objectbytes)10301031///////////////////////////////////////////////////////////////////////////////1032// //1033// Memorypool //1034// //1035// A structure for memory allocation. (It is written by J. Shewchuk) //1036// //1037// firstblock is the first block of items. nowblock is the block from which //1038// items are currently being allocated. nextitem points to the next slab //1039// of free memory for an item. deaditemstack is the head of a linked list //1040// (stack) of deallocated items that can be recycled. unallocateditems is //1041// the number of items that remain to be allocated from nowblock. //1042// //1043// Traversal is the process of walking through the entire list of items, and //1044// is separate from allocation. Note that a traversal will visit items on //1045// the "deaditemstack" stack as well as live items. pathblock points to //1046// the block currently being traversed. pathitem points to the next item //1047// to be traversed. pathitemsleft is the number of items that remain to //1048// be traversed in pathblock. //1049// //1050///////////////////////////////////////////////////////////////////////////////10511052class memorypool {10531054public:10551056void **firstblock, **nowblock;1057void *nextitem;1058void *deaditemstack;1059void **pathblock;1060void *pathitem;1061int alignbytes;1062int itembytes, itemwords;1063int itemsperblock;1064long items, maxitems;1065int unallocateditems;1066int pathitemsleft;10671068memorypool();1069memorypool(int, int, int, int);1070~memorypool();10711072void poolinit(int, int, int, int);1073void restart();1074void *alloc();1075void dealloc(void*);1076void traversalinit();1077void *traverse();1078};10791080///////////////////////////////////////////////////////////////////////////////1081// //1082// badface //1083// //1084// Despite of its name, a 'badface' can be used to represent one of the //1085// following objects: //1086// - a face of a tetrahedron which is (possibly) non-Delaunay; //1087// - an encroached subsegment or subface; //1088// - a bad-quality tetrahedron, i.e, has too large radius-edge ratio; //1089// - a sliver, i.e., has good radius-edge ratio but nearly zero volume; //1090// - a recently flipped face (saved for undoing the flip later). //1091// //1092///////////////////////////////////////////////////////////////////////////////10931094class badface {1095public:1096triface tt;1097face ss;1098REAL key, cent[6]; // circumcenter or cos(dihedral angles) at 6 edges.1099point forg, fdest, fapex, foppo, noppo;1100badface *nextitem;1101badface() : key(0), forg(0), fdest(0), fapex(0), foppo(0), noppo(0),1102nextitem(0) {}1103};11041105///////////////////////////////////////////////////////////////////////////////1106// //1107// insertvertexflags //1108// //1109// A collection of flags that pass to the routine insertvertex(). //1110// //1111///////////////////////////////////////////////////////////////////////////////11121113class insertvertexflags {11141115public:11161117int iloc; // input/output.1118int bowywat, lawson;1119int splitbdflag, validflag, respectbdflag;1120int rejflag, chkencflag, cdtflag;1121int assignmeshsize;1122int sloc, sbowywat;11231124// Used by Delaunay refinement.1125int refineflag; // 0, 1, 2, 31126triface refinetet;1127face refinesh;1128int smlenflag; // for useinsertradius.1129REAL smlen; // for useinsertradius.1130point parentpt;11311132insertvertexflags() {1133iloc = bowywat = lawson = 0;1134splitbdflag = validflag = respectbdflag = 0;1135rejflag = chkencflag = cdtflag = 0;1136assignmeshsize = 0;1137sloc = sbowywat = 0;11381139refineflag = 0;1140refinetet.tet = NULL;1141refinesh.sh = NULL;1142smlenflag = 0;1143smlen = 0.0;1144}1145};11461147///////////////////////////////////////////////////////////////////////////////1148// //1149// flipconstraints //1150// //1151// A structure of a collection of data (options and parameters) which pass //1152// to the edge flip function flipnm(). //1153// //1154///////////////////////////////////////////////////////////////////////////////11551156class flipconstraints {11571158public:11591160// Elementary flip flags.1161int enqflag; // (= flipflag)1162int chkencflag;11631164// Control flags1165int unflip; // Undo the performed flips.1166int collectnewtets; // Collect the new tets created by flips.1167int collectencsegflag;11681169// Optimization flags.1170int remove_ndelaunay_edge; // Remove a non-Delaunay edge.1171REAL bak_tetprism_vol; // The value to be minimized.1172REAL tetprism_vol_sum;1173int remove_large_angle; // Remove a large dihedral angle at edge.1174REAL cosdihed_in; // The input cosine of the dihedral angle (> 0).1175REAL cosdihed_out; // The improved cosine of the dihedral angle.11761177// Boundary recovery flags.1178int checkflipeligibility;1179point seg[2]; // A constraining edge to be recovered.1180point fac[3]; // A constraining face to be recovered.1181point remvert; // A vertex to be removed.118211831184flipconstraints() {1185enqflag = 0;1186chkencflag = 0;11871188unflip = 0;1189collectnewtets = 0;1190collectencsegflag = 0;11911192remove_ndelaunay_edge = 0;1193bak_tetprism_vol = 0.0;1194tetprism_vol_sum = 0.0;1195remove_large_angle = 0;1196cosdihed_in = 0.0;1197cosdihed_out = 0.0;11981199checkflipeligibility = 0;1200seg[0] = NULL;1201fac[0] = NULL;1202remvert = NULL;1203}1204};12051206///////////////////////////////////////////////////////////////////////////////1207// //1208// optparameters //1209// //1210// Optimization options and parameters. //1211// //1212///////////////////////////////////////////////////////////////////////////////12131214class optparameters {12151216public:12171218// The one of goals of optimization.1219int max_min_volume; // Maximize the minimum volume.1220int min_max_aspectratio; // Minimize the maximum aspect ratio.1221int min_max_dihedangle; // Minimize the maximum dihedral angle.12221223// The initial and improved value.1224REAL initval, imprval;12251226int numofsearchdirs;1227REAL searchstep;1228int maxiter; // Maximum smoothing iterations (disabled by -1).1229int smthiter; // Performed iterations.123012311232optparameters() {1233max_min_volume = 0;1234min_max_aspectratio = 0;1235min_max_dihedangle = 0;12361237initval = imprval = 0.0;12381239numofsearchdirs = 10;1240searchstep = 0.01;1241maxiter = -1; // Unlimited smoothing iterations.1242smthiter = 0;12431244}1245};124612471248///////////////////////////////////////////////////////////////////////////////1249// //1250// Labels (enumeration declarations) used by TetGen. //1251// //1252///////////////////////////////////////////////////////////////////////////////12531254// Labels that signify the type of a vertex.1255enum verttype {UNUSEDVERTEX, DUPLICATEDVERTEX, RIDGEVERTEX, ACUTEVERTEX,1256FACETVERTEX, VOLVERTEX, FREESEGVERTEX, FREEFACETVERTEX,1257FREEVOLVERTEX, NREGULARVERTEX, DEADVERTEX};12581259// Labels that signify the result of triangle-triangle intersection test.1260enum interresult {DISJOINT, INTERSECT, SHAREVERT, SHAREEDGE, SHAREFACE,1261TOUCHEDGE, TOUCHFACE, ACROSSVERT, ACROSSEDGE, ACROSSFACE};12621263// Labels that signify the result of point location.1264enum locateresult {TGUNKNOWN, OUTSIDE, INTETRAHEDRON, ONFACE, ONEDGE, ONVERTEX,1265ENCVERTEX, ENCSEGMENT, ENCSUBFACE, NEARVERTEX, NONREGULAR,1266INSTAR, BADELEMENT};12671268///////////////////////////////////////////////////////////////////////////////1269// //1270// Variables of TetGen //1271// //1272///////////////////////////////////////////////////////////////////////////////12731274// Pointer to the input data (a set of nodes, a PLC, or a mesh).1275tetgenio *in, *addin;12761277// Pointer to the switches and parameters.1278tetgenbehavior *b;12791280// Pointer to a background mesh (contains size specification map).1281tetgenmesh *bgm;12821283// Memorypools to store mesh elements (points, tetrahedra, subfaces, and1284// segments) and extra pointers between tetrahedra, subfaces, and segments.1285memorypool *tetrahedrons, *subfaces, *subsegs, *points;1286memorypool *tet2subpool, *tet2segpool;12871288// Memorypools to store bad-quality (or encroached) elements.1289memorypool *badtetrahedrons, *badsubfacs, *badsubsegs;12901291// A memorypool to store faces to be flipped.1292memorypool *flippool;1293arraypool *unflipqueue;1294badface *flipstack;12951296// Arrays used for point insertion (the Bowyer-Watson algorithm).1297arraypool *cavetetlist, *cavebdrylist, *caveoldtetlist;1298arraypool *cavetetshlist, *cavetetseglist, *cavetetvertlist;1299arraypool *caveencshlist, *caveencseglist;1300arraypool *caveshlist, *caveshbdlist, *cavesegshlist;13011302// Stacks used for CDT construction and boundary recovery.1303arraypool *subsegstack, *subfacstack, *subvertstack;13041305// Arrays of encroached segments and subfaces (for mesh refinement).1306arraypool *encseglist, *encshlist;13071308// The map between facets to their vertices (for mesh refinement).1309int *idx2facetlist;1310point *facetverticeslist;13111312// The map between segments to their endpoints (for mesh refinement).1313point *segmentendpointslist;13141315// The infinite vertex.1316point dummypoint;1317// The recently visited tetrahedron, subface.1318triface recenttet;1319face recentsh;13201321// PI is the ratio of a circle's circumference to its diameter.1322static REAL PI;13231324// Array (size = numberoftetrahedra * 6) for storing high-order nodes of1325// tetrahedra (only used when -o2 switch is selected).1326point *highordertable;13271328// Various variables.1329int numpointattrib; // Number of point attributes.1330int numelemattrib; // Number of tetrahedron attributes.1331int sizeoftensor; // Number of REALs per metric tensor.1332int pointmtrindex; // Index to find the metric tensor of a point.1333int pointparamindex; // Index to find the u,v coordinates of a point.1334int point2simindex; // Index to find a simplex adjacent to a point.1335int pointmarkindex; // Index to find boundary marker of a point.1336int pointinsradiusindex; // Index to find the insertion radius of a point.1337int elemattribindex; // Index to find attributes of a tetrahedron.1338int volumeboundindex; // Index to find volume bound of a tetrahedron.1339int elemmarkerindex; // Index to find marker of a tetrahedron.1340int shmarkindex; // Index to find boundary marker of a subface.1341int areaboundindex; // Index to find area bound of a subface.1342int checksubsegflag; // Are there segments in the tetrahedralization yet?1343int checksubfaceflag; // Are there subfaces in the tetrahedralization yet?1344int checkconstraints; // Are there variant (node, seg, facet) constraints?1345int nonconvex; // Is current mesh non-convex?1346int autofliplinklevel; // The increase of link levels, default is 1.1347int useinsertradius; // Save the insertion radius for Steiner points.1348long samples; // Number of random samples for point location.1349unsigned long randomseed; // Current random number seed.1350REAL cosmaxdihed, cosmindihed; // The cosine values of max/min dihedral.1351REAL cossmtdihed; // The cosine value of a bad dihedral to be smoothed.1352REAL cosslidihed; // The cosine value of the max dihedral of a sliver.1353REAL minfaceang, minfacetdihed; // The minimum input (dihedral) angles.1354REAL tetprism_vol_sum; // The total volume of tetrahedral-prisms (in 4D).1355REAL longest; // The longest possible edge length.1356REAL minedgelength; // = longest * b->epsion.1357REAL xmax, xmin, ymax, ymin, zmax, zmin; // Bounding box of points.13581359// Counters.1360long insegments; // Number of input segments.1361long hullsize; // Number of exterior boundary faces.1362long meshedges; // Number of mesh edges.1363long meshhulledges; // Number of boundary mesh edges.1364long steinerleft; // Number of Steiner points not yet used.1365long dupverts; // Are there duplicated vertices?1366long unuverts; // Are there unused vertices?1367long nonregularcount; // Are there non-regular vertices?1368long st_segref_count, st_facref_count, st_volref_count; // Steiner points.1369long fillregioncount, cavitycount, cavityexpcount;1370long flip14count, flip26count, flipn2ncount;1371long flip23count, flip32count, flip44count, flip41count;1372long flip31count, flip22count;1373unsigned long totalworkmemory; // Total memory used by working arrays.137413751376///////////////////////////////////////////////////////////////////////////////1377// //1378// Mesh manipulation primitives //1379// //1380///////////////////////////////////////////////////////////////////////////////13811382// Fast lookup tables for mesh manipulation primitives.1383static int bondtbl[12][12], fsymtbl[12][12];1384static int esymtbl[12], enexttbl[12], eprevtbl[12];1385static int enextesymtbl[12], eprevesymtbl[12];1386static int eorgoppotbl[12], edestoppotbl[12];1387static int facepivot1[12], facepivot2[12][12];1388static int orgpivot[12], destpivot[12], apexpivot[12], oppopivot[12];1389static int tsbondtbl[12][6], stbondtbl[12][6];1390static int tspivottbl[12][6], stpivottbl[12][6];1391static int ver2edge[12], edge2ver[6], epivot[12];1392static int sorgpivot [6], sdestpivot[6], sapexpivot[6];1393static int snextpivot[6];13941395void inittables();13961397// Primitives for tetrahedra.1398inline tetrahedron encode(triface& t);1399inline tetrahedron encode2(tetrahedron* ptr, int ver);1400inline void decode(tetrahedron ptr, triface& t);1401inline void bond(triface& t1, triface& t2);1402inline void dissolve(triface& t);1403inline void esym(triface& t1, triface& t2);1404inline void esymself(triface& t);1405inline void enext(triface& t1, triface& t2);1406inline void enextself(triface& t);1407inline void eprev(triface& t1, triface& t2);1408inline void eprevself(triface& t);1409inline void enextesym(triface& t1, triface& t2);1410inline void enextesymself(triface& t);1411inline void eprevesym(triface& t1, triface& t2);1412inline void eprevesymself(triface& t);1413inline void eorgoppo(triface& t1, triface& t2);1414inline void eorgoppoself(triface& t);1415inline void edestoppo(triface& t1, triface& t2);1416inline void edestoppoself(triface& t);1417inline void fsym(triface& t1, triface& t2);1418inline void fsymself(triface& t);1419inline void fnext(triface& t1, triface& t2);1420inline void fnextself(triface& t);1421inline point org (triface& t);1422inline point dest(triface& t);1423inline point apex(triface& t);1424inline point oppo(triface& t);1425inline void setorg (triface& t, point p);1426inline void setdest(triface& t, point p);1427inline void setapex(triface& t, point p);1428inline void setoppo(triface& t, point p);1429inline REAL elemattribute(tetrahedron* ptr, int attnum);1430inline void setelemattribute(tetrahedron* ptr, int attnum, REAL value);1431inline REAL volumebound(tetrahedron* ptr);1432inline void setvolumebound(tetrahedron* ptr, REAL value);1433inline int elemindex(tetrahedron* ptr);1434inline void setelemindex(tetrahedron* ptr, int value);1435inline int elemmarker(tetrahedron* ptr);1436inline void setelemmarker(tetrahedron* ptr, int value);1437inline void infect(triface& t);1438inline void uninfect(triface& t);1439inline bool infected(triface& t);1440inline void marktest(triface& t);1441inline void unmarktest(triface& t);1442inline bool marktested(triface& t);1443inline void markface(triface& t);1444inline void unmarkface(triface& t);1445inline bool facemarked(triface& t);1446inline void markedge(triface& t);1447inline void unmarkedge(triface& t);1448inline bool edgemarked(triface& t);1449inline void marktest2(triface& t);1450inline void unmarktest2(triface& t);1451inline bool marktest2ed(triface& t);1452inline int elemcounter(triface& t);1453inline void setelemcounter(triface& t, int value);1454inline void increaseelemcounter(triface& t);1455inline void decreaseelemcounter(triface& t);1456inline bool ishulltet(triface& t);1457inline bool isdeadtet(triface& t);14581459// Primitives for subfaces and subsegments.1460inline void sdecode(shellface sptr, face& s);1461inline shellface sencode(face& s);1462inline shellface sencode2(shellface *sh, int shver);1463inline void spivot(face& s1, face& s2);1464inline void spivotself(face& s);1465inline void sbond(face& s1, face& s2);1466inline void sbond1(face& s1, face& s2);1467inline void sdissolve(face& s);1468inline point sorg(face& s);1469inline point sdest(face& s);1470inline point sapex(face& s);1471inline void setsorg(face& s, point pointptr);1472inline void setsdest(face& s, point pointptr);1473inline void setsapex(face& s, point pointptr);1474inline void sesym(face& s1, face& s2);1475inline void sesymself(face& s);1476inline void senext(face& s1, face& s2);1477inline void senextself(face& s);1478inline void senext2(face& s1, face& s2);1479inline void senext2self(face& s);1480inline REAL areabound(face& s);1481inline void setareabound(face& s, REAL value);1482inline int shellmark(face& s);1483inline void setshellmark(face& s, int value);1484inline void sinfect(face& s);1485inline void suninfect(face& s);1486inline bool sinfected(face& s);1487inline void smarktest(face& s);1488inline void sunmarktest(face& s);1489inline bool smarktested(face& s);1490inline void smarktest2(face& s);1491inline void sunmarktest2(face& s);1492inline bool smarktest2ed(face& s);1493inline void smarktest3(face& s);1494inline void sunmarktest3(face& s);1495inline bool smarktest3ed(face& s);1496inline void setfacetindex(face& f, int value);1497inline int getfacetindex(face& f);14981499// Primitives for interacting tetrahedra and subfaces.1500inline void tsbond(triface& t, face& s);1501inline void tsdissolve(triface& t);1502inline void stdissolve(face& s);1503inline void tspivot(triface& t, face& s);1504inline void stpivot(face& s, triface& t);15051506// Primitives for interacting tetrahedra and segments.1507inline void tssbond1(triface& t, face& seg);1508inline void sstbond1(face& s, triface& t);1509inline void tssdissolve1(triface& t);1510inline void sstdissolve1(face& s);1511inline void tsspivot1(triface& t, face& s);1512inline void sstpivot1(face& s, triface& t);15131514// Primitives for interacting subfaces and segments.1515inline void ssbond(face& s, face& edge);1516inline void ssbond1(face& s, face& edge);1517inline void ssdissolve(face& s);1518inline void sspivot(face& s, face& edge);15191520// Primitives for points.1521inline int pointmark(point pt);1522inline void setpointmark(point pt, int value);1523inline enum verttype pointtype(point pt);1524inline void setpointtype(point pt, enum verttype value);1525inline int pointgeomtag(point pt);1526inline void setpointgeomtag(point pt, int value);1527inline REAL pointgeomuv(point pt, int i);1528inline void setpointgeomuv(point pt, int i, REAL value);1529inline void pinfect(point pt);1530inline void puninfect(point pt);1531inline bool pinfected(point pt);1532inline void pmarktest(point pt);1533inline void punmarktest(point pt);1534inline bool pmarktested(point pt);1535inline void pmarktest2(point pt);1536inline void punmarktest2(point pt);1537inline bool pmarktest2ed(point pt);1538inline void pmarktest3(point pt);1539inline void punmarktest3(point pt);1540inline bool pmarktest3ed(point pt);1541inline tetrahedron point2tet(point pt);1542inline void setpoint2tet(point pt, tetrahedron value);1543inline shellface point2sh(point pt);1544inline void setpoint2sh(point pt, shellface value);1545inline point point2ppt(point pt);1546inline void setpoint2ppt(point pt, point value);1547inline tetrahedron point2bgmtet(point pt);1548inline void setpoint2bgmtet(point pt, tetrahedron value);1549inline void setpointinsradius(point pt, REAL value);1550inline REAL getpointinsradius(point pt);1551inline bool issteinerpoint(point pt);15521553// Advanced primitives.1554inline void point2tetorg(point pt, triface& t);1555inline void point2shorg(point pa, face& s);1556inline point farsorg(face& seg);1557inline point farsdest(face& seg);15581559///////////////////////////////////////////////////////////////////////////////1560// //1561// Memory management //1562// //1563///////////////////////////////////////////////////////////////////////////////15641565void tetrahedrondealloc(tetrahedron*);1566tetrahedron *tetrahedrontraverse();1567tetrahedron *alltetrahedrontraverse();1568void shellfacedealloc(memorypool*, shellface*);1569shellface *shellfacetraverse(memorypool*);1570void pointdealloc(point);1571point pointtraverse();15721573void makeindex2pointmap(point*&);1574void makepoint2submap(memorypool*, int*&, face*&);1575void maketetrahedron(triface*);1576void makeshellface(memorypool*, face*);1577void makepoint(point*, enum verttype);15781579void initializepools();15801581///////////////////////////////////////////////////////////////////////////////1582// //1583// Advanced geometric predicates and calculations //1584// //1585// TetGen uses a simplified symbolic perturbation scheme from Edelsbrunner, //1586// et al [*]. Hence the point-in-sphere test never returns a zero. The idea //1587// is to perturb the weights of vertices in the fourth dimension. TetGen //1588// uses the indices of the vertices decide the amount of perturbation. It is //1589// implemented in the routine insphere_s().1590// //1591// The routine tri_edge_test() determines whether or not a triangle and an //1592// edge intersect in 3D. If they intersect, their intersection type is also //1593// reported. This test is a combination of n 3D orientation tests (n is bet- //1594// ween 3 and 9). It uses the robust orient3d() test to make the branch dec- //1595// isions. The routine tri_tri_test() determines whether or not two triang- //1596// les intersect in 3D. It also uses the robust orient3d() test. //1597// //1598// There are a number of routines to calculate geometrical quantities, e.g., //1599// circumcenters, angles, dihedral angles, face normals, face areas, etc. //1600// They are so far done by the default floating-point arithmetics which are //1601// non-robust. They should be improved in the future. //1602// //1603///////////////////////////////////////////////////////////////////////////////16041605// Symbolic perturbations (robust)1606REAL insphere_s(REAL*, REAL*, REAL*, REAL*, REAL*);1607REAL orient4d_s(REAL*, REAL*, REAL*, REAL*, REAL*,1608REAL, REAL, REAL, REAL, REAL);16091610// Triangle-edge intersection test (robust)1611int tri_edge_2d(point, point, point, point, point, point, int, int*, int*);1612int tri_edge_tail(point, point, point, point, point, point, REAL, REAL, int,1613int*, int*);1614int tri_edge_test(point, point, point, point, point, point, int, int*, int*);16151616// Triangle-triangle intersection test (robust)1617int tri_edge_inter_tail(point, point, point, point, point, REAL, REAL);1618int tri_tri_inter(point, point, point, point, point, point);16191620// Linear algebra functions1621inline REAL dot(REAL* v1, REAL* v2);1622inline void cross(REAL* v1, REAL* v2, REAL* n);1623bool lu_decmp(REAL lu[4][4], int n, int* ps, REAL* d, int N);1624void lu_solve(REAL lu[4][4], int n, int* ps, REAL* b, int N);16251626// An embedded 2-dimensional geometric predicate (non-robust)1627REAL incircle3d(point pa, point pb, point pc, point pd);16281629// Geometric calculations (non-robust)1630REAL orient3dfast(REAL *pa, REAL *pb, REAL *pc, REAL *pd);1631inline REAL norm2(REAL x, REAL y, REAL z);1632inline REAL distance(REAL* p1, REAL* p2);1633void facenormal(point pa, point pb, point pc, REAL *n, int pivot, REAL *lav);1634REAL shortdistance(REAL* p, REAL* e1, REAL* e2);1635REAL triarea(REAL* pa, REAL* pb, REAL* pc);1636REAL interiorangle(REAL* o, REAL* p1, REAL* p2, REAL* n);1637void projpt2edge(REAL* p, REAL* e1, REAL* e2, REAL* prj);1638void projpt2face(REAL* p, REAL* f1, REAL* f2, REAL* f3, REAL* prj);1639bool tetalldihedral(point, point, point, point, REAL*, REAL*, REAL*);1640void tetallnormal(point, point, point, point, REAL N[4][3], REAL* volume);1641REAL tetaspectratio(point, point, point, point);1642bool circumsphere(REAL*, REAL*, REAL*, REAL*, REAL* cent, REAL* radius);1643bool orthosphere(REAL*,REAL*,REAL*,REAL*,REAL,REAL,REAL,REAL,REAL*,REAL*);1644void planelineint(REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*);1645int linelineint(REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*, REAL*);1646REAL tetprismvol(REAL* pa, REAL* pb, REAL* pc, REAL* pd);1647bool calculateabovepoint(arraypool*, point*, point*, point*);1648void calculateabovepoint4(point, point, point, point);16491650// PLC error reports.1651void report_overlapping_facets(face*, face*, REAL dihedang = 0.0);1652int report_selfint_edge(point, point, face* sedge, triface* searchtet,1653enum interresult);1654int report_selfint_face(point, point, point, face* sface, triface* iedge,1655int intflag, int* types, int* poss);16561657///////////////////////////////////////////////////////////////////////////////1658// //1659// Local mesh transformations //1660// //1661// A local transformation replaces a small set of tetrahedra with another //1662// set of tetrahedra which fills the same space and the same boundaries. //1663// In 3D, the most simplest local transformations are the elementary flips //1664// performed within the convex hull of five vertices: 2-to-3, 3-to-2, 1-to-4,//1665// and 4-to-1 flips, where the numbers indicate the number of tetrahedra //1666// before and after each flip. The 1-to-4 and 4-to-1 flip involve inserting //1667// or deleting a vertex, respectively. //1668// There are complex local transformations which can be decomposed as a //1669// combination of elementary flips. For example,a 4-to-4 flip which replaces //1670// two coplanar edges can be regarded by a 2-to-3 flip and a 3-to-2 flip. //1671// Note that the first 2-to-3 flip will temporarily create a degenerate tet- //1672// rahedron which is removed immediately by the followed 3-to-2 flip. More //1673// generally, a n-to-m flip, where n > 3, m = (n - 2) * 2, which removes an //1674// edge can be done by first performing a sequence of (n - 3) 2-to-3 flips //1675// followed by a 3-to-2 flip. //1676// //1677// The routines flip23(), flip32(), and flip41() perform the three element- //1678// ray flips. The flip14() is available inside the routine insertpoint(). //1679// //1680// The routines flipnm() and flipnm_post() implement a generalized edge flip //1681// algorithm which uses a combination of elementary flips. //1682// //1683// The routine insertpoint() implements a variant of Bowyer-Watson's cavity //1684// algorithm to insert a vertex. It works for arbitrary tetrahedralization, //1685// either Delaunay, or constrained Delaunay, or non-Delaunay. //1686// //1687///////////////////////////////////////////////////////////////////////////////16881689// The elementary flips.1690void flip23(triface*, int, flipconstraints* fc);1691void flip32(triface*, int, flipconstraints* fc);1692void flip41(triface*, int, flipconstraints* fc);16931694// A generalized edge flip.1695int flipnm(triface*, int n, int level, int, flipconstraints* fc);1696int flipnm_post(triface*, int n, int nn, int, flipconstraints* fc);16971698// Point insertion.1699int insertpoint(point, triface*, face*, face*, insertvertexflags*);1700void insertpoint_abort(face*, insertvertexflags*);17011702///////////////////////////////////////////////////////////////////////////////1703// //1704// Delaunay tetrahedralization //1705// //1706// The routine incrementaldelaunay() implemented two incremental algorithms //1707// for constructing Delaunay tetrahedralizations (DTs): the Bowyer-Watson //1708// (B-W) algorithm and the incremental flip algorithm of Edelsbrunner and //1709// Shah, "Incremental topological flipping works for regular triangulation," //1710// Algorithmica, 15:233-241, 1996. //1711// //1712// The routine incrementalflip() implements the flip algorithm of [Edelsbru- //1713// nner and Shah, 1996]. It flips a queue of locally non-Delaunay faces (in //1714// an arbitrary order). The success is guaranteed when the Delaunay tetrah- //1715// edralization is constructed incrementally by adding one vertex at a time. //1716// //1717// The routine locate() finds a tetrahedron contains a new point in current //1718// DT. It uses a simple stochastic walk algorithm: starting from an arbitr- //1719// ary tetrahedron in DT, it finds the destination by visit one tetrahedron //1720// at a time, randomly chooses a tetrahedron if there are more than one //1721// choices. This algorithm terminates due to Edelsbrunner's acyclic theorem. //1722// Choose a good starting tetrahedron is crucial to the speed of the walk. //1723// TetGen originally uses the "jump-and-walk" algorithm of Muecke, E.P., //1724// Saias, I., and Zhu, B. "Fast Randomized Point Location Without Preproces- //1725// sing." In Proceedings of the 12th ACM Symposium on Computational Geometry,//1726// 274-283, 1996. It first randomly samples several tetrahedra in the DT //1727// and then choosing the closet one to start walking. //1728// The above algorithm slows download dramatically as the number of points //1729// grows -- reported in Amenta, N., Choi, S. and Rote, G., "Incremental //1730// construction con {BRIO}," In Proceedings of 19th ACM Symposium on //1731// Computational Geometry, 211-219, 2003. On the other hand, Liu and //1732// Snoeyink showed that the point location can be made in constant time if //1733// the points are pre-sorted so that the nearby points in space have nearby //1734// indices, then adding the points in this order. They sorted the points //1735// along the 3D Hilbert curve. //1736// //1737// The routine hilbert_sort3() sorts a set of 3D points along the 3D Hilbert //1738// curve. It recursively splits a point set according to the Hilbert indices //1739// mapped to the subboxes of the bounding box of the point set. //1740// The Hilbert indices is calculated by Butz's algorithm in 1971. A nice //1741// exposition of this algorithm can be found in the paper of Hamilton, C., //1742// "Compact Hilbert Indices", Technical Report CS-2006-07, Computer Science, //1743// Dalhousie University, 2006 (the Section 2). My implementation also refer- //1744// enced Steven Witham's implementation of "Hilbert walk" (hopefully, it is //1745// still available at: http://www.tiac.net/~sw/2008/10/Hilbert/). //1746// //1747// TetGen sorts the points using the method in the paper of Boissonnat,J.-D.,//1748// Devillers, O. and Hornus, S. "Incremental Construction of the Delaunay //1749// Triangulation and the Delaunay Graph in Medium Dimension," In Proceedings //1750// of the 25th ACM Symposium on Computational Geometry, 2009. //1751// It first randomly sorts the points into subgroups using the Biased Rand-//1752// omized Insertion Ordering (BRIO) of Amenta et al 2003, then sorts the //1753// points in each subgroup along the 3D Hilbert curve. Inserting points in //1754// this order ensures a randomized "sprinkling" of the points over the //1755// domain, while sorting of each subset ensures locality. //1756// //1757///////////////////////////////////////////////////////////////////////////////17581759void transfernodes();17601761// Point sorting.1762int transgc[8][3][8], tsb1mod3[8];1763void hilbert_init(int n);1764int hilbert_split(point* vertexarray, int arraysize, int gc0, int gc1,1765REAL, REAL, REAL, REAL, REAL, REAL);1766void hilbert_sort3(point* vertexarray, int arraysize, int e, int d,1767REAL, REAL, REAL, REAL, REAL, REAL, int depth);1768void brio_multiscale_sort(point*,int,int threshold,REAL ratio,int* depth);17691770// Point location.1771unsigned long randomnation(unsigned int choices);1772void randomsample(point searchpt, triface *searchtet);1773enum locateresult locate(point searchpt, triface *searchtet,1774int chkencflag = 0);17751776// Incremental flips.1777void flippush(badface*&, triface*);1778int incrementalflip(point newpt, int, flipconstraints *fc);17791780// Incremental Delaunay construction.1781void initialdelaunay(point pa, point pb, point pc, point pd);1782void incrementaldelaunay(clock_t&);17831784///////////////////////////////////////////////////////////////////////////////1785// //1786// Surface triangulation //1787// //1788///////////////////////////////////////////////////////////////////////////////17891790void flipshpush(face*);1791void flip22(face*, int, int);1792void flip31(face*, int);1793long lawsonflip();1794int sinsertvertex(point newpt, face*, face*, int iloc, int bowywat, int);1795int sremovevertex(point delpt, face*, face*, int lawson);17961797enum locateresult slocate(point, face*, int, int, int);1798enum interresult sscoutsegment(face*, point, int, int, int);1799void scarveholes(int, REAL*);1800int triangulate(int, arraypool*, arraypool*, int, REAL*);18011802void unifysegments();1803void identifyinputedges(point*);1804void mergefacets();1805void meshsurface();18061807void interecursive(shellface** subfacearray, int arraysize, int axis,1808REAL, REAL, REAL, REAL, REAL, REAL, int* internum);1809void detectinterfaces();18101811///////////////////////////////////////////////////////////////////////////////1812// //1813// Constrained Delaunay tetrahedralization //1814// //1815// A constrained Delaunay tetrahedralization (CDT) is a variation of a Dela- //1816// unay tetrahedralization (DT) that is constrained to respect the boundary //1817// of a 3D PLC (domain). In a CDT of a 3D PLC, every vertex or edge of the //1818// PLC is also a vertex or an edge of the CDT, every polygon of the PLC is a //1819// union of triangles of the CDT. A crucial difference between a CDT and a //1820// DT is that triangles in the PLC's polygons are not required to be locally //1821// Delaunay, which frees the CDT to better respect the PLC's polygons. CDTs //1822// have optimal properties similar to those of DTs. //1823// //1824// Steiner Points and Steiner CDTs. It is known that even a simple 3D polyh- //1825// edron may not have a tetrahedralization which only uses its own vertices. //1826// Some extra points, so-called "Steiner points" are needed in order to form //1827// a tetrahedralization of such polyhedron. It is true for tetrahedralizing //1828// a 3D PLC as well. A Steiner CDT of a 3D PLC is a CDT containing Steiner //1829// points. The CDT algorithms of TetGen in general create Steiner CDTs. //1830// Almost all of the Steiner points are added in the edges of the PLC. They //1831// guarantee the existence of a CDT of the modified PLC. //1832// //1833// The routine constraineddelaunay() starts from a DT of the vertices of a //1834// PLC and creates a (Steiner) CDT of the PLC (including Steiner points). It //1835// is constructed by two steps, (1) segment recovery and (2) facet (polygon) //1836// recovery. Each step is accomplished by its own algorithm. //1837// //1838// The routine delaunizesegments() implements the segment recovery algorithm //1839// of Si, H. and Gaertner, K. "Meshing Piecewise Linear Complexes by Constr- //1840// ained Delaunay Tetrahedralizations," In Proceedings of the 14th Internat- //1841// ional Meshing Roundtable, 147--163, 2005. It adds Steiner points into //1842// non-Delaunay segments until all subsegments appear together in a DT. The //1843// running time of this algorithm is proportional to the number of added //1844// Steiner points. //1845// //1846// There are two incremental facet recovery algorithms: the cavity re-trian- //1847// gulation algorithm of Si, H. and Gaertner, K. "3D Boundary Recovery by //1848// Constrained Delaunay Tetrahedralization," International Journal for Numer-//1849// ical Methods in Engineering, 85:1341-1364, 2011, and the flip algorithm //1850// of Shewchuk, J. "Updating and Constructing Constrained Delaunay and //1851// Constrained Regular Triangulations by Flips." In Proceedings of the 19th //1852// ACM Symposium on Computational Geometry, 86-95, 2003. //1853// //1854// It is guaranteed in theory, no Steiner point is needed in both algorithms //1855// However, a facet with non-coplanar vertices might cause the additions of //1856// Steiner points. It is discussed in the paper of Si, H., and Shewchuk, J.,//1857// "Incrementally Constructing and Updating Constrained Delaunay //1858// Tetrahedralizations with Finite Precision Coordinates." In Proceedings of //1859// the 21th International Meshing Roundtable, 2012. //1860// //1861// Our implementation of the facet recovery algorithms recover a "missing //1862// region" at a time. Each missing region is a subset of connected interiors //1863// of a polygon. The routine formcavity() creates the cavity of crossing //1864// tetrahedra of the missing region. //1865// //1866// The cavity re-triangulation algorithm is implemented by three subroutines,//1867// delaunizecavity(), fillcavity(), and carvecavity(). Since it may fail due //1868// to non-coplanar vertices, the subroutine restorecavity() is used to rest- //1869// ore the original cavity. //1870// //1871// The routine flipinsertfacet() implements the flip algorithm. The subrout- //1872// ine flipcertify() is used to maintain the priority queue of flips. //1873// //1874// The routine refineregion() is called when the facet recovery algorithm //1875// fail to recover a missing region. It inserts Steiner points to refine the //1876// missing region. In order to avoid inserting Steiner points very close to //1877// existing segments. The classical encroachment rules of the Delaunay //1878// refinement algorithm are used to choose the Steiner points. //1879// //1880// The routine constrainedfacets() does the facet recovery by using either //1881// the cavity re-triangulation algorithm (default) or the flip algorithm. It //1882// results a CDT of the (modified) PLC (including Steiner points). //1883// //1884///////////////////////////////////////////////////////////////////////////////18851886void makesegmentendpointsmap();18871888enum interresult finddirection(triface* searchtet, point endpt);1889enum interresult scoutsegment(point, point, face*, triface*, point*,1890arraypool*);1891int getsteinerptonsegment(face* seg, point refpt, point steinpt);1892void delaunizesegments();18931894int scoutsubface(face* searchsh,triface* searchtet,int shflag);1895void formregion(face*, arraypool*, arraypool*, arraypool*);1896int scoutcrossedge(triface& crosstet, arraypool*, arraypool*);1897bool formcavity(triface*, arraypool*, arraypool*, arraypool*, arraypool*,1898arraypool*, arraypool*);1899// Facet recovery by cavity re-triangulation [Si and Gaertner 2011].1900void delaunizecavity(arraypool*, arraypool*, arraypool*, arraypool*,1901arraypool*, arraypool*);1902bool fillcavity(arraypool*, arraypool*, arraypool*, arraypool*,1903arraypool*, arraypool*, triface* crossedge);1904void carvecavity(arraypool*, arraypool*, arraypool*);1905void restorecavity(arraypool*, arraypool*, arraypool*, arraypool*);1906// Facet recovery by flips [Shewchuk 2003].1907void flipcertify(triface *chkface, badface **pqueue, point, point, point);1908void flipinsertfacet(arraypool*, arraypool*, arraypool*, arraypool*);19091910int insertpoint_cdt(point, triface*, face*, face*, insertvertexflags*,1911arraypool*, arraypool*, arraypool*, arraypool*,1912arraypool*, arraypool*);1913void refineregion(face&, arraypool*, arraypool*, arraypool*, arraypool*,1914arraypool*, arraypool*);1915void constrainedfacets();19161917void constraineddelaunay(clock_t&);19181919///////////////////////////////////////////////////////////////////////////////1920// //1921// Constrained tetrahedralizations. //1922// //1923///////////////////////////////////////////////////////////////////////////////19241925int checkflipeligibility(int fliptype, point, point, point, point, point,1926int level, int edgepivot, flipconstraints* fc);19271928int removeedgebyflips(triface*, flipconstraints*);1929int removefacebyflips(triface*, flipconstraints*);19301931int recoveredgebyflips(point, point, face*, triface*, int fullsearch);1932int add_steinerpt_in_schoenhardtpoly(triface*, int, int chkencflag);1933int add_steinerpt_in_segment(face*, int searchlevel);1934int addsteiner4recoversegment(face*, int);1935int recoversegments(arraypool*, int fullsearch, int steinerflag);19361937int recoverfacebyflips(point, point, point, face*, triface*);1938int recoversubfaces(arraypool*, int steinerflag);19391940int getvertexstar(int, point searchpt, arraypool*, arraypool*, arraypool*);1941int getedge(point, point, triface*);1942int reduceedgesatvertex(point startpt, arraypool* endptlist);1943int removevertexbyflips(point steinerpt);19441945int suppressbdrysteinerpoint(point steinerpt);1946int suppresssteinerpoints();19471948void recoverboundary(clock_t&);19491950///////////////////////////////////////////////////////////////////////////////1951// //1952// Mesh reconstruction //1953// //1954///////////////////////////////////////////////////////////////////////////////19551956void carveholes();19571958void reconstructmesh();19591960int scoutpoint(point, triface*, int randflag);1961REAL getpointmeshsize(point, triface*, int iloc);1962void interpolatemeshsize();19631964void insertconstrainedpoints(point *insertarray, int arylen, int rejflag);1965void insertconstrainedpoints(tetgenio *addio);19661967void collectremovepoints(arraypool *remptlist);1968void meshcoarsening();19691970///////////////////////////////////////////////////////////////////////////////1971// //1972// Mesh refinement //1973// //1974// The purpose of mesh refinement is to obtain a tetrahedral mesh with well- //1975// -shaped tetrahedra and appropriate mesh size. It is necessary to insert //1976// new Steiner points to achieve this property. The questions are (1) how to //1977// choose the Steiner points? and (2) how to insert them? //1978// //1979// Delaunay refinement is a technique first developed by Chew [1989] and //1980// Ruppert [1993, 1995] to generate quality triangular meshes in the plane. //1981// It provides guarantee on the smallest angle of the triangles. Rupper's //1982// algorithm guarantees that the mesh is size-optimal (to within a constant //1983// factor) among all meshes with the same quality. //1984// Shewchuk generalized Ruppert's algorithm into 3D in his PhD thesis //1985// [Shewchuk 1997]. A short version of his algorithm appears in "Tetrahedral //1986// Mesh Generation by Delaunay Refinement," In Proceedings of the 14th ACM //1987// Symposium on Computational Geometry, 86-95, 1998. It guarantees that all //1988// tetrahedra of the output mesh have a "radius-edge ratio" (equivalent to //1989// the minimal face angle) bounded. However, it does not remove slivers, a //1990// type of very flat tetrahedra which can have no small face angles but have //1991// very small (and large) dihedral angles. Moreover, it may not terminate if //1992// the input PLC contains "sharp features", e.g., two edges (or two facets) //1993// meet at an acute angle (or dihedral angle). //1994// //1995// TetGen uses the basic Delaunay refinement scheme to insert Steiner points.//1996// While it always maintains a constrained Delaunay mesh. The algorithm is //1997// described in Si, H., "Adaptive Constrained Delaunay Mesh Generation," //1998// International Journal for Numerical Methods in Engineering, 75:856-880. //1999// This algorithm always terminates and sharp features are easily preserved. //2000// The mesh has good quality (same as Shewchuk's Delaunay refinement algori- //2001// thm) in the bulk of the mesh domain. Moreover, it supports the generation //2002// of adaptive mesh according to a (isotropic) mesh sizing function. //2003// //2004///////////////////////////////////////////////////////////////////////////////20052006void makefacetverticesmap();2007int segsegadjacent(face *, face *);2008int segfacetadjacent(face *checkseg, face *checksh);2009int facetfacetadjacent(face *, face *);2010void save_segmentpoint_insradius(point segpt, point parentpt, REAL r);2011void save_facetpoint_insradius(point facpt, point parentpt, REAL r);2012void enqueuesubface(memorypool*, face*);2013void enqueuetetrahedron(triface*);20142015int checkseg4encroach(point pa, point pb, point checkpt);2016int checkseg4split(face *chkseg, point&, int&);2017int splitsegment(face *splitseg, point encpt, REAL, point, point, int, int);2018void repairencsegs(int chkencflag);20192020int checkfac4encroach(point, point, point, point checkpt, REAL*, REAL*);2021int checkfac4split(face *chkfac, point& encpt, int& qflag, REAL *ccent);2022int splitsubface(face *splitfac, point, point, int qflag, REAL *ccent, int);2023void repairencfacs(int chkencflag);20242025int checktet4split(triface *chktet, int& qflag, REAL *ccent);2026int splittetrahedron(triface* splittet,int qflag,REAL *ccent, int);2027void repairbadtets(int chkencflag);20282029void delaunayrefinement();20302031///////////////////////////////////////////////////////////////////////////////2032// //2033// Mesh optimization //2034// //2035///////////////////////////////////////////////////////////////////////////////20362037long lawsonflip3d(flipconstraints *fc);2038void recoverdelaunay();20392040int gettetrahedron(point, point, point, point, triface *);2041long improvequalitybyflips();20422043int smoothpoint(point smtpt, arraypool*, int ccw, optparameters *opm);2044long improvequalitybysmoothing(optparameters *opm);20452046int splitsliver(triface *, REAL, int);2047long removeslivers(int);20482049void optimizemesh();20502051///////////////////////////////////////////////////////////////////////////////2052// //2053// Mesh check and statistics //2054// //2055///////////////////////////////////////////////////////////////////////////////20562057// Mesh validations.2058int checkmesh(int topoflag);2059int checkshells();2060int checksegments();2061int checkdelaunay(int perturb = 1);2062int checkregular(int);2063int checkconforming(int);20642065// Mesh statistics.2066void printfcomma(unsigned long n);2067void qualitystatistics();2068void memorystatistics();2069void statistics();20702071///////////////////////////////////////////////////////////////////////////////2072// //2073// Mesh output //2074// //2075///////////////////////////////////////////////////////////////////////////////20762077void jettisonnodes();2078void highorder();2079void indexelements();2080void numberedges();2081void outnodes(tetgenio*);2082void outmetrics(tetgenio*);2083void outelements(tetgenio*);2084void outfaces(tetgenio*);2085void outhullfaces(tetgenio*);2086void outsubfaces(tetgenio*);2087void outedges(tetgenio*);2088void outsubsegments(tetgenio*);2089void outneighbors(tetgenio*);2090void outvoronoi(tetgenio*);2091void outsmesh(char*);2092void outmesh2medit(char*);2093void outmesh2vtk(char*);2094209520962097///////////////////////////////////////////////////////////////////////////////2098// //2099// Constructor & destructor //2100// //2101///////////////////////////////////////////////////////////////////////////////21022103void initializetetgenmesh()2104{2105in = addin = NULL;2106b = NULL;2107bgm = NULL;21082109tetrahedrons = subfaces = subsegs = points = NULL;2110badtetrahedrons = badsubfacs = badsubsegs = NULL;2111tet2segpool = tet2subpool = NULL;2112flippool = NULL;21132114dummypoint = NULL;2115flipstack = NULL;2116unflipqueue = NULL;21172118cavetetlist = cavebdrylist = caveoldtetlist = NULL;2119cavetetshlist = cavetetseglist = cavetetvertlist = NULL;2120caveencshlist = caveencseglist = NULL;2121caveshlist = caveshbdlist = cavesegshlist = NULL;21222123subsegstack = subfacstack = subvertstack = NULL;2124encseglist = encshlist = NULL;2125idx2facetlist = NULL;2126facetverticeslist = NULL;2127segmentendpointslist = NULL;21282129highordertable = NULL;21302131numpointattrib = numelemattrib = 0;2132sizeoftensor = 0;2133pointmtrindex = 0;2134pointparamindex = 0;2135pointmarkindex = 0;2136point2simindex = 0;2137pointinsradiusindex = 0;2138elemattribindex = 0;2139volumeboundindex = 0;2140shmarkindex = 0;2141areaboundindex = 0;2142checksubsegflag = 0;2143checksubfaceflag = 0;2144checkconstraints = 0;2145nonconvex = 0;2146autofliplinklevel = 1;2147useinsertradius = 0;2148samples = 0l;2149randomseed = 1l;2150minfaceang = minfacetdihed = PI;2151tetprism_vol_sum = 0.0;2152longest = minedgelength = 0.0;2153xmax = xmin = ymax = ymin = zmax = zmin = 0.0;21542155insegments = 0l;2156hullsize = 0l;2157meshedges = meshhulledges = 0l;2158steinerleft = -1;2159dupverts = 0l;2160unuverts = 0l;2161nonregularcount = 0l;2162st_segref_count = st_facref_count = st_volref_count = 0l;2163fillregioncount = cavitycount = cavityexpcount = 0l;2164flip14count = flip26count = flipn2ncount = 0l;2165flip23count = flip32count = flip44count = flip41count = 0l;2166flip22count = flip31count = 0l;2167totalworkmemory = 0l;216821692170} // tetgenmesh()21712172void freememory()2173{2174if (bgm != NULL) {2175delete bgm;2176}21772178if (points != (memorypool *) NULL) {2179delete points;2180delete [] dummypoint;2181}2182if (tetrahedrons != (memorypool *) NULL) {2183delete tetrahedrons;2184}2185if (subfaces != (memorypool *) NULL) {2186delete subfaces;2187delete subsegs;2188}2189if (tet2segpool != NULL) {2190delete tet2segpool;2191delete tet2subpool;2192}21932194if (badtetrahedrons) {2195delete badtetrahedrons;2196}2197if (badsubfacs) {2198delete badsubfacs;2199}2200if (badsubsegs) {2201delete badsubsegs;2202}2203if (encseglist) {2204delete encseglist;2205}2206if (encshlist) {2207delete encshlist;2208}22092210if (flippool != NULL) {2211delete flippool;2212delete unflipqueue;2213}22142215if (cavetetlist != NULL) {2216delete cavetetlist;2217delete cavebdrylist;2218delete caveoldtetlist;2219delete cavetetvertlist;2220}22212222if (caveshlist != NULL) {2223delete caveshlist;2224delete caveshbdlist;2225delete cavesegshlist;2226delete cavetetshlist;2227delete cavetetseglist;2228delete caveencshlist;2229delete caveencseglist;2230}22312232if (subsegstack != NULL) {2233delete subsegstack;2234delete subfacstack;2235delete subvertstack;2236}22372238if (idx2facetlist != NULL) {2239delete [] idx2facetlist;2240delete [] facetverticeslist;2241}22422243if (segmentendpointslist != NULL) {2244delete [] segmentendpointslist;2245}22462247if (highordertable != NULL) {2248delete [] highordertable;2249}22502251initializetetgenmesh();2252}22532254tetgenmesh()2255{2256initializetetgenmesh();2257}22582259~tetgenmesh()2260{2261freememory();2262} // ~tetgenmesh()22632264}; // End of class tetgenmesh.22652266///////////////////////////////////////////////////////////////////////////////2267// //2268// tetrahedralize() Interface for using TetGen's library to generate //2269// Delaunay tetrahedralizations, constrained Delaunay //2270// tetrahedralizations, quality tetrahedral meshes. //2271// //2272// 'in' is an object of 'tetgenio' which contains a PLC you want to tetrahed-//2273// ralize or a previously generated tetrahedral mesh you want to refine. It //2274// must not be a NULL. 'out' is another object of 'tetgenio' for storing the //2275// generated tetrahedral mesh. It can be a NULL. If so, the output will be //2276// saved to file(s). If 'bgmin' != NULL, it contains a background mesh which //2277// defines a mesh size function. //2278// //2279///////////////////////////////////////////////////////////////////////////////22802281void tetrahedralize(tetgenbehavior *b, tetgenio *in, tetgenio *out,2282tetgenio *addin = NULL, tetgenio *bgmin = NULL);22832284#ifdef TETLIBRARY2285void tetrahedralize(char *switches, tetgenio *in, tetgenio *out,2286tetgenio *addin = NULL, tetgenio *bgmin = NULL);22872288extern "C"2289#ifdef WIN322290__declspec(dllexport)2291#endif2292void delegate_tetrahedralize(int bs, tetgenbehavior *b, char *switches,2293tetgenio *in, tetgenio *out, tetgenio *addin = NULL, tetgenio *bgmin = NULL);22942295#endif // #ifdef TETLIBRARY22962297///////////////////////////////////////////////////////////////////////////////2298// //2299// terminatetetgen() Terminate TetGen with a given exit code. //2300// //2301///////////////////////////////////////////////////////////////////////////////23022303inline void terminatetetgen(tetgenmesh *m, int x)2304{2305#ifdef TETLIBRARY2306throw x;2307#else2308switch (x) {2309case 1: // Out of memory.2310printf("Error: Out of memory.\n");2311break;2312case 2: // Encounter an internal error.2313printf("Please report this bug to [email protected]. Include\n");2314printf(" the message above, your input data set, and the exact\n");2315printf(" command line you used to run this program, thank you.\n");2316break;2317case 3:2318printf("A self-intersection was detected. Program stopped.\n");2319printf("Hint: use -d option to detect all self-intersections.\n");2320break;2321case 4:2322printf("A very small input feature size was detected. Program stopped.\n");2323if (m) {2324printf("Hint: use -T option to set a smaller tolerance. Current is %g\n",2325m->b->epsilon);2326}2327break;2328case 5:2329printf("Two very close input facets were detected. Program stopped.\n");2330printf("Hint: use -Y option to avoid adding Steiner points in boundary.\n");2331break;2332case 10:2333printf("An input error was detected. Program stopped.\n");2334break;2335} // switch (x)2336exit(x);2337#endif // #ifdef TETLIBRARY2338}23392340///////////////////////////////////////////////////////////////////////////////2341// //2342// Primitives for tetrahedra //2343// //2344///////////////////////////////////////////////////////////////////////////////23452346// encode() compress a handle into a single pointer. It relies on the2347// assumption that all addresses of tetrahedra are aligned to sixteen-2348// byte boundaries, so that the last four significant bits are zero.23492350inline tetgenmesh::tetrahedron tetgenmesh::encode(triface& t) {2351return (tetrahedron) ((uintptr_t) (t).tet | (uintptr_t) (t).ver);2352}23532354inline tetgenmesh::tetrahedron tetgenmesh::encode2(tetrahedron* ptr, int ver) {2355return (tetrahedron) ((uintptr_t) (ptr) | (uintptr_t) (ver));2356}23572358// decode() converts a pointer to a handle. The version is extracted from2359// the four least significant bits of the pointer.23602361inline void tetgenmesh::decode(tetrahedron ptr, triface& t) {2362(t).ver = (int) ((uintptr_t) (ptr) & (uintptr_t) 15);2363(t).tet = (tetrahedron *) ((uintptr_t) (ptr) ^ (uintptr_t) (t).ver);2364}23652366// bond() connects two tetrahedra together. (t1,v1) and (t2,v2) must2367// refer to the same face and the same edge.23682369inline void tetgenmesh::bond(triface& t1, triface& t2) {2370t1.tet[t1.ver & 3] = encode2(t2.tet, bondtbl[t1.ver][t2.ver]);2371t2.tet[t2.ver & 3] = encode2(t1.tet, bondtbl[t2.ver][t1.ver]);2372}237323742375// dissolve() a bond (from one side).23762377inline void tetgenmesh::dissolve(triface& t) {2378t.tet[t.ver & 3] = NULL;2379}23802381// enext() finds the next edge (counterclockwise) in the same face.23822383inline void tetgenmesh::enext(triface& t1, triface& t2) {2384t2.tet = t1.tet;2385t2.ver = enexttbl[t1.ver];2386}23872388inline void tetgenmesh::enextself(triface& t) {2389t.ver = enexttbl[t.ver];2390}23912392// eprev() finds the next edge (clockwise) in the same face.23932394inline void tetgenmesh::eprev(triface& t1, triface& t2) {2395t2.tet = t1.tet;2396t2.ver = eprevtbl[t1.ver];2397}23982399inline void tetgenmesh::eprevself(triface& t) {2400t.ver = eprevtbl[t.ver];2401}24022403// esym() finds the reversed edge. It is in the other face of the2404// same tetrahedron.24052406inline void tetgenmesh::esym(triface& t1, triface& t2) {2407(t2).tet = (t1).tet;2408(t2).ver = esymtbl[(t1).ver];2409}24102411inline void tetgenmesh::esymself(triface& t) {2412(t).ver = esymtbl[(t).ver];2413}24142415// enextesym() finds the reversed edge of the next edge. It is in the other2416// face of the same tetrahedron. It is the combination esym() * enext().24172418inline void tetgenmesh::enextesym(triface& t1, triface& t2) {2419t2.tet = t1.tet;2420t2.ver = enextesymtbl[t1.ver];2421}24222423inline void tetgenmesh::enextesymself(triface& t) {2424t.ver = enextesymtbl[t.ver];2425}24262427// eprevesym() finds the reversed edge of the previous edge.24282429inline void tetgenmesh::eprevesym(triface& t1, triface& t2) {2430t2.tet = t1.tet;2431t2.ver = eprevesymtbl[t1.ver];2432}24332434inline void tetgenmesh::eprevesymself(triface& t) {2435t.ver = eprevesymtbl[t.ver];2436}24372438// eorgoppo() Finds the opposite face of the origin of the current edge.2439// Return the opposite edge of the current edge.24402441inline void tetgenmesh::eorgoppo(triface& t1, triface& t2) {2442t2.tet = t1.tet;2443t2.ver = eorgoppotbl[t1.ver];2444}24452446inline void tetgenmesh::eorgoppoself(triface& t) {2447t.ver = eorgoppotbl[t.ver];2448}24492450// edestoppo() Finds the opposite face of the destination of the current2451// edge. Return the opposite edge of the current edge.24522453inline void tetgenmesh::edestoppo(triface& t1, triface& t2) {2454t2.tet = t1.tet;2455t2.ver = edestoppotbl[t1.ver];2456}24572458inline void tetgenmesh::edestoppoself(triface& t) {2459t.ver = edestoppotbl[t.ver];2460}24612462// fsym() finds the adjacent tetrahedron at the same face and the same edge.24632464inline void tetgenmesh::fsym(triface& t1, triface& t2) {2465decode((t1).tet[(t1).ver & 3], t2);2466t2.ver = fsymtbl[t1.ver][t2.ver];2467}246824692470#define fsymself(t) \2471t1ver = (t).ver; \2472decode((t).tet[(t).ver & 3], (t));\2473(t).ver = fsymtbl[t1ver][(t).ver]24742475// fnext() finds the next face while rotating about an edge according to2476// a right-hand rule. The face is in the adjacent tetrahedron. It is2477// the combination: fsym() * esym().24782479inline void tetgenmesh::fnext(triface& t1, triface& t2) {2480decode(t1.tet[facepivot1[t1.ver]], t2);2481t2.ver = facepivot2[t1.ver][t2.ver];2482}248324842485#define fnextself(t) \2486t1ver = (t).ver; \2487decode((t).tet[facepivot1[(t).ver]], (t)); \2488(t).ver = facepivot2[t1ver][(t).ver]248924902491// The following primtives get or set the origin, destination, face apex,2492// or face opposite of an ordered tetrahedron.24932494inline tetgenmesh::point tetgenmesh::org(triface& t) {2495return (point) (t).tet[orgpivot[(t).ver]];2496}24972498inline tetgenmesh::point tetgenmesh:: dest(triface& t) {2499return (point) (t).tet[destpivot[(t).ver]];2500}25012502inline tetgenmesh::point tetgenmesh:: apex(triface& t) {2503return (point) (t).tet[apexpivot[(t).ver]];2504}25052506inline tetgenmesh::point tetgenmesh:: oppo(triface& t) {2507return (point) (t).tet[oppopivot[(t).ver]];2508}25092510inline void tetgenmesh:: setorg(triface& t, point p) {2511(t).tet[orgpivot[(t).ver]] = (tetrahedron) (p);2512}25132514inline void tetgenmesh:: setdest(triface& t, point p) {2515(t).tet[destpivot[(t).ver]] = (tetrahedron) (p);2516}25172518inline void tetgenmesh:: setapex(triface& t, point p) {2519(t).tet[apexpivot[(t).ver]] = (tetrahedron) (p);2520}25212522inline void tetgenmesh:: setoppo(triface& t, point p) {2523(t).tet[oppopivot[(t).ver]] = (tetrahedron) (p);2524}25252526#define setvertices(t, torg, tdest, tapex, toppo) \2527(t).tet[orgpivot[(t).ver]] = (tetrahedron) (torg);\2528(t).tet[destpivot[(t).ver]] = (tetrahedron) (tdest); \2529(t).tet[apexpivot[(t).ver]] = (tetrahedron) (tapex); \2530(t).tet[oppopivot[(t).ver]] = (tetrahedron) (toppo)25312532// Check or set a tetrahedron's attributes.25332534inline REAL tetgenmesh::elemattribute(tetrahedron* ptr, int attnum) {2535return ((REAL *) (ptr))[elemattribindex + attnum];2536}25372538inline void tetgenmesh::setelemattribute(tetrahedron* ptr, int attnum,2539REAL value) {2540((REAL *) (ptr))[elemattribindex + attnum] = value;2541}25422543// Check or set a tetrahedron's maximum volume bound.25442545inline REAL tetgenmesh::volumebound(tetrahedron* ptr) {2546return ((REAL *) (ptr))[volumeboundindex];2547}25482549inline void tetgenmesh::setvolumebound(tetrahedron* ptr, REAL value) {2550((REAL *) (ptr))[volumeboundindex] = value;2551}25522553// Get or set a tetrahedron's index (only used for output).2554// These two routines use the reserved slot ptr[10].25552556inline int tetgenmesh::elemindex(tetrahedron* ptr) {2557int *iptr = (int *) &(ptr[10]);2558return iptr[0];2559}25602561inline void tetgenmesh::setelemindex(tetrahedron* ptr, int value) {2562int *iptr = (int *) &(ptr[10]);2563iptr[0] = value;2564}25652566// Get or set a tetrahedron's marker.2567// Set 'value = 0' cleans all the face/edge flags.25682569inline int tetgenmesh::elemmarker(tetrahedron* ptr) {2570return ((int *) (ptr))[elemmarkerindex];2571}25722573inline void tetgenmesh::setelemmarker(tetrahedron* ptr, int value) {2574((int *) (ptr))[elemmarkerindex] = value;2575}25762577// infect(), infected(), uninfect() -- primitives to flag or unflag a2578// tetrahedron. The last bit of the element marker is flagged (1)2579// or unflagged (0).25802581inline void tetgenmesh::infect(triface& t) {2582((int *) (t.tet))[elemmarkerindex] |= 1;2583}25842585inline void tetgenmesh::uninfect(triface& t) {2586((int *) (t.tet))[elemmarkerindex] &= ~1;2587}25882589inline bool tetgenmesh::infected(triface& t) {2590return (((int *) (t.tet))[elemmarkerindex] & 1) != 0;2591}25922593// marktest(), marktested(), unmarktest() -- primitives to flag or unflag a2594// tetrahedron. Use the second lowerest bit of the element marker.25952596inline void tetgenmesh::marktest(triface& t) {2597((int *) (t.tet))[elemmarkerindex] |= 2;2598}25992600inline void tetgenmesh::unmarktest(triface& t) {2601((int *) (t.tet))[elemmarkerindex] &= ~2;2602}26032604inline bool tetgenmesh::marktested(triface& t) {2605return (((int *) (t.tet))[elemmarkerindex] & 2) != 0;2606}26072608// markface(), unmarkface(), facemarked() -- primitives to flag or unflag a2609// face of a tetrahedron. From the last 3rd to 6th bits are used for2610// face markers, e.g., the last third bit corresponds to loc = 0.26112612inline void tetgenmesh::markface(triface& t) {2613((int *) (t.tet))[elemmarkerindex] |= (4 << (t.ver & 3));2614}26152616inline void tetgenmesh::unmarkface(triface& t) {2617((int *) (t.tet))[elemmarkerindex] &= ~(4 << (t.ver & 3));2618}26192620inline bool tetgenmesh::facemarked(triface& t) {2621return (((int *) (t.tet))[elemmarkerindex] & (4 << (t.ver & 3))) != 0;2622}26232624// markedge(), unmarkedge(), edgemarked() -- primitives to flag or unflag an2625// edge of a tetrahedron. From the last 7th to 12th bits are used for2626// edge markers, e.g., the last 7th bit corresponds to the 0th edge, etc.2627// Remark: The last 7th bit is marked by 2^6 = 64.26282629inline void tetgenmesh::markedge(triface& t) {2630((int *) (t.tet))[elemmarkerindex] |= (int) (64 << ver2edge[(t).ver]);2631}26322633inline void tetgenmesh::unmarkedge(triface& t) {2634((int *) (t.tet))[elemmarkerindex] &= ~(int) (64 << ver2edge[(t).ver]);2635}26362637inline bool tetgenmesh::edgemarked(triface& t) {2638return (((int *) (t.tet))[elemmarkerindex] &2639(int) (64 << ver2edge[(t).ver])) != 0;2640}26412642// marktest2(), unmarktest2(), marktest2ed() -- primitives to flag and unflag2643// a tetrahedron. The 13th bit (2^12 = 4096) is used for this flag.26442645inline void tetgenmesh::marktest2(triface& t) {2646((int *) (t.tet))[elemmarkerindex] |= (int) (4096);2647}26482649inline void tetgenmesh::unmarktest2(triface& t) {2650((int *) (t.tet))[elemmarkerindex] &= ~(int) (4096);2651}26522653inline bool tetgenmesh::marktest2ed(triface& t) {2654return (((int *) (t.tet))[elemmarkerindex] & (int) (4096)) != 0;2655}26562657// elemcounter(), setelemcounter() -- primitives to read or set a (small)2658// integer counter in this tet. It is saved from the 16th bit. On 32 bit2659// system, the range of the counter is [0, 2^15 = 32768].26602661inline int tetgenmesh::elemcounter(triface& t) {2662return (((int *) (t.tet))[elemmarkerindex]) >> 16;2663}26642665inline void tetgenmesh::setelemcounter(triface& t, int value) {2666int c = ((int *) (t.tet))[elemmarkerindex];2667// Clear the old counter while keep the other flags.2668c &= 65535; // sum_{i=0^15} 2^i2669c |= (value << 16);2670((int *) (t.tet))[elemmarkerindex] = c;2671}26722673inline void tetgenmesh::increaseelemcounter(triface& t) {2674int c = elemcounter(t);2675setelemcounter(t, c + 1);2676}26772678inline void tetgenmesh::decreaseelemcounter(triface& t) {2679int c = elemcounter(t);2680setelemcounter(t, c - 1);2681}26822683// ishulltet() tests if t is a hull tetrahedron.26842685inline bool tetgenmesh::ishulltet(triface& t) {2686return (point) (t).tet[7] == dummypoint;2687}26882689// isdeadtet() tests if t is a tetrahedron is dead.26902691inline bool tetgenmesh::isdeadtet(triface& t) {2692return ((t.tet == NULL) || (t.tet[4] == NULL));2693}26942695///////////////////////////////////////////////////////////////////////////////2696// //2697// Primitives for subfaces and subsegments //2698// //2699///////////////////////////////////////////////////////////////////////////////27002701// Each subface contains three pointers to its neighboring subfaces, with2702// edge versions. To save memory, both information are kept in a single2703// pointer. To make this possible, all subfaces are aligned to eight-byte2704// boundaries, so that the last three bits of each pointer are zeros. An2705// edge version (in the range 0 to 5) is compressed into the last three2706// bits of each pointer by 'sencode()'. 'sdecode()' decodes a pointer,2707// extracting an edge version and a pointer to the beginning of a subface.27082709inline void tetgenmesh::sdecode(shellface sptr, face& s) {2710s.shver = (int) ((uintptr_t) (sptr) & (uintptr_t) 7);2711s.sh = (shellface *) ((uintptr_t) (sptr) ^ (uintptr_t) (s.shver));2712}27132714inline tetgenmesh::shellface tetgenmesh::sencode(face& s) {2715return (shellface) ((uintptr_t) s.sh | (uintptr_t) s.shver);2716}27172718inline tetgenmesh::shellface tetgenmesh::sencode2(shellface *sh, int shver) {2719return (shellface) ((uintptr_t) sh | (uintptr_t) shver);2720}27212722// sbond() bonds two subfaces (s1) and (s2) together. s1 and s2 must refer2723// to the same edge. No requirement is needed on their orientations.27242725inline void tetgenmesh::sbond(face& s1, face& s2)2726{2727s1.sh[s1.shver >> 1] = sencode(s2);2728s2.sh[s2.shver >> 1] = sencode(s1);2729}27302731// sbond1() bonds s1 <== s2, i.e., after bonding, s1 is pointing to s2,2732// but s2 is not pointing to s1. s1 and s2 must refer to the same edge.2733// No requirement is needed on their orientations.27342735inline void tetgenmesh::sbond1(face& s1, face& s2)2736{2737s1.sh[s1.shver >> 1] = sencode(s2);2738}27392740// Dissolve a subface bond (from one side). Note that the other subface2741// will still think it's connected to this subface.27422743inline void tetgenmesh::sdissolve(face& s)2744{2745s.sh[s.shver >> 1] = NULL;2746}27472748// spivot() finds the adjacent subface (s2) for a given subface (s1).2749// s1 and s2 share at the same edge.27502751inline void tetgenmesh::spivot(face& s1, face& s2)2752{2753shellface sptr = s1.sh[s1.shver >> 1];2754sdecode(sptr, s2);2755}27562757inline void tetgenmesh::spivotself(face& s)2758{2759shellface sptr = s.sh[s.shver >> 1];2760sdecode(sptr, s);2761}27622763// These primitives determine or set the origin, destination, or apex2764// of a subface with respect to the edge version.27652766inline tetgenmesh::point tetgenmesh::sorg(face& s)2767{2768return (point) s.sh[sorgpivot[s.shver]];2769}27702771inline tetgenmesh::point tetgenmesh::sdest(face& s)2772{2773return (point) s.sh[sdestpivot[s.shver]];2774}27752776inline tetgenmesh::point tetgenmesh::sapex(face& s)2777{2778return (point) s.sh[sapexpivot[s.shver]];2779}27802781inline void tetgenmesh::setsorg(face& s, point pointptr)2782{2783s.sh[sorgpivot[s.shver]] = (shellface) pointptr;2784}27852786inline void tetgenmesh::setsdest(face& s, point pointptr)2787{2788s.sh[sdestpivot[s.shver]] = (shellface) pointptr;2789}27902791inline void tetgenmesh::setsapex(face& s, point pointptr)2792{2793s.sh[sapexpivot[s.shver]] = (shellface) pointptr;2794}27952796#define setshvertices(s, pa, pb, pc)\2797setsorg(s, pa);\2798setsdest(s, pb);\2799setsapex(s, pc)28002801// sesym() reserves the direction of the lead edge.28022803inline void tetgenmesh::sesym(face& s1, face& s2)2804{2805s2.sh = s1.sh;2806s2.shver = (s1.shver ^ 1); // Inverse the last bit.2807}28082809inline void tetgenmesh::sesymself(face& s)2810{2811s.shver ^= 1;2812}28132814// senext() finds the next edge (counterclockwise) in the same orientation2815// of this face.28162817inline void tetgenmesh::senext(face& s1, face& s2)2818{2819s2.sh = s1.sh;2820s2.shver = snextpivot[s1.shver];2821}28222823inline void tetgenmesh::senextself(face& s)2824{2825s.shver = snextpivot[s.shver];2826}28272828inline void tetgenmesh::senext2(face& s1, face& s2)2829{2830s2.sh = s1.sh;2831s2.shver = snextpivot[snextpivot[s1.shver]];2832}28332834inline void tetgenmesh::senext2self(face& s)2835{2836s.shver = snextpivot[snextpivot[s.shver]];2837}283828392840// Check or set a subface's maximum area bound.28412842inline REAL tetgenmesh::areabound(face& s)2843{2844return ((REAL *) (s.sh))[areaboundindex];2845}28462847inline void tetgenmesh::setareabound(face& s, REAL value)2848{2849((REAL *) (s.sh))[areaboundindex] = value;2850}28512852// These two primitives read or set a shell marker. Shell markers are used2853// to hold user boundary information.28542855inline int tetgenmesh::shellmark(face& s)2856{2857return ((int *) (s.sh))[shmarkindex];2858}28592860inline void tetgenmesh::setshellmark(face& s, int value)2861{2862((int *) (s.sh))[shmarkindex] = value;2863}2864286528662867// sinfect(), sinfected(), suninfect() -- primitives to flag or unflag a2868// subface. The last bit of ((int *) ((s).sh))[shmarkindex+1] is flagged.28692870inline void tetgenmesh::sinfect(face& s)2871{2872((int *) ((s).sh))[shmarkindex+1] =2873(((int *) ((s).sh))[shmarkindex+1] | (int) 1);2874}28752876inline void tetgenmesh::suninfect(face& s)2877{2878((int *) ((s).sh))[shmarkindex+1] =2879(((int *) ((s).sh))[shmarkindex+1] & ~(int) 1);2880}28812882// Test a subface for viral infection.28832884inline bool tetgenmesh::sinfected(face& s)2885{2886return (((int *) ((s).sh))[shmarkindex+1] & (int) 1) != 0;2887}28882889// smarktest(), smarktested(), sunmarktest() -- primitives to flag or unflag2890// a subface. The last 2nd bit of the integer is flagged.28912892inline void tetgenmesh::smarktest(face& s)2893{2894((int *) ((s).sh))[shmarkindex+1] =2895(((int *)((s).sh))[shmarkindex+1] | (int) 2);2896}28972898inline void tetgenmesh::sunmarktest(face& s)2899{2900((int *) ((s).sh))[shmarkindex+1] =2901(((int *)((s).sh))[shmarkindex+1] & ~(int)2);2902}29032904inline bool tetgenmesh::smarktested(face& s)2905{2906return ((((int *) ((s).sh))[shmarkindex+1] & (int) 2) != 0);2907}29082909// smarktest2(), smarktest2ed(), sunmarktest2() -- primitives to flag or2910// unflag a subface. The last 3rd bit of the integer is flagged.29112912inline void tetgenmesh::smarktest2(face& s)2913{2914((int *) ((s).sh))[shmarkindex+1] =2915(((int *)((s).sh))[shmarkindex+1] | (int) 4);2916}29172918inline void tetgenmesh::sunmarktest2(face& s)2919{2920((int *) ((s).sh))[shmarkindex+1] =2921(((int *)((s).sh))[shmarkindex+1] & ~(int)4);2922}29232924inline bool tetgenmesh::smarktest2ed(face& s)2925{2926return ((((int *) ((s).sh))[shmarkindex+1] & (int) 4) != 0);2927}29282929// The last 4th bit of ((int *) ((s).sh))[shmarkindex+1] is flagged.29302931inline void tetgenmesh::smarktest3(face& s)2932{2933((int *) ((s).sh))[shmarkindex+1] =2934(((int *)((s).sh))[shmarkindex+1] | (int) 8);2935}29362937inline void tetgenmesh::sunmarktest3(face& s)2938{2939((int *) ((s).sh))[shmarkindex+1] =2940(((int *)((s).sh))[shmarkindex+1] & ~(int)8);2941}29422943inline bool tetgenmesh::smarktest3ed(face& s)2944{2945return ((((int *) ((s).sh))[shmarkindex+1] & (int) 8) != 0);2946}294729482949// Each facet has a unique index (automatically indexed). Starting from '0'.2950// We save this index in the same field of the shell type.29512952inline void tetgenmesh::setfacetindex(face& s, int value)2953{2954((int *) (s.sh))[shmarkindex + 2] = value;2955}29562957inline int tetgenmesh::getfacetindex(face& s)2958{2959return ((int *) (s.sh))[shmarkindex + 2];2960}29612962///////////////////////////////////////////////////////////////////////////////2963// //2964// Primitives for interacting between tetrahedra and subfaces //2965// //2966///////////////////////////////////////////////////////////////////////////////29672968// tsbond() bond a tetrahedron (t) and a subface (s) together.2969// Note that t and s must be the same face and the same edge. Moreover,2970// t and s have the same orientation.2971// Since the edge number in t and in s can be any number in {0,1,2}. We bond2972// the edge in s which corresponds to t's 0th edge, and vice versa.29732974inline void tetgenmesh::tsbond(triface& t, face& s)2975{2976if ((t).tet[9] == NULL) {2977// Allocate space for this tet.2978(t).tet[9] = (tetrahedron) tet2subpool->alloc();2979// Initialize.2980for (int i = 0; i < 4; i++) {2981((shellface *) (t).tet[9])[i] = NULL;2982}2983}2984// Bond t <== s.2985((shellface *) (t).tet[9])[(t).ver & 3] =2986sencode2((s).sh, tsbondtbl[t.ver][s.shver]);2987// Bond s <== t.2988s.sh[9 + ((s).shver & 1)] =2989(shellface) encode2((t).tet, stbondtbl[t.ver][s.shver]);2990}29912992// tspivot() finds a subface (s) abutting on the given tetrahdera (t).2993// Return s.sh = NULL if there is no subface at t. Otherwise, return2994// the subface s, and s and t must be at the same edge with the same2995// orientation.29962997inline void tetgenmesh::tspivot(triface& t, face& s)2998{2999if ((t).tet[9] == NULL) {3000(s).sh = NULL;3001return;3002}3003// Get the attached subface s.3004sdecode(((shellface *) (t).tet[9])[(t).ver & 3], (s));3005(s).shver = tspivottbl[t.ver][s.shver];3006}30073008// Quickly check if the handle (t, v) is a subface.3009#define issubface(t) \3010((t).tet[9] && ((t).tet[9])[(t).ver & 3])30113012// stpivot() finds a tetrahedron (t) abutting a given subface (s).3013// Return the t (if it exists) with the same edge and the same3014// orientation of s.30153016inline void tetgenmesh::stpivot(face& s, triface& t)3017{3018decode((tetrahedron) s.sh[9 + (s.shver & 1)], t);3019if ((t).tet == NULL) {3020return;3021}3022(t).ver = stpivottbl[t.ver][s.shver];3023}30243025// Quickly check if this subface is attached to a tetrahedron.30263027#define isshtet(s) \3028((s).sh[9 + ((s).shver & 1)])30293030// tsdissolve() dissolve a bond (from the tetrahedron side).30313032inline void tetgenmesh::tsdissolve(triface& t)3033{3034if ((t).tet[9] != NULL) {3035((shellface *) (t).tet[9])[(t).ver & 3] = NULL;3036}3037}30383039// stdissolve() dissolve a bond (from the subface side).30403041inline void tetgenmesh::stdissolve(face& s)3042{3043(s).sh[9] = NULL;3044(s).sh[10] = NULL;3045}30463047///////////////////////////////////////////////////////////////////////////////3048// //3049// Primitives for interacting between subfaces and segments //3050// //3051///////////////////////////////////////////////////////////////////////////////30523053// ssbond() bond a subface to a subsegment.30543055inline void tetgenmesh::ssbond(face& s, face& edge)3056{3057s.sh[6 + (s.shver >> 1)] = sencode(edge);3058edge.sh[0] = sencode(s);3059}30603061inline void tetgenmesh::ssbond1(face& s, face& edge)3062{3063s.sh[6 + (s.shver >> 1)] = sencode(edge);3064//edge.sh[0] = sencode(s);3065}30663067// ssdisolve() dissolve a bond (from the subface side)30683069inline void tetgenmesh::ssdissolve(face& s)3070{3071s.sh[6 + (s.shver >> 1)] = NULL;3072}30733074// sspivot() finds a subsegment abutting a subface.30753076inline void tetgenmesh::sspivot(face& s, face& edge)3077{3078sdecode((shellface) s.sh[6 + (s.shver >> 1)], edge);3079}30803081// Quickly check if the edge is a subsegment.30823083#define isshsubseg(s) \3084((s).sh[6 + ((s).shver >> 1)])30853086///////////////////////////////////////////////////////////////////////////////3087// //3088// Primitives for interacting between tetrahedra and segments //3089// //3090///////////////////////////////////////////////////////////////////////////////30913092inline void tetgenmesh::tssbond1(triface& t, face& s)3093{3094if ((t).tet[8] == NULL) {3095// Allocate space for this tet.3096(t).tet[8] = (tetrahedron) tet2segpool->alloc();3097// Initialization.3098for (int i = 0; i < 6; i++) {3099((shellface *) (t).tet[8])[i] = NULL;3100}3101}3102((shellface *) (t).tet[8])[ver2edge[(t).ver]] = sencode((s));3103}31043105inline void tetgenmesh::sstbond1(face& s, triface& t)3106{3107((tetrahedron *) (s).sh)[9] = encode(t);3108}31093110inline void tetgenmesh::tssdissolve1(triface& t)3111{3112if ((t).tet[8] != NULL) {3113((shellface *) (t).tet[8])[ver2edge[(t).ver]] = NULL;3114}3115}31163117inline void tetgenmesh::sstdissolve1(face& s)3118{3119((tetrahedron *) (s).sh)[9] = NULL;3120}31213122inline void tetgenmesh::tsspivot1(triface& t, face& s)3123{3124if ((t).tet[8] != NULL) {3125sdecode(((shellface *) (t).tet[8])[ver2edge[(t).ver]], s);3126} else {3127(s).sh = NULL;3128}3129}31303131// Quickly check whether 't' is a segment or not.31323133#define issubseg(t) \3134((t).tet[8] && ((t).tet[8])[ver2edge[(t).ver]])31353136inline void tetgenmesh::sstpivot1(face& s, triface& t)3137{3138decode((tetrahedron) s.sh[9], t);3139}31403141///////////////////////////////////////////////////////////////////////////////3142// //3143// Primitives for points //3144// //3145///////////////////////////////////////////////////////////////////////////////31463147inline int tetgenmesh::pointmark(point pt) {3148return ((int *) (pt))[pointmarkindex];3149}31503151inline void tetgenmesh::setpointmark(point pt, int value) {3152((int *) (pt))[pointmarkindex] = value;3153}315431553156// These two primitives set and read the type of the point.31573158inline enum tetgenmesh::verttype tetgenmesh::pointtype(point pt) {3159return (enum verttype) (((int *) (pt))[pointmarkindex + 1] >> (int) 8);3160}31613162inline void tetgenmesh::setpointtype(point pt, enum verttype value) {3163((int *) (pt))[pointmarkindex + 1] =3164((int) value << 8) + (((int *) (pt))[pointmarkindex + 1] & (int) 255);3165}31663167// Read and set the geometry tag of the point (used by -s option).31683169inline int tetgenmesh::pointgeomtag(point pt) {3170return ((int *) (pt))[pointmarkindex + 2];3171}31723173inline void tetgenmesh::setpointgeomtag(point pt, int value) {3174((int *) (pt))[pointmarkindex + 2] = value;3175}31763177// Read and set the u,v coordinates of the point (used by -s option).31783179inline REAL tetgenmesh::pointgeomuv(point pt, int i) {3180return pt[pointparamindex + i];3181}31823183inline void tetgenmesh::setpointgeomuv(point pt, int i, REAL value) {3184pt[pointparamindex + i] = value;3185}31863187// pinfect(), puninfect(), pinfected() -- primitives to flag or unflag3188// a point. The last bit of the integer '[pointindex+1]' is flagged.31893190inline void tetgenmesh::pinfect(point pt) {3191((int *) (pt))[pointmarkindex + 1] |= (int) 1;3192}31933194inline void tetgenmesh::puninfect(point pt) {3195((int *) (pt))[pointmarkindex + 1] &= ~(int) 1;3196}31973198inline bool tetgenmesh::pinfected(point pt) {3199return (((int *) (pt))[pointmarkindex + 1] & (int) 1) != 0;3200}32013202// pmarktest(), punmarktest(), pmarktested() -- more primitives to3203// flag or unflag a point.32043205inline void tetgenmesh::pmarktest(point pt) {3206((int *) (pt))[pointmarkindex + 1] |= (int) 2;3207}32083209inline void tetgenmesh::punmarktest(point pt) {3210((int *) (pt))[pointmarkindex + 1] &= ~(int) 2;3211}32123213inline bool tetgenmesh::pmarktested(point pt) {3214return (((int *) (pt))[pointmarkindex + 1] & (int) 2) != 0;3215}32163217inline void tetgenmesh::pmarktest2(point pt) {3218((int *) (pt))[pointmarkindex + 1] |= (int) 4;3219}32203221inline void tetgenmesh::punmarktest2(point pt) {3222((int *) (pt))[pointmarkindex + 1] &= ~(int) 4;3223}32243225inline bool tetgenmesh::pmarktest2ed(point pt) {3226return (((int *) (pt))[pointmarkindex + 1] & (int) 4) != 0;3227}32283229inline void tetgenmesh::pmarktest3(point pt) {3230((int *) (pt))[pointmarkindex + 1] |= (int) 8;3231}32323233inline void tetgenmesh::punmarktest3(point pt) {3234((int *) (pt))[pointmarkindex + 1] &= ~(int) 8;3235}32363237inline bool tetgenmesh::pmarktest3ed(point pt) {3238return (((int *) (pt))[pointmarkindex + 1] & (int) 8) != 0;3239}32403241// These following primitives set and read a pointer to a tetrahedron3242// a subface/subsegment, a point, or a tet of background mesh.32433244inline tetgenmesh::tetrahedron tetgenmesh::point2tet(point pt) {3245return ((tetrahedron *) (pt))[point2simindex];3246}32473248inline void tetgenmesh::setpoint2tet(point pt, tetrahedron value) {3249((tetrahedron *) (pt))[point2simindex] = value;3250}32513252inline tetgenmesh::point tetgenmesh::point2ppt(point pt) {3253return (point) ((tetrahedron *) (pt))[point2simindex + 1];3254}32553256inline void tetgenmesh::setpoint2ppt(point pt, point value) {3257((tetrahedron *) (pt))[point2simindex + 1] = (tetrahedron) value;3258}32593260inline tetgenmesh::shellface tetgenmesh::point2sh(point pt) {3261return (shellface) ((tetrahedron *) (pt))[point2simindex + 2];3262}32633264inline void tetgenmesh::setpoint2sh(point pt, shellface value) {3265((tetrahedron *) (pt))[point2simindex + 2] = (tetrahedron) value;3266}326732683269inline tetgenmesh::tetrahedron tetgenmesh::point2bgmtet(point pt) {3270return ((tetrahedron *) (pt))[point2simindex + 3];3271}32723273inline void tetgenmesh::setpoint2bgmtet(point pt, tetrahedron value) {3274((tetrahedron *) (pt))[point2simindex + 3] = value;3275}327632773278// The primitives for saving and getting the insertion radius.3279inline void tetgenmesh::setpointinsradius(point pt, REAL value)3280{3281pt[pointinsradiusindex] = value;3282}32833284inline REAL tetgenmesh::getpointinsradius(point pt)3285{3286return pt[pointinsradiusindex];3287}32883289inline bool tetgenmesh::issteinerpoint(point pt) {3290return (pointtype(pt) == FREESEGVERTEX) || (pointtype(pt) == FREEFACETVERTEX)3291|| (pointtype(pt) == FREEVOLVERTEX);3292}32933294// point2tetorg() Get the tetrahedron whose origin is the point.32953296inline void tetgenmesh::point2tetorg(point pa, triface& searchtet)3297{3298decode(point2tet(pa), searchtet);3299if ((point) searchtet.tet[4] == pa) {3300searchtet.ver = 11;3301} else if ((point) searchtet.tet[5] == pa) {3302searchtet.ver = 3;3303} else if ((point) searchtet.tet[6] == pa) {3304searchtet.ver = 7;3305} else {3306searchtet.ver = 0;3307}3308}33093310// point2shorg() Get the subface/segment whose origin is the point.33113312inline void tetgenmesh::point2shorg(point pa, face& searchsh)3313{3314sdecode(point2sh(pa), searchsh);3315if ((point) searchsh.sh[3] == pa) {3316searchsh.shver = 0;3317} else if ((point) searchsh.sh[4] == pa) {3318searchsh.shver = (searchsh.sh[5] != NULL ? 2 : 1);3319} else {3320searchsh.shver = 4;3321}3322}33233324// farsorg() Return the origin of the subsegment.3325// farsdest() Return the destination of the subsegment.33263327inline tetgenmesh::point tetgenmesh::farsorg(face& s)3328{3329face travesh, neighsh;33303331travesh = s;3332while (1) {3333senext2(travesh, neighsh);3334spivotself(neighsh);3335if (neighsh.sh == NULL) break;3336if (sorg(neighsh) != sorg(travesh)) sesymself(neighsh);3337senext2(neighsh, travesh);3338}3339return sorg(travesh);3340}33413342inline tetgenmesh::point tetgenmesh::farsdest(face& s)3343{3344face travesh, neighsh;33453346travesh = s;3347while (1) {3348senext(travesh, neighsh);3349spivotself(neighsh);3350if (neighsh.sh == NULL) break;3351if (sdest(neighsh) != sdest(travesh)) sesymself(neighsh);3352senext(neighsh, travesh);3353}3354return sdest(travesh);3355}33563357///////////////////////////////////////////////////////////////////////////////3358// //3359// Linear algebra operators. //3360// //3361///////////////////////////////////////////////////////////////////////////////33623363// dot() returns the dot product: v1 dot v2.3364inline REAL tetgenmesh::dot(REAL* v1, REAL* v2)3365{3366return v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2];3367}33683369// cross() computes the cross product: n = v1 cross v2.3370inline void tetgenmesh::cross(REAL* v1, REAL* v2, REAL* n)3371{3372n[0] = v1[1] * v2[2] - v2[1] * v1[2];3373n[1] = -(v1[0] * v2[2] - v2[0] * v1[2]);3374n[2] = v1[0] * v2[1] - v2[0] * v1[1];3375}33763377// distance() computes the Euclidean distance between two points.3378inline REAL tetgenmesh::distance(REAL* p1, REAL* p2)3379{3380return sqrt((p2[0] - p1[0]) * (p2[0] - p1[0]) +3381(p2[1] - p1[1]) * (p2[1] - p1[1]) +3382(p2[2] - p1[2]) * (p2[2] - p1[2]));3383}33843385inline REAL tetgenmesh::norm2(REAL x, REAL y, REAL z)3386{3387return (x) * (x) + (y) * (y) + (z) * (z);3388}338933903391#endif // #ifndef tetgenH3392339333943395