API Reference

Python module to interface with wrapped TetGen C++ code

class tetgen.TetGen(*args)

Bases: object

Input, clean, and tetrahedralize surface meshes using TetGen

Parameters

args (str, pyvista.PolyData or (np.ndarray, np.ndarray)) – Either a pyvista surface mesh or a nx3 vertex array and nx3 face array.

Examples

Tetrahedralize a sphere using pyvista >>> import pyvista >>> import tetgen >>> sphere = pyvista.Sphere(theta_resolution=10, phi_resolution=10) >>> tgen = tetgen.TetGen(sphere) >>> nodes, elem = tgen.tetrahedralize() >>> tgen.grid.plot(show_edges=True)

Tetrahedralize a cube using numpy arrays >>> import numpy as np >>> import tetgen >>> v = np.array([[0, 0, 0], [1, 0, 0],

[1, 1, 0], [0, 1, 0], [0, 0, 1], [1, 0, 1], [1, 1, 1], [0, 1, 1],])

>>> f = np.vstack([[0, 1, 2], [2, 3, 0],
                   [0, 1, 5], [5, 4, 0],
                   [1, 2, 6], [6, 5, 1],
                   [2, 3, 7], [7, 6, 2],
                   [3, 0, 4], [4, 7, 3],
                   [4, 5, 6], [6, 7, 4]])
>>> tgen = tetgen.TetGen(v, f)
>>> nodes, elems = tgen.tetrahedralize()
grid

Returns a pyvista.UnstructuredGrid

make_manifold(verbose=False)

Reconstruct a manifold clean surface from input mesh. Updates mesh in-place.

Requires pymeshfix

Parameters

verbose (bool, optional) – Controls output printing. Default False.

mesh

Return the surface mesh

plot(**kwargs)

Displays input mesh

See help(pyvista.plot()) for available arguments.

tetrahedralize(psc=0, refine=0, quality=1, nobisect=True, coarsen=0, metric=0, weighted=0, brio_hilbert=1, incrflip=0, flipinsert=0, varvolume=0, fixedvolume=0, noexact=0, nostaticfilter=0, insertaddpoints=0, regionattrib=0, cdtrefine=0, diagnose=0, convex=0, zeroindex=0, facesout=0, edgesout=0, neighout=0, voroout=0, meditview=0, vtkview=0, nobound=0, nonodewritten=1, noelewritten=1, nofacewritten=1, noiterationnum=0, nomergefacet=0, nomergevertex=0, nojettison=0, docheck=0, quiet=0, verbose=0, vertexperblock=4092, tetrahedraperblock=8188, shellfaceperblock=4092, nobisect_nomerge=1, supsteiner_level=2, addsteiner_algo=1, coarsen_param=0, weighted_param=0, fliplinklevel=-1, flipstarsize=-1, fliplinklevelinc=1, reflevel=3, optscheme=7, optlevel=2, delmaxfliplevel=1, order=2, reversetetori=0, steinerleft=10000, no_sort=0, hilbert_order=52, hilbert_limit=8, brio_threshold=64, brio_ratio=0.125, facet_separate_ang_tol=179.9, facet_overlap_ang_tol=0.001, facet_small_ang_tol=15.0, maxvolume=-1.0, minratio=2.0, mindihedral=0.0, optmaxdihedral=165.0, optminsmtdihed=179.0, optminslidihed=179.0, epsilon=1e-08, coarsen_percent=1.0)

Generates tetrahedrals interior to the surface mesh described by the vertex and face arrays already loaded. Returns nodes and elements belonging to the all tetrahedral mesh.

The tetrahedral generator uses the C++ library TetGen and can be configured by either using a string of ‘switches’ or by changing the underlying behavior using optional inputs.

Should the user desire more control over the mesh tetrahedralization or wish to control the tetrahedralization in a more pythonic manner, use the optional inputs rather than inputting switches.

