Information theoretic estimators

JIDT Estimators (CPU)

Provide JIDT estimators.

class idtxl.estimators_jidt.JidtDiscrete(settings)[source]

Abstract class for implementation of discrete JIDT-estimators.

Abstract class for implementation of plug-in JIDT-estimators for discrete data. Child classes implement estimators for mutual information (MI), conditional mutual information (CMI), actice information storage (AIS), and transfer entropy (TE). See parent class for references.

Set common estimation parameters for discrete JIDT-estimators. For usage of these estimators see documentation for the child classes.

Args:
settingsdict [optional]

set estimator parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • discretise_method : str [optional] - if and how to discretise incoming continuous data, can be ‘max_ent’ for maximum entropy binning, ‘equal’ for equal size bins, and ‘none’ if no binning is required (default=’none’)

Note:

Discrete JIDT estimators require the data’s alphabet size for instantiation. Hence, opposed to the Kraskov and Gaussian estimators, the JAVA class is added to the object instance, while for Kraskov/ Gaussian estimators an instance of that class is added (because for the latter, objects can be instantiated independent of data properties).

estimate_surrogates_analytic(n_perm=200, **data)[source]

Return estimate of the analytical surrogate distribution.

This method must be implemented because this class’ is_analytic_null_estimator() method returns true.

Args:
n_permsint [optional]

number of permutations (default=200)

datanumpy arrays

realisations of random variables required for the calculation (varies between estimators, e.g. 2 variables for MI, 3 for CMI). Formatted as per the estimate method for this estimator.

Returns:
float | numpy array

n_perm surrogates of the average MI/CMI/TE over all samples under the null hypothesis of no relationship between var1 and var2 (in the context of conditional)

abstract get_analytic_distribution(**data)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
datanumpy arrays

realisations of random variables required for the calculation (varies between estimators, e.g. 2 variables for MI, 3 for CMI). Formatted as per the estimate method for this estimator.

Returns:
Java object

JIDT calculator that was used here

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

class idtxl.estimators_jidt.JidtDiscreteAIS(settings)[source]

Calculate AIS with JIDT’s discrete-variable implementation.

Calculate the active information storage (AIS) for one process. Call JIDT via jpype and use the discrete estimator. See parent class for references.

Results are returned in bits.

Args:
settingsdict

set estimator parameters:

  • history : int - number of samples in the target’s past used as embedding (>= 0)

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • discretise_method : str [optional] - if and how to discretise incoming continuous data, can be ‘max_ent’ for maximum entropy binning, ‘equal’ for equal size bins, and ‘none’ if no binning is required (default=’none’)

  • n_discrete_bins : int [optional] - number of discrete bins/ levels or the base of each dimension of the discrete variables (default=2). If set, this parameter overwrites/sets alph. (>= 2)

  • alph : int [optional] - number of discrete bins/levels for var1 (default=2 , or the value set for n_discrete_bins). (>= 2)

estimate(process, return_calc=False)[source]

Estimate active information storage.

Args:
processnumpy array

realisations as either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

return_calcboolean

return the calculator used here as well as the numeric calculated value(s)

Returns:
float | numpy array

average AIS over all samples or local AIS for individual samples if ‘local_values’=True

Java object

JIDT calculator that was used here. Only returned if return_calc was set.

Raises:
ex.JidtOutOfMemoryError

Raised when JIDT object cannot be instantiated due to mem error

get_analytic_distribution(process)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
processnumpy array

realisations as either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

Returns:
Java object

JIDT calculator that was used here

class idtxl.estimators_jidt.JidtDiscreteCMI(settings=None)[source]

Calculate CMI with JIDT’s implementation for discrete variables.

Calculate the conditional mutual information between two variables given the third. Call JIDT via jpype and use the discrete estimator. See parent class for references.

Results are returned in bits.

Args:
settingsdict [optional]

sets estimation parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • discretise_method : str [optional] - if and how to discretise incoming continuous data, can be ‘max_ent’ for maximum entropy binning, ‘equal’ for equal size bins, and ‘none’ if no binning is required (default=’none’)

  • n_discrete_bins : int [optional] - number of discrete bins/ levels or the base of each dimension of the discrete variables (default=2). If set, this parameter overwrites/sets alph1, alph2 and alphc

  • alph1 : int [optional] - number of discrete bins/levels for var1 (default=2, or the value set for n_discrete_bins)

  • alph2 : int [optional] - number of discrete bins/levels for var2 (default=2, or the value set for n_discrete_bins)

  • alphc : int [optional] - number of discrete bins/levels for conditional (default=2, or the value set for n_discrete_bins)

