From f9d53301c33204fcffa75674229809392da8ac5c Mon Sep 17 00:00:00 2001 From: acorreia61201 Date: Mon, 3 Jun 2024 18:57:16 +0000 Subject: [PATCH] Added reference_time to methods as kwarg; revised methods to be more consistent with LIGO detector class --- pycbc/detector.py | 25 ++++++++++++++++++------- 1 file changed, 18 insertions(+), 7 deletions(-) diff --git a/pycbc/detector.py b/pycbc/detector.py index b8f1ad0d75f..cf05edf7e62 100644 --- a/pycbc/detector.py +++ b/pycbc/detector.py @@ -683,7 +683,7 @@ def __init__(self, orbits=ESAOrbits(), use_gpu=False): # initialize whether to use gpu; FLR has handles if this cannot be done self.use_gpu = use_gpu - def project_wf(self, hp, hc, lamb, beta, t0=10000., use_gpu=None): + def get_links(self, hp, hc, lamb, beta, t0=10000., use_gpu=None, reference_time=None): """ Project a radiation frame (SSB) waveform to the LISA constellation. This *does not* perform TDI; this only produces a detector frame @@ -711,6 +711,10 @@ def project_wf(self, hp, hc, lamb, beta, t0=10000., use_gpu=None): Flag whether to use GPU support. Default to class input. CuPy is required if use_gpu == True; an ImportError will be raised if CuPy could not be imported. + + reference_time : float (optional) + The GPS start time of the GW signal in the SSB frame. Default to + value in input signals hp and hc. Returns ------- @@ -724,9 +728,12 @@ def project_wf(self, hp, hc, lamb, beta, t0=10000., use_gpu=None): n = len(hp) # number of points # configure the orbits file according to wf times + if reference_time is not None: + hp.start_time = reference_time + self.sample_times = hp.sample_times.numpy() + ### FIXME: This currently doesn't work. There seems to be a bug in ### LISAanalysistools that breaks the code when specifying a time array. - # self.sample_times = hp.sample_times.numpy() # self.orbits.configure(t_arr=self.sample_times) # format wf to FLR standard hp + i*hc as ndarray @@ -748,7 +755,7 @@ def project_wf(self, hp, hc, lamb, beta, t0=10000., use_gpu=None): if self.response_init is None: # initialize the class ### TDI is set to '2nd generation' here to match default value in get_tdi() - # see FIXME in get_tdi() + ### see FIXME in get_tdi() response_init = pyResponseTDI(1/dt, n, orbits=self.orbits, use_gpu=use_gpu, tdi='2nd generation') self.response_init = response_init else: @@ -769,10 +776,10 @@ def project_wf(self, hp, hc, lamb, beta, t0=10000., use_gpu=None): return wf_proj - def get_tdi(self, hp, hc, lamb, beta, tdi='2nd generation', tdi_chan='XYZ', tdi_orbits=None, + def project_wave(self, hp, hc, lamb, beta, reference_time=None, tdi='2nd generation', tdi_chan='XYZ', tdi_orbits=None, use_gpu=None, remove_garbage=True): """ - Calculate the TDI variables given a waveform and LISA orbit data. + Evaluate the TDI observables Parameters ---------- @@ -787,7 +794,11 @@ def get_tdi(self, hp, hc, lamb, beta, tdi='2nd generation', tdi_chan='XYZ', tdi_ beta : float The ecleptic longitude in the SSB frame. - + + reference_time : float (optional) + The GPS start time of the GW signal in the SSB frame. Default to + value in input signals hp and hc. + tdi : string (optional) TDI channel configuration. Accepts "1st generation" or "2nd generation" as inputs. Default "2nd generation". @@ -846,7 +857,7 @@ def get_tdi(self, hp, hc, lamb, beta, tdi='2nd generation', tdi_chan='XYZ', tdi_ self.response_init._init_tdi_delays() # generate the Doppler time series - links = self.project_wf(hp, hc, lamb, beta, t0=self.t0, use_gpu=use_gpu) + links = self.get_links(hp, hc, lamb, beta, t0=self.t0, use_gpu=use_gpu, reference_time=reference_time) # generate the TDI channels tdi_obs = self.response_init.get_tdi_delays()