Parameter Estimation (param_est) package
Submodules
simple_pe.param_est.metric module
- class simple_pe.param_est.metric.Metric(x, dx_directions, mismatch, f_low, psd, approximant='IMRPhenomD', tolerance=0.01, prob=None, snr=None)[source]
Bases:
object
A class to store the parameter space metric at a point x in parameter space, where the variations are given by dx_directions
- Parameters:
x – dictionary with parameter values for initial point
dx_directions – a list of directions to vary
mismatch – the mismatch value to use when calculating metric
f_low – low frequency cutoff
psd – the power spectrum to use in calculating the match
approximant – the approximant generator to use
tolerance – tolerance for scaling vectors
prob – probability to enclose within ellipse
snr – snr of signal, used in scaling mismatch
- calc_metric_error()[source]
We are looking for a metric and corresponding basis for which the basis is orthogonal and whose normalization is given by the desired mismatch This function checks for the largest error
- calculate_metric()[source]
Calculate the metric using the existing basis vectors, stored as self.dxs.
- generate_ellipse(npts=100, projected=False, scale=1.0)[source]
Generate an ellipse of points of the stored mismatch. Scale the radius by a factor `scale’ If the metric is projected, and we know the snr, then scale the mismatch appropriately for the number of projected dimensions. :param npts: number of points in the ellipse :param projected: use the projected metric if True, else use metric :param scale: scaling to apply to the mismatch :return ellipse_dict: SimplePESamples with ellipse of points
- generate_match_grid(npts=10, projected=False, scale=1.0)[source]
Generate a grid of points with extent governed by stored mismatch. If the metric is projected, scale to the projected number of dimensions Apply an overall scaling of scale
- Parameters:
npts – number of points in each dimension of the grid
projected – use the projected metric if True, else use metric
scale – a factor by which to scale the extent of the grid
- Return ellipse_dict:
SimplePESamples with grid of points and match at each point
- generate_posterior_grid(npts=10, projected=False, scale=None)[source]
Generate a grid of points with extent governed by requested mismatch
- Parameters:
npts – number of points in each dimension of the grid
projected – use the projected metric if True, else use metric
scale – the scale to apply to the greatest extent of the grid
- Return ellipse_dict:
SimplePESamples with grid of points and match at each point
- generate_samples(npts=100000, mins=None, maxs=None)[source]
Generate a given number of samples
- Parameters:
npts – number of points to generate
- Return phys_samples:
SimplePESamples with samples
- iteratively_update_metric(max_iter=20, verbose=True)[source]
A method to re-calculate the metric gij based on the matches obtained for the eigenvectors of the original metric :param max_iter: maximum number of iterations :param verbose: print information messages during update
- project_metric(projected_directions)[source]
Project out the unwanted directions of the metric :param projected_directions: list of parameters that we want to keep
- simple_pe.param_est.metric.average_mismatch(x, dx, scaling, f_low, psd, approximant='IMRPhenomD', verbose=False)[source]
This function calculates the average match for steps of +dx and -dx It also takes care of times when one of the steps moves beyond the edge of the physical parameter space
:param x : dictionary with the initial point :param dx: dictionary with parameter variations :param scaling: the scaling to apply to dx :param f_low: low frequency cutoff :param psd: the power spectrum to use in calculating the match :param approximant: the approximant generator to use :param verbose: print debugging information :return m: The average match from steps of +/- scaling * dx
- simple_pe.param_est.metric.find_metric_and_eigendirections(x, dx_directions, snr, f_low, psd, approximant='IMRPhenomD', tolerance=0.05, max_iter=20)[source]
Calculate the eigendirections in parameter space, normalized to enclose a 90% confidence region at the requested SNR
- Parameters:
x – dictionary with parameter values for initial point
dx_directions – list of parameters for which to calculate waveform variations
snr – the observed SNR
f_low – low frequency cutoff
psd – the power spectrum to use in calculating the match
approximant – the approximant to use
tolerance – the allowed error in the metric is (tolerance * mismatch)
max_iter – the maximum number of iterations
- Return g:
the mismatch metric in physical space
- simple_pe.param_est.metric.scale_dx(x, dx, desired_mismatch, f_low, psd, approximant='IMRPhenomD', tolerance=0.01)[source]
This function scales the input vectors so that the mismatch between a waveform at point x and one at x + v[i] is equal to the specified mismatch, up to the specified tolerance.
- Parameters:
x – dictionary with parameter values for initial point
dx – a dictionary with the initial change in parameters
desired_mismatch – the desired mismatch (1 - match)
f_low – low frequency cutoff
psd – the power spectrum to use in calculating the match
approximant – the approximant generator to use
tolerance – the maximum fractional error in the mismatch
- Return scale:
The required scaling of dx to achieve the desired mismatch
simple_pe.param_est.pe module
- class simple_pe.param_est.pe.SimplePESamples(*args, logger_warn='warn', autoscale=True)[source]
Bases:
SamplesDict
Class for holding Simple PE Samples, and generating PE distributions
- add_fixed(name, value)[source]
generate an additional parameter called ‘name’ with constant ‘value’
- Parameters:
name – the name of the parameter
value – its value
- calculate_hm_prec_probs(hm_snr=None, prec_snr=None, snr_2pol=None, overlaps=None)[source]
Calculate the precession SNR
- Parameters:
hm_snr – dictionary of measured SNRs in higher modes
prec_snr – measured precession SNR
snr_2pol – the SNR in the second polarization
overlaps – dictionary of the measured overlaps between modes
- calculate_rho_2nd_pol(a_net, net_snr)[source]
Calculate the SNR in the second polarization
- Parameters:
a_net – network sensitivity to x polarization (in DP frame)
net_snr – the network SNR
- calculate_rho_lm(psd, f_low, net_snr, modes, interp_directions, interp_points=5, approximant='IMRPhenomXPHM', alpha_lm_grid=None)[source]
Calculate the higher mode SNRs
- Parameters:
psd – the PSD to use
f_low – low frequency cutoff
net_snr – the network SNR
modes – modes for which to calculate SNR
interp_directions – directions to interpolate
interp_points – number of points to interpolate alpha_lm
approximant – waveform approximant
- calculate_rho_p(psd, f_low, net_snr, interp_directions, interp_points=5, approximant='IMRPhenomXP', beta_22_grid=None)[source]
Calculate the precession SNR
- Parameters:
psd – the PSD to use
f_low – low frequency cutoff
net_snr – the network SNR
interp_directions – directions to interpolate
interp_points – number of points to interpolate alpha_lm
approximant – waveform approximant
- generate_all_posterior_samples(function=None, **kwargs)[source]
Convert samples stored in the SamplesDict according to a conversion function
Parameters
- function: func, optional
function to use when converting posterior samples. Must take a dictionary as input and return a dictionary of converted posterior samples. Default simple_pe.param_est.convert
- **kwargs: dict, optional
All additional kwargs passed to function
- generate_chi_p(chi_p_dist='uniform', overwrite=False)[source]
generate chi_p points with the desired distribution and include in the existing samples dict
- Parameters:
chi_p_dist – the distribution to use for chi_p. Currently supports ‘uniform’ and ‘isotropic_on_sky’
overwrite – if True, then overwrite existing values, otherwise don’t
- generate_distance(fiducial_distance, fiducial_sigma, psd, f_low, interp_directions, interp_points=5, approximant='IMRPhenomXPHM', overwrite=False, sigma_22_grid=None)[source]
generate distance points using the existing theta_JN samples and fiducial distance. Interpolate sensitivity over the parameter space
- Parameters:
fiducial_distance – distance for a fiducial set of parameters
fiducial_sigma – the range for a fiducial set of parameters
psd – the PSD to use
f_low – low frequency cutoff
interp_directions – directions to interpolate
interp_points – number of points to interpolate alpha_lm
approximant – waveform approximant
overwrite – if True, then overwrite existing values,
otherwise don’t
- generate_prec_spin(overwrite=False)[source]
Generate component spins from chi_eff/chi_align/spin_z and chi_p
- Parameters:
overwrite – if True, then overwrite existing values, otherwise don’t
- generate_spin_z(overwrite=False)[source]
Generate z-component spins from chi_eff or chi_align
- Parameters:
overwrite – if True, then overwrite existing values, otherwise don’t
- generate_theta_jn(theta_dist='uniform', snr_left=0.0, snr_right=0.0, overwrite=False)[source]
generate theta JN points with the desired distribution and include in the SimplePESamples
- Parameters:
theta_dist – the distribution to use for theta. Currently, supports ‘uniform’, ‘left_circ’, ‘right_circ’, ‘left_right’
snr_left – left snr
snr_right – right snr
overwrite – if True, then overwrite existing values,
otherwise don’t
- jitter_distance(net_snr, response_sigma=0.0)[source]
jitter distance values based upon existing distances jitter due to SNR and variation of network response
- Parameters:
net_snr – the network SNR
response_sigma – standard deviation of network response (over sky)
- property neff
- trim_unphysical(maxs=None, mins=None, set_to_bounds=False)[source]
Trim unphysical points from SimplePESamples
- Parameters:
maxs – the maximum permitted values of the physical parameters
mins – the minimum physical values of the physical parameters
set_to_bounds – move points that lie outside physical space to the boundary of allowed space
- Return physical_samples:
SamplesDict with points outside the param max and min given
- simple_pe.param_est.pe.calculate_interpolated_snrs(samples, psd, f_low, dominant_snr, modes, response_sigma, fiducial_sigma, dist_interp_dirs, hm_interp_dirs, prec_interp_dirs, interp_points, approximant, localization_method, **kwargs)[source]
Wrapper function to calculate the SNR in the (l,m) multipoles, the SNR in the second polarisation and the SNR in precession.
Parameters
- samples: simple_pe.param_est.pe.SimplePESamples
table of posterior distributions
- psd: pycbc.types.frequencyseries
frequency series containing the PSD
- f_low: float
low frequency cut-off to use for SNR calculations
- dominant_snr: float
SNR in the dominant 22 multipole
- modes: list
list of higher order multipoles that you wish to calculate the SNR for
- alpha_net: float
network sensitivity to x polarization (in DP frame) used to calculate the SNR in the second
- response_sigma: float
standard deviation of network response over sky region
- fiducial_distance: float
distance at which a face on signal would give the observed dominant SNR
- fiducial_sigma: float
distance at which a face on signal would give SNR=8 at (using params for fiducial_distance)
- dist_interp_dirs: list
directions to interpolate the distance
- hm_interp_dirs: list
directions to interpolate the higher multipole SNR calculation
- prec_interp_dirs: list
directions to interpolate the precession SNR calculation
- interp_points: int
number of points to interpolate the SNRs
- approximant: str
approximant to use when calculating the SNRs
- localization_method: str
method to use when localizing the event. Must either be ‘average’ or ‘fullsky’
- simple_pe.param_est.pe.component_mass_prior_weight(samples, dx_directions)[source]
Re-weight points to uniform in mass ratio rather than symmetric mass ratio. Since the transformation is singular at equal mass we truncate at close to equal mass (eta = 0.2499)
Parameters
- samples: SimplePESamples
set of samples to re-weight
- dx_directions: list
directions used in generating samples
Returns
reweighted_samples: SimplePESamples
- simple_pe.param_est.pe.interpolate_alpha_lm(param_max, param_min, fixed_pars, psd, f_low, grid_points, modes, approximant)[source]
generate interpolating functions for the amplitudes of the lm multipoles
- Parameters:
param_max – A dictionary containing the maximum value of each parameter
param_min – A dictionary containing the maximum value of each parameter
fixed_pars – A dictionary containing values of fixed parameters
psd – the PSD to use
f_low – low frequency cutoff
grid_points – number of points to interpolate alpha_33 and beta
modes – waveform modes to calculate
approximant – waveform approximant
- Return alpha:
dictionary of alpha[lm] values interpolated across the grid
- Return pts:
set of points used in each direction
- simple_pe.param_est.pe.interpolate_opening(param_max, param_min, fixed_pars, psd, f_low, grid_points, approximant)[source]
generate interpolating functions for the amplitudes of the opening angle
- Parameters:
param_max – A dictionary containing the maximum value of each parameter
param_min – A dictionary containing the maximum value of each parameter
fixed_pars – the fixed parameters needed to generate the waveform
psd – the psd to use in calculating mean frequency, used for opening angle
f_low – the low frequency cutoff to use
grid_points – number of points to interpolate opening angle
approximant – the waveform approximant to use
- Return opening:
array of opening angle values interpolated across the grid
- Return pts:
set of points used in each direction
- simple_pe.param_est.pe.interpolate_sigma(param_max, param_min, fixed_pars, psd, f_low, grid_points, approximant)[source]
generate interpolating function for sigma
- Parameters:
param_max – A dictionary containing the maximum value of each parameter
param_min – A dictionary containing the maximum value of each parameter
fixed_pars – A dictionary containing values of fixed parameters
psd – the PSD to use
f_low – low frequency cutoff
grid_points – number of points to interpolate alpha_33 and beta
approximant – waveform approximant
- Return alpha:
dictionary of alpha[lm] values interpolated across the grid
- Return pts:
set of points used in each direction
- simple_pe.param_est.pe.isotropic_spin_prior_weight(samples, dx_directions)[source]
Re-weight points to prior proportional to chi_p (1 - chi_p), rather than uniform in chi_p or chi_p2. This prior approximately matches the one that arises from using uniform priors on spin magnitudes and orientations.
Parameters
- samples: SimplePESamples
set of samples to reweight
- dx_directions: list
directions used in generating samples
Returns
reweighted_samples: SimplePESamples
- simple_pe.param_est.pe.reweight_based_on_observed_snrs(samples, **kwargs)[source]
Resample a table of posterior distributions based on the observed SNR in the higher multipoles, precession and second polarisation.
Parameters
- samples: simple_pe.param_est.pe.SimplePESamples
table containing posterior distributions
- **kwargs: dict, optional
all kwargs passed to the samples.calculate_hm_prec_probs function
simple_pe.param_est.filter module
- simple_pe.param_est.filter.filter_grid_components(event_info, grid_data, strain_f, f_low, f_high, delta_f, psds, approximant, t_start, t_end, ifos, prec_snr=True, hm_snr=True)[source]
- simple_pe.param_est.filter.find_peak_snr(ifos, data, psds, t_start, t_end, x, dx_directions, f_low, approximant='IMRPhenomD', method='scipy', harm2=False, bounds=None, initial_mismatch=0.03, final_mismatch=0.001, tolerance=0.01, verbose=False, _net_snr=None)[source]
A function to find the maximum SNR. Either calculate a metric at the point x in dx_directions and walk to peak, or use scipy optimization regime
- Parameters:
ifos – list of ifos to use
data – a dictionary of data from the given ifos
psds – a dictionary of power spectra for the ifos
t_start – start time to consider SNR peak
t_end – end time to consider SNR peak
x – dictionary with parameter values for initial point
dx_directions – list of parameters for which to calculate waveform
variations :param f_low: low frequency cutoff :param approximant: the approximant to use :param method: how to find the maximum (either ‘scipy’ or ‘metric’) :param harm2: use SNR from second harmonic :param bounds: give initial bounds for the range of parameters to investigate :param initial_mismatch: the mismatch for calculating the metric :param final_mismatch: the mismatch required to stop iteration :param tolerance: the allowed error in the metric is (tolerance * mismatch) :param verbose: if True then print info :return x_prime: the point in the grid with the highest snr :return snr_peak: the SNR squared at this point
- simple_pe.param_est.filter.matched_filter_network(ifos, data, psds, t_start, t_end, h, f_low, dominant_mode=0)[source]
Find the maximum SNR in the network for a waveform h within a given time range
- Parameters:
ifos – list of ifos
data – a dictionary containing data from the ifos
psds – a dictionary containing psds from the given ifos
t_start – start time to consider SNR peak
t_end – end time to consider SNR peak
h – waveform (either a time series or dictionary of time series)
f_low – low frequency cutoff
dominant_mode – the dominant waveform mode
(if a dictionary was passed) :return snr: the network snr :return smax: the max snr in each ifo :return tmax: return the time of max snr in each ifo
simple_pe.param_est.result module
- class simple_pe.param_est.result.Result(f_low=None, psd=None, approximant=None, snr_threshold=4.0, localization_file=None, data_from_matched_filter={})[source]
Bases:
GWSingleAnalysisRead
- property alpha_net
- property distance_face_on
- property f_net
- generate_all_posterior_samples(function=None, **kwargs)[source]
Generate all posterior samples via the conversion module
Parameters
- **kwargs: dict
all kwargs passed to the conversion module
- generate_metric(metric_directions, template_parameters=None, dominant_snr=None, tolerance=0.01, max_iter=10)[source]
- generate_samples_from_aligned_spin_template_parameters(metric_directions, prec_interp_dirs, hm_interp_dirs, dist_interp_dirs, modes=['33'], alpha_net=None, interp_points=7, template_parameters=None, dominant_snr=None, reweight_to_isotropic_spin_prior=True, reweight_to_component_mass_prior=True, localization_method='fullsky', metric=None, neff=5000, nsamples=None)[source]
- generate_samples_from_metric(*args, npts=100000, metric=None, mins=None, maxs=None, **kwargs)[source]
- property left_snr
- property metric
- property overlaps
- property prec_snr
- property response_sigma
- property right_snr
- property samples_dict
- property sigma
- property snrs
- property template_parameters