TetGen.jl

TetGen is a mesh generator written in C++. It generates Delaunay tetrahedralizations, constrained Delaunay tetrahedralizations, and quality tetrahedral meshes. TetGen.jl provides a Julia interface to TetGen.

Mesh based API

This API uses instances of types from GeometryBasics.jl to describe input and output of TetGen.

TetGen.tetrahedralizeFunction
tetrahedralize(mesh; ...)
tetrahedralize(mesh, command; marker, holes)

Tetrahedralize a mesh of polygons with optional facet markers. Returns a mesh of tetrahdra.

source
TetGen.tetrahedralizeFunction
tetrahedralize(mesh)
tetrahedralize(mesh, command)

Tetrahedralize a domain described by a mesh of triangles. Returns a mesh of tetrahdra.

source
TetGen.voronoiMethod
voronoi(points)

Create voronoi diagram of point set.

Returns a mesh of triangles.

source

Raw API

This API is closer to TetGen's C++ API in the sense that input and output are described using arrays of integers and floats, without conversion to any other higher level data structure.

TetGen.RawFacetType
struct RawFacet{T}

A complex facet as part to the input to TetGen.

  • polygonlist::Vector{Vector{Int32}}: Polygons given as arrays of indices which point into the pointlist array describing the input points.
  • holelist::Matrix: Array of points given by their coordinates marking polygons describing holes in the facet.
source
TetGen.RawTetGenIOType
mutable struct RawTetGenIO{T}

A structure for transferring data into and out of TetGen's internal representation.

The input of TetGen is either a 3D point set, or a 3D piecewise linear complex (PLC), or a tetrahedral mesh. Depending on the input object and the specified options, the output of TetGen is either a Delaunay (or weighted Delaunay) tetrahedralization, or a constrained (Delaunay) tetrahedralization, or a quality tetrahedral mesh.

A piecewise linear complex (PLC) represents a 3D polyhedral domain with possibly internal boundaries(subdomains). It is introduced in [Miller et al, 1996]. Basically it is a set of "cells", i.e., vertices, edges, polygons, and polyhedra, and the intersection of any two of its cells is the union of other cells of it.

The 'RawTetGenIO' structure is a collection of arrays of data, i.e., points, facets, tetrahedra, and so forth. All data are compatible to the representation in C++ and can be used without copying.

  • pointlist::Matrix: 'pointlist': Array of point coordinates with size(pointlist,1)==3.
  • pointattributelist::Matrix: 'pointattributelist': Array of point attributes. The number of attributes per point is determined by size(pointattributelist,1)
  • pointmtrlist::Matrix: 'pointmtrlist': An array of metric tensors at points.
  • pointmarkerlist::Vector{Int32}: 'pointmarkerlist': An array of point markers; one integer per point.
  • tetrahedronlist::Matrix{Int32}: 'tetrahedronlist': An array of tetrahedron corners represented by indices of points in pointlist. Unless option -o2 is given, one has size(tetrahedronlist,1)==4, i.e. each column describes the four corners of a tetrahedron. Otherwise, size(tetrahedronlist,1)==10 and the 4 corners are followed by 6 edge midpoints.
  • tetrahedronattributelist::Matrix: 'tetrahedronattributelist': An array of tetrahedron attributes.
  • tetrahedronvolumelist::Vector: 'tetrahedronvolumelist': An array of constraints, i.e. tetrahedron's volume; Input only. This can be used for triggering local refinement.
  • neighborlist::Matrix{Int32}: 'neighborlist': An array of tetrahedron neighbors; 4 ints per element. Output only.
  • facetlist::Array{RawFacet{T}, 1} where T: 'facetlist': An array of facets. Each entry is a structure of facet.
  • facetmarkerlist::Vector{Int32}: 'facetmarkerlist': An array of facet markers; one int per facet.
  • holelist::Matrix: 'holelist': An array of holes (in volume). Each hole is given by a point which lies strictly inside it.
  • regionlist::Matrix: 'regionlist': An array of regions (subdomains). Each region is given by a seed (point) which lies strictly inside it. For each column, the point coordinates ade at indices [1], [2] and [3], followed by the regional attribute at index [4], followed by the maximum volume at index [5].

    ote that each regional attribute is used only if you select the 'A' switch, and each volume constraint is used only if you select the 'a' switch (with no number following).

  • facetconstraintlist::Matrix: 'facetconstraintlist': An array of facet constraints. Each constraint specifies a maximum area bound on the subfaces of that facet. Each column contains a facet marker at index [1] and its maximum area bound at index [2]. Note: the facet marker is actually an integer.
  • segmentconstraintlist::Matrix: 'segmentconstraintlist': An array of segment constraints. Each constraint specifies a maximum length bound on the subsegments of that segment. Each columb consists of the two endpoints of the segment at index [1] and [2], and the maximum length bound at index [3]. Note the segment endpoints are actually integers.
  • trifacelist::Matrix{Int32}: 'trifacelist': An array of face (triangle) corners.
  • trifacemarkerlist::Vector{Int32}: 'trifacemarkerlist': An array of face markers; one int per face.
  • edgelist::Matrix{Int32}: 'edgelist': An array of edge endpoints.
  • edgemarkerlist::Array{Int32}: 'edgemarkerlist': An array of edge markers.