estimate(var1, var2, conditional=None, return_calc=False)[source]

Estimate conditional mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

var2numpy array

realisations of the second variable (similar to var1)

conditionalnumpy array [optional]

realisations of the conditioning variable (similar to var), if no conditional is provided, return MI between var1 and var2

return_calcboolean

return the calculator used here as well as the numeric calculated value(s)

Returns:
float | numpy array

average CMI over all samples or local CMI for individual samples if ‘local_values’=True

Java object

JIDT calculator that was used here. Only returned if return_calc was set.

Raises:
ex.JidtOutOfMemoryError

Raised when JIDT object cannot be instantiated due to mem error

get_analytic_distribution(var1, var2, conditional=None)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

var2numpy array

realisations of the second variable (similar to var1)

conditionalnumpy array [optional]

realisations of the conditioning variable (similar to var), if no conditional is provided, return MI between var1 and var2

Returns:
Java object

JIDT calculator that was used here

class idtxl.estimators_jidt.JidtDiscreteMI(settings=None)[source]

Calculate MI with JIDT’s discrete-variable implementation.

Calculate the mutual information (MI) between two variables. Call JIDT via jpype and use the discrete estimator. See parent class for references.

Results are returned in bits.

Args:
settingsdict [optional]

sets estimation parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • discretise_method : str [optional] - if and how to discretise incoming continuous data, can be ‘max_ent’ for maximum entropy binning, ‘equal’ for equal size bins, and ‘none’ if no binning is required (default=’none’)

  • n_discrete_bins : int [optional] - number of discrete bins/ levels or the base of each dimension of the discrete variables (default=2). If set, this parameter overwrites/sets alph1 and alph2

  • alph1 : int [optional] - number of discrete bins/levels for var1 (default=2, or the value set for n_discrete_bins)

  • alph2 : int [optional] - number of discrete bins/levels for var2 (default=2, or the value set for n_discrete_bins)

  • lag_mi : int [optional] - time difference in samples to calculate the lagged MI between processes (default=0)

estimate(var1, var2, return_calc=False)[source]

Estimate mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

var2numpy array

realisations of the second variable (similar to var1)

return_calcboolean

return the calculator used here as well as the numeric calculated value(s)

Returns:
float | numpy array

average MI over all samples or local MI for individual samples if ‘local_values’=True

Java object

JIDT calculator that was used here. Only returned if return_calc was set.

Raises:
ex.JidtOutOfMemoryError

Raised when JIDT object cannot be instantiated due to mem error

get_analytic_distribution(var1, var2)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

var2numpy array

realisations of the second variable (similar to var1)

Returns:
Java object

JIDT calculator that was used here

class idtxl.estimators_jidt.JidtDiscreteTE(settings)[source]

Calculate TE with JIDT’s implementation for discrete variables.

Calculate the transfer entropy between two time series processes. Call JIDT via jpype and use the discrete estimator. Transfer entropy is defined as the conditional mutual information between the source’s past state and the target’s current value, conditional on the target’s past. See parent class for references.

Results are returned in bits.

Args:
settingsdict

sets estimation parameters:

  • history_target : int - number of samples in the target’s past used as embedding. (>= 0)

  • history_source : int [optional] - number of samples in the source’s past used as embedding (default=same as the target history). (>= 1)

  • tau_source : int [optional] - source’s embedding delay (default=1). (>= 1)

  • tau_target : int [optional] - target’s embedding delay (default=1). (>= 1)

  • source_target_delay : int [optional] - information transfer delay between source and target (default=1) (>= 0)

  • discretise_method : str [optional] - if and how to discretise incoming continuous data, can be ‘max_ent’ for maximum entropy binning, ‘equal’ for equal size bins, and ‘none’ if no binning is required (default=’none’)

  • n_discrete_bins : int [optional] - number of discrete bins/ levels or the base of each dimension of the discrete variables (default=2). If set, this parameter overwrites/sets alph1 and alph2. (>= 2)

  • alph1 : int [optional] - number of discrete bins/levels for source (default=2, or the value set for n_discrete_bins). (>= 2)

  • alph2 : int [optional] - number of discrete bins/levels for target (default=2, or the value set for n_discrete_bins). (>= 2)

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

estimate(source, target, return_calc=False)[source]

Estimate transfer entropy from a source to a target variable.

Args:
sourcenumpy array

realisations of source variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

targetnumpy array

realisations of target variable (similar to var1)

return_calcboolean

return the calculator used here as well as the numeric calculated value(s)

Returns:
float | numpy array

average TE over all samples or local TE for individual samples if ‘local_values’=True

