cedalion.geometry package

Subpackages

Submodules

cedalion.geometry.ellipsoid module

Ellipsoid head model for estimating landmark positions from head measurements.

cedalion.geometry.ellipsoid.get_landmarks_for_headsize(
circumference: Annotated[Quantity, '[length]'],
nz_cz_iz: Annotated[Quantity, '[length]'],
lpa_cz_rpa: Annotated[Quantity, '[length]'],
) Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')][source]

Estimate 10-20 landmark positions for a given head size on an ellipsoidal model.

Fits a tri-axial ellipsoid whose arc lengths match the provided head circumference and two arc measurements, then returns the positions of Cz, LPA, RPA, Nz, and Iz in the fitted ellipsoid coordinate system.

Parameters:
  • circumference – Head circumference (e.g. 560 mm).

  • nz_cz_iz – Arc length from Nz to Iz passing through Cz (anterior–posterior arc, scaled by 1.2/2 internally to match the 10-20 convention).

  • lpa_cz_rpa – Arc length from LPA to RPA passing through Cz (lateral arc, scaled by 1.2/2 internally).

Returns:

LabeledPoints with five landmarks (Cz, LPA, RPA, Nz, Iz) in the "ellipsoid" CRS, with units of mm.

cedalion.geometry.landmarks module

Module for constructing the 10-10-system on the scalp surface.

cedalion.geometry.landmarks.normalize_landmarks_labels(
geo3d: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
) Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')][source]

Normalize landmark labels to canonical names.

Maps commonly used alternative landmark names and capitalizations to canonical names.

When multiple labels normalize to the same canonical name: - If the canonical name already exists, alternative forms are not altered. - If multiple alternatives exist without the canonical form, a ValueError is

raised and the user must resolve the ambiguity.

Parameters:

geo3d – LabeledPoints with potentially non-canonical landmark names.

Returns:

LabeledPoints with normalized landmark labels.

class cedalion.geometry.landmarks.LandmarksBuilder1010(
scalp_surface: Surface,
landmarks: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
)[source]

Bases: object

Construct the 10-10-system on scalp surface based on Oostenveld and Praamstra [OP01].

scalp_surface[source]

a triangle-mesh representing the scalp

Type:

Surface

landmarks_mm[source]

positions of all 10-10 landmarks in mm

Type:

LabeledPoints

vtk_mesh[source]

the scalp surface as a VTK mesh

Type:

vtk.vtkPolyData

lines[source]

points along the lines connecting the landmarks

Type:

List[np.ndarray]

build()[source]

Construct the 10-10-system on the scalp surface.

plot()[source]

Plot scalp surface with landmarks.

cedalion.geometry.landmarks.order_ref_points_6(
landmarks: DataArray,
twoPoints: str,
) DataArray[source]

Reorder a set of six landmarks based on spatial relationships and give labels.

Parameters:
  • landmarks (xr.DataArray) – coordinates for six landmark points

  • twoPoints (str) – two reference points (‘Nz’ or ‘Iz’) for orientation.

Returns:

the landmarks ordered as “Nz”, “Iz”, “RPA”, “LPA”, “Cz”

Return type:

xr.DataArray

cedalion.geometry.meshing module

Mesh decimation, voxel-to-vertex mapping, and scalar upscaling utilities.

cedalion.geometry.meshing.upscale_scalars(
highres_mesh: TrimeshSurface,
lowres_mesh: TrimeshSurface,
lowres_scalars: ndarray,
)[source]

Upscale a scalar function from a low-resolution mesh to a higher-resolution one.

Uses barycentric interpolation: for each high-res vertex, the closest point on the low-res mesh is found and the scalar value is interpolated from the enclosing face’s three vertices. The low-res mesh must be a spatial subset of the high-res mesh (e.g. produced by decimate_mesh()).

Parameters:
  • highres_mesh – Target surface at higher resolution.

  • lowres_mesh – Source surface at lower resolution.

  • lowres_scalars – Scalar values defined at each vertex of lowres_mesh, shape (n_lowres_vertices,).

