Simulation grids (simqso.sqgrids)

Introduction

simqso.sqgrids provides classes for managing “points” (quasars) in a simulation. A collection of points is referred to as a grid, although they are not required to be uniformly distributed.

simqso.sqgrids.QsoSimObjects is the base container class for QSO grids. Objects can be distributed in a regular grid (simqso.sqgrids.QsoSimGrid), or not, e.g. according to a luminosity function (simqso.sqgrids.QsoSimPoints).

Each point is defined by a set of variables (simqso.sqgrids.QsoSimVar). The simplest descriptive variables are absolute magnitude (simqso.sqgrids.AbsMagVar), apparent magnitude (simqso.sqgrids.AppMagVar), and redshift (simqso.sqgrids.RedshiftVar). Variables that inherit from simqso.sqgrids.SpectralFeatureVar define a function simqso.sqgrids.SpectralFeatureVar.add_to_spec() that knows how to apply the variable to a spectrum. Variables can be multidimensional (simqso.sqgrids.MultiDimVar), e.g., simqso.sqgrids.BrokenPowerLawContinuumVar has an extra dimension to account for the set of power-law slopes defining the full continuum.

Each variable must know how to be randomly sampled to produce a given realization of the simulation. Thus they are instantiated with a simqso.sqgrids.Sampler instance, which defines how to produce n samples of the variable. Basic samplers include simqso.sqgrids.UniformSampler which uniformly distributes values between two bounds, and simqso.sqgrids.ConstSampler which sets all values to be a single constant. Many samplers inherit from simqso.sqgrids.CdfSampler, which means values are drawn from a cumulative distribution function. For example, simqso.sqgrids.GaussianSampler, simqso.sqgrids.ExponentialSampler, and simqso.sqgrids.LogNormalSampler all generate points sampled from the cdf of the eponymous function.

A simple tutorial for working with samplers, variables, and grids is provided in https://github.com/imcgreer/simqso/blob/master/examples/GridExamples.ipynb.

Reference/API

class simqso.sqgrids.Sampler(low, high)[source]

Bases: object

Base class for sampling one-dimensional values within a given bound.

Subclasses must define the sample() function.

Parameters:
low,high : float

Lower and upper bounds for the sampler.

sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.FixedSampler(vals)[source]

Bases: simqso.sqgrids.Sampler

Use a fixed set of values as the sample.

>>> from simqso.sqgrids import FixedSampler
>>> s = FixedSampler([1.,2.,3.])
>>> s(3)
[1.0, 2.0, 3.0]
sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.NullSampler[source]

Bases: simqso.sqgrids.Sampler

Special container for variables which are not sampled.

sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.IndexSampler[source]

Bases: simqso.sqgrids.Sampler

Special container for variables which need an index into the grid.

sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.RandomSubSampler(n)[source]

Bases: simqso.sqgrids.Sampler

sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.ConstSampler(*val)[source]

Bases: simqso.sqgrids.Sampler

Returns a constant for all samples.

>>> from simqso.sqgrids import ConstSampler
>>> s = ConstSampler(17)
>>> s(3)
array([17, 17, 17])
sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.UniformSampler(low, high)[source]

Bases: simqso.sqgrids.Sampler

Returns values uniformly sampled between low and high, inclusive.

>>> from simqso.sqgrids import UniformSampler
>>> s = UniformSampler(0,1)
>>> s(3)
array([ 0. ,  0.5,  1. ])
sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.CdfSampler(low, high)[source]

Bases: simqso.sqgrids.Sampler

Returns values sampled from a cumulative distribution function, within the bounds passed during instantiation.

Subclasses must implement the cdf(x) and ppf(x) functions.

Parameters:
low,high : float

Lower and upper bounds for the sampler.

sample(n, **kwargs)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.PowerLawSampler(low, high, a)[source]

Bases: simqso.sqgrids.CdfSampler

