cluster_slug: Bayesian Inference of Star Cluster Properties

The slugpy.cluster_slug module computes posterior probabilities for the mass, age, and extinction of star clusters from a set of input photometry. It is implemented as a wrapper around bayesphot: Bayesian Inference for Stochastic Stellar Populations, so for details on how the calculation is performed see the bayesphot documentation.

Getting the Default Library

The cluster_slug module requires a pre-computed library of slug simulations to use as a “training set” for its calculations. Due to its size, the default library is not included in the slug git repository. Instead, it is provided for download from the SLUG data products website. Download the two files clusterslug_mw_cluster_phot.fits and clusterslug_mw_cluster_prop.fits and save them in the cluster_slug directory of the main respository. If you do not do so, and do not provide your own library when you attempt to use cluster_slug, you will be prompted to download the default library.

Basic Usage

For an example of how to use cluster_slug, see the file cluster_slug/cluster_slug_example.py in the repository. All funtionality is provided through the cluster_slug class. The basic steps are as follows:

  1. Import the library and instantiate an sfr_slug object (see Full Documentation of slugpy.cluster_slug for full details):

    from slugpy.cluster_slug import cluster_slug
    cs = cluster_slug(photsystem=photsystem)
    

This creates a cluster_slug object, using the default simulation library. If you have another library of simulations you’d rather use, you can use the libname keyword to the cluster_slug constructor to select it. The optional argument photsystem specifies the photometric system you will be using for your data. Possible values are L_nu (flux per unit frequency, in erg/s/Hz), L_lambda (flux per unit wavelength, in erg/s/Angstrom), AB (AB magnitudes), STMAG (ST magnitudes), and Vega (Vega magnitudes). If left unspecified, the photometric system will be whatever the library was written in; the default library is in the L_nu system. Finally, if you have already read a library into memory using read_cluster, you can set the keyword lib in the cluster_slug constructor to specify that library should be used.

  1. Specify your filter(s), for example:

    cs.add_filters(['WFC3_UVIS_F336W', 'WFC3_UVIS_F438W', 'WFC3_UVIS_F555W',
                    'WFC3_UVIS_F814W', 'WFC3_UVIS_F657N'])
    

The add_filter method takes as an argument a string or list of strings specifying which filters were used for the observations you’re going to analyze. You can have more than one set of filters active at a time (just by calling add_filters more than once), and then specify which set of filters you’re using for any given calculation.

  1. Specify your priors, for example:

    # Set priors to be flat in log T and A_V, but vary with log M as
    # p(log M) ~ 1/M
    def priorfunc(physprop):
       # Note: physprop is an array of shape (N, 3) where physprop[:,0] =
       # log M, physprop[:,1] = log T, physprop[:,2] = A_V
       return 1.0/exp(physprop[:,0])
    cs.priors = prorfunc
    

The priors property specifies the assumed prior probability distribution on the physical properties of star clusters. It can be either None (in which case all simulations in the library are given equal prior probability), an array with as many elements as there are simulations in the library giving the prior for each one, or a callable that takes a vector of physical properties as input and returns the prior for it.

  1. Generate a marginal posterior probability distribuiton via:

    logm, pdf = cs.mpdf(idx, phot, photerr = photerr)
    

The first argument idx is an index for which posterior distribution should be computed – a value of 0 generates the posterior in log mass, a value of 1 generates the posterion on log age, and a value of generates the posterior in A_V. The second argument phot is an array giving the photometric values in the filters specified in step 2; make sure you’re using the same photometric system you used in step 1. For the array phot, the trailing dimension must match the number of filters, and the marginal posterior-finding exercise is repeated over every value in the leading dimensions. If you have added two or more filter sets, you need to specify which one you want to use via the filters keyword. The optional argument photerr can be used to provide errors on the photometric values. The shape rules on it are the same as on phot, and the two leading dimensions of the two arrays will be broadcast together using normal broadcasting rules.

The cluster_slug.mpdf method returns a tuple of two quantities. The first is a grid of values for log M, log T, or A_V, depending on the value of idx. The second is the posterior probability distribution at each value of of the grid. Posteriors are normalized to have unit integral. If the input consisted of multiple sets of photometric values, the output will contains marginal posterior probabilities for each input. The output grid will be created automatically be default, but all aspects of it (shape, size, placement of grid points) can be controlled by keywords – see Full Documentation of slugpy.cluster_slug.

Using cluster_slug in Parallel

The cluster_slug module has full support for threaded computation using the python multiprocessing module. This allows efficient use of multiple cores on a shared memory machine, without the need for every project to read a large simulation library or store it in memory. See Parallelism in bayesphot for full details on the recommended paradigm for parallel computing. The full list of thread-safe cluster_slug methods is:

  • cluster_slug.logL
  • cluster_slug.mpdf
  • cluster_slug.mcmc
  • cluster_slug.bestmatch
  • cluster_slug.make_approx_phot
  • cluster_slug.make_approx_phys
  • cluster_slug.squeeze_rep
  • cluster_slug.mpdf_approx

Making Your Own Library

You can generate your own library by running slug; you might want to do this, for example, to have a library that works at different metallicity or for a different set of stellar tracks. An example parameter file (the one that was used to generate the default clusterslug_mw library) is included in the cluster_slug directory. This file uses slug’s capability to pick the output time and the cluster mass from specified PDFs.

One subtle post-processing step you should take once you’ve generated your library is to read it in using slugpy – The Python Helper Library and then write the photometry back out using the slugpy.write_cluster_phot routine with the format set to fits2. This uses an alternative FITS format that is faster to search when you want to load only a few filters out of a large library. For large data sets, this can reduce cluster_slug load times by an order of magnitude. (To be precise: the default format for FITS outputs to put all filters into a single binary table HDU, while the fits2 format puts each filter in its own HDU. This puts all data for a single filter into a contiguous block, rather than all the data for a single cluster into a contiguous block, and is therefore faster to load when one wants to load the data filter by filter.)

Variable Mode IMF

If your library was run with variable IMF parameters, these can also be used in cluster_slug. When creating a cluster_slug object, you can pass the array vp_list as an argument. This list should have an element for each variable parameter in your library. Each element should then be either True or False depending on whether you wish to include this parameter in the analysis. For example, for a library with four variable parameters you could have:

vp_list=[True,False,True,True]

Full Documentation of slugpy.cluster_slug

class slugpy.cluster_slug.cluster_slug(libname=None, filters=None, photsystem=None, lib=None, bw_phys=0.1, bw_phot=None, ktype='gaussian', priors=None, sample_density=None, pobs=None, reltol=0.01, abstol=1e-08, leafsize=16, use_nebular=True, use_extinction=True, use_phot_as_phys=[], thread_safe=True, pruning=False, caching='none', vp_list=[])[source]

A class that can be used to estimate the PDF of star cluster physical properties from a set of input photometry in various bands.

Properties
priors : array, shape (N) | callable | None
prior probability on each data point; interpretation depends on the type passed; array, shape (N): values are interpreted as the prior probability of each data point; callable: the callable must take as an argument an array of shape (N, nphys), and return an array of shape (N) giving the prior probability at each data point; None: all data points have equal prior probability
abstol : float
absolute error tolerance for kernel density estimation
reltol : float
relative error tolerance for kernel density estimation
thread_safe : bool
if True, the computation routines will run in thread-safe mode, allowing use with multiprocessing; this incurs a small performance penalty
Methods
physprop() :
returns the list of physical properties available for inference
filters() :
returns list of filters available in the library
filtersets() :
return a list of the currently-loaded filter sets
filter_units() :
returns units for available filters
add_filters() :
adds a set of filters for use in parameter estimation
logL() :
compute log likelihood at a particular set of physical and photometric parameters
mpdf() :
computer marginal posterior probability distribution for one or more physical properties from a set of photometric measurements
mcmc() :
do MCMC estimation of the posterior PDF on a set of photometric measurments
bestmatch() :
find the simulations in the library that are the closest matches to the input photometry
make_approx_phot() :
given a set of physical properties, return a set of points that can be used for fast approximation of the corresponding photometric properties
make_approx_phys() :
given a set of photometric properties, return a set of points that can be used for fast approximation of the corresponding physical properties
draw_sample():
draw a random sample of cluster physical and photometric properties
draw_phot():
given a set of physical properties, return a randomly-selected set of photometric properties
__init__(libname=None, filters=None, photsystem=None, lib=None, bw_phys=0.1, bw_phot=None, ktype='gaussian', priors=None, sample_density=None, pobs=None, reltol=0.01, abstol=1e-08, leafsize=16, use_nebular=True, use_extinction=True, use_phot_as_phys=[], thread_safe=True, pruning=False, caching='none', vp_list=[])[source]