Returns:

NumPy array of scalar values interpolated onto highres_mesh vertices, shape (n_highres_vertices,).

cedalion.geometry.meshing.decimate_mesh(
surface: TrimeshSurface,
nvertex_target: int,
vertex_quality=None,
selected=False,
selection_threshold=0.5,
)[source]

Decimate a triangulated surface mesh to a target vertex count.

Uses PyMeshLab’s quadric edge collapse algorithm. When vertex_quality is provided, it is used as a per-vertex quality scalar to guide the decimation. Optionally, only vertices below a quality threshold are decimated when selected=True.

Parameters:
  • surface – Input TrimeshSurface.

  • nvertex_target – Desired number of vertices in the output mesh.

  • vertex_quality – Optional per-vertex quality scalar array of shape (n_vertices,) used to weight the decimation.

  • selected – If True, only decimate vertices whose quality is below selection_threshold.

  • selection_threshold – Quality threshold for vertex selection when selected=True.

Returns:

New TrimeshSurface with approximately nvertex_target vertices. Per-vertex coordinates stored in surface.vertex_coords are transferred via nearest-neighbour lookup.

cedalion.geometry.meshing.map_voxels_to_vertices(
surface: TrimeshSurface,
cell_coords,
)[source]

Map voxel centres to their nearest surface vertex.

Projects each voxel coordinate onto the surface mesh and then finds the nearest vertex. Processed in chunks to limit memory use.

Parameters:
  • surface – Target TrimeshSurface.

  • cell_coords – Array of voxel-centre coordinates, shape (N, 3).

Returns:

  • voxel2vertex_indices (np.ndarray[int], shape (N,)): index of the nearest vertex for each voxel.

  • voxel_count (np.ndarray[int], shape (n_vertices,)): number of voxels mapped to each vertex.

Return type:

Tuple (voxel2vertex_indices, voxel_count)

cedalion.geometry.meshing.parcel_aware_voxels_to_vertices_map(
surface: TrimeshSurface,
cell_coords,
skip_parcels=('Background+FreeSurfer_Defined_Medial_Wall_LH', 'Background+FreeSurfer_Defined_Medial_Wall_RH'),
voxel_stealing=False,
)[source]

Map voxel centres to surface vertices respecting parcel boundaries.

Each voxel is mapped only to vertices within the same cortical parcel, preventing leakage across parcel boundaries. Optionally, vertices that would otherwise receive no voxel can “steal” the nearest voxel from a neighbouring vertex (voxel_stealing).

Parameters:
  • surface – Target TrimeshSurface whose vertices carry a "parcel" coordinate.

  • cell_coords – xr.DataArray of voxel-centre coordinates with a "parcel" coordinate that groups voxels by parcel.

  • skip_parcels – Parcel labels to ignore (typically medial-wall parcels).

  • voxel_stealing – If True, reassign voxels to ensure every parcel vertex gets at least one voxel via linear assignment.

Returns:

  • voxel2vertex_indices (xr.DataArray[int]): nearest vertex index for each voxel (-1 for skipped parcels).

  • voxel_count (np.ndarray[int], shape (n_vertices,)): number of voxels mapped to each vertex.

Return type:

Tuple (voxel2vertex_indices, voxel_count)

cedalion.geometry.registration module

Registrating optodes to scalp surfaces.

class cedalion.geometry.registration.SpringICPResult(
initial_positions: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
nominal_distances: dict,
spring_errors: DataArray,
landmark_errors: DataArray,
snap_displacement_per_iter: ndarray,
n_iterations: int,
converged: bool,
)[source]

Bases: object

Quality-control details returned by register_optodes_spring_icp().

All distances and displacements are in the same units as the scalp surface.

initial_positions: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')][source]

Positions after the initial landmark alignment, before spring relaxation.

nominal_distances: dict[source]