Returns values sampled from a power law distribution with index a.

Unlike scipy.stats.powerlaw, allows a<0, but then requires low>0 in that case.

Examples

>>> from simqso.sqgrids import PowerLawSampler
>>> s = PowerLawSampler(1,2,-2)
>>> s(3)
array([ 1.4537,  1.1208,  1.1691])
class simqso.sqgrids.GaussianSampler(mean, sigma, low=-inf, high=inf)[source]

Bases: simqso.sqgrids.CdfSampler

Returns values sampled from a Gaussian distibution N(mean,sigma).

Examples

>>> from simqso.sqgrids import GaussianSampler
>>> s = GaussianSampler(50.,10.)
>>> s(3)
array([ 50.07  ,  42.0223,  58.9512])
class simqso.sqgrids.ExponentialSampler(scale, low=0, high=inf)[source]

Bases: simqso.sqgrids.CdfSampler

Returns values sampled from an exponential distibution with a given scale parameter.

Examples

>>> from simqso.sqgrids import ExponentialSampler
>>> s = ExponentialSampler(0.1)
>>> s(3)
array([ 0.08072409,  0.45771082,  0.03769428])
class simqso.sqgrids.LinearTrendWithAsymScatterSampler(coeffs, pts, low=-inf, high=inf)[source]

Bases: simqso.sqgrids.Sampler

Returns values sampled from a set of linear trends that define the Gaussian mean and sigma at each point x.

Must be calibrated with a set of input points that define where to sample the linear trends.

class simqso.sqgrids.BaldwinEffectSampler(coeffs, absMag, x=None, low=-inf, high=inf)[source]

Bases: simqso.sqgrids.LinearTrendWithAsymScatterSampler

Uses LinearTrendWithAsymScatterSampler to implement the Baldwin Effect, by sampling from mean, upper, and lower log-linear trends as a function of absolute magnitude.

sample(n=None, ii=None)[source]

Return a set of n values obtained from the sampler.

class simqso.sqgrids.QsoSimVar(sampler, name=None, seed=None, meta=None)[source]

Bases: object

Base class for variables used to define points within simulation grid. Each variable must have a name and a Sampler instance for generating values of the variable.

Parameters:
sampler : simqso.sqgrids.Sampler instance
name : str

Unique name for variable.

seed : int

Seed to apply to RNG before sampling the variable. If None there is no call to random.seed().

resample(*args, **kwargs)[source]

Update the samplers of any dependent variables and then resample.

set_seed(seed, overwrite=False)[source]

Update the random seed used when sampling the variable. If overwrite is False an exisiting seed will be preserved.

updateMeta(meta, axPfx)[source]

Update the meta-data dictionary associated with the variable.

class simqso.sqgrids.MultiDimVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

Special case of QsoSimVar that handles multi-dimensional variables. The last dimension must be a sequence of Sampler instances, which can be nested in as many outer dimensions as necessary.

resample(*args, **kwargs)[source]

Update the samplers of any dependent variables and then resample.

class simqso.sqgrids.SpectralFeatureVar[source]

Bases: object

Mix-in class to define variables that act on spectra.

Subclasses must define the render() function.

add_to_spec(spec, par, **kwargs)[source]

Applies the variable to an input spectrum.

Parameters:
spec : simqso.sqbase.Spectrum instance
par : sampled values of the variable that are passed to render()
class simqso.sqgrids.AppMagVar(sampler, obsBand, **kwargs)[source]

Bases: simqso.sqgrids.QsoSimVar

An apparent magnitude variable, defined in an observed bandpass band.

updateMeta(meta, axPfx)[source]

Update the meta-data dictionary associated with the variable.

class simqso.sqgrids.AbsMagVar(sampler, restWave=None, **kwargs)[source]

Bases: simqso.sqgrids.QsoSimVar

An absolute magnitude variable, defined at rest-frame wavelength restWave in Angstroms.

