NIRSimple API
preprocessing
The preprocessing module contains functions to preprocess fNIRS signals.
- preprocessing.intensities_to_od_changes(intensities, refs=None)
Converts intensities into optical density changes. Changes are relative to the average intensity or a reference intensity for each channel.
Optical density from light intensity: OD = log10(I_0/I_t)
Optical density changes relative to average transmitted intensity: delta_OD = log10(I_0/I_t) - log10(I_0/I_average) delta_OD = -log10(I_t/I_average)
- Parameters:
intensities (array) – numpy array of absolute intensities, must have the correct shape (channels, data points).
refs (list of floats) – List of reference intensities to use instead of averages, length must be equal to the number of channels.
- Returns:
delta_od – numpy array of optical density changes, relative to average intensities or a reference for each channel, of shape (channels, data points).
- Return type:
array
- preprocessing.od_to_od_changes(optical_densities, refs=None)
Convert optical densities into optical density changes, relative to the average optical density or a reference optical density for each channel.
Optical density changes relative to average optical density: delta_OD = OD - OD_average
- Parameters:
optical_densities (array) – numpy array of optical densities, must have the correct shape (channels, data points).
refs (list of floats) – List of reference optical densities to use instead of averages, length must be equal to the number of channels.
- Returns:
delta_od – numpy array of optical density changes, relative to average optical densities or a reference for each channel, of shape (channels, data points).
- Return type:
array
- preprocessing.get_dpf(wavelength, age)
Get the differential pathlength factor (DPF) from wavelength and age using the general equation from Scholkmann & Wolf, 2013.
General DPF equation: DPF = (223.3 + 0.05624*age**0.8493 - 5.723e-7*wavelength**3 + 0.001245*wavelength**2 - 0.9025*wavelength)
- Parameters:
wavelength (float) – Wavelength in nm.
age (float) – Participant’s age in years.
- Returns:
dpf – Differential pathlength factor (DPF).
- Return type:
float
- preprocessing.mbll(delta_od, ch_names, ch_wls, ch_dpfs, ch_distances, unit, table='wray')
Apply the modified Beer-Lambert law (from Delpy et al., 1988) to optical density changes in order to obtain concentration changes in oxygenated hemoglobin (HbO) and deoxygenated hemoglobin (HbR).
Modified Beer-Lambert law: Dod_wl = e_HbO_wl*Dc_HbO*l*DPF_wl + e_HbR_wl*Dc_HbR*l*DPF_wl
With two different wavelengths we obtain a set of linear equations solved with matrices: [[Dod_1] = [[e_HbO_1*l*DPF_1, e_HbR_1*l*DPF_1] . [[Dc_HbO] [Dod_2]] [e_HbO_2*l*DPF_2, e_HbR_2*l*DPF_2]] [Dc_HbR]]
Equivalent to: [[Dc_HbO] = [[e_HbO_1, e_HbR_1] -1 . [[Dod_1/(l*DPF_1)] [DC_HbR]] [e_HbO_2, e_HbR_2]] [Dod_2/(l*DPF_2)]]
- Parameters:
delta_od (array) – numpy array of optical density changes, relative to average intensities for each channel, of shape (channels, data points).
ch_names (list of strings) – List of channel names.
ch_wls (list of integers) – List of channel wavelengths (in nm).
ch_dpfs (list of floats) – List of channel differential pathlength factors (DPF) (or partial pathlength factors (PPF)).
ch_distances (list of floats) – List of channel source-detector distances.
unit (string) – Unit for ch_distances (‘cm’ or ‘mm’).
table (string) – Table to use as molar extinction coefficients. ‘wray’: data from S. Wray et al., 1988 ‘cope’: data from M. Cope, 1991 ‘gratzer’: data from W.B. Gratzer and K. Kollias compiled by S. Prahl ‘moaveni’: data from M.K. Moaveni and J.M. Schmitt compiled by S. Prahl ‘takatani’: data from S. Takatani and M.D. Graham compiled by S. Prahl
- Returns:
delta_c (array) – numpy array of hemoglobin concentration changes in [moles/liter] or [M] for each channel, of shape (channels, data points).
new_ch_names (list of strings) – New list of channel names.
new_ch_types (list of strings) – New list of channel types (‘hbo’ for oxygenated hemoglobin and ‘hbr’ for deoxygenated hemoglobin).
processing
The processing module contains functions to process fNIRS signals.
- processing.cbsi(delta_c, ch_names, ch_types)
Apply correlation based signal improvement (from Cui at al., 2010) to hemoglobin concentration changes.
Correlation based signal improvement (CBSI): x_0 = (1/2)*(x-alpha*x) y_0 = -(1/alpha)*x_0
- Parameters:
delta_c (array) – numpy array of hemoglobin concentration changes in [moles/liter] or [M] for each channel, of shape (channels, data points).
ch_names (list of strings) – List of channel names.
ch_types (list of strings) – List of channel types (‘hbo’ for oxygenated hemoglobin and ‘hbr’ for deoxygenated hemoglobin).
- Returns:
delta_c_0 (array) – numpy array of corrected activation signals in [moles/liter] or [M] for each channel, of shape (channels, data points).
new_ch_names (list of strings) – New list of channel names.
new_ch_types (list of strings) – New list of channel types (‘hbo’ for oxygenated hemoglobin and ‘hbr’ for deoxygenated hemoglobin).