Python API

Module eos

Basic Classes

EOS provides its basic functionality via the main eos module.

class eos.Analysis(priors, likelihood, external_likelihood=[], global_options={}, manual_constraints={}, fixed_parameters={}, parameters=None)

Represents a statistical analysis.

Describes a Bayesian analysis in terms of a set of parameters, a log(likelihood), and a set containing one or more log(prior)s.

Parameters:
  • global_options (dict, optional) – The options as (key, value) pairs that shall be forwarded to all theory predictions.

  • priors (iterable) – The priors for this analysis as a list of prior descriptions. See below for what consitutes a valid prior description.

  • likelihood (iterable) – The likelihood as a list of individual constraints from the internal data base of experimental and theoretical constraints; cf. the complete list of constraints.

  • external_likelihood (list or iterable of eos.LogLikelihoodBlock.) – The external likelihood blocks as a list or iterable of objects returned by eos.LogLikelihoodBlock.External().

  • manual_constraints (dict, optional) – Additional manually-specified constraints that shall be added to the log(likelihood).

  • fixed_parameters (dict, optional) – Values of parameters that are set when the analysis is defined.

  • parameters (eos.Parameters or None, optional) – The optional set of parameters that shall be used for this analysis. Defaults to None which means that a new instance of eos.Parameters is created.

Each prior description is a dictionary with the following mandatory elements:

  • type (str) – The type specification of the prior. Must be one of uniform, flat, or gaussian.

  • parameter (str) – The name of the parameter for which the prior shall apply.

  • min (float) – The lower boundary of the prior’s support.

  • max (float) – The upper boundary of the prior’s support.

A uniform or flat prior does not require any further description. A gaussian prior requires in addition providing the following two elements:

  • central (float) – The median value of the parameter.

  • sigma (float, or list, tuple of float) – The width of the 68% probability interval. If a list or tuple of two numbers is provided, the prior will by a asymmetric but continuous. The two values are then taken to be the distance to the lower and upper end of the 68% probability interval.

clone()

Returns an independent instance of eos.Analysis.

goodness_of_fit()

Returns a GoodnessOfFit object that summarizes the quality of the fit for the current parameter point.

log_likelihood(p, *args)

Adapter for use with external sampling software (e.g. dynesty) to aid when sampling from the log(likelihood).

Parameters:
  • p (iterable) – Parameter point, with the elements in the same order as in eos.Analysis.varied_parameters.

  • args (optional) – Dummy parameter (ignored)

log_pdf(u, *args)

Adapter for use with external optimization software (e.g. pypmc) to aid when optimizing the log(posterior).

Parameters:
  • u (iterable) – Parameter point in u space, with the elements in the same order as in eos.Analysis.varied_parameters.

  • args (optional) – Dummy parameter (ignored)

negative_log_pdf(u, *args)

Adapter for use with external optimization software (e.g. scipy.optimize.minimize) to aid when optimizing the log(posterior).

Parameters:
  • u (iterable) – Parameter point in u space, with the elements in the same order as in eos.Analysis.varied_parameters.

  • args (optional) – Dummy parameter (ignored)

optimize(start_point=None, rng=<module 'numpy.random.mtrand' from '/opt/venv/lib/python3.10/site-packages/numpy/random/mtrand.cpython-310-x86_64-linux-gnu.so'>, **kwargs)

Optimize the log(posterior) and returns a best-fit-point summary.

Parameters:
  • start_point (iterable, optional) – Parameter point from which to start the optimization, with the elements in the same order as in eos.Analysis.varied_parameters. If set to “random”, optimization starts at the random point in the space of the priors. If not specified, optimization starts at the current parameter point.

  • rng – Optional random number generator

  • **kwargs – Are passed to scipy.optimize.minimize

parameters

The set of parameters used for this analysis.

reset_parameters()

Resets the analysis parameters to their default values.

sample(N=1000, stride=5, pre_N=150, preruns=3, cov_scale=0.1, observables=None, start_point=None, rng=<module 'numpy.random.mtrand' from '/opt/venv/lib/python3.10/site-packages/numpy/random/mtrand.cpython-310-x86_64-linux-gnu.so'>, return_uspace=False)

Return samples of the parameters, log(weights), and optionally posterior-predictive samples for a sequence of observables.

Obtains random samples of the log(posterior) using an adaptive Markov Chain Monte Carlo with PyPMC. A prerun with adaptations is carried out first and its samples are discarded.

Parameters:
  • N – Number of samples that shall be returned

  • stride – Stride, i.e., the number by which the actual amount of samples shall be thinned to return N samples.

  • pre_N – Number of samples in each prerun.

  • preruns – Number of preruns.

  • cov_scale – Scale factor for the initial guess of the covariance matrix.

  • observables (list-like, optional) – Observables for which posterior-predictive samples shall be obtained.

  • start_point (list-like, optional) – Optional starting point for the chain

  • rng – Optional random number generator (must be compatible with the requirements of pypmc.sampler.markov_chain.MarkovChain)

Returns:

A tuple of the parameters as array of size N, the logarithmic weights as array of size N, and optionally the posterior-predictive samples of the observables as array of size N x len(observables).

Note

This method requiries the PyPMC python module, which can be installed from PyPI.

sample_nested(bound='multi', nlive=250, dlogz=1.0, maxiter=None, seed=10, print_progress=True)

Return samples of the parameters.

Obtains random samples of log(likelihood) using dynamic nested sampling with dynesty.

Parameters:
  • bound (str, optional) – The option for bounding the target distribution. For valid values, see the dynesty documentation. Defaults to ‘multi’.

  • nlive (int, optional) – The number of live points.

  • dlogz (float, optional) – The relative tolerance for the remaining evidence. Defaults to 1.0.

  • maxiter (int, optional) – The maximum number of iterations. Iterations may stop earlier if the termination condition is reached.

  • seed ({None, int, array_like[ints], SeedSequence}, optional) – The seed used to initialize the Mersenne Twister pseudo-random number generator.

Note

This method requires the dynesty python module, which can be installed from PyPI.

sample_pmc(log_proposal, step_N=1000, steps=10, final_N=5000, rng=<module 'numpy.random.mtrand' from '/opt/venv/lib/python3.10/site-packages/numpy/random/mtrand.cpython-310-x86_64-linux-gnu.so'>, return_final_only=True, final_perplexity_threshold=1.0, weight_threshold=1e-10, pmc_iterations=1, pmc_rel_tol=1e-10, pmc_abs_tol=1e-05, pmc_lookback=1)

Return samples of the parameters and log(weights), and a mixture density adapted to the posterior.

Obtains random samples of the log(posterior) using adaptive importance sampling following the Population Monte Carlo approach with PyPMC.

Parameters:
  • log_proposal (pypmc.density.mixture.MixtureDensity) – Initial gaussian mixture density that shall be adapted to the posterior density.

  • step_N (int) – Number of samples that shall be drawn in each adaptation step.

  • steps (int) – Number of adaptation steps.

  • final_N (int) – Number of samples that shall be drawn after all adaptation steps.

  • rng – Optional random number generator (must be compatible with the requirements of pypmc.sampler.importance_sampler.ImportanceSampler)

  • return_final_only – If set to True, only returns the samples and weights of the final sampling step, after all adaptations have finished.

  • final_perplexity_threshold – Adaptations are stopped if the perplexity of the last adaptation step is above this threshold value.

  • weight_threshold – Mixture components with a weight smaller than this threshold are pruned.

  • pmc_iterations – (advanced) Maximum number of update of the PMC, changing this value may make the update unstable.

  • pmc_rel_tol – (advanced) Relative tolerance of the PMC. If two consecutive values of the current density log-likelihood are relatively smaller than this value, the convergence is declared.

  • pmc_abs_tol – (advanced) Absolute tolerance of the PMC. If two consecutive values of the current density log-likelihood are smaller than this value, the convergence is declared.

  • pmc_lookback – (advanced) Use reweighted samples from the previous update steps when adjusting the mixture density. The parameter determines the number of update steps to “look back”. The default value of 1 disables this feature, a value of 0 means that all previous steps are used.

Returns:

A tuple of the parameters as array of length N = step_N * steps + final_N, the (linear) weights as array of length N, the posterior values as array of length N, and the final proposal function as pypmc.density.mixture.MixtureDensity.

This method should be called after obtaining approximate samples of the log(posterior) by other means, e.g., by using eos.Analysis.sample(). A possible (incomplete) example could look as follows:

from pypmc.mix_adapt.r_value import make_r_gaussmix
chains = []
for i in range(10):
    # run Markov Chains for your problem
    chain, _ = analysis.sample(...)
    chains.append(chain)

# please consult the pypmc documentation for details on the call below
proposal_density = make_r_gaussmix(chains, K_g=3, critical_r=1.1)

# adapt the proposal to the posterior and obtain high-quality samples
analysis.sample_pmc(proposal_density, ...)

Note

This method requires the PyPMC python module, which can be installed from PyPI.

class eos.BestFitPoint(analysis, point)

Represents the best-fit point of a Bayesian analysis undertaken with the Analysis class.

class eos.Constraints(prefix=None, name=None, suffix=None)
insert((_Constraints)arg1, (QualifiedName)arg2, (str)arg3) ConstraintEntry
class eos.GoodnessOfFit

Represents the goodness of fit characteristics of the log(posterior).

total_chi_square((GoodnessOfFit)arg1) float :

Returns the total \(\chi^2\) value of the log(likelihood). Only (multivariate) gaussian likelihoods are considered for this result.

total_degrees_of_freedom((GoodnessOfFit)arg1) int :

Returns the total number of degrees of freedom in the log(posterior).

class eos.Kinematics

Represents the set of kinematic variables relevant to an observable.

Initialize a new set of kinematic variables. The inital set of variables and their initial set of values can be provided through keyword arguments, e.g. using

k = eos.Kinematics(q2=0.4, k2=0.0)                      # default keyword arguments
k = eos.Kinematics({'q2': 0.4, 'cos(theta_l)': -1.0})   # use a dictionary if variable names are not
                                                        # valid python identifiers
declare((Kinematics)self, (str)name, (float)value) KinematicVariable :

Declares a new kinematic variable.

Parameters:
  • name (str) – The name of the new kinematic variable.

  • value (float) – The initial value for the new kinematic variable.

class eos.LogLikelihood

Represents the log(likelihood) of a Bayesian analysis undertaken with the Analysis class.

add((LogLikelihood)arg1, (Constraint)arg2) None

add( (LogLikelihood)arg1, (LogLikelihoodBlock)arg2) -> None

evaluate((LogLikelihood)arg1) float
observable_cache((LogLikelihood)arg1) ObservableCache
class eos.LogPrior

Represents a Bayesian prior on the log scale.

New LogPrior objects can only be created using the capitalized static methods: LogPrior.Uniform(), LogPrior.Gaussian(), and LogPrior.Scale().

static CurtailedGauss((_Parameters)arg1, (str)parameters, (float)name, (float)range, (float)lower, (float)central, (float)upper) LogPrior :

Returns a new (curtailed) Gaussian prior as a LogPrior.

The prior’s support is provided by the pair of min and max parameters, with the approximate 68% probability interval [lower, upper] and the mode provided by the parameter central.

Parameters:
  • parameters (eos.Parameters) – The parameters to which this LogPrior is bound.

  • name (str) – The name of the parameter for which the LogPrior is defined.

  • min (float) – The minimum value that the parameter is allowed to take.

  • max (float) – The maximum value that the parameter is allowed to take.

  • lower (float) – The lower boundary of the 68% probability interval.

  • central (float) – The mode and median of the prior.

  • upper (float) – The upper boundary of the 68% probability interval.

static Flat((_Parameters)parameters, (str)name, (float)min, (float)max) LogPrior :

Alias for LogPrior.Uniform().

static Gaussian((_Parameters)parameters, (QualifiedName)name, (float)mu, (float)sigma) LogPrior :

Returns a new Gaussian prior as a LogPrior.

The priors support is infinite. The mode is provided by the parameter mu and the standard deviation by the parameter sigma.

Parameters:
  • parameters (eos.Parameters) – The parameters to which this LogPrior is bound.

  • name (str) – The name of the parameter for which the LogPrior is defined.

  • mu (float) – The mode of the prior.

  • sigma (float, strictly positive) – The standard deviation of the prior.

Poisson((_Parameters)parameters, (str)name, (float)k) LogPrior :

Returns a new Poisson prior as a LogPrior.

The priors support is infinite. The mode is provided by the parameter k.

Parameters:
  • parameters (eos.Parameters) – The parameters to which this LogPrior is bound.

  • name (str) – The name of the parameter for which the LogPrior is defined.

  • k (float, strictly positive) – The mode of the prior.

static Scale((_Parameters)parameters, (str)name, (float)min, (float)max, (float)mu_0, (float)scale) LogPrior :

Returns a new Scale prior as a LogPrior.

The prior’s support is provided by the pair of min and max parameters, which should coincide with [mu_0 / lambda, mu_0 * lambda]. The PDF is chosen such that a renormalization scale is varied in this range and with central value mu_0 such that \(\ln x / \mu_0\) is uniformly distributed in the interval \([-\ln \lambda, +\ln \lambda]\).

Parameters:
  • parameters (eos.Parameters) – The parameters to which this LogPrior is bound.

  • name (str) – The name of the parameter for which the LogPrior is defined.

  • min (float) – The minimum value that the parameter is allowed to take.

  • max (float) – The maximum value that the parameter is allowed to take.

  • mu_0 (float, strictly positive) – The central value of the parameter.

  • lambda (float, strictly positive) – The scale factor.

static Uniform((_Parameters)parameters, (str)name, (float)min, (float)max) LogPrior :

Returns a new uniform prior as a LogPrior.

The prior’s support is provided by the range parameter.

Parameters:
  • parameters (eos.Parameters) – The parameters to which this LogPrior is bound.

  • name (str) – The name of the parameter for which the LogPrior is defined.

  • min (float) – The minimum value that the parameter is allowed to take.

  • max (float) – The maximum value that the parameter is allowed to take.

compute_cdf((LogPrior)arg1) None :

Returns the cumulative probabilities \(p\) assigned to each parameter via its Parameter.evaluate_generator() method.

evaluate((LogPrior)arg1) float :

Returns the logarithm of the prior’s probability density at the current parameter values.

sample((LogPrior)arg1) None :

Sets its parameters’ values corresponding to the cumulative propability \(p\) assigned to each parameter via its Parameter.set_generator() method.

varied_parameters((object)arg1) object
class eos.Observable

Represents an observable or pseudo observable known to EOS.

New observable objects are created using the make static method. See also the complete list of observables.

evaluate((Observable)self) float :

Evaluates the observable for the present values of its bound set of parameters and set of kinematic variables.

Returns:

The value of the observable.

Return type:

float

kinematics((Observable)arg1) Kinematics :

Returns the set of kinematic variables bound to this observable.

static make((QualifiedName)name, (_Parameters)parameters, (Kinematics)kinematics, (Options)options) Observable :

Makes a new Observable object.

Parameters:
Returns:

The new observable object.

Return type:

eos.Observable

name((Observable)arg1) QualifiedName :

Returns the name of the observable.

options((Observable)arg1) Options :

Returns the set of options used when creating the observable.

parameters((Observable)arg1) _Parameters :

Returns the set of parameters bound to this observable.

class eos.Observables(prefix=None, name=None, suffix=None, showall=False)

Represents the complete list of observables known to EOS.

Objects of this class are visualized as tables in Jupyter notebooks for easy overview. Filters can be applied through keyword arguments to the initialization.

Parameters:
  • prefix (str) – Only show observables whose qualified names contain the provided prefix in their prefix part.

  • name (str) – Only show observables whose qualified names contain the provided name in their name part.

  • suffix (str) – Only show observables whose qualified names contain the provided suffix in their suffix part.

See also the complete list of observables in this documentation.

insert((_Observables)arg1, (QualifiedName)name, (str)latex, (Unit)unit, (Options)options, (str)expression) None :

Insert a new observable to EOS by parsing the input string.

