pycbc.distributions package

Submodules

pycbc.distributions.angular module

This modules provides classes for evaluating angular distributions.

class pycbc.distributions.angular.CosAngle(**params)[source]

Bases: SinAngle

A cosine distribution. This is the same thing as a sine distribution, but with the domain shifted to [-pi/2, pi/2]. See SinAngle for more details.

Parameters:

**params – The keyword arguments should provide the names of parameters and (optionally) their corresponding bounds, as either boundaries.Bounds instances or tuples. The bounds must be in [-PI/2, PI/2].

name = 'cos_angle'
class pycbc.distributions.angular.SinAngle(**params)[source]

Bases: UniformAngle

A sine distribution; the pdf of each parameter theta is given by:

..math::

p(theta) = frac{sin theta}{costheta_0 - costheta_1}, theta_0 leq theta < theta_1,

and 0 otherwise. Here, \(\theta_0, \theta_1\) are the bounds of the parameter.

The domain of this distribution is [0, pi]. This is accomplished by putting hard boundaries at [0, pi]. Bounds may be provided to further limit the range for which the pdf has support. As with UniformAngle, these are initialized in radians.

Parameters:

**params – The keyword arguments should provide the names of parameters and (optionally) their corresponding bounds, as either boundaries.Bounds instances or tuples. The bounds must be in [0,PI]. These are converted to radians for storage. None may also be passed; in that case, the domain bounds will be used.

name = 'sin_angle'
class pycbc.distributions.angular.UniformAngle(cyclic_domain=False, **params)[source]

Bases: Uniform

A uniform distribution in which the dependent variable is between [0,2pi).

The domain of the distribution may optionally be made cyclic using the cyclic_domain parameter.

Bounds may be provided to limit the range for which the pdf has support. If provided, the parameter bounds are in radians.

Parameters:
  • cyclic_domain ({False, bool}) – If True, cyclic bounds on [0, 2pi) are applied to all values when evaluating the pdf. This is done prior to any additional bounds specified for a parameter are applied. Default is False.

  • **params – The keyword arguments should provide the names of parameters and (optionally) their corresponding bounds, as either boundaries.Bounds instances or tuples. The bounds must be in [0,2PI). These are converted to radians for storage. None may also be passed; in that case, the domain bounds will be used.

Notes

For more information, see Uniform.

apply_boundary_conditions(**kwargs)[source]

Maps values to be in [0, 2pi) (the domain) first, before applying any additional boundary conditions.

Parameters:

**kwargs – The keyword args should be the name of a parameter and value to apply its boundary conditions to. The arguments need not include all of the parameters in self.

Returns:

A dictionary of the parameter names and the conditioned values.

Return type:

dict

property domain

Returns the domain of the distribution.

classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file.

The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file. By default, only the name of the distribution (uniform_angle) needs to be specified. This will results in a uniform prior on [0, 2pi). To make the domain cyclic, add cyclic_domain =. To specify boundaries that are not [0, 2pi), add (min|max)-var arguments, where var is the name of the variable.

For example, this will initialize a variable called theta with a uniform distribution on [0, 2pi) without cyclic boundaries:

[{section}-theta]
name = uniform_angle

This will make the domain cyclic on [0, 2pi):

[{section}-theta]
name = uniform_angle
cyclic_domain =
Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance from the pycbc.inference.prior module.

Return type:

UniformAngle

name = 'uniform_angle'
class pycbc.distributions.angular.UniformSolidAngle(polar_angle=None, azimuthal_angle=None, polar_bounds=None, azimuthal_bounds=None, azimuthal_cyclic_domain=False)[source]

Bases: BoundedDist

A distribution that is uniform in the solid angle of a sphere. The names of the two angluar parameters can be specified on initalization.

Parameters:
  • polar_angle ({'theta', str}) – The name of the polar angle.

  • azimuthal_angle ({'phi', str}) – The name of the azimuthal angle.

  • polar_bounds ({None, tuple}) – Limit the polar angle to the given bounds. If None provided, the polar angle will vary from 0 (the north pole) to pi (the south pole). The bounds should be specified as factors of pi. For example, to limit the distribution to the northern hemisphere, set polar_bounds=(0,0.5).

  • azimuthal_bounds ({None, tuple}) – Limit the azimuthal angle to the given bounds. If None provided, the azimuthal angle will vary from 0 to 2pi. The bounds should be specified as factors of pi. For example, to limit the distribution to the one hemisphere, set azimuthal_bounds=(0,1).

  • azimuthal_cyclic_domain ({False, bool}) – Make the domain of the azimuthal angle be cyclic; i.e., azimuthal values are constrained to be in [0, 2pi) using cyclic boundaries prior to applying any other boundary conditions and prior to evaluating the pdf. Default is False.

apply_boundary_conditions(**kwargs)[source]

Maps the given values to be within the domain of the azimuthal and polar angles, before applying any other boundary conditions.

Parameters:

