Skip to content
Snippets Groups Projects
Commit 6f693479 authored by Gregory Ashton's avatar Gregory Ashton
Browse files

Merge branch 'fixing-GW150914-advanced-example' into 'master'

Fixing gw150914 advanced example

Closes #422

See merge request !626
parents 67571bf5 51003802
No related branches found
No related tags found
1 merge request!626Fixing gw150914 advanced example
Pipeline #85373 passed
...@@ -14,3 +14,4 @@ ra = Uniform(name='ra', minimum=0, maximum=2 * np.pi, boundary='periodic') ...@@ -14,3 +14,4 @@ ra = Uniform(name='ra', minimum=0, maximum=2 * np.pi, boundary='periodic')
theta_jn = Sine(name='theta_jn', boundary='reflective') theta_jn = Sine(name='theta_jn', boundary='reflective')
psi = Uniform(name='psi', minimum=0, maximum=np.pi, boundary='periodic') psi = Uniform(name='psi', minimum=0, maximum=np.pi, boundary='periodic')
phase = Uniform(name='phase', minimum=0, maximum=2 * np.pi, boundary='periodic') phase = Uniform(name='phase', minimum=0, maximum=2 * np.pi, boundary='periodic')
geocent_time = Uniform(name='geocent_time', minimum =1126259460.4, maximum=1126259464.4)
...@@ -12,66 +12,83 @@ List of events in BILBY dict: run -> help(bilby.gw.utils.get_event_time(event)) ...@@ -12,66 +12,83 @@ List of events in BILBY dict: run -> help(bilby.gw.utils.get_event_time(event))
""" """
from __future__ import division, print_function from __future__ import division, print_function
import bilby import bilby
from gwpy.timeseries import TimeSeries
outdir = 'outdir' outdir = 'outdir'
label = 'GW150914' label = 'GW150914'
logger = bilby.core.utils.logger
time_of_event = bilby.gw.utils.get_event_time(label) time_of_event = bilby.gw.utils.get_event_time(label)
bilby.core.utils.setup_logger(outdir=outdir, label=label) bilby.core.utils.setup_logger(outdir=outdir, label=label)
# GET DATA FROM INTERFEROMETER # GET DATA FROM INTERFEROMETER
interferometer_names = ['H1', 'L1'] # include 'V1' for appropriate O2 events interferometer_names = ['H1', 'L1'] # include 'V1' for appropriate O2 events
duration = 4 # length of data segment containing the signal duration = 4 # length of data segment containing the signal
roll_off = 0.2 # smoothness of transition from no signal post_trigger_duration = 2 # time between trigger time and end of segment
# to max signal in a Tukey Window. end_time = time_of_event + post_trigger_duration
psd_offset = -1024 # PSD is estimated using data from start_time = end_time - duration
# 'center_time + psd_offset' to 'center_time + psd_offset + psd_duration'.
# This determines the time window used to fetch open data. roll_off = 0.4 # smoothness in a tukey window, default is 0.4s
psd_duration = 100 # This determines the time window used to fetch open data
psd_duration = 32 * duration
psd_start_time = start_time - psd_duration
psd_end_time = start_time
filter_freq = None # low pass filter frequency to cut signal content above filter_freq = None # low pass filter frequency to cut signal content above
# Nyquist frequency. The condition is 2 * filter_freq >= sampling_frequency # Nyquist frequency. The condition is 2 * filter_freq >= sampling_frequency
# Keyword args are passed to 'gwpy.timeseries.TimeSeries.fetch_open_data()'
# sample_rate: most data are stored by LOSC at 4096 Hz, ifo_list = bilby.gw.detector.InterferometerList([])
# there may be event-related data releases with a 16384 Hz rate. for det in interferometer_names:
kwargs = {"sample_rate": 4096} logger.info("Downloading analysis data for ifo {}".format(det))
# For O2 events a "tag" is required to download the data. ifo = bilby.gw.detector.get_empty_interferometer(det)
# CLN = clean data; C00 or C01 = raw data data = TimeSeries.fetch_open_data(det, start_time, end_time)
# kwargs = {"tag": 'CLN'} ifo.set_strain_data_from_gwpy_timeseries(data)
# For some events can specify channels: source data stream for LOSC data. # Additional arguments you might need to pass to TimeSeries.fetch_open_data:
# kwargs = {"channel": {'H1': 'H1:DCS-CALIB_STRAIN_C02', # - sample_rate = 4096, most data are stored by LOSC at this frequency
# 'L1': 'L1:DCS-CALIB_STRAIN_C02', # there may be event-related data releases with a 16384Hz rate.
# 'V1': 'V1:FAKE_h_16384Hz_4R'}} # - tag = 'CLN' for clean data; C00/C01 for raw data (different releases)
# note that for O2 events a "tag" is required to download the data.
interferometers = bilby.gw.detector.get_event_data( # - channel = {'H1': 'H1:DCS-CALIB_STRAIN_C02',
label, interferometer_names=interferometer_names, # 'L1': 'L1:DCS-CALIB_STRAIN_C02',
duration=duration, roll_off=roll_off, psd_offset=psd_offset, # 'V1': 'V1:FAKE_h_16384Hz_4R'}}
psd_duration=psd_duration, cache=True, outdir=outdir, # for some events can specify channels: source data stream for LOSC data.
filter_freq=filter_freq, **kwargs) logger.info("Downloading psd data for ifo {}".format(det))
psd_data = TimeSeries.fetch_open_data(det, psd_start_time, psd_end_time)
psd_alpha = 2 * roll_off / duration # shape parameter of tukey window
psd = psd_data.psd(
fftlength=duration,
overlap=0,
window=("tukey", psd_alpha),
method="median")
ifo.power_spectral_density = bilby.gw.detector.PowerSpectralDensity(
frequency_array=psd.frequencies.value, psd_array=psd.value)
ifo_list.append(ifo)
logger.info("Saving data plots to {}".format(outdir))
bilby.core.utils.check_directory_exists_and_if_not_mkdir(outdir)
ifo_list.plot_data(outdir=outdir, label=label)
# CHOOSE PRIOR FILE # CHOOSE PRIOR FILE
# DEFAULT PRIOR FILES: GW150914.prior, binary_black_holes.prior, # DEFAULT PRIOR FILES: GW150914.prior, binary_black_holes.prior,
# binary_neutron_stars.prior (if bns, use BNSPriorDict) # binary_neutron_stars.prior (if bns, use BNSPriorDict)
# Needs to specify path for any other prior file. # Needs to specify path if you want to use any other prior file.
prior = bilby.gw.prior.BBHPriorDict(filename='GW150914.prior') prior = bilby.gw.prior.BBHPriorDict(filename='GW150914.prior')
# GENERATE WAVEFORM # GENERATE WAVEFORM
duration = None # duration and sampling frequency will be overwritten sampling_frequency = 4096. # same at which the data is stored
# to match the ones in the interferometers. conversion = bilby.gw.conversion.convert_to_lal_binary_black_hole_parameters
sampling_frequency = kwargs["sample_rate"] # same at which the data is stored
start_time = 0 # set the starting time of the time array
# reference_frequency: either low freq. limit or most sensitive freq.
# OVERVIEW OF APPROXIMANTS: # OVERVIEW OF APPROXIMANTS:
# https://www.lsc-group.phys.uwm.edu/ligovirgo/cbcnote/Waveforms/Overview # https://www.lsc-group.phys.uwm.edu/ligovirgo/cbcnote/Waveforms/Overview
waveform_arguments = {'waveform_approximant': 'IMRPhenomPv2', waveform_arguments = {
'reference_frequency': 50} 'waveform_approximant': 'IMRPhenomPv2',
source_model = bilby.gw.source.lal_binary_black_hole 'reference_frequency': 50 # most sensitive frequency
}
waveform_generator = bilby.gw.WaveformGenerator( waveform_generator = bilby.gw.WaveformGenerator(
duration=duration, sampling_frequency=sampling_frequency, parameter_conversion=conversion,
start_time=start_time, frequency_domain_source_model=bilby.gw.source.lal_binary_black_hole,
frequency_domain_source_model=source_model,
waveform_arguments=waveform_arguments) waveform_arguments=waveform_arguments)
# CHOOSE LIKELIHOOD FUNCTION # CHOOSE LIKELIHOOD FUNCTION
...@@ -80,26 +97,27 @@ waveform_generator = bilby.gw.WaveformGenerator( ...@@ -80,26 +97,27 @@ waveform_generator = bilby.gw.WaveformGenerator(
# Phase marginalisation is done analytically using a Bessel function. # Phase marginalisation is done analytically using a Bessel function.
# If prior given, used in the distance and phase marginalization. # If prior given, used in the distance and phase marginalization.
likelihood = bilby.gw.likelihood.GravitationalWaveTransient( likelihood = bilby.gw.likelihood.GravitationalWaveTransient(
interferometers, waveform_generator, time_marginalization=False, interferometers=ifo_list, waveform_generator=waveform_generator,
distance_marginalization=False, phase_marginalization=False) priors=prior, time_marginalization=False, distance_marginalization=False,
phase_marginalization=False)
# RUN SAMPLER # RUN SAMPLER
# Can use log_likelihood_ratio, rather than just the log_likelihood. # Can use log_likelihood_ratio, rather than just the log_likelihood.
# If using simulated data, pass a dictionary of injection parameters.
# A function can be specified in 'conversion_function' and applied # A function can be specified in 'conversion_function' and applied
# to posterior to generate additional parameters e.g. source-frame masses. # to posterior to generate additional parameters e.g. source-frame masses.
# Implemented Samplers: # Implemented Samplers:
# LIST OF AVAILABLE SAMPLERS: Run -> bilby.sampler.implemented_samplers # LIST OF AVAILABLE SAMPLERS: Run -> bilby.sampler.implemented_samplers
npoints = 500 # number of live points for the MCMC (dynesty) npoints = 512 # number of live points for the nested sampler
n_steps = 100 # min number of steps before proposing a new live point,
# defaults `ndim * 10`
sampler = 'dynesty' sampler = 'dynesty'
# Different samplers can have different additional kwargs, # Different samplers can have different additional kwargs,
# visit https://lscsoft.docs.ligo.org/bilby/samplers.html for details. # visit https://lscsoft.docs.ligo.org/bilby/samplers.html for details.
result = bilby.run_sampler( result = bilby.run_sampler(
likelihood, prior, outdir=outdir, label=label, likelihood, prior, outdir=outdir, label=label,
sampler=sampler, npoints=npoints, use_ratio=False, sampler=sampler, nlive=npoints, use_ratio=False,
injection_parameters=None, walks=n_steps, n_check_point=10000, check_point_plot=True,
conversion_function=bilby.gw.conversion.generate_all_bbh_parameters) conversion_function=bilby.gw.conversion.generate_all_bbh_parameters)
result.plot_corner() result.plot_corner()
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