Code reference

module-level methods

dclab.new_dataset(data, identifier=None)[source]

Initialize a new RT-DC dataset

Parameters:
  • data

    can be one of the following:

    • dict
    • .tdms file
    • .rtdc file
    • subclass of RTDCBase (will create a hierarchy child)
  • identifier (str) – A unique identifier for this dataset. If set to None an identifier is generated.
Returns:

dataset – A new dataset instance

Return type:

subclass of dclab.rtdc_dataset.RTDCBase

global definitions

These definitionas are used throughout the dclab/Shape-In/Shape-Out ecosystem.

configuration

Valid configuration sections and keys are described in: Analysis metadata and Experiment metadata.

dclab.dfn.CFG_ANALYSIS

User-editable configuration for data analysis.

dclab.dfn.CFG_METADATA

Measurement-specific metadata.

dclab.dfn.config_funcs

Dictionary of dictionaries containing functions to convert input data to the predefined data type

dclab.dfn.config_keys

Dictionary with sections as keys and configuration parameter names as values

dclab.dfn.config_types

Dictionary of dictionaries containing the data type of each configuration parameter

features

Features are discussed in more detail in: Features.

dclab.dfn.FEATURES_SCALAR

Scalar features

dclab.dfn.FEATURES_NON_SCALAR

Non-scalar features

dclab.dfn.feature_names

List of valid feature names

dclab.dfn.feature_labels

List of human-readable labels for each valid feature

dclab.dfn.feature_name2label

Dictionary that maps feature names to feature labels

dclab.dfn.scalar_feature_names

List of valid scalar feature names

parse functions

dclab.parse_funcs.fbool(value)[source]

boolean

dclab.parse_funcs.fint(value)[source]

integer

dclab.parse_funcs.fintlist(alist)[source]

A list of integers

dclab.parse_funcs.lcstr(astr)[source]

lower-case string

dclab.parse_funcs.func_types = {<function fbool>: <class 'bool'>, <function fint>: <class 'int'>, <function fintlist>: <class 'list'>, <function lcstr>: <class 'str'>}

maps functions to their expected output types

RT-DC dataset manipulation

Base class

class dclab.rtdc_dataset.RTDCBase(identifier=None)[source]

RT-DC measurement base class

Notes

Besides the filter arrays for each data feature, there is a manual boolean filter array RTDCBase.filter.manual that can be edited by the user - a boolean value of False means that the event is excluded from all computations.

apply_filter(force=[])[source]

Compute the filters for the dataset

get_downsampled_scatter(xax='area_um', yax='deform', downsample=0, xscale='linear', yscale='linear', remove_invalid=False, ret_mask=False)[source]

Downsampling by removing points at dense locations

Parameters:
  • xax (str) – Identifier for x axis (e.g. “area_um”, “aspect”, “deform”)
  • yax (str) – Identifier for y axis
  • downsample (int) –

    Number of points to draw in the down-sampled plot. This number is either

    • >=1: exactly downsample to this number by randomly adding
      or removing points
    • 0 : do not perform downsampling
  • xscale (str) – If set to “log”, take the logarithm of the x-values before performing downsampling. This is useful when data are are displayed on a log-scale. Defaults to “linear”.
  • yscale (str) – See xscale.
  • remove_invalid (bool) – Remove nan and inf values before downsampling; if set to True, the actual number of samples returned might be smaller than downsample due to infinite or nan values (e.g. due to logarithmic scaling).
  • ret_mask (bool) – If set to True, returns a boolean array of length len(self) where True values identify the filtered data.
Returns:

  • xnew, xnew (1d ndarray of lenght N) – Filtered data; N is either identical to downsample or smaller (if remove_invalid==True)
  • mask (1d boolean array of length len(RTDCBase)) – Array for identifying the downsampled data points

get_kde_contour(xax='area_um', yax='deform', xacc=None, yacc=None, kde_type='histogram', kde_kwargs={}, xscale='linear', yscale='linear')[source]

Evaluate the kernel density estimate for contour plots

Parameters:
  • xax (str) – Identifier for X axis (e.g. “area_um”, “aspect”, “deform”)
  • yax (str) – Identifier for Y axis
  • xacc (float) – Contour accuracy in x direction
  • yacc (float) – Contour accuracy in y direction
  • kde_type (str) – The KDE method to use
  • kde_kwargs (dict) – Additional keyword arguments to the KDE method
  • xscale (str) – If set to “log”, take the logarithm of the x-values before computing the KDE. This is useful when data are are displayed on a log-scale. Defaults to “linear”.
  • yscale (str) – See xscale.
Returns:

X, Y, Z – The kernel density Z evaluated on a rectangular grid (X,Y).

Return type:

coordinates

get_kde_scatter(xax='area_um', yax='deform', positions=None, kde_type='histogram', kde_kwargs={}, xscale='linear', yscale='linear')[source]

Evaluate the kernel density estimate for scatter plots