**kwargs – The keyword args must include values for both the azimuthal and polar angle, using the names they were initilialized with. For example, if polar_angle=’theta’ and azimuthal_angle=`phi, then the keyword args must be theta={val1}, phi={val2}.

Returns:

A dictionary of the parameter names and the conditioned values.

Return type:

dict

property azimuthal_angle

The name of the azimuthal angle.

Type:

str

property bounds

The bounds on each angle. The keys are the names of the polar and azimuthal angles, the values are the minimum and maximum of each, in radians. For example, if the distribution was initialized with polar_angle=’theta’, polar_bounds=(0,0.5) then the bounds will have ‘theta’: 0, 1.5707963267948966 as an entry.

Type:

dict

classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file.

The section must have the names of the polar and azimuthal angles in the tag part of the section header. For example:

[prior-theta+phi]
name = uniform_solidangle

If nothing else is provided, the default names and bounds of the polar and azimuthal angles will be used. To specify a different name for each angle, set the polar-angle and azimuthal-angle attributes. For example:

[prior-foo+bar]
name = uniform_solidangle
polar-angle = foo
azimuthal-angle = bar

Note that the names of the variable args in the tag part of the section name must match the names of the polar and azimuthal angles.

Bounds may also be specified for each angle, as factors of pi. For example:

[prior-theta+phi]
polar-angle = theta
azimuthal-angle = phi
min-theta = 0
max-theta = 0.5

This will return a distribution that is uniform in the upper hemisphere.

By default, the domain of the azimuthal angle is [0, 2pi). To make this domain cyclic, add azimuthal_cyclic_domain =.

Parameters:
  • cp (ConfigParser instance) – The config file.

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

  • variable_args (str) – The names of the parameters for this distribution, separated by VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance from the pycbc.inference.prior module.

Return type:

UniformSolidAngle

name = 'uniform_solidangle'
property polar_angle

The name of the polar angle.

Type:

str

pycbc.distributions.arbitrary module

This modules provides classes for evaluating arbitrary distributions from a file.

class pycbc.distributions.arbitrary.Arbitrary(bounds=None, bandwidth='scott', **kwargs)[source]

Bases: BoundedDist

A distribution constructed from a set of parameter values using a kde. Bounds may be optionally provided to limit the range.

Parameters:
  • bounds (dict, optional) – Independent bounds on one or more parameters may be provided to limit the range of the kde.

  • bandwidth (str, optional) – Set the bandwidth method for the KDE. See scipy.stats.gaussian_kde() for details. Default is “scott”.

  • **params – The keyword arguments should provide the names of the parameters and a list of their parameter values. If multiple parameters are provided, a single kde will be produced with dimension equal to the number of parameters.

classmethod from_config(cp, section, variable_args)[source]

Raises a NotImplementedError; to load from a config file, use FromFile.

static get_kde_from_arrays(*arrays)[source]

Constructs a KDE from the given arrays.

*arrays :

Each argument should be a 1D numpy array to construct the kde from. The resulting KDE will have dimension given by the number of parameters.

property kde
name = 'arbitrary'
property params

The list of parameter names.

Type:

list of strings

rvs(size=1, param=None)[source]

Gives a set of random values drawn from the kde.

Parameters:
  • size ({1, int}) – The number of values to generate; default is 1.

  • param ({None, string}) – If provided, will just return values for the given parameter. Otherwise, returns random values for each parameter.

Returns:

The random values in a numpy structured array. If a param was specified, the array will only have an element corresponding to the given parameter. Otherwise, the array will have an element for each parameter in self’s params.

Return type:

structured array

set_bandwidth(set_bw='scott')[source]
class pycbc.distributions.arbitrary.FromFile(filename=None, datagroup=None, **params)[source]

Bases: Arbitrary

A distribution that reads the values of the parameter(s) from an hdf file, computes the kde to construct the pdf, and draws random variables from it.

Parameters:
  • filename (str) – The path to an hdf file containing the values of the parameters that want to be used to construct the distribution. Each parameter should be a separate dataset in the hdf file, and all datasets should have the same size. For example, to give a prior for mass1 and mass2 from file f, f[‘mass1’] and f[‘mass2’] contain the n values for each parameter.

  • datagroup (str, optional) – The name of the group to look in for the samples. For example, if datagroup = 'samples', then parameter param will be retrived from f['samples'][param]. If none provided (the default) the data sets will be assumed to be in the top level directory of the file.

  • **params – The keyword arguments should provide the names of the parameters to be read from the file and (optionally) their bounds. If no parameters are provided, it will use all the parameters found in the file. To provide bounds, specify e.g. mass1=[10,100]. Otherwise, mass1=None.

norm

The normalization of the multi-dimensional pdf.

Type:

float

lognorm

The log of the normalization.

Type:

float

kde

The kde obtained from the values in the file.

property filename

The path to the file containing values for the parameter(s).

Type:

str

classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file.

The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

The file to construct the distribution from must be provided by setting filename. Boundary arguments can be provided in the same way as described in get_param_bounds_from_config.

[{section}-{tag}]
name = fromfile
filename = ra_prior.hdf
min-ra = 0
max-ra = 6.28
Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by prior.VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance from the pycbc.inference.prior module.

Return type:

BoundedDist

get_arrays_from_file(params_file, params=None)[source]

Reads the values of one or more parameters from an hdf file and returns as a dictionary.

Parameters:
  • params_file (str) – The hdf file that contains the values of the parameters.

  • params ({None, list}) – If provided, will just retrieve the given parameter names.

Returns:

A dictionary of the parameters mapping param_name -> array.

Return type:

dict

name = 'fromfile'

pycbc.distributions.bounded module

This modules provides classes for evaluating distributions with bounds.

class pycbc.distributions.bounded.BoundedDist(**params)[source]

Bases: object

A generic class for storing common properties of distributions in which each parameter has a minimum and maximum value.

Parameters:

**params – The keyword arguments should provide the names of parameters and their corresponding bounds, as either tuples or a boundaries.Bounds instance.

apply_boundary_conditions(**kwargs)[source]

Applies any boundary conditions to the given values (e.g., applying cyclic conditions, and/or reflecting values off of boundaries). This is done by running apply_conditions of each bounds in self on the corresponding value. See boundaries.Bounds.apply_conditions for details.

Parameters:

**kwargs – The keyword args should be the name of a parameter and value to apply its boundary conditions to. The arguments need not include all of the parameters in self. Any unrecognized arguments are ignored.

Returns:

A dictionary of the parameter names and the conditioned values.

Return type:

dict

property bounds

A dictionary of the parameter names and their bounds.

Type:

dict

cdfinv(**kwds)[source]

Return the inverse cdf to map the unit interval to parameter bounds. You must provide a keyword for every parameter.

classmethod from_config(cp, section, variable_args, bounds_required=False)[source]

Returns a distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by prior.VARARGS_DELIM. These must appear in the “tag” part of the section header.

  • bounds_required ({False, bool}) – If True, raise a ValueError if a min and max are not provided for every parameter. Otherwise, the prior will be initialized with the parameter set to None. Even if bounds are not required, a ValueError will be raised if only one bound is provided; i.e., either both bounds need to provided or no bounds.

Returns:

A distribution instance from the pycbc.distribution subpackage.

Return type:

BoundedDist

logpdf(**kwargs)[source]

Returns the log of the pdf at the given values. The keyword arguments must contain all of parameters in self’s params. Unrecognized arguments are ignored. Any boundary conditions are applied to the values before the pdf is evaluated.

property params

The list of parameter names.

Type:

list of strings

pdf(**kwargs)[source]

Returns the pdf at the given values. The keyword arguments must contain all of parameters in self’s params. Unrecognized arguments are ignored. Any boundary conditions are applied to the values before the pdf is evaluated.

rvs(size=1, **kwds)[source]

Draw random value

pycbc.distributions.bounded.bounded_from_config(cls, cp, section, variable_args, bounds_required=False, additional_opts=None)[source]

Returns a bounded distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Parameters:
  • cls (pycbc.prior class) – The class to initialize with.

  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by prior.VARARGS_DELIM. These must appear in the “tag” part of the section header.

  • bounds_required ({False, bool}) – If True, raise a ValueError if a min and max are not provided for every parameter. Otherwise, the prior will be initialized with the parameter set to None. Even if bounds are not required, a ValueError will be raised if only one bound is provided; i.e., either both bounds need to provided or no bounds.

  • additional_opts ({None, dict}) – Provide additional options to be passed to the distribution class; should be a dictionary specifying option -> value. If an option is provided that also exists in the config file, the value provided will be used instead of being read from the file.

Returns:

An instance of the given class.

Return type:

cls

pycbc.distributions.bounded.get_param_bounds_from_config(cp, section, tag, param)[source]

Gets bounds for the given parameter from a section in a config file.

Minimum and maximum values for bounds are specified by adding min-{param} and max-{param} options, where {param} is the name of the parameter. The types of boundary (open, closed, or reflected) to create may also be specified by adding options btype-min-{param} and btype-max-{param}. Cyclic conditions can be adding option cyclic-{param}. If no btype arguments are provided, the left bound will be closed and the right open.

For example, the following will create right-open bounds for parameter foo:

[{section}-{tag}]
min-foo = -1
max-foo = 1

This would make the boundaries cyclic:

[{section}-{tag}]
min-foo = -1
max-foo = 1
cyclic-foo =

For more details on boundary types and their meaning, see boundaries.Bounds.

If the parameter is not found in the section will just return None (in this case, all btype and cyclic arguments are ignored for that parameter). If bounds are specified, both a minimum and maximum must be provided, else a Value or Type Error will be raised.

Parameters:
  • cp (ConfigParser instance) – The config file.

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

  • tag (str) – Any tag in the section name. The full section name searched for in the config file is {section}(-{tag}).

  • param (str) – The name of the parameter to retrieve bounds for.

Returns:

bounds – If bounds were provided, a boundaries.Bounds instance representing the bounds. Otherwise, None.

Return type:

{Bounds instance | None}

pycbc.distributions.constraints module

This modules provides classes for evaluating multi-dimensional constraints.

class pycbc.distributions.constraints.Constraint(constraint_arg, static_args=None, transforms=None, **kwargs)[source]

Bases: object

Creates a constraint that evaluates to True if parameters obey the constraint and False if they do not.

name = 'custom'
class pycbc.distributions.constraints.SupernovaeConvexHull(constraint_arg, transforms=None, **kwargs)[source]

Bases: Constraint

Pre defined constraint for core-collapse waveforms that checks whether a given set of coefficients lie within the convex hull of the coefficients of the principal component basis vectors.

name = 'supernovae_convex_hull'
required_parameters = ['coeff_0', 'coeff_1']

pycbc.distributions.external module

This modules provides classes for evaluating PDF, logPDF, CDF and inverse CDF from external arbitrary distributions, and drawing samples from them.

class pycbc.distributions.external.DistributionFunctionFromFile(params=None, file_path=None, column_index=None, **kwargs)[source]

Bases: External

Evaluating PDF, logPDF, CDF and inverse CDF from the external

density function.

Instances of this class can be called like a distribution in the .ini file, when used with pycbc.distributions.external.External. Please see the example in the External class.

Parameters:
  • parameter ({'file_path', 'column_index'}) – The path of the external density function’s .txt file, and the column index of the density distribution. By default, the first column should be the values of a certain parameter, such as “mass”, other columns should be the corresponding density values (as a function of that parameter). If you add the name of the parameter in the first row, please add the ‘#’ at the beginning.

  • **kwargs – All other keyword args are passed to scipy.integrate.quad to control the numerical accuracy of the inverse CDF. If not be provided, will use the default values in self.__init__.

Notes

This class is different from pycbc.distributions.arbitrary.FromFile, which needs samples from the hdf file to construct the PDF by using KDE. This class reads in any continuous functions of the parameter.

name = 'external_func_fromfile'
class pycbc.distributions.external.External(params=None, logpdf=None, rvs=None, cdfinv=None, **kwds)[source]

Bases: object

Distribution defined by external cdfinv and logpdf functions

To add to an inference configuration file:

[prior-param1+param2]
name = external
module = custom_mod
logpdf = custom_function_name
cdfinv = custom_function_name2

Or call DistributionFunctionFromFile in the .ini file:

[prior-param]
name = external_func_fromfile
module = pycbc.distributions.external
file_path = path
column_index = index
logpdf = _logpdf
cdfinv = _cdfinv
Parameters:
  • params (list) – list of parameter names

  • custom_mod (module) – module from which logpdf and cdfinv functions can be imported

  • logpdf (function) – function which returns the logpdf

  • cdfinv (function) – function which applies the invcdf

Examples

To instantate by hand and example of function format. You must provide the logpdf function, and you may either provide the rvs or cdfinv function. If the cdfinv is provided, but not the rvs, the random values will be calculated using the cdfinv function.

>>> import numpy
>>> params = ['x', 'y']
>>> def logpdf(x=None, y=None):
...     p = numpy.ones(len(x))
...     return p
>>>
>>> def cdfinv(**kwds):
...     return kwds
>>> e = External(['x', 'y'], logpdf, cdfinv=cdfinv)
>>> e.rvs(size=10)
apply_boundary_conditions(**params)[source]
classmethod from_config(cp, section, variable_args)[source]
name = 'external'
rvs(size=1, **kwds)[source]

Draw random value

pycbc.distributions.fixedsamples module

This modules provides classes for evaluating distributions based on a fixed set of points

class pycbc.distributions.fixedsamples.FixedSamples(params, samples)[source]

Bases: object

A distribution consisting of a collection of a large number of fixed points. Only these values can be drawn from, so the number of points may need to be large to properly reflect the paramter space. This distribution is intended to aid in using nested samplers for semi-abitrary or complicated distributions where it is possible to provide or draw samples but less straightforward to provide an analytic invcdf. This class numerically approximates the invcdf for 1 or 2 dimensional distributions (but no higher).

Parameters:
  • params – This of parameters this distribution should use

  • samples (dict of arrays or FieldArray) – Sampled points of the distribution. May contain transformed parameters which are different from the original distribution. If so, an inverse mapping is provided to associate points with other parameters provided.

apply_boundary_conditions(**params)[source]

Apply boundary conditions (none here)

cdfinv(**original)[source]

Map unit cube to parameters in the space

classmethod from_config(cp, section, tag)[source]

Return instance based on config file

Return a new instance based on the config file. This will draw from a single distribution section provided in the config file and apply a single transformation section if desired. If a transformation is applied, an inverse mapping is also provided for use in the config file.

name = 'fixed_samples'
rvs(size=1, **kwds)[source]

Draw random value

pycbc.distributions.gaussian module

This modules provides classes for evaluating Gaussian distributions.

class pycbc.distributions.gaussian.Gaussian(**params)[source]

Bases: BoundedDist

A Gaussian distribution on the given parameters; the parameters are independent of each other.

Bounds can be provided on each parameter, in which case the distribution will be a truncated Gaussian distribution. The PDF of a truncated Gaussian distribution is given by:

\[p(x|a, b, \mu,\sigma) = \frac{1}{\sqrt{2 \pi \sigma^2}}\frac{e^{- \frac{\left( x - \mu \right)^2}{2 \sigma^2}}}{\Phi(b|\mu, \sigma) - \Phi(a|\mu, \sigma)},\]

where \(\mu\) is the mean, \(\sigma^2\) is the variance, \(a,b\) are the bounds, and \(\Phi\) is the cumulative distribution of an unbounded normal distribution, given by:

\[\Phi(x|\mu, \sigma) = \frac{1}{2}\left[1 + \mathrm{erf}\left(\frac{x-\mu}{\sigma \sqrt{2}}\right)\right].\]

Note that if \([a,b) = [-\infty, \infty)\), this reduces to a standard Gaussian distribution.

Instances of this class can be called like a function. By default, logpdf will be called, but this can be changed by setting the class’s __call__ method to its pdf method.

Parameters:

**params – The keyword arguments should provide the names of parameters and (optionally) some bounds, as either a tuple or a boundaries.Bounds instance. The mean and variance of each parameter can be provided by additional keyword arguments that have _mean and _var adding to the parameter name. For example, foo=(-2,10), foo_mean=3, foo_var=2 would create a truncated Gaussian with mean 3 and variance 2, bounded between \([-2, 10)\). If no mean or variance is provided, the distribution will have 0 mean and unit variance. If None is provided for the bounds, the distribution will be a normal, unbounded Gaussian (equivalent to setting the bounds to [-inf, inf)).

Examples

Create an unbounded Gaussian distribution with zero mean and unit variance: >>> dist = distributions.Gaussian(mass1=None)

Create a bounded Gaussian distribution on \([1,10)\) with a mean of 3 and a variance of 2: >>> dist = distributions.Gaussian(mass1=(1,10), mass1_mean=3, mass1_var=2)

Create a bounded Gaussian distribution with the same parameters, but with cyclic boundary conditions: >>> dist = distributions.Gaussian(mass1=Bounds(1,10, cyclic=True), mass1_mean=3, mass1_var=2)

cdf(param, value)[source]

Returns the CDF of the given parameter value.

classmethod from_config(cp, section, variable_args)[source]

Returns a Gaussian distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Boundary arguments should be provided in the same way as described in get_param_bounds_from_config. In addition, the mean and variance of each parameter can be specified by setting {param}_mean and {param}_var, respectively. For example, the following would create a truncated Gaussian distribution between 0 and 6.28 for a parameter called phi with mean 3.14 and variance 0.5 that is cyclic:

[{section}-{tag}]
min-phi = 0
max-phi = 6.28
phi_mean = 3.14
phi_var = 0.5
cyclic =
Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by prior.VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance from the pycbc.inference.prior module.

Return type:

Gaussian

property mean
name = 'gaussian'
property var

pycbc.distributions.joint module

This module provides classes to describe joint distributions

class pycbc.distributions.joint.JointDistribution(variable_args, *distributions, **kwargs)[source]

Bases: object

Callable class that calculates the joint distribution built from a set of distributions.

Parameters:
  • variable_args (list) – A list of strings that contain the names of the variable parameters and the order they are expected when the class is called.

  • *distributions – The rest of the arguments must be instances of distributions describing the individual distributions on the variable parameters. A single distribution may contain multiple parameters. The set of all params across the distributions (retrieved from the distributions’ params attribute) must be the same as the set of variable_args provided.

  • **kwargs – Valid keyword arguments include: constraints : a list of functions that accept a dict of parameters with the parameter name as the key. If the constraint is satisfied the function should return True, if the constraint is violated, then the function should return False. n_test_samples : number of random draws used to fix pdf normalization factor after applying constraints.

variable_args

The parameters expected when the evaluator is called.

Type:

tuple

distributions

The distributions for the parameters.

Type:

list

constraints

A list of functions to test if parameter values obey multi-dimensional constraints.

Type:

list

Examples

An example of creating a joint distribution with constraint that total mass must be below 30.

>>> from pycbc.distributions import Uniform, JointDistribution
>>> def mtotal_lt_30(params):
    ...    return params["mass1"] + params["mass2"] < 30
>>> mass_lim = (2, 50)
>>> uniform_prior = Uniform(mass1=mass_lim, mass2=mass_lim)
>>> prior_eval = JointDistribution(["mass1", "mass2"], uniform_prior,
    ...                               constraints=[mtotal_lt_30])
>>> print(prior_eval(mass1=20, mass2=1))
apply_boundary_conditions(**params)[source]

Applies each distributions’ boundary conditions to the given list of parameters, returning a new list with the conditions applied.

Parameters:

**params – Keyword arguments should give the parameters to apply the conditions to.

Returns:

A dictionary of the parameters after each distribution’s apply_boundary_conditions function has been applied.

Return type:

dict

property bounds

Get the dict of boundaries

cdfinv(**original)[source]

Apply the inverse cdf to the array of values [0, 1]. Every variable parameter must be given as a keyword argument.

contains(params)[source]
Evaluates whether the given parameters satisfy the boundary

conditions, boundaries, and constraints. This method is different from within_constraints, that method only check the constraints.

Parameters:

params (dict, FieldArray, numpy.record, or numpy.ndarray) – The parameter values to evaluate.

Returns:

If params was an array, or if params a dictionary and one or more of the parameters are arrays, will return an array of booleans. Otherwise, a boolean.

Return type:

(array of) bool

property cyclic

Get list of which parameters are cyclic

name = 'joint'
rvs(size=1)[source]

Rejection samples the parameter space.

property well_reflected

Get list of which parameters are well reflected

within_constraints(params)[source]

Evaluates whether the given parameters satisfy the constraints.

Parameters:

params (dict, FieldArray, numpy.record, or numpy.ndarray) – The parameter values to evaluate.

Returns:

If params was an array, or if params a dictionary and one or more of the parameters are arrays, will return an array of booleans. Otherwise, a boolean.

Return type:

(array of) bool

pycbc.distributions.mass module

This modules provides classes for evaluating distributions for mchirp and q (i.e., mass ratio) from uniform component mass.

class pycbc.distributions.mass.MchirpfromUniformMass1Mass2(dim=2, **params)[source]

Bases: UniformPowerLaw

A distribution for chirp mass from uniform component mass + constraints given by chirp mass. This is a special case for UniformPowerLaw with index 1. For more details see UniformPowerLaw.

The parameters (i.e. **params) are independent of each other. Instances of this class can be called like a function. By default, logpdf will be called, but this can be changed by setting the class’s __call__ method to its pdf method.

Derivation for the probability density function:

\[P(m_1,m_2)dm_1dm_2 = P(\mathcal{M}_c,q)d\mathcal{M}_cdq\]

Where \(\mathcal{M}_c\) is chirp mass and \(q\) is mass ratio, \(m_1\) and \(m_2\) are component masses. The jacobian to transform chirp mass and mass ratio to component masses is

\[\frac{\partial(m_1,m_2)}{\partial(\mathcal{M}_c,q)} = \ \mathcal{M}_c \left(\frac{1+q}{q^3}\right)^{2/5}\]

(https://github.com/gwastro/pycbc/blob/master/pycbc/transforms.py#L416.)

Because \(P(m_1,m_2) = const\), then

\[P(\mathcal{M}_c,q) = P(\mathcal{M}_c)P(q)\propto \mathcal{M}_c \left(\frac{1+q}{q^3}\right)^{2/5}`.\]

