Grid
- class zdm.grid.Grid(survey, state, zDMgrid, zvals, dmvals, smear_mask, wdist=None, prev_grid=None)[source]
Bases:
object2D grid for computing FRB detection rates as a function of z and DM.
The Grid class is the core computational object in zdm. It builds a 2D probability distribution representing expected FRB detection rates across the redshift-DM plane for a given survey and parameter set.
Assumptions: - Each z bin represents FRBs originating at that redshift (not integrated) - Linear uniform spacing in both z and DM - DM includes cosmic + host contributions convolved together
- rates
2D array of expected FRB rates per (z, DM) bin.
- Type:
ndarray
- zvals
Redshift bin centers.
- Type:
ndarray
- dmvals
DM bin centers in pc/cm^3.
- Type:
ndarray
- state
Parameter state used for grid calculation.
- Type:
- survey
Associated survey object.
- Type:
Methods Summary
EF([alpha, bandwidth])Calculates the fluence--energy conversion factors as a function of redshift This does NOT account for the central frequency
GenMCFRB(Emax_boost)Generates a single FRB according to the grid distributions
GenMCSample(N[, Poisson])Generate a MC sample of FRB events
build_sz()calc_dV([reINIT])Calculates volume per steradian probed by a survey.
calc_pdv([beam_b, beam_o])Calculates the rate per cell. Assumed model: a power-law between Emin and Emax (erg) with slope gamma. Efficiencies: list of efficiency response to DM So-far: does NOT include time x solid-angle factor.
multiplies the rate per cell with the appropriate pdm plot
calc_thresholds(F0, eff_table[, bandwidth, ...])Sets the effective survey threshold on the zdm grid
check_grid([TOLERANCE])Check that the grid values are behaving as expected
chk_upd_param(param, vparams[, update])Check to see whether a parameter is differs from that in self.state
get_dm_coeffs(DMlist)Returns indices and coefficients for interpolating between DM values
get_p_zgdm(DMs)Calcuates the probability of redshift given a DM We already have our grid of observed DM values.
Function asking the grid to return the p(w) distribution.
Returns rates, multiplied by the relevant constant, and accounting for any DM preference via a DM mask
get_z_coeffs(zlist)Returns indices and coefficients for interpolating between z values
initMC()Initialises the MC sample, if it has not been doen already This uses a great deal of RAM - hence, do not do this lightly!
Set the luminsoity function for FRB energetics
load_grid(gridfile, zfile, dmfile)parse_grid(zDMgrid, zvals, dmvals)Scales volumetric rate by SFR
smear_dm(smear)Smears DM using the supplied array.
update(vparams[, ALL, prev_grid])Update the grid based on a set of input parameters
Methods Documentation
- EF(alpha=0, bandwidth=1000000000.0)[source]
Calculates the fluence–energy conversion factors as a function of redshift This does NOT account for the central frequency
- GenMCFRB(Emax_boost)[source]
Generates a single FRB according to the grid distributions
Samples beam position b, FRB DM, z, s=SNR/SNRth, and w Currently: no interpolation included. This should be implemented in s,DM, and z.
- NOTE: currently, the actual FRB widths are not part of ‘grid’
only the relative probabilities of any given width. Hence, this routine only returns the integer of the width bin not the width itelf.
- Parameters:
pwb (optional) – probability(width,beam)
Emax_boost (float, optional) – Allow for larger energies than Emax The factor is logarithmic, i.e. Emax_boost = 2. allows for 10**2 higher energies than Emax
- Returns:
FRBparams=[MCz,MCDM,MCb,j,MCs], pwb values These are:
MCz: redshift MCDM: dispersion measure (extragalactic) MCb: beam value j: MCs: SNR/SNRth value of FRB MCw: width value of FRB
[MCz, MCDM, MCb, j, MCs, MCw]
- Return type:
- GenMCSample(N, Poisson=False)[source]
Generate a MC sample of FRB events
If Poisson=True, then interpret N as a Poisson expectation value Otherwise, generate precisely N FRBs
Generated values are [MCz, MCDM, MCb, MCs, MCw] NOTE: the routine GenMCFRB does not know ‘w’, merely
which w bin it generates.
- calc_dV(reINIT=False)[source]
Calculates volume per steradian probed by a survey.
Does this only in the z-dimension (for obvious reasons!)
- calc_pdv(beam_b=None, beam_o=None)[source]
Calculates the rate per cell. Assumed model: a power-law between Emin and Emax (erg)
with slope gamma.
Efficiencies: list of efficiency response to DM So-far: does NOT include time x solid-angle factor
NOW: this includes a solid-angle and beam factor if initialised
This will recalculate beam factors if they are passed, however during iteration this is not recalculated
- calc_thresholds(F0: float, eff_table, bandwidth=1000000000.0, nuRef=1300000000.0, weights=None)[source]
Sets the effective survey threshold on the zdm grid
- Parameters:
F0 (float) – base survey threshold
eff_table ([type]) – table of efficiencies corresponding to DM-values. 1, 2, or 3 dimensions!
bandwidth ([type], optional) – [description]. Defaults to 1e9.
nuObs ([float], optional) – survey frequency (affects sensitivity via alpha - only for alpha method) Defaults to 1.3e9.
nuRef ([float], optional) – reference frequency we are calculating thresholds at Defaults to 1.3e9.
weights ([type], optional) – [description]. Defaults to None.
- Raises:
ValueError – [description]
ValueError – [description]
- check_grid(TOLERANCE=1e-06)[source]
Check that the grid values are behaving as expected
- TOLERANCE: defines the max relative difference in expected
and found values that will be tolerated
- chk_upd_param(param: str, vparams: dict, update=False)[source]
Check to see whether a parameter is differs from that in self.state
- get_dm_coeffs(DMlist)[source]
Returns indices and coefficients for interpolating between DM values
dmlist: np.ndarray of dispersion measures (extragalactic!)
- get_p_zgdm(DMs: ndarray)[source]
Calcuates the probability of redshift given a DM We already have our grid of observed DM values. Just take slices!
- Parameters:
DMs (np.ndarray) – array of DM values
- Returns:
array of priors for the DMs
- Return type:
np.ndarray
- get_pw_dist()[source]
Function asking the grid to return the p(w) distribution.
This will be an “all-burst” distribution in case of a repeater inherited class.
Note that a grid does not actually know what a “width” means: it is simply an abstract category of FRBs corresponding to a particular efficiency and fraction of the population.
Args: None
- Returns:
Rate per width bin Wzs (np.ndarray: Nw x Nz): Rate as a function of w and z Wdms (np.ndarray: Nw x Ndm): Rate as a function of w and DM
- Return type:
Wtots (np.ndarray)
- get_rates()[source]
Returns rates, multiplied by the relevant constant, and accounting for any DM preference via a DM mask
- get_z_coeffs(zlist)[source]
Returns indices and coefficients for interpolating between z values
zlist: np.ndarray of dispersion measures (extragalactic!)
- initMC()[source]
Initialises the MC sample, if it has not been doen already This uses a great deal of RAM - hence, do not do this lightly!
- smear_dm(smear: ndarray)[source]
Smears DM using the supplied array. Example use: DMX contribution
smear_grid is created in place
- Parameters:
smear (np.ndarray) – Smearing array
- update(vparams: dict, ALL=False, prev_grid=None)[source]
Update the grid based on a set of input parameters
Hierarchy: Each indent corresponds to one ‘level’. This is used in the program control below to dictate how far each tree should proceed in calculation. Direct variable inputs are always listed first We see that sfr evolution and dm smearing lie just before the pdv step Hence, we deal with these first, and calc rates as a final step regardless of what else has changed.
- Parameters:
- calc_rates:
- calc_pdv
Emin Emax gamma H0 calc_thresholds
F0 alpha bandwidth
- set_evolution
sfr_n H0
- smear_grid
grid mask
dmx_params (lmean, lsigma)
dV zdm_grid
H0
Note that the grid is independent of the constant C (trivial dependence)
- Parameters:
vparams (dict) – [description]
- __init__(survey, state, zDMgrid, zvals, dmvals, smear_mask, wdist=None, prev_grid=None)[source]
Initialize the Grid for a survey and parameter state.
- Parameters:
survey (survey.Survey) – Survey object with telescope properties and FRB data.
state (parameters.State) – Parameter state defining the model. Note: grids share the same State object, so modifications affect all grids.
zDMgrid (ndarray) – 2D array of p(DM|z) probabilities, shape (nz, ndm).
zvals (ndarray) – Redshift bin centers. Bins span [z - dz/2, z + dz/2].
dmvals (ndarray) – DM bin centers in pc/cm^3. Bins span [DM - dDM/2, DM + dDM/2].
smear_mask (ndarray) – 1D convolution kernel for host DM smearing.
wdist (bool, optional) – If True, include width distribution effects.
prev_grid (Grid, optional) – Another Grid with same z/DM values but different survey. Allows reusing pre-computed cosmological quantities.
- get_dm_coeffs(DMlist)[source]
Returns indices and coefficients for interpolating between DM values
dmlist: np.ndarray of dispersion measures (extragalactic!)
- get_z_coeffs(zlist)[source]
Returns indices and coefficients for interpolating between z values
zlist: np.ndarray of dispersion measures (extragalactic!)
- check_grid(TOLERANCE=1e-06)[source]
Check that the grid values are behaving as expected
- TOLERANCE: defines the max relative difference in expected
and found values that will be tolerated
- calc_dV(reINIT=False)[source]
Calculates volume per steradian probed by a survey.
Does this only in the z-dimension (for obvious reasons!)
- EF(alpha=0, bandwidth=1000000000.0)[source]
Calculates the fluence–energy conversion factors as a function of redshift This does NOT account for the central frequency
- calc_pdv(beam_b=None, beam_o=None)[source]
Calculates the rate per cell. Assumed model: a power-law between Emin and Emax (erg)
with slope gamma.
Efficiencies: list of efficiency response to DM So-far: does NOT include time x solid-angle factor
NOW: this includes a solid-angle and beam factor if initialised
This will recalculate beam factors if they are passed, however during iteration this is not recalculated
- get_pw_dist()[source]
Function asking the grid to return the p(w) distribution.
This will be an “all-burst” distribution in case of a repeater inherited class.
Note that a grid does not actually know what a “width” means: it is simply an abstract category of FRBs corresponding to a particular efficiency and fraction of the population.
Args: None
- Returns:
Rate per width bin Wzs (np.ndarray: Nw x Nz): Rate as a function of w and z Wdms (np.ndarray: Nw x Ndm): Rate as a function of w and DM
- Return type:
Wtots (np.ndarray)
- get_rates()[source]
Returns rates, multiplied by the relevant constant, and accounting for any DM preference via a DM mask
- calc_thresholds(F0: float, eff_table, bandwidth=1000000000.0, nuRef=1300000000.0, weights=None)[source]
Sets the effective survey threshold on the zdm grid
- Parameters:
F0 (float) – base survey threshold
eff_table ([type]) – table of efficiencies corresponding to DM-values. 1, 2, or 3 dimensions!
bandwidth ([type], optional) – [description]. Defaults to 1e9.
nuObs ([float], optional) – survey frequency (affects sensitivity via alpha - only for alpha method) Defaults to 1.3e9.
nuRef ([float], optional) – reference frequency we are calculating thresholds at Defaults to 1.3e9.
weights ([type], optional) – [description]. Defaults to None.
- Raises:
ValueError – [description]
ValueError – [description]
- smear_dm(smear: ndarray)[source]
Smears DM using the supplied array. Example use: DMX contribution
smear_grid is created in place
- Parameters:
smear (np.ndarray) – Smearing array
- get_p_zgdm(DMs: ndarray)[source]
Calcuates the probability of redshift given a DM We already have our grid of observed DM values. Just take slices!
- Parameters:
DMs (np.ndarray) – array of DM values
- Returns:
array of priors for the DMs
- Return type:
np.ndarray
- GenMCSample(N, Poisson=False)[source]
Generate a MC sample of FRB events
If Poisson=True, then interpret N as a Poisson expectation value Otherwise, generate precisely N FRBs
Generated values are [MCz, MCDM, MCb, MCs, MCw] NOTE: the routine GenMCFRB does not know ‘w’, merely
which w bin it generates.
- initMC()[source]
Initialises the MC sample, if it has not been doen already This uses a great deal of RAM - hence, do not do this lightly!
- GenMCFRB(Emax_boost)[source]
Generates a single FRB according to the grid distributions
Samples beam position b, FRB DM, z, s=SNR/SNRth, and w Currently: no interpolation included. This should be implemented in s,DM, and z.
- NOTE: currently, the actual FRB widths are not part of ‘grid’
only the relative probabilities of any given width. Hence, this routine only returns the integer of the width bin not the width itelf.
- Parameters:
pwb (optional) – probability(width,beam)
Emax_boost (float, optional) – Allow for larger energies than Emax The factor is logarithmic, i.e. Emax_boost = 2. allows for 10**2 higher energies than Emax
- Returns:
FRBparams=[MCz,MCDM,MCb,j,MCs], pwb values These are:
MCz: redshift MCDM: dispersion measure (extragalactic) MCb: beam value j: MCs: SNR/SNRth value of FRB MCw: width value of FRB
[MCz, MCDM, MCb, j, MCs, MCw]
- Return type:
- update(vparams: dict, ALL=False, prev_grid=None)[source]
Update the grid based on a set of input parameters
Hierarchy: Each indent corresponds to one ‘level’. This is used in the program control below to dictate how far each tree should proceed in calculation. Direct variable inputs are always listed first We see that sfr evolution and dm smearing lie just before the pdv step Hence, we deal with these first, and calc rates as a final step regardless of what else has changed.
- Parameters:
- calc_rates:
- calc_pdv
Emin Emax gamma H0 calc_thresholds
F0 alpha bandwidth
- set_evolution
sfr_n H0
- smear_grid
grid mask
dmx_params (lmean, lsigma)
dV zdm_grid
H0
Note that the grid is independent of the constant C (trivial dependence)
- Parameters:
vparams (dict) – [description]