Module description¶
Fit lifetime decays
-
class
lifefit.tcspc.
Anisotropy
(VV, VH, HV, HH)[source]¶ Bases:
object
Create an Anisotropy object with four polarization resolved lifetime decays
Parameters: - VV (ndarray) – vertical excitation - vertical emission
- VH (ndarray) – vertical excitation - horizontal emission
- HV (ndarray) – horizontal excitation - vertical emission
- HH (ndarray) – horizontal excitation - horizontal emission
Example
>>> lf.tcspc.Anisotropy(decay['VV'], decay['VH'], decay['HV'],decay['HH'])
-
static
G_factor
(HV, HH)[source]¶ Compute G-factor to correct for differences in transmittion effiency of the horizontal and vertical polarized light
Parameters: - HV (ndarray) – horizontal excitation - vertical emission
- HH (ndarray) – horizontal excitation - horizontal emission
Returns: G (float) – G-factor
Notes
The G-factor is defined as follows:
\[G = \frac{\int HV}{\int HH}\]
-
static
aniso_decay
(VV, VH, G)[source]¶ Parameters: - VV (ndarray) – vertical excitation - vertical emission
- VH (ndarray) – vertical excitation - horizontal emission
- G (float) – G-factor
Returns: r (ndarray) – anisotropy decay
Notes
The anisotropy decay is calculated from the parallel and perpendicular lifetime decays as follows:
\[r(t) = \frac{I_\text{VV} - GI_\text{VH}}{I_\text{VV} + 2GI_\text{VH}}\]
-
export
(filename)[source]¶ Export the data and the fit parameters to a json file
Parameters: filename (str) –
-
static
hindered_rotation
(time, r0, tau_r, r_inf)[source]¶ Hindered rotation in-a-cone model
Parameters: - time (array_like) – time bins
- r0 (float) – fundamental anisotropy
- tau_r (float) – rotational correlation time
- r_inf (float) – residual anisotropy at time->inf
Returns: ndarray – hindered rotation anisotropy decay
-
static
local_global_rotation
(time, r0, tau_rloc, r_inf, tau_rglob)[source]¶ Local-global rotation in-a-cone model
Parameters: - time (array_like) – time bins
- r0 (float) – fundamental anisotropy
- tau_rloc (float) – local rotational correlation time
- r_inf (float) – residual anisotropy at time->inf
- tau_rglob (float) – global rotational correlation time
Returns: ndarray – local-global rotation anisotropy decay
-
static
one_rotation
(time, r0, tau)[source]¶ Single rotator model
Parameters: - time (array_like) – time bins
- r0 (float) – fundamental anisotropy
- tau_r (float) – rotational correlation time
Returns: ndarray – two-rotation anisotropy decay
-
rotation_fit
(p0=[0.4, 1], model='one_rotation', manual_interval=None, bounds=(0, inf), verbose=True, ns_before_VVmax=1, signal_percentage=0.01)[source]¶ Fit rotation model to anisotropy decay.
Parameters: - p0 (array_like) – start values of the chosen anisotropy fit model
- model (str) – one of the following anisotropy models: {‘one_rotation’, ‘two_rotations’, ‘hindered_rotation’, ‘local_global_rotation’}
- manual_interval (2-tuple of float, optional) –
- bounds (2-tuple of float or array_like) – lower and upper bounds for each parameter in p0. Can be either a tuple of two scalars (same bound for all parameters) or a tuple of array_like with the same length as p0. To deactivate parameter bounds set: bounds=(-np.inf, np.inf)
- verbose (bool) – print anisotropy fit result
- ns_before_VVmax (float, optional) – how many nanoseconds before the maximum of the VV decay should the search for r0 start
- signal_percentage (float, optional) – percentage of the VV decay serving as a treshold to define the end of the anisotropy fit interval
Example
>>> obj.rotation_fit(p0=[0.4, 1, 10, 1], model='two_rotations')
-
static
two_rotations
(time, r0, b, tau_r1, tau_r2)[source]¶ Two-rotator model
Parameters: - time (array_like) – time bins
- r0 (float) – fundamental anisotropy
- b (float) – amplitude of second decay
- tau_r1 (float) – first rotational correlation time
- tau_r2 (float) – second rotational correlation time
Returns: ndarray – two-rotation anisotropy decay
-
class
lifefit.tcspc.
Lifetime
(fluor_decay, fluor_ns_per_chan, irf_decay=None, gauss_sigma=None, gauss_amp=None, shift_time=False)[source]¶ Bases:
object
Create lifetime class
Parameters: - fluor_decay (ndarray) – n x 2 array containing numbered channels and intensity counts of the fluorescence decay
- fluor_ns_per_chan (float) – nanoseconds per channel
- irf_decay (ndarray, optional) – n x 2 array containing numbered channels and intensity counts for instrument reponse function (IRF). If None, then IRF is approximated by a Gaussian
Variables: - ns_per_chan (float) – nanoseconds per channel
- fluor (ndarray) – n x 4 array containing time, channel number, intensity counts and associated Poissonian weights of the fluorescence decay
- irf (ndarray) – n x 3 array containing time, channel number and intensity counts of the IRF
- irf_type (str) – type of IRF: {‘Gaussian’, ‘experimental’}
- fit_param (ndarray) –
- fit_param_std (ndarray) –
Example
>>> fluor, fluor_nsperchan = lf.tcspc.read_decay(pathToFluorDecay) >>> irf, irf_nsperchan = lf.tcspc.read_decay(pathToIRF) >>> lf.tcspc.Lifetime(fluor, fluor_nsperchan, irf)
-
static
average_lifetime
(a, tau_val, tau_std)[source]¶ Calculate average lifetime according to [1]
Parameters: - a (array_like) – weighting factors of tau
- tau_val (array_like) – fluorescence lifetimes
- tau_std (array_like) – standard deviation of the fluorescence lifetimes
Returns: av_lt (tuple) – average lifetime and associated standard deviation
References
[1] - Lakowicz, Principles of Fluorescence, 3rd ed., Springer, 2010.
-
static
convolution
(irf, sgl_exp)[source]¶ Compute convolution of irf with a single exponential decay
Parameters: - irf (array_like) – intensity counts of the instrument reponse function (experimental of Gaussian shaped)
- sgl_exp (array_like) – single-exponential decay
Returns: convolved (ndarray) – convoluted signal of IRF and exponential decay
-
static
exp_decay
(time, tau)[source]¶ Create a single-exponential decay
Parameters: - time (array_like) – time bins
- tau (float) – fluorescence lifetime
Returns: sgl_exp (array_like) – single-exponential decay
-
classmethod
from_filenames
(fluor_file, irf_file=None, fileformat='HORIBA', gauss_sigma=None, gauss_amp=None, shift_time=False)[source]¶ Alternative constructor for the Lifetime class by reading in filename for the fluorophore and IRF decay
Parameters: - fluor_file (str) – filename of the fluorophore decay
- irf_file (str) – filename of the IRF decay
- fileformat (str, optional) – currently implemented formats: {‘HORIBA’}
Example
>>> lf.tcspc.Lifetime.from_filenames(pathToFluorDecay, pathToIRFDecay)
-
static
gauss_irf
(time, mu, sigma=0.01, A=10000)[source]¶ Calculate a Gaussian-shaped instrument response function (IRF)
Parameters: - time (ndarray) – time bins
- mu (float) – mean of the Gaussian distribution
- sigma (float, optional) – standard deviation of the Gaussian distribution
- A (float, optional) – amplitude of the Gaussian distribution
Returns: irf (ndarray) – Gaussian shaped instrument response function (IRF)
-
nnls_convol_irfexp
(x_data, p0)[source]¶ Solve non-negative least squares for series of IRF-convolved single-exponential decays. First, the IRF is shifted, then convolved with each exponential decay individually (decays 1,…,n), merged into an m x n array (=A) and finally plugged into scipy.optimize.nnls(A, experimental y-data) to compute argmin_x || Ax - y ||_2. This optimizes the relative weight of the exponential decays whereas the curve_fit function optimizes the decay parameters (tau1, taus2, etc.)
Parameters: - x_data (array_like) – array of the independent variable
- p0 (array_like) – start values for the fit model
Returns: - A (ndarray) – matrix containing irf-convoluted single-exponential decays in the first n columns and ones in the last column (background counts)
- x (ndarray) – vector that minimizes || Ax - y ||_2
- y (ndarray) – fit vector computed as y = Ax
-
reconvolution_fit
(tau0=[1], tau_bounds=(0, inf), irf_shift=0, sigma=None, verbose=True)[source]¶ Fit the experimental lifetime decay to a series of exponentials via interative reconvolution with the instrument reponse function (IRF).
Parameters: - tau0 (int or array_like) – start value(s) of the fluorescence lifetime(s)
- tau_bounds (2-tuple of float or 2-tuple of array_like, optional) – lower and upper bounds for each parameter in tau0. Can be either a tuple of two scalars (same bound for all parameters) or a tuple of array_like with the same length as tau0. To deactivate parameter bounds set: bounds=(-np.inf, np.inf)
- irf_shift (int, optional) – shift of the IRF on the time axis (in channel units)
- sigma (array_like , optional) – uncertainty of the decay (same length as y_data)
- verbose (bool, optional) – print lifetime fit result
Example
>>> obj.reconvolution_fit([1,5])
-
lifefit.tcspc.
fit
(fun, x_data, y_data, p0, bounds=([0, 0, 0], [inf, inf, inf]), sigma=None)[source]¶ Wrapper for the curve_fit function of the scipy.optimize module The curve_fit optimizes the decay parameters (tau1, tau2, etc.) while the nnls weights the exponential decays.
Parameters: - fun (callable) – The model function f(x,…) taking x values as a first argument followed by the function parameters
- x_data (array_like) – array of the independent variable
- y_data (array_like) – array of the dependent variable
- p0 (array_like) – start values for the fit model
- bounds (2-tuple of float or 2-tuple of array_like, optional) – lower and upper bounds for each parameter in p0. Can be either a tuple of two scalars (same bound for all parameters) or a tuple of array_like with the same length as p0. To deactivate parameter bounds set: bounds=(-np.inf, np.inf)
- sigma (array_like, optional) – uncertainty of the decay (same length as y_data)
Returns: - p (ndarray) – optimized fit parameters
- p_std (ndarray) – standard deviation of optimized fit parameters
-
lifefit.tcspc.
parseCmd
()[source]¶ Parse the command line to get the experimental decay and instrument reponse function (IRF) file.
Returns: - fluor_file (str) – filename of the fluorescence decay
- irf_file (str) – filename of the IRF (if None then the IRF is approximated by a Gaussian)
-
lifefit.tcspc.
parse_file
(decay_file, fileformat='Horiba')[source]¶ Parse the decay file
Parameters: - decay_file (StringIO) –
- fileformat (str, optional) – currently implemented formats: {‘HORIBA’}
Returns: - decay_data (ndarray) – n x 2 decay containing numbered channels and intensity counts for instrument reponse function (IRF)
- ns_per_chan (float)
-
lifefit.tcspc.
read_decay
(filepath_or_buffer, fileformat='Horiba')[source]¶ Read TCSPC decay file from HORIBA or another data format
Parameters: - filepath_or_buffer (str, os.PathLike, StringIO) – filename of the decay or StringIO object
- fileformat (str, optional) – currently implemented formats: {‘HORIBA’}
Returns: - decay_data (ndarray) – n x 2 decay containing numbered channels and intensity counts for instrument response function (IRF)
- ns_per_chan (float)