Therefore,

\[P(\mathcal{M}_c) \propto \mathcal{M}_c\]

and

\[P(q) \propto \left(\frac{1+q}{q^3}\right)^{2/5}\]

Examples

Generate 10000 random numbers from this distribution in [5,100]

>>> from pycbc import distributions as dist
>>> minmc = 5, maxmc = 100, size = 10000
>>> mc = dist.MchirpfromUniformMass1Mass2(value=(minmc,maxmc)).rvs(size)

The settings in the configuration file for pycbc_inference should be

[variable_params]
mchirp =
[prior-mchirp]
name = mchirp_from_uniform_mass1_mass2
min-mchirp = 10
max-mchirp = 80
Parameters:

**params – The keyword arguments should provide the names of parameters and their corresponding bounds, as either tuples or a boundaries.Bounds instance.

name = 'mchirp_from_uniform_mass1_mass2'
class pycbc.distributions.mass.QfromUniformMass1Mass2(**params)[source]

Bases: BoundedDist

A distribution for mass ratio (i.e., q) from uniform component mass + constraints given by q.

The parameters (i.e. **params) are independent of each other. Instances of this class can be called like a function. By default, logpdf will be called, but this can be changed by setting the class’s __call__ method to its pdf method.

For mathematical derivation see the documentation above in the class MchirpfromUniformMass1Mass2.