source
TetGen.facetlist!Method
facetlist!(
    tio::RawTetGenIO{T},
    facets::AbstractMatrix
) -> RawTetGenIO
Set list of input facets from AbstractMatrix desribing polygons of the same
size (e.g. triangles)
source
TetGen.facetlist!Method
facetlist!(
    tio::RawTetGenIO{T},
    facets::Vector
) -> RawTetGenIO
Set list of input facets from a vector of polygons of different size
source
TetGen.surfacemeshMethod
surfacemesh(tgio::RawTetGenIO) -> GeometryBasics.Mesh

Create GeometryBasics.Mesh from the triface list (for quick visualization purposes using Makie's wireframe).

source
TetGen.tetrahedralizeMethod
tetrahedralize(
    input::RawTetGenIO{Float64},
    flags::String
) -> RawTetGenIO{Float64}

Tetrahedralize input.

  flags: -pYrq_Aa_miO_S_T_XMwcdzfenvgkJBNEFICQVh 
    -p  Tetrahedralizes a piecewise linear complex (PLC).
    -Y  Preserves the input surface mesh (does not modify it).
    -r  Reconstructs a previously generated mesh.
    -q  Refines mesh (to improve mesh quality).
    -R  Mesh coarsening (to reduce the mesh elements).
    -A  Assigns attributes to tetrahedra in different regions.
    -a  Applies a maximum tetrahedron volume constraint.
    -m  Applies a mesh sizing function.
    -i  Inserts a list of additional points.
    -O  Specifies the level of mesh optimization.
    -S  Specifies maximum number of added points.
    -T  Sets a tolerance for coplanar test (default 1e-8).
    -X  Suppresses use of exact arithmetic.
    -M  No merge of coplanar facets or very close vertices.
    -w  Generates weighted Delaunay (regular) triangulation.
    -c  Retains the convex hull of the PLC.
    -d  Detects self-intersections of facets of the PLC.
    -z  Numbers all output items starting from zero.
    -f  Outputs all faces to .face file.
    -e  Outputs all edges to .edge file.
    -n  Outputs tetrahedra neighbors to .neigh file.
    -v  Outputs Voronoi diagram to files.
    -g  Outputs mesh to .mesh file for viewing by Medit.
    -k  Outputs mesh to .vtk file for viewing by Paraview.
    -J  No jettison of unused vertices from output .node file.
    -B  Suppresses output of boundary information.
    -N  Suppresses output of .node file.
    -E  Suppresses output of .ele file.
    -F  Suppresses output of .face and .edge file.
    -I  Suppresses mesh iteration numbers.
    -C  Checks the consistency of the final mesh.
    -Q  Quiet:  No terminal output except errors.
    -V  Verbose:  Detailed information, more terminal output.
    -h  Help:  A brief instruction for using TetGen.
source
TetGen.volumemeshMethod
volumemesh(
    tgio::RawTetGenIO
) -> GeometryBasics.Mesh{_A, T, GeometryBasics.Ngon{_A, T, N, Point}, V} where {_A, T<:Real, N, Point<:GeometryBasics.AbstractPoint{_A, T}, V<:(GeometryBasics.FaceView{GeometryBasics.Ngon{_A, T, N, Point}, _A1, GeometryBasics.TriangleFace{Int32}, Vector{_A1}, Vector{GeometryBasics.TriangleFace{Int32}}} where _A1<:GeometryBasics.AbstractPoint)}

Create GeometryBasics.Mesh of all tetrahedron faces (for quick visualization purposes using Makie's wireframe).

source
TetGen.tetunsuitable!Function
tetunsuitable!(unsuitable; check_signature)

Set tetunsuitable function called from C wrapper. Setting this function is valid only for one subsequent call to tetrahedralize. The function to be passed has the signature

unsuitable(p1::Vector{Float64},p2::Vector{Float64},p3::Vector{Float64},p4::Vector{Float64})::Bool

where p1...p4 are 3-Vectors describing the corners of a tetrahedron, and the return value is true if its volume is seen as too large.

source

TetGen C++ code