Java object

JIDT calculator that was used here. Only returned if return_calc was set.

Raises:
ex.JidtOutOfMemoryError

Raised when JIDT object cannot be instantiated due to mem error

get_analytic_distribution(source, target)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
sourcenumpy array

realisations of source variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations], array type can be float (requires discretisation) or int

targetnumpy array

realisations of target variable (similar to var1)

Returns:
Java object

JIDT calculator that was used here

class idtxl.estimators_jidt.JidtEstimator(settings=None)[source]

Abstract class for implementation of JIDT estimators.

Abstract class for implementation of JIDT estimators, child classes implement estimators for mutual information (MI), conditional mutual information (CMI), active information storage (AIS), transfer entropy (TE) using the Kraskov-Grassberger-Stoegbauer estimator for continuous data, plug-in estimators for discrete data, and Gaussian estimators for continuous Gaussian data.

References:

  • Lizier, Joseph T. (2014). JIDT: an information-theoretic toolkit for studying the dynamics of complex systems. Front Robot AI, 1(11).

  • Kraskov, A., Stoegbauer, H., & Grassberger, P. (2004). Estimating mutual information. Phys Rev E, 69(6), 066138.

  • Lizier, Joseph T., Mikhail Prokopenko, and Albert Y. Zomaya. (2012). Local measures of information storage in complex distributed computation. Inform Sci, 208, 39-54.

  • Schreiber, T. (2000). Measuring information transfer. Phys Rev Lett, 85(2), 461.

Set common estimation parameters for JIDT estimators. For usage of these estimators see documentation for the child classes.

Args:
settingsdict [optional]

set estimator parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

is_parallel()[source]

Indicate if estimator supports parallel estimation over chunks.

Return true if the supports parallel estimation over chunks, where a chunk is one independent data set.

Returns:

bool

class idtxl.estimators_jidt.JidtGaussian(CalcClass, settings)[source]

Abstract class for implementation of JIDT Gaussian-estimators.

Abstract class for implementation of JIDT Gaussian-estimators, child classes implement estimators for mutual information (MI), conditional mutual information (CMI), actice information storage (AIS), transfer entropy (TE) using JIDT’s Gaussian estimator for continuous data. See parent class for references.

Set common estimation parameters for JIDT Kraskov-estimators. For usage of these estimators see documentation for the child classes.

Results are returned in nats.

Args:
CalcClassJAVA class

JAVA class returned by jpype.JPackage

settingsdict [optional]

set estimator parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

estimate_surrogates_analytic(n_perm=200, **data)[source]

Estimate the surrogate distribution analytically. This method must be implemented because this class’ is_analytic_null_estimator() method returns true

Args:
n_permsint

number of permutations (default=200)

datanumpy arrays

realisations of random variables required for the calculation (varies between estimators, e.g. 2 variables for MI, 3 for CMI). Formatted as per estimate_parallel for this estimator.

Returns:
float | numpy array

n_perm surrogates of the average MI/CMI/TE over all samples under the null hypothesis of no relationship between var1 and var2 (in the context of conditional)

get_analytic_distribution(**data)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
datanumpy arrays

realisations of random variables required for the calculation (varies between estimators, e.g. 2 variables for MI, 3 for CMI). Formatted as per the estimate method for this estimator.

Returns:
Java object

JIDT calculator that was used here

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

class idtxl.estimators_jidt.JidtGaussianAIS(settings)[source]

Calculate active information storage with JIDT’s Gaussian implementation.

Calculate active information storage (AIS) for some process using JIDT’s implementation of the Gaussian estimator. AIS is defined as the mutual information between the processes’ past state and current value.

The past state needs to be defined in the settings dictionary, where a past state is defined as a uniform embedding with parameters history and tau. The history describes the number of samples taken from a processes’ past, tau describes the embedding delay, i.e., the spacing between every two samples from the processes’ past.

See parent class for references.Results are returned in nats.

Args:
settingsdict

sets estimation parameters:

  • history : int - number of samples in the processes’ past used as embedding

  • tau : int [optional] - the processes’ embedding delay (default=1)

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the AIS estimator to save computation time. The Theiler window ignores trial boundaries. The AIS estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(process)[source]

Estimate active information storage.

Args:
processnumpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

Returns:
float | numpy array

average AIS over all samples or local AIS for individual samples if ‘local_values’=True

class idtxl.estimators_jidt.JidtGaussianCMI(settings=None)[source]

Calculate conditional mutual infor with JIDT’s Gaussian implementation.

Computes the differential conditional mutual information of two multivariate sets of observations, conditioned on another, assuming that the probability distribution function for these observations is a multivariate Gaussian distribution. Call JIDT via jpype and use ConditionalMutualInfoCalculatorMultiVariateGaussian estimator. If no conditional is given (is None), the function returns the mutual information between var1 and var2.

See parent class for references. Results are returned in nats.

Args:
settingsdict [optional]

sets estimation parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the CMI estimator to save computation time. The Theiler window ignores trial boundaries. The CMI estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(var1, var2, conditional=None)[source]

Estimate conditional mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of the second variable (similar to var1)

conditionalnumpy array [optional]

realisations of the conditioning variable (similar to var), if no conditional is provided, return MI between var1 and var2

Returns:
float | numpy array

average CMI over all samples or local CMI for individual samples if ‘local_values’=True

get_analytic_distribution(var1, var2, conditional=None)[source]

Return a JIDT AnalyticNullDistribution object.

Required so that our estimate_surrogates_analytic method can use the common_estimate_surrogates_analytic() method, where data is formatted as per the estimate method for this estimator.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of the second variable (similar to var1)

conditionalnumpy array [optional]

realisations of the conditioning variable (similar to var), if no conditional is provided, return MI between var1 and var2

Returns:
Java object

JIDT calculator that was used here

class idtxl.estimators_jidt.JidtGaussianMI(settings=None)[source]

Calculate mutual information with JIDT’s Gaussian implementation.

Calculate the mutual information between two variables. Call JIDT via jpype and use the Gaussian estimator. See parent class for references.

Results are returned in nats.

Args:
settingsdict [optional]

sets estimation parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • lag_mi : int [optional] - time difference in samples to calculate the lagged MI between processes (default=0)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the MI estimator to save computation time. The Theiler window ignores trial boundaries. The MI estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(var1, var2)[source]

Estimate mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of the second variable (similar to var1)

Returns:
float | numpy array

average MI over all samples or local MI for individual samples if ‘local_values’=True

class idtxl.estimators_jidt.JidtGaussianTE(settings)[source]

Calculate transfer entropy with JIDT’s Gaussian implementation.

Calculate transfer entropy between a source and a target variable using JIDT’s implementation of the Gaussian estimator. Transfer entropy is defined as the conditional mutual information between the source’s past state and the target’s current value, conditional on the target’s past.

Past states need to be defined in the settings dictionary, where a past state is defined as a uniform embedding with parameters history and tau. The history describes the number of samples taken from a variable’s past, tau descrices the embedding delay, i.e., the spacing between every two samples from the processes’ past.

See parent class for references. Results are returned in nats.

Args:
settingsdict

sets estimation parameters:

  • history_target : int - number of samples in the target’s past used as embedding

  • history_source : int [optional] - number of samples in the source’s past used as embedding (default=same as the target history)

  • tau_source : int [optional] - source’s embedding delay (default=1)

  • tau_target : int [optional] - target’s embedding delay (default=1)

  • source_target_delay : int [optional] - information transfer delay between source and target (default=1)

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the CMI estimator to save computation time. The Theiler window ignores trial boundaries. The CMI estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(source, target)[source]

Estimate transfer entropy from a source to a target variable.

Args:
sourcenumpy array

realisations of source variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of target variable (similar to var1)

Returns:
float | numpy array

average TE over all samples or local TE for individual samples if ‘local_values’=True

class idtxl.estimators_jidt.JidtKraskov(CalcClass, settings=None)[source]

Abstract class for implementation of JIDT Kraskov-estimators.

Abstract class for implementation of JIDT Kraskov-estimators, child classes implement estimators for mutual information (MI), conditional mutual information (CMI), actice information storage (AIS), transfer entropy (TE) using the Kraskov-Grassberger-Stoegbauer estimator for continuous data. See parent class for references.

Set common estimation parameters for JIDT Kraskov-estimators. For usage of these estimators see documentation for the child classes.

Args:
CalcClassJAVA class

JAVA class returned by jpype.JPackage

settingsdict [optional]

set estimator parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • num_threads : int | str [optional] - number of threads used for estimation (default=’USE_ALL’, note that this uses all available threads on the current machine)

  • algorithm_num : int [optional] - which Kraskov algorithm (1 or 2) to use (default=1). Only applied at this method for TE and AIS (is already applied for MI/CMI). Note that default algorithm of 1 here is different to the default ALG_NUM argument for the JIDT AIS KSG estimator.

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

class idtxl.estimators_jidt.JidtKraskovAIS(settings)[source]

Calculate active information storage with JIDT’s Kraskov implementation.