Parameters:
  • xax (str) – Identifier for X axis (e.g. “area_um”, “aspect”, “deform”)
  • yax (str) – Identifier for Y axis
  • positions (list of two 1d ndarrays or ndarray of shape (2, N)) – The positions where the KDE will be computed. Note that the KDE estimate is computed from the the points that are set in self.filter.all.
  • kde_type (str) – The KDE method to use
  • kde_kwargs (dict) – Additional keyword arguments to the KDE method
  • xscale (str) – If set to “log”, take the logarithm of the x-values before computing the KDE. This is useful when data are are displayed on a log-scale. Defaults to “linear”.
  • yscale (str) – See xscale.
Returns:

density – The kernel density evaluated for the filtered data points.

Return type:

1d ndarray

static get_kde_spacing(a, scale='linear', method=<function bin_width_doane>, method_kw={}, feat='undefined', ret_scaled=False)[source]

Convenience function for computing the contour spacing

Parameters:
  • a (ndarray) – feature data
  • scale (str) – how the data should be scaled (“log” or “linear”)
  • method (callable) – KDE method to use (see kde_methods submodule)
  • method_kw (dict) – keyword arguments to method
  • feat (str) – feature name for debugging
  • ret_scale (bol) – whether or not to return the scaled array of a
polygon_filter_add(filt)[source]

Associate a Polygon Filter with this instance

Parameters:filt (int or instance of PolygonFilter) – The polygon filter to add
polygon_filter_rm(filt)[source]

Remove a polygon filter from this instance

Parameters:filt (int or instance of PolygonFilter) – The polygon filter to remove
reset_filter()[source]

Reset the current filter

config = None

Configuration of the measurement

export = None

Export functionalities; instance of dclab.rtdc_dataset.export.Export.

features

All available features

features_innate

All features excluding ancillary features

features_loaded

All features including ancillary that have been computed

Notes

Features that are computationally cheap to compute are always included. They are defined in dclab.rtdc_dataset.ancillaries.FEATURES_RAPID.

filter = None

Filtering functionalities; instance of dclab.rtdc_dataset.filter.Filter.

format = None

Dataset format (derived from class name)

hash

Reproducible dataset hash (defined by derived classes)

identifier

Unique (unreproducible) identifier

logs = None

Dictionary of log files. Each log file is a list of strings (one string per line).

title = None

Title of the measurement

Dictionary format

class dclab.rtdc_dataset.RTDC_Dict(ddict, *args, **kwargs)[source]

Dictionary-based RT-DC dataset

Parameters:
  • ddict (dict) – Dictionary with keys from dclab.definitions.feature_names (e.g. “area_cvx”, “deform”, “image”) with which the class will be instantiated. The configuration is set to the default configuration of dclab.
  • *args – Arguments for RTDCBase
  • **kwargs – Keyword arguments for RTDCBase

HDF5 (.rtdc) format

class dclab.rtdc_dataset.RTDC_HDF5(h5path, *args, **kwargs)[source]

HDF5 file format for RT-DC measurements

Parameters:
  • h5path (str or pathlib.Path) – Path to a ‘.tdms’ measurement file.
  • *args – Arguments for RTDCBase
  • **kwargs – Keyword arguments for RTDCBase
path

Path to the experimental HDF5 (.rtdc) file

Type:pathlib.Path
static parse_config(h5path)[source]

Parse the RT-DC configuration of an hdf5 file

dclab.rtdc_dataset.fmt_hdf5.MIN_DCLAB_EXPORT_VERSION = '0.3.3.dev2'

rtdc files exported with dclab prior to this version are not supported

Hierarchy format

class dclab.rtdc_dataset.RTDC_Hierarchy(hparent, apply_filter=True, *args, **kwargs)[source]

Hierarchy dataset (filtered from RTDCBase)

A few words on hierarchies: The idea is that a subclass of RTDCBase can use the filtered data of another subclass of RTDCBase and interpret these data as unfiltered events. This comes in handy e.g. when the percentage of different subpopulations need to be distinguished without the noise in the original data.

Children in hierarchies always update their data according to the filtered event data from their parent when apply_filter is called. This makes it easier to save and load hierarchy children with e.g. Shape-Out and it makes the handling of hierarchies more intuitive (when the parent changes, the child changes as well).

Parameters:
  • hparent (instance of RTDCBase) – The hierarchy parent
  • apply_filter (bool) – Whether to apply the filter durint instantiation; If set to False, apply_filter must be called manually.
  • *args – Arguments for RTDCBase
  • **kwargs – Keyword arguments for RTDCBase
hparent

Hierarchy parent of this instance

Type:RTDCBase

TDMS format

class dclab.rtdc_dataset.RTDC_TDMS(tdms_path, *args, **kwargs)[source]

TDMS file format for RT-DC measurements

Parameters:
  • tdms_path (str or pathlib.Path) – Path to a ‘.tdms’ measurement file.
  • *args – Arguments for RTDCBase
  • **kwargs – Keyword arguments for RTDCBase
path

Path to the experimental dataset (main .tdms file)

Type:pathlib.Path
dclab.rtdc_dataset.fmt_tdms.get_project_name_from_path(path, append_mx=False)[source]

Get the project name from a path.

For a path “/home/peter/hans/HLC12398/online/M1_13.tdms” or For a path “/home/peter/hans/HLC12398/online/data/M1_13.tdms” or without the “.tdms” file, this will return always “HLC12398”.

Parameters:
  • path (str) – path to tdms file
  • append_mx (bool) – append measurement number, e.g. “M1”
dclab.rtdc_dataset.fmt_tdms.get_tdms_files(directory)[source]

Recursively find projects based on ‘.tdms’ file endings

Searches the directory recursively and return a sorted list of all found ‘.tdms’ project files, except fluorescence data trace files which end with _traces.tdms.

config

class dclab.rtdc_dataset.config.Configuration(files=[], cfg={})[source]

Configuration class for RT-DC datasets

This class has a dictionary-like interface to access and set configuration values, e.g.

cfg = load_from_file("/path/to/config.txt")
# access the channel width
cfg["setup"]["channel width"]
# modify the channel width
cfg["setup"]["channel width"] = 30
Parameters:
  • files (list of files) – The config files with which to initialize the configuration
  • cfg (dict-like) – The dictionary with which to initialize the configuration
copy()[source]

Return copy of current configuration

keys()[source]

Return the configuration keys (sections)

save(filename)[source]

Save the configuration to a file

tostring(sections=None)[source]

Convert the configuration to its string representation

The optional argument sections allows to export only specific sections of the configuration, i.e. sections=dclab.dfn.CFG_METADATA will only export configuration data from the original measurement and no filtering data.

update(newcfg)[source]

Update current config with a dictionary

dclab.rtdc_dataset.config.load_from_file(cfg_file)[source]

Load the configuration from a file

Parameters:cfg_file (str) – Path to configuration file
Returns:cfg – Dictionary with configuration parameters
Return type:CaseInsensitiveDict

export

exception dclab.rtdc_dataset.export.NoImageWarning[source]
class dclab.rtdc_dataset.export.Export(rtdc_ds)[source]

Export functionalities for RT-DC datasets

avi(path, filtered=True, override=False)[source]

Exports filtered event images to an avi file

Parameters:
  • path (str) – Path to a .avi file. The ending .avi is added automatically.
  • filtered (bool) – If set to True, only the filtered data (index in ds.filter.all) are used.
  • override (bool) – If set to True, an existing file path will be overridden. If set to False, raises OSError if path exists.

Notes

Raises OSError if current dataset does not contain image data

fcs(path, features, meta_data={}, filtered=True, override=False)[source]

Export the data of an RT-DC dataset to an .fcs file

Parameters:
  • mm (instance of dclab.RTDCBase) – The dataset that will be exported.
  • path (str) – Path to an .fcs file. The ending .fcs is added automatically.
  • features (list of str) – The features in the resulting .fcs file. These are strings that are defined in dclab.definitions.scalar_feature_names, e.g. “area_cvx”, “deform”, “frame”, “fl1_max”, “aspect”.
  • meta_data (dict) – User-defined, optional key-value pairs that are stored in the primary TEXT segment of the FCS file; the version of dclab is stored there by default
  • filtered (bool) – If set to True, only the filtered data (index in ds.filter.all) are used.
  • override (bool) – If set to True, an existing file path will be overridden. If set to False, raises OSError if path exists.

Notes

Due to incompatibility with the .fcs file format, all events with NaN-valued features are not exported.

hdf5(path, features, filtered=True, override=False, compression='gzip')[source]

Export the data of the current instance to an HDF5 file

Parameters:
  • path (str) – Path to an .rtdc file. The ending .rtdc is added automatically.
  • features (list of str) – The features in the resulting .rtdc file. These are strings that are defined in dclab.definitions.feature_names, e.g. “area_cvx”, “deform”, “frame”, “fl1_max”, “image”.
  • filtered (bool) – If set to True, only the filtered data (index in ds.filter.all) are used.
  • override (bool) – If set to True, an existing file path will be overridden. If set to False, raises OSError if path exists.
  • compression (str or None) – Compression method for “contour”, “image”, and “trace” data as well as logs; one of [None, “lzf”, “gzip”, “szip”].
tsv(path, features, meta_data={'dclab version': '0.21.2'}, filtered=True, override=False)[source]

Export the data of the current instance to a .tsv file

Parameters:
  • path (str) – Path to a .tsv file. The ending .tsv is added automatically.
  • features (list of str) – The features in the resulting .tsv file. These are strings that are defined in dclab.definitions.scalar_feature_names, e.g. “area_cvx”, “deform”, “frame”, “fl1_max”, “aspect”.
  • meta_data (dict) – User-defined, optional key-value pairs that are stored at the beginning of the tsv file - one key-value pair is stored per line which starts with a hash. The version of dclab is stored there by default.
  • filtered (bool) – If set to True, only the filtered data (index in ds.filter.all) are used.
  • override (bool) – If set to True, an existing file path will be overridden. If set to False, raises OSError if path exists.

filter

class dclab.rtdc_dataset.filter.Filter(rtdc_ds)[source]

Boolean filter arrays for RT-DC measurements

Parameters:rtdc_ds (instance of RTDCBase) – The RT-DC dataset the filter applies to
reset()[source]

Reset all filters

update(rtdc_ds, force=[])[source]

Update the filters according to rtdc_ds.config[“filtering”]

