API Reference
This section provides complete API documentation for the einit library.
einit – Ellipsoid ICP initialization
- Provides register_ellipsoid(src, dst) that:
Centers src & dst by their centroids.
Computes ellipsoid matrices; eigen-decomposes to get principal axes.
Searches all 8 diagonal ±1 reflections for best alignment.
Returns a 4×4 homogeneous transform for OpenCV.
- einit.einit.register_ellipsoid(src_points, dst_points, src_features=None, dst_features=None, feature_weight=0.0, max_correspondence_distance=None, min_inlier_fraction=0.5, leafsize=16, positive_only=False)[source]
Compute initial transformation between 3D point clouds using ellipsoid analysis.
This function computes an initial rigid transformation that aligns the source point cloud with the destination point cloud by analyzing their ellipsoids of inertia. The algorithm uses KD-tree correspondence recovery to handle point clouds with different orderings, partial overlaps, and outliers.
- Parameters:
src_points (array_like, shape (N, 3)) – Source point cloud as N×3 array of 3D coordinates.
dst_points (array_like, shape (M, 3)) – Destination point cloud as M×3 array of 3D coordinates. N and M can be different (partial overlap).
src_features (array_like, shape (N, k), optional) – Per-point feature vectors for the source cloud (e.g. RGB colour, LiDAR intensity, surface normals). 1-D arrays of shape (N,) are accepted and treated as (N, 1). If provided, dst_features must also be given.
dst_features (array_like, shape (M, k), optional) – Per-point feature vectors for the destination cloud. Must have the same number of columns as src_features.
feature_weight (float, default 0.0) –
Controls how strongly features influence the alignment. 0.0 (default) gives the original geometry-only behaviour. Typical useful range: 0.1–1.0.
Features enter the algorithm in two places:
Ellipsoid step – the spatial covariance is augmented by the spatial-feature cross-covariance:
E_aug = P^T P + (feature_weight / N) * (P^T F) (P^T F)^T
This biases the principal axes toward spatial directions where features vary most, breaking eigenvalue degeneracy for symmetric shapes (spheres, cubes) where geometry alone gives arbitrary axes.
KD-tree step – correspondences are searched in the augmented space [x, y, z, w·f₁, …, w·fₖ], so that feature similarity guides which destination point is selected as nearest neighbour. Inlier filtering and error scoring remain in spatial (coordinate) units.
Features are normalised to unit standard deviation (per column, estimated from dst_features) before scaling, so feature_weight is dimensionless and comparable across datasets.
max_correspondence_distance (float, optional) – Maximum distance for valid point correspondences. Points farther than this distance from their nearest neighbors are considered outliers. If None (default), automatically estimated as 3× the median nearest- neighbor distance within the destination point cloud. Always interpreted in spatial (coordinate) units, even when features are used.
min_inlier_fraction (float, default 0.5) – Minimum fraction of source points that must have valid correspondences within max_correspondence_distance. Transformations with fewer inliers are rejected. Must be between 0 and 1.
leafsize (int, default 16) – KD-tree leaf size parameter. Affects search performance vs memory usage. Smaller values may improve accuracy for small point clouds but increase build time. Typical range: 8-32.
positive_only (bool, default False) – If True, only search proper rotations (determinant +1) by considering only sign combinations with an even number of negative values. This prevents reflections and ensures chirality preservation. Recommended when point distributions are spatially biased (e.g., bounding box overlap).
- Returns:
T – Homogeneous transformation matrix that transforms src_points to align with dst_points. Apply as: dst_aligned = (src @ T[:3,:3].T) + T[:3,3]
- Return type:
ndarray, shape (4, 4)
- Raises:
ValueError – If input arrays don’t have shape (N, 3) or (M, 3), if features are supplied inconsistently, or if feature_weight is negative.
Examples
Basic usage:
>>> import numpy as np >>> from einit import register_ellipsoid >>> src = np.random.randn(100, 3) >>> dst = np.random.randn(80, 3) # Different size OK >>> T = register_ellipsoid(src, dst) >>> T.shape (4, 4)
With custom parameters:
>>> T = register_ellipsoid( ... src, dst, ... max_correspondence_distance=0.1, ... min_inlier_fraction=0.7, ... leafsize=8, ... positive_only=True ... )
With per-point features (e.g. LiDAR intensity or RGB colour):
>>> src_intensity = np.random.rand(100, 1) # scalar reflectance per point >>> dst_intensity = np.random.rand(80, 1) >>> T = register_ellipsoid( ... src, dst, ... src_features=src_intensity, ... dst_features=dst_intensity, ... feature_weight=0.3, ... )
Notes
The algorithm is permutation-invariant: point ordering in the input arrays does not affect the result. It handles partial overlaps, noise, and outliers through KD-tree correspondence recovery and distance-based filtering.
Time complexity is O(N + M log M) where N and M are the number of points in the source and destination clouds respectively.