Multiple detectors, software injection (circular polarisation)#
Here we compare lalpulsar_parameter_estimation_nested
with cwinpy
in the case of
simulated Gaussian noise from multiple detectors (H1, L1 and V1) containing a software injected signal with close-to circular polarisation. The parameters being
estimated are \(h_0\), \(\phi_0\), \(\psi\) and \(\cos{\iota}\), all with uniform priors.
The script for this comparison, using the dynesty nested sampling algorithm, is shown at the bottom of the page. It produces the following comparison data:
Method |
\(\ln{(Z)}\) |
\(\ln{(Z)}\) noise |
\(\ln{}\) Odds |
---|---|---|---|
|
488839.112 |
488320.478 |
518.634±0.109 |
|
488302.601 |
487788.855 |
513.746±0.190 |
|
488302.351 |
513.497 |
Method |
\(h_0\) |
\(\phi_0\) (rad) |
\(\psi\) (rad) |
\(\cos{\iota}\) |
---|---|---|---|---|
|
1.27±0.16×10-25 |
0.96±0.44 |
0.75±0.44 |
0.82±0.11 |
90% credible intervals |
[1.06, 1.56]×10-25 |
[0.23, 1.63] |
[0.08, 1.48] |
[0.64, 0.98] |
|
1.27±0.16×10-25 |
0.97±0.44 |
0.74±0.44 |
0.81±0.11 |
90% credible intervals |
[1.05, 1.55]×10-25 |
[0.24, 1.62] |
[0.09, 1.47] |
[0.64, 0.98] |
Method |
\(h_0\) |
\(\phi_0\) (rad) |
\(\psi\) (rad) |
\(\cos{\iota}\) |
\(\ln{(L)}\) max |
---|---|---|---|---|---|
|
1.25×10-25 |
1.11 |
0.60 |
0.83 |
488852.25 |
|
1.22×10-25 |
1.33 |
0.38 |
0.84 |
488315.67 |
#!/usr/bin/env python
"""
Compare cwinpy with lalpulsar_parameter_estimation_nested for noise-only
data for a multiple detectors (H1, L1 and V1) containing a software injection
with close-to-circular polarisation.
"""
import os
import subprocess as sp
import h5py
import matplotlib
import numpy as np
from astropy.utils.data import download_file
from bilby.core.prior import Uniform
from comparitors import comparisons
from lalinference import LALInferenceHDF5PosteriorSamplesDatasetName
from lalinference.io import read_samples
from matplotlib import pyplot as plt
from cwinpy import HeterodynedData
from cwinpy.pe import pe
from cwinpy.plot import Plot
matplotlib.use("Agg")
# URL for ephemeris files
DOWNLOAD_URL = "https://git.ligo.org/lscsoft/lalsuite/raw/master/lalpulsar/lib/{}"
# create a fake pulsar parameter file
parcontent = """\
PSRJ J0123+3456
RAJ 01:23:45.6789
DECJ 34:56:54.321
F0 567.89
F1 -1.2e-12
PEPOCH 56789
H0 1.1e-25
COSIOTA 0.95
PSI 0.5
PHI0 1.2
"""
injection_parameters = {}
injection_parameters["h0"] = 1.1e-25
injection_parameters["phi0"] = 1.2
injection_parameters["psi"] = 0.5
injection_parameters["cosiota"] = 0.95
label = "multi_detector_software_injection_circular"
outdir = "outputs"
if not os.path.isdir(outdir):
os.makedirs(outdir)
# add content to the par file
parfile = os.path.join(outdir, "{}.par".format(label))
with open(parfile, "w") as fp:
fp.write(parcontent)
# create some fake heterodyned data
detectors = ["H1", "L1", "V1"] # the detector to use
asd = 1e-24 # noise amplitude spectral density
times = np.linspace(1000000000.0, 1000086340.0, 1440) # times
het = {}
hetfiles = []
for detector in detectors:
het[detector] = HeterodynedData(
times=times,
par=parfile,
injpar=parfile,
inject=True,
fakeasd=asd,
detector=detector,
)
# output the data
hetfile = os.path.join(outdir, "{}_{}_data.txt".format(label, detector))
het[detector].write(hetfile)
hetfiles.append(hetfile)
# create priors
phi0range = [0.0, np.pi]
psirange = [0.0, np.pi / 2.0]
cosiotarange = [-1.0, 1.0]
h0range = [0.0, 1e-23]
# set prior for lalpulsar_parameter_estimation_nested
priorfile = os.path.join(outdir, "{}_prior.txt".format(label))
priorcontent = """H0 uniform {} {}
PHI0 uniform {} {}
PSI uniform {} {}
COSIOTA uniform {} {}
"""
with open(priorfile, "w") as fp:
fp.write(priorcontent.format(*(h0range + phi0range + psirange + cosiotarange)))
# set prior for bilby
priors = {}
priors["h0"] = Uniform(h0range[0], h0range[1], "h0", latex_label=r"$h_0$")
priors["phi0"] = Uniform(
phi0range[0], phi0range[1], "phi0", latex_label=r"$\phi_0$", unit="rad"
)
priors["psi"] = Uniform(
psirange[0], psirange[1], "psi", latex_label=r"$\psi$", unit="rad"
)
priors["cosiota"] = Uniform(
cosiotarange[0], cosiotarange[1], "cosiota", latex_label=r"$\cos{\iota}$"
)
# run lalpulsar_parameter_estimation_nested
try:
execpath = os.environ["CONDA_PREFIX"]
except KeyError:
raise KeyError(
"Please work in a conda environment with lalsuite and cwinpy installed"
)
execpath = os.path.join(execpath, "bin")
lppen = os.path.join(execpath, "lalpulsar_parameter_estimation_nested")
n2p = os.path.join(execpath, "lalinference_nest2pos")
Nlive = 1000 # number of nested sampling live points
Nmcmcinitial = 0 # set to 0 so that prior samples are not resampled
outfile = os.path.join(outdir, "{}_nest.hdf".format(label))
# set ephemeris files
efile = download_file(DOWNLOAD_URL.format("earth00-40-DE405.dat.gz"), cache=True)
sfile = download_file(DOWNLOAD_URL.format("sun00-40-DE405.dat.gz"), cache=True)
tfile = download_file(DOWNLOAD_URL.format("te405_2000-2040.dat.gz"), cache=True)
# set the command line arguments
runcmd = " ".join(
[
lppen,
"--verbose",
"--input-files",
",".join(hetfiles),
"--detectors",
",".join(detectors),
"--par-file",
parfile,
"--prior-file",
priorfile,
"--Nlive",
"{}".format(Nlive),
"--Nmcmcinitial",
"{}".format(Nmcmcinitial),
"--outfile",
outfile,
"--ephem-earth",
efile,
"--ephem-sun",
sfile,
"--ephem-timecorr",
tfile,
]
)
with sp.Popen(
runcmd,
stdout=sp.PIPE,
stderr=sp.PIPE,
shell=True,
bufsize=1,
universal_newlines=True,
) as p:
for line in p.stderr:
print(line, end="")
# convert nested samples to posterior samples
outpost = os.path.join(outdir, "{}_post.hdf".format(label))
runcmd = " ".join([n2p, "-p", outpost, outfile])
with sp.Popen(
runcmd,
stdout=sp.PIPE,
stderr=sp.PIPE,
shell=True,
bufsize=1,
universal_newlines=True,
) as p:
for line in p.stdout:
print(line, end="")
# get posterior samples
post = read_samples(outpost, tablename=LALInferenceHDF5PosteriorSamplesDatasetName)
lp = len(post["H0"])
postsamples = np.zeros((lp, len(priors)))
for i, p in enumerate(priors.keys()):
postsamples[:, i] = post[p.upper()]
# get evidence
hdf = h5py.File(outpost, "r")
a = hdf["lalinference"]["lalinference_nest"]
evsig = a.attrs["log_evidence"]
evnoise = a.attrs["log_noise_evidence"]
hdf.close()
# run bilby via the pe interface
runner = pe(
data_file=hetfiles,
par_file=parfile,
prior=priors,
detector=detectors,
outdir=outdir,
label=label,
)
result = runner.result
# evaluate the likelihood on a grid
gridpoints = 50
grid_size = dict()
for p in priors.keys():
grid_size[p] = np.linspace(
np.min(result.posterior[p]), np.max(result.posterior[p]), gridpoints
)
grunner = pe(
data_file=hetfiles,
par_file=parfile,
prior=priors,
detector=detectors,
outdir=outdir,
label=label,
grid=True,
grid_kwargs={"grid_size": grid_size},
)
grid = grunner.grid
# output comparisons
comparisons(label, outdir, grid, priors, cred=0.9)
# create results plot
allresults = {
"lalpulsar_parameter_estimation_nested": outpost,
"cwinpy_pe": result,
"cwinpy_pe (grid)": grid,
}
colors = {
key: plt.rcParams["axes.prop_cycle"].by_key()["color"][i]
for i, key in enumerate(allresults.keys())
}
plot = Plot(
results=allresults,
parameters=list(priors.keys()),
plottype="corner",
pulsar=parfile,
)
plot.plot(
bins=50,
smooth=0.9,
quantiles=[0.16, 0.84],
levels=(1 - np.exp(-0.5), 1 - np.exp(-2), 1 - np.exp(-9 / 2.0)),
fill_contours=True,
colors=colors,
)
plot.savefig(os.path.join(outdir, "{}_corner.png".format(label)), dpi=150)