Initialize a cluster_slug object.

Parameters
libname : string
name of the SLUG model to load; if left as None, the default is $SLUG_DIR/cluster_slug/modp020_chabrier_MW
lib : object
a library read by the read_cluster function; if specified this overrides the libname option; the library must contain both physical properties and photometry, and must include filter data; if this is not None, then the photsystem keyword is ignored
filters : iterable of stringlike
list of filter names to be used for inferenence
photsystem : None or string
If photsystem is None, the library will be left in whatever photometric system was used to write it. Alternately, if it is a string, the data will be converted to the specified photometric system. Allowable values are ‘L_nu’, ‘L_lambda’, ‘AB’, ‘STMAG’, and ‘Vega’, corresponding to the options defined in the SLUG code. Once this is set, any subsequent photometric data input are assumed to be in the same photometric system.
bw_phys : ‘auto’ | float | array, shape (2) | array, shape (3)
bandwidth for the physical quantities in the kernel density estimation; if set to ‘auto’, the bandwidth will be estimated automatically; if set to a scalar quantity, this will be used for all physical quantities; if set to an array, the array must have 2 elements if use_extinction is False, or 3 if it is True
bw_phot : None | ‘auto’ | float | array
bandwidth for the photometric quantities; if set to None, defaults to 0.25 mag / 0.1 dex; if set to ‘auto’, bandwidth is estimated automatically; if set to a float, this bandwidth is used for all photometric dimensions; if set to an array, the array must have the same number of dimensions as len(filters)
ktype : string
type of kernel to be used in densty estimation; allowed values are ‘gaussian’ (default), ‘epanechnikov’, and ‘tophat’; only Gaussian can be used with error bars
priors : array, shape (N) | callable | None
prior probability on each data point; interpretation depends on the type passed; array, shape (N): values are interpreted as the prior probability of each data point; callable: the callable must take as an argument an array of shape (N, nphys), and return an array of shape (N) giving the prior probability at each data point; None: all data points have equal prior probability
pobs : array, shape (N) | callable | None
probability of being observed on each data point; interpretation depends on the type passed; array, shape (N): values are interpreted as the prior probability of each data point; callable: the callable must take as an argument an array of shape (N, nphys), and return an array of shape (N) giving the prior probability at each data point; None: all data points have equal prior probability
sample_density : array, shape (N) | callable | ‘auto’ | None
the density of the data samples at each data point; this need not match the prior density; interpretation depends on the type passed; array, shape (N): values are interpreted as the density of data sampling at each sample point; callable: the callable must take as an argument an array of shape (N, nphys), and return an array of shape (N) giving the sampling density at each point; ‘auto’: the sample density will be computed directly from the data set; note that this can be quite slow for large data sets, so it is preferable to specify this analytically if it is known; None: data are assumed to be uniformly sampled, or to be sampled as the default library is if libname is also None
reltol : float
relative error tolerance; errors on all returned probabilities p will satisfy either abs(p_est - p_true) <= reltol * p_est OR abs(p_est - p_true) <= abstol, where p_est is the returned estimate and p_true is the true value
abstol : float
absolute error tolerance; see above
leafsize : int
number of data points in each leaf of the KD tree
use_nebular : bool
if True, photometry including nebular emission will be used if available; if not, nebular emission will be omitted
use_extinction : bool
if True, photometry including extinction will be used; if not, it will be omitted, and in this case no results making use of the A_V dimension will be available
use_phot_as_phys : listlike
a list of photometric filters to be treated as “physical properties”, so that inference can be carried out on them from other photometric bands; for example, setting this to [ “Lbol”, “QH0” ] would make it possible to derive the marginal posterior probability of bolometric luminosity or hydrogen-ionizing luminosity frm the observed photometry
thread_safe : bool
if True, cluster_slug will make extra copies of internals as needed to ensure thread safety when the computation routines (logL, mpdf, mcmc, bestmatch, make_approx_phot, make_approx_phys, mpdf_approx) are used with multiprocessing; this incurs a minor performance penalty, and can be disabled by setting to False if the code will not be run with the multiprocessing module
pruning : bool
if True, the underlying bayesphot objects will be pruned of clusters for which pobs or priors are zero, speeding up evaluations; the current implementation is limited in that pruning for priors is only done when the cluster_slug object is instantiated, and pruning for pobs is only done when each filter set is added, so the list of pruned clusters cannot be modified later
caching : ‘aggressive’ | ‘lazy’ | ‘none’

