cedalion.imagereco package
Submodules
cedalion.imagereco.forward_model module
- class cedalion.imagereco.forward_model.ForwardModel(head_model, geo3d, measurement_list)
Bases:
object
Forward model for simulating light transport in the head.
…
- Parameters:
head_model (
TwoSurfaceHeadModel
) – TwoSurfaceHeadModel Head model containing voxel projections to brain and scalp surfaces.optode_pos – cdt.LabeledPointCloud Optode positions.
optode_dir – xr.DataArray Optode orientations (directions of light beams).
tissue_properties – xr.DataArray Tissue properties for each tissue type.
volume – xr.DataArray Voxelated head volume from segmentation masks.
unitinmm – float Unit of head model, optodes expressed in mm.
measurement_list (
DataFrame
) – pd.DataFrame List of measurements of experiment with source, detector, channel and wavelength.
- compute_fluence(nphoton)
Compute fluence for each channel and wavelength from photon simulation.
- compute_sensitivity(fluence_all, fluence_at_optodes)
Compute sensitivity matrix from fluence.
- compute_fluence_mcx(nphoton=100000000.0)
Compute fluence for each channel and wavelength using MCX package.
- Parameters:
nphoton (
int
) – int Number of photons to simulate.- Returns:
- xr.DataArray
Fluence in each voxel for each channel and wavelength.
References
(Fang and Boas [FB09]) Qianqian Fang and David A. Boas, “Monte Carlo Simulation of Photon Migration in 3D Turbid Media Accelerated by Graphics Processing Units,” Optics Express, vol.17, issue 22, pp. 20178-20190 (2009).
(Yu et al. [YNPKF18]) Leiming Yu, Fanny Nina-Paravecino, David Kaeli, Qianqian Fang, “Scalable and massively parallel Monte Carlo photon transport simulations for heterogeneous computing platforms,” J. Biomed. Opt. 23(1), 010504 (2018).
(Yan and Fang [YF20]) Shijie Yan and Qianqian Fang* (2020), “Hybrid mesh and voxel based Monte Carlo algorithm for accurate and efficient photon transport modeling in complex bio-tissues,” Biomed. Opt. Express, 11(11) pp. 6262-6270. https://www.osapublishing.org/boe/abstract.cfm?uri=boe-11-11-6262
- compute_fluence_nirfaster(meshingparam=None)
Compute fluence for each channel and wavelength using NIRFASTer package.
- Parameters:
meshingparam – ff.utils.MeshingParam Parameters to be used by the CGAL mesher. Note:they should all be double
Returns: xr.DataArray
Fluence in each voxel for each channel and wavelength.
References
(Dehghani et al. [DEY+09]) Dehghani, Hamid, et al. “Near infrared optical tomography using NIRFAST: Algorithm for numerical model and image reconstruction.” Communications in numerical methods in engineering 25.6 (2009): 711-732.
- compute_sensitivity(fluence_all, fluence_at_optodes)
Compute sensitivity matrix from fluence.
- Parameters:
fluence_all – xr.DataArray Fluence in each voxel for each wavelength.
fluence_at_optodes – xr.DataArray Fluence at all optode positions for each wavelength.
- Returns:
- xr.DataArray
Sensitivity matrix for each channel, vertex and wavelength.
- compute_sensitivity_all(fluence_all, fluence_at_optodes)
Compute sensitivity matrix from fluence.
- Parameters:
fluence_all – xr.DataArray Fluence in each voxel for each wavelength.
fluence_at_optodes – xr.DataArray Fluence at all optode positions for each wavelength.
Returns: xr.DataArray
Sensitivity matrix for each channel, vertex and wavelength.
- static compute_stacked_sensitivity(sensitivity)
Compute stacked HbO and HbR sensitivity matrices from fluence.
- Parameters:
sensitivity (
DataArray
) – xr.DataArray Sensitivity matrix for each vertex and wavelength.- Returns:
- xr.DataArray
Stacked sensitivity matrix for each channel and vertex.
- class cedalion.imagereco.forward_model.TwoSurfaceHeadModel(segmentation_masks, brain, scalp, landmarks, t_ijk2ras, t_ras2ijk, voxel_to_vertex_brain, voxel_to_vertex_scalp)
Bases:
object
Head Model class to represent a segmented head.
Its main functions are reduced to work on voxel projections to scalp and cortex surfaces.
- segmentation_masks
xr.DataArray Segmentation masks of the head for each tissue type.
- brain
cdc.Surface Surface of the brain.
- scalp
cdc.Surface Surface of the scalp.
- landmarks
cdt.LabeledPointCloud Anatomical landmarks in RAS space.
- t_ijk2ras
cdt.AffineTransform Affine transformation from ijk to RAS space.
- t_ras2ijk
cdt.AffineTransform Affine transformation from RAS to ijk space.
- voxel_to_vertex_brain
scipy.sparse.spmatrix Mapping from voxel to brain vertices.
- voxel_to_vertex_scalp
scipy.sparse.spmatrix Mapping from voxel to scalp vertices.
- crs
str Coordinate reference system of the head model.
- 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.
- apply_transform(transform)
Apply a coordinate transformation to the head model.
- save(foldername)
Save the head model to a folder.
- load(foldername)
Load the head model from a folder.
- align_and_snap_to_scalp(points)
Align and snap optodes or points to the scalp surface.
- align_and_snap_to_scalp(points)
Align and snap optodes or points to the scalp surface.
- Parameters:
points (
DataArray
) – Points to be aligned and snapped to the scalp surface.- Return type:
DataArray
- Returns:
Points aligned and snapped to the scalp surface.
- apply_transform(transform)
Apply a coordinate transformation to the head model.
- Parameters:
transform (
DataArray
) – Affine transformation matrix (4x4) to be applied.- Return type:
- Returns:
Transformed head model.
- property crs
Coordinate reference system of the head model.
- classmethod from_segmentation(segmentation_dir, mask_files={'csf': 'csf.nii', 'gm': 'gm.nii', 'scalp': 'scalp.nii', 'skull': 'skull.nii', 'wm': 'wm.nii'}, landmarks_ras_file=None, brain_seg_types=['gm', 'wm'], scalp_seg_types=['scalp'], smoothing=0.5, brain_face_count=180000, scalp_face_count=60000, fill_holes=False)
Constructor from binary masks as gained from segmented MRI scans.
- Parameters:
segmentation_dir (str) – Folder containing the segmentation masks in NIFTI format.
mask_files (dict[str, str]) – Dictionary mapping segmentation types to NIFTI filenames.
landmarks_ras_file (Optional[str]) – Filename of the landmarks in RAS space.
brain_seg_types (list[str]) – List of segmentation types to be included in the brain surface.
scalp_seg_types (list[str]) – List of segmentation types to be included in the scalp surface.
smoothing (float) – Smoothing factor for the brain and scalp surfaces.
brain_face_count (Optional[int]) – Number of faces for the brain surface.
- Return type:
-
landmarks:
Annotated
[DataArray
]
- classmethod load(foldername)
Load the head model from a folder.
- Parameters:
foldername (
str
) – Folder to load the head model from.- Returns:
Loaded head model.
- save(foldername)
Save the head model to a folder.
- Parameters:
foldername (
str
) – Folder to save the head model into.- Returns:
None
-
segmentation_masks:
DataArray
-
t_ijk2ras:
DataArray
-
t_ras2ijk:
DataArray
-
voxel_to_vertex_brain:
spmatrix
-
voxel_to_vertex_scalp:
spmatrix
cedalion.imagereco.solver module
- cedalion.imagereco.solver.pseudo_inverse_stacked(Adot, alpha=0.01)
cedalion.imagereco.tissue_properties module
- class cedalion.imagereco.tissue_properties.TissueType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)
Bases:
Enum
- CSF = 4
- DM = 3
- GM = 5
- OTHER = 7
- SKIN = 1
- SKULL = 2
- WM = 6
- cedalion.imagereco.tissue_properties.get_tissue_properties(segmentation_masks)
Return tissue properties for the given segmentation mask.
- Return type:
ndarray
cedalion.imagereco.utils module
- cedalion.imagereco.utils.create_mock_activation_below_point(head_model, point, time_length, sampling_rate, spatial_size, vmax)
- cedalion.imagereco.utils.map_segmentation_mask_to_surface(segmentation_mask, transform_vox2ras, surface)
Find for each voxel the closest vertex on the surface.
- cedalion.imagereco.utils.normal_hrf(t, t_peak, t_std, vmax)