Parameters:
  • name (eos.QualifiedName) – The name of the new observable.

  • latex (std::string) – The latex representation of the new observable.

  • unit (eos.Unit) – The unit of the new observable.

  • options (eos.Options) – The set of options relevant to this new observable. These options apply to all the observables in the expression.

  • expression (std::string) – The expression to be parsed.

sections((object)arg1) object
class eos.ObservableCache

Provides a cache for the efficient evaluation of observables.

add((ObservableCache)arg1, (Observable)observable) int :

Add an existing observable to the cache.

Parameters:

observable (eos.Observable) – The observable to add to the cache.

Returns:

An internal handle to the cached observable. The observable’s value can be retrieved using cache[handle].

Return type:

int

parameters((ObservableCache)arg1) _Parameters :

Retrieve the set of parameters bound to this cache.

update((ObservableCache)arg1) None :

Update the cache for the current parameter point.

class eos.Options

Represents the set of options provided to an observable.

Options are pairs of (key, value) pairs. The list of valid keys and their respective valid options are specific to each observable. The initialization accepts keyword arguments, e.g.:

o = eos.Options(model='WET')                   # default keyword arguments
o = eos.Options({'form-factors': 'BSZ2015'})   # use a dictionary if option keys are not
                                               # valid python identifiers
declare((Options)arg1, (str)arg2, (str)arg3) None
class eos.Parameter

Represents a single real-valued scalar parameter in EOS.

Users cannot directly create new objects of this class. However, new named sets of parameters can be created, via the eos.Parameters class, from which the parameter of interest can be extracted, inspected, and altered.

central((Parameter)arg1) float
evaluate((Parameter)arg1) float :

Return the current value of a parameter.

evaluate_generator((Parameter)arg1) float :

Return the current generator value of a parameter.

latex((Parameter)arg1) str :

Returns the LaTeX representation of the parameter.

max((Parameter)arg1) float
min((Parameter)arg1) float
name((Parameter)arg1) str :

Returns the name of the parameter.

set((Parameter)arg1, (float)arg2) None :

Set the value of a parameter.

Parameters:

value (float) – The value to set the parameter to.

set_generator((Parameter)arg1, (float)arg2) None :

Set the generator value of a parameter.

Parameters:

value (float) – The value to set the parameter’s generator to.

set_max((Parameter)arg1, (float)arg2) None
set_min((Parameter)arg1, (float)arg2) None
unit((Parameter)arg1) Unit
class eos.Parameters(*args, **kwargs)

Represents the set of parameters known to EOS.

Objects of this class are visualized as tables in Jupyter notebooks for easy overview. Filters can be applied through keyword arguments to the initialization. An independent instance of the default parameters is obtained by each call to the constructor.

Parameters:
  • prefix (str) – Only show parameters whose qualified names contain the provided prefix in their prefix part.

  • name (str) – Only show parameters whose qualified names contain the provided name in their name part.

  • suffix (str) – Only show parameters whose qualified names contain the provided suffix in their suffix part.

See also the complete list of parameters.

static Defaults() _Parameters
static FromWCxf(w)

Imports Wilson coefficients into an eos.Parameters object from a wilson.Wilson object.

Parameters:

w (wilson.Wilson) – WET Wilson coefficients in the EOS basis.

by_id((_Parameters)arg1, (int)arg2) Parameter
static declare((QualifiedName)name, (str)latex, (Unit)unit, (float)value, (float)min, (float)max) int :

Declare a new parameter as part of the default parameter set.

Parameters:
  • name (eos.QualifiedName) – The name of the parameter to declare.

  • latex (str) – The LaTeX representation for the parameter.

  • unit (eos.Unit) – The unit in which the parameter is expressed.

  • value (float) – The initial value for the parameter.

  • min (float) – The minimum value that the parameter can attain.

  • max (float) – The maximum value that the parameter can attain.

declare_and_insert((_Parameters)arg1, (QualifiedName)name, (str)latex, (Unit)unit, (float)value, (float)min, (float)max) Parameter :

Declare a new parameter as part of the default parameter set and insert it into this parameter set.

Parameters:
  • name (eos.QualifiedName) – The name of the parameter to declare.

  • latex (str) – The LaTeX representation for the parameter.

  • unit (eos.Unit) – The unit in which the parameter is expressed.

  • value (float) – The initial value for the parameter.

  • min (float) – The minimum value that the parameter can attain.

  • max (float) – The maximum value that the parameter can attain.

Returns:

The newly inserted parameter.

Return type:

eos.Parameter

dump(file, **kwargs)

Dumps an eos.Parameters object to a YAML file.

Parameters:
  • file (str) – Name of the file to which the parameters shall be written.

  • names (iterable of str (optional)) – Names of the parameters that shall be converted.

has((_Parameters)arg1, (QualifiedName)arg2) bool
override_from_file((_Parameters)arg1, (str)arg2) None
static redirect((QualifiedName)name, (int)id) None :

Redirect a parameter name to a different parameter id in the default set of parameters.

The internal mapping of the parameter name will be redirected to the new id. If the the parameter’s previous id is not already aliased, it will become inaccessible. This is useful for example to alias a parameter name to a different parameter object.

Parameters:
  • name – The name of the parameter to be redirected.

  • id (eos.Parameter::Id) – The id of the parameter to which the name shall be redirected.

sections((object)arg1) object
set((_Parameters)arg1, (QualifiedName)arg2, (float)arg3) None :

Set the value of a parameter.

Parameters:
  • name (str) – The name of the parameter to set.

  • value (float) – The value to set the parameter to.

to_yaml(names=None)

Converts an eos.Parameters object to a YAML representation.

Parameters:

names – Names of the parameters that shall be converted.

class eos.QualifiedName

Represent a qualified (i.e. complete and syntactically correct) name.

EOS uses qualified names when naming any observable or constraint. The composition is approximately:

PREFIX::NAME@SUFFIX;OPTIONS
full((QualifiedName)arg1) str :

Returns the full name, i.e., the concatenation of all parts.

name_part((QualifiedName)arg1) qnpName :

Returns the name part of the name, i.e., the part following the ‘::’ and preceeding any optional ‘@’.

options_part((QualifiedName)arg1) Options :

Returns the optional options part of the name, i.e., the part following the optional ‘;’.

prefix_part((QualifiedName)arg1) qnpPrefix :

Returns the prefix part of the name, i.e., the part preceeding the ‘::’.

suffix_part((QualifiedName)arg1) qnpSuffix :

Returns the optional suffix part of the name, i.e., the part following the optional ‘@’.

class eos.SignalPDF
evaluate((_SignalPDF)self) float :

Evaluates the (unnormalized) PDF for the present values of the sets of parameters and kinematic variables that it is bound to.

Returns:

The value of the PDF.

Return type:

float

kinematics((_SignalPDF)arg1) Kinematics :

Returns the set of kinematic variables bound to this PDF.

log_pdf(x)

Adapter for use with external optimization software (e.g. pypmc) to aid when sampling from the log(PDF).

Parameters:

x (iterable of float) – Phase space point, with the elements corresponding to the sorted list of variables in lexicographical order; inspect self.variables.

static make(name, parameters, kinematics, options)

Makes a new SignalPDF object.

Parameters:
Returns:

The new PDF object.

Return type:

eos.SignalPDF

name((_SignalPDF)arg1) QualifiedName :

Returns the name of the PDF.

normalization((_SignalPDF)arg1) float :

Evaluates the normalization of the PDF.

To speed up sampling from the PDF, the evaluate returns values of an unnormalized function proportional to the actual PDF. To ensure that the integral over the PDF is normalized to 1, the values returned by evaluate need to be divided by the return value of this method.

options((_SignalPDF)arg1) Options :

Returns the set of options used when creating the PDF.

parameters((_SignalPDF)arg1) _Parameters :

Returns the set of parameters bound to this PDF.

sample_mcmc(N, stride, pre_N, preruns, cov_scale=0.1, start_point=None, rng=<module 'numpy.random.mtrand' from '/opt/venv/lib/python3.10/site-packages/numpy/random/mtrand.cpython-310-x86_64-linux-gnu.so'>)

Return samples of the kinematic variables and the log(PDF).

Obtains random samples of the log(PDF) using an adaptive Markov Chain Monte Carlo with PyPMC. A prerun with adaptations is carried out first and its samples are discarded.