strategy for caching subsets of the data with some dimensions marginalised out; behavior is as follows:

‘agressive’
on construction, store sorted data for fast calculation of 1D PDFs of variables by themselves, and 1D PDFs of all physical variables marginalised over all other physical variables; this significantly increases the memory footprint and construction time, but greatly speeds up subsequent evaluation of any of these quantities, and is generally the best choice for prudction work in parallel
‘lazy’
sorted data sets for fast computation are created and cached as needed; this will make the first computation of any marginal PDF slower, but speed up all subsequent ones, without imposing any extra time at initial construction; this mode is generally best for interactive work, but is not thread_safe = True; memory cost depends on how many different marginal PDF combinations are calculated, but is always less than aggressive
‘none’
no caching is performed automatically; the user may still manually cache data by calling the make_cache method
vp_list : list
A list with an element for each of the variable parameters in the data. An element is set to True if we wish to use that parameter here, or False if we do not.
Returns
Nothing
Raises
IOError, if the library cannot be found
__weakref__

list of weak references to the object (if defined)

add_filters(filters, bandwidth=None, pobs=None)[source]

Add a set of filters to use for cluster property estimation

Parameters
filters : iterable of stringlike
list of filter names to be used for inferenence
bandwidth : None | ‘auto’ | float | array
bandwidth for the photometric quantities; if set to None, the bandwidth is unchanged for an existing filter set, and for a newly-created one the default physical and photometric bandwidths are used; if set to ‘auto’, bandwidth is estimated automatically; if set to a float, this bandwidth is used for all physical and photometric dimensions; if set to an array, the array must have the same number of entries as nphys+len(filters)
pobs : array, shape (N) | callable | ‘equal’ | None
the probability that a particular object would be observed, which is used, like prior, to weight the library; interpretation depends on type. ‘equal’ means all objects are equally likely to be observed, array is an array giving the observation probability of each object in the library, and callable means must be a function that takes an array containing the photometry, of shape (N, nhpot), as an argument, and returns an array of shape (N) giving the probability of observation for that object. Finally, None leaves the observational probability unchanged
Returns
nothing
bestmatch(phot, photerr=None, nmatch=1, bandwidth_units=False, filters=None)[source]

Searches through the simulation library and returns the closest matches to an input set of photometry.

Parameters:
phot : arraylike, shape (nfilter) or (…, nfilter)
array giving the photometric values; for a multidimensional array, the operation is vectorized over the leading dimensions
photerr : arraylike, shape (nfilter) or (…, nfilter)
array giving photometric errors, which must have the same shape as phot; if this is not None, then distances will be measured in units of the photometric error if bandwidth_units is False, or in units of the bandwidth added in quadrature with the errors if it is True
nmatch : int
number of matches to return; returned matches will be ordered by distance from the input
bandwidth_units : bool
if False, distances are computed based on the logarithmic difference in luminosity; if True, they are measured in units of the bandwidth
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
matches : array, shape (…, nmatch, nphys + nfilter)
best matches to the input photometry; shape in the leading dimensions will be the same as for phot, and if nmatch == 1 then that dimension will be omitted; in the final dimension, the first 3 elements give log M, log T, and A_V, while the last nfilter give the photometric values; if created with use_extinct = False, the A_V dimension is omitted
dist : array, shape (…, nmatch)
distances between the matches and the input photometry
bestmatch_phys(phys, nmatch=1, bandwidth_units=False, filters=None)[source]

Searches through the simulation library and returns the closest matches to an input set of photometry.

Parameters:
phys : arraylike, shape (nphys) or (…, nphys)
array giving the physical values; for a multidimensional array, the operation is vectorized over the leading dimensions
nmatch : int
number of matches to return; returned matches will be ordered by distance from the input
bandwidth_units : bool
if False, distances are computed based on the logarithmic difference in physical properties; if True, they are measured in units of the bandwidth
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
matches : array, shape (…, nmatch, nphys + nfilter)
best matches to the input properties; shape in the leading dimensions will be the same as for phot, and if nmatch == 1 then that dimension will be omitted
dist : array, shape (…, nmatch)
distances between the matches and the input physical properties
clear_cache(margindims=None)[source]