Calculate active information storage (AIS) for some process using JIDT’s implementation of the Kraskov type 1 estimator. AIS is defined as the mutual information between the processes’ past state and current value.

The past state needs to be defined in the settings dictionary, where a past state is defined as a uniform embedding with parameters history and tau. The history describes the number of samples taken from a processes’ past, tau describes the embedding delay, i.e., the spacing between every two samples from the processes’ past.

See parent class for references. Results are returned in nats.

Args:
settingsdict

sets estimation parameters:

  • history : int - number of samples in the processes’ past used as embedding

  • tau : int [optional] - the processes’ embedding delay (default=1)

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • num_threads : int | str [optional] - number of threads used for estimation (default=’USE_ALL’, note that this uses all available threads on the current machine)

  • algorithm_num : int [optional] - which Kraskov algorithm (1 or 2) to use (default=1)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the AIS estimator to save computation time. The Theiler window ignores trial boundaries. The AIS estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(process)[source]

Estimate active information storage.

Args:
processnumpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

Returns:
float | numpy array

average AIS over all samples or local AIS for individual samples if ‘local_values’=True

class idtxl.estimators_jidt.JidtKraskovCMI(settings=None)[source]

Calculate conditional mutual inform with JIDT’s Kraskov implementation.

Calculate the conditional mutual information (CMI) between three variables. Call JIDT via jpype and use the Kraskov 1 estimator. If no conditional is given (is None), the function returns the mutual information between var1 and var2. See parent class for references.

Results are returned in nats.

Args:
settingsdict [optional]

set estimator parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • num_threads : int | str [optional] - number of threads used for estimation (default=’USE_ALL’, note that this uses all available threads on the current machine)

  • algorithm_num : int [optional] - which Kraskov algorithm (1 or 2) to use (default=1)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the CMI estimator to save computation time. The Theiler window ignores trial boundaries. The CMI estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(var1, var2, conditional=None)[source]

Estimate conditional mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of the second variable (similar to var1)

conditionalnumpy array [optional]

realisations of the conditioning variable (similar to var), if no conditional is provided, return MI between var1 and var2

Returns:
float | numpy array

average CMI over all samples or local CMI for individual samples if ‘local_values’=True

class idtxl.estimators_jidt.JidtKraskovMI(settings=None)[source]

Calculate mutual information with JIDT’s Kraskov implementation.

Calculate the mutual information between two variables. Call JIDT via jpype and use the Kraskov 1 estimator. See parent class for references.

Results are returned in nats.

Args:
settingsdict [optional]

sets estimation parameters:

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • num_threads : int | str [optional] - number of threads used for estimation (default=’USE_ALL’, note that this uses all available threads on the current machine)

  • algorithm_num : int [optional] - which Kraskov algorithm (1 or 2) to use (default=1)

  • lag_mi : int [optional] - time difference in samples to calculate the lagged MI between processes (default=0)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the MI estimator to save computation time. The Theiler window ignores trial boundaries. The MI estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(var1, var2)[source]

Estimate mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of the second variable (similar to var1)

Returns:
float | numpy array

average MI over all samples or local MI for individual samples if ‘local_values’=True

class idtxl.estimators_jidt.JidtKraskovTE(settings)[source]

Calculate transfer entropy with JIDT’s Kraskov implementation.

Calculate transfer entropy between a source and a target variable using JIDT’s implementation of the Kraskov type 1 estimator. Transfer entropy is defined as the conditional mutual information between the source’s past state and the target’s current value, conditional on the target’s past.

Past states need to be defined in the settings dictionary, where a past state is defined as a uniform embedding with parameters history and tau. The history describes the number of samples taken from a variable’s past, tau descrices the embedding delay, i.e., the spacing between every two samples from the processes’ past.

See parent class for references. Results are returned in nats.

Args:
settingsdict

sets estimation parameters:

  • history_target : int - number of samples in the target’s past used as embedding

  • history_source : int [optional] - number of samples in the source’s past used as embedding (default=same as the target history)

  • tau_source : int [optional] - source’s embedding delay (default=1)

  • tau_target : int [optional] - target’s embedding delay (default=1)

  • source_target_delay : int [optional] - information transfer delay between source and target (default=1)

  • debug : bool [optional] - return debug information when calling JIDT (default=False)

  • local_values : bool [optional] - return local TE instead of average TE (default=False)

  • algorithm_num : int [optional] - which Kraskov algorithm (1 or 2) to use (default=1)

Note:

Some technical details: JIDT normalises over realisations, IDTxl normalises over raw data once, outside the CMI estimator to save computation time. The Theiler window ignores trial boundaries. The CMI estimator does add noise to the data as a default. To make analysis runs replicable set noise_level to 0.

