spateo.tdr.models.mesh#

Submodules#

Package Contents#

Functions#

construct_cells(pc, cell_size[, geometry, xyz_scale, ...])

Reconstructing cells from point clouds.

construct_surface(→ Tuple[pyvista.PolyData, ...)

Surface mesh reconstruction based on 3D point cloud model.

clean_mesh(→ pyvista.PolyData)

Removes unused points and degenerate cells.

fix_mesh(→ pyvista.PolyData)

Repair the mesh where it was extracted and subtle holes along complex parts of the mesh.

uniform_larger_pc(→ pyvista.PolyData)

Generates a uniform point cloud with a larger number of points.

uniform_mesh(→ pyvista.PolyData)

Generate a uniformly meshed surface using voronoi clustering.

alpha_shape_mesh(→ pyvista.PolyData)

Computes a triangle mesh from a point cloud based on the alpha shape algorithm.

ball_pivoting_mesh(pc[, radii])

Computes a triangle mesh from an oriented point cloud based on the ball pivoting algorithm.

marching_cube_mesh(pc[, levelset, mc_scale_factor])

Computes a triangle mesh from a point cloud based on the marching cube algorithm.

poisson_mesh(→ pyvista.PolyData)

Computes a triangle mesh from an oriented point cloud based on the screened poisson reconstruction.

pv_mesh(→ pyvista.PolyData)

Generate a 3D tetrahedral mesh from a scattered points and extract surface mesh of the 3D tetrahedral mesh.

spateo.tdr.models.mesh.construct_cells(pc: pyvista.PolyData, cell_size: numpy.ndarray, geometry: Literal[cube, sphere, ellipsoid] = 'cube', xyz_scale: tuple = (1, 1, 1), n_scale: tuple = (1, 1), factor: float = 0.5)[source]#

Reconstructing cells from point clouds.

Parameters
pc

A point cloud object, including pc.point_data["obs_index"].

geometry

The geometry of generating cells. Available geometry are:

  • geometry = 'cube'

  • geometry = 'sphere'

  • geometry = 'ellipsoid'

cell_size

A numpy.ndarray object including the relative radius/length size of each cell.

xyz_scale

The scale factor for the x-axis, y-axis and z-axis.

n_scale

The squareness parameter in the x-y plane adn z axis. Only works if geometry = 'ellipsoid'.

factor

Scale factor applied to scaling array.

Returns

A cells mesh including ds_glyph.point_data[“cell_size”], ds_glyph.point_data[“cell_centroid”] and the data contained in the pc.

Return type

ds_glyph

spateo.tdr.models.mesh.construct_surface(pc: pyvista.PolyData, key_added: str = 'groups', label: str = 'surface', color: Optional[str] = 'gainsboro', alpha: Union[float, int] = 1.0, uniform_pc: bool = False, uniform_pc_alpha: Union[float, int] = 0, cs_method: Literal[pyvista, alpha_shape, ball_pivoting, poisson, marching_cube] = 'marching_cube', cs_args: Optional[dict] = None, nsub: Optional[int] = 3, nclus: int = 20000, smooth: Optional[int] = 1000, scale_distance: Union[float, int, list, tuple] = None, scale_factor: Union[float, int, list, tuple] = None) Tuple[pyvista.PolyData, pyvista.PolyData][source]#

Surface mesh reconstruction based on 3D point cloud model.

Parameters
pc

A point cloud model.

key_added

The key under which to add the labels.

label

The label of reconstructed surface mesh model.

color

Color to use for plotting mesh. The default color is 'gainsboro'.

alpha

The opacity of the color to use for plotting mesh. The default alpha is 0.8.

uniform_pc

Generates a uniform point cloud with a larger number of points.

uniform_pc_alpha

Specify alpha (or distance) value to control output of this filter.

cs_method

The methods of generating a surface mesh. Available cs_method are:

  • 'pyvista': Generate a 3D tetrahedral mesh based on pyvista.

  • 'alpha_shape': Computes a triangle mesh on the alpha shape algorithm.

  • 'ball_pivoting': Computes a triangle mesh based on the Ball Pivoting algorithm.

  • 'poisson': Computes a triangle mesh based on thee Screened Poisson Reconstruction.

  • 'marching_cube': Computes a triangle mesh based on the marching cube algorithm.

cs_args

Parameters for various surface reconstruction methods. Available cs_args are: * 'pyvista': {‘alpha’: 0} * 'alpha_shape': {‘alpha’: 2.0} * 'ball_pivoting': {‘radii’: [1]} * 'poisson': {‘depth’: 8, ‘width’=0, ‘scale’=1.1, ‘linear_fit’: False, ‘density_threshold’: 0.01} * 'marching_cube': {‘levelset’: 0, ‘mc_scale_factor’: 1}

nsub

Number of subdivisions. Each subdivision creates 4 new triangles, so the number of resulting triangles is nface*4**nsub where nface is the current number of faces.

nclus

Number of voronoi clustering.

smooth

Number of iterations for Laplacian smoothing.

scale_distance

The distance by which the model is scaled. If scale_distance is float, the model is scaled same distance along the xyz axis; when the scale factor is list, the model is scaled along the xyz axis at different distance. If scale_distance is None, there will be no scaling based on distance.

scale_factor

The scale by which the model is scaled. If scale factor is float, the model is scaled along the xyz axis at the same scale; when the scale factor is list, the model is scaled along the xyz axis at different scales. If scale_factor is None, there will be no scaling based on scale factor.

Returns

A reconstructed surface mesh, which contains the following properties:

uniform_surf.cell_data[key_added], the label array; uniform_surf.cell_data[f'{key_added}_rgba'], the rgba colors of the label array.

inside_pc: A point cloud, which contains the following properties:

inside_pc.point_data['obs_index'], the obs_index of each coordinate in the original adata. inside_pc.point_data[key_added], the groupby information. inside_pc.point_data[f'{key_added}_rgba'], the rgba colors of the groupby information.

Return type

uniform_surf

spateo.tdr.models.mesh.clean_mesh(mesh: pyvista.PolyData) pyvista.PolyData[source]#

Removes unused points and degenerate cells.

spateo.tdr.models.mesh.fix_mesh(mesh: pyvista.PolyData) pyvista.PolyData[source]#

Repair the mesh where it was extracted and subtle holes along complex parts of the mesh.

spateo.tdr.models.mesh.uniform_larger_pc(pc: pyvista.PolyData, alpha: Union[float, int] = 0, nsub: Optional[int] = 5, nclus: int = 20000) pyvista.PolyData[source]#

Generates a uniform point cloud with a larger number of points. If the number of points in the original point cloud is too small or the distribution of the original point cloud is not uniform, making it difficult to construct the surface, this method can be used for preprocessing.

Parameters
pc

A point cloud model.

alpha

Specify alpha (or distance) value to control output of this filter. For a non-zero alpha value, only edges or triangles contained within a sphere centered at mesh vertices will be output. Otherwise, only triangles will be output.

nsub

Number of subdivisions. Each subdivision creates 4 new triangles, so the number of resulting triangles is nface*4**nsub where nface is the current number of faces.

nclus

Number of voronoi clustering.

Returns

A uniform point cloud with a larger number of points.

Return type

new_pc

spateo.tdr.models.mesh.uniform_mesh(mesh: pyvista.PolyData, nsub: Optional[int] = 3, nclus: int = 20000) pyvista.PolyData[source]#

Generate a uniformly meshed surface using voronoi clustering.

Parameters
mesh

A mesh model.

nsub

Number of subdivisions. Each subdivision creates 4 new triangles, so the number of resulting triangles is nface*4**nsub where nface is the current number of faces.

nclus

Number of voronoi clustering.

Returns

A uniform mesh model.

Return type

new_mesh

spateo.tdr.models.mesh.alpha_shape_mesh(pc: pyvista.PolyData, alpha: float = 2.0) pyvista.PolyData[source]#

Computes a triangle mesh from a point cloud based on the alpha shape algorithm. Algorithm Overview:

For each real number α, define the concept of a generalized disk of radius 1/α as follows:

If α = 0, it is a closed half-plane; If α > 0, it is a closed disk of radius 1/α; If α < 0, it is the closure of the complement of a disk of radius −1/α.

Then an edge of the alpha-shape is drawn between two members of the finite point set whenever there exists a generalized disk of radius 1/α containing none of the point set and which has the property that the two points lie on its boundary. If α = 0, then the alpha-shape associated with the finite point set is its ordinary convex hull.

Parameters
pc

A point cloud model.

alpha

Parameter to control the shape. With decreasing alpha value the shape shrinks and creates cavities. A very big value will give a shape close to the convex hull.

Returns

A mesh model.

spateo.tdr.models.mesh.ball_pivoting_mesh(pc: pyvista.PolyData, radii: List[float] = None)[source]#

Computes a triangle mesh from an oriented point cloud based on the ball pivoting algorithm. Algorithm Overview:

The main assumption this algorithm is based on is the following: Given three vertices, and a ball of radius r, the three vertices form a triangle if the ball is getting “caught” and settle between the points, without containing any other point. The algorithm stimulates a virtual ball of radius r. Each iteration consists of two steps:

  • Seed triangle - The ball rolls over the point cloud until it gets “caught” between three vertices and

    settles between in them. Choosing the right r promises no other point is contained in the formed triangle. This triangle is called “Seed triangle”.

  • Expanding triangle - The ball pivots from each edge in the seed triangle, looking for a third point. It

    pivots until it gets “caught” in the triangle formed by the edge and the third point. A new triangle is formed, and the algorithm tries to expand from it. This process continues until the ball can’t find any point to expand to.

At this point, the algorithm looks for a new seed triangle, and the process described above starts all over.

Useful Notes:
  1. The point cloud is “dense enough”;

  2. The chosen r size should be “slightly” larger than the average space between points.

Parameters
pc

A point cloud model.

radii

The radii of the ball that are used for the surface reconstruction. This is a list of multiple radii that will create multiple balls of different radii at the same time.

Returns

A mesh model.

spateo.tdr.models.mesh.marching_cube_mesh(pc: pyvista.PolyData, levelset: Union[int, float] = 0, mc_scale_factor: Union[int, float] = 1.0)[source]#

Computes a triangle mesh from a point cloud based on the marching cube algorithm. Algorithm Overview:

The algorithm proceeds through the scalar field, taking eight neighbor locations at a time (thus forming an imaginary cube), then determining the polygon(s) needed to represent the part of the iso-surface that passes through this cube. The individual polygons are then fused into the desired surface.

Parameters
pc

A point cloud model.

levelset

The levelset of iso-surface. It is recommended to set levelset to 0 or 0.5.

mc_scale_factor

The scale of the model. The scaled model is used to construct the mesh model.

Returns

A mesh model.

spateo.tdr.models.mesh.poisson_mesh(pc: pyvista.PolyData, depth: int = 8, width: float = 0, scale: float = 1.1, linear_fit: bool = False, density_threshold: Optional[float] = None) pyvista.PolyData[source]#

Computes a triangle mesh from an oriented point cloud based on the screened poisson reconstruction.

Parameters
pc

A point cloud model.

depth

Maximum depth of the tree that will be used for surface reconstruction. Running at depth d corresponds to solving on a grid whose resolution is no larger than 2^d x 2^d x 2^d.

Note that since the reconstructor adapts the octree to the sampling density, the specified reconstruction depth is only an upper bound.

The depth that defines the depth of the octree used for the surface reconstruction and hence implies the resolution of the resulting triangle mesh. A higher depth value means a mesh with more details.

width

Specifies the target width of the finest level octree cells. This parameter is ignored if depth is specified.

scale

Specifies the ratio between the diameter of the cube used for reconstruction and the diameter of the samples’ bounding cube.

linear_fit

If true, the reconstructor will use linear interpolation to estimate the positions of iso-vertices.

density_threshold

The threshold of the low density.

Returns

A mesh model.

spateo.tdr.models.mesh.pv_mesh(pc: pyvista.PolyData, alpha: float = 2.0) pyvista.PolyData[source]#

Generate a 3D tetrahedral mesh from a scattered points and extract surface mesh of the 3D tetrahedral mesh.

Parameters
pc

A point cloud model.

alpha

Distance value to control output of this filter. For a non-zero alpha value, only vertices, edges, faces, or tetrahedron contained within the circumspect (of radius alpha) will be output. Otherwise, only tetrahedron will be output.

Returns

A mesh model.