This method deletes from the cache

Parameters
margindims : listlike of integers
list of marginalised dimensions that should be removed from the cache, in the same format as make_cache; if left as None, the cache is completely emptied
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns
Nothing
del_filters(filters)[source]

Remove a set of filters, freeing the memory associated with them. Note that this does not delete the underlying library data, just the data for the KD tree used internally.

Parameters
filters : iterable of stringlike
list of filter names
Returns
Nothing
Raises
KeyError if the input set of filters is not loaded
draw_phot(physprop, physidx=None, photerr=None, nsample=1, filters=None)[source]

Returns a randomly-drawn sample of clusters for a given filter set

Parameters:
physprop : arraylike
physical properties to be used; the final dimension of the input must have len(physidx) indices, or nphys indicates if physidx is None; if the input is a multidimensional array, the operation is vectorized over the leading dimensions physical properties
physidx : arraylike
indices of the physical quantities being constrained; if left as None, all physical properties are set, and physprop must have a trailing dimension of size equal to nphys; otherwise this must be an arraylike of <= nphys positive integers, each unique and in the range [0, nphys), specying which physical dimensions are constrained
photerr : arraylike, shape (nphot)
photometric errors to apply to the output photometry; these are added in quadrature with the kernel density estimation bandwidth
nsample : int
number of random samples to draw for each set of physical properties; must be positive
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
samples : array, shape (nsample, nphys+nphot)
random sample drawn from the kernel density object; for the final dimension in the output, the first nphys elements are the physical quantities, the next nphot are the photometric quantities
draw_sample(photerr=None, nsample=1, filters=None)[source]

Returns a randomly-drawn sample of clusters for a given filter set

Parameters:
photerr : arraylike, shape (nphot)
photometric errors to apply to the output photometry; these are added in quadrature with the kernel density estimation bandwidth
nsample : int
number of random samples to draw for each set of physical properties; must be positive
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
samples : array, shape (nsample, nphys+nphot)
random sample drawn from the kernel density object; for the final dimension in the output, the first nphys elements are the physical quantities, the next nphot are the photometric quantities
filter_units()[source]

Returns list of all available filter units

Parameters:
None
Returns:
units : list of strings
list of available filter units
filters()[source]

Returns list of all available filters

Parameters:
None
Returns:
filters : list of strings
list of available filter names
filtersets()[source]

Returns list of all currently-loaded filter sets

Parameters:
None
Returns:
filtersets : list of list of strings
list of currently-loaded filter sets
load_data(filter_name, bandwidth=None, force_reload=False)[source]

Loads photometric data for the specified filter into memory

Parameters:
filter_name : string
name of filter to load
bandwidth : float
default bandwidth for this filter
force_reload : bool
if True, reinitialize the data even if has already been stored
Returns:
None
Raises:
ValueError, if filter_name is not one of the available filters
logL(physprop, photprop, photerr=None, filters=None, margindim=None)[source]

This function returns the natural log of the likelihood function evaluated at a particular log mass, log age, extinction, and set of log luminosities

Parameters:
physprop : arraylike, shape (nhpys) or (…, nphys)
array giving values of the log M, log T, and A_V; for a multidimensional array, the operation is vectorized over the leading dimensions; if created with use_extinct = False, the A_V dimension should be omitted. Will also include variable any variable parameters VPx if they are requested.
photprop : arraylike, shape (nfilter) or (…, nfilter)
array giving the photometric values; for a multidimensional array, the operation is vectorized over the leading dimensions
photerr : arraylike, shape (nfilter) or (…, nfilter)
array giving photometric errors; for a multidimensional array, the operation is vectorized over the leading dimensions
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
margindim : int | arraylike of ints | None
The index or indices of the physical or photometric properties to be maginalized over, numbered from 0 - nphys-1 for physical properties and from nphys - nfilter + nfilter - 1 for photometric properties. If this keyword is set, then physprop and/or photprop should have fewer than nphys or nphot elements due to the omission of marginalised dimensions. If all physical or photometric dimensions are marginalised out, that corresponding argument for physprop or photprop should be set to None
Returns:
logL : float or arraylike
natural log of the likelihood function
make_approx_phot(phys, squeeze=True, filter_ignore=None, filters=None)[source]