estimate(source, target)[source]

Estimate transfer entropy from a source to a target variable.

Args:
sourcenumpy array

realisations of source variable, either a 2D numpy array where array dimensions represent [realisations x variable dimension] or a 1D array representing [realisations]

var2numpy array

realisations of target variable (similar to var1)

Returns:
float | numpy array

average TE over all samples or local TE for individual samples if ‘local_values’=True

idtxl.estimators_jidt.common_estimate_surrogates_analytic(estimator, n_perm=200, **data)[source]

Estimate the surrogate distribution analytically for JidtEstimator.

Estimate the surrogate distribution analytically for a JidtEstimator which is_analytic_null_estimator(), by sampling estimates at random p-values in the analytic distribution.

Args:
estimatora JidtEstimator object, which returns True to a call to

its is_analytic_null_estimator() method

n_permsint

number of permutations (default=200)

datanumpy arrays

realisations of random variables required for the calculation (varies between estimators, e.g. 2 variables for MI, 3 for CMI)

Returns:
float | numpy array

n_perm surrogates of the average MI/CMI/TE over all samples under the null hypothesis of no relationship between var1 and var2 (in the context of conditional)

Python Estimators (CPU)

class idtxl.estimators_python.PythonKraskovCMI(settings)[source]

Estimate conditional mutual information using Kraskov’s first estimator.

Args:
settingsdict [optional]

set estimator parameters:

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • base : float - base of returned values (default=np=e)

  • normalise : bool [optional] - z-standardise data (default=False)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • rng_seed : int | None [optional] - random seed if noise level > 0

  • num_threads : int | str [optional] - number of threads used for estimation (default=’USE_ALL’, note that this uses all available threads on the current machine)

  • knn_finder : str [optional] - knn algorithm to use, can be ‘scipy_kdtree’ (default), ‘sklearn_kdtree’, or ‘sklearn_balltree’

estimate(var1: numpy.ndarray, var2: numpy.ndarray, conditional=None)[source]

Estimate conditional mutual information between var1 and var2, given conditional.

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

is_parallel()[source]

Indicate if estimator supports parallel estimation over chunks.

Return true if the supports parallel estimation over chunks, where a chunk is one independent data set.

Returns:

bool

OpenCL Estimators (GPU)

class idtxl.estimators_opencl.OpenCLKraskov(settings=None)[source]

Abstract class for implementation of OpenCL estimators.

Abstract class for implementation of OpenCL estimators, child classes implement estimators for mutual information (MI) and conditional mutual information (CMI) using the Kraskov-Grassberger-Stoegbauer estimator for continuous data.

References:

  • Kraskov, A., Stoegbauer, H., & Grassberger, P. (2004). Estimating mutual information. Phys Rev E, 69(6), 066138.

  • Lizier, Joseph T., Mikhail Prokopenko, and Albert Y. Zomaya. (2012). Local measures of information storage in complex distributed computation. Inform Sci, 208, 39-54.

  • Schreiber, T. (2000). Measuring information transfer. Phys Rev Lett, 85(2), 461.

Estimators can be used to perform multiple, independent searches in parallel. Each of these parallel searches is called a ‘chunk’. To search multiple chunks, provide point sets as 2D arrays, where the first dimension represents samples or points, and the second dimension represents the points’ dimensions. Concatenate chunk data in the first dimension and pass the number of chunks to the estimators. Chunks must be of equal size.

Set common estimation parameters for OpenCL estimators. For usage of these estimators see documentation for the child classes.

Args:
settingsdict [optional]

set estimator parameters:

  • gpuid : int [optional] - device ID used for estimation (if more than one device is available on the current platform) (default=0)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • padding : bool [optional] - pad data to a length that is a multiple of 1024, workaround for a

  • debug : bool [optional] - calculate intermediate results, i.e. neighbour counts from range searches and KNN distances, print debug output to console (default=False)

  • return_counts : bool [optional] - return intermediate results, i.e. neighbour counts from range searches and KNN distances (default=False)

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

is_parallel()[source]

Indicate if estimator supports parallel estimation over chunks.

Return true if the supports parallel estimation over chunks, where a chunk is one independent data set.

Returns:

bool

class idtxl.estimators_opencl.OpenCLKraskovCMI(settings=None)[source]

Calculate conditional mutual inform with OpenCL Kraskov implementation.

Calculate the conditional mutual information (CMI) between three variables using OpenCL GPU-code. If no conditional is given (is None), the function returns the mutual information between var1 and var2. See parent class for references.

Results are returned in nats.

Args:
settingsdict [optional]