Parameters:

**params – The keyword arguments should provide the names of parameters and their corresponding bounds, as either tuples or a boundaries.Bounds instance.

Examples

Generate 10000 random numbers from this distribution in [1,8]

>>> from pycbc import distributions as dist
>>> minq = 1, maxq = 8, size = 10000
>>> q = dist.QfromUniformMass1Mass2(value=(minq,maxq)).rvs(size)
classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Example:

[variable_params]
q =
[prior-q]
name = q_from_uniform_mass1_mass2
min-q = 1
max-q = 8
Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

  • QfromUniformMass1Mass2 – A distribution instance from the pycbc.distributions.bounded

  • module.

property lognorm

The log of the normalization.

Type:

float

name = 'q_from_uniform_mass1_mass2'
property norm

The normalization of the multi-dimensional pdf.

Type:

float

rvs(size=1, param=None)[source]

Gives a set of random values drawn from this distribution.

Parameters:
  • size ({1, int}) – The number of values to generate; default is 1.

  • param ({None, string}) – If provided, will just return values for the given parameter. Otherwise, returns random values for each parameter.

Returns:

The random values in a numpy structured array. If a param was specified, the array will only have an element corresponding to the given parameter. Otherwise, the array will have an element for each parameter in self’s params.

Return type:

structured array

