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 | numpy.ndarray) – Either a pyvista surface mesh or a n x 3 vertex array and n x 3 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, attr, triface_markers = 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, attr, triface_markers = tgen.tetrahedralize()
add_hole(point_in_hole: tuple[float, float, float])#

Add a hole to the mesh.

Parameters:

point_in_hole (tuple, list, np.array of float) – A single point inside the hole, specified as (x, y, z).

Examples

Create a sphere as a hole in a cube in PyVista

>>> import pyvista as pv
>>> import tetgen
>>> cube = pv.Cube().triangulate()
>>> sphere = pv.Sphere(theta_resolution=16, phi_resolution=16, radius=0.25)
>>> mesh = pv.merge([sphere, cube])
>>> tgen = tetgen.TetGen(mesh)
>>> tgen.add_hole([0.0, 0.0, 0.0])
>>> nodes, elem = tgen.tetrahedralize(switches="pzq1.4")
>>> grid = tgen.grid
>>> grid
UnstructuredGrid (0x28233f6d420)
N Cells:    3533
N Points:   781
X Bounds:   -5.000e-01, 5.000e-01
Y Bounds:   -5.000e-01, 5.000e-01
Z Bounds:   -5.000e-01, 5.000e-01
N Arrays:   0
>>> grid.slice(normal="z").plot(show_edges=True, cpos="xy")
add_region(id: int, point_in_region: tuple[float, float, float], max_vol: float = 0.0)#

Add a region to the mesh.

Parameters:
  • id (int) – Unique identifier for the region.

  • point_in_region (tuple, list, np.array of float) – A single point inside the region, specified as (x, y, z).

  • max_vol (float, default: 0.0) – Maximum volume for the region.

Examples

Create a sphere in a cube in PyVista and mesh the region in the sphere at a higher density.

>>> import pyvista as pv
>>> import tetgen
>>> cube = pv.Cube().triangulate()
>>> sphere = pv.Sphere(theta_resolution=16, phi_resolution=16, radius=0.25)
>>> mesh = pv.merge([sphere, cube])
>>> tgen = tetgen.TetGen(mesh)
>>> tgen.add_region(1, (0.0, 0.0, 0.0), 8e-6)  # sphere
>>> tgen.add_region(2, [0.99, 0.0, 0.0], 4e-4)  # cube
>>> nodes, elem, attrib = tgen.tetrahedralize(switches="pzq1.4Aa")
>>> grid = tgen.grid
>>> grid
UnstructuredGrid (0x7cc05412bb80)
  N Cells:    23768
  N Points:   3964
  X Bounds:   -5.000e-01, 5.000e-01
  Y Bounds:   -5.000e-01, 5.000e-01
  Z Bounds:   -5.000e-01, 5.000e-01
  N Arrays:   1

Supply the region info to the grid and then plot a slice through the grid.

>>> grid["regions"] = attrib.ravel()
>>> grid.slice().plot(show_edges=True, cpos="zy")
property edge_markers: ndarray[tuple[Any, ...], dtype[int32]]#

Return the (n, ) array of edge markers denoting if an edge is internal.

Interior edges are 0 and exterior edges are -1.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a sphere and create a mask of internal edges, is_internal.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> is_internal = tet.edge_markers == 0
>>> is_internal[:10]
array([ True,  True, False, False,  True, False,  True,  True, False,
       False])
property edges: ndarray[tuple[Any, ...], dtype[int32]]#

Return the (n, 2) array of edges composing the tetrahedralized grid.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a sphere and the first 10 edges composing the tetrahedralized grid.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> tet.edges[:10]
array([[46, 47],
       [46, 40],
       [38, 40],
       [38, 47],
       [38, 46],
       [40, 47],
       [11,  4],
       [11,  2],
       [ 3,  2],
       [ 3,  4]], dtype=int32)
property elem: ndarray[tuple[Any, ...], dtype[float64]]#