set estimator parameters:

  • gpuid : int [optional] - device ID used for estimation (if more than one device is available on the current platform) (default=0)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • debug : bool [optional] - return intermediate results, i.e. neighbour counts from range searches and KNN distances (default=False)

  • return_counts : bool [optional] - return intermediate results, i.e. neighbour counts from range searches and KNN distances (default=False)

estimate(var1, var2, conditional=None, n_chunks=1)[source]

Estimate conditional mutual information.

If conditional is None, the mutual information between var1 and var2 is calculated.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [(realisations * n_chunks) x variable dimension] or a 1D array representing [realisations], array type should be int32

var2numpy array

realisations of the second variable (similar to var1)

conditionalnumpy array

realisations of conditioning variable (similar to var1)

n_chunksint

number of data chunks, no. data points has to be the same for each chunk

Returns:
float | numpy array

average CMI over all samples or local CMI for individual samples if ‘local_values’=True

numpy arrays

distances and neighborhood counts for var1 and var2 if debug=True and return_counts=True

class idtxl.estimators_opencl.OpenCLKraskovMI(settings=None)[source]

Calculate mutual information with OpenCL Kraskov implementation.

Calculate the mutual information (MI) between two variables using OpenCL GPU-code. See parent class for references.

Results are returned in nats.

Args:
settingsdict [optional]

set estimator parameters:

  • gpuid : int [optional] - device ID used for estimation (if more than one device is available on the current platform) (default=0)

  • kraskov_k : int [optional] - no. nearest neighbours for KNN search (default=4)

  • normalise : bool [optional] - z-standardise data (default=False)

  • theiler_t : int [optional] - no. next temporal neighbours ignored in KNN and range searches (default=0)

  • noise_level : float [optional] - random noise added to the data (default=1e-8)

  • debug : bool [optional] - return intermediate results, i.e. neighbour counts from range searches and KNN distances (default=False)

  • return_counts : bool [optional] - return intermediate results, i.e. neighbour counts from range searches and KNN distances (default=False)

  • lag_mi : int [optional] - time difference in samples to calculate the lagged MI between processes (default=0)

estimate(var1, var2, n_chunks=1)[source]

Estimate mutual information.

Args:
var1numpy array

realisations of first variable, either a 2D numpy array where array dimensions represent [(realisations * n_chunks) x variable dimension] or a 1D array representing [realisations], array type should be int32

var2numpy array

realisations of the second variable (similar to var1)

n_chunksint

number of data chunks, no. data points has to be the same for each chunk

Returns:
float | numpy array

average MI over all samples or local MI for individual samples if ‘local_values’=True

numpy arrays

distances and neighborhood counts for var1 and var2 if debug=True and return_counts=True

MPI Estimators (CPU)

class idtxl.estimators_mpi.MPIEstimator(est, settings)[source]

MPI Wrapper for arbitrary Estimator implementations

Make sure to have an “if __name__==’__main__’:” guard in your main script to avoid infinite recursion!

To use MPI, add MPI=True to the Estimator settings dictionary and optionally provide max_workers

Call using mpiexec:
>>> mpiexec -n 1 -usize <max workers + 1> python <python script>
or, if MPI does not support spawning new workers (i.e. MPI version < 2)
>>> mpiexec -n <max workers + 1> python -m mpi4py.futures <python script>
Call using slurm:
>>> srun -n $SLURM_NTASKS --mpi=pmi2 python -m mpi4py.futures <python script>
estimate(*, n_chunks=1, **data)[source]

Distributes the given chunks of a task to Estimators on worker ranks using MPI.

Needs to be called with kwargs only.

Args:
n_chunksint [optional]

Number of chunks to split the data into, default=1.

datadict[str, Sequence]

Dictionary of random variable realizations

Returns:
numpy array

Estimates of information-theoretic quantities as np.double values

estimate_surrogates_analytic(**data)[source]

Forward analytic estimation to the base Estimator.

Analytic estimation is assumed to have shorter runtime and is thus performed on rank 0 alone for now.

is_analytic_null_estimator()[source]

Test if the base Estimator is an analytic null estimator.

is_parallel()[source]

Indicate if estimator supports parallel estimation over chunks.

Return true if the supports parallel estimation over chunks, where a chunk is one independent data set.

Returns:

bool

PID Estimators

Partical information decomposition for discrete random variables.

This module provides an estimator for partial information decomposition as proposed in

Bertschinger, N., Rauh, J., Olbrich, E., Jost, J., & Ay, N. (2014). Quantifying Unique Information. Entropy, 16(4), 2161–2183. http://doi.org/10.3390/e16042161

