Mesh#

This module contains meshing-related classes and functions. Standalone mesh-tools (functions) are also available as mesh-methods.

Core

Mesh(points, cells[, cell_type])

A mesh with points, cells and optional a specified cell type.

MeshContainer(meshes[, merge, decimals])

A container which operates on a list of meshes with identical dimensions.

Geometries

mesh.Line([a, b, n])

A line shaped 1d-mesh with lines between a and b with n points.

Rectangle([a, b, n])

A rectangular 2d-mesh with quads between a and b with n points per axis.

Cube([a, b, n])

A cube shaped 3d-mesh with hexahedrons between a and b with n points per axis.

Grid(*xi[, indexing])

A grid shaped 3d-mesh with hexahedrons.

Circle([radius, centerpoint, n, sections, ...])

A circular shaped 2d-mesh with quads and n points on the circumferential edge of a 45-degree section.

mesh.Triangle([a, b, c, n, decimals])

A triangular shaped 2d-mesh with quads and n points at the edges of the three sub-quadrilaterals.

mesh.RectangleArbitraryOrderQuad([a, b, order])

A rectangular 2d-mesh with arbitrary order quads between a and b with n points per axis.

Tools

mesh.expand(points, cells, cell_type[, n, z])

Expand a 1d-Line to a 2d-Quad or a 2d-Quad to a 3d-Hexahedron Mesh.

mesh.rotate(points, cells, cell_type, ...[, ...])

Rotate a Mesh.

mesh.revolve(points, cells, cell_type[, n, ...])

Revolve a 2d-Quad to a 3d-Hexahedron Mesh.

mesh.sweep(points, cells, cell_type[, decimals])

Merge duplicated points and update cells of a Mesh.

mesh.mirror(points, cells, cell_type[, ...])

Mirror points by plane normal and ensure positive cell volumes for tria, tetra, quad and hexahedron cell types.

mesh.concatenate(meshes)

Join a sequence of meshes with identical cell types.

mesh.runouts(points, cells, cell_type[, ...])

Add simple rubber-runouts for realistic rubber-metal structures.

mesh.triangulate(points, cells, cell_type[, ...])

Triangulate a quad or a hex mesh.

mesh.convert(points, cells, cell_type[, ...])

Convert mesh to a given order (only order=0 and order=2 from order=1 are supported).

mesh.collect_edges(points, cells, cell_type)

"Collect all unique edges, calculate and return midpoints on edges as well as the additional cells array.

mesh.collect_faces(points, cells, cell_type)

"Collect all unique faces, calculate and return midpoints on faces as well as the additional cells array.

mesh.collect_volumes(points, cells, cell_type)

"Collect all volumes, calculate and return midpoints on volumes as well as the additional cells array.

mesh.add_midpoints_edges(points, cells, ...)

"Add midpoints on edges for given points and cells and update cell_type accordingly.

mesh.add_midpoints_faces(points, cells, ...)

"Add midpoints on faces for given points and cells and update cell_type accordingly.

mesh.add_midpoints_volumes(points, cells, ...)

"Add midpoints on volumes for given points and cells and update cell_type accordingly.

mesh.flip(points, cells, cell_type[, mask])

Ensure positive cell volumes for tria, tetra, quad and hexahedron cell types.

mesh.fill_between(mesh, other_mesh[, n])

Fill a 2d-Quad Mesh between two 1d-Line Meshes, embedded in 2d-space, or a 3d-Hexahedron Mesh between two 2d-Quad Meshes, embedded in 3d-space, by expansion.

mesh.dual(points, cells, cell_type[, ...])

Create a new dual mesh with given points per cell.

Detailed API Reference

class felupe.Mesh(points, cells, cell_type=None)[source]#

A mesh with points, cells and optional a specified cell type.

Parameters:
  • points (ndarray) – Point coordinates.

  • cells (ndarray) – Point-connectivity of cells.

  • cell_type (str or None, optional) – An optional string in VTK-convention that specifies the cell type (default is None). Necessary when a mesh is saved to a file.

points#

Point coordinates.

Type:

ndarray

cells#

Point-connectivity of cells.

Type:

ndarray

cell_type#

A string in VTK-convention that specifies the cell type.

Type:

str or None

add_midpoints_edges(cells, cell_type, cell_type_new=None)#

“Add midpoints on edges for given points and cells and update cell_type accordingly.

add_midpoints_faces(cells, cell_type, cell_type_new=None)#

“Add midpoints on faces for given points and cells and update cell_type accordingly.

add_midpoints_volumes(cells, cell_type, cell_type_new=None)#

“Add midpoints on volumes for given points and cells and update cell_type accordingly.

add_runouts(cells, cell_type, values=[0.1, 0.1], centerpoint=[0, 0, 0], axis=0, exponent=5, mask=slice(None, None, None), normalize=False)#

Add simple rubber-runouts for realistic rubber-metal structures.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • values (list or ndarray, optional) – Relative amount of runouts (per coordinate) perpendicular to the axis (default is 10% per coordinate, i.e. [0.1, 0.1]).

  • centerpoint (list or ndarray, optional) – Center-point coordinates (default is [0, 0, 0]).

  • axis (int or None, optional) – Axis (default is 0).

  • exponent (int, optional) – Positive exponent to control the shape of the runout. The higher the exponent, the steeper the transition (default is 5).

  • mask (list or None, optional) – List of points to be considered (default is None).

  • normalize (bool, optional) – Normalize the runouts to create indents, i.e. maintain the original shape at the ends (default is False).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

as_meshio(**kwargs)[source]#

Export the mesh as meshio.Mesh.

collect_edges(cells, cell_type)#

“Collect all unique edges, calculate and return midpoints on edges as well as the additional cells array.

collect_faces(cells, cell_type)#

“Collect all unique faces, calculate and return midpoints on faces as well as the additional cells array.

collect_volumes(cells, cell_type)#

“Collect all volumes, calculate and return midpoints on volumes as well as the additional cells array.

convert(cells, cell_type, order=0, calc_points=False, calc_midfaces=False, calc_midvolumes=False)#

Convert mesh to a given order (only order=0 and order=2 from order=1 are supported).

copy(points=None, cells=None, cell_type=None)#

Return a deepcopy.

disconnect(points_per_cell=None, calc_points=True)[source]#

Return a new instance of a Mesh with disconnected cells. Optionally, the points-per-cell may be specified (must be lower or equal the number of points- per-cell of the original Mesh). If the Mesh is to be used as a dual Mesh, then the point-coordinates do not have to be re-created because they are not used.

dual(cells, cell_type, points_per_cell=None, disconnect=True, calc_points=False, offset=0, npoints=None)#

Create a new dual mesh with given points per cell. The point coordinates are not used in a dual mesh and hence, by default they are all zero.

expand(cells, cell_type, n=11, z=1)#

Expand a 1d-Line to a 2d-Quad or a 2d-Quad to a 3d-Hexahedron Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • n (int, optional) – Number of n-point repetitions or (n-1)-cell repetitions, default is 11.

  • z (float or ndarray, optional) – Total expand dimension as float (edge length in expand direction is z / n), default is 1. Optionally, if an array is passed these entries are taken as expansion and n is ignored.

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

fill_between(other_mesh, n=11)#

Fill a 2d-Quad Mesh between two 1d-Line Meshes, embedded in 2d-space, or a 3d-Hexahedron Mesh between two 2d-Quad Meshes, embedded in 3d-space, by expansion. Both meshes must have equal number of points and cells. The cells-array is taken from the first mesh.

Parameters:
  • mesh (felupe.Mesh) – The base line- or quad-mesh.

  • other_mesh (felupe.Mesh) – The other line- or quad-mesh.

  • n (int or ndarray) – Number of n-point repetitions or (n-1)-cell repetitions, (default is 11). If an array is given, then its values are used for the relative positions in a reference configuration (-1, 1) between the two meshes.

Returns:

mesh – The expanded mesh.

Return type:

felupe.Mesh

flip(cells, cell_type, mask=None)#

Ensure positive cell volumes for tria, tetra, quad and hexahedron cell types.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • mask (list or ndarray, optional) – Boolean mask for selected cells to flip.