updateMeta(meta, axPfx)[source]

Update the meta-data dictionary associated with the variable.

class simqso.sqgrids.RedshiftVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

A redshift variable.

class simqso.sqgrids.ContinuumVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar, simqso.sqgrids.SpectralFeatureVar

Base class for variables that define the quasar spectral continuum.

class simqso.sqgrids.BrokenPowerLawContinuumVar(samplers, breakPts, **kwargs)[source]

Bases: simqso.sqgrids.ContinuumVar, simqso.sqgrids.MultiDimVar

Representation of a quasar continuum as a series of broken power laws.

Parameters:
samplers : sequence of simqso.sqgrids.Sampler instances

Each sampler instance defines the power law spectral index at a given section of the continuum, as alpha_nu where f_nu = nu^alpha_nu.

breakPts : sequence of floats

Break wavelengths in Angstroms.

Examples

>>> from simqso.sqgrids import BrokenPowerLawContinuumVar,GaussianSampler
>>> v = BrokenPowerLawContinuumVar([GaussianSampler(-1.5,0.3),
                                    GaussianSampler(-0.5,0.3)],
                                   [1215.7])
>>> v(3)
array([[-1.801, -1.217],
       [-1.56 , -0.594],
       [-1.605, -0.248]])
render(wave, z, slopes, fluxNorm=None, assocvals=None)[source]

Renders the broken power law continuum at redshift z given the set of sampled slopes. Aribrarily normalized unless the fluxNorm parameter is supplied.

Parameters:
fluxNorm : dict
wavelength : float

rest-frame wavelength in Angstroms at which to normalize spectrum.

M_AB : float

absolute AB magnitude at wavelength

DM : function to return distance modulus, as in DM(z)

updateMeta(meta, axPfx)[source]

Update the meta-data dictionary associated with the variable.

class simqso.sqgrids.EmissionFeatureVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar, simqso.sqgrids.SpectralFeatureVar

Base class for variables that define quasar spectral emission features.

class simqso.sqgrids.GaussianEmissionLineVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.EmissionFeatureVar, simqso.sqgrids.MultiDimVar

A single Gaussian emission line. Must be instantiated with three samplers for the profile, namely (wavelength, equivalent width, sigma). All parameters are given in the rest-frame and in Angstroms.

Examples

>>> from simqso.sqgrids import GaussianEmissionLineVar
>>> v = GaussianEmissionLineVar([GaussianSampler(1215.7,0.1),GaussianSampler(100.,10.),GaussianSampler(10.,1.)])
>>> v(3)
array([[ 1215.645,   113.125,     9.099],
       [ 1215.987,   109.654,     9.312],
       [ 1215.74 ,   101.765,    10.822]])
class simqso.sqgrids.GaussianLineEqWidthVar(sampler, name, wave0, width0, log=False, **kwargs)[source]

Bases: simqso.sqgrids.EmissionFeatureVar

this is an arguably kludgy way of making it possible to include line EW as a variable in grids, by reducing the line to a single parameter

Parameters:
sampler : simqso.sqgrids.Sampler instance

Sampler for generating equivalent width values.

name : str

Name of emission line.

wave0,width0 : float

Fixed Gaussian parameters for the rest-frame wavelength and sigma in Angstroms. Only the equivalent width is sampled.

class simqso.sqgrids.GaussianEmissionLinesTemplateVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.EmissionFeatureVar, simqso.sqgrids.MultiDimVar

A multidimensional variable representing a template of Gaussian-profile emission lines.

class simqso.sqgrids.BossDr9EmissionLineTemplateVar(samplers, lineNames=None, **kwargs)[source]

Bases: simqso.sqgrids.GaussianEmissionLinesTemplateVar

Subclass of GaussianEmissionLinesTemplateVar that obtains log-linear trends for the emission lines from the BOSS DR9 model (Ross et al. 2013). TODO: this should really start with the file