Parameters:
  • N – Number of samples that shall be returned

  • stride – Stride, i.e., the number by which the actual amount of samples shall be thinned to return N samples.

  • pre_N – Number of samples in each prerun.

  • preruns – Number of preruns.

  • cov_scale – Scale factor for the initial guess of the covariance matrix.

  • start_point (list-like, optional) – Optional starting point for the chain

  • rng – Optional random number generator (must be compatible with the requirements of pypmc.sampler.markov_chain.MarkovChain)

Returns:

A tuple of the kinematic variables as array of size N and the log(PDF) as array of size N.

Note

This method requires the PyPMC python module, which can be installed from PyPI.

class eos.SignalPDFs(prefix=None, name=None, suffix=None, showall=False)

Represents the complete list of probability density functions (PDFs) known to EOS.

Objects of this class are visualized as tables in Jupyter notebooks for easy overview. Filters can be applied through keyword arguments to the initialization.

Parameters:
  • prefix (str) – Only show observables whose qualified names contain the provided prefix in their prefix part.

  • name (str) – Only show observables whose qualified names contain the provided name in their name part.

  • suffix (str) – Only show observables whose qualified names contain the provided suffix in their suffix part.

See also the complete list of signal PDFs in this documentation.

sections((object)arg1) object
class eos.Unit

Represents the unit of the observables.

Seven possible entries are currently implemented in EOS: - Undefined - None - GeV - GeV2 - InverseGeV2 - InverseGeV4 - InversePicoSecond

static GeV() Unit
static GeV2() Unit
static InverseGeV2() Unit
static InverseGeV4() Unit
static InversePicoSecond() Unit
static Undefined() Unit
static Unity() Unit
latex((Unit)arg1) str

Module eos.tasks

Common Tasks

In addition to the basic classes, EOS provides functions to carry out often-repeated tasks, such as finding the mode of a posterior, sampling from a posterior, and similar. Tasks require a description of the statistical analysis (see eos.Analysis) by means of an analysis file. For a documentation of its format, see eos.AnalysisFile.

Tasks store their results in a hierarchy for directories below a user-provided base directory. This ensures that tasks can readily use results produced by another task.

eos.corner_plot(analysis_file: str, posterior: str, base_directory: str = './', format: str = 'pdf', distribution: str = 'posterior', begin: int = 0, end: int = None)

Generates a corner plot of the 1-D and 2-D marginalized posteriors.

The input files are expected in EOS_BASE_DIRECTORY/POSTERIOR/samples. The output files will be stored in EOS_BASE_DIRECTORY/POSTERIOR/plots.

Parameters:
  • analysis_file (str or eos.AnalysisFile) – The name of the analysis file that describes the named posterior, or an object of class eos.AnalysisFile.

  • posterior (str) – The name of the posterior.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • format (str or list of str, optional) – The file extension of the data files. Can also be a list of file extensions.

  • distribution (str, optional) – The distribution to plot. Can be either ‘posterior’ or the name of a prediction.

  • begin (int, optional) – The index of the first parameter to plot. Defaults to 0.

  • end (int or NoneType, optional) – The index beyond the last parameter to plot. Defaults to None.

eos.find_clusters(posterior: str, base_directory: str = './', threshold: float = 2.0, K_g: int = 1, analysis_file: str = None)

Finds clusters among posterior MCMC samples, grouped by Gelman-Rubin R value, and creates a Gaussian mixture density.

Finding clusters and creating a Gaussian mixture density is a necessary intermediate step before using the sample-pmc subcommand. The input files are expected in EOS_BASE_DIRECTORY/POSTERIOR/mcmc-*. All MCMC input files present will be used in the clustering. The output files will be stored in EOS_BASE_DIRECTORY/POSTERIOR/clusters.

Parameters:
  • posterior (str) – The name of the posterior.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • threshold (float > 1.0, optional) – The R value threshold. If two sample subsets have an R value larger than this threshold, they will be treated as two distinct clusters. Defaults to 2.0.

  • K_g (int >= 1, optional) – The number of mixture components per cluster. Default to 1.

eos.find_mode(analysis_file: str, posterior: str, base_directory: str = './', optimizations: int = 3, start_point: list = None, chain: int = None, importance_samples: bool = None, seed: int = None, label: str = 'default')

Finds the mode of the named posterior.

The optimization process can be initialized either with a random point, a provided parameter point, or by extracting the point with the largest posterior from among previously obtained MCMC or importance samples. The optimization can be iterated to increase the accuracy of the result.

Parameters:
  • analysis_file (str or eos.AnalysisFile) – The name of the analysis file that describes the named posterior, or an object of class eos.AnalysisFile.

  • posterior (str) – The name of the posterior PDF from which to draw the samples.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • optimizations (int, optional) – The number of calls to the optimization algorithm.

  • start_point (numpy.ndarray, optional) – If provided, the optimization will start at this point.

  • chain (int, optional) – If provided, the optimization will start at the point with the largest posterior of previously computed MCMC samples. The samples are expected to be stored in the mcmc-XXXX subdirectories of the base_directory, where ‘XXXX’ is the provided int.

  • importance_samples (bool, optional) – If set to True, the optimization will start at the point with the largest posterior of previously computed importance samples.

  • seed (int, optional) – If provided, the optimization will start from a random point and the random number generator will be seeded with this value.

eos.mixture_product(posterior: str, posteriors: list, base_directory: str = './', analysis_file: str = None)

Compute the cartesian product of the densities listed in posteriors. Note that this product is not commutative.

The input densities are read from EOS_BASE_DIRECTORY/POSTERIOR_i/pmc, where POSTERIOR_i is listed in posteriors. The output density will be stored in EOS_BASE_DIRECTORY/POSTERIOR/product.

Parameters:
  • posterior (str) – The name of the posterior.

  • posteriors (iterable of str) – The list of names of the posteriors whose mixture densities will be concatenated.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

eos.predict_observables(analysis_file: str, posterior: str, prediction: str, base_directory: str = './', begin: int = 0, end: int = None)

Predicts a set of observables based on previously obtained importance samples.

The input files are expected in EOS_BASE_DIRECTORY/POSTERIOR/samples. The output files will be stored in EOS_BASE_DIRECTORY/POSTERIOR/pred-PREDICTION.

Parameters:
  • analysis_file (str or eos.AnalysisFile) – The name of the analysis file that describes the named posterior, or an object of class eos.AnalysisFile.

  • posterior (str) – The name of the posterior.

  • prediction (str) – The name of the set of observables to predict.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • begin (int) – The index of the first sample to use for the predictions. Defaults to 0.

  • end – The index beyond the last sample to use for the predictions. Defaults to None.

eos.sample_mcmc(analysis_file: str, posterior: str, chain: int, base_directory: str = './', pre_N: int = 150, preruns: int = 3, N: int = 1000, stride: int = 5, cov_scale: float = 0.1, start_point: list = None)

Samples from a named posterior PDF using Markov Chain Monte Carlo (MCMC) methods.

The output file will be stored in EOS_BASE_DIRECTORY/POSTERIOR/mcmc-CHAIN.

Parameters:
  • analysis_file (str or eos.AnalysisFile) – The name of the analysis file that describes the named posterior, or an object of class eos.AnalysisFile.

  • posterior (str) – The name of the posterior PDF from which to draw the samples.

  • chain (int >= 0) – The index assigned to the Markov chain. This value is used to seed the RNG for a reproducible analysis.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • pre_N (int, optional) – The number of samples to be used for an adaptation in each prerun steps. These samples will be discarded.

  • preruns (int, optional) – The number of prerun steps, which are used to adapt the MCMC proposal to the posterior.

  • N (int, optional) – The number of samples to be stored in the output file. Defaults to 1000.

  • stride (int, optional) – The ratio of samples drawn over samples stored. For every S samples, S - 1 will be discarded. Defaults to 5.

  • cov_scale (float, optional) – Scale factor for the initial guess of the covariance matrix.

  • start_point (list-like, optional) – Optional starting point for the chain

eos.sample_nested(analysis_file: str, posterior: str, base_directory: str = './', bound: str = 'multi', nlive: int = 250, dlogz: float = 1.0, maxiter: int = None, seed: int = 10)