Returns:

  • points (ndarray) – Point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

mirror(cells, cell_type, normal=[1, 0, 0], centerpoint=[0, 0, 0], axis=None)#

Mirror points by plane normal and ensure positive cell volumes for tria, tetra, quad and hexahedron cell types.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • normal (list or ndarray, optional) – Mirror-plane normal vector (default is [1, 0, 0]).

  • centerpoint (list or ndarray, optional) – Center-point coordinates on the mirror plane (default is [0, 0, 0]).

  • axis (int or None, optional) – Mirror axis (default is None).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

plot(*args, **kwargs)[source]#

Plot the mesh.

See also

felupe.Scene.plot

Plot method of a scene.

revolve(cells, cell_type, n=11, phi=180, axis=0)#

Revolve a 2d-Quad to a 3d-Hexahedron Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • n (int, optional) – Number of n-point revolutions (or (n-1) cell revolutions), default is 11.

  • phi (float or ndarray, optional) – Revolution angle in degree (default is 180).

  • axis (int, optional) – Revolution axis (default is 0).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

rotate(cells, cell_type, angle_deg, axis, center=None)#

Rotate a Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • angle_deg (int) – Rotation angle in degree.

  • axis (int) – Rotation axis.

  • center (list or ndarray or None, optional) – Center point coordinates (default is None).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

save(filename='mesh.vtk', **kwargs)[source]#

Export the mesh as VTK file. For XDMF-export please ensure to have h5py (as an optional dependancy of meshio) installed.

Parameters:

filename (str, optional) – The filename of the mesh (default is mesh.vtk).

screenshot(*args, filename='field.png', transparent_background=None, scale=None, **kwargs)[source]#

Take a screenshot of the mesh.

See also

pyvista.Plotter.screenshot

Take a screenshot of a PyVista plotter.

sweep(cells, cell_type, decimals=None)#

Merge duplicated points and update cells of a Mesh.

WARNING: This function re-sorts points.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • decimals (int or None, optional) – Number of decimals for point coordinate comparison (default is None).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

translate(cells, cell_type, move, axis)#

Translate (move) a Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • move (float) – Translation along given axis.

  • axis (int) – Translation axis.

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

triangulate(cells, cell_type, mode=3)#

Triangulate a quad or a hex mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • mode (int, optional) – Choose a mode how to convert hexahedrons to tets [1] (default is 3).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

References

[1] Dompierre, J., Labbé, P., Vallet, M. G., & Camarero, R. (1999). How to Subdivide Pyramids, Prisms, and Hexahedra into Tetrahedra. IMR, 99, 195.

update(points=None, cells=None, cell_type=None)#

Update the cell and dimension attributes with a given cell array.

view(point_data=None, cell_data=None, cell_type=None)[source]#

View the mesh with optional given dicts of point- and cell-data items.

Parameters:
  • point_data (dict or None, optional) – Additional point-data dict (default is None).

  • cell_data (dict or None, optional) – Additional cell-data dict (default is None).

  • cell_type (pyvista.CellType or None, optional) – Cell-type of PyVista (default is None).

Returns:

A object which provides visualization methods for felupe.Mesh.

Return type:

felupe.ViewMesh

See also

felupe.ViewMesh

Visualization methods for felupe.Mesh.

class felupe.MeshContainer(meshes, merge=False, decimals=None)[source]#

A container which operates on a list of meshes with identical dimensions.

Parameters:

meshes ([felupe.Mesh, ...]) – A list with meshes.

dim#

The (identical) dimension of all underlying meshes.

Type:

int

points#

Point coordinates.

Type:

ndarray

meshes#

A list with meshes.

Type:

[felupe.Mesh, …]

append(mesh)[source]#

Append a Mesh to the list of meshes.

as_meshio(combined=True, **kwargs)[source]#

Export a (combined) mesh object as meshio.Mesh.

cells()[source]#

Return a list of tuples with cell-types and cell-connectivities.

copy()[source]#

Return a deepcopy of the mesh container.

merge_duplicate_points(decimals=None)[source]#

Merge duplicate points and update meshes.

pop(index)[source]#

