cedalion package
Subpackages
- cedalion.data package
- cedalion.dataclasses package
- Submodules
- cedalion.dataclasses.accessors module
- cedalion.dataclasses.geometry module
PointType
Surface
Voxels
TrimeshSurface
TrimeshSurface.mesh
TrimeshSurface.crs
TrimeshSurface.units
TrimeshSurface.mesh
TrimeshSurface.vertices
TrimeshSurface.nvertices
TrimeshSurface.nfaces
TrimeshSurface.apply_transform()
TrimeshSurface.decimate()
TrimeshSurface.smooth()
TrimeshSurface.get_vertex_normals()
TrimeshSurface.fix_vertex_normals()
TrimeshSurface.from_vtksurface()
VTKSurface
SimpleMesh
PycortexSurface
PycortexSurface.mesh
PycortexSurface.vertices
PycortexSurface.nvertices
PycortexSurface.nfaces
PycortexSurface.apply_transform()
PycortexSurface.decimate()
PycortexSurface.get_vertex_normals()
PycortexSurface.ppts
PycortexSurface.connected
PycortexSurface.adj
PycortexSurface.face_normals
PycortexSurface.vertex_normals
PycortexSurface.face_areas
PycortexSurface.cotangent_weights
PycortexSurface.laplace_operator
PycortexSurface.avg_edge_length
PycortexSurface.surface_gradient()
PycortexSurface.geodesic_distance()
PycortexSurface.geodesic_path()
PycortexSurface.from_trimeshsurface()
PycortexSurface.from_vtksurface()
affine_transform_from_numpy()
- cedalion.dataclasses.recording module
Recording
Recording.timeseries
Recording.masks
Recording.geo3d
Recording.geo2d
Recording.stim
Recording.aux_ts
Recording.aux_obj
Recording.head_model
Recording.meta_data
Recording.timeseries
Recording.masks
Recording.geo3d
Recording.geo2d
Recording.stim
Recording.aux_ts
Recording.aux_obj
Recording.head_model
Recording.meta_data
Recording.get_timeseries()
Recording.set_timeseries()
Recording.get_mask()
Recording.set_mask()
Recording.get_timeseries_type()
Recording.source_labels
Recording.detector_labels
Recording.wavelengths
Recording.trial_types
- cedalion.dataclasses.schemas module
- Module contents
- cedalion.geometry package
- cedalion.imagereco package
- Submodules
- cedalion.imagereco.forward_model module
TwoSurfaceHeadModel
TwoSurfaceHeadModel.segmentation_masks
TwoSurfaceHeadModel.brain
TwoSurfaceHeadModel.scalp
TwoSurfaceHeadModel.landmarks
TwoSurfaceHeadModel.t_ijk2ras
TwoSurfaceHeadModel.t_ras2ijk
TwoSurfaceHeadModel.voxel_to_vertex_brain
TwoSurfaceHeadModel.voxel_to_vertex_scalp
TwoSurfaceHeadModel.crs
TwoSurfaceHeadModel.apply_transform()
TwoSurfaceHeadModel.save()
TwoSurfaceHeadModel.load()
TwoSurfaceHeadModel.align_and_snap_to_scalp()
TwoSurfaceHeadModel.segmentation_masks
TwoSurfaceHeadModel.brain
TwoSurfaceHeadModel.scalp
TwoSurfaceHeadModel.landmarks
TwoSurfaceHeadModel.t_ijk2ras
TwoSurfaceHeadModel.t_ras2ijk
TwoSurfaceHeadModel.voxel_to_vertex_brain
TwoSurfaceHeadModel.voxel_to_vertex_scalp
TwoSurfaceHeadModel.from_segmentation()
TwoSurfaceHeadModel.from_surfaces()
TwoSurfaceHeadModel.crs
TwoSurfaceHeadModel.apply_transform()
TwoSurfaceHeadModel.save()
TwoSurfaceHeadModel.load()
TwoSurfaceHeadModel.align_and_snap_to_scalp()
TwoSurfaceHeadModel.snap_to_scalp_voxels()
ForwardModel
apply_inv_sensitivity()
- cedalion.imagereco.solver module
- cedalion.imagereco.tissue_properties module
- cedalion.imagereco.utils module
- Module contents
- cedalion.io package
- Submodules
- cedalion.io.anatomy module
- cedalion.io.bids module
- cedalion.io.forward_model module
- cedalion.io.nirs module
- cedalion.io.photogrammetry module
- cedalion.io.probe_geometry module
- cedalion.io.snirf module
DataType
DataType.CW_AMPLITUDE
DataType.CW_FLUORESCENCE_AMPLITUDE
DataType.FD_AC_AMPLITUDE
DataType.FD_PHASE
DataType.FD_FLUORESCENCE_AMPLITUDE
DataType.FD_FLUORESCENCE_PHASE
DataType.TDG_AMPLITUDE
DataType.TDG_FLUORESCENCE_AMPLITUDE
DataType.TDM_AMPLITUDE
DataType.TDM_FLUORESCENCE_AMPLITUDE
DataType.DCS_G2
DataType.DCS_BFI
DataType.PROCESSED
DataTypeLabel
DataTypeLabel.DOD
DataTypeLabel.DMEAN
DataTypeLabel.DVAR
DataTypeLabel.DSKEW
DataTypeLabel.MUA
DataTypeLabel.MUSP
DataTypeLabel.HBO
DataTypeLabel.HBR
DataTypeLabel.HBT
DataTypeLabel.H2O
DataTypeLabel.LIPID
DataTypeLabel.BFI
DataTypeLabel.HRF_DOD
DataTypeLabel.HRF_DMEAN
DataTypeLabel.HRF_DVAR
DataTypeLabel.HRF_DSKEW
DataTypeLabel.HRF_HBO
DataTypeLabel.HRF_HBR
DataTypeLabel.HRF_HBT
DataTypeLabel.HRF_BFI
DataTypeLabel.RAW_SATORI
DataTypeLabel.RAW_NIRX
parse_data_type()
parse_data_type_label()
reduce_ndim_sourceLabels()
labels_and_positions()
geometry_from_probe()
measurement_list_to_dataframe()
meta_data_tags_to_dict()
stim_to_dataframe()
read_aux()
add_number_to_name()
read_data_elements()
read_nirs_element()
read_snirf()
denormalize_measurement_list()
measurement_list_from_stacked()
write_snirf()
- Module contents
- cedalion.models package
- cedalion.sigdecomp package
- cedalion.sigproc package
- cedalion.sim package
- cedalion.vis package
Submodules
cedalion.datasets module
Cedalion datasets and utility functions.
- cedalion.datasets.get_ninja_uhd_cap_probe()[source]
Load the fullhead Ninja NIRS ultra HD cap probe.
- cedalion.datasets.get_colin27_segmentation(downsampled=False)[source]
Retrieves the Colin27 segmentation dataset, based on Holmes et al. [HHC+98].
- cedalion.datasets.get_colin27_headmodel()[source]
Retrieves the Colin27 headmodel, based on Holmes et al. [HHC+98].
- cedalion.datasets.get_fingertapping() Recording [source]
Retrieves a finger tapping recording in BIDS format.
Data is from Luke and McAlpine [LM21a]
- cedalion.datasets.get_fingertappingDOT() Recording [source]
Retrieves a finger tapping DOT example dataset from the IBS Lab.
- cedalion.datasets.get_precomputed_fluence(
- dataset: str,
- head_model: str,
Precomputed forward model results for examples and documentation.
- Parameters:
dataset – “fingertapping” or “fingertappingDOT”
head_model – “colin27” or “icbm152”
- Returns:
fluence_all, fluence_at_optodes
cedalion.errors module
Cedalion-specific exceptions.
cedalion.nirs module
Functions for preliminary processing of near-infrared spectroscopy (NIRS) data.
- cedalion.nirs.get_extinction_coefficients(spectrum: str, wavelengths: ArrayLike)[source]
Provide a matrix of extinction coefficients from tabulated data.
- Parameters:
spectrum –
The type of spectrum to use. Currently supported options are:
”prahl”: Extinction coefficients based on the Prahl absorption spectrum (Prahl1998).
wavelengths – An array-like object containing the wavelengths at which to calculate the extinction coefficients.
- Returns:
- A matrix of extinction coefficients with dimensions “chromo”
(chromophore, e.g. HbO/HbR) and “wavelength” (e.g. 750, 850, …) at which the coefficients for each chromophore are given in units of “mm^-1 / M”.
- Return type:
xr.DataArray
References
(Prahl 1998) - taken from Homer2/3, Copyright 2004 - 2006 - The General Hospital Corporation and President and Fellows of Harvard University.
“These values for the molar extinction coefficient e in [cm-1/(moles/liter)] were compiled by Scott Prahl (prahl@ece.ogi.edu) using data from W. B. Gratzer, Med. Res. Council Labs, Holly Hill, London N. Kollias, Wellman Laboratories, Harvard Medical School, Boston To convert this data to absorbance A, multiply by the molar concentration and the pathlength. For example, if x is the number of grams per liter and a 1 cm cuvette is being used, then the absorbance is given by
[(1/cm)/(moles/liter)] (x) [g/liter] (1) [cm]
- A = —————————————————
66,500 [g/mole]
using 66,500 as the gram molecular weight of hemoglobin. To convert this data to absorption coefficient in (cm-1), multiply by the molar concentration and 2.303, µa = (2.303) e (x g/liter)/(66,500 g Hb/mole) where x is the number of grams per liter. A typical value of x for whole blood is x=150 g Hb/liter.”
- cedalion.nirs.channel_distances(amplitudes: cdt.NDTimeSeries, geo3d: cdt.LabeledPointCloud)[source]
Calculate distances between channels.
- Parameters:
amplitudes – A DataArray representing the amplitudes with dimensions (channel, *).
geo3d (xr.DataArray) – A DataArray containing the 3D coordinates of the channels with dimensions (channel, pos).
- Returns:
- A DataArray containing the calculated distances between
source and detector channels. The resulting DataArray has the dimension ‘channel’.
- Return type:
dists (xr.DataArray)
- cedalion.nirs.int2od(amplitudes: cdt.NDTimeSeries)[source]
Calculate optical density from intensity amplitude data.
- Parameters:
amplitudes (xr.DataArray, (time, channel, *)) – amplitude data.
- Returns:
(xr.DataArray, (time, channel,*): The optical density data.
- Return type:
od
- cedalion.nirs.od2conc(
- od: cdt.NDTimeSeries,
- geo3d: cdt.LabeledPointCloud,
- dpf: DataArray,
- spectrum: str = 'prahl',
Calculate concentration changes from optical density data.
- Parameters:
od (xr.DataArray, (channel, wavelength, *)) – The optical density data array
geo3d (xr.DataArray) – The 3D coordinates of the optodes.
dpf (xr.DataArray, (wavelength, *)) – The differential pathlength factor data
spectrum (str, optional) – The type of spectrum to use for calculating extinction coefficients. Defaults to “prahl”.
- Returns:
A data array containing concentration changes by channel.
- Return type:
conc (xr.DataArray, (channel, *))
- cedalion.nirs.conc2od(
- conc: cdt.NDTimeSeries,
- geo3d: cdt.LabeledPointCloud,
- dpf: DataArray,
- spectrum: str = 'prahl',
Calculate optical density data from concentration changes.
- Parameters:
conc (xr.DataArray, (channel, *)) – The concentration changes by channel.
geo3d (xr.DataArray) – The 3D coordinates of the optodes.
dpf (xr.DataArray, (wavelength, *)) – The differential pathlength factor data.
spectrum (str, optional) – The type of spectrum to use for calculating extinction coefficients. Defaults to “prahl”.
- Returns:
- A data array containing
optical density data.
- Return type:
od (xr.DataArray, (channel, wavelength, *))
- cedalion.nirs.beer_lambert(
- amplitudes: cdt.NDTimeSeries,
- geo3d: cdt.LabeledPointCloud,
- dpf: DataArray,
- spectrum: str = 'prahl',
Calculate concentration changes from amplitude using the modified BL law.
- Parameters:
amplitudes (xr.DataArray, (channel, wavelength, *)) – The input data array containing the raw intensities.
geo3d (xr.DataArray) – The 3D coordinates of the optodes.
dpf (xr.DataArray, (wavelength,*)) – The differential pathlength factors
spectrum (str, optional) – The type of spectrum to use for calculating extinction coefficients. Defaults to “prahl”.
- Returns:
- A data array containing
concentration changes according to the mBLL.
- Return type:
conc (xr.DataArray, (channel, *))
- cedalion.nirs.split_long_short_channels(
- ts: cdt.NDTimeSeries,
- geo3d: cdt.LabeledPointCloud,
- distance_threshold: cdt.QLength = <Quantity(1.5,
- 'centimeter')>,
Split a time series into two based on channel distances.
- Parameters:
ts (cdt.NDTimeSeries) – Time series to split.
geo3d (cdt.LabeledPointCloud) – 3D coordinates of the channels.
distance_threshold (Quantity) – Distance threshold for splitting the channels.
- Returns:
time series with channel distances >= distance_threshold ts_short : time series with channel distances < distance_threshold
- Return type:
ts_long
cedalion.physunits module
Builds on pint_xarray’s unit registry.
cedalion.plots module
Plotting functions for visualization of montages, meshes, etc.
- cedalion.plots.plot_montage3D(
- amp: DataArray,
- geo3d: DataArray,
Plots a 3D visualization of a montage.
- Parameters:
amp (xr.DataArray) – Time series data array.
geo3d (xr.DataArray) – Landmark coordinates.
- cedalion.plots.plot3d(
- brain_mesh,
- scalp_mesh,
- geo3d,
- timeseries,
- poly_lines=[],
- brain_scalars=None,
- plotter=None,
Plots a 3D visualization of brain and scalp meshes.
- Parameters:
brain_mesh (TrimeshSurface) – The brain mesh as a TrimeshSurface object.
scalp_mesh (TrimeshSurface) – The scalp mesh as a TrimeshSurface object.
geo3d (xarray.Dataset) – Dataset containing 3-dimentional point centers.
timeseries – Time series data array.
poly_lines – List of lists of points to be plotted as polylines.
brain_scalars – Scalars to be used for coloring the brain mesh.
plotter (pv.Plotter, optional) – An existing PyVista plotter instance to use for plotting. If None, a new PyVista plotter instance is created. Default: None.
- Initial Contributors:
Eike Middell | middell@tu-berlin.de | 2024
- cedalion.plots.plot_surface(
- plotter: pv.Plotter,
- surface: cdc.Surface,
- color: pv.ColorLike | None = None,
- opacity: float = 1.0,
- pick_landmarks: bool = False,
- **kwargs,
Plots a surface mesh with optional landmark picking in a PyVista plotter.
- Parameters:
plotter – A PyVista plotter instance used for rendering the surface.
surface – The surface object to be plotted.
color – Color of the mesh.
opacity – Opacity of the mesh, ranging from 0 (transparent) to 1 (opaque). Default is 1.0.
pick_landmarks – If True, enables interactive picking of landmarks on the surface. Default is False.
**kwargs – Additional keyword arguments are passed to pv.add_mesh.
- Returns:
If pick_landmarks is True, returns a function that when called, provides the current picked points and their labels. This function prints warnings if some labels are missing or are repeated.
- Return type:
function
- Initial Contributors:
Eike Middell | middell@tu-berlin.de | 2024
Masha Iudina | mashayudi@gmail.com | 2024
- cedalion.plots.plot_labeled_points(
- plotter: pv.Plotter,
- points: cdt.LabeledPointCloud,
- color: pv.ColorLike = None,
- show_labels: bool = False,
- ppoints: bool = None,
- labels: list[str] | None = None,
Plots a labeled point cloud with optional interaction for picking points.
This function visualizes a point cloud where each point can have a label. Points can be interactively picked if enabled. Picked point is indicated by increased radius.
- Parameters:
plotter – A PyVista plotter instance used for rendering the points.
points – A labeled point cloud data structure containing points and optional labels.
color – Override color for all points. If None, colors are assigned based on point types.
show_labels – If True, labels are displayed next to the points.
ppoints – A list to store indices of picked points, enables picking if not None.
labels – List of labels to show if show_labels is True. If None and show_labels is True, the labels from points are used.
- Initial Contributors:
Eike Middell | middell@tu-berlin.de | 2024
- cedalion.plots.plot_vector_field(
- plotter: Plotter,
- points: cdt.LabeledPointCloud,
- vectors: DataArray,
- ppoints=None,
Plots a vector field on a PyVista plotter.
- Parameters:
plotter (pv.Plotter) – A PyVista plotter instance used for rendering the vector field.
points (cdt.LabeledPointCloud) – A labeled point cloud data structure containing point coordinates.
vectors (xr.DataArray) – A data array containing the vector field.
ppoints (list, optional) – A list to store indices of picked points, enables picking if not None. Default is None.
- class cedalion.plots.OptodeSelector(surface, points, normals=None, plotter=None, labels=None)[source]
Bases:
object
A class for visualizing point clouds with interactive features in PyVista.
This class provides functionality to visualize and interact with labeled point clouds using a PyVista plotter. It allows points to be dynamically added or removed by picking them directly from the plot interface.
- labels[source]
Labels corresponding to the points, displayed if provided.
- Type:
list of str, optional
- actors[source]
List of PyVista actor objects representing the points in the visualization.
- Type:
list
- color[source]
Default color for points if not specified by point type.
- Type:
str or tuple, optional
- Initial Contributors:
Masha Iudina | mashayudi@gmail.com | 2024
- cedalion.plots.plot_stim_markers(
- ax,
- stim: DataFrame,
- fmt: dict[str, dict] | None = None,
- y: float = 0.03,
Add stimulus indicators to an Axes.
For each trial a Rectangle is plotted in x from onset to onset+duration. The height of the rectangle is specified in axes coordinates. In the default setting a small bar at bottom of the axes is drawn. By setting y to 1. the stimulus marker covers the full height of the axes.
- Parameters:
ax – the matplotlib axes to operate on
stim – a stimulas data frame
fmt – for each trial_type a dictioniary of keyword arguments can be provided. These kwargs are passed to matplotlib.patches.Rectangle to format the stimulus indicator.
y – the height of the Rectangle in axes coordinates.
- Initial Contributors:
Eike Middell | middell@tu-berlin.de | 2024
- cedalion.plots.scalp_plot(
- ts: cdt.NDTimeSeries,
- geo3d: cdt.LabeledPointCloud,
- metric: xr.DataArray | ArrayLike,
- ax,
- title: str | None = None,
- vmin: float | None = None,
- vmax: float | None = None,
- cmap: str | matplotlib.colors.Colormap = 'bwr',
- norm: object | None = None,
- bad_color: ColorType = [0.7, 0.7, 0.7],
- min_dist: Quantity | None = None,
- min_metric: float | None = None,
- channel_lw: float = 2.0,
- optode_size: float = 36.0,
- optode_labels: bool = False,
- cb_label: str | None = None,
- cb_ticks_labels: list[float, str] | None = None,
- add_colorbar: bool = True,
- zorder: str | None = None,
Creates a 2D plot of the head with channels coloured according to a given metric.
- Parameters:
ts – a NDTimeSeries to provide channel definitions
geo3d – a LabeledPointCloud to provide the probe geometry
metric ((
DataArray
, (channel,) | ArrayLike)) – the scalar metric to be plotted for each channel. If provided as a DataArray it needs a channel dimension. If provided as a plain array or list it must have the same length as ts.channel and the matching is done by position.ax – the matplotlib.Axes object into which to draw
title – the axes title
vmin – the minimum value of the metric
vmax – the maximum value of the metric
cmap – the name of the colormap
bad_color – the color to use when the metric contains NaNs
min_dist – if provided channels below this distance threshold are not drawn
min_metric – if provided channels below this metric threshold are toned down
channel_lw – channel line width
optode_size – optode marker size
optode_labels – if True draw optode labels instead of markers
cb_label – colorbar label
zorder – ‘ascending’ or ‘descending’ or None. Controls whether channels with high or low metric values are plotted on top.
- Initial Contributors:
Laura Carlton | lcarlton@bu.edu | 2024
Eike Middell | middell@tu-berlin.de | 2024
- cedalion.plots.scalp_plot_gif(
- data_ts: cdt.NDTimeSeries,
- geo3d: cdt.LabeledPointCloud,
- filename: str,
- time_range: tuple = None,
- cmap: str | Colormap = 'seismic',
- scl=None,
- fps: int = 10,
- optode_size: float = 6,
- optode_labels: bool = False,
- str_title: str = '',
Generate a GIF of scalp topographies over time from time-series data.
- Parameters:
data_ts – xarray.DataArray A 2D DataArray with dimensions (channel, time). Must include coordinate labels for ‘source’ and ‘detector’ in the ‘channel’ dimension.
geo3d – cedalion.core.LabeledPointCloud 3D geometry object defining optode locations for projecting onto the scalp surface.
filename – str Full path to the output GIF file without file extension.
time_range – tuple, optional Provides (start_time, stop_time, step_time) in quantity ‘s’ for generating animation.
cmap – string, optional A matplotlib colormap name or a Colormap object. Default is ‘seismic’.
scl – tuple of (float, float), optional Tuple defining the (vmin, vmax) for the color scale. If None, the color scale is set to ± the maximum absolute value of the data.
fps – int, optional Frames per second for the output GIF. Default is 10.
optode_size – float, optional Size of optode markers on the plot. Default is 6.
optode_labels – bool, optional Whether to show text labels for optodes instead of markers. Default is False.
str_title – str, optional Extra string to append to the title of each frame.
Returns: None
The function saves a GIF file to the specified location.
Initial Contributors: - David Boas | dboas@bu.edu | 2025 - Alexander von Lühmann | vonluehmann@tu-berlin.de | 2025
- cedalion.plots.image_recon(
- X: cdt.NDTimeSeries,
- head: TwoSurfaceHeadModel,
- cmap: str | Colormap = 'seismic',
- clim=None,
- view_type: str = 'hbo_brain',
- view_position: str = 'superior',
- p0=None,
- title_str: str = None,
- off_screen: bool = False,
- plotshape=(1, 1),
- iax=(0, 0),
- show_scalar_bar: bool = False,
- wdw_size: tuple = (1024, 768),
Render a single frame of brain or scalp activity on a specified view.
This function creates (or reuses) a PyVista plotter, applies a custom colormap, sets the camera view according to the given view_position, adds the surface mesh with the scalar data (extracted from X), and returns the plotter, the mesh, and a text label.
- Parameters:
X – cdt.NDTimeSeries (or similar) Scalar data for the current frame. Expected to have a boolean attribute is_brain indicating brain vs. non-brain vertices, and HbO / HbR chromophore dimension
head – TwoSurfaceHeadModel A head model containing attributes such as head.brain and head.scalp.
cmap – str or matplotlib.colors.Colormap, default ‘seismic’ The colormap to use.
clim – tuple, optional Color limits. If None, they are computed from the data.
view_type – str, default ‘hbo_brain’ Indicates whether to plot brain (‘hbo_brain’ or ‘hbr_brain’) or scalp (‘hbo_scalp’ or ‘hbr_scalp’) data.
view_position – str, default ‘superior’ The view direction. Options are: ‘superior’, ‘anterior’, ‘posterior’,’left’, ‘right’, and ‘scale_bar’.
p0 – PyVista Plotter instance, optional If provided the mesh is added to this plotter; else a new plotter is created.
title_str – str, optional Title to use on the scalar bar.
off_screen – bool, default False Whether to use off-screen rendering.
plotshape – tuple, default (1, 1) The subplot grid shape.
iax – tuple, default (0, 0) The target subplot index (row, col).
show_scalar_bar – bool, optional Flag to control scalar bar visibility
wdw_size – tuple, default (1024, 768) The window size for the plotter (the plot resolution)
- Returns:
p0: the PyVista Plotter instance.
surf: the wrapped surface mesh (a pyvista mesh).
surf_label: a text actor (e.g., the scalar bar label).
- Return type:
A tuple (p0, surf, surf_label) where
Initial Contributors: - David Boas | dboas@bu.edu | 2025 - Laura Carlton | lcarlton@bu.edu | 2025 - Alexander von Lühmann | vonluehmann@tu-berlin.de | 2025
- cedalion.plots.image_recon_view(
- X_ts: cdt.NDTimeSeries,
- head: TwoSurfaceHeadModel,
- cmap: str | Colormap = 'seismic',
- clim=None,
- view_type: str = 'hbo_brain',
- view_position: str = 'superior',
- title_str: str = None,
- filename: str = None,
- SAVE: bool = False,
- time_range: tuple = None,
- fps: int = 6,
- geo3d_plot: cdt.LabeledPointCloud = None,
- wdw_size: tuple = (1024, 768),
Generate a single-view visualization of head activity.
For static data (2D: vertex × channel) the function can display (or save) a single frame. For time series data (3D: vertex × channel × time) the function can create an animated GIF by looping over the specified frame indices.
- Parameters:
X_ts – xarray.DataArray or NDTimeSeries Activity data. If 2D, a single static frame is plotted; if 3D, a time series is used. Expected to have a boolean attribute is_brain indicating brain vs. non-brain vertices, and HbO / HbR chromophore dimension
head – TwoSurfaceHeadModel The head mesh data to plot activity on.
cmap – str or matplotlib.colors.Colormap, default ‘seismic’ The colormap to use.
view_position – str, default ‘superior’ The view to render.
clim – tuple, optional Color limits. If None, they are computed from the data.
view_type – str, default ‘hbo_brain’ Indicates whether to plot brain (‘hbo_brain’ or ‘hbr_brain’) or scalp (‘hbo_scalp’ or ‘hbr_scalp’) data.
view_position – str, default ‘superior’ The view direction. Options are: ‘superior’, ‘anterior’, ‘posterior’,’left’, ‘right’, and ‘scale_bar’.
title_str – str, optional Title to use on the scalar bar.
filename – str, optional The output filename (without extension) for saving the image/GIF.
SAVE – bool, default False If True, the resulting still image is saved, otherwise only shown. Rendered gifs are always saved.
time_range – tuple, optional Provides (start_time, stop_time, step_time) in quantity ‘s’ for generating animation.
fps – int, default 6 Frames per second for the GIF.
geo3d_plot – cdt.LabeledPointCloud, optional A 3D point cloud for plotting labeled points (e.g. optodes) on the mesh.
wdw_size – tuple, default (1024, 768) The window size for the plotter (the plot resolution)
Returns: Nothing
Initial Contributors: - David Boas | dboas@bu.edu | 2025 - Laura Carlton | lcarlton@bu.edu | 2025 - Alexander von Lühmann | vonluehmann@tu-berlin.de | 2025
- cedalion.plots.image_recon_multi_view(
- X_ts: cdt.NDTimeSeries,
- head: TwoSurfaceHeadModel,
- cmap: str | Colormap = 'seismic',
- clim=None,
- view_type: str = 'hbo_brain',
- title_str: str = None,
- filename: str = None,
- SAVE: bool = True,
- time_range: tuple = None,
- fps: int = 6,
- geo3d_plot: cdt.LabeledPointCloud = None,
- wdw_size: tuple = (1024, 768),
Generate a multi-view (2×3 grid) visualization of head activity across different views.
For static data (2D: vertex × channel) the function can display (or save) a single frame. For time series data (3D: vertex × channel × time) the function creates an animated GIF where each frame updates all views.
- Parameters:
X_ts – xarray.DataArray or NDTimeSeries Activity data. If 2D, a single static frame is plotted; if 3D, a time series is used. Expected to have a boolean attribute is_brain indicating brain vs. non-brain vertices, and HbO / HbR chromophore dimension
head – TwoSurfaceHeadModel The head mesh data to plot activity on.
cmap – str or matplotlib.colors.Colormap, default ‘seismic’ The colormap to use.
view_position – str, default ‘superior’ The view to render.
clim – tuple, optional Color limits. If None, they are computed from the data.
view_type – str, default ‘hbo_brain’ Indicates whether to plot brain (‘hbo_brain’ or ‘hbr_brain’) or scalp (‘hbo_scalp’ or ‘hbr_scalp’) data.
title_str – str, optional Title to use on the scalar bar.
filename – str, optional The output filename (without extension) for saving the image/GIF.
SAVE – bool, default False If True, the resulting still image is saved, otherwise only shown. Rendered gifs are always saved.
time_range – tuple, optional Provides (start_time, stop_time, step_time) in quantity ‘s’ for generating animation.
fps – int, default 6 Frames per second for the GIF.
geo3d_plot – cdt.LabeledPointCloud, optional A 3D point cloud for plotting labeled points (e.g. optodes) on the mesh.
wdw_size – tuple, default (1024, 768) The window size for the plotter (the plot resolution)
Returns: Nothing
Initial Contributors: - David Boas | dboas@bu.edu | 2025 - Laura Carlton | lcarlton@bu.edu | 2025 - Alexander von Lühmann | vonluehmann@tu-berlin.de | 2025
cedalion.tasks module
cedalion.typing module
Type aliases for Cedalion dataclasses.
Cedalion relies as much as possible on generic data types (like xarray DataArrays). We then use type aliases and annotations to augment these data types with additional information about the data they carry. For DataArrays there is a basic mechanism to specify and validate data schemas that specify dimension and coordinate names. This way we can distinguish between time series DataArrays (NDTimeSeries) and DataArrays representing points in space (LabeledPointCloud). By using these aliases in type hints we indicate to user which kind of DataArray is expected.
Parameters with physical units are represented by cedalion.Quantity. Aliases are defined to indicate the dimensionality of quantities.
- cedalion.typing.LabeledPointCloud[source]
DataArrays representing labeled points in space.
alias of
Annotated
[DataArray
,DataArraySchema
(dims=(‘label’,), coords=((‘label’, (‘label’, ‘type’)),))]
- cedalion.typing.NDTimeSeries[source]
DataArrays representing time series.
alias of
Annotated
[DataArray
,DataArraySchema
(dims=(‘time’,), coords=((‘time’, (‘time’, ‘samples’)),))]
cedalion.utils module
Utility functions.
cedalion.validators module
Validators for common data structures.
cedalion.vtktutils module
cedalion.xrutils module
Utility functions for xarray objects.
- cedalion.xrutils.pinv(array: DataArray) DataArray [source]
Calculate the pseudoinverse of a 2D xr.DataArray.
- FIXME: handles unitless and quantified DataArrays but not
DataArrays with units in their attrs.
- Parameters:
array (xr.DataArray) – Input array
- Returns:
Pseudoinverse of the input array
- Return type:
array_inv (xr.DataArray)
- cedalion.xrutils.norm(
- array: DataArray,
- dim: str,
Calculate the vector norm along a given dimension.
Extends the behavior of numpy.linalg.norm to xarray DataArrays.
- Parameters:
array (xr.DataArray) – Input array
dim (str) – Dimension along which to calculate the norm
- Returns:
Array with the norm along the specified dimension
- Return type:
normed (xr.DataArray)
- cedalion.xrutils.mask(
- array: DataArray,
- initval: bool,
Create a boolean mask array with the same shape as the input array.
- cedalion.xrutils.apply_mask(
- data_array: DataArray,
- mask: DataArray,
- operator: str,
- dim_collapse: str,
Apply a boolean mask to a DataArray according to the defined “operator”.
- Parameters:
data_array – NDTimeSeries, input time series data xarray
mask – input boolean mask array with a subset of dimensions matching data_array
operator – operators to apply to the mask and data_array “nan”: inserts NaNs in the data_array where mask is False “drop”: drops value in the data_array where mask is False
dim_collapse – Mask dimension to collapse to, merging boolean masks along all other dimensions. Can be skipped with “none”. Example: collapsing to “channel” dimension will drop or nan a channel if it is “False” along any other dimensions
- Returns:
Input data_array with applied mask masked_elements: List of elements in data_array that were masked (e.g.
dropped or set to NaN)
- Return type:
masked_data_array
- cedalion.xrutils.convolve(
- data_array: DataArray,
- kernel: ndarray,
- dim: str,
Convolve a DataArray along a given dimension “dim” with a “kernel”.
- cedalion.xrutils.other_dim(data_array: DataArray, *dims: str) str [source]
Get the dimension name not listed in *dims.
Checks that there is only one more dimension than given in dims and returns its name.
- Parameters:
data_array – an xr.DataArray
*dims – names of dimensions
- Returns:
The name of the dimension of data_array.
- cedalion.xrutils.coords_from_other(
- source: DataArray,
- dims: list[str] = None,
- **aux_coords,
Create a dictionary of coordinates from source for matching dims in target.
- Parameters:
source – the DataArray to copy the coordinates from.
dims – a list of dimensions names. If specified, copy only coords for those dims.
aux_coords – additional key-value pairs to add to the resulting coords dict.
- Returns:
A dictionary that can be passed to DataArray.assign_coords.