Fast Response Function
- class fastlisaresponse.response.pyResponseTDI(sampling_frequency, num_pts, order=25, tdi='1st generation', orbits: ~lisatools.detector.Orbits | None = <class 'lisatools.detector.EqualArmlengthOrbits'>, tdi_orbits: ~lisatools.detector.Orbits | None = None, tdi_chan='XYZ', use_gpu=False)
Bases:
object
Class container for fast LISA response function generation.
The class computes the generic time-domain response function for LISA. It takes LISA constellation orbital information as input and properly determines the response for these orbits numerically. This includes both the projection of the gravitational waves onto the LISA constellation arms and combinations of projections into TDI observables. The methods and maths used can be found [here](https://arxiv.org/abs/2204.06633).
This class is also GPU-accelerated, which is very helpful for Bayesian inference methods.
- Parameters:
sampling_frequency (double) – The sampling rate in Hz.
num_pts (int) – Number of points to produce for the final output template.
order (int, optional) – Order of Lagrangian interpolation technique. Lower orders will be faster. The user must make sure the order is sufficient for the waveform being used. (default: 25)
tdi (str or list, optional) – TDI setup. Currently, the stock options are
'1st generation'
and'2nd generation'
. Or the user can provide a list of tdi_combinations of the form{"link": 12, "links_for_delay": [21, 13, 31], "sign": 1, "type": "delay"}
.'link'
(int) the link index (12, 21, 13, 31, 23, 32) for the projection (\(y_{ij}\)).'links_for_delay'
(list) are the link indexes as a list used for delays applied to the link projections.'sign'
is the sign in front of the contribution to the TDI observable. It takes the value of +1 or -1.type
is either"delay"
or"advance"
. It is optional and defaults to"delay"
. (default:"1st generation"
)orbits (
Orbits
, optional) – Orbits class from LISA Analysis Tools. Works with LISA Orbits outputs: lisa-simulation.pages.in2p3.fr/orbits/. (default:EqualArmlengthOrbits
)tdi_chan (str, optional) – Which TDI channel combination to return. Choices are
'XYZ'
,AET
, orAE
. (default:'XYZ'
)tdi_orbits (
Orbits
, optional) –Set if different orbits from projection. Orbits class from LISA Analysis Tools. Works with LISA Orbits outputs: lisa-simulation.pages.in2p3.fr/orbits/. (default:
EqualArmlengthOrbits
)use_gpu (bool, optional) – If True, run code on the GPU. (default:
False
)
- A_in
Array containing y values for linear spline of A during Lagrangian interpolation.
- Type:
xp.ndarray
- buffer_integer
Self-determined buffer necesary for the given value for
order
.- Type:
int
- channels_no_delays
Carrier of link index and sign information for arms that do not get delayed during TDI computation.
- Type:
2D np.ndarray
- deps
The spacing between Epsilon values in the interpolant for the A quantity in Lagrangian interpolation. Hard coded to 1/(
num_A
- 1).- Type:
double
- dt
Inverse of the sampling_frequency.
- Type:
double
- E_in
Array containing y values for linear spline of E during Lagrangian interpolation.
- Type:
xp.ndarray
- half_order
Half of
order
adjusted to beint
.- Type:
int
- link_inds
Link indexes for delays in TDI.
- Type:
xp.ndarray
- link_space_craft_0_in
Link indexes for receiver on each arm of the LISA constellation.
- Type:
xp.ndarray
- link_space_craft_1_in
Link indexes for emitter on each arm of the LISA constellation.
- Type:
xp.ndarray
- nlinks
The number of links in the constellation. Typically 6.
- Type:
int
- num_A
Number of points to use for A spline values used in the Lagrangian interpolation. This is hard coded to 1001.
- Type:
int
- num_channels
- Type:
int
- num_pts
Number of points to produce for the final output template.
- Type:
int
- order
Order of Lagrangian interpolation technique.
- Type:
int
- sampling_frequency
The sampling rate in Hz.
- Type:
double
- tdi
TDI setup.
- Type:
str or list
- tdi_buffer
The buffer necessary for all information needed at early times for the TDI computation. This is set to 200.
- Type:
int
- use_gpu
If True, run on GPU.
- Type:
bool
- xp
Either Numpy or Cupy.
- Type:
obj
- property response_gen: callable
CPU/GPU function for generating the projections.
- property tdi_gen: callable
CPU/GPU function for generating tdi.
- property response_orbits: Orbits
Response function orbits.
- property tdi_orbits: Orbits
TDI function orbits.
- property citation
Get citations for use of this code
- property tdi_combinations: List
TDI Combination setup
- property y_gw
Projections along the arms
- get_projections(input_in, lam, beta, t0=10000.0)
Compute projections of GW signal on to LISA constellation
- Parameters:
input_in (xp.ndarray) – Input complex time-domain signal. It should be of the form: \(h_+ + ih_x\). If using the GPU for the response, this should be a CuPy array.
lam (double) – Ecliptic Longitude in radians.
beta (double) – Ecliptic Latitude in radians.
t0 (double, optional) – Time at which to the waveform. Because of the delays and interpolation towards earlier times, the beginning of the waveform is garbage.
t0
tells the waveform generator where to start the waveform compraed tot=0
.
- Raises:
ValueError – If
t0
is not large enough.
- property XYZ
Return links as an array
- get_tdi_delays(y_gw=None)
Get TDI combinations from projections.
This functions generates the TDI combinations from the projections computed with
get_projections
. It can return XYZ, AET, or AE depending on what was input fortdi_chan
into__init__
.- Parameters:
y_gw (xp.ndarray, optional) – Projections along each link. Must be a 2D
numpy
orcupy
array with shape:(nlinks, num_pts)
. The links must be entered in the proper order in the code: 21, 12, 31, 13, 32, 23. (Default: None)- Returns:
(X,Y,Z) or (A,E,T) or (A,E)
- Return type:
tuple
- Raises:
ValueError – If
tdi_chan
is not one of the options.
Response Function Wrapper
- class fastlisaresponse.response.ResponseWrapper(waveform_gen, Tobs, dt, index_lambda, index_beta, t0=10000.0, flip_hx=False, remove_sky_coords=False, is_ecliptic_latitude=True, use_gpu=False, remove_garbage=True, n_overide=None, orbits: ~lisatools.detector.Orbits | None = <class 'lisatools.detector.EqualArmlengthOrbits'>, **kwargs)
Bases:
object
Wrapper to produce LISA TDI from TD waveforms
This class takes a waveform generator that produces \(h_+ \pm ih_x\). (
flip_hx
is used if the waveform produces \(h_+ - ih_x\)). It takes the complex waveform in the SSB frame and produces the TDI channels according to settings chosen forpyResponseTDI
.The waveform generator must have
kwargs
withT
for the observation time in years anddt
for the time step in seconds.- Parameters:
waveform_gen (obj) – Function or class (with a
__call__
function) that takes parameters and produces \(h_+ \pm h_x\).Tobs (double) – Observation time in years.
dt (double) – Time between time samples in seconds. The inverse of the sampling frequency.
index_lambda (int) – The user will input parameters. The code will read these in with the
*args
formalism producing a list.index_lambda
tells the class the index of the ecliptic longitude within this list of parameters.index_beta (int) – The user will input parameters. The code will read these in with the
*args
formalism producing a list.index_beta
tells the class the index of the ecliptic latitude (or ecliptic polar angle) within this list of parameters.t0 (double, optional) – Start of returned waveform in seconds leaving ample time for garbage at the beginning of the waveform. It also removed the same amount from the end. (Default: 10000.0)
flip_hx (bool, optional) – If True,
waveform_gen
produces \(h_+ - ih_x\).pyResponseTDI
takes \(h_+ + ih_x\), so this setting will multiply the cross polarization term out of the waveform generator by -1. (Default:False
)remove_sky_coords (bool, optional) – If True, remove the sky coordinates from the
*args
list. This should be set to True if the waveform generator does not take in the sky information. (Default:False
)is_ecliptic_latitude (bool, optional) – If True, the latitudinal sky coordinate is the ecliptic latitude. If False, thes latitudinal sky coordinate is the polar angle. In this case, the code will convert it with \(eta=\pi / 2 - \Theta\). (Default:
True
)use_gpu (bool, optional) – If True, use GPU. (Default:
False
)remove_garbage (bool or str, optional) – If True, it removes everything before
t0
and after the end time -t0
. Ifstr
, it must be"zero"
. If"zero"
, it will not remove the points, but set them to zero. This is ideal for PE. (Default:True
)n_overide (int, optional) – If not
None
, this will override the determination of the number of points,n
, fromint(T/dt)
to then_overide
. This is used if there is an issue matching points between the waveform generator and the response model.orbits (
Orbits
, optional) –Orbits class from LISA Analysis Tools. Works with LISA Orbits outputs: lisa-simulation.pages.in2p3.fr/orbits/. (default:
EqualArmlengthOrbits
)**kwargs (dict, optional) – Keyword arguments passed to
pyResponseTDI
.
- property citation
Get citations for use of this code
- __call__(*args, **kwargs)
Run the waveform and response generation
- Parameters:
*args (list) – Arguments to the waveform generator. This must include the sky coordinates.
**kwargs (dict) – kwargs necessary for the waveform generator.
- Returns:
TDI Channels.
- Return type:
list
Utilities
- fastlisaresponse.utils.utility.get_overlap(sig1, sig2, phase_maximize=False, use_gpu=False)
Calculate the mismatch across TDI channels
Calculates the overlap between two sets of TDI observables in the time domain. The overlap is complex allowing for the addition of overlap over all channels. It can be phase maximized as well.
This function has GPU capabilities.
- Parameters:
sig1 (list or xp.ndarray) – TDI observables for first signal. Must be
list
ofxp.ndarray
or a singlexp.ndarray
. Must have same length assig2
in terms of number of channels and length of the indivudal channels.sig2 (list or xp.ndarray) – TDI observables for second signal. Must be
list
ofxp.ndarray
or a singlexp.ndarray
. Must have same length assig1
in terms of number of channels and length of the individual channels.phase_maximize (bool, optional) – If
True
, maximize over the phase in the overlap. This is equivalent to getting the magnitude of the phasor that is the complex overlap. (Defaut:False
)use_gpu (bool, optional) – If
True
, use the GPU. This setsxp=cupy
. IfFalse, use the CPU and set ``xp=numpy
.
- Returns:
Overlap as a real value.
- Return type:
double