Skip to content
Snippets Groups Projects
Commit 15241888 authored by Moritz Huebner's avatar Moritz Huebner
Browse files

Merge branch 'gw-improvements' into 'master'

Improvements to bilby.gw

See merge request !948
parents da290af1 570ae9d8
No related branches found
No related tags found
1 merge request!948Improvements to bilby.gw
Pipeline #218624 passed with warnings
......@@ -156,8 +156,12 @@ class GravitationalWaveTransient(Likelihood):
priors['geocent_time'] = float(self.interferometers.start_time)
if self.jitter_time:
priors['time_jitter'] = Uniform(
minimum=- self._delta_tc / 2, maximum=self._delta_tc / 2,
boundary='periodic')
minimum=- self._delta_tc / 2,
maximum=self._delta_tc / 2,
boundary='periodic',
name="time_jitter",
latex_label="$t_j$"
)
self._marginalized_parameters.append('geocent_time')
elif self.jitter_time:
logger.debug(
......@@ -1353,7 +1357,8 @@ class ROQGravitationalWaveTransient(GravitationalWaveTransient):
np.vdot(np.abs(h_plus_quadratic + h_cross_quadratic)**2,
self.weights[interferometer.name + '_quadratic'])
complex_matched_filter_snr = d_inner_h / (optimal_snr_squared**0.5)
with np.errstate(invalid="ignore"):
complex_matched_filter_snr = d_inner_h / (optimal_snr_squared**0.5)
d_inner_h_squared_tc_array = None
return self._CalculatedSNRs(
......
......@@ -20,6 +20,7 @@ from .conversion import (
chirp_mass_and_mass_ratio_to_total_mass,
total_mass_and_mass_ratio_to_component_masses)
from .cosmology import get_cosmology
from .source import PARAMETER_SETS
DEFAULT_PRIOR_DIR = os.path.join(os.path.dirname(__file__), 'prior_files')
......@@ -369,11 +370,21 @@ class UniformInComponentsMassRatio(Prior):
def rescale(self, val):
self.test_valid_for_rescaling(val)
return self.icdf(val)
resc = self.icdf(val)
if resc.ndim == 0:
return resc.item()
else:
return resc
def prob(self, val):
in_prior = (val >= self.minimum) & (val <= self.maximum)
return (1. + val)**(2. / 5.) / (val**(6. / 5.)) / self.norm * in_prior
with np.errstate(invalid="ignore"):
prob = (1. + val)**(2. / 5.) / (val**(6. / 5.)) / self.norm * in_prior
return prob
def ln_prob(self, val):
with np.errstate(divide="ignore"):
return np.log(self.prob(val))
class AlignedSpin(Interped):
......@@ -488,7 +499,8 @@ class ConditionalChiInPlane(ConditionalBasePrior):
)
def ln_prob(self, val, **required_variables):
return np.log(self.prob(val, **required_variables))
with np.errstate(divide="ignore"):
return np.log(self.prob(val, **required_variables))
def cdf(self, val, **required_variables):
r"""
......@@ -539,9 +551,11 @@ class ConditionalChiInPlane(ConditionalBasePrior):
return chi_aligned * ((self._reference_maximum / chi_aligned) ** (2 * val) - 1) ** 0.5
def _condition_function(self, reference_params, **kwargs):
return dict(minimum=0, maximum=(
self._reference_maximum ** 2 - kwargs[self._required_variables[0]] ** 2
) ** 0.5)
with np.errstate(invalid="ignore"):
maximum = np.sqrt(
self._reference_maximum ** 2 - kwargs[self._required_variables[0]] ** 2
)
return dict(minimum=0, maximum=maximum)
def __repr__(self):
return Prior.__repr__(self)
......@@ -633,6 +647,43 @@ class CBCPriorDict(ConditionalPriorDict):
logger.warning("Unable to determine minimum component mass")
return None
def is_nonempty_intersection(self, pset):
""" Check if keys in self exist in the PARAMETER_SETS pset """
if len(PARAMETER_SETS[pset].intersection(self.non_fixed_keys)) > 0:
return True
else:
return False
@property
def spin(self):
""" Return true if priors include any spin parameters """
return self.is_nonempty_intersection("spin")
@property
def precession(self):
""" Return true if priors include any precession parameters """
return self.is_nonempty_intersection("precession_only")
@property
def intrinsic(self):
""" Return true if priors include any intrinsic parameters """
return self.is_nonempty_intersection("intrinsic")
@property
def extrinsic(self):
""" Return true if priors include any extrinsic parameters """
return self.is_nonempty_intersection("extrinsic")
@property
def mass(self):
""" Return true if priors include any mass parameters """
return self.is_nonempty_intersection("mass")
@property
def phase(self):
""" Return true if priors include phase parameters """
return self.is_nonempty_intersection("phase")
class BBHPriorDict(CBCPriorDict):
def __init__(self, dictionary=None, filename=None, aligned_spin=False,
......@@ -827,6 +878,11 @@ class BNSPriorDict(CBCPriorDict):
redundant = True
return redundant
@property
def tidal(self):
""" Return true if priors include phase parameters """
return self.is_nonempty_intersection("tidal")
Prior._default_latex_labels = {
'mass_1': '$m_1$',
......@@ -853,12 +909,16 @@ Prior._default_latex_labels = {
'psi': '$\psi$',
'phase': '$\phi$',
'geocent_time': '$t_c$',
'time_jitter': '$t_j$',
'lambda_1': '$\\Lambda_1$',
'lambda_2': '$\\Lambda_2$',
'lambda_tilde': '$\\tilde{\\Lambda}$',
'delta_lambda_tilde': '$\\delta\\tilde{\\Lambda}$',
'chi_1': '$\\chi_1$',
'chi_2': '$\\chi_2$'}
'chi_2': '$\\chi_2$',
'chi_1_in_plane': '$\\chi_{1, \perp}$',
'chi_2_in_plane': '$\\chi_{2, \perp}$',
}
class CalibrationPriorDict(PriorDict):
......
......@@ -590,3 +590,38 @@ def supernova_pca_model(
pc_coeff4 * pc4 + pc_coeff5 * pc5)
return {'plus': h_plus, 'cross': h_cross}
precession_only = {
"tilt_1", "tilt_2", "phi_12", "phi_jl", "chi_1_in_plane", "chi_2_in_plane",
}
spin = {
"a_1", "a_2", "tilt_1", "tilt_2", "phi_12", "phi_jl", "chi_1", "chi_2",
"chi_1_in_plane", "chi_2_in_plane",
}
mass = {
"chirp_mass", "mass_ratio", "total_mass", "mass_1", "mass_2",
"symmetric_mass_ratio",
}
primary_spin_and_q = {
"a_1", "chi_1", "mass_ratio"
}
tidal = {
"lambda_1", "lambda_2", "lambda_tilde", "delta_lambda_tilde"
}
phase = {
"phase", "delta_phase",
}
extrinsic = {
"azimuth", "zenith", "luminosity_distance", "psi", "theta_jn",
"cos_theta_jn", "geocent_time", "time_jitter", "ra", "dec",
"H1_time", "L1_time", "V1_time",
}
PARAMETER_SETS = dict(
spin=spin, mass=mass, phase=phase, extrinsic=extrinsic,
tidal=tidal, primary_spin_and_q=primary_spin_and_q,
intrinsic=spin.union(mass).union(phase).union(tidal),
precession_only=precession_only,
)
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment