Meshing routinesο
plane2mesh()
and volume2mesh()
are high-level functions for
generating a mesh. For finer control, use the classes Mesher2D
for
image data and Mesher3D
for volume data.
Both classes derive from Mesher
. Initiating a Mesher
instance
will create the appropriate meshing class.
Classes
|
Utility class to generate a mesh from image data. |
|
Utility class to generate a triangular mesh from 2D image data. |
|
Utility class to generate a tetrahedral mesh from 3D (volumetric) image data. |
Functions
|
Generate a triangular mesh from a 2D segmented image. |
|
Generate a tetrahedral mesh from a 3D segmented image. |
Referenceο
- nanomesh.plane2mesh(image: np.ndarray | Plane, *, level: float = None, max_edge_dist: int = 5, opts: str = 'q30a10', plot: bool = False) MeshContainer [source]ο
Generate a triangular mesh from a 2D segmented image.
- Parameters
image ((i,j) numpy.ndarray or Plane) β Input image to mesh.
level (float, optional) β Level to generate contours at from image
max_edge_dist (int, optional) β Maximum distance between neighbouring pixels in contours.
opts (str, optional) β Options passed to
triangulate()
. For more info, see: https://rufat.be/triangle/API.html#triangle.triangulate
- Returns
mesh β Triangulated 2D mesh with domain labels.
- Return type
- nanomesh.volume2mesh(image: np.ndarray | Volume, *, level: float = None, **kwargs) MeshContainer [source]ο
Generate a tetrahedral mesh from a 3D segmented image.
- Parameters
image ((i,j,k) numpy.ndarray or Volume) β Input image to mesh.
level (float, optional) β Contour value to search for isosurfaces (i.e. the threshold value). By default takes the average of the min and max value. Can be ignored if a binary image is passed as image.
**kwargs β Optional keyword arguments passed to
tetrahedralize()
- Returns
volume_mesh β Instance of
MeshContainer
- Return type
- class nanomesh.Mesher(image: Union[ndarray, Image])[source]ο
Bases:
object
Utility class to generate a mesh from image data.
Depending on the number of the image data, the appropriate subclass will be chosen if possible.
- Parameters
image (np.array) β N-dimensional numpy array containing image data.
- imageο
Reference to image data
- Type
- image_origο
Keep reference to original image data
- Type
Methods:
- class nanomesh.Mesher2D(image: Union[ndarray, Image])[source]ο
Bases:
Mesher
Utility class to generate a triangular mesh from 2D image data.
Depending on the number of the image data, the appropriate subclass will be chosen if possible.
- Parameters
image (np.array) β N-dimensional numpy array containing image data.
- imageο
Reference to image data
- Type
- image_origο
Keep reference to original image data
- Type
Attributes:
Return bbox attribute.
Return bbox from image shape.
Methods:
generate_contour
([level,Β precision,Β ...])Generate contours using marching cubes algorithm.
pad_contour
(**kwargs)Pad the contour using
image2mesh._mesher2d.pad()
.plot_contour
([ax,Β cmap])Plot contours on image.
triangulate
([opts])Triangulate contours.
- property bbox: ndarrayο
Return bbox attribute.
If not explicity set, returns
Mesher2D.image_bbox
.- Sequence:
x0, y0 x1, y0 x1, y1 x0, y0
- Returns
bbox β Bounding box set for output mesh.
- Return type
np.ndarray
- generate_contour(level: Optional[float] = None, precision: int = 1, max_edge_dist: int = 5, group_regions: bool = True)[source]ο
Generate contours using marching cubes algorithm.
Contours are approximated by a polygon, where the maximum distance between points is decided by max_edge_dist.
- Parameters
level (float, optional) β Contour value to search for isosurfaces (i.e. the threshold value). By default takes the average of the min and max value. Can be ignored if a binary image is passed to
Mesher2D
.precision (int, optional) β Maximum distance from original points in polygon approximation routine.
max_edge_dist (int, optional) β Divide long edges so that maximum distance between points does not exceed this value.
group_regions (bool, optional) β If True, assign the same label to all features If False, label regions sequentially
- property image_bbox: ndarrayο
Return bbox from image shape.
- Returns
bbox β Coordinates of bounding box contour.
- Return type
(4,2) np.array
- pad_contour(**kwargs)[source]ο
Pad the contour using
image2mesh._mesher2d.pad()
.- Parameters
mesh (LineMesh) β The mesh to pad.
side (str) β Side to pad, must be one of left, right, top, bottom.
width (int) β Width of the padded area.
label (int, optional) β The label to assign to the padded area. If not defined, generates the next unique label based on the existing ones.
name (str, optional) β Name of the added region. Note that in case of conflicts, the label takes presedence over the name.
- Returns
new_mesh β Padded line mesh.
- Return type
- Raises
ValueError β When the value of side is invalid.
- plot_contour(ax: Optional[Axes] = None, cmap: Optional[str] = None, **kwargs)[source]ο
Plot contours on image.
- Parameters
ax (matplotlib.axes.Axes) β Axes to use for plotting.
cmap (str) β Matplotlib color map for
matplotlib.pyplot.imshow()
**kwargs β These parameters are passed to
plotting.linemeshplot()
- Returns
ax
- Return type
- triangulate(opts='pAq30a10', **kwargs) MeshContainer [source]ο
Triangulate contours.
- Mandatory switches for opts:
e: ensure edges get returned
p: planar straight line graph
A: assign regional attribute to each triangle
If missing, these will be added.
- Parameters
opts (str, optional) β Options passed to
triangulate()
. For more info, see: https://rufat.be/triangle/API.html#triangle.triangulate- Returns
mesh β Triangulated 2D mesh with domain labels
- Return type
- class nanomesh.Mesher3D(image: Union[ndarray, Image])[source]ο
Bases:
Mesher
Utility class to generate a tetrahedral mesh from 3D (volumetric) image data.
Depending on the number of the image data, the appropriate subclass will be chosen if possible.
- Parameters
image (np.array) β N-dimensional numpy array containing image data.
- imageο
Reference to image data
- Type
- image_origο
Keep reference to original image data
- Type
Methods:
generate_contour
([level,Β group_regions])Generate contours using marching cubes algorithm (
skimage.measure.marching_cubes()
).pad_contour
(**kwargs)Pad the contour using
image2mesh._mesher3d.pad()
.plot_contour
(**kwargs)Wrapper for
pyvista.plot()
.set_region_markers
(region_markers)Sets custom region markers for tetrahedralization.
tetrahedralize
(**kwargs)Tetrahedralize surface contour mesh using:func:tetrahedralize.
- generate_contour(level: Optional[float] = None, group_regions: bool = True)[source]ο
Generate contours using marching cubes algorithm (
skimage.measure.marching_cubes()
).Also generates an envelope around the entire data volume corresponding to the bounding box. The bounding box equals the dimensions of the data volume.
By default, the 0-value after segmentation will map to the βbackgroundβ, and the 1-value to feature.
- Parameters
level (float, optional) β Contour value to search for isosurfaces (i.e. the threshold value). By default takes the average of the min and max value. Can be ignored if a binary image is passed to
Mesher3D
.group_regions (bool, optional) β If True, assign the same label to all features If False, label regions sequentially
- pad_contour(**kwargs)[source]ο
Pad the contour using
image2mesh._mesher3d.pad()
.Note that tetgen will assign a different label for physically separate regions, even when they are given the same label/name.
- Parameters
mesh (TriangleMesh) β The mesh to pad.
side (str) β Side to pad, must be one of left, right, top, bottom, back, front.
width (int) β Width of the padded area.
label (int, optional) β The label to assign to the padded area. If not defined, generates the next unique label based on the existing ones.
name (str, optional) β Name of the added region. Note that in case of conflicts, the label takes presedence over the name.
- Returns
new_mesh β Padded contour triangle mesh.
- Return type
- Raises
ValueError β When the value of side is invalid.
- plot_contour(**kwargs)[source]ο
Wrapper for
pyvista.plot()
.- Parameters
**kwargs β These parameters are passed to
pyvista.plot()
- set_region_markers(region_markers: RegionMarkerList)[source]ο
Sets custom region markers for tetrahedralization.
This overwrites any contours generated by .generate_contour().
- Parameters
region_markers (RegionMarkerList) β List of RegionMarker objects.
- tetrahedralize(**kwargs) MeshContainer [source]ο
Tetrahedralize surface contour mesh using:func:tetrahedralize.
- Parameters
mesh (TriangleMesh) β Input contour mesh
opts (str, optional) β
Command-line options passed to tetgen.
More info: http://wias-berlin.de/software/tetgen/1.5/doc/manual/manual005.html
Some useful flags:
-A: Assigns attributes to tetrahedra in different regions.
-p: Tetrahedralizes a piecewise linear complex (PLC).
-q: Refines mesh (to improve mesh quality).
-a: Applies a maximum tetrahedron volume constraint.
Can be passed as a raw string, opts=β-pAq1.2β, or dict, `opts=dict(βpβ= True, βAβ= True, βqβ=1.2).
default_opts (dict, optional) β Dictionary with default options. These will be merged with opts.
- Returns
Tetrahedralized mesh.
- Return type