Returns an object that can be used for a fast approximation of the PDF of photometric properties that corresponds to a set of physical properties. The PDF produced by summing over the points returned is guaranteed to account for at least 1-reltol of the marginal photometric probability, and to represent the shape of the PDF in photometric space within a local accuracy of reltol as well.

Parameters:
phys : arraylike, shape (nphys) or (N, nphys)
the set or sets of physical properties for which the approximation is to be generated
squeeze : bool
if True, the representation returned will be squeezed to minimize the number of points included, using reltol as the error tolerance
filter_ignore : None or listlike of bool
if None, the kernel density representation returned covers all filters; otherwise this must be a listlike of bool, one entry per filter, with a value of False indicating that filter should be excluded from the values returned; suppressing filters can allow for more efficient representations
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
x : array, shape (M, nphot), or a list of such arrays
an array containing the list of points to be used for the approximation
wgts : array, shape (M), or a list of such arrays
an array containing the weights of the points
Notes:
if the requested relative tolerance cannot be reached for numerical reasons (usually because the input point is too far from the library to allow accurate computation), x and wgts will be return as None, and a warning will be issued
make_approx_phys(phot, photerr=None, squeeze=True, phys_ignore=None, filters=None, tol=None)[source]

Returns an object that can be used for a fast approximation of the PDF of physical properties that corresponds to a set of photometric properties. The PDF produced by summing over the points returned is guaranteed to account for at least 1-reltol of the marginal photometric probability, and to represent the shape of the PDF in photometric space within a local accuracy of reltol as well.

Parameters:
phot : arraylike, shape (nfilter) or (N, nfilter)
the set or sets of photometric properties for which the approximation is to be generated
photerr : arraylike, shape (nfilter) or (N, nfilter)
array giving photometric errors; the number of elements in the output lists will be the size that results from broadcasting together the leading dimensions of phot and photerr
squeeze : bool
if True, the representation returned will be squeezed to minimize the number of points included, using reltol as the error tolerance
phys_ignore : None or listlike of bool
if None, the kernel density representation returned covers all physical properties; otherwise this must be a listlike of bool, one entry per physical dimension, with a value of False indicating that dimension should be excluded from the values returned; suppressing dimensions can allow for more efficient representations
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
tol : float
if set, this tolerance overrides the value of reltol
Returns:
x : array, shape (M, nphys), or a list of such arrays
an array containing the list of points to be used for the approximation, where nphys is the number of physical dimensions being returned
wgts : array, shape (M), or a list of such arrays
an array containing the weights of the points
Notes:
if the requested relative tolerance cannot be reached for numerical reasons (usually because the input point is too far from the library to allow accurate computation), x and wgts will be return as None, and a warning will be issued
make_cache(margindims, filters=None)[source]

This method builds a cache to do faster calculation of PDFs where certain dimensions are marginalised out. If such caches exist, they are used automatically by all the computation methods.

Parameters
margindims : listlike of integers
list of dimensions to be marginalised out; physical dimensions go from 0 - nphys-1, photometric dimensions from nphys to nphys + nphot - 1; note that the indexing depends on the filter set specified by filters
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns
Nothing
mcmc(photprop, photerr=None, mc_walkers=100, mc_steps=500, mc_burn_in=50, filters=None)[source]

This function returns a sample of MCMC walkers for cluster physical properties

Parameters:
photprop : arraylike, shape (nfilter) or (…, nfilter)
array giving the photometric values; for a multidimensional array, the operation is vectorized over the leading dimensions
photerr : arraylike, shape (nfilter) or (…, nfilter)
array giving photometric errors; for a multidimensional array, the operation is vectorized over the leading dimensions
mc_walkers : int
number of walkers to use in the MCMC
mc_steps : int
number of steps in the MCMC
mc_burn_in : int
number of steps to consider “burn-in” and discard
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns
samples : array
array of sample points returned by the MCMC
mpdf(idx, photprop, photerr=None, ngrid=128, qmin=None, qmax=None, grid=None, norm=True, filters=None)[source]

Returns the marginal probability for one or mode physical quantities for one or more input sets of photometric properties. Output quantities are computed on a grid of values, in the same style as meshgrid