Parameters:
  • rtdc_ds (dclab.rtdc_dataset.core.RTDCBase) – The measurement to which the filter is applied
  • force (list) – A list of feature names that must be refiltered with min/max values.

Notes

This function is called when ds.apply_filter is called.

low-level functionalities

downsampling

Content-based downsampling of ndarrays

dclab.downsampling.downsample_rand(a, samples, remove_invalid=False, ret_idx=False)[source]

Downsampling by randomly removing points

Parameters:
  • a (1d ndarray) – The input array to downsample
  • samples (int) – The desired number of samples
  • remove_invalid (bool) – Remove nan and inf values before downsampling
  • ret_idx (bool) – Also return a boolean array that corresponds to the downsampled indices in a.
Returns:

  • dsa (1d ndarray of size samples) – The pseudo-randomly downsampled array a
  • idx (1d boolean array with same shape as a) – Only returned if ret_idx is True. A boolean array such that a[idx] == dsa

dclab.downsampling.norm(a, ref1, ref2)[source]

Normalize a with min/max values of ref1, using all elements of ref1 where the ref1 and ref2 are not nan or inf

dclab.downsampling.valid(a, b)[source]

Check whether a and b are not inf or nan

features

dclab.features.contour.get_contour(mask)[source]

Compute the image contour from a mask

The contour is computed in a very inefficient way using scikit-image and a conversion of float coordinates to pixel coordinates.

Parameters:mask (binary ndarray of shape (M,N) or (K,M,N)) – The mask outlining the pixel positions of the event. If a 3d array is given, then K indexes the individual contours.
Returns:cont – A 2D array that holds the contour of an event (in pixels) e.g. obtained using mm.contour where mm is an instance of RTDCBase. The first and second columns of cont correspond to the x- and y-coordinates of the contour.
Return type:ndarray or list of K ndarrays of shape (J,2)
dclab.features.bright.get_bright(mask, image, ret_data='avg, sd')[source]

Compute avg and/or std of the event brightness

The event brightness is defined by the gray-scale values of the image data within the event mask area.

Parameters:
  • mask (ndarray or list of ndarrays of shape (M,N) and dtype bool) – The mask values, True where the event is located in image.
  • image (ndarray or list of ndarrays of shape (M,N)) – A 2D array that holds the image in form of grayscale values of an event.
  • ret_data (str) – A comma-separated list of metrices to compute - “avg”: compute the average - “sd”: compute the standard deviation Selected metrics are returned in alphabetical order.
Returns:

  • bright_avg (float or ndarray of size N) – Average image data within the contour
  • bright_std (float or ndarray of size N) – Standard deviation of image data within the contour

dclab.features.emodulus.get_emodulus(area_um, deform, medium='CellCarrier', channel_width=20.0, flow_rate=0.16, px_um=0.34, temperature=23.0, copy=True)[source]

Compute apparent Young’s modulus using a look-up table

Parameters:
  • area_um (float or ndarray) – Apparent (2D image) area [µm²] of the event(s)
  • deform (float or ndarray) – The deformation (1-circularity) of the event(s)
  • medium (str or float) – The medium to compute the viscosity for. If a string is given, the viscosity is computed. If a float is given, this value is used as the viscosity in mPa*s (Note that temperature must be set to None in this case).
  • channel_width (float) – The channel width [µm]
  • flow_rate (float) – Flow rate [µl/s]
  • px_um (float) – The detector pixel size [µm] used for pixelation correction. Set to zero to disable.
  • temperature (float, ndarray, or None) – Temperature [°C] of the event(s)
  • copy (bool) – Copy input arrays. If set to false, input arrays are overridden.
Returns:

elasticity – Apparent Young’s modulus in kPa

Return type:

float or ndarray

Notes

  • The look-up table used was computed with finite elements methods according to [Mokbel2017].
  • The computation of the Young’s modulus takes into account corrections for the viscosity (medium, channel width, flow rate, and temperature) [Mietke2015] and corrections for pixelation of the area and the deformation which are computed from a (pixelated) image [Herold2017].

See also

dclab.features.emodulus_viscosity.get_viscosity()
compute viscosity for known media
dclab.features.emodulus_viscosity.get_viscosity(medium='CellCarrier', channel_width=20.0, flow_rate=0.16, temperature=23.0)[source]

Returns the viscosity for RT-DC-specific media

Parameters:
  • medium (str) – The medium to compute the viscosity for; Valid values are defined in KNOWN_MEDIA.
  • channel_width (float) – The channel width in µm
  • flow_rate (float) – Flow rate in µl/s
  • temperature (float or ndarray) – Temperature in °C
Returns:

viscosity – Viscosity in mPa*s

Return type:

float or ndarray

Notes

  • CellCarrier and CellCarrier B media are optimized for RT-DC measurements.
  • Values for the viscosity of water are computed using equation (15) from [Kestin_1978].
  • A TemperatureOutOfRangeWarning is issued if the input temperature range exceeds the temperature ranges given by [Herold2017] and [Kestin_1978].
dclab.features.emodulus_viscosity.KNOWN_MEDIA = ['CellCarrier', 'CellCarrierB', 'water']

Media for which computation of viscosity is defined