Samples from a likelihood associated with a named posterior using dynamic nested sampling.

The output will be stored in EOS_BASE_DIRECTORY/POSTERIOR/nested. In addition, an ImportanceSamples object is exported to EOS_BASE_DIRECTORY/POSTERIOR/samples.

Parameters:
  • analysis_file (str or eos.AnalysisFile) – The name of the analysis file that describes the named posterior, or an object of class eos.AnalysisFile.

  • posterior (str) – The name of the posterior.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • bound (str, optional) – The option for bounding the target distribution. For valid values, see the dynesty documentation. Defaults to ‘multi’.

  • nlive (int, optional) – The number of live points.

  • dlogz (float, optional) – Relative tolerance for the remaining evidence. Defaults to 5%.

  • maxiter (int, optional) – The maximum number of iterations. Iterations may stop earlier if the termination condition is reached.

  • seed (int, optional) – The seed used to initialize the Mersenne Twister pseudo-random number generator.

eos.sample_pmc(analysis_file: str, posterior: str, base_directory: str = './', step_N: int = 500, steps: int = 10, final_N: int = 5000, perplexity_threshold: float = 1.0, weight_threshold: float = 1e-10, sigma_test_stat: list = None, initial_proposal: str = 'clusters', pmc_iterations: int = 1, pmc_rel_tol: float = 1e-10, pmc_abs_tol: float = 1e-05, pmc_lookback: int = 1)

Samples from a named posterior using the Population Monte Carlo (PMC) methods.

The results of the find-cluster command are expected in EOS_BASE_DIRECTORY/POSTERIOR/clusters. The output will be stored in EOS_BASE_DIRECTORY/POSTERIOR/pmc. In addition, an ImportanceSamples object is exported to EOS_BASE_DIRECTORY/POSTERIOR/samples.

Parameters:
  • analysis_file (str or eos.AnalysisFile) – The name of the analysis file that describes the named posterior, or an object of class eos.AnalysisFile.

  • posterior (str) – The name of the posterior.

  • base_directory (str, optional) – The base directory for the storage of data files. Can also be set via the EOS_BASE_DIRECTORY environment variable.

  • step_N (int > 0, optional) – The number of samples to be used in each adaptation step. These samples will be discarded. Defaults to 500.

  • steps (int > 0, optional) – The number of adaptation steps, which are used to adapt the PMC proposal to the posterior. Defaults to 10.

  • final_N (int > 0, optional) – The number of samples to be stored in the output file. Defaults to 5000,

  • perplexity_threshold (0.0 < float <= 1.0, optional) – The threshold for the perplexity in the last step after which further adaptation steps are to be skipped. Defaults to 1.0.

  • weight_threshold (0.0 < float <= 1.0, optional.) – Mixture components with a weight smaller than this threshold are pruned.

  • sigma_test_stat (list or iterable) – If provided, the inverse CDF of -2*log(PDF) will be evaluated, using the provided values as the respective significance.

  • initial_proposal (str, optional) –

    Specify where the initial proposal should be taken from:

    • clusters: use the proposal obtained using find-clusters (default)

    • product: use the proposal obtained from mixture_product

    • pmc: continue sampling from the previous sample-pmc results.

  • pmc_iterations (int > 0, optional, advanced) – Maximum number of update of the PMC, changing this value may make the update unstable.

  • pmc_rel_tol (float > 0.0, optional, advanced) – Relative tolerance of the PMC. If two consecutive values of the current density log-likelihood are relatively smaller than this value, the convergence is declared.

  • pmc_abs_tol (float > 0.0, optional, advanced) – Absolute tolerance of the PMC. If two consecutive values of the current density log-likelihood are smaller than this value, the convergence is declared.

  • pmc_lookback (int >= 0, optional) – Use reweighted samples from the previous update steps when adjusting the mixture density. The parameter determines the number of update steps to “look back”. The default value of 1 disables this feature, a value of 0 means that all previous steps are used.

Module eos.data

EOS provides access to save and load the various (intermediate) results of analyses via the eos.data module.

class eos.data.ImportanceSamples(path)
static create(path, parameters, samples, weights, posterior_values=None)

Write a new ImportanceSamples object to disk.

Parameters:
  • path (str) – Path to the storage location, which will be created as a directory.

  • parameters (list or iterable of eos.Parameter) – Parameter descriptions as a 1D array of shape (P, ).

  • samples (2D numpy array) – Samples as a 2D array of shape (N, P).

  • weights (1D numpy array, optional) – Weights on a linear scale as a 2D array of shape (N, 1).

class eos.data.MarkovChain(path)
static create(path, parameters, samples, usamples, weights=None)

Write a new MarkovChain object to disk.

Parameters:
  • path (str) – Path to the storage location, which will be created as a directory.

  • parameters (list or iterable of eos.Parameter) – Parameter descriptions as a 1D array of shape (N, ).

  • samples (2D numpy array) – Samples in parameter space as a 2D array of shape (N, P).

  • usamples (2D numpy array) – Samples in u space as a 2D array of shape (N, P).

  • weights (2D numpy array, optional) – Weights on a linear scale as a 2D array of shape (N, 1).

class eos.data.MixtureDensity(path)
static cartesian_product(densities)

Construct the cartesian product of a list of mixture densities. The means and covariances of the components are copied and concatenated and the weights are multiplied. Note that this product is not commutative, the order of densities affects the order of parameters in the output density. Only Gaussian mixtures are supported.

Parameters:

densities (iterable of pypmc.density.MixtureDensity) – List of mixture densities.

static create(path, density)

Write a new MixtureDensity object to disk.

Parameters:
  • path (str) – Path to the storage location, which will be created as a directory.

  • density (pypmc.density.MixtureDensity) – Mixture density.

density()

Return a pypmc.density.MixtureDensity.

class eos.data.Mode(path)
static create(path, parameters, mode, pvalue, local_pvalues)

Write a new Mode object to disk.

Parameters:
  • path (str) – Path to the storage location, which will be created as a directory.

  • parameters (list or iterable of eos.Parameter) – Parameter descriptions as a 1D array of shape (N, ).

  • mode (numpy.ndarray) – The mode to be stored.

class eos.data.PMCSampler(path)
static create(path, parameters, proposal, sigma_test_stat=None, samples=None, weights=None)

Write a new PMCSampler object to disk.

Parameters:
  • path (str) – Path to the storage location, which will be created as a directory.

  • parameters (list or iterable of eos.Parameter) – Parameter descriptions as a 1D array of shape (P, ).

  • sigma_test_stat (list or iterable) – (optional) If provided, the inverse CDF of -2*log(PDF) will be evaluated, using the provided values as the respective significance.

  • samples (2D numpy array, optional) – Samples as a 2D array of shape (N, P). Needed to generate the test statistic.

  • weights (1D numpy array, optional) – Weights on a linear scale as a 2D array of shape (N, 1). Needed to generate the test statistic.

density()

Return a pypmc.density.MixtureDensity.

class eos.data.Prediction(path)
static create(path, observables, samples, weights)

Write a new Prediction object to disk.

Parameters:
  • path (str) – Path to the storage location, which will be created as a directory.

  • observables (list or iterable of eos.Observable) – Observables as a 1D array of shape (O, ).

  • samples (2D numpy array) – Samples as a 2D array of shape (N, O).

  • weights (1D numpy array) – Weights on a linear scale as a 1D array of shape (N, ).

Module eos.plot

EOS provides a plotting framework based on Matplotlib. Plots can readily be created from a Python script, from within a Jupyter notebook, or in the command line using the eos-plot script. For all of these cases a description of the plot is required in the format described below. For the command-line script eos-plot, the Python dictionary describing the plots must be provided as a YAML file.

Note

Import eos.plot before you do something like import matplotlib.pyplot as plt, because the eos.plot module sets its default plot style and a matplotlib backend. All options (except the backend) can be overwritten by updating matplotlib.rcParams[...]; see also the matplotlib documentation. Note that the default settings use LaTeX to create labels and math expressions, so for this to work latex needs to be available on your system.

class eos.plot.Plotter(description, output=None)