class idtxl.estimators_pid.SydneyPID(settings)[source]

Estimate partial information decomposition of discrete variables.

Fast implementation of the BROJA partial information decomposition (PID) estimator for discrete data (Bertschinger, 2014). The estimator does not require JAVA or GPU modules to run.

The estimator finds shared information, unique information and synergistic information between the two inputs s1 and s2 with respect to the output t.

Improved version with larger initial swaps and checking for convergence of both the unique information from sources 1 and 2. The function counts the empirical observations, calculates probabilities and the initial CMI, then does the vitrualised swaps until it has converged, and finally calculates the PID. The virtualised swaps stage contains two loops. An inner loop which actually does the virtualised swapping, keeping the changes if the CMI decreases; and an outer loop which decreases the size of the probability mass increment the virtualised swapping utilises.

References

  • Bertschinger, N., Rauh, J., Olbrich, E., Jost, J., & Ay, N. (2014). Quantifying unique information. Entropy, 16(4), 2161–2183. http://doi.org/10.3390/e16042161

Args:
settingsdict

estimation parameters

  • alph_s1 : int - alphabet size of s1

  • alph_s2 : int - alphabet size of s2

  • alph_t : int - alphabet size of t

  • max_unsuc_swaps_row_parm : int - soft limit for virtualised swaps based on the number of unsuccessful swaps attempted in a row. If there are too many unsuccessful swaps in a row, then it will break the inner swap loop; the outer loop decrements the size of the probability mass increment and then attemps virtualised swaps again with the smaller probability increment. The exact number of unsuccessful swaps allowed before breaking is the total number of possible swaps (given our alphabet sizes) times the control parameter max_unsuc_swaps_row_parm, e.g., if the parameter is set to 3, this gives a high degree of confidence that nearly (if not) all of the possible swaps have been attempted before this soft limit breaks the swap loop.

  • num_reps : int - number of times the outer loop will halve the size of the probability increment used for the virtualised swaps. This is in direct correspondence with the number of times the empirical data was replicated in your original implementation.

  • max_iters : int - provides a hard upper bound on the number of times it will attempt to perform virtualised swaps in the inner loop. However, this hard limit is (practically) never used as it should always hit the soft limit defined above (parameter may be removed in the future).

  • verbose : bool [optional] - print output to console (default=False)

estimate(s1, s2, t)[source]
Args:
s1numpy array

1D array containing realizations of a discrete random variable

s2numpy array

1D array containing realizations of a discrete random variable

tnumpy array

1D array containing realizations of a discrete random variable

Returns:
dict

estimated decomposition, contains the joint distribution, unique, shared, and synergistic information

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

is_parallel()[source]

Indicate if estimator supports parallel estimation over chunks.

Return true if the supports parallel estimation over chunks, where a chunk is one independent data set.

Returns:

bool

class idtxl.estimators_pid.TartuPID(settings)[source]

Estimate partial information decomposition for two inputs and one output

Implementation of the partial information decomposition (PID) estimator for discrete data. The estimator finds shared information, unique information and synergistic information between the two inputs s1 and s2 with respect to the output t.

The algorithm uses exponential cone programming and requires the Python package for ECOS: Embedded Cone Solver (https://pypi.python.org/pypi/ecos).

References:

  • Makkeh, A., Theis, D.O., & Vicente, R. (2017). Bivariate Partial Information Decomposition: The Optimization Perspective. Entropy, 19(10), 530.

  • Makkeh, A., Theis, D.O., & Vicente, R. (2018). BROJA-2PID: A cone programming based Partial Information Decomposition estimator. Entropy, 20(271), https://github.com/Abzinger/BROJA_2PID.

Args:
settingsdict

estimation parameters (with default parameters)

  • verbose : bool [optional] - print output to console (default=False)

  • cone_solver : str [optional] - which cone solver to use (default=’ECOS’)

  • solver_args : dict [optional] - solver arguments (default={})

estimate(s1, s2, t)[source]
Args:
s1numpy array

1D array containing realizations of a discrete random variable

s2numpy array

1D array containing realizations of a discrete random variable

tnumpy array

1D array containing realizations of a discrete random variable

Returns:
dict

estimated decomposition, solver used, numerical error

is_analytic_null_estimator()[source]

Indicate if estimator supports analytic surrogates.

Return true if the estimator implements estimate_surrogates_analytic() where data is formatted as per the estimate method for this estimator.

Returns:

bool

is_parallel()[source]

Indicate if estimator supports parallel estimation over chunks.

Return true if the supports parallel estimation over chunks, where a chunk is one independent data set.

Returns:

bool