dclab.features.fl_crosstalk.correct_crosstalk(fl1, fl2, fl3, fl_channel, ct21=0, ct31=0, ct12=0, ct32=0, ct13=0, ct23=0)[source]

Perform crosstalk correction

Parameters:
  • fli (int, float, or np.ndarray) – Measured fluorescence signals
  • fl_channel (int (1, 2, or 3)) – The channel number for which the crosstalk-corrected signal should be computed
  • cij (float) – Spill (crosstalk or bleed-through) from channel i to channel j This spill is computed from the fluorescence signal of e.g. single-stained positive control cells; It is defined by the ratio of the fluorescence signals of the two channels, i.e cij = flj / fli.

See also

get_compensation_matrix()
compute the inverse crosstalk matrix

Notes

If there are only two channels (e.g. fl1 and fl2), then the crosstalk to and from the other channel (ct31, ct32, ct13, ct23) should be set to zero.

dclab.features.fl_crosstalk.get_compensation_matrix(ct21, ct31, ct12, ct32, ct13, ct23)[source]

Compute crosstalk inversion matrix

The spillover matrix is

| c11 c12 c13 |
| c21 c22 c23 |
| c31 c32 c33 |

The diagonal elements are set to 1, i.e.

ct11 = c22 = c33 = 1

Parameters:cij (float) – Spill from channel i to channel j
Returns:inv – Compensation matrix (inverted spillover matrix)
Return type:np.ndarray
dclab.features.inert_ratio.get_inert_ratio_cvx(cont)[source]

Compute the inertia ratio of the convex hull of a contour

The inertia ratio is computed from the central second order of moments along x (mu20) and y (mu02) via sqrt(mu20/mu02).

Parameters:cont (ndarray or list of ndarrays of shape (N,2)) – A 2D array that holds the contour of an event (in pixels) e.g. obtained using mm.contour where mm is an instance of RTDCBase. The first and second columns of cont correspond to the x- and y-coordinates of the contour.
Returns:inert_ratio_cvx – The inertia ratio of the contour’s convex hull
Return type:float or ndarray of size N

Notes

The contour moments mu20 and mu02 are computed the same way they are computed in OpenCV’s moments.cpp.

See also

get_inert_ratio_raw()
Compute inertia ratio of a raw contour

References

dclab.features.inert_ratio.get_inert_ratio_raw(cont)[source]

Compute the inertia ratio of a contour

The inertia ratio is computed from the central second order of moments along x (mu20) and y (mu02) via sqrt(mu20/mu02).

Parameters:cont (ndarray or list of ndarrays of shape (N,2)) – A 2D array that holds the contour of an event (in pixels) e.g. obtained using mm.contour where mm is an instance of RTDCBase. The first and second columns of cont correspond to the x- and y-coordinates of the contour.
Returns:inert_ratio_raw – The inertia ratio of the contour
Return type:float or ndarray of size N

Notes

The contour moments mu20 and mu02 are computed the same way they are computed in OpenCV’s moments.cpp.

See also

get_inert_ratio_cvx()
Compute inertia ratio of the convex hull of a contour

References

dclab.features.volume.get_volume(cont, pos_x, pos_y, pix)[source]

Calculate the volume of a polygon revolved around an axis

