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 an x 3
vertex array andn 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 = 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()
- property 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.
- property mesh#
Return the surface mesh
- plot(**kwargs)#
Displays input mesh
See help(
pyvista.plot()
) for available arguments.
- tetrahedralize(plc=True, psc=0.0, refine=0.0, quality=True, nobisect=False, cdt=0.0, cdtrefine=7.0, coarsen=0.0, weighted=0.0, brio_hilbert=1.0, flipinsert=0.0, metric=0.0, varvolume=0.0, fixedvolume=0.0, regionattrib=0.0, insertaddpoints=0.0, diagnose=0.0, convex=0.0, nomergefacet=0.0, nomergevertex=0.0, noexact=0.0, nostaticfilter=0.0, zeroindex=0.0, facesout=0.0, edgesout=0.0, neighout=0.0, voroout=0.0, meditview=0.0, vtkview=0.0, vtksurfview=0.0, nobound=0.0, nonodewritten=0.0, noelewritten=0.0, nofacewritten=0.0, noiterationnum=0.0, nojettison=0.0, docheck=0.0, quiet=0.0, nowarning=0.0, verbose=0.0, vertexperblock=4092.0, tetrahedraperblock=8188.0, shellfaceperblock=2044.0, supsteiner_level=2.0, addsteiner_algo=1.0, coarsen_param=0.0, weighted_param=0.0, fliplinklevel=-1.0, flipstarsize=-1.0, fliplinklevelinc=1.0, opt_max_flip_level=3.0, opt_scheme=7.0, opt_iterations=3.0, smooth_cirterion=1.0, smooth_maxiter=7.0, delmaxfliplevel=1.0, order=1.0, reversetetori=0.0, steinerleft=100000.0, unflip_queue_limit=1000.0, no_sort=0.0, hilbert_order=52.0, hilbert_limit=8.0, brio_threshold=64.0, brio_ratio=0.125, epsilon=1e-08, facet_separate_ang_tol=179.9, collinear_ang_tol=179.9, facet_small_ang_tol=15.0, maxvolume=-1.0, maxvolume_length=-1.0, minratio=2.0, opt_max_asp_ratio=1000.0, opt_max_edge_ratio=100.0, mindihedral=0.0, optmaxdihedral=177.0, metric_scale=1.0, smooth_alpha=0.3, coarsen_percent=1.0, elem_growth_ratio=0.0, refine_progress_ratio=0.333, switches=None)#
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:
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.0Testing 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.0Testing 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.
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)
Using the switches option:
>>> node, elem = 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'
psc
'-s'
refine
'-r'
quality
'-q'
nobisect
'-Y'
cdt
'-D'
cdtrefine
'-D#'
coarsen
'-R'
weighted
'-w'
brio_hilbert
'-b'
flipinsert
'-L'
metric
'-m'
varvolume
'-a'
fixedvolume
'-a'
regionattrib
'-A'
insertaddpoints
'-i'
diagnose
'-d'
convex
'-c'
nomergefacet
'-M'
nomergevertex
'-M'
noexact
'-X'
nostaticfilter
'-X'
zeroindex
'-z'
facesout
'-f'
edgesout
'-e'
neighout
'-n'
voroout
'-v'
meditview
'-g'
vtkview
'-k'
vtksurfview
'-k'
nobound
'-B'
nonodewritten
'-N'
noelewritten
'-E'
nofacewritten
'-F'
noiterationnum
'-I'
nojettison
'-J'
docheck
'-C'
quiet
'-Q'
nowarning
'-W'
verbose
'-V'
Parameters of TetGen.
Option
Switch
Default
vertexperblock
'-x'
tetrahedraperblock
'-x'
shellfaceperblock
'-x'
supsteiner_level
'-Y/'
addsteiner_algo
'-Y//'
coarsen_param
'-R'
weighted_param
'-w'
opt_max_flip_level
'-O'
opt_scheme
'-O/#'
opt_iterations
'-O//#'
smooth_cirterion
'-s'
smooth_maxiter
'-s'
order
'-o'
reversetetori
'-o/'
steinerleft
'-S'
unflip_queue_limit
'-U#'
hilbert_order
'-b///'
hilbert_limit
'-b//'
brio_threshold
'-b'
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.
- 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".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.
Examples
>>> 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.