Return an (n, 4) or (n, 10) array of elements composing the grid.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a sphere and return the first five elements.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> tet.elem[:5]
array([[ 81,   9,   1, 107],
       [ 58,  50,   0,  98],
       [  0, 104,  58,  98],
       [ 98, 105,  90,  87],
       [ 82,  91,  92, 102]], dtype=int32)
property face2tet: ndarray[tuple[Any, ...], dtype[_ScalarT]]#

Return the mapping between each triface and the tetrahedral elements.

property grid: UnstructuredGrid#

Return a pyvista.UnstructuredGrid of the tetrahedralized surface.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a pyvista.PolyData surface mesh into a pyvista.UnstructuredGrid.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> grid = tet.grid
>>> grid
UnstructuredGrid (...)
  N Cells:    367
  N Points:   110
  X Bounds:   -4.924e-01, 4.924e-01
  Y Bounds:   -4.683e-01, 4.683e-01
  Z Bounds:   -5.000e-01, 5.000e-01
  N Arrays:   0
make_manifold(verbose: bool = False) None#

Reconstruct a manifold clean surface from input mesh.

Updates mesh in-place.

Requires pymeshfix.

Parameters:

verbose (bool, default: False) – Enable verbose output.

Examples

Create a mesh and ensure it’s manfold.

>>> import pyvista
>>> import tetgen
>>> sphere = pyvista.Sphere(theta_resolution=10, phi_resolution=10)
>>> tgen = tetgen.TetGen(sphere)
>>> tgen.make_manifold()
property mesh: PolyData#

Return the input surface mesh.

Returns:

Input surface mesh.

Return type:

pyvista.PolyData

Examples

Generate a tetgen.TetGen and return a pyvista.PolyData.

>>> import pyvista
>>> import tetgen
>>> sphere = pyvista.Sphere(theta_resolution=10, phi_resolution=10)
>>> tgen = tetgen.TetGen(sphere)
>>> tgen.mesh
PolyData (0x7fa3c97138e0)
  N Cells:    160
  N Points:   82
  N Strips:   0
  X Bounds:   -4.924e-01, 4.924e-01
  Y Bounds:   -4.683e-01, 4.683e-01
  Z Bounds:   -5.000e-01, 5.000e-01
  N Arrays:   0
property node: ndarray[tuple[Any, ...], dtype[float64]]#

Return an (n, 3) array of nodes composing the tetrahedralized surface.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a sphere and return the first three nodes.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> tet.node[:3]
array([[ 0.        ,  0.        ,  0.5       ],
       [ 0.        ,  0.        , -0.5       ],
       [ 0.17101008,  0.        ,  0.46984631]])
plot(**kwargs: Any) Any#

Display the input mesh.

See pyvista.plot() for available arguments.

Examples

Plot the input mesh.