Produces publication-quality plots

Plots can contain EOS observables, EOS constraints, and Analysis results. See Plot description format for documentation of how to create a plot.

Parameters:
  • description (dict) – Description of the plot and its contents, see Plot description format.

  • output (string, optional) – Name of the output file. The file format is automatically determined based on the file’s extension.

class Band(plotter, item)

Plots a shaded band

class BasePlot(plotter, item)

Base class for any of the plots supported by Plotter

class Constraint(plotter, item)

Plots constraints from the EOS library of experimental and theoretical likelihoods

class Constraint2D(plotter, item)

Plot 2D contours for two correlated observables from a single constraint

class ConstraintOverview(plotter, item)

Plots overview of several constraints from the EOS library of experimental and theoretical likelihoods

class Contours2D(plotter, item)

Plots 2D contours of a pair of parameters based on pre-existing random samples

class ErrorBar(plotter, item)

Plots a single errors bar

class Expression(plotter, item)

Plots a given expression

class ExternalLikelihood2D(plotter, item)

Plots contours of a user-provided function

class Graph(plotter, item)

Plots the graph of a function

class Histogram1D(plotter, item)

Plots a 1D histogram of pre-existing random samples

class Histogram2D(plotter, item)

Plots a 2D histogram of pre-existing random samples

class KernelDensityEstimate1D(plotter, item)

Plots a 1D Kernel Density Estimate (KDE) of pre-existing random samples

class KernelDensityEstimate2D(plotter, item)

Plots contours of a 2D Kernel Density Estimate (KDE) of pre-existing random samples

class Observable(plotter, item)

Plots a single EOS observable w/o uncertainties as a function of one kinematic variable or one parameter

class Point(plotter, item)

Plots a single point

class SignalPDF(plotter, item)

Plots a single EOS signal PDF w/o uncertainties as a function of one kinemtic variable

class Uncertainty(plotter, item)

Plots an uncertainty band as a function of one kinematic variable

This routine expects the uncertainty propagation to have produced an EOS data file

class UncertaintyBinned(plotter, item)

Plots one or more uncertainty band integrated over one kinematic variable

This routine expects the uncertainty propagation to have produced an data file

class UncertaintyOverview(plotter, item)

Plots an overview of uncertainty estimates

This routine expects the uncertainty propagation to have produced an HDF5 file

class Watermark(plotter, item)

Inserts an EOS watermark into the plots

property next_color

Returns the next available color

property next_z_order

Returns the next available z-order value, incremented for each plot w/o pre-defined z-order value

plot()

Produces the plot

Returns:

  • fig (matplotlib.figure.Figure) - the created figure

  • ax (matplotlib.axes.Axes) - the created axes

plot_contents()

Plots the contents specified in the instructions provided to Plotter

setup_plot()

Setting up the plot based on the provided instruction

Plot description format

The input must be formatted as a dictionary containing the keys plot and contents. The plot key must be mapped to a dictionary; it describes the layout of the plot, including axis labels, positioning of the key, and similar settings. The contents key must be mapped to a list; it describes the contents of the plot, expressed in terms of independent plot items.

plot_desc = {
    'plot': {
        'x': { ... },       # description of the x axis
        'y': { ... },       # description of the y axis
        'legend': { ... },  # description of the legend
        ...                 # further layouting options
    },
    'contents': [
        { ... }, # first plot item
        { ... }, # second plot item
    ]
}
eos.plot.Plotter(plot_desc, FILENAME).plot()

In the above, FILENAME is an optional argument naming the file into which the plot shall be placed. The format is automatically determined based on the file name extension.

Plot Layouting

By default plots lack any axis labels and units, and any legend.

An axis’ description is provided through the following key/value pairs, which can apply equally to the x and y axis:

  • label (str, may contain LaTeX commands) – The axis’ label.

  • unit (str, may contain LaTeX commands) – The axis’ unit, which will be appended to the axis’ label in square brackets.

  • range (list or tuple of two float) – The tuple of [minimal, maximal] values, which will be displayed along the axis.

  • scale (str) – Can be either linear or log.

  • scaling_factor (number) – The axis’ scale by which all tick coordinates will be divided. The scale will not be appended to the axis’ label automatically. Using this argument is not recommended as it prevents automatic axis tick formatting and providing the argument below is required.

  • format (str, Python 3 format string) – The axis’ tick label format can be provided and is only in use when scale is used, to avoid a bad string representation of the axis ticks. For example, the user might need to determine the necessary number of digits manually. Due to a Matplotlib peculiarity, the format string must always format the variable x. See Matplotlib format strings for details.

The legend description presently only includes options for its location:

  • location (str, valid Matplotlib legend location) – The legend’s location within the plot.

Further layouting options are:

  • title (str) – The plot’s title.

  • size (tuple of two numbers) – The plot’s size in x and y directions provided in centimeters.

  • axes (str, equal) – Enforces equal scaling of the plot’s x and y axes, if set.

  • grid (str, either major, minor, or both) – Enables the plot’s gridline, if set.

An example illustrating plot layouting follows:

plot_args = {
    'plot': {
        'x': { 'label': r'$q^2$', 'unit': r'$\textnormal{GeV}^2$', 'range': [0.0, 11.60] },
        'y': { 'label': r'$d\mathcal{B}/dq^2$',                    'range': [0.0,  5e-3] },
        'legend': { 'location': 'upper center' },
        'size': [10, 5]
    },
    'contents': [
        ...
    ]
}

Plot Contents

Each item in a plot’s contents is represented by a dictionary. Each features the mandatory type key and more type-specific (mandatory or optional) keys.

  • type (str, mandatory) – The type of the plot item, from one of the following recognized item types:

    value

    description

    band

    Plots a shaded band

    constraint

    Plots constraints from the EOS library of experimental and theoretical likelihoods

    constraint2D

    Plot 2D contours for two correlated observables from a single constraint

    constraint-overview

    Plots overview of several constraints from the EOS library of experimental and theoretical likelihoods

    graph

    Plots the graph of a function

    expression

    Plots a given expression

    errorbar

    Plots a single errors bar

    histogram

    Plots a 1D histogram of pre-existing random samples

    histogram2D

    Plots a 2D histogram of pre-existing random samples

    kde

    Plots a 1D Kernel Density Estimate (KDE) of pre-existing random samples

    kde2D

    Plots contours of a 2D Kernel Density Estimate (KDE) of pre-existing random samples

    likelihood2D

    Plots contours of a user-provided function

    observable

    Plots a single EOS observable w/o uncertainties as a function of one kinematic variable or one parameter

    point

    Plots a single point

    signal-pdf

    Plots a single EOS signal PDF w/o uncertainties as a function of one kinemtic variable

    uncertainty

    Plots an uncertainty band as a function of one kinematic variable

    uncertainty-binned

    Plots one or more uncertainty band integrated over one kinematic variable

    watermark

    Inserts an EOS watermark into the plots

All item types recognize the following optional keys:

  • name (str, optional) – The name of the plot item, for convenience when reporting warnings and errors.

  • alpha (float, between 0.0 and 1.0) – The opacity of the plot item expressed as an alpha value. 0.0 means completely transparent, 1.0 means completely opaque.

  • color (str, containing any valid Matplotlib color specification) – The color of the plot item. Defaults to one of the colors in the Matplotlib default color cycler.

  • label (str, may contain LaTeX commands) – The label that appears in the plot’s legend for this plot item.

  • style (str, containing any valid Matplotlib style specification) – The style of the plot item.

Plotting Constraints

Contents items of type constraints are used to display one of the built-in experimental or theoretical constraints. The following keys are mandatory:

  • constraints (QualifiedName or iterable thereof) – The name or the list of names of the constraints that will be plotted. Must identify at least one of the constraints known to EOS; see the complete list of constraints.

  • variable (str) – The name of the kinematic variable to which the x axis will be mapped.

When plotting multivariate constraints, the following key is also mandatory:

  • observable (QualifiedName) – The name of the observable whose constraints will be plotted. Must identify one of the observables known to EOS; see the complete list of observables. This is only mandatory in multivariate constraints, since these can constrain more than one observable simultaneously.