class simqso.sqgrids.FeTemplateVar(feGrid=None, **kwargs)[source]

Bases: simqso.sqgrids.EmissionFeatureVar

Variable used to store an iron emission template, and then render it at an input redshift.

Since the template is fixed it uses a simqso.sqgrids.NullSampler instance internally.

class simqso.sqgrids.DustBlackbodyVar(*args, **kwargs)[source]

Bases: simqso.sqgrids.ContinuumVar, simqso.sqgrids.MultiDimVar

Variable used to represent a warm dust component as a single blackbody, treated as a continuum.

class simqso.sqgrids.SightlineVar(forest, losMap=None, **kwargs)[source]

Bases: simqso.sqgrids.QsoSimVar

Variable used to associate quasars with lines-of-sight.

Since the spectra are precomputed a simqso.sqgrids.IndexSampler instance is used internally to map the sightlines to individual spectra.

class simqso.sqgrids.HIAbsorptionVar(forest, losMap=None, **kwargs)[source]

Bases: simqso.sqgrids.SightlineVar, simqso.sqgrids.SpectralFeatureVar

Variable used to store IGM HI absorption spectra.

Since the spectra are precomputed a simqso.sqgrids.IndexSampler instance is used internally to map the forest sightlines to individual spectra.

add_to_spec(spec, sightLine, advance=True, **kwargs)[source]

Applies the variable to an input spectrum.

Parameters:
spec : simqso.sqbase.Spectrum instance
par : sampled values of the variable that are passed to render()
class simqso.sqgrids.DustExtinctionVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar, simqso.sqgrids.SpectralFeatureVar

Base class for dust extinction features. Dust curves are provided in the rest frame and convolved with input spectra.

add_to_spec(spec, ebv, **kwargs)[source]

Applies the variable to an input spectrum.

Parameters:
spec : simqso.sqbase.Spectrum instance
par : sampled values of the variable that are passed to render()
class simqso.sqgrids.SMCDustVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.DustExtinctionVar

SMC dust extinction curve from XXX.

class simqso.sqgrids.CalzettiDustVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.DustExtinctionVar

Calzetti XXX dust extinction curve for starburst galaxies.

class simqso.sqgrids.BlackHoleMassVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

A black hole mass variable, in units of log(Msun).

class simqso.sqgrids.EddingtonRatioVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

A dimensionless Eddington ratio variable, as lambda_edd = L/L_edd.

class simqso.sqgrids.AbsMagFromAppMagVar(appMag, z, kcorr, cosmo, restWave=None, **kwargs)[source]

Bases: simqso.sqgrids.AbsMagVar

A variable that provides a conversion from apparent magnitude to absolute magnitude.

Internally uses a simqso.sqgrids.FixedSampler instance after converting to absMag.

Parameters:
appMag : ndarray

Apparent magnitudes (usually from an AppMagVar).

m2M : function

Conversion from apparent to absolute mag, as m2M(z) = K(z) + DM(z)

restWave : float

Rest wavelength in Angstroms for the absolute magnitudes.

class simqso.sqgrids.AbsMagFromBHMassEddRatioVar(logBhMass, logEddRatio, restWave=None, **kwargs)[source]

Bases: simqso.sqgrids.AbsMagVar

A variable that provides a conversion from black hole mass and Eddington ratio to absolute magnitude.

Internally uses a simqso.sqgrids.FixedSampler instance after converting to absMag.

TODO: uses a fixed BC estimate, should be an input.

Parameters:
logBhMass : ndarray

Log of black hole mass in Msun. (e.g., from an BlackHoleMassVar).

logEddRatio : ndarray

Log of dimensionless Eddington ratio (e.g., from an EddingtonRatioVar).

restWave : float

Rest wavelength in Angstroms for the absolute magnitudes.

class simqso.sqgrids.TimeVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

Variable to associate instantaneous quasar spectrum with a epoch.

