cedalion.geometry package
Subpackages
- cedalion.geometry.photogrammetry package
- Submodules
- cedalion.geometry.photogrammetry.processors module
minEnclosingCircle()
ScanProcessor
pca()
ColoredStickerProcessorDetails
ColoredStickerProcessorDetails.cluster_coords
ColoredStickerProcessorDetails.cluster_circles
ColoredStickerProcessorDetails.cluster_colors
ColoredStickerProcessorDetails.vertex_colors
ColoredStickerProcessorDetails.vertex_hue
ColoredStickerProcessorDetails.vertex_value
ColoredStickerProcessorDetails.cfg_colors
ColoredStickerProcessorDetails.vertex_sat
ColoredStickerProcessorDetails.plot_cluster_circles()
ColoredStickerProcessorDetails.plot_vertex_colors1()
ColoredStickerProcessorDetails.plot_vertex_colors()
ColoredStickerProcessor
geo3d_from_scan()
- Module contents
Submodules
cedalion.geometry.landmarks module
Module for constructing the 10-10-system on the scalp surface.
- class cedalion.geometry.landmarks.LandmarksBuilder1010(
- scalp_surface: Surface,
- landmarks: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')],
Bases:
object
Construct the 10-10-system on scalp surface based on Oostenveld and Praamstra [OP01].
- cedalion.geometry.landmarks.order_ref_points_6(
- landmarks: DataArray,
- twoPoints: str,
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.registration module
Registrating optodes to scalp surfaces.
- 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')],
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 (LabeledPointCloud) – Target point cloud.
coords_trafo (LabeledPointCloud) – 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')],
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 (LabeledPointCloud) – Target point cloud.
coords_trafo (LabeledPointCloud) – 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')],
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 (LabeledPointCloud) – Target point cloud.
coords_trafo (LabeledPointCloud) – 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,
Iterative Closest Point algorithm for registration.
- Parameters:
surface (Surface) – Surface mesh to which to register the points.
landmarks (LabeledPointCloud) – Landmarks to use for registration.
geo3d (LabeledPointCloud) – 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,
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')],
Projects 3D coordinates onto a 2D plane using a simple scalp projection.
- Parameters:
geo3d (LabeledPointCloud) – 3D coordinates of points to project. Requires the landmarks Nz, LPA, and RPA.
- Returns:
A LabeledPointCloud containing the 2D coordinates of the projected points.
cedalion.geometry.segmentation module
Funtionality 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,
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,
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 cell coordinates from a 3D volume.
- Parameters:
volume (np.ndarray) – 3D volume.
flat (bool, optional) – If True, return coordinates as a flat array, by default False.
Returns
-------
xr.DataArray – Cell coordinates in voxel space.
- 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,
Postprocessing of the segmented SPM12 MRI segmentation files.
- Parameters:
segmentation_dir (str) – Directory where the segmented files are stored.
mask_files (dict[str, str], optional) – Dictionary containing the filenames of the segmented tissues.
isSmooth (bool, optional) – Smooth the segmented tissues using Gaussian filter.
fixCSF (bool, optional) – Fix the CSF continuity.
removeDisconnected (bool, optional) – Remove disconnected voxels.
labelUnassigned (bool, optional) – Label empty voxels to the nearest tissue type.
removeAir (bool, optional) – Remove air cavities.
subtractTissues (bool, optional) – Subtract tissues from each others
Returns:
- mask_filesdict
Dictionary containing the filenames of the postprocessed masks.
References:
This whole postprocessing is based on the following references: Huang et al. [HDS+13] Harmening and Miklody [HM22]
cedalion.geometry.utils module
Utility functions for geometric calculations.
- cedalion.geometry.utils.m_trans(t: ndarray) ndarray [source]
Calculate the affine transformation matrix for a tranlation t.
- cedalion.geometry.utils.m_scale3(s: ndarray) ndarray [source]
Calculate the affine transformation matrix for scaling s.
Apply different scaling factors for each dimension.
- cedalion.geometry.utils.m_scale1(s: ndarray) ndarray [source]
Calculate the affine transformation matrix for scaling s.
Apply one scaling factor for all dimensions.
- cedalion.geometry.utils.m_rot(angles: ndarray) ndarray [source]
Calculate the affine transformation matrix for a 3D rotation.
R = Rz(alpha)Ry(beta)Rx(gamma)
https://en.wikipedia.org/wiki/Rotation_matrix#General_rotations
- cedalion.geometry.utils.cart2sph(
- x: ndarray,
- y: ndarray,
- z: ndarray,
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.
Module contents
Tools for geometric calculations.