pycbc.distributions.power_law module

This modules provides classes for evaluating distributions where the probability density function is a power law.

class pycbc.distributions.power_law.UniformPowerLaw(dim=None, **params)[source]

Bases: BoundedDist

For a uniform distribution in power law. The parameters are independent of each other. Instances of this class can be called like a function. By default, logpdf will be called, but this can be changed by setting the class’s __call__ method to its pdf method.

The cumulative distribution function (CDF) will be the ratio of volumes:

\[F(r) = \frac{V(r)}{V(R)}\]

Where \(R\) is the radius of the sphere. So we can write our probability density function (PDF) as:

\[f(r) = c r^n\]

For generality we use \(n\) for the dimension of the volume element, eg. \(n=2\) for a 3-dimensional sphere. And use \(c\) as a general constant.

So now we calculate the CDF in general for this type of PDF:

\[F(r) = \int f(r) dr = \int c r^n dr = \frac{1}{n + 1} c r^{n + 1} + k\]

Now with the definition of the CDF at radius \(r_{l}\) is equal to 0 and at radius \(r_{h}\) is equal to 1 we find that the constant from integration from this system of equations:

\[1 = \frac{1}{n + 1} c ((r_{h})^{n + 1} - (r_{l})^{n + 1}) + k\]

Can see that \(c = (n + 1) / ((r_{h})^{n + 1} - (r_{l})^{n + 1}))\). And \(k\) is:

