Extending GBGPU Waveforms
The base class gbgpu.gbgpu.GBGPUBase
can be inherited in order to build other waveforms in the model of FastGB. To do this, there are required methods that need to be added to the new waveform class. They are described in the Abstract Base Class gbgpu.gbgpu.GBGPUBase
. After the base class, waveform models that have already been extended beyond the base are described.
GBGPUBase
base class
- class gbgpu.gbgpu.GBGPUBase(orbits: Orbits = None, force_backend=None)
Bases:
GBGPUParallelModule
,ABC
Generate Galactic Binary Waveforms
This class generates galactic binary waveforms in the frequency domain, in the form of LISA TDI channels X, A, and E. It generates waveforms in batches. It can also provide injection signals and calculate likelihoods in batches. These batches are run on GPUs or CPUs. When CPUs are used, all available threads are leveraged with OpenMP. To adjust the available threads, use
OMP_NUM_THREADS
environmental variable orgbgpu.utils.set_omp_num_threads()
.This class can generate waveforms for two different types of GB sources:
Circular Galactic binaries
Circular Galactic binaries with an eccentric third body (inherited)
- Parameters:
force_backend (str, optional) – Change to general backend in use. Example options are
'cpu'
and'gpu'
. Default isNone
.
- get_basis_tensors
Cython function.
- Type:
obj
- GenWave
Cython function.
- Type:
obj
- GenWaveThird
Cython function.
- Type:
obj
- unpack_data_1
Cython function.
- Type:
obj
- XYZ
Cython function.
- Type:
obj
- num_bin
Number of binaries in the current calculation.
- Type:
int
- N_max
Maximum points in a waveform based on maximum harmonic mode considered.
- Type:
int
- start_inds
Start indices into data stream array. q - N/2.
- Type:
list of 1D int xp.ndarray
- df
Fourier bin spacing.
- Type:
double
- X_out, A_out, E_out
X, A, or E channel TDI templates. Each array is a 2D complex array of shape (number of points, number of binaries) that is flattened. These can be accessed in python with the properties
X
,A
,E
.- Type:
1D complex xp.ndarrays
- N
Last N value used.
- Type:
int
- d_d
<d|d> term in the likelihood.
- Type:
double
- classmethod supported_backends()
List of backends supported by a parallel module by order of preference.
- property get_ll_func
get_ll c func.
- property fill_global_func
fill_global c func.
- property orbits: Orbits
Orbits class.
- property citation
Get citations for this class
- run_wave(amp, f0, fdot, fddot, phi0, iota, psi, lam, beta, *args, N=None, T=126232599.05418241, dt=10.0, oversample=1, tdi2=False)
Create waveforms in batches.
This call creates the TDI templates in batches.
The parameters and code below are based on an implementation of Fast GB in the LISA Data Challenges’
ldc
package.This class can be inherited to build fast waveforms for systems with additional astrophysical effects.
- Parameters:
amp (double or 1D double np.ndarray) – Amplitude parameter.
f0 (double or 1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
fdot (double or 1D double np.ndarray) – Initial time derivative of the frequency given as Hz/s.
fddot (double or 1D double np.ndarray) – Initial second derivative with respect to time of the frequency given in Hz/s^2.
phi0 (double or 1D double np.ndarray) – Initial phase angle of gravitational wave given in radians.
iota (double or 1D double np.ndarray) – Inclination of the Galactic binary orbit given in radians.
psi (double or 1D double np.ndarray) – Polarization angle of the Galactic binary orbit in radians.
lam (double or 1D double np.ndarray) – Ecliptic longitutude of the source given in radians.
beta (double or 1D double np.ndarray) – Ecliptic Latitude of the source given in radians. This is converted to the spherical polar angle.
*args (tuple, optional) – Flexible parameter to allow for a flexible number of argmuments when inherited by other classes. If running a circular Galactic binarys,
args = ()
. Iflen(args) != 0
, then the inheriting class must have aprepare_additional_args
method.N (int, optional) – Number of points in waveform. This should be determined by the initial frequency,
f0
. Default isNone
. IfNone
, will usegbgpu.utils.utility.get_N()
function to determine properN
.T (double, optional) – Observation time in seconds. Default is
4 * YRSID_SI
.dt (double, optional) – Observation cadence in seconds. Default is
10.0
seconds.oversample (int, optional) – Oversampling factor compared to the determined
N
value. Final N will beoversample * N
. This is only used if N is not provided. Default is1
.tdi2 (bool, optional) – If
True
, produce the TDI channels for TDI 2nd-generation. IfFalse
, produce TDI 1st-generation. Technically, the current TDI computation is not valid for generic LISA orbits, which are dealth with with 2nd-generation TDI, only those with an “equal-arm length” condition. Default isFalse
.Raises – ValueError: Length of
*args
is not 0 or 5.
- property X_out
X channel.
- property Y_out
Y channel.
- property Z_out
Z channel.
- property A_out
A channel.
- property E_out
E channel.
- property T_out
T channel.
- property A
return A channel reshaped based on number of binaries
- property E
return E channel reshaped based on number of binaries
- property T
return T channel reshaped based on number of binaries
- property X
return X channel reshaped based on number of binaries
- property Y
return Y channel reshaped based on number of binaries
- property Z
return Z channel reshaped based on number of binaries
- property freqs
Return frequencies associated with each signal
- get_ll(params, data, psd, phase_marginalize=False, start_freq_ind=0, data_index=None, noise_index=None, **kwargs)
Get batched log likelihood
Generate the individual log likelihood for a batched set of Galactic binaries. This is also GPU/CPU agnostic.
- Parameters:
params (2D double np.ndarrays) – Parameters of all binaries to be calculated. The shape is
(number of parameters, number of binaries)
.data (length 2 list of 1D or 2D complex128 xp.ndarrays) – List of arrays representing the data stream. These should be CuPy arrays if running on the GPU, NumPy arrays if running on a CPU. The list should be [A channel, E channel]. Should be 1D if only one data stream is analyzed. If 2D, shape is
(number of data streams, data_length)
. If 2D, user must also providedata_index
kwarg.psd (length 2 list of 1D or 2D double xp.ndarrays) – List of arrays representing the power spectral density (PSD) in the noise. These should be CuPy arrays if running on the GPU, NumPy arrays if running on a CPU. The list should be [A channel, E channel]. Should be 1D if only one PSD is analyzed. If 2D, shape is
(number of PSDs, data_length)
. If 2D, user must also providenoise_index
kwarg.phase_marginalize (bool, optional) – If True, marginalize over the initial phase. Default is False.
start_freq_ind (int, optional) – Starting index into the frequency-domain data stream for the first entry of
data
/psd
. This is used if a subset of a full data stream is presented for the computation. If providing mutliple data streams indata
, this single start index value will apply to all of them.data_index (1D xp.int32 array, optional) – If providing 2D
data
, need to providedata_index
to indicate the data stream associated with each waveform for which the log-Likelihood is being computed. For example, if you have 100 binaries with 5 different data streams,data_index
will be a length-100 xp.int32 array with values 0 to 4, indicating the specific data stream to use for each source. IfNone
, this will be filled with zeros and only analyzed with the first data stream given. Default isNone
.noise_index (1D xp.int32 array, optional) – If providing 2D
psd
, need to providenoise_index
to indicate the PSD associated with each waveform for which the log-Likelihood is being computed. For example, if you have 100 binaries with 5 different PSDs,noise_index
will be a length-100 xp.int32 array with values 0 to 4, indicating the specific PSD to use for each source. IfNone
, this will be filled with zeros and only analyzed with the first PSD given. Default isNone
.**kwargs (dict, optional) – Passes keyword arguments to the
run_wave()
method.
- Raises:
TypeError – If data arrays are NumPy/CuPy while template arrays are CuPy/NumPy.
- Returns:
Log likelihood values associated with each binary.
- Return type:
1D double np.ndarray
- fill_global_template(group_index, templates, A, E, start_inds, N=None, start_freq_ind=0)
Fill many global templates with waveforms
This method takes already generated waveforms (
A, E, start_inds
) and their associated grouping index (group_index
) and fills buffer tempalte arrays (templates
).This method combines waveforms that have already been created. When a user does not have the waveforms in hand, they should use the
generate_global_template()
method.- Parameters:
group_index (1D double int32 xp.ndarray) – Index indicating to which template each individual binary belongs.
templates (3D complex128 xp.ndarray) – Buffer array for template output to filled in place. The shape is
(number of templates, 2, data_length)
. The2
is for theA
andE
TDI channels in that order.A (1D or 2D complex128 xp.ndarray) – TDI A channel template values for each individual binary. The shape if 2D is
(number of binaries, N)''. In 1D, the array should be arranged so that it resembles ``(number of binaries, N).transpose().flatten()
. After running waveforms, this is howself.A_out
is arranged.E (1D 2D complex128 xp.ndarray) – TDI E channel template values for each individual binary. The shape if 2D is
(number of binaries, N)''. In 1D, the array should be arranged so that it resembles ``(number of binaries, N).transpose().flatten()
. After running waveforms, this is howself.E_out
is arranged.start_inds (1D int32 xp.ndarray) – The start indices of each binary waveform in the full Fourier transform:
int(f0/T) - N/2
.N (int, optional) – The length of the A and E channels for each individual binary. When
A
andE
are 1D,N
must be given. Default isNone
.start_freq_ind (int, optional) – Starting index into the frequency-domain data stream for the first entry of
templates
. This is used if a subset of a full data stream is presented for the computation.
- Raises:
TypeError – If data arrays are NumPy/CuPy while tempalte arrays are CuPy/NumPy.
ValueError – Inputs are not correctly provided.
- generate_global_template(params, group_index, templates, start_freq_ind=0, **kwargs)
Generate global templates from binary parameters
Generate waveforms in batches and then combine them into global fit templates. This method wraps
fill_global_template()
by building the waveforms first.- Parameters:
params (2D double np.ndarrays) – Parameters of all binaries to be calculated. The shape is
(number of parameters, number of binaries)
.group_index (1D double int32 xp.ndarray) – Index indicating to which template each individual binary belongs.
templates (3D complex128 xp.ndarray) – Buffer array for template output to filled in place. The shape is
(number of templates, 2, data_length)
. The2
is for theA
andE
TDI channels in that order.start_freq_ind (int, optional) – Starting index into the frequency-domain data stream for the first entry of
templates
. This is used if a subset of a full data stream is presented for the computation.**kwargs (dict, optional) – Passes keyword arguments to
run_wave()
function above.
- inject_signal(*args, fmax=None, T=126232599.05418241, dt=10.0, **kwargs)
Inject a single signal
Provides the injection of a single signal into a data stream with frequencies spanning from 0.0 to fmax with 1/T spacing (from Fourier transform).
- Parameters:
*args (list, tuple, or 1D double np.array) – Arguments to provide to
run_wave()
to build the TDI templates for injection.fmax (double, optional) – Maximum frequency to use in data stream. If
None
, will use1/(2 * dt)
. Default isNone
.T (double, optional) – Observation time in seconds. Default is
4 * YRSID_SI
.dt (double, optional) – Observation cadence in seconds. Default is
10.0
seconds.**kwargs (dict, optional) – Passes kwargs to
run_wave()
.
- Returns:
- NumPy arrays for the A channel and
E channel:
(A channel, E channel)
. Need to convert to CuPy if working on GPU.
- Return type:
Tuple of 1D np.ndarrays
- information_matrix(params, eps=1e-09, parameter_transforms={}, inds=None, N=1024, psd_func=None, psd_kwargs={}, easy_central_difference=False, return_gpu=False, **kwargs)
Get the information matrix for a batch.
This function computes the Information matrix for a batch of Galactic binaries. It uses a 2nd order calculation for the derivative if
easy_central_difference
isFalse
:..math:: frac{dh}{dlambda_i} = frac{-h(lambda_i + 2epsilon) + h(lambda_i - 2epsilon) + 8(h(lambda_i + epsilon) - h(lambda_i - epsilon))}{12epsilson}
Otherwise, it will just calculate the derivate with a first-order central difference.
This function maps all parameter values to 1. For example, if the square root of the diagonal of the associated covariance matrix is 1e-7 for the frequency parameter, then the standard deviation in the frequency is
1e-7 * f0
. To properly use with covariance values not on the diagonal, they will have to be multipled by the parameters: \(C_{ij} \vec{\theta}_j\).- Parameters:
params (2D double np.ndarrays) – Parameters of all binaries to be calculated. The shape is
(number of parameters, number of binaries)
. Seerun_wave()
for more information on the adjustable number of parameters when calculating for a third body.eps (double, optional) – Step to take when calculating the derivative. The step is relative difference. Default is
1e-9
.parameter_transforms (dict, optional) – Dictionary containing the parameter transform functions. The keys in the dict should be the index associated with the parameter. The items should be the actual transform function. Default is no transforms (
{}
).inds (1D int np.ndarray, optional) – Numpy array with the indexes of the parameters to test in the Information matrix. Default is
None
. When it is not given, it defaults to all parameters.N (int, optional) – Number of points for the waveform. Same as the
N
parameter inrun_wave()
. We recommend using higherN
in the Information Matrix computation because of the numerical derivatives. Default is1024
.psd_func (object, optional) – Function to compute the PSD for the A and E channels. Must take on argument: the frequencies as an xp.ndarray. When
None
, it attemps to use the sensitivity functions from LISA Analysis Tools.psd_kwargs (dict, optional) – Keyword arguments for the TDI noise generator. Default is
None
.easy_central_difference (bool, optional) – If
True
, compute the derivatives with a first-order central difference computation. IfFalse
, use the higher order derivative that computes two more waveforms during the derivative calculation. Default isFalse
.return_gpu (False, optional) – If
True
and backend is GPU-based, return Information matrices in cupy array. Default isFalse
.
- Returns:
Information Matrices for all binaries with shape:
(number of binaries, number of parameters, number of parameters)
.- Return type:
3D xp.ndarray
- Raises:
ValueError – Step size issues.
ModuleNotFoundError – LISA Analysis Tools package not available. Occurs when NOT providing
psd_func
kwarg.
- abstract add_to_argS(argS, f0, fdot, fddot, xi, *args)
Update
argS
in FastGB formalism for third-body effectargS
is an effective phase that goes intokdotP
in the construction of the slow part of the waveform.kdotP
is then included directly in the transfer function. Seegbgpu.gbgpu.GBGPU._construct_slow_part()
for the use of argS in the larger code.- Parameters:
argS (3D double xp.ndarray) – Special phase evaluation that goes into
kdotP
. Shape is(num binaries, 3 spacecraft, N)
.f0 (1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
fdot (1D double np.ndarray) – Initial time derivative of the frequency given as Hz/s.
fddot (1D double np.ndarray) – Initial second derivative with respect to time of the frequency given in Hz/s^2.
xi (3D double xp.ndarray) – Time at each spacecraft. The shape is
(num binaries, 3 spacecraft, N)
.T (double) – Observation time in seconds.
*args (tuple) – Args returned from
prepare_additional_args()
.
- Returns:
Updated
argS
with third-body effect- Return type:
3D double xp.ndarray
- abstract prepare_additional_args(*args)
Prepare the arguments special to this class
This function must take in the extra
args
input intoGBGPU.run_wave()
and transform them as needed to input into the rest of the code. If using GPUs, this is where the parameters are copied to GPUs.- Parameters:
*args (tuple) – Any additional args to be dealt with.
- Returns:
New args. In the rest of the code this is
add_args
.- Return type:
Tuple
- abstract special_get_N(amp, f0, T, *args, oversample=1)
Determine proper sampling rate in time domain for slow-part.
- Parameters:
amp (double or 1D double np.ndarray) – Amplitude parameter.
f0 (double or 1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
T (double) – Observation time in seconds.
*args (tuple) – Args input for beyond-GBGPU functionality.
oversample (int, optional) – Oversampling factor compared to the determined
N
value. Final N will beoversample * N
. This is only used if N is not provided. Default is1
.
- Returns:
Number of time-domain points recommended for each binary.
- Return type:
1D int32 xp.ndarray
- abstract shift_frequency(fi, xi, *args)
Shift the evolution of the frequency in the slow part
- Parameters:
fi (3D double xp.ndarray) – Instantaneous frequencies of the wave before applying third-body effect at each spacecraft as a function of time. The shape is
(num binaries, 3 spacecraft, N)
.xi (3D double xp.ndarray) – Time at each spacecraft. The shape is
(num binaries, 3 spacecraft, N)
.*args (tuple) – Args returned from
prepare_additional_args()
.
- Returns:
Updated frequencies with third-body effect.
- Return type:
3D double xp.ndarray
Third-body inclusion
gbgpu.gbgpu.GBGPUBase
has been extended to include a third-body in long orbit around the inner binary. This waveform was first built in arXiv:1806.00500. It was more recently used and adapted into this code base for # TODO: add arxiv. Please cite both papers, as well as the base FastGB papers, if you use this waveform.
Third-body waveform
- class gbgpu.thirdbody.GBGPUThirdBody(orbits: Orbits = None, force_backend=None)
Bases:
GBGPUBase
Build the effect of a third body into Galactic binary waveforms.
The third-body components are originally by Travis Robson for the paper arXiv:1806.00500.
- Parameters:
force_backend (str, optional) – Change to general backend in use. Example options are
'cpu'
and'gpu'
. Default isNone
.
- property citation
Get citations for this class
- prepare_additional_args(A2, varpi, e2, P2, T2)
Prepare the arguments special to this class
- Parameters:
A2 (double or 1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (double or 1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (double or 1D double np.ndarray) – Eccentricity of the third body orbit.
P2 (double or 1D double np.ndarray) – Period of the third body orbit in Years.
T2 (double or 1D double np.ndarray) – Time of pericenter passage of the third body in Years. This parameter is effectively a constant of integration.
- Returns:
- (A2, varpi, e2, n2, T2) adjusted for GPU usage if necessary.
(\(n_2=\frac{2\pi}{P_2}\) is the angular frequency of the orbit.)
- Return type:
Tuple
- special_get_N(amp, f0, T, A2, varpi, e2, P2, T2, oversample=1)
Determine proper sampling rate in time domain for slow-part.
- Parameters:
amp (double or 1D double np.ndarray) – Amplitude parameter.
f0 (double or 1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
T (double) – Observation time in seconds.
A2 (double or 1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above. (Not needed in this function)
varpi (double or 1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above. (Not needed in this function)
e2 (double or 1D double np.ndarray) – Eccentricity of the third body orbit. (Not needed in this function)
P2 (double or 1D double np.ndarray) – Period of the third body orbit in Years.
T2 (double or 1D double np.ndarray) – Time of pericenter passage of the third body in Years. This parameter is effectively a constant of integration. (Not needed in this function)
- Returns:
Number of time-domain points recommended for each binary.
- Return type:
1D int32 xp.ndarray
- Raises:
AssertionError – Shapes of inputs are wrong.
- shift_frequency(fi, xi, A2, varpi, e2, n2, T2)
Shift the evolution of the frequency in the slow part
- Parameters:
fi (3D double xp.ndarray) – Instantaneous frequencies of the wave before applying third-body effect at each spacecraft as a function of time. The shape is
(num binaries, 3 spacecraft, N)
.xi (3D double xp.ndarray) – Time at each spacecraft. The shape is
(num binaries, 3 spacecraft, N)
.A2 (1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (1D double np.ndarray) – Eccentricity of the third body orbit.
n2 (1D double np.ndarray) – Angular frequency of the third body orbit in per seconds.
T2 (1D double np.ndarray) – Time of pericenter passage of the third body in seconds. This parameter is effectively a constant of integration.
- Returns:
Updated frequencies with third-body effect.
- Return type:
3D double xp.ndarray
- add_to_argS(argS, f0, fdot, fddot, xi, A2, varpi, e2, n2, T2)
Update
argS
in FastGB formalism for third-body effectargS
is an effective phase that goes intokdotP
in the construction of the slow part of the waveform.kdotP
is then included directly in the transfer function. Seegbgpu.gbgpu.GBGPU._construct_slow_part()
for the use of argS in the larger code.- Parameters:
argS (3D double xp.ndarray) – Special phase evaluation that goes into
kdotP
. Shape is(num binaries, 3 spacecraft, N)
.f0 (1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
fdot (1D double np.ndarray) – Initial time derivative of the frequency given as Hz/s.
fddot (1D double np.ndarray) – Initial second derivative with respect to time of the frequency given in Hz/s^2.
xi (3D double xp.ndarray) – Time at each spacecraft. The shape is
(num binaries, 3 spacecraft, N)
.T (double) – Observation time in seconds.
A2 (1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (1D double np.ndarray) – Eccentricity of the third body orbit.
n2 (1D double np.ndarray) – Angular frequency of the third body orbit in per seconds.
T2 (1D double np.ndarray) – Time of pericenter passage of the third body in seconds. This parameter is effectively a constant of integration.
- Returns:
Updated
argS
with third-body effect- Return type:
3D double xp.ndarray
- get_u(l, e)
Invert Kepler’s Equation to get eccentric anomaly
Invert Kepler’s equation (\(l = u - e \sin{u}\)) using Mikkola’s method (1987) referenced in Tessmer & Gopakumar 2007.
- Parameters:
l (1D double xp.ndarray) – Mean anomaly in radians.
e (1D double xp.ndarray) – Eccentricity.
- Returns:
Eccentric anomaly
- Return type:
3D double xp.ndarray
- get_phi(t, T, e, n)
Get phi value for Line-of-sight velocity
See arXiv:1806.00500.
- Parameters:
t (3D double xp.ndarray) – Time values to evaluate \(\bar{\phi}\).
T (1D double xp.ndarray) – Time of periastron passage (
T2
) in seconds.e (1D double xp.ndarray) – Eccentricity.
n (1D double xp.ndarray) – Angular frequency of third-body orbit (
n2
) in per seconds.
- Returns:
Phi values for line-of-sight velocities.
- Return type:
3D double xp.ndarray
- get_vLOS(t, A2, varpi, e2, n2, T2)
Calculate the line-of-site velocity
See equation 13 in arXiv:1806.00500.
- Parameters:
t (3D double xp.ndarray) – Time values to evaluate \(\bar{\phi}\).
A2 (1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (1D double np.ndarray) – Eccentricity of the third body orbit.
n2 (1D double np.ndarray) – Angular frequency of the third body orbit in per seconds.
T2 (1D double np.ndarray) – Time of pericenter passage of the third body in seconds. This parameter is effectively a constant of integration.
- Returns:
LOS velocity.
- Return type:
3D double xp.ndarray
- parab_step_ET(f0, fdot, fddot, A2, varpi, e2, n2, T2, t0, t0_old)
Determine phase difference caused by third-body
Takes a step in the integration of the orbit. In this setup, the calculations can all be done in parallel because we are just inverted Kepler’s equation rather than integrating an ODE where a serial operation is required. TODO: Check this again. Was checked in past.
- Parameters:
f0 (1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
fdot (1D double np.ndarray) – Initial time derivative of the frequency given as Hz/s.
fddot (1D double np.ndarray) – Initial second derivative with respect to time of the frequency given in Hz/s^2.
A2 (1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (1D double np.ndarray) – Eccentricity of the third body orbit.
n2 (1D double np.ndarray) – Angular frequency of the third body orbit in per seconds.
T2 (1D double np.ndarray) – Time of pericenter passage of the third body in seconds. This parameter is effectively a constant of integration.
t0 (3D double xp.ndarray) – Time values at end of step. Shape is
(num binaries, 3 spacecraft, N - 1)
.t0_old (3D double xp.ndarray) – Time values at start of step. Shape is
(num binaries, 3 spacecraft, N - 1)
.
- Returns:
Phase shifts due to third-body effect.
- Return type:
3D double xp.ndarray
- get_aLOS(A2, varpi, e2, P2, T2, t, eps=1e-09)
Get line-of-sight acceleration
Use central difference with LOS velocity to get LOS acceleration.
- Parameters:
A2 (1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (1D double np.ndarray) – Eccentricity of the third body orbit.
n2 (1D double np.ndarray) – Angular frequency of the third body orbit in per seconds.
T2 (1D double np.ndarray) – Time of pericenter passage of the third body in seconds. This parameter is effectively a constant of integration.
t (3D double xp.ndarray) – Time values to evaluate. Shape is
(num binaries, 3 spacecraft, N)
.
- Returns:
LOS acceleration
- Return type:
3D double xp.ndarray
- get_f_derivatives(f0, fdot, fddot, A2, varpi, e2, P2, T2, t=None, eps=50000.0)
Get instantaneous frequency derivatives in third-body waveform
Computes the instantaneous frequency and frequency derivatives by calculating the effect of the third-body over the course of the orbit.
Central difference is used for both first and second derivatives.
- Parameters:
f0 (1D double np.ndarray) – Initial frequency of gravitational wave in Hz.
fdot (1D double np.ndarray) – Initial time derivative of the frequency given as Hz/s.
fddot (1D double np.ndarray) – Initial second derivative with respect to time of the frequency given in Hz/s^2.
A2 (1D double np.ndarray) – Special amplitude parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
varpi (1D double np.ndarray) – Special angular frequency parameter related to the line-of-site velocity for the third body orbit as defined in the paper given in the description above.
e2 (1D double np.ndarray) – Eccentricity of the third body orbit.
n2 (1D double np.ndarray) – Angular frequency of the third body orbit in per seconds.
T2 (1D double np.ndarray) – Time of pericenter passage of the third body in seconds. This parameter is effectively a constant of integration.
t (3D double xp.ndarray, optional) – Time values for derivative calculation. Shape is
(num binaries, anything, 3)
. The 3 here represents the times for each derivative computation. The derivatives are to be calucated at index 1. The step before (after) the time of the derivative (for central differencing) is at index 0 (1). In other words, the final dimension oft
should be[t_deriv - eps, t_deriv, t_deriv + eps]
. Default isNone
. IfNone
, will fill[t_deriv - eps, t_deriv, t_deriv + eps]
witht_deriv = 0.0
andeps
will be the kwarg value.eps (double) – Step size for central differencing. Only used if
t
is not provided.
- Returns:
Phase shifts due to third-body effect.
- Return type:
3D double xp.ndarray
Third-body utility functions
- gbgpu.thirdbody.third_body_factors(M, mc, P2, e2, iota, Omega2, omega2, phi2, lam, beta, third_mass_unit='Mjup', third_period_unit='yrs')
Get
A2,varpi,T2
from third-body parametersGet all the third-body factors that go into the waveform computation:
A2
,varpi
, andT2
.- Parameters:
M (double or double np.ndarray) – Total mass of inner Galactic binary in Solar Masses.
mc (double or double np.ndarray) – Mass of third body in units of
third_mass_unit
kwarg.P2 (double or double np.ndarray) – Orbital period of third body in units of
third_period_unit
kwarg.e2 (double or double np.ndarray) – Orbital eccentricity of third body.
iota (double or double np.ndarray) – Orbital inclination of third body in radians. This orbital inclination is one of the three euler angles describing the rotation of the third-body orbital frame to ecliptic frame.
Omega2 (double or double np.ndarray) –
The other two of three euler angles describing the third-body orbital frame rotation from the ecliptic frame. See the Figure 1 in arXiv:1806.00500.
omega2 (double or double np.ndarray) –
The other two of three euler angles describing the third-body orbital frame rotation from the ecliptic frame. See the Figure 1 in arXiv:1806.00500.
phi2 (double or double np.ndarray) – Orbital phase in radians.
lam (double or double np.ndarray) – Ecliptic longitude in radians.
beta (double or double np.ndarray) – Ecliptic latitude in radians.
third_mass_unit (str, optional) – Mass unit for third-body mass. Options are
"Mjup"
for Jupiter Masses or"MSUN
” for Solar Masses. Default is"Mjup"
.third_period_unit (str, optional) – Time unit for third-body period. Options are
"sec"
for seconds or"yrs
” for years. Default is"yrs"
.
- Returns:
(A2, varpi, T2)
associated with the input parameters.- Return type:
Tuple
- gbgpu.thirdbody.get_T2(P2, e2, phi2, third_period_unit='yrs')
Get
T2
from third-body parameters- Parameters:
P2 (double or double np.ndarray) – Orbital period of third body in units of
third_period_unit
kwarg.e2 (double or double np.ndarray) – Orbital eccentricity of third body.
phi2 (double or double np.ndarray) – Orbital phase in radians.
third_period_unit (str, optional) – Time unit for third-body period. Options are
"sec"
for seconds or"yrs
” for years. Default is"yrs"
.
- Returns:
T2
associated with the input parameters.- Return type:
double xp.ndarray