Nominal channel distances {(src_label, det_label): float} used as spring rest lengths. Measured in scalp units from the phase-1 aligned positions unless supplied by the caller.

spring_errors: DataArray[source]

Per-channel actual_distance - nominal_distance at convergence. Dim channel; auxiliary coords source and detector carry the individual optode labels.

landmark_errors: DataArray[source]

Per-anchor Euclidean distance between the final optode position and the target landmark position at convergence. Dim label. Empty when no optode labels overlap with landmarks_scalp.

snap_displacement_per_iter: ndarray[source]

Maximum surface-projection displacement (over all optodes) per iteration. Shape (n_iterations,). Use as a convergence diagnostic.

n_iterations: int[source]

Number of iterations actually performed.

converged: bool[source]

True when the convergence criterion was met before n_iter.

cedalion.geometry.registration.register_identity(
coords_target: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
coords_trafo: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
)[source]

Affine transformation that just adapts units.

Parameters:
  • coords_target (LabeledPoints) – Target point cloud.

  • coords_trafo (LabeledPoints) – Source point cloud.

Returns:

Affine transformation between the two point clouds.

Return type:

cdt.AffineTransform

cedalion.geometry.registration.register_trans_rot(
coords_target: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
coords_trafo: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
)[source]

Finds affine transformation between coords_target and coords_trafo.

Uses translation and roatation rotation. Requires at least 3 common labels between the two point clouds.

Parameters:
  • coords_target (LabeledPoints) – Target point cloud.

  • coords_trafo (LabeledPoints) – Source point cloud.

Returns:

Affine transformation between the two point clouds.

Return type:

cdt.AffineTransform

cedalion.geometry.registration.register_general_affine(
coords_target: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
coords_trafo: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
)[source]

Finds affine transformation between coords_target and coords_trafo.

This method fits all 12 parameters of a 3D affine transformation without constraints. The transform thus contains scaling, rotation, translation, shear and mirroring, i.e. it can transform between left and right handed coordinate systems.

Parameters:
  • coords_target (LabeledPoints) – Target point cloud.

  • coords_trafo (LabeledPoints) – Source point cloud.

Returns:

Affine transformation between the two point clouds.

Return type:

cdt.AffineTransform

cedalion.geometry.registration.register_trans_rot_isoscale(
coords_target: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
coords_trafo: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
)[source]

Finds affine transformation between coords_target and coords_trafo.

Uses translation, rotation and isotropic scaling. Requires at least 3 common labels between the two point clouds.

Parameters:
  • coords_target (LabeledPoints) – Target point cloud.

  • coords_trafo (LabeledPoints) – Source point cloud.

Returns:

Affine transformation between the two point clouds.

Return type:

cdt.AffineTransform

cedalion.geometry.registration.register_trans_rot_scale(
coords_target: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
coords_trafo: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
)[source]

Finds affine transformation between coords_target and coords_trafo.

Uses translation, rotation and scaling. Requires at least 3 common labels between the two point clouds.

Parameters:
  • coords_target (LabeledPoints) – Target point cloud.

  • coords_trafo (LabeledPoints) – Source point cloud.

Returns:

Affine transformation between the two point clouds.

Return type:

cdt.AffineTransform

cedalion.geometry.registration.gen_xform_from_pts(p1: ndarray, p2: ndarray) ndarray[source]

Calculate the affine transformation matrix T that transforms p1 to p2.

Parameters:
  • p1 (np.ndarray) – Source points (p x m) where p is the number of points and m is the number of dimensions.

  • p2 (np.ndarray) – Target points (p x m) where p is the number of points and m is the number of dimensions.

Returns:

Affine transformation matrix T.

cedalion.geometry.registration.register_icp(
surface: Surface,
landmarks: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
geo3d: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
niterations=1000,
random_sample_fraction=0.5,
)[source]

Iterative Closest Point algorithm for registration.