\[k = - \frac{r_{l}^{n + 1}}{(r_{h})^{n + 1} - (r_{l})^{n + 1}}\]

Can see that \(c= \frac{n + 1}{R^{n + 1}}\). So can see that the CDF is:

\[F(r) = \frac{1}{(r_{h})^{n + 1} - (r_{l})^{n + 1}} r^{n + 1} - \frac{r_{l}^{n + 1}}{(r_{h})^{n + 1} - (r_{l})^{n + 1}}\]

And the PDF is the derivative of the CDF:

\[f(r) = \frac{(n + 1)}{(r_{h})^{n + 1} - (r_{l})^{n + 1}} (r)^n\]

Now we use the probabilty integral transform method to get sampling on uniform numbers from a continuous random variable. To do this we find the inverse of the CDF evaluated for uniform numbers:

\[F(r) = u = \frac{1}{(r_{h})^{n + 1} - (r_{l})^{n + 1}} r^{n + 1} - \frac{r_{l}^{n + 1}}{(r_{h})^{n + 1} - (r_{l})^{n + 1}}\]

And find \(F^{-1}(u)\) gives:

\[u = \frac{1}{n + 1} \frac{(r_{h})^{n + 1} - (r_{l})^{n + 1}} - \frac{r_{l}^{n + 1}}{(r_{h})^{n + 1} - (r_{l})^{n + 1}}\]

And solving for \(r\) gives:

\[r = ( ((r_{h})^{n + 1} - (r_{l})^{n + 1}) u + (r_{l})^{n + 1})^{\frac{1}{n + 1}}\]

Therefore the radius can be sampled by taking the n-th root of uniform numbers and multiplying by the radius offset by the lower bound radius.

**params :

The keyword arguments should provide the names of parameters and their corresponding bounds, as either tuples or a boundaries.Bounds instance.

dim

The dimension of volume space. In the notation above dim is \(n+1\). For a 3-dimensional sphere this is 3.

Type:

int

classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by prior.VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance from the pycbc.inference.prior module.

Return type:

Uniform

property lognorm

The log of the normalization.

Type:

float

name = 'uniform_power_law'
property norm

The normalization of the multi-dimensional pdf.

Type:

float

class pycbc.distributions.power_law.UniformRadius(dim=3, **params)[source]

Bases: UniformPowerLaw

For a uniform distribution in volume using spherical coordinates, this is the distriubtion to use for the radius.

For more details see UniformPowerLaw.

name = 'uniform_radius'

pycbc.distributions.qnm module

class pycbc.distributions.qnm.UniformF0Tau(f0=None, tau=None, final_mass=None, final_spin=None, rdfreq='f0', damping_time='tau', norm_tolerance=0.001, norm_seed=0)[source]

Bases: Uniform

A distribution uniform in QNM frequency and damping time.

Constraints may be placed to exclude frequencies and damping times corresponding to specific masses and spins.

To ensure a properly normalized pdf that accounts for the constraints on final mass and spin, a renormalization factor is calculated upon initialization. This is calculated numerically: f0 and tau are drawn randomly, then the norm is scaled by the fraction of points that yield final masses and spins within the constraints. The norm_tolerance keyword arguments sets the error on the estimate of the norm from this numerical method. If this value is too large, such that no points are found in the allowed region, a ValueError is raised.

Parameters:
  • f0 (tuple or boundaries.Bounds) – The range of QNM frequencies (in Hz).

  • tau (tuple or boundaries.Bounds) – The range of QNM damping times (in s).

  • final_mass (tuple or boundaries.Bounds, optional) – The range of final masses to allow. Default is [0,inf).

  • final_spin (tuple or boundaries.Bounds, optional) – The range final spins to allow. Must be in [-0.996, 0.996], which is the default.

  • rdfreq (str, optional) – Use the given string as the name for the f0 parameter. Default is ‘f0’.

  • damping_time (str, optional) – Use the given string as the name for the tau parameter. Default is ‘tau’.

  • norm_tolerance (float, optional) – The tolerance on the estimate of the normalization. Default is 1e-3.

  • norm_seed (int, optional) – Seed to use for the random number generator when estimating the norm. Default is 0. After the norm is estimated, the random number generator is set back to the state it was in upon initialization.

Examples

Create a distribution:

>>> dist = UniformF0Tau(f0=(10., 2048.), tau=(1e-4,1e-2))

Check that all random samples drawn from the distribution yield final masses > 1:

>>> from pycbc import conversions
>>> samples = dist.rvs(size=1000)
>>> (conversions.final_mass_from_f0_tau(samples['f0'],
        samples['tau']) > 1.).all()
True

Create a distribution with tighter bounds on final mass and spin:

>>> dist = UniformF0Tau(f0=(10., 2048.), tau=(1e-4,1e-2),
        final_mass=(20., 200.), final_spin=(0,0.996))

Check that all random samples drawn from the distribution are in the final mass and spin constraints:

>>> samples = dist.rvs(size=1000)
>>> (conversions.final_mass_from_f0_tau(samples['f0'],
        samples['tau']) >= 20.).all()
