cedalion.dot.head_model
Two-surface head model and standard atlas loading for DOT forward modelling.
Functions
|
Load the inflated cortex surface for a standard atlas head model. |
|
Create a TwoSurfaceHeadmodel for common atlases. |
Classes
|
Head Model class to represent a segmented head. |
- class cedalion.dot.head_model.TwoSurfaceHeadModel(
- segmentation_masks: xr.DataArray,
- brain: cdc.Surface,
- scalp: cdc.Surface,
- landmarks: cdt.LabeledPoints,
- t_ijk2ras: cdt.AffineTransform,
- t_ras2ijk: cdt.AffineTransform,
- voxel_to_vertex_brain: scipy.sparse.spmatrix,
- voxel_to_vertex_scalp: scipy.sparse.spmatrix,
Bases:
objectHead Model class to represent a segmented head.
Its main functions are reduced to work on voxel projections to scalp and cortex surfaces.
- segmentation_masks[source]
xr.DataArray Segmentation masks of the head for each tissue type.
- Type:
xarray.core.dataarray.DataArray
- landmarks[source]
cdt.LabeledPoints Anatomical landmarks in RAS space.
- Type:
xarray.core.dataarray.DataArray
- t_ijk2ras[source]
cdt.AffineTransform Affine transformation from ijk to RAS space.
- Type:
xarray.core.dataarray.DataArray
- t_ras2ijk[source]
cdt.AffineTransform Affine transformation from RAS to ijk space.
- Type:
xarray.core.dataarray.DataArray
- voxel_to_vertex_brain[source]
scipy.sparse.spmatrix Mapping from voxel to brain vertices.
- Type:
scipy.sparse._matrix.spmatrix
- voxel_to_vertex_scalp[source]
scipy.sparse.spmatrix Mapping from voxel to scalp vertices.
- Type:
scipy.sparse._matrix.spmatrix
- from_segmentation(cls, segmentation_dir, mask_files, landmarks_ras_file,
brain_seg_types, scalp_seg_types, smoothing, brain_face_count, scalp_face_count): Construct instance from segmentation masks in NIfTI format.
- landmarks: Annotated[DataArray, DataArraySchema(dims='label', coords='label', 'label', 'type')][source]
- classmethod from_segmentation(
- segmentation_dir: str,
- mask_files: dict[str, str] = {'csf': 'csf.nii', 'gm': 'gm.nii', 'scalp': 'scalp.nii', 'skull': 'skull.nii', 'wm': 'wm.nii'},
- landmarks_ras_file: str | None = None,
- brain_seg_types: list[str] = ['gm', 'wm'],
- scalp_seg_types: list[str] = ['scalp'],
- smoothing: float = 0.0,
- brain_face_count: int | None = 180000,
- scalp_face_count: int | None = 60000,
- fill_holes: bool = True,
Constructor from binary masks as gained from segmented MRI scans.
- Parameters:
segmentation_dir – Folder containing the segmentation masks in NIFTI format.
mask_files – Dictionary mapping segmentation types to NIFTI filenames.
landmarks_ras_file – Filename of the landmarks in RAS space.
brain_seg_types – List of segmentation types to be included in the brain surface.
scalp_seg_types – List of segmentation types to be included in the scalp surface.
smoothing – Smoothing factor for the brain and scalp surfaces.
brain_face_count – Number of faces for the brain surface.
scalp_face_count – Number of faces for the scalp surface.
fill_holes – Whether to fill holes in the segmentation masks.
- classmethod from_surfaces(
- segmentation_dir: str,
- mask_files: dict[str, str] = {'csf': 'csf.nii', 'gm': 'gm.nii', 'scalp': 'scalp.nii', 'skull': 'skull.nii', 'wm': 'wm.nii'},
- brain_surface_file: str = None,
- scalp_surface_file: str = None,
- landmarks_ras_file: str | None = None,
- brain_seg_types: list[str] = ['gm', 'wm'],
- scalp_seg_types: list[str] = ['scalp'],
- smoothing: float = 0.0,
- brain_face_count: int | None = None,
- scalp_face_count: int | None = None,
- fill_holes: bool = False,
- coordinates_file: Path | str | None = None,
- voxel_to_vertex_mapping_file_brain: Path | None = None,
Constructor from seg.masks, brain and head surfaces as gained from MRI scans.
- Parameters:
segmentation_dir – Folder containing the segmentation masks in NIFTI format.
mask_files – Dictionary mapping segmentation types to NIFTI filenames.
brain_surface_file – Path to the brain surface.
scalp_surface_file – Path to the scalp surface.
landmarks_ras_file – Filename of the landmarks in RAS space.
brain_seg_types – List of segmentation types to be included in the brain surface.
scalp_seg_types – List of segmentation types to be included in the scalp surface.
smoothing – Smoothing factor for the brain and scalp surfaces.
brain_face_count – Number of faces for the brain surface.
scalp_face_count – Number of faces for the scalp surface.
fill_holes – Whether to fill holes in the segmentation masks.
coordinates_file – path to file with additional vertex coordinates.
voxel_to_vertex_mapping_file_brain – precomputed voxel to vertex mapping
- Returns:
An instance of the TwoSurfaceHeadModel class.
- Return type:
- apply_transform(
- transform: cdt.AffineTransform,
Apply a coordinate transformation to the head model.
- Parameters:
transform – Affine transformation matrix (4x4) to be applied.
- Returns:
Transformed head model.
- copy(copy_segmasks=False)[source]
Create a copy of this head model.
- Parameters:
copy_segmasks – if False, the new instance gets a reference to the existing segmentation masks instead of a full copy.
- save(foldername: str)[source]
Save the head model to a folder.
- Parameters:
foldername (str) – Folder to save the head model into.
- Returns:
None
- classmethod load(foldername: str)[source]
Load the head model from a folder.
- Parameters:
foldername (str) – Folder to load the head model from.
- Returns:
Loaded head model.
- Return type:
- align_and_snap_to_scalp(
- points: cdt.LabeledPoints,
- mode: str = 'general',
Align and snap optodes or points to the scalp surface.
- Parameters:
points (cdt.LabeledPoints) – Points to be aligned and snapped to the scalp surface.
mode – method to derive the affine transform. Could be either ‘trans_rot_isoscale’ or ‘general’. See cedalion.geometry.registraion for details.
- Returns:
Points aligned and snapped to the scalp surface.
- Return type:
cdt.LabeledPoints
- align_and_relax_to_scalp(
- points: cdt.LabeledPoints,
- channels: list[tuple[str, str]] | cdt.NDTimeSeries | pd.DataFrame,
- 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',
Align and project optodes onto the scalp using spring-relaxation ICP.
Unlike
align_and_snap_to_scalp(), which projects every optode independently to its nearest scalp point, this method preserves the probe geometry by coupling source–detector pairs through Hooke’s-law springs. The algorithm alternates between a spring-force step — attracting each pair toward its nominal channel distance — and an ICP projection step that keeps every optode on the scalp surface. Anatomical landmarks shared between the probe and the head model are enforced as strong anchor springs. The approach is based on the spring-relaxation registration introduced in AtlasViewer (Aasted et al. [AYucelC+15]).Two-phase registration:
Initial alignment — a global transform (translation + rotation + optional scale) fitted to matched landmark pairs brings
pointsinto the scalp coordinate system.Spring-relaxation ICP — optodes are iteratively refined: channel springs resist geometry distortion; anchor springs enforce anatomical positions; every optode is projected back onto the scalp mesh after each force step.
- Parameters:
points – Probe coordinates (sources, detectors and landmarks).
channels – Channel definitions — accepted forms are a
NDTimeSeriesDataArray (coordssourceanddetector), a measurement-listpandas.DataFrame(columnssource,detector,wavelength), or a list of(source_label, detector_label)tuples.nominal_distances – Optional pre-specified rest lengths for each channel spring, keyed by
(src_label, det_label). When None, distances are measured from the initially aligned positions.n_iter – Maximum number of spring-relaxation iterations.
k_spring – Spring constant for channel springs (Hooke’s law).
k_anchor – Spring constant for landmark-anchor springs. Should exceed
k_springto enforce anatomical constraints strongly.step_size – Fraction of the net force vector applied per iteration. Reduce if the relaxation is numerically unstable.
convergence_tol – Early-stop threshold: halt when the maximum surface-projection displacement across all optodes falls below this value (in scalp units, typically mm).
initial_align_mode – Global alignment method applied before relaxation. One of
"general"(12-DOF full affine),"trans_rot_isoscale"(7-DOF),"trans_rot"(6-DOF), or"identity"(unit conversion only).
- Returns:
(registered_points, details)whereregistered_pointsis aLabeledPointsDataArray with the final optode positions on the scalp, anddetailsis aSpringICPResultcontaining convergence diagnostics and per-channel quality metrics.- Return type:
tuple
References
Paper & Code: Aasted et al. [AYucelC+15]
See also
align_and_snap_to_scalp()for a faster but geometry-unaware alternative.
- snap_to_scalp_voxels(points: cdt.LabeledPoints) cdt.LabeledPoints[source]
Snap optodes or points to the closest scalp voxel.
- Parameters:
points (cdt.LabeledPoints) – Points to be snapped to the closest scalp voxel.
- Returns:
- Points aligned and snapped to the closest scalp
voxel.
- Return type:
cdt.LabeledPoints
- scale_to_landmarks(
- target_landmarks: cdt.LabeledPoints,
- mode='general',
Scale the head model to match a set of anatomical landmarks.
Computes a general affine transform that maps the head model’s RAS landmarks to
target_landmarksand applies it to both surfaces and the stored affine transforms.- Parameters:
target_landmarks – Target landmark positions (e.g. from a digitizer) in any CRS. Must contain the same label subset as the model’s landmarks.
mode – method to derive the affine transform. Could be either ‘trans_rot_isoscale’ or ‘general’. See cedalion.geometry.registraion for details.
- Returns:
New
TwoSurfaceHeadModelscaled and aligned totarget_landmarks.
- scale_to_headsize(
- circumference: cdt.QLength,
- nz_cz_iz: cdt.QLength,
- lpa_cz_rpa: cdt.QLength,
Scale the head model to match anthropometric head-size measurements.
Fits a tri-axial ellipsoid model to the three measurements, derives corresponding landmark positions, and delegates to
scale_to_landmarks().- Parameters:
circumference – Head circumference.
nz_cz_iz – Nasion–Cz–Inion arc length.
lpa_cz_rpa – Left preauricular–Cz–right preauricular arc length.
- Returns:
New
TwoSurfaceHeadModelscaled to the specified head size.
- get_brain_mni152_coords() DataArray[source]
Return MNI152 vertex coordinates of the brain surface.
- Returns:
xr.DataArray with dimensions
(label, mni152)containing themni152_r,mni152_a, andmni152_svertex coordinates.- Raises:
ValueError – If the brain surface does not have
mni152_r/a/svertex coordinates (i.e. not a standard atlas head model).
- assign_parcels_via_mni_coords(
- coordinate_label: str,
- label_mapping: dict[int, str],
- voxel_label_niftii: Path,
- voxel_label_crs: str = 'mni152',
- background_label='Background',
- mni_eps=5.0,
Assign parcel labels to brain vertices using a volumetric atlas in MNI space.
For each brain vertex, the nearest non-background voxel within
mni_epsin MNI152 space is found and its string label is assigned to the veretx. Vertices with no labeled neighbor within the search radius receivebackground_label.- Parameters:
coordinate_label – Name of the new vertex coordinate.
label_mapping – A dictionary that maps numeric voxel labels in the NIfTI file to string parcel names.
voxel_label_niftii – Path to the NIfTI atlas file whose voxels carry numeric region labels.
voxel_label_crs – Coordinate reference system of the NIfTI file, either
"mni152"or"mni305". Defaults to"mni152".background_label – String label used for unlabeled vertices and, if present in
label_mapping, to identify background voxels. Defaults to"Background".mni_eps – Radius of the ball search around each vertex. Defaults to
5.0. Increase this value to find non-background voxels for each vertex.
- Returns:
A copy of the head model with parcel labels assigned as a new coordinate to brain vertices.
- cedalion.dot.head_model.get_standard_headmodel(model: str) TwoSurfaceHeadModel[source]
Create a TwoSurfaceHeadmodel for common atlases.
Once created, this function caches a head model.
- Parameters:
model – either colin27 or icbm152
- Returns:
The loaded head model with 1010-landmarks and parcel labels assigned.
- cedalion.dot.head_model.get_inflated_cortex_surface(model: str) TrimeshSurface[source]
Load the inflated cortex surface for a standard atlas head model.
The result is cached after the first call.
- Parameters:
model – Atlas name — one of
"colin27"or"icbm152".- Returns:
TrimeshSurfacein the"inflated"CRS.- Raises:
ValueError – If
modelis not one of the supported atlases.