Parameters:
  • surface (Surface) – Surface mesh to which to register the points.

  • landmarks (LabeledPoints) – Landmarks to use for registration.

  • geo3d (LabeledPoints) – Points to register to the surface.

  • niterations (int) – Number of iterations for the ICP algorithm (default 1000).

  • random_sample_fraction (float) – Fraction of points to use in each iteration (default 0.5).

Returns:

Tuple containing the losses and transformations

Return type:

Tuple[np.ndarray, np.ndarray]

cedalion.geometry.registration.icp_with_full_transform(
opt_centers: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
montage_points: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
max_iterations: int = 50,
tolerance: float = 500.0,
)[source]

Perform Iterative Closest Point algorithm with full transformation capabilities.

Parameters:
  • opt_centers – Source point cloud for alignment.

  • montage_points – Target reference point cloud.

  • max_iterations – Maximum number of iterations for convergence.

  • tolerance – Tolerance for convergence check.

Returns:

Transformed source points as a numpy array with their coordinates

updated to reflect the best alignment.

np.ndarray: Transformation parameters array consisting of

[tx, ty, tz, rx, ry, rz, sx, sy, sz], where ‘t’ stands for translation components, ‘r’ for rotation components (in radians), and ‘s’ for scaling components.

np.ndarray: Indices of the target points that correspond to each source point as

per the nearest neighbor search.

Return type:

np.ndarray

cedalion.geometry.registration.find_spread_points(points_xr: DataArray) ndarray[source]

Selects three points that are spread apart from each other in the dataset.

Parameters:

points_xr – An xarray DataArray containing the points from which to select.

Returns:

Indices of the initial, farthest, and median-distanced points from the initial point as determined by their positions in the original dataset.

cedalion.geometry.registration.simple_scalp_projection(
geo3d: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
) Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')][source]

Projects 3D coordinates onto a 2D plane using a simple scalp projection.

Parameters:

geo3d (LabeledPoints) – 3D coordinates of points to project. Requires the landmarks Nz, LPA, and RPA.

Returns:

A LabeledPoints containing the 2D coordinates of the projected points.