>>> import pyvista
>>> import tetgen
>>> sphere = pyvista.Sphere(theta_resolution=10, phi_resolution=10)
>>> tgen = tetgen.TetGen(sphere)
>>> tgen.plot()
tetrahedralize(plc: bool = True, psc: bool = False, refine: bool = False, quality: bool = True, nobisect: bool = False, cdt: bool = False, coarsen: bool = False, weighted: bool = False, brio_hilbert: bool = True, flipinsert: bool = False, metric: bool = False, varvolume: bool = False, fixedvolume: bool = False, regionattrib: bool = False, insertaddpoints: bool = False, diagnose: bool = False, convex: bool = False, nomergefacet: bool = False, nomergevertex: bool = False, noexact: bool = False, nostaticfilter: bool = False, zeroindex: bool = False, facesout: bool = False, edgesout: bool = False, neighout: bool = False, voroout: bool = False, meditview: bool = False, vtkview: bool = False, vtksurfview: bool = False, nobound: bool = False, nonodewritten: bool = False, noelewritten: bool = False, nofacewritten: bool = False, noiterationnum: bool = False, nojettison: bool = False, docheck: bool = False, quiet: bool = False, nowarning: bool = False, verbose: bool = 0.0, cdtrefine: float = 7.0, vertexperblock: int = 4092, tetrahedraperblock: int = 8188, shellfaceperblock: int = 2044, supsteiner_level: int = 2, addsteiner_algo: int = 1, coarsen_param: int = 0, weighted_param: int = 0, fliplinklevel: int = -1, flipstarsize: int = -1, fliplinklevelinc: int = 1, opt_max_flip_level: int = 3, opt_scheme: int = 7, opt_iterations: int = 3, smooth_cirterion: int = 1, smooth_maxiter: int = 7, delmaxfliplevel: int = 1, order: int = 1, reversetetori: int = 0, steinerleft: int = 100000, unflip_queue_limit: int = 1000, no_sort: bool = False, hilbert_order: int = 52, hilbert_limit: int = 8, brio_threshold: int = 64, brio_ratio: float = 0.125, epsilon: float = 1e-08, facet_separate_ang_tol: float = 179.9, collinear_ang_tol: float = 179.9, facet_small_ang_tol: float = 15.0, maxvolume: float = -1.0, maxvolume_length: float = -1.0, minratio: float = 2.0, opt_max_asp_ratio: float = 1000.0, opt_max_edge_ratio: float = 100.0, mindihedral: float = 0.0, optmaxdihedral: float = 177.0, metric_scale: float = 1.0, smooth_alpha: float = 0.3, coarsen_percent: float = 1.0, elem_growth_ratio: float = 0.0, refine_progress_ratio: float = 0.333, switches: str | None = None, bgmeshfilename: str = '', bgmesh: UnstructuredGrid | None = None)#

Generate tetrahedrals interior to the surface mesh.

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:
  • quality (bool, optional) – Enables/disables mesh improvement. Enabled by default. Disable this to speed up mesh generation while sacrificing quality. Default True.

  • minratio (double, default: 2.0) –

    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.

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

  • mindihedral (double, default: 0.0) –

    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.

    Testing has shown that 10.0 is a reasonable input.

  • verbose (bool, default: False) – 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 is no output

  • nobisect (bool, default: False) –

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

    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, default: 100000) –

    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.

    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.

  • order (int, default: 1) – Controls whether TetGen creates linear tetrahedrals or quadradic tetrahedrals. Set order to 2 to output quadradic tetrahedrals.

  • bgmeshfilename (str, optional) – Filename of the background mesh.

  • regionattrib (bool, default: False) – Return region attributes.

  • bgmesh (pv.UnstructuredGrid) – Background mesh to be processed. Must be composed of only linear tetra. Cannot specify both bgmeshfilename and bgmesh.

Returns:

  • nodes (np.ndarray[np.float64]) – Array of nodes representing the tetrahedral mesh.

  • elems (np.ndarray[np.int32]) – Array of elements representing the tetrahedral mesh.

  • attr (np.ndarray[np.int64] | None) – Region attributes. None unless regionattrib=True or the tetgen flag A is passed.

  • triface_markers (np.ndarray[np.int32]) – Marker for each face in TetGen.triface_list.

Examples

The following switches “pq1.1/10Y” would be:

>>> nodes, elems = tgen.tetrahedralize(
...     nobisect=True, quality=True, minratio=1.1, mindihedral=10
... )

Using the switches option:

>>> nodes, elems = tgen.tetrahedralize(switches="pq1.1/10Y")

Notes

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

Switches of TetGen:

Option

Switch

Default

plc

'-p'

False

psc

'-s'

False

refine

'-r'

False

quality

'-q'

False

nobisect

'-Y'

False

cdt

'-D'

False

coarsen

'-R'

False

weighted

'-w'

False

brio_hilbert

'-b'

True

flipinsert

'-L'

False

metric

'-m'

False

varvolume

'-a'

False

fixedvolume

'-a'

False

regionattrib

'-A'

False

insertaddpoints

'-i'

False

diagnose

'-d'

False

convex

'-c'

False

nomergefacet

'-M'

False

nomergevertex

'-M'

False

noexact

'-X'

False

nostaticfilter

'-X'