True
>>> (conversions.final_mass_from_f0_tau(samples['f0'],
        samples['tau']) < 200.).all()
True
>>> (conversions.final_spin_from_f0_tau(samples['f0'],
        samples['tau']) >= 0.).all()
True
>>> (conversions.final_spin_from_f0_tau(samples['f0'],
        samples['tau']) < 0.996).all()
True
classmethod from_config(cp, section, variable_args)[source]

Initialize this class from a config file.

Bounds on f0, tau, final_mass and final_spin should be specified by providing min-{param} and max-{param}. If the f0 or tau param should be renamed, rdfreq and damping_time should be provided; these must match variable_args. If rdfreq and damping_time are not provided, variable_args are expected to be f0 and tau.

Only min/max-f0 and min/max-tau need to be provided.

Example:

[{section}-f0+tau]
name = uniform_f0_tau
min-f0 = 10
max-f0 = 2048
min-tau = 0.0001
max-tau = 0.010
min-final_mass = 10
Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – WorkflowConfigParser instance to read.

  • section (str) – The name of the section to read.

  • variable_args (str) – The name of the variable args. These should be separated by pycbc.VARARGS_DELIM.

Returns:

This class initialized with the parameters provided in the config file.

Return type:

UniformF0Tau

name = 'uniform_f0_tau'
rvs(size=1)[source]

Draw random samples from this distribution.

Parameters:

size (int, optional) – The number of draws to do. Default is 1.

Returns:

A structured array of the random draws.

Return type:

array

pycbc.distributions.sky_location module

This modules provides classes for evaluating sky distributions in right ascension and declination.

class pycbc.distributions.sky_location.FisherSky(**params)[source]

Bases: object

A distribution that returns a random angle drawn from an approximate Von_Mises-Fisher distribution. Assumes that the Fisher concentration parameter is large, so that we can draw the samples from a simple rotationally-invariant distribution centered at the North Pole (which factors as a uniform distribution for the right ascension, and a Rayleigh distribution for the declination, as described in Fabrycky and Winn 2009 ApJ 696 1230) and then rotate the samples to be centered around the specified mean position. As in UniformSky, the declination varies from π/2 to -π/2 and the right ascension varies from 0 to 2π.

Parameters:
  • mean_ra (float) – RA of the center of the distribution.

  • mean_dec (float) – Declination of the center of the distribution.

  • sigma (float) – Spread of the distribution. For the precise interpretation, see Eq 8 of Briggs et al 1999 ApJS 122 503. This should be smaller than about 20 deg for the approximation to be valid.

  • angle_unit (str) – Unit for the angle parameters: either “deg” or “rad”.

classmethod from_config(cp, section, variable_args)[source]
name = 'fisher_sky'
property params
rvs(size)[source]
class pycbc.distributions.sky_location.UniformSky(polar_angle=None, azimuthal_angle=None, polar_bounds=None, azimuthal_bounds=None, azimuthal_cyclic_domain=False)[source]

Bases: UniformSolidAngle

A distribution that is uniform on the sky. This is the same as UniformSolidAngle, except that the polar angle varies from pi/2 (the north pole) to -pi/2 (the south pole) instead of 0 to pi. Also, the default names are “dec” (declination) for the polar angle and “ra” (right ascension) for the azimuthal angle, instead of “theta” and “phi”.

name = 'uniform_sky'

pycbc.distributions.spins module

This modules provides spin distributions of CBCs.

class pycbc.distributions.spins.IndependentChiPChiEff(mass1=None, mass2=None, chi_eff=None, chi_a=None, xi_bounds=None, nsamples=None, seed=None)[source]

Bases: Arbitrary

A distribution such that \(\chi_{\mathrm{eff}}\) and \(\chi_p\) are uniform and independent of each other.

To ensure constraints are applied correctly, this distribution produces all three components of both spins as well as the component masses.

Parameters:
  • mass1 (BoundedDist, Bounds, or tuple) – The distribution or bounds to use for mass1. Must be either a BoundedDist giving the distribution on mass1, or bounds (as either a Bounds instance or a tuple) giving the minimum and maximum values to use for mass1. If the latter, a Uniform distribution will be used.

  • mass2 (BoundedDist, Bounds, or tuple) – The distribution or bounds to use for mass2. Syntax is the same as mass1.

  • chi_eff (BoundedDist, Bounds, or tuple; optional) – The distribution or bounds to use for \(chi_eff\). Syntax is the same as mass1, except that None may also be passed. In that case, (-1, 1) will be used for the bounds. Default is None.

  • chi_a (BoundedDist, Bounds, or tuple; optional) – The distribution or bounds to use for \(chi_a\). Syntax is the same as mass1, except that None may also be passed. In that case, (-1, 1) will be used for the bounds. Default is None.

  • xi_bounds (Bounds or tuple, optional) – The bounds to use for \(\xi_1\) and \(\xi_2\). Must be \(\in (0, 1)\). If None (the default), will be (0, 1).

  • nsamples (int, optional) – The number of samples to use for the internal kde. The larger the number of samples, the more accurate the pdf will be, but the longer it will take to evaluate. Default is 10000.

  • seed (int, optional) – Seed value to use for the number generator for the kde. The current random state of numpy will be saved prior to setting the seed. After the samples are generated, the state will be set back to what it was. If None provided, will use 0.

apply_boundary_conditions(**kwargs)[source]

Applies any boundary conditions to the given values (e.g., applying cyclic conditions, and/or reflecting values off of boundaries). This is done by running apply_conditions of each bounds in self on the corresponding value. See boundaries.Bounds.apply_conditions for details.

Parameters:

**kwargs – The keyword args should be the name of a parameter and value to apply its boundary conditions to. The arguments need not include all of the parameters in self. Any unrecognized arguments are ignored.

Returns:

A dictionary of the parameter names and the conditioned values.

Return type:

dict

classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by prior.VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance.

Return type:

IndependentChiPChiEff

name = 'independent_chip_chieff'
rvs(size=1, **kwargs)[source]

Returns random values for all of the parameters.

pycbc.distributions.uniform module

This modules provides classes for evaluating uniform distributions.