cedalion.geometry.registration.register_optodes_spring_icp(
scalp: TrimeshSurface,
geo3d: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
channels: list[tuple[str, str]] | Annotated[DataArray, DataArraySchema(dims='time', coords='time', 'time', 'samples')] | DataFrame,
landmarks_scalp: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
nominal_distances: dict[tuple[str, str], float] | None = None,
n_iter: int = 400,
k_spring: float = 1.0,
k_anchor: float = 10.0,
step_size: float = 0.1,
convergence_tol: float = 0.01,
initial_align_mode: str = 'general',
) tuple[Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')], SpringICPResult][source]

Register optode coordinates onto a scalp mesh via spring relaxation and ICP.

Two-phase registration:

  1. Initial alignment — a global transform (translation + rotation + optional scale) is fitted to the matched landmark pairs and applied to geo3d to bring them into the scalp coordinate system.

  2. Spring-relaxation ICP — optodes are iteratively refined:

    • Hooke’s-law springs between each channel-forming source–detector pair, with rest length equal to the nominal inter-optode distance, resist distortion of channel geometry.

    • Strong anchor springs pull any coordinate in geo3d whose label appears in landmarks_scalp toward its known anatomical position on the scalp.

    • After each force step every optode is projected onto the nearest point on the scalp triangulated surface (ICP step).

Parameters:
  • scalp – Triangulated scalp surface in the target coordinate system.

  • geo3d – Probe coordinates (sources, detectors and landmarks) in the probe coordinate system.

  • channels(source_label, detector_label) pairs defining which optode pairs form channels and therefore get a spring.

  • landmarks_scalp – Anatomical landmark positions in the scalp CRS (same CRS as scalp). Any label that also appears in optodes is used as a strong anchor during the relaxation phase.

  • nominal_distances – Optional pre-specified rest lengths for each channel spring, keyed by (src_label, det_label). When None, distances are measured from the phase-1 aligned positions.

  • n_iter – Maximum number of relaxation iterations.

  • k_spring – Spring constant for channel springs. Force magnitude scales linearly with extension from the nominal distance.

  • k_anchor – Spring constant for landmark-anchor springs. Should exceed k_spring to enforce anatomical constraints strongly.

  • step_size – Fraction of the net force vector added to each position per iteration. Reduce if the relaxation is unstable.

  • convergence_tol – Stop early when the maximum surface-projection displacement across all optodes is below this value (in scalp units).

  • initial_align_mode – Phase-1 transform type. One of "trans_rot_isoscale" (7 DOF, default), "trans_rot" (6 DOF), "general" (12 DOF full affine), or "identity" (only adapt units).

Returns:

SpringICPResult with the final optode positions and per-iteration / per-channel quality-control metrics.

Raises:
  • CRSMismatchError – If landmarks_scalp is not in the same CRS as scalp.

  • ValueError – If initial_align_mode is not recognised, or if a channel label is absent from optodes.

cedalion.geometry.segmentation module

Functionality to work with segmented MRI scans.

cedalion.geometry.segmentation.voxels_from_segmentation(
segmentation_mask: DataArray,
segmentation_types: List[str],
isovalue=0.9,
fill_holes_in_mask=False,
) Voxels[source]

Generate voxels from a segmentation mask.

Parameters:
  • segmentation_mask – xr.DataArray Segmentation mask.

  • segmentation_types – List[str] List of segmentation types.

  • isovalue – float, optional Isovalue for marching cubes, by default 0.9.

  • fill_holes_in_mask – bool, optional Fill holes in the mask, by default False.

Returns:

cdc.Voxels

Voxels in voxel space.

cedalion.geometry.segmentation.surface_from_segmentation(
segmentation_mask: DataArray,
segmentation_types: List[str],
isovalue=0.9,
fill_holes_in_mask=False,
) Surface[source]

Create a surface from a segmentation mask.

Parameters:
  • segmentation_mask (xr.DataArray) – Segmentation mask with dimensions segmentation type, i, j, k.

  • segmentation_types (List[str]) – A list of segmentation types to include in the surface.

  • isovalue (Float) – The isovalue to use for the marching cubes algorithm.

  • fill_holes_in_mask (Bool) – Whether to fill holes in the mask before creating the surface.

Returns:

A cedalion.Surface object.

cedalion.geometry.segmentation.cell_coordinates(volume, flat: bool = False)[source]

Generate voxel-centre coordinates from a 3D volume.

Parameters:
  • volume – 3-D array whose shape defines the grid dimensions.

  • flat – If True, return a flat (N, 3) array with label, i, j, k coordinates; otherwise return a (ni, nj, nk, 3) grid array.

Returns:

xr.DataArray of voxel-centre coordinates in voxel (ijk) space, quantified as dimensionless pint units.

cedalion.geometry.segmentation.segmentation_postprocessing(
segmentation_dir: str,
mask_files: dict[str, str] = {'air': 'c6.nii', 'bone': 'c4.nii', 'csf': 'c3.nii', 'gray': 'c1.nii', 'skin': 'c5.nii', 'white': 'c2.nii'},
isSmooth: bool = True,
fixCSF: bool = True,
removeDisconnected: bool = True,
labelUnassigned: bool = True,
removeAir: bool = True,
subtractTissues: bool = True,
) dict[source]

Postprocess SPM12 tissue-probability maps into binary segmentation masks.

Applies a configurable sequence of steps: Gaussian smoothing, CSF continuity repair, removal of disconnected components, labelling of unassigned voxels, removal of air cavities, and tissue subtraction. Results are saved as NIfTI files alongside the originals.

Parameters:
  • segmentation_dir – Directory containing the SPM12 output files.

  • mask_files – Mapping from tissue name to filename within segmentation_dir. Defaults to the standard SPM12 output names (c1.niic6.nii).

  • isSmooth – Apply Gaussian smoothing to each tissue probability map before thresholding.

  • fixCSF – Repair CSF continuity between brain tissue and bone.

  • removeDisconnected – Remove small disconnected components from each mask.

  • labelUnassigned – Assign unlabelled voxels to the nearest tissue type via Gaussian distance blurring.

  • removeAir – Remove air cavities outside the head.

  • subtractTissues – Enforce non-overlapping masks by subtracting inner tissues from outer ones (inside-out order).

Returns:

Dictionary mapping tissue names to the filenames of the saved postprocessed NIfTI masks.

References

Huang et al. [HDS+13], Harmening and Miklody [HM22]

cedalion.geometry.segmentation.binaryMaskGenerator(d)[source]

Convert multi-channel probability maps to mutually exclusive binary masks.

For each voxel the tissue with the highest probability wins; ties go to the first tissue (index 0 = empty).

Parameters:

d – Array of shape (n_tissues, i, j, k) with tissue probability values.

Returns:

Binary mask array of shape (n_tissues + 1, i, j, k) where index 0 is the “empty/unassigned” class and indices 1…n_tissues correspond to the input tissues.

cedalion.geometry.segmentation.sizeOfObject(img, conn=None)[source]

Return connected-component sizes for a binary 3-D image, sorted descending.

Parameters:
  • img – Binary 3-D (or 2-D) NumPy array.

  • conn – Connectivity (number of neighbours to consider). Defaults to 8 for 2-D and 26 for 3-D images.

Returns:

  • size_descend (list[int]): component sizes sorted in descending order.

  • ind (list[int]): original component indices sorted by descending size.

Return type:

Tuple (size_descend, ind)

cedalion.geometry.utils module

Utility functions for geometric calculations.

cedalion.geometry.utils.m_trans(t: ndarray) ndarray[source]

Return a 4×4 homogeneous translation matrix.

Parameters:

t – Translation vector [tx, ty, tz].

Returns:

4×4 NumPy array encoding the translation as a homogeneous transform.

cedalion.geometry.utils.m_scale3(s: ndarray) ndarray[source]

Return a 4×4 homogeneous anisotropic scaling matrix.

Parameters:

s – Scale factors [sx, sy, sz] for each axis independently.

Returns:

4×4 NumPy array encoding the anisotropic scaling as a homogeneous transform.

cedalion.geometry.utils.m_scale1(s: ndarray) ndarray[source]

Return a 4×4 homogeneous isotropic scaling matrix.

Parameters:

s – Array whose first element is the uniform scale factor applied to all axes.

Returns:

4×4 NumPy array encoding the isotropic scaling as a homogeneous transform.

cedalion.geometry.utils.m_rot(angles: ndarray) ndarray[source]

Return a 4×4 homogeneous rotation matrix R = Rz(α)·Ry(β)·Rx(γ).

See https://en.wikipedia.org/wiki/Rotation_matrix#General_rotations.

Parameters:

angles – Euler angles [alpha, beta, gamma] in radians, corresponding to rotations about Z, Y, and X axes respectively.

Returns:

4×4 NumPy array encoding the combined rotation as a homogeneous transform.

cedalion.geometry.utils.cart2sph(
x: ndarray,
y: ndarray,
z: ndarray,
) tuple[ndarray, ndarray, ndarray][source]

Convert 3D cartesian into spherical coordinates.

Parameters:
  • x – cartesian x coordinates

  • y – cartesian y coordinates

  • z – cartesian z coordinates

Returns:

The spherical coordinates azimuth, elevation and radius as np.ndarrays.

cedalion.geometry.utils.pol2cart(
theta: ndarray,
rho: ndarray,
) tuple[ndarray, ndarray][source]

Convert 2D polar into 2D cartesian coordinates.

Parameters:
  • theta – polar theta/angle coordinates

  • rho – polar rho/radius coordinates

Returns:

The cartesian coordinates x and y as np.ndarrays.

Module contents

Tools for geometric calculations.