spateo.tdr.models.models_individual.mesh_methods#

This code contains methods to reconstruct the mesh model based on the 3D point cloud.
  1. 3D Delaunay

  2. marching cube algorithm

  3. alpha shape algorithm

  4. ball pivoting algorithm

  5. poisson reconstruction

Module Contents#

Functions#

pv_mesh(→ pyvista.PolyData)

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

rigid_transform(→ numpy.ndarray)

Compute optimal transformation based on the two sets of points and apply the transformation to other points.

marching_cube_mesh(pc[, levelset, mc_scale_factor, ...])

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

_pv2o3d(pc)

Convert a point cloud in PyVista to Open3D.

_o3d2pv(→ pyvista.PolyData)

Convert a triangle mesh in Open3D to PyVista.

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.

poisson_mesh(→ pyvista.PolyData)

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

spateo.tdr.models.models_individual.mesh_methods.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.

spateo.tdr.models.models_individual.mesh_methods.rigid_transform(coords: numpy.ndarray, coords_refA: numpy.ndarray, coords_refB: numpy.ndarray) numpy.ndarray[source]#

Compute optimal transformation based on the two sets of points and apply the transformation to other points.

Parameters:
coords

Coordinate matrix needed to be transformed.

coords_refA

Referential coordinate matrix before transformation.

coords_refB

Referential coordinate matrix after transformation.

Returns:

The coordinate matrix after transformation

spateo.tdr.models.models_individual.mesh_methods.marching_cube_mesh(pc: pyvista.PolyData, levelset: int | float = 0, mc_scale_factor: int | float = 1.0, dist_sample_num: int | None = None)[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.

dist_sample_num

The down-sampling number when calculating the scaling factor using the minimum distance. Set to 100 for computation efficiency.

Returns:

A mesh model.

spateo.tdr.models.models_individual.mesh_methods._pv2o3d(pc: pyvista.PolyData)[source]#

Convert a point cloud in PyVista to Open3D.

spateo.tdr.models.models_individual.mesh_methods._o3d2pv(trimesh) pyvista.PolyData[source]#

Convert a triangle mesh in Open3D to PyVista.

spateo.tdr.models.models_individual.mesh_methods.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.models_individual.mesh_methods.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.models_individual.mesh_methods.poisson_mesh(pc: pyvista.PolyData, depth: int = 8, width: float = 0, scale: float = 1.1, linear_fit: bool = False, density_threshold: float | None = 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.