False

zeroindex

'-z'

False

facesout

'-f'

False

edgesout

'-e'

False

neighout

'-n'

False

voroout

'-v'

False

meditview

'-g'

False

vtkview

'-k'

False

vtksurfview

'-k'

False

nobound

'-B'

False

nonodewritten

'-N'

False

noelewritten

'-E'

False

nofacewritten

'-F'

False

noiterationnum

'-I'

False

nojettison

'-J'

False

docheck

'-C'

False

quiet

'-Q'

False

nowarning

'-W'

False

verbose

'-V'

False

Parameters of TetGen:

Option

Switch

Default

cdtrefine

'-D#'

7

vertexperblock

'-x'

4092

tetrahedraperblock

'-x'

8188

shellfaceperblock

'-x'

2044

supsteiner_level

'-Y/'

2

addsteiner_algo

'-Y//'

1

coarsen_param

'-R'

0

weighted_param

'-w'

0

opt_max_flip_level

'-O'

3

opt_scheme

'-O/#'

7

opt_iterations

'-O//#'

3

smooth_cirterion

'-s'

1

smooth_maxiter

'-s'

7

order

'-o'

1

reversetetori

'-o/'

0

steinerleft

'-S'

100000

unflip_queue_limit

'-U#'

1000

hilbert_order

'-b///'

52

hilbert_limit

'-b//'

8

brio_threshold

'-b'

64

brio_ratio

'-b/'

0.125.

epsilon

'-T'

1.0e-8.

facet_separate_ang_tol

'-p'

179.9.

collinear_ang_tol

'-p/'

179.9.

facet_small_ang_tol

'-p//'

15.0.

maxvolume

'-a'

-1.0.

maxvolume_length

'-a'

-1.0.

minratio

'-q'

0.0.

mindihedral

'-q'

5.0.

optmaxdihedral

'-o/#'

177.0.

metric_scale

'-m#'

1.0.

smooth_alpha

'-s'

0.3.

coarsen_percent

'-R1/#'

1.0.

elem_growth_ratio

'-r#'

0.0.

refine_progress_ratio

'-r/#'

0.333.

property triface_markers: ndarray[tuple[Any, ...], dtype[int32]]#

Return the (n, ) array of markers for each triangular face in TetGen.trifaces.

Interior faces are denoted with 0 and exterior faces are marked as -1.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a sphere and return the array of the triangular faces that compose the tetrahedralized grid.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> is_interior = tet.triface_markers == 0
>>> is_interior[:10]
array([ True,  True,  True, False,  True,  True,  True, False,  True,
        True])
property trifaces: ndarray[tuple[Any, ...], dtype[int32]]#

Return the (n, 3) array of triangular faces composing the tetrahedral mesh.

The indices of these faces correspond to a node in TetGen.node.

This attribute is only available after running TetGen.tetrahedralize().

Examples

Tetrahedralize a sphere and return the array of the triangular faces that compose the tetrahedralized grid.

>>> import tetgen
>>> import pyvista as pv
>>> sphere = pv.Sphere(theta_resolution=10, phi_resolution=10)
>>> tet = tetgen.TetGen(sphere)
>>> tet.tetrahedralize(switches="pq1.1/10YQ")
>>> tet.trifaces
array([[107,   1,   9],
       [107,  81,   1],
       [  9,  81, 107],
       ...,
       [ 15,   6, 109],
       [ 15,   7,   6],
       [ 15,   6,  14]], shape=(814, 3), dtype=int32)
write(filename: str | Path, binary: bool = True) None#

Write an unstructured grid to disk.

Parameters:
  • filename (str | pathlib.Path) –

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

    • ".vtk" will use the vtk legacy writer

    • ".vtu" will select the VTK XML writer

  • binary (bool, default: True) – Writes as a binary file by default. Set to False to write ASCII.

Examples

Write to a VTK file.

>>> tgen.write("grid.vtk", binary=True)

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.

You can use utilities like meshio to convert to other formats in order to import into FEA software.