Pop an item of the list of meshes.

class felupe.mesh.Line(a=0, b=1, n=2)[source]#

Bases: Mesh

A line shaped 1d-mesh with lines between a and b with n points.

class felupe.Rectangle(a=(0, 0), b=(1, 1), n=(2, 2))[source]#

Bases: Mesh

A rectangular 2d-mesh with quads between a and b with n points per axis.

class felupe.Cube(a=(0, 0, 0), b=(1, 1, 1), n=(2, 2, 2))[source]#

Bases: Mesh

A cube shaped 3d-mesh with hexahedrons between a and b with n points per axis.

class felupe.Grid(*xi, indexing='ij', **kwargs)[source]#

Bases: Mesh

A grid shaped 3d-mesh with hexahedrons. Basically a wrapper for numpy.meshgrid() with default indexing="ij".

class felupe.Circle(radius=1, centerpoint=[0, 0], n=2, sections=[0, 90, 180, 270], value=0.15, exponent=2, decimals=10)[source]#

Bases: Mesh

A circular shaped 2d-mesh with quads and n points on the circumferential edge of a 45-degree section. 90-degree sections are placed at given angles in degree.

class felupe.mesh.Triangle(a=(0, 0), b=(1, 0), c=(0, 1), n=2, decimals=10)[source]#

Bases: Mesh

A triangular shaped 2d-mesh with quads and n points at the edges of the three sub-quadrilaterals.

class felupe.mesh.RectangleArbitraryOrderQuad(a=(0, 0), b=(1, 1), order=2)[source]#

Bases: Mesh

A rectangular 2d-mesh with arbitrary order quads between a and b with n points per axis.

class felupe.mesh.CubeArbitraryOrderHexahedron(a=(0, 0, 0), b=(1, 1, 1), order=2)[source]#

Bases: Mesh

A cube shaped 3d-mesh with arbitrary order hexahedrons between a and b with n points per axis.

felupe.mesh.add_midpoints_edges(points, cells, cell_type, cell_type_new=None)[source]#

“Add midpoints on edges for given points and cells and update cell_type accordingly.

felupe.mesh.add_midpoints_faces(points, cells, cell_type, cell_type_new=None)[source]#

“Add midpoints on faces for given points and cells and update cell_type accordingly.

felupe.mesh.add_midpoints_volumes(points, cells, cell_type, cell_type_new=None)[source]#

“Add midpoints on volumes for given points and cells and update cell_type accordingly.

felupe.mesh.collect_edges(points, cells, cell_type)[source]#

“Collect all unique edges, calculate and return midpoints on edges as well as the additional cells array.

felupe.mesh.collect_faces(points, cells, cell_type)[source]#

“Collect all unique faces, calculate and return midpoints on faces as well as the additional cells array.

felupe.mesh.collect_volumes(points, cells, cell_type)[source]#

“Collect all volumes, calculate and return midpoints on volumes as well as the additional cells array.

felupe.mesh.concatenate(meshes)[source]#

Join a sequence of meshes with identical cell types.

felupe.mesh.convert(points, cells, cell_type, order=0, calc_points=False, calc_midfaces=False, calc_midvolumes=False)[source]#

Convert mesh to a given order (only order=0 and order=2 from order=1 are supported).

felupe.mesh.dual(points, cells, cell_type, points_per_cell=None, disconnect=True, calc_points=False, offset=0, npoints=None)[source]#

Create a new dual mesh with given points per cell. The point coordinates are not used in a dual mesh and hence, by default they are all zero.

felupe.mesh.expand(points, cells, cell_type, n=11, z=1)[source]#

Expand a 1d-Line to a 2d-Quad or a 2d-Quad to a 3d-Hexahedron Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • n (int, optional) – Number of n-point repetitions or (n-1)-cell repetitions, default is 11.

  • z (float or ndarray, optional) – Total expand dimension as float (edge length in expand direction is z / n), default is 1. Optionally, if an array is passed these entries are taken as expansion and n is ignored.

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.fill_between(mesh, other_mesh, n=11)[source]#