class simqso.sqgrids.DrwTimeSeriesVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.MultiDimVar

Variable that associates damped random walk (DRW) parameters \(\tau\) and \(\sigma\) with a quasar.

class simqso.sqgrids.SynMagVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

Container for synthetic magnitudes.

class simqso.sqgrids.SynFluxVar(sampler, name=None, seed=None, meta=None)[source]

Bases: simqso.sqgrids.QsoSimVar

Container for synthetic fluxes.

class simqso.sqgrids.QsoSimObjects(qsoVars=[], cosmo=None, units=None, seed=None)[source]

Bases: object

A collection of simulated quasar objects. Objects are defined by a set of variables (QsoSimVar). The values for the variables are maintained internally as an astropy.table.Table, which can be saved and restored.

Parameters:
qsoVars : list of QsoSimVar instances

Set of variables used to initialize the simulation grid.

cosmo : astropy.cosmology.FLRW instance

Cosmology used for the simulation.

units : str

One of “flux” or “luminosity”, XXX should be handled internally…

seed : int

Seed for RNG. Note the seed is only applied at initialization, individual variables can be seeded with QsoSimVar.__init___().

addVar(var, noVals=False)[source]

Add a variable to the simulation.

addVars(newVars, noVals=False)[source]

Add a list of variables to the simulation.

getVars(varType=<class 'simqso.sqgrids.QsoSimVar'>)[source]

Return all variables that are instances of varType. If varType is a string, return the variable with name varType.

read(gridFile, clean=False, extname=None)[source]

Read a simulation grid from a file.

write(outFn=None, simPars=None, outputDir='.', extname=None, overwrite=False)[source]

Write a simulation grid to a FITS file as a binary table, storing meta-data in the header.

class simqso.sqgrids.QsoSimPoints(qsoVars, n=None, **kwargs)[source]

Bases: simqso.sqgrids.QsoSimObjects

Simulation grid represented as a list of points.

Parameters:
qsoVars : list of QsoSimVar instances

Set of variables used to initialize the simulation grid.

n : int

Number of points in the grid. Not required (None) if the input variables already know how to sample the correct number of points (e.g., if they all use a FixedSampler).

class simqso.sqgrids.QsoSimGrid(*args, **kwargs)[source]

Bases: simqso.sqgrids.QsoSimObjects

Simulation grid represented as a uniform grid. Within each grid cell nPerBin objects are randomly sampled to fill the cell.

Parameters:
qsoVars : list of QsoSimVar instances

Set of variables used to initialize the simulation grid.

nBins : tuple

Number of bins along each grid axis (i.e., each variable).

nPerBin : int

Number of objects within each grid cell.

simqso.sqgrids.generateQlfPoints(qlf, mRange, zRange, kcorr, **kwargs)[source]

Generate a QsoSimPoints grid fed by AppMagVar and RedshiftVar instances which are sampled from an input luminosity function.

Parameters:
qlf : lumfun.LuminosityFunction instance

Representation of QLF to sample from.

mRange : sequence (m_min,m_max)

Range of apparent magnitudes to sample within [e.g., (17,22)].

zRange : sequence (z_min,z_max)

Range of redshifts to sample within [e.g., (2.2,3.5)].

kcorr : callable

K-correcton defined by either an sqbase.SimKCorr object or an equivalent callable that accepts as arguments (m,z,inverse=False) where m is apparent (absolute) magnitude if inverse is False (True), and z is redshift.

qlfseed : int

Seed for RNG applied before sampling from QLF.

gridseed : int

Seed for RNG applied after sampling from QLF, but before sampling additional variables added to the grid.

zin : sequence

Optional input redshifts. If supplied only the apparent magnitudes are sampled from the QLF.

kwargs : dict

Additional parameters passed to LuminosityFunction.sample_from_fluxrange().

Inheritance diagram of simqso.sqgrids