The volume estimation assumes rotational symmetry. Green`s theorem and the Gaussian divergence theorem allow to formulate the volume as a line integral.

Parameters:
  • cont (ndarray or list of ndarrays of shape (N,2)) – A 2D array that holds the contour of an event [px] e.g. obtained using mm.contour where mm is an instance of RTDCBase. The first and second columns of cont correspond to the x- and y-coordinates of the contour.
  • pos_x (float or ndarray of length N) – The x coordinate(s) of the centroid of the event(s) [µm] e.g. obtained using mm.pos_x
  • pos_y (float or ndarray of length N) – The y coordinate(s) of the centroid of the event(s) [µm] e.g. obtained using mm.pos_y
  • px_um (float) – The detector pixel size in µm. e.g. obtained using: mm.config[“image”][“pix size”]
Returns:

volume – volume in um^3

Return type:

float or ndarray

Notes

The computation of the volume is based on a full rotation of the upper and the lower halves of the contour from which the average is then used.

The volume is computed radially from the the center position given by (pos_x, pos_y). For sufficiently smooth contours, such as densely sampled ellipses, the center position does not play an important role. For contours that are given on a coarse grid, as is the case for RT-DC, the center position must be given.

References

isoelastics

Isoelastics management

exception dclab.isoelastics.IsoelasticsEmodulusMeaninglessWarning[source]
class dclab.isoelastics.Isoelastics(paths=[])[source]
add(isoel, col1, col2, channel_width, flow_rate, viscosity, method)[source]

Add isoelastics

Parameters:
  • isoel (list of ndarrays) – Each list item resembles one isoelastic line stored as an array of shape (N,3). The last column contains the emodulus data.
  • col1 (str) – Name of the first feature of all isoelastics (e.g. isoel[0][:,0])
  • col2 (str) – Name of the second feature of all isoelastics (e.g. isoel[0][:,1])
  • channel_width (float) – Channel width in µm
  • flow_rate (float) – Flow rate through the channel in µl/s
  • viscosity (float) – Viscosity of the medium in mPa*s
  • method (str) – The method used to compute the isoelastics (must be one of VALID_METHODS).

Notes

The following isoelastics are automatically added for user convenience: - isoelastics with col1 and col2 interchanged - isoelastics for circularity if deformation was given

static add_px_err(isoel, col1, col2, px_um, inplace=False)[source]

Undo pixelation correction

Isoelasticity lines are already corrected for pixelation effects as described in

Mapping of Deformation to Apparent Young’s Modulus in Real-Time Deformability Cytometry Christoph Herold, arXiv:1704.00572 [cond-mat.soft] (2017) https://arxiv.org/abs/1704.00572.

If the isoealsticity lines are displayed with deformation data that are not corrected, then the lines must be “un”-corrected, i.e. the pixelation error must be added to the lines to match the experimental data.

Parameters:
  • isoel (list of 2d ndarrays of shape (N, 3)) – Each item in the list corresponds to one isoelasticity line. The first column is defined by col1, the second by col2, and the third column is the emodulus.
  • col2 (col1,) – Define the fist to columns of each isoelasticity line. One of [“area_um”, “circ”, “deform”]
  • px_um (float) – Pixel size [µm]
static check_col12(col1, col2)[source]
static convert(isoel, col1, col2, channel_width_in, channel_width_out, flow_rate_in, flow_rate_out, viscosity_in, viscosity_out, inplace=False)[source]

Convert isoelastics in area_um-deform space

Parameters:
  • isoel (list of 2d ndarrays of shape (N, 3)) – Each item in the list corresponds to one isoelasticity line. The first column is defined by col1, the second by col2, and the third column is the emodulus.
  • col2 (col1,) – Define the fist to columns of each isoelasticity line. One of [“area_um”, “circ”, “deform”]
  • channel_width_in (float) – Original channel width [µm]
  • channel_width_out (float) – Target channel width [µm]
  • flow_rate_in (float) – Original flow rate [µl/s]
  • flow_rate_in – Target flow rate [µl/s]
  • viscosity_in (float) – Original viscosity [mPa*s]
  • viscosity_out (float) – Target viscosity [mPa*s]

Notes

If only the positions of the isoelastics are of interest and not the value of the elastic modulus, then it is sufficient to supply values for the channel width and set the values for flow rate and viscosity to a constant (e.g. 1).

See also

dclab.features.emodulus.convert()
conversion method used
get(col1, col2, method, channel_width, flow_rate=None, viscosity=None, add_px_err=False, px_um=None)[source]

Get isoelastics

Parameters:
  • col1 (str) – Name of the first feature of all isoelastics (e.g. isoel[0][:,0])
  • col2 (str) – Name of the second feature of all isoelastics (e.g. isoel[0][:,1])
  • method (str) – The method used to compute the isoelastics (must be one of VALID_METHODS).
  • channel_width (float) – Channel width in µm
  • flow_rate (float or None) – Flow rate through the channel in µl/s. If set to None, the flow rate of the imported data will be used (only do this if you do not need the correct values for elastic moduli).
  • viscosity (float or None) – Viscosity of the medium in mPa*s. If set to None, the flow rate of the imported data will be used (only do this if you do not need the correct values for elastic moduli).
  • add_px_err (bool) – If True, add pixelation errors according to C. Herold (2017), https://arxiv.org/abs/1704.00572
  • px_um (float) – Pixel size [µm], used for pixelation error computation

See also

dclab.features.emodulus.convert()
conversion in-between channel sizes and viscosities
dclab.features.emodulus.corrpix_deform_delta()
pixelation error that is applied to the deformation data
get_with_rtdcbase(col1, col2, method, dataset, viscosity=None, add_px_err=False)[source]

Convenience method that extracts the metadata from RTDCBase

Parameters:
  • col1 (str) – Name of the first feature of all isoelastics (e.g. isoel[0][:,0])
  • col2 (str) – Name of the second feature of all isoelastics (e.g. isoel[0][:,1])
  • method (str) – The method used to compute the isoelastics (must be one of VALID_METHODS).
  • dataset (dclab.rtdc_dataset.RTDCBase) – The dataset from which to obtain the metadata.
  • viscosity (float, None, or False) – Viscosity of the medium in mPa*s. If set to None, the viscosity is computed from the meta data (medium, flow rate, channel width, temperature) in the [setup] config section. If this is not possible, the flow rate of the imported data is used and a warning will be issued.
  • add_px_err (bool) – If True, add pixelation errors according to C. Herold (2017), https://arxiv.org/abs/1704.00572
load_data(path)[source]

Load isoelastics from a text file

The text file is loaded with numpy.loadtxt and must have three columns, representing the two data columns and the elastic modulus with units defined in definitions.py. The file header must have a section defining meta data of the content like so:

# […] # # - column 1: area_um # - column 2: deform # - column 3: emodulus # - channel width [um]: 20 # - flow rate [ul/s]: 0.04 # - viscosity [mPa*s]: 15 # - method: analytical # # […]
Parameters:path (str) – Path to a isoelastics text file
class dclab.isoelastics.IsoelasticsDict[source]
dclab.isoelastics.get_default()[source]

Return default isoelasticity lines

kde_contours

dclab.kde_contours.find_contours_level(density, x, y, level, closed=False)[source]

Find iso-valued density contours for a given level value

Parameters:
  • density (2d ndarray of shape (M, N)) – Kernel density estimate for which to compute the contours
  • x (2d ndarray of shape (M, N) or 1d ndarray of size M) – X-values corresponding to density
  • y (2d ndarray of shape (M, N) or 1d ndarray of size M) – Y-values corresponding to density
  • level (float between 0 and 1) – Value along which to find contours in density relative to its maximum
Returns:

contours – Contours found for the given level value

Return type:

list of ndarrays of shape (P, 2)

See also

skimage.measure.find_contours()
Contour finding algorithm used
dclab.kde_contours.get_quantile_levels(density, x, y, xp, yp, q, normalize=True)[source]

Compute density levels for given quantiles by interpolation

For a given 2D density, compute the density levels at which the resulting contours contain the fraction 1-q of all data points. E.g. for a measurement of 1000 events, all contours at the level corresponding to a quantile of q=0.95 (95th percentile) contain 50 events (5%).

Parameters:
  • density (2d ndarray of shape (M, N)) – Kernel density estimate for which to compute the contours
  • x (2d ndarray of shape (M, N) or 1d ndarray of size M) – X-values corresponding to density
  • y (2d ndarray of shape (M, N) or 1d ndarray of size M) – Y-values corresponding to density
  • xp (1d ndarray of size D) – Event x-data from which to compute the quantile
  • yp (1d ndarray of size D) – Event y-data from which to compute the quantile
  • q (array_like or float between 0 and 1) – Quantile along which to find contours in density relative to its maximum
  • normalize (bool) – Whether output levels should be normalized to the maximum of density
Returns:

level – Contours level(s) corresponding to the given quantile

Return type:

np.ndarray or float

Notes

NaN-values events in xp and yp are ignored.

kde_methods

Kernel Density Estimation methods

dclab.kde_methods.bin_num_doane(a)[source]

Compute number of bins based on Doane’s formula

dclab.kde_methods.bin_width_doane(a)[source]

Compute contour spacing based on Doane’s formula

References

Notes

Doane’s formula is actually designed for histograms. This function is kept here for backwards-compatibility reasons. It is highly recommended to use bin_width_percentile() instead.

dclab.kde_methods.bin_width_percentile(a)[source]

Compute contour spacing based on data percentiles

The 10th and the 90th percentile of the input data are taken. The spacing then computes to the difference between those two percentiles divided by 23.

Notes

The Freedman–Diaconis rule uses the interquartile range and normalizes to the third root of len(a). Such things do not work very well for RT-DC data, because len(a) is huge. Here we use just the top an bottom 10th percentiles with a fixed normalization.

dclab.kde_methods.get_bad_vals(x, y)[source]
dclab.kde_methods.ignore_nan_inf(kde_method)[source]

Ignores nans and infs from the input data

Invalid positions in the resulting density are set to nan.

dclab.kde_methods.kde_gauss(events_x, events_y, xout=None, yout=None, *args, **kwargs)[source]

Gaussian Kernel Density Estimation

Parameters:
  • events_y (events_x,) – The input points for kernel density estimation. Input is flattened automatically.
  • yout (xout,) – The coordinates at which the KDE should be computed. If set to none, input coordinates are used.
Returns:

density – The KDE for the points in (xout, yout)

Return type:

ndarray, same shape as xout

See also

scipy.stats.gaussian_kde

Notes

This is a wrapped version that ignores nan and inf values.

dclab.kde_methods.kde_histogram(events_x, events_y, xout=None, yout=None, *args, **kwargs)[source]

Histogram-based Kernel Density Estimation

Parameters:
  • events_y (events_x,) – The input points for kernel density estimation. Input is flattened automatically.
  • yout (xout,) – The coordinates at which the KDE should be computed. If set to none, input coordinates are used.
  • bins (tuple (binsx, binsy)) – The number of bins to use for the histogram.
Returns:

density – The KDE for the points in (xout, yout)

Return type:

ndarray, same shape as xout

See also

numpy.histogram2d scipy.interpolate.RectBivariateSpline

Notes

This is a wrapped version that ignores nan and inf values.

dclab.kde_methods.kde_multivariate(events_x, events_y, xout=None, yout=None, *args, **kwargs)[source]

Multivariate Kernel Density Estimation

Parameters:
  • events_y (events_x,) – The input points for kernel density estimation. Input is flattened automatically.
  • bw (tuple (bwx, bwy) or None) – The bandwith for kernel density estimation.
  • yout (xout,) – The coordinates at which the KDE should be computed. If set to none, input coordinates are used.
Returns:

density – The KDE for the points in (xout, yout)

Return type:

ndarray, same shape as xout

See also

statsmodels.nonparametric.kernel_density.KDEMultivariate

Notes

This is a wrapped version that ignores nan and inf values.

dclab.kde_methods.kde_none(events_x, events_y, xout=None, yout=None)[source]

No Kernel Density Estimation

Parameters:
  • events_y (events_x,) – The input points for kernel density estimation. Input is flattened automatically.
  • yout (xout,) – The coordinates at which the KDE should be computed. If set to none, input coordinates are used.
Returns:

density – The KDE for the points in (xout, yout)

Return type:

ndarray, same shape as xout

Notes

This method is a convenience method that always returns ones in the shape that the other methods in this module produce.

polygon_filter

exception dclab.polygon_filter.FilterIdExistsWarning[source]
exception dclab.polygon_filter.PolygonFilterError[source]
class dclab.polygon_filter.PolygonFilter(axes=None, points=None, inverted=False, name=None, filename=None, fileid=0, unique_id=None)[source]

An object for filtering RTDC data based on a polygonial area

Parameters:
  • axes (tuple of str) – The axes/features on which the polygon is defined. The first axis is the x-axis. Example: (“area_um”, “deform”).
  • points (array-like object of shape (N,2)) – The N coordinates (x,y) of the polygon. The exact order is important.
  • inverted (bool) – Invert the polygon filter. This parameter is overridden if filename is given.
  • name (str) – A name for the polygon (optional).
  • filename (str) – A path to a .poly file as create by this classes’ save method. If filename is given, all other parameters are ignored.
  • fileid (int) – Which filter to import from the file (starting at 0).
  • unique_id (int) – An integer defining the unique id of the new instance.

Notes

The minimal arguments to this class are either filename OR (axes, points). If filename is set, all parameters are taken from the given .poly file.

static clear_all_filters()[source]

Remove all filters and reset instance counter

copy(invert=False)[source]

Return a copy of the current instance

Parameters:invert (bool) – The copy will be inverted w.r.t. the original
filter(datax, datay)[source]

Filter a set of datax and datay according to self.points

static get_instance_from_id(unique_id)[source]

Get an instance of the PolygonFilter using a unique id

static import_all(path)[source]

Import all polygons from a .poly file.

Returns a list of the imported polygon filters

static instace_exists(unique_id)[source]

Determine whether an instance with this unique id exists

static point_in_poly(p, poly)[source]

Determine whether a point is within a polygon area

Uses the ray casting algorithm.

Parameters:
  • p (float) – Coordinates of the point
  • poly (array_like of shape (N, 2)) – Polygon (PolygonFilter.points)
Returns:

insideTrue, if point is inside.

Return type:

bool

Notes

If p lies on a side of the polygon, it is defined as

  • “inside” if it is on the top or right
  • “outside” if it is on the lower or left
static remove(unique_id)[source]

Remove a polygon filter from PolygonFilter.instances

save(polyfile, ret_fobj=False)[source]

Save all data to a text file (appends data if file exists).

Polyfile can be either a path to a file or a file object that was opened with the write “w” parameter. By using the file object, multiple instances of this class can write their data.

If ret_fobj is True, then the file object will not be closed and returned.

static save_all(polyfile)[source]

Save all polygon filters

static unique_id_exists(pid)[source]

Whether or not a filter with this unique id exists

hash

Hash of axes, points, and inverted

instances = [<dclab.polygon_filter.PolygonFilter object>]
points
dclab.polygon_filter.get_polygon_filter_names()[source]

Get the names of all polygon filters in the order of creation

statistics

Statistics computation for RT-DC dataset instances

exception dclab.statistics.BadMethodWarning[source]
class dclab.statistics.Statistics(name, method, req_feature=False)[source]

A helper class for computing statistics

All statistical methods are registered in the dictionary Statistics.available_methods.

get_feature(ds, feat)[source]

Return filtered feature data

The features are filtered according to the user-defined filters, using the information in ds.filter.all. In addition, all nan and inf values are purged.

Parameters:
available_methods = {'%-gated': <dclab.statistics.Statistics object>, 'Events': <dclab.statistics.Statistics object>, 'Flow rate': <dclab.statistics.Statistics object>, 'Mean': <dclab.statistics.Statistics object>, 'Median': <dclab.statistics.Statistics object>, 'Mode': <dclab.statistics.Statistics object>, 'SD': <dclab.statistics.Statistics object>}
dclab.statistics.flow_rate(ds)[source]

Return the flow rate of an RT-DC dataset

dclab.statistics.get_statistics(ds, methods=None, features=None)[source]

Compute statistics for an RT-DC dataset

Parameters:
  • ds (dclab.rtdc_dataset.RTDCBase) – The dataset for which to compute the statistics.
  • methods (list of str or None) – The methods wih which to compute the statistics. The list of available methods is given with dclab.statistics.Statistics.available_methods.keys() If set to None, statistics for all methods are computed.
  • features (list of str) – Feature name identifiers are defined in dclab.definitions.scalar_feature_names. If set to None, statistics for all axes are computed.
Returns:

  • header (list of str) – The header (feature + method names) of the computed statistics.
  • values (list of float) – The computed statistics.

dclab.statistics.mode(data)[source]

Compute an intelligent value for the mode

The most common value in experimental is not very useful if there are a lot of digits after the comma. This method approaches this issue by rounding to bin size that is determined by the Freedman–Diaconis rule.

Parameters:data (1d ndarray) – The data for which the mode should be computed.
Returns:mode – The mode computed with the Freedman-Diaconis rule.
Return type:float