Fill a 2d-Quad Mesh between two 1d-Line Meshes, embedded in 2d-space, or a 3d-Hexahedron Mesh between two 2d-Quad Meshes, embedded in 3d-space, by expansion. Both meshes must have equal number of points and cells. The cells-array is taken from the first mesh.

Parameters:
  • mesh (felupe.Mesh) – The base line- or quad-mesh.

  • other_mesh (felupe.Mesh) – The other line- or quad-mesh.

  • n (int or ndarray) – Number of n-point repetitions or (n-1)-cell repetitions, (default is 11). If an array is given, then its values are used for the relative positions in a reference configuration (-1, 1) between the two meshes.

Returns:

mesh – The expanded mesh.

Return type:

felupe.Mesh

felupe.mesh.flip(points, cells, cell_type, mask=None)[source]#

Ensure positive cell volumes for tria, tetra, quad and hexahedron cell types.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • mask (list or ndarray, optional) – Boolean mask for selected cells to flip.

Returns:

  • points (ndarray) – Point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.mirror(points, cells, cell_type, normal=[1, 0, 0], centerpoint=[0, 0, 0], axis=None)[source]#

Mirror points by plane normal and ensure positive cell volumes for tria, tetra, quad and hexahedron cell types.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • normal (list or ndarray, optional) – Mirror-plane normal vector (default is [1, 0, 0]).

  • centerpoint (list or ndarray, optional) – Center-point coordinates on the mirror plane (default is [0, 0, 0]).

  • axis (int or None, optional) – Mirror axis (default is None).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.revolve(points, cells, cell_type, n=11, phi=180, axis=0)[source]#

Revolve a 2d-Quad to a 3d-Hexahedron Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • n (int, optional) – Number of n-point revolutions (or (n-1) cell revolutions), default is 11.

  • phi (float or ndarray, optional) – Revolution angle in degree (default is 180).

  • axis (int, optional) – Revolution axis (default is 0).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.rotate(points, cells, cell_type, angle_deg, axis, center=None)[source]#

Rotate a Mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • angle_deg (int) – Rotation angle in degree.

  • axis (int) – Rotation axis.

  • center (list or ndarray or None, optional) – Center point coordinates (default is None).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.runouts(points, cells, cell_type, values=[0.1, 0.1], centerpoint=[0, 0, 0], axis=0, exponent=5, mask=slice(None, None, None), normalize=False)[source]#

Add simple rubber-runouts for realistic rubber-metal structures.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • values (list or ndarray, optional) – Relative amount of runouts (per coordinate) perpendicular to the axis (default is 10% per coordinate, i.e. [0.1, 0.1]).

  • centerpoint (list or ndarray, optional) – Center-point coordinates (default is [0, 0, 0]).

  • axis (int or None, optional) – Axis (default is 0).

  • exponent (int, optional) – Positive exponent to control the shape of the runout. The higher the exponent, the steeper the transition (default is 5).

  • mask (list or None, optional) – List of points to be considered (default is None).

  • normalize (bool, optional) – Normalize the runouts to create indents, i.e. maintain the original shape at the ends (default is False).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.sweep(points, cells, cell_type, decimals=None)[source]#

Merge duplicated points and update cells of a Mesh.

WARNING: This function re-sorts points.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • decimals (int or None, optional) – Number of decimals for point coordinate comparison (default is None).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (list or ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

felupe.mesh.triangulate(points, cells, cell_type, mode=3)[source]#

Triangulate a quad or a hex mesh.

Parameters:
  • points (list or ndarray) – Original point coordinates.

  • cells (list or ndarray) – Original point-connectivity of cells.

  • cell_type (str) – A string in VTK-convention that specifies the cell type.

  • mode (int, optional) – Choose a mode how to convert hexahedrons to tets [1] (default is 3).

Returns:

  • points (ndarray) – Modified point coordinates.

  • cells (ndarray) – Modified point-connectivity of cells.

  • cell_type (str or None) – A string in VTK-convention that specifies the cell type.

References

[1] Dompierre, J., Labbé, P., Vallet, M. G., & Camarero, R. (1999). How to Subdivide Pyramids, Prisms, and Hexahedra into Tetrahedra. IMR, 99, 195.