class pycbc.distributions.uniform.Uniform(**params)[source]

Bases: BoundedDist

A uniform distribution on the given parameters. The parameters are independent of each other. Instances of this class can be called like a function. By default, logpdf will be called, but this can be changed by setting the class’s __call__ method to its pdf method.

Parameters:

**params – The keyword arguments should provide the names of parameters and their corresponding bounds, as either tuples or a boundaries.Bounds instance.

Examples

Create a 2 dimensional uniform distribution:

>>> from pycbc import distributions
>>> dist = distributions.Uniform(mass1=(10.,50.), mass2=(10.,50.))

Get the log of the pdf at a particular value:

>>> dist.logpdf(mass1=25., mass2=10.)
    -7.3777589082278725

Do the same by calling the distribution:

>>> dist(mass1=25., mass2=10.)
    -7.3777589082278725

Generate some random values:

>>> dist.rvs(size=3)
    array([(36.90885758394699, 51.294212757995254),
           (39.109058546060346, 13.36220145743631),
           (34.49594465315212, 47.531953033719454)],
          dtype=[('mass1', '<f8'), ('mass2', '<f8')])

Initialize a uniform distribution using a boundaries.Bounds instance, with cyclic bounds:

>>> dist = distributions.Uniform(phi=Bounds(10, 50, cyclic=True))

Apply boundary conditions to a value:

>>> dist.apply_boundary_conditions(phi=60.)
    {'mass1': array(20.0)}

The boundary conditions are applied to the value before evaluating the pdf; note that the following returns a non-zero pdf. If the bounds were not cyclic, the following would return 0:

>>> dist.pdf(phi=60.)
    0.025
classmethod from_config(cp, section, variable_args)[source]

Returns a distribution based on a configuration file. The parameters for the distribution are retrieved from the section titled “[section-variable_args]” in the config file.

Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – A parsed configuration file that contains the distribution options.

  • section (str) – Name of the section in the configuration file.

  • variable_args (str) – The names of the parameters for this distribution, separated by VARARGS_DELIM. These must appear in the “tag” part of the section header.

Returns:

A distribution instance from the pycbc.inference.prior module.

Return type:

Uniform

property lognorm

The log of the normalization

Type:

float

name = 'uniform'
property norm

The normalization of the multi-dimensional pdf.

Type:

float

pycbc.distributions.uniform_log module

This modules provides classes for evaluating distributions whose logarithm are uniform.

class pycbc.distributions.uniform_log.UniformLog10(**params)[source]

Bases: Uniform

A uniform distribution on the log base 10 of the given parameters. The parameters are independent of each other. Instances of this class can be called like a function. By default, logpdf will be called.

Parameters:

**params – The keyword arguments should provide the names of parameters and their corresponding bounds, as either tuples or a boundaries.Bounds instance.

name = 'uniform_log10'

pycbc.distributions.utils module

This module provides functions for drawing samples from a standalone .ini file in a Python script, rather than in the command line.

pycbc.distributions.utils.draw_samples_from_config(path, num=1, seed=150914)[source]

Generate sampling points from a standalone .ini file.

Parameters:
  • path (str) – The path to the .ini file.

  • num (int) – The number of samples.

  • seed (int) – The random seed for sampling.

Returns:

samples – The parameter values and names of sample(s).

Return type:

pycbc.io.record.FieldArray

Examples

Draw a sample from the distribution defined in the .ini file:

>>> import numpy as np
>>> from pycbc.distributions.utils import draw_samples_from_config
>>> # A path to the .ini file.
>>> CONFIG_PATH = "./pycbc_bbh_prior.ini"
>>> random_seed = np.random.randint(low=0, high=2**32-1)
>>> sample = draw_samples_from_config(
>>>          path=CONFIG_PATH, num=1, seed=random_seed)
>>> # Print all parameters.
>>> print(sample.fieldnames)
>>> print(sample)
>>> # Print a certain parameter, for example 'mass1'.
>>> print(sample[0]['mass1'])
pycbc.distributions.utils.prior_from_config(cp, prior_section='prior')[source]

Loads a prior distribution from the given config file.

Parameters:
  • cp (pycbc.workflow.WorkflowConfigParser) – The config file to read.

  • sections (list of str, optional) – The sections to retrieve the prior from. If None (the default), will look in sections starting with ‘prior’.

Returns:

The prior distribution.

Return type:

distributions.JointDistribution

Module contents

This modules provides classes and functions for drawing and calculating the probability density function of distributions.

pycbc.distributions.read_constraints_from_config(cp, transforms=None, static_args=None, constraint_section='constraint')[source]

Loads parameter constraints from a configuration file.

Parameters:
  • cp (WorkflowConfigParser) – An open config parser to read from.

  • transforms (list, optional) – List of transforms to apply to parameters before applying constraints.

  • static_args (dict, optional) – Dictionary of static parameters and their values to be applied to constraints.

  • constraint_section (str, optional) – The section to get the constraints from. Default is ‘constraint’.

Returns:

List of Constraint objects. Empty if no constraints were provided.

Return type:

list

pycbc.distributions.read_distributions_from_config(cp, section='prior')[source]

Returns a list of PyCBC distribution instances for a section in the given configuration file.

Parameters:
  • cp (WorflowConfigParser) – An open config file to read.

  • section ({"prior", string}) – Prefix on section names from which to retrieve the distributions.

Returns:

A list of the parsed distributions.

Return type:

list

pycbc.distributions.read_params_from_config(cp, prior_section='prior', vargs_section='variable_args', sargs_section='static_args')[source]

Loads static and variable parameters from a configuration file.

Parameters:
  • cp (WorkflowConfigParser) – An open config parser to read from.

  • prior_section (str, optional) – Check that priors exist in the given section. Default is ‘prior.’

  • vargs_section (str, optional) – The section to get the parameters that will be varied/need priors defined for them. Default is ‘variable_args’.

  • sargs_section (str, optional) – The section to get the parameters that will remain fixed. Default is ‘static_args’.

Returns:

  • variable_args (list) – The names of the parameters to vary in the PE run.

  • static_args (dict) – Dictionary of names -> values giving the parameters to keep fixed.