Parameters:
idx : int or listlike containing ints
index of the physical quantity whose PDF is to be computed; 0 = log M, 1 = log T, 2 = A_V, (2 or 3)+x = VPx; if this is an iterable, the joint distribution of the indicated quantities is returned
photprop : arraylike, shape (nfilter) or (…, nfilter)
array giving the photometric values; for a multidimensional array, the operation is vectorized over the leading dimensions
photerr : arraylike, shape (nfilter) or (…, nfilter)
array giving photometric errors; for a multidimensional array, the operation is vectorized over the leading dimensions
ngrid : int or listlike containing ints
number of points in each dimension of the output grid; if this is an iterable, it must have the same number of elements as idx
qmin : float or listlike
minimum value in the output grid in each quantity; if left as None, defaults to the minimum value in the library; if this is an iterable, it must contain the same number of elements as idx
qmax : float or listlike
maximum value in the output grid in each quantity; if left as None, defaults to the maximum value in the library; if this is an iterable, it must contain the same number of elements as idx
grid : listlike of arrays
set of values defining the grid on which the PDF is to be evaluated, in the same format used by meshgrid
norm : bool
if True, returned pdf’s will be normalized to integrate to 1
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
grid_out : array
array of values at which the PDF is evaluated; contents are the same as returned by meshgrid
pdf : array
array of marginal posterior probabilities at each point of the output grid, for each input cluster; the leading dimensions match the leading dimensions produced by broadcasting the leading dimensions of photprop and photerr together, while the trailing dimensions match the dimensions of the output grid
mpdf_approx(x, wgts, dims='phys', dims_return=None, ngrid=64, qmin='all', qmax='all', grid=None, norm=True, filters=None)[source]

Returns the marginal posterior PDF computed from a kernel density approximation returned by make_approx_phys or make_approx_phot. Outputs are computed on a grid of values, in the same style as meshgrid.

Parameters:
x : array, shape (M, ndim), or a list of such arrays
array of points retured by make_approx_phot or make_approx_phys
wgts : array, shape (M) or a list of such arrays
array of weights returned by make_approx_phot or make_approx_phys
dims : ‘phys’ | ‘phot’ | arraylike of ints
dimensions covered by x and wgts; the strings ‘phys’ or ‘phot’ indicate that they cover all physical or photometric dimensions, and correspond to the defaults returned by make_approx_phys and make_approx_phot, respectively; if dims is an array of ints, these specify the dimensions covered by x and wgts, where the physical dimensions are numbered 0, 1, … nphys-1, and the photometric ones are nphys, nphys+1, … nphys+nphot-1
dims_return : None or arraylike of ints
if None, the output PDF has the same dimensions as specified in dims; if not, then dims_return must be a subset of dims, and a marginal PDF in certain dimensions will be generated
ngrid : int or listlike containing ints
number of points in each dimension of the output grid; if this is an iterable, it must have the same number of elements as idx
qmin : float | listlike | ‘zoom’ | ‘all’
minimum value in the output grid in each quantity; if this a float, it is applied to each dimension; if it is an iterable, it must contain the same number of elements as the number of dimensions being returned, as gives the minimum in each dimension; if it is ‘zoom’ or ‘all’, the minimum is chosen automatically, with ‘zoom’ focusing on a region encompassing the probability maximum, and ‘all’ encompassing all the points in the representation
qmax : float | listlike | ‘zoom’ | ‘all’
same as qmin, but for the maximum of the output grid
grid : listlike of arrays
set of values defining the grid on which the PDF is to be evaluated, in the same format used by meshgrid
norm : bool
if True, returned pdf’s will be normalized to integrate to 1
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
grid_out : array
array of values at which the PDF is evaluated; contents are the same as returned by meshgrid
pdf : array
array of marginal posterior probabilities at each point of the output grid, for each input cluster; the leading dimensions match the leading dimensions produced by broadcasting the leading dimensions of photprop and photerr together, while the trailing dimensions match the dimensions of the output grid
mpdf_gen(fixeddim, fixedprop, margindim, ngrid=128, qmin=None, qmax=None, grid=None, norm=True, filters=None)[source]

Returns the marginal probability for one or more physical or photometric properties, keeping other properties fixed and marginalizing over other quantities. This is the most general marginal PDF routine provided.

Parameters:
fixeddim : int | arraylike of ints | None
The index or indices of the physical or photometric properties to be held fixed; physical properties are numbered 0 … nphys-1, and photometric ones are numbered nphys … nphys + nphot - 1. This can also be set to None, in which case no properties are held fixed.
fixedprop : array | None
The values of the properties being held fixed; the size of the final dimension must be equal to the number of elements in fixeddim, and if fixeddim is None, this must be too
margindim : int | arraylike of ints | None
The index or indices of the physical or photometric properties to be maginalized over, numbered in the same way as with fixeddim; if set to None, no marginalization is performed
ngrid : int or listlike containing ints
number of points in each dimension of the output grid; if this is an iterable, it must have nphys + nphot - len(fixeddim) - len(margindim) elements
qmin : float | arraylike
minimum value in the output grid in each quantity; if left as None, defaults to the minimum value in the library; if this is an iterable, it must contain a number of elements equal to nphys + nphot - len(fixeddim) - len(margindim)
qmax : float | arraylike
maximum value in the output grid in each quantity; if left as None, defaults to the maximum value in the library; if this is an iterable, it must have the same number of elements as qmin
grid : listlike of arrays
set of values defining the grid on which the PDF is to be evaluated, in the same format used by meshgrid
norm : bool
if True, returned pdf’s will be normalized to integrate to 1
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
grid_out : array
array of values at which the PDF is evaluated; contents are the same as returned by meshgrid
pdf : array
array of marginal posterior probabilities at each point of the output grid, for each input set of properties; the leading dimensions match the leading dimensions produced by broadcasting the leading dimensions of photprop and photerr together, while the trailing dimensions match the dimensions of the output grid
mpdf_phot(idx, physprop, ngrid=128, qmin=None, qmax=None, grid=None, norm=True, filters=None)[source]

Returns the marginal probability for one or more photometric quantities corresponding to an input set or distribution of physical properties. Output quantities are computed on a grid of values, in the same style as meshgrid.

Parameters:
idx : int or listlike containing ints
index of the photometric quantity whose PDF is to be computed, starting at 0; indices correspond to the order of elements in the filters argument; if this is an iterable, the joint distribution of the indicated quantities is returned
physprop : arraylike, shape (nphys) or (…, nphys)
physical properties to be used; if this is an array of nphys elements, these give the physical properties; if it is a multidimensional array, the operation is vectoried over the leading dimensions physical properties – the function must take an array of (nphys) elements as an input, and return a floating point value representing the PDF evaluated at that set of physical properties as an output
ngrid : int or listlike containing ints
number of points in each dimension of the output grid; if this is an iterable, it must have the same number of elements as idx
qmin : float or listlike
minimum value in the output grid in each quantity; if left as None, defaults to the minimum value in the library; if this is an iterable, it must contain the same number of elements as idx
qmax : float or listlike
maximum value in the output grid in each quantity; if left as None, defaults to the maximum value in the library; if this is an iterable, it must contain the same number of elements as idx
grid : listlike of arrays
set of values defining the grid on which the PDF is to be evaluated, in the same format used by meshgrid
norm : bool
if True, returned pdf’s will be normalized to integrate to 1
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
grid_out : array
array of values at which the PDF is evaluated; contents are the same as returned by meshgrid
pdf : array
array of marginal posterior probabilities at each point of the output grid, for each input set of properties; the leading dimensions match the leading dimensions produced by broadcasting the leading dimensions of photprop and photerr together, while the trailing dimensions match the dimensions of the output grid
physprop()[source]

Returns list of all physical properties available for inference

Parameters:
None
Returns:
physprop : list of string
list of available physical properties
squeeze_rep(x, wgts, dims=None, filters=None)[source]

Takes an input array of positions and weights that form a kernel density representation and approximates them using fewer points, using an error tolerance of reltol

Parameters:
x : array, shape (N, ndim)
an array of points forming a kernel density representation; on exit, x will be resized to (M, ndim) with M <= N
wgts : array, shape (N)
an array of weights for the kernel density representation; on exit, wgts will be resized to (M), with M <= N
dims : array, shape (ndim)
array specifying which dimensions in the kernel density representation the coordinates in x correspond to; if left as None, they are assumed to correspond to the first ndim dimensions in the data set
filters : listlike of strings
list of photometric filters to use; if left as None, and only 1 set of photometric filters has been defined for the cluster_slug object, that set will be used by default
Returns:
Nothing