Parameters
  • facet_overlap_ang_tol (double, optional) – Threshold angle at which TetGen will consider to faces overlapping. Raising this will require a higher quality mesh input and may cause tetrahedralize to fail. Default 0.001.

  • quality (bool, optional) – Enables/disables mesh improvement. Enabled by default. Disable this to speed up mesh generation while sacrificing quality. Default True.

  • minratio (double, optional.) –

    Maximum allowable radius-edge ratio. Must be greater than 1.0 the closer to 1.0, the higher the quality of the mesh. Be sure to raise steinerleft to allow for the addition of points to improve the quality of the mesh. Avoid overly restrictive requirements, otherwise, meshing will appear to hang. Default 2.0

    Testing has showed that 1.1 is a reasonable input for a high quality mesh.

  • mindihedral (double, optional) –

    Minimum allowable dihedral angle. The larger this number, the higher the quality of the resulting mesh. Be sure to raise steinerleft to allow for the addition of points to improve the quality of the mesh. Avoid overly restrictive requirements, otherwise, meshing will appear to hang. Default 0.0

    Testing has shown that 10 is a reasonable input

  • verbose (int, optional) – Controls the underlying TetGen library to output text to console. Users using iPython will not see this output. Setting to 1 enables some information about the mesh generation while setting verbose to 2 enables more debug output. Default 0, or no output.

  • nobisect (bool, optional) –

    Controls if Steiner points are added to the input surface mesh. When enabled, the surface mesh will be modified. Default False.

    Testing has shown that if your input surface mesh is already well shaped, disabling this setting will improve meshing speed and mesh quality.

  • steinerleft (int, optional) –

    Steiner points are points added to the original surface mesh to create a valid tetrahedral mesh. Settings this to -1 will allow tetgen to create an unlimited number of steiner points, but the program will likely hang if this is used in combination with narrow quality requirements. Default 100000.

    The first type of Steiner points are used in creating an initial tetrahedralization of PLC. These Steiner points are mandatory in order to create a valid tetrahedralization

    The second type of Steiner points are used in creating quality tetra- hedral meshes of PLCs. These Steiner points are optional, while they may be necessary in order to improve the mesh quality or to conform the size of mesh elements.

  • double (optmaxdihedral, optional) – Setting unreachable using switches. Controls the optimial maximum dihedral. Settings closer, but not exceeding, 180 degrees results in a lower quality mesh. Should be between 135 and 180 degrees. Default 165.0

  • order (int optional) – Controls whether TetGen creates linear tetrahedrals or quadradic tetrahedrals. Set order to 2 to output quadradic tetrahedrals. Default 2.

Examples

The following switches “pq1.1/10Y” would be: >>> node, elem = tgen.Tetrahedralize(nobisect=True, quality=True,

minratio=1.1, mindihedral=10)

Notes

There are many other options and the TetGen documentation contains descritpions only for the switches of the original C++ program. This is the relationship between tetgen switches and python optinal inputs:

PYTHON OPTION TETGEN SWITCH int psc; // -s int refine; // -r int quality; // -q int nobisect; // -Y int coarsen; // -R int weighted; // -w int brio_hilbert; // -b int incrflip; // -l int flipinsert; // -L int metric; // -m int varvolume; // -a int fixedvolume; // -a int regionattrib; // -A int cdtrefine; // -D int insertaddpoints; // -i int diagnose; // -d int convex; // -c int nomergefacet; // -M int nomergevertex; // -M int noexact; // -X int nostaticfilter; // -X int zeroindex; // -z int voroout; // -v int meditview; // -g int vtkview; // -k int nobound; // -B int noiterationnum; // -I int nojettison; // -J int docheck; // -C int quiet; // -Q int verbose; // -V

PYTHON OPTION TETGEN SWITCH int vertexperblock; // ‘-x’, 4092. int tetrahedraperblock; // ‘-x’, 8188. int shellfaceperblock; // ‘-x’, 2044. int nobisect_nomerge; // ‘-Y’, 1. int supsteiner_level; // ‘-Y/’, 2. int addsteiner_algo; // ‘-Y//’, 1. int coarsen_param; // ‘-R’, 0. int weighted_param; // ‘-w’, 0. int fliplinklevel; // -1. int flipstarsize; // -1. int fliplinklevelinc; // 1. int reflevel; // ‘-D’, 3. int optlevel; // ‘-O’, 2. int optscheme; // ‘-O’, 7. int delmaxfliplevel; // 1. int order; // ‘-o’, 1. int reversetetori; // ‘-o/’, 0. int steinerleft; // ‘-S’, 0. int no_sort; // 0. int hilbert_order; // ‘-b///’, 52. int hilbert_limit; // ‘-b//’ 8. int brio_threshold; // ‘-b’ 64. REAL brio_ratio; // ‘-b/’ 0.125. REAL facet_separate_ang_tol; // ‘-p’, 179.9. REAL facet_overlap_ang_tol; // ‘-p/’, 0.1. REAL facet_small_ang_tol; // ‘-p//’, 15.0. REAL maxvolume; // ‘-a’, -1.0. REAL minratio; // ‘-q’, 0.0. REAL mindihedral; // ‘-q’, 5.0. REAL optmaxdihedral; // 165.0. REAL optminsmtdihed; // 179.0. REAL optminslidihed; // 179.0. REAL epsilon; // ‘-T’, 1.0e-8. REAL coarsen_percent; // -R1/#, 1.0.

write(filename, binary=False)

Writes an unstructured grid to disk.

Parameters
  • filename (str) –

    Filename of grid to be written. The file extension will select the type of writer to use.

    • ”.vtk” will use the vtk legacy writer, while

    • ”.vtu” will select the VTK XML writer.

    • ”.cdb” will write an ANSYS APDL archive file.

  • binary (bool, optional) – Writes as a binary file by default. Set to False to write ASCII. Ignored when output is a cdb.

Notes

Binary files write much faster than ASCII, but binary files written on one system may not be readable on other systems. Binary can be used only with the legacy writer.