The following keys are optional:

  • xrange (list of int) – The interval in which the observable is plotted in the case of a multivariate constraint.

  • rescale-by-width (bool) – Rescales binned constraints by the inverse of the bin width. This is often required to compare theory (integrated) predictions and experimental (averaged) measurements. Defaults to false.

Example:

plot_args = {
    'plot': { ... },
    'contents': [
        {
            'label': r'Belle 2015 $\ell=e,\, q=d$',
            'type': 'constraint',
            'color': 'C0',
            'constraints': 'B^0->D^+e^-nu::BRs@Belle:2015A',
            'observable': 'B->Dlnu::BR',
            'variable': 'q2',
            'rescale-by-width': False
        }
    ]
}

Plotting a Function Graph

The graph of a function can be easily plotted by providing the coordinates (x, f(x)) of the function. The coordinates are provided using a data object containing:

  • xvalues (array-like of float) – The values on the x axis at which the function has been evaluated.

  • yvalues (array-like of float) – The values oof the function at the points provided by xvalues.

Example:

xvalues = np.linspace(0, 5, 20)
plot_args = {
    'plot': {
        'x': { 'label': r'$q^2$' },
        'y': { 'label': r'$q^4$' }
    },
    'contents': [
        {
            'label': r'$\ell=\mu$',
            'type': 'graph',
            'data': { 'yvalues': xvalues**2, 'xvalues': xvalues },
            'range': [0, 5],
        },
    ]
}

Plotting Histograms

Contents items of type histogram are used to display samples of a probability density, be it a prior, a posterior, or a signal PDF. The following key is mandatory:

  • data (dict, see below) – The data on probability density that will be histogramed.

Within the data object, the following keys are understood.

  • samples (list of float) – The samples that will be histogramed. Mandatory.

  • weights or log_weights (list of float, optional) – The weights of the samples, on a linear or logarithmic scale. Defaults to uniform weights.

Example:

analysis = ... # eos.Analysis object as discussed in the example notebook `inference.ipynb`
parameter_samples, _, = analysis.sample(N=5000, pre_N=1000)
plot_args = {
    'plot': {
        'x': { 'label': r'$|V_{cb}|$' },
        'y': { 'label': r'$d\mathcal{B}/dq^2$' }
    },
    'contents': [
        {
            'type': 'histogram',
            'data': { 'samples': parameter_samples[:, 0] },
        },
    ]
}

Plotting 2D Histograms

Contents items of type histogram2D are used to display samples of a two-dimension probability density, be it a prior, a posterior, or a signal PDF. The following key is mandatory:

  • data (dict, see below) – The data on probability density that will be histogramed.

Within the data object, the following keys are understood.

  • samples (list of float with shape (N, 2)) – The samples that will be histogramed. Mandatory.

  • weights or log_weights (list of float with length N, optional) – The weights of the samples, on a linear or logarithmic scale, respectively. Defaults to uniform weights.

Additional optional keys are:

  • bins (int or list of int, optional) – The number of bins for both dimensions together or in each dimension seperately. Defaults to 100 for both dimensions.

Example:

analysis = ... # eos.Analysis object as discussed in the example notebook `inference.ipynb`
dstarlnu_kinematics = ... # create kineamtics and SignalPDF as discussed in the example notebook `simulation.ipynb`
dstarlnu_pdf        = ...
dstarlnu_samples, _ = dstarlnu_pdf.sample_mcmc(N=50000, stride=5, pre_N=1000, preruns=3, rng=rng)
plot_args = {
    'plot': {
        'x': { 'label': r'$q^2$', 'unit': r'$        extnormal{GeV}^2$', 'range': [ 0.0, 10.50] },
        'y': { 'label': r'$cos(      heta_\ell)$',                     'range': [-1.0,  +1.0] },
    },
    'contents': [
        {
            'label': r'samples ($\ell=\mu$)',
            'type': 'histogram2D',
            'data': {
                'samples': dstarlnu_samples[:, (0, 1)]
            },
            'bins': 40
        },
    ]
}
eos.plot.Plotter(plot_args).plot()

Plotting Kernel Density Estimates

Contents items of type kde are used to display a kernel density estimate (a smooth histogram) of samples of a probability density, be it a prior, a posterior, or a signal PDF.

The following key is mandatory:

  • data (dict, see below) – The data on the probability density that will be histogramed.

    Within the data object, the following keys are understood.

    • samples (list of float) – The samples that will be histogramed. Mandatory.

    • weights or log_weights (list of float, optional) – The weights of the samples, on a linear or logarithmic scale. Defaults to uniform weights.

The following keys are optional:

  • bandwidth (float) – The factor by which the automatically determined kernel bandwidth is scaled. See the SciPy documentation for gaussian_kde, bw_method='silverman'. Defaults to 1.

  • range (tuple of two float) – The minimum and maximum value of the x coordinate for which the smooth histogram is plotted.

  • level (float) – The probability level of the contour. Defaults to 68.27.

Example:

analysis = ... # eos.Analysis object as discussed in the example notebook `inference.ipynb`
parameter_samples, _, = analysis.sample(N=5000, pre_N=1000)
plot_args = {
    'plot': {
        'x': { 'label': r'$|V_{cb}|$', 'range':[38e-3, 45e-3] }
    },
    'contents': [
        {
            'type': 'kde', 'color': C0, 'label': 'posterior', 'bandwidth': 1,
            'range': [40e-3, 45e-3],
            'data': { 'samples': parameter_samples[:, 0] }
        }
    ]
}

Plotting 2D Kernel Density Estimates

Contents items of type kde2D are used to display contours of a two-dimensional kernel density estimate (a 2D smooth histogram) of samples of a probability density, be it a prior, a posterior, or a signal PDF.

The following key is mandatory:

  • data (dict, see below) – The data on the probability density that will be histogramed.

    Within the data object, the following keys are understood.

    • samples (list of float with shape (N, 2)) – The samples that will be histogramed. Mandatory.

    • weights or log_weights (list of float, optional) – The weights of the samples, on a linear or logarithmic scale. Defaults to uniform weights.

The following keys are optional:

  • bandwidth (float) – The factor by which the automatically determined kernel bandwidth is scaled. See the SciPy documentation for gaussian_kde, bw_method='silverman'. Defaults to 1.

  • contours (a list containing 'lines', 'areas', or both) – The setting for the illustration of the contours. If 'lines' is provided, the contour lines are drawn. If 'areas' is provided, the contour areas are filled. Defaults to ['lines'].

  • levels (list of float) – The probability levels of the contours. Defaults to [68, 95, 99].

Example:

analysis = ... # eos.Analysis object as discussed in the example notebook `inference.ipynb`
parameter_samples, _, = analysis.sample(N=5000, pre_N=1000)
plot_args = {
    'plot': {
        'x': { 'label': r'$|V_{cb}|$', 'range': [38e-3, 47e-3] },
        'y': { 'label': r'$f_+(0)$',   'range': [0.6, 0.75] },
    },
    'contents': [
        {
            'type': 'kde2D', 'color': 'C1', 'label': 'posterior',
            'levels': [68, 95], 'contours': ['lines','areas'], 'bandwidth':3,
            'data': { 'samples': parameter_samples[:, (0,1)] }
        }
    ]
}
eos.plot.Plotter(plot_args).plot()

Plotting a user-provided Likelihood

The following key is mandatory:

  • likelihood (evaluable, see below) – The 2D probability density that will be plotted.

The following keys are optional:

  • contours (a list containing 'lines', 'areas', or both) – The setting for the illustration of the contours. If 'lines' is provided, the contour lines are drawn. If 'areas' is provided, the contour areas are filled. Defaults to ['lines'].

  • levels (list of float) – The probability levels of the contours. Defaults to [68, 95, 99].

  • xrange (list of float) – The x-axis range where the likelihood is evaluated, defaults to the plot ranges if provided.

  • yrange (list of float) – The y-axis range where the likelihood is evaluated, defaults to the plot ranges if provided.

Plotting Observables

Contents items of type observable are used to display one of the built-in observables. The following keys are mandatory:

  • observable (QualifiedName) – The name of the observable that will be plotted. Must identify one of the observables known to EOS; see the complete list of observables.

  • range (list or tuple of two float) –The tuple of [minimal, maximal] values of the specified kinematic variable for which the observable will be evaluated.

  • variable (str) – The name of the kinematic or parameter variable to which the x axis will be mapped; see the complete list of parameters.

Example:

plot_args = {
    'plot': { ... },
    'contents': [
        {
            'label': r'$\ell=\mu$',
            'type': 'observable',
            'observable': 'B->Dlnu::dBR/dq2;l=mu',
            'variable': 'q2',
            'range': [0.02, 11.60],
        },
    ]
}

Plotting Points

Content items of type point are used to display a single data point manually. The following keys are mandatory:

  • x (number) – The point’s x coordinate.

  • y (number) – The point’s y coordinate.

Beside the common set of optional keys, this item type recognizes the following optional keys:

  • marker (str, a valid Matplotlib marker style) – The point’s marker style.

  • markersize (number, a valid Matplotlib marker size) – The point’s marker size in pts.

Example:

plot_args = {
    'plot': { ... },
    'contents': [
        {
            'label': r'LCSR (Bharucha 2012)',
            'type': 'point',
            'color': 'C0',
            'x': 0,
            'y': 0.261,
            'marker': 'o',
            'markersize': 12
        }
    ]
}

Plotting Uncertainty Bands

Contents items of type uncertainty are used to display uncertainty bands at 68% probability for one of the built-in observables. This item type requires that either a predictive distribution of the observables has been previously produced.

Exactly one of the following keys is mandatory:

  • data (dict, see below) – The data on predictive distribution of the observable whose uncertainty band will be plotted.

  • data-file (str, path to an existing data file of type eos.data.Prediction) – The path to a data file that was generated with the eos-analysis command-line client.

For data object, the following keys are mandatory:

  • xvalues (array-like of float) – The values on the x axis at which the observable has been evaluated.

  • samples (list of tuples of float) – The list of samples of the predictive distribution. Each tuple of samples corresponds to an evaluation of the observables at the kinematic configuration corresponding to the xvalues entry with the same index.

  • weights (array-like of float, optional) – The weights of the samples, on a linear scale. Defaults to uniform weights.

The following keys are optional:

  • band (a list containing 'lines', 'areas', or both) – The setting for the illustration of the band. If 'outer' is provided, the band’s outer lines are drawn. If 'median' is provided, the band’s median line is drawn. If 'area' is provided, the band’s areas are filled. Defaults to ['area', 'outer', 'median'].

Example:

analysis = ... # eos.Analysis object as discussed in the example notebook `predictions.ipynb`
mu_q2values = numpy.unique(numpy.concatenate((numpy.linspace(0.02,1.00,20), numpy.linspace(1.00,11.60,20))))
mu_obs      = [eos.Observable.make('B->Dlnu::dBR/dq2', prior.parameters,
                                   eos.Kinematics(q2=q2), eos.Options({'form-factors':'BSZ2015','l':'mu'}))
               for q2 in mu_q2values]
_, _, mu_samples = analysis.sample(N=5000, pre_N=1000, observables=mu_obs)
plot_args = {
    'plot': {
        'x': { 'label': r'$q^2$' },
        'y': { 'label': r'$d\mathcal{B}/dq^2$' }
    },
    'contents': [
        {
            'label': r'$\ell=\mu$',
            'type': 'uncertainty',
            'data': { 'samples': mu_samples, 'xvalues': mu_q2values },
            'range': [0.02, 11.60],
        },
    ]
}

Module eos.figure

EOS also provides a new plotting framework with a different interface, which will replace the existing framework as soon as it is feature complete.

Note

Import eos.figure before you do something like import matplotlib.pyplot as plt, because the eos.figure module modified the default plot style and the matplotlib backend. All options (except the backend) can be overwritten by updating matplotlib.rcParams[...]; see also the matplotlib documentation. Note that the default settings use LaTeX to create labels and math expressions, so for this to work latex needs to be available on your system.

Figures

The eos.figure submodule provides means to produce publication-quality figures.

A figure contains one or more plots, each of which is a set of axis that can be thought of as a ‘canvas’ to draw on. A plot further contains items for display, for example a constraint in the form of data points with error bars or an observable as a function of a kinematic variable in the form or a line graph.

eos.figure.make(description: dict) Figure

This function takes a figure description and handles the contruction of a publication-grade figure in a convenient way. The description contains the arguments used for initialization of the objects in the eos.figure module.

A figure can contain one or more plots, where each plot containts one or more items such as a data series from the evaluation of an observable or a kernel density estimate to visualize a set of statistical samples.

The description can be either a dictionary (figure description) describing a single plot or a list of dictionaries describing a grid of plots.

figure description (dict)

description

plot description (dict)

items

list of item description

plot description (dict)

plot_type

valid plot type (str)

one or several keyword arguments

item description (dict)

item_type

valid item type (str)

one or several keyword arguments

Parameters:

description (dict or list of dict) – A figure description as documented

Returns:

EOS figure

Return type:

eos.figure.Figure

Plot types

Below we list the available plot types. Each key plot_type corresponds to a class which is instantiated in the background. The user can supply keyword arguments which are documented under eos.figure.plot.

plot_type

class

description

2D-plot

eos.figure.plot.TwoDimPlot

Represents a 2d plot with x and y axes.

empty

eos.figure.plot.EmptyPlot

Can be used in a grid as empty space instead of a plot.

The classes including the arguments the user can pass are documented in detail below.

class eos.figure.plot.EmptyPlot(ax, contents: list, **kwargs)

Can be used in a grid as empty space instead of a plot.

draw()

Draw nothing

class eos.figure.plot.Plot(ax, items: list, title=None, legend: dict | None = None)

Represents a plot along a single set of axes.

Parameters:
  • ax (matplotlib.pyplot.Axes or similar) – Set of axes

  • items (list[dict] or iterable[dict]) – List of content items displayed in this plot

  • title (str or NoneType (optional)) – Title of the plot

  • legend (dict) – Keyword arguments to pass to matplotlib.pyplot.legend

draw()

Draw the plot including all items

class eos.figure.plot.TwoDimPlot(ax, items: list, x_range: tuple[float, float] = (None, None), y_range: tuple[float, float] = (None, None), x_label: str | None = None, y_label: str | None = None, **kwargs)

Represents a 2d plot with x and y axes.

Parameters:
  • ax (matplotlib.pyplot.Axes or similar) – Set of axes

  • items (list[dict] or iterable[dict]) – List of content items displayed in this plot

  • x_range (tuple of float (optional)) – Range of the x axis (automatic if value is [None, None])

  • y_range (tuple of float (optional)) – Range of the y axis (automatic if value is [None, None])

  • x_label (str (optional)) – Label of the x axis

  • y_label (str (optional)) – Label of the y axis

  • **kwargs – Additional keyword arguments to pass to the base class eos.figure.plot.Plot

Item types

Below we list the available item types. Each key item_type corresponds to a class which is instantiated in the background. The user can supply keyword arguments which are documented under eos.figure.item.

item_type

class

description

observable

eos.figure.item.Observable

Represents an observable as a function of a variable.

The classes including the arguments the user can pass are documented in detail below.

class eos.figure.item.Observable(observable, variable: str, x_values: list, fixed_kinematics: dict | None = None, fixed_parameters: dict | None = None, fixed_parameters_from_file: str | None = None, options: dict | None = None, label: str | None = None, **kwargs)

Represents an observable as a function of a variable. The variable can be either a kinematic variable or a parameter.

Parameters:
  • observable (str) – EOS qualified name of an observable

  • variable (str) – Name of a kinematic variable or of a parameter

  • x_values (list of float numbers) – Values to be used as sampling points

  • fixed_kineamtics (dict) – Names and values of fixed kinematic variables

  • fixed_parameters (dict) – Names and values of fixed parameters

  • fixed_parameters_from_file (str) – Path to file that contains names and values of fixed parameters in the YAML format

  • options (dict) – Names and values of options passed to the EOS observable

  • label (str) – Label used in the legend

  • **kwargs – Additional keyword arguments to pass to matplotlib.axes.Axes.plot

draw(ax)

Draw a line plot of the observable

prepare()

Evaluate the observable at the sample points