Python API

Module eos

EOS provides its basic functionality via the main eos module.

class eos.Analysis(priors, likelihood, global_options={}, manual_constraints={}, fixed_parameters={})

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.

  • 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.

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 goodness-of-fit summary for the current parameter point.

log_pdf(x, *args)

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

Parameters
  • x (iterable) – Parameter point, with the elements in the same order as in eos.Analysis.varied_parameters, rescaled so that every element is in the interval [-1, +1].

  • args (optional) – Dummy parameter (ignored)

negative_log_pdf(x, *args)

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

Parameters
  • x (iterable) – Parameter point, with the elements in the same order as in eos.Analysis.varied_parameters, rescaled so that every element is in the interval [-1, +1].

  • args (optional) – Dummy parameter (ignored)

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

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

Parameters
  • start_point – 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.

  • start_point – iterable, optional

  • rng – Optional random number generator

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.6/site-packages/numpy/random/mtrand.cpython-36m-x86_64-linux-gnu.so'>)

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_pmc(log_proposal, step_N=1000, steps=10, final_N=5000, rng=<module 'numpy.random.mtrand' from '/opt/venv/lib/python3.6/site-packages/numpy/random/mtrand.cpython-36m-x86_64-linux-gnu.so'>, return_final_only=True, final_perplexity_threshold=1.0)

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 Popoulation 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.ImportancSampler)

  • 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 perpexlity of the last adaptation step is above this threshold value.

Returns

A tuple of the parameters as array of length N = pre_N * steps + final_N, the (linear) weights 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.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
evaluate((LogLikelihood)arg1) float
observable_cache((LogLikelihood)arg1) ObservableCache
class eos.LogPosterior
add((LogPosterior)arg1, (LogPrior)arg2, (bool)arg3) bool
evaluate((LogPosterior)arg1) float
log_likelihood((LogPosterior)arg1) LogLikelihood
log_priors((object)arg1) object
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 Flat((_Parameters)parameters, (str)name, (ParameterRange)range) LogPrior :

Alias for LogPrior.Uniform().

static Gauss((_Parameters)parameters, (str)name, (ParameterRange)range, (float)lower, (float)central, (float)upper) LogPrior :

Returns a new Gaussian prior as a LogPrior.

The prior’s support is provided by the range parameter, with the 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.

  • range (tuple of two floating point numbers) – The range [min, max] for the values 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 Scale((_Parameters)parameters, (str)name, (ParameterRange)range, (float)mu_0, (float)scale) LogPrior :

Returns a new Scale prior as a LogPrior.

The prior’s support is provided by the range parameter, 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.

  • range (tuple of two floating point numbers) – The range [min, max] for the values 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, (ParameterRange)range) 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.

  • range (tuple of two floating point numbers) – The range [min, max] for the values that the parameter is allowed to take.

inverse_cdf((LogPrior)arg1, (float)p) float :

Returns the parameter value corresponding to the cumulative propability \(p\).

Parameters

p (float, [0.0, 1.0]) – The cumulative propability.

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.

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.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 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.

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.

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
name_part((QualifiedName)arg1) qnpName :

Returns the name part of the name, i.e., the part following the ‘::’ and preceeding any 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
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

sample_mcmc(N, stride, pre_N, preruns, cov_scale=0.1, start_point=None, rng=<module 'numpy.random.mtrand' from '/opt/venv/lib/python3.6/site-packages/numpy/random/mtrand.cpython-36m-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.

Module eos.data

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

class eos.data.MarkovChain(path)
static create(path, parameters, samples, 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 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 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.PMCSampler(path)
static create(path, parameters, samples, weights, proposal, compute_test_stat=False, sigma_test=array([0.5, 1.0, 1.5, 2.0, 2.5, 3.0, 3.5, 4.0, 4.5, 5.0]))

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 (N, ).

  • 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).

  • compute_test_stat (bool) – (optional) If true, the weighted samples will be used to numerically estimate the inverse CDF of -2 * log(PDF) of the mixture.

  • sigma_test (list or iterable) – (optional) If compute_test_stat is true, the inverse CDF of -2*log(PDF) will be evaluated for the significance values contained in sigma_test.

density()

Return a pypmc.density.MixtureDensity.

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 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 kinemtic 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 HDF5 file

class UncertaintyBinned(plotter, item)

Plots one or more uncertainty band integrated over one kinematic variable

This routine expects the uncertainty propagation to have produces an HDF5 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

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 (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

    contours2D

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

    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

    observable

    Plots a single EOS observable w/o uncertainties as a function of one kinemtic 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

    uncertainty-overview

    Plots an overview of uncertainty estimates

    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 styale 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.

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 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, optional) – The number of bins in each dimension. Defaults to 100.

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.

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 [0.68, 0.95, 0.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 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.

Exactly one of the following keys is mandatory, to specify either a kinematic variable or a parameter to which the x coordinate will be mapped:

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

  • kinematic (str) – Alias for variable.

  • parameter (str) – The name of the parameter 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 (tuple 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.

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],
        },
    ]
}