Utilities

Mode Filtering

class few.utils.modeselector.ModeSelector(l_arr, m_arr, n_arr, sensitivity_fn=None, force_backend=None)

Bases: ParallelModuleBase

Filter teukolsky amplitudes based on power contribution.

This module takes teukolsky modes, combines them with their associated ylms, and determines the power contribution from each mode. It then filters the modes bases on the fractional accuracy on the total power (eps) parameter. Additionally, if a sensitivity curve is provided, the mode power is also weighted according to the PSD of the sensitivity.

The mode filtering is a major contributing factor to the speed of these waveforms as it removes large numbers of useles modes from the final summation calculation.

Be careful as this is built based on the construction that input mode arrays will in order of \(m=0\), \(m>0\), and then \(m<0\).

Parameters:
  • l_arr (ndarray) – The l-mode indices for each mode index.

  • m_arr (ndarray) – The m-mode indices for each mode index.

  • n_arr (ndarray) – The n-mode indices for each mode index.

  • sensitivity_fn (object | None) – Sensitivity curve function that takes a frequency (Hz) array as input and returns the Power Spectral Density (PSD) of the sensitivity curve. Default is None. If this is not none, this sennsitivity is used to weight the mode values when determining which modes to keep. Note: if the sensitivity function is provided, and GPUs are used, then this function must accept CuPy arrays as input.

  • **kwargs – Optional keyword arguments for the base class: few.utils.baseclasses.ParallelModuleBase.

  • force_backend (str | Backend | None)

m0mask

Mask for m=0 modes.

Type:

array

num_m_zero_up

Number of modes with \(m\geq0\).

Type:

int

num_m_1_up

Number of modes with \(m\geq1\).

Type:

int

num_m0

Number of modes with \(m=0\).

Type:

int

sensitivity_fn

sensitivity generating function for power-weighting.

Type:

object

classmethod supported_backends()

List of backends supported by a parallel module by order of preference.

property is_predictive

Whether this mode selector should be used before or after amplitude generation

__call__(teuk_modes, ylms, modeinds, fund_freq_args=None, eps=1e-05)

Call to sort and filer teukolsky modes.

This is the call function that takes the teukolsky modes, ylms, mode indices and fractional accuracy of the total power and returns filtered teukolsky modes and ylms.

Parameters:
  • teuk_modes (ndarray) – Complex teukolsky amplitudes from the amplitude modules. Shape: (number of trajectory points, number of modes).

  • ylms (ndarray) – Array of ylm values for each mode, including m<0. Shape is (num of m==0,) + (num of m>0,) + (num of m<0). Number of m<0 and m>0 is the same, but they are ordered as (m==0) first then m>0 then m<0.

  • modeinds (list[ndarray]) – List containing the mode index arrays. If in an equatorial model, need \((l,m,n)\) arrays. If generic, \((l,m,k,n)\) arrays. e.g. [l_arr, m_arr, n_arr].

  • fund_freq_args (tuple | None) – Args necessary to determine fundamental frequencies along trajectory. The tuple will represent \((M, a, e, p, \cos\iota)\) where the large black hole mass (\(M\)) and spin (\(a\)) are scalar and the other three quantities are self.xp.ndarrays. This must be provided if sensitivity weighting is used. Default is None.

  • eps (float) – Fractional accuracy of the total power used to determine the contributing modes. Lowering this value will calculate more modes slower the waveform down, but generally improving accuracy. Increasing this value removes modes from consideration and can have a considerable affect on the speed of the waveform, albeit at the cost of some accuracy (usually an acceptable loss). Default that gives good mismatch qualities is 1e-5.

Return type:

ndarray

static CPU_ONLY()

List of supported backend for CPU only class

Return type:

list[str]

List of supported backends for CPU-recommended class with GPU support

Return type:

list[str]

static GPU_ONLY()

List of supported backends for GPU-only class

Return type:

list[str]

List of supported backends for GPU-recommended class with CPU support

Return type:

list[str]

adapt_backend_kwargs(kwargs=None)

Adapt a set of keyword arguments to add/set ‘force_backend’ to current backend

Parameters:

kwargs (dict | None)

Return type:

dict

property backend: Backend

Access the underlying backend.

property backend_name: str

Return the name of current backend

build_with_same_backend(module_class, args=None, kwargs=None)

Build an instance of module_class with same backend as current object.

Parameters:
  • module_class (type[ParallelModuleDerivate]) – class of the object to be built, must derive from ParallelModuleBase

  • args (list, optional) – positional arguments for module_class constructor

  • kwargs (dict, optional) – keyword arguments for module_class constructor (the ‘force_backend’ argument will be ignored and replaced)

Return type:

ParallelModuleDerivate

classmethod citation()

Return the module references as a printable BibTeX string.

Return type:

str

classmethod module_references()

Method implemented by each class to define its list of references

Return type:

List[REFERENCE | str]

property xp: ModuleType

Return the module providing ndarray capabilities

class few.utils.modeselector.NeuralModeSelector(l_arr, m_arr, n_arr, threshold=0.5, mode_selector_location=None, return_type='tuples', keep_inds=None, maximise_over_theta=False, force_backend=None)

Bases: ParallelModuleBase

Filter teukolsky amplitudes based on power contribution.

This module uses a combination of a pre-computed mask and a feed-forward neural network to predict the mode content of the waveform given its parameters. Therefore, the results will vary compared to manually computing the mode selection, but it should be very accurate especially for the stronger modes.

The mode filtering is a major contributing factor to the speed of these waveforms as it removes large numbers of useles modes from the final summation calculation.

Parameters:
  • l_arr (1D int self.xp.ndarray) – The l-mode indices for each mode index.

  • m_arr (1D int self.xp.ndarray) – The m-mode indices for each mode index.

  • n_arr (1D int self.xp.ndarray) – The n-mode indices for each mode index.

  • threshold (double) – The network threshold value for mode retention. Decrease to keep more modes, minimising missed modes but slowing down the waveform computation. Defaults to 0.5 (the optimal value for accuracy).

  • **kwargs (dict, optional) – Keyword arguments for the base classes: few.utils.baseclasses.ParallelModuleBase. Default is {}.

  • force_backend (str | Backend | None)

static CPU_ONLY()

List of supported backend for CPU only class

Return type:

list[str]

List of supported backends for CPU-recommended class with GPU support

Return type:

list[str]

static GPU_ONLY()

List of supported backends for GPU-only class

Return type:

list[str]

List of supported backends for GPU-recommended class with CPU support

Return type:

list[str]

adapt_backend_kwargs(kwargs=None)

Adapt a set of keyword arguments to add/set ‘force_backend’ to current backend

Parameters:

kwargs (dict | None)

Return type:

dict

property backend: Backend

Access the underlying backend.

property backend_name: str

Return the name of current backend

build_with_same_backend(module_class, args=None, kwargs=None)

Build an instance of module_class with same backend as current object.

Parameters:
  • module_class (type[ParallelModuleDerivate]) – class of the object to be built, must derive from ParallelModuleBase

  • args (list, optional) – positional arguments for module_class constructor

  • kwargs (dict, optional) – keyword arguments for module_class constructor (the ‘force_backend’ argument will be ignored and replaced)

Return type:

ParallelModuleDerivate

classmethod citation()

Return the module references as a printable BibTeX string.

Return type:

str

classmethod module_references()

Method implemented by each class to define its list of references

Return type:

List[REFERENCE | str]

property xp: ModuleType

Return the module providing ndarray capabilities

classmethod supported_backends()

List of backends supported by a parallel module by order of preference.

property is_predictive

Whether this mode selector should be used before or after amplitude generation

__call__(M, mu, a, p0, e0, xI, theta, phi, T, eps)

Call to predict the mode content of the waveform.

This is the call function that takes the waveform parameters, applies a precomputed mask and then evaluates the remaining mode content with a neural network classifier.

Parameters:
  • M (double) – Mass of larger black hole in solar masses.

  • mu (double) – Mass of compact object in solar masses.

  • p0 (double) – Initial semilatus rectum (Must be greater than the separatrix at the the given e0 and x0). See documentation for more information on \(p_0<10\).

  • e0 (double) – Initial eccentricity.

  • theta (double) – Polar viewing angle.

  • phi (double) – Azimuthal viewing angle.

  • T (double) – Duration of waveform in years.

  • eps (double) – Mode selection threshold power.

(-2) Spin-Weighted Spherical Harmonics

class few.utils.ylm.GetYlms(assume_positive_m=False, **kwargs)

Bases: ParallelModuleBase

(-2) Spin-weighted Spherical Harmonics

The class generates (-2) spin-weighted spherical harmonics, \(Y_{lm}(\Theta,\phi)\). Important Note: this class also applies the parity operator (\(-1^l\)) to modes with \(m<0\).

Parameters:
  • assume_positive_m (bool) – Set true if only providing \(m\geq0\), it will return twice the number of requested modes with the seconds half as modes with \(m<0\). Warning: It will also duplicate the \(m=0\) modes. Default is False.

  • **kwargs (dict | None) – Optional keyword arguments for the base class: few.utils.baseclasses.ParallelModuleBase.

classmethod supported_backends()

List of backends supported by a parallel module by order of preference.

__call__(l_in, m_in, theta, phi)

Call method for Ylms.

This returns ylms based on requested \((l,m)\) values and viewing angles.

Parameters:
  • l_in (ndarray) – \(l\) values requested.

  • m_in (ndarray) – \(m\) values requested.

  • theta (float) – Polar viewing angle.

  • phi (float) – Azimuthal viewing angle.

Returns:

Complex array of Ylm values.

Return type:

ndarray

static CPU_ONLY()

List of supported backend for CPU only class

Return type:

list[str]

List of supported backends for CPU-recommended class with GPU support

Return type:

list[str]

static GPU_ONLY()

List of supported backends for GPU-only class

Return type:

list[str]

List of supported backends for GPU-recommended class with CPU support

Return type:

list[str]

adapt_backend_kwargs(kwargs=None)

Adapt a set of keyword arguments to add/set ‘force_backend’ to current backend

Parameters:

kwargs (dict | None)

Return type:

dict

property backend: Backend

Access the underlying backend.

property backend_name: str

Return the name of current backend

build_with_same_backend(module_class, args=None, kwargs=None)

Build an instance of module_class with same backend as current object.

Parameters:
  • module_class (type[ParallelModuleDerivate]) – class of the object to be built, must derive from ParallelModuleBase

  • args (list, optional) – positional arguments for module_class constructor

  • kwargs (dict, optional) – keyword arguments for module_class constructor (the ‘force_backend’ argument will be ignored and replaced)

Return type:

ParallelModuleDerivate

classmethod citation()

Return the module references as a printable BibTeX string.

Return type:

str

classmethod module_references()

Method implemented by each class to define its list of references

Return type:

List[REFERENCE | str]

property xp: ModuleType

Return the module providing ndarray capabilities

Frequency Domain Utilities

few.utils.fdutils.get_convolution(a, b)

Calculate the convolution of two arrays.

Parameters:
  • a (ndarray) – First array to convolve.

  • b (ndarray) – First array to convolve.

Returns:

convolution of the two arrays.

Return type:

ndarray

few.utils.fdutils.get_fft_td_windowed(signal, window, dt)

Calculate the Fast Fourier Transform of a windowed time domain signal.

Parameters:
  • signal (list) – A length-2 list containing the signals plus and cross polarizations.

  • window (ndarray) – Array to be applied in time domain to each signal.

  • dt (float) – Time sampling interval of the signal and window.

Returns:

Fast Fourier Transform of the windowed time domain signals.

Return type:

list[ndarray]

few.utils.fdutils.get_fd_windowed(signal, window=None, window_in_fd=False)

Calculate the convolution of a frequency domain signal with a window in time domain.

Parameters:
  • signal (list) – A length-2 list containing the signals plus and cross polarizations in frequency domain.

  • window (bool | None) – Array of the time domain window. If None, do not apply window. This is added for flexibility. Default is None.

  • window_in_fd (bool) – If True, window is given in the frequency domain. If False, window is given in the time domain. Default is False.

Returns:

convolution of a frequency domain signal with a time domain window.

Return type:

list[ndarray]

class few.utils.fdutils.GetFDWaveformFromFD(waveform_generator, positive_frequency_mask, dt, non_zero_mask=None, window=None, window_in_fd=False)

Bases: object

Generic frequency domain class

This class allows to obtain the frequency domain signal given the frequency domain waveform class from the FEW package.

Parameters:
  • waveform_generator (object) – FEW waveform class.

  • positive_frequency_mask (ndarray) – boolean array to indicate where the frequencies are positive.

  • dt (float) – time sampling interval of the signal and window.

  • non_zero_mask (ndarray | None) – boolean array to indicate where the waveform needs to be set to zero.

  • window (ndarray | None) – Array of the time domain window. If None, do not apply window. This is added for flexibility. Default is None.

  • window_in_fd (bool | None) – If True, window is given in the frequency domain. If False, window is given in the time domain. Default is False.

__call__(*args, **kwargs)

Run the waveform generator.

Parameters:
  • *args – Arguments passed to waveform generator.

  • **kwargs – Keyword arguments passed to waveform generator.

Returns:

FD Waveform as [h+, hx]

Return type:

list[ndarray]

class few.utils.fdutils.GetFDWaveformFromTD(waveform_generator, positive_frequency_mask, dt, non_zero_mask=None, window=None)

Bases: object

Generic time domain class

This class allows to obtain the frequency domain signal given the time domain waveform class from the FEW package.

Parameters:
  • waveform_generator (object) – FEW waveform class.

  • positive_frequency_mask (ndarray) – boolean array to indicate where the frequencies are positive.

  • dt (float) – time sampling interval of the signal and window.

  • non_zero_mask (ndarray | None) – boolean array to indicate where the waveform needs to be set to zero.

  • window (ndarray | None) – Array of the time domain window. If None, do not apply window. This is added for flexibility. Default is None.

__call__(*args, **kwargs)

Run the waveform generator.

Parameters:
  • *args – Arguments passed to waveform generator.

  • **kwargs – Keyword arguments passed to waveform generator.

Returns:

FD Waveform as [h+, hx] (fft from TD)

Return type:

list[ndarray]

Analysis and Other Tools

few.utils.utility.get_overlap(time_series_1, time_series_2, use_gpu=False)

Calculate the overlap.

Takes two time series and finds which one is shorter in length. It then shortens the longer time series if necessary. Then it performs a normalized correlation calulation on the two time series to give the overlap. The overlap of \(a(t)\) and \(b(t)\), \(\gamma_{a,b}\), is given by,

\[\gamma_{a,b} = <a,b>/(<a,a><b,b>)^{(1/2)},\]

where \(<a,b>\) is the inner product of the two time series.

Parameters:
  • time_series_1 (ndarray) – Strain time series 1.

  • time_series_2 (ndarray) – Strain time series 2.

  • use_gpu (bool) – If True use cupy. If False, use numpy. Default is False.

Return type:

float

few.utils.utility.get_mismatch(time_series_1, time_series_2, use_gpu=False)

Calculate the mismatch.

The mismatch is 1 - overlap. Therefore, see documentation for few.utils.utility.overlap() for information on the overlap calculation.

Parameters:
  • time_series_1 (ndarray) – Strain time series 1.

  • time_series_2 (ndarray) – Strain time series 2.

  • use_gpu (bool) – If True use cupy. If False, use numpy. Default is False.

Return type:

float

few.utils.utility.ELQ_to_pex(a, E, Lz, Q)

Convert from Kerr constants of motion to orbital elements.

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • E (float | ndarray) – Values of energy, \(E\).

  • Lz (float | ndarray) – Values of angular momentum, \(L_z\).

  • Q (float | ndarray) – Values of the Carter constant, \(Q\).

Returns:

Tuple of (OmegaPhi, OmegaTheta, OmegaR). These are 1D arrays or scalar values depending on inputs.

Return type:

tuple[float | ndarray]

few.utils.utility.get_fundamental_frequencies(a, p, e, x, use_gpu=False)

Get dimensionless fundamental frequencies.

Determines fundamental frequencies in generic Kerr from Schmidt 2002.

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • p (float | ndarray) – Values of separation, \(p\).

  • e (float | ndarray) – Values of eccentricity, \(e\).

  • x (float | ndarray) – Values of cosine of the inclination, \(x=\cos{I}\). Please note this is different from \(Y=\cos{\iota}\).

  • use_gpu (bool)

Returns:

Tuple of (OmegaPhi, OmegaTheta, OmegaR). These are 1D arrays or scalar values depending on inputs.

Return type:

tuple[float | ndarray]

few.utils.utility.get_fundamental_frequencies_spin_corrections(a, p, e, x)

Get the leading-order correction term to the fundamental frequencies due to the spin of the secondary compact object.

Currently only supported for equatorial orbits.

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • p (float | ndarray) – Values of separation, \(p\).

  • e (float | ndarray) – Values of eccentricity, \(e\).

  • x (float | ndarray) – Values of cosine of the inclination, \(x=\cos{I}\). Please note this is different from \(Y=\cos{\iota}\).

Returns:

Tuple of (OmegaPhi, OmegaTheta, OmegaR). These are 1D arrays or scalar values depending on inputs.

Return type:

tuple[float | ndarray]

few.utils.utility.get_kerr_geo_constants_of_motion(a, p, e, x)

Get Kerr constants of motion.

Determines the constants of motion: \((E, L, Q)\) associated with a geodesic orbit in the generic Kerr spacetime.

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • p (float | ndarray) – Values of separation, \(p\).

  • e (float | ndarray) – Values of eccentricity, \(e\).

  • x (float | ndarray) – Values of cosine of the inclination, \(x=\cos{I}\). Please note this is different from \(Y=\cos{\iota}\).

Returns:

Tuple of (E, L, Q). These are 1D arrays or scalar values depending on inputs.

Return type:

tuple

few.utils.utility.get_separatrix(a, e, x, tol=1e-13, use_gpu=False)

Get separatrix in generic Kerr.

Determines separatrix in generic Kerr from Stein & Warburton 2020.

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • e (float | ndarray) – Values of eccentricity, \(e\).

  • x (float | ndarray) – Values of cosine of the inclination, \(x=\cos{I}\). Please note this is different from \(Y=\cos{\iota}\).

  • tol (float) – Tolerance for root-finding. Default is 1e-13.

  • use_gpu (bool)

Returns:

Separatrix value with shape based on input shapes.

Return type:

float | ndarray

few.utils.utility.get_at_t(traj_module, traj_args, bounds, t_out, index_of_interest, traj_kwargs=None, xtol=2e-12, rtol=8.881784197001252e-16)

Root finding wrapper using Brent’s method.

This function uses scipy’s brentq routine to find root.

Parameters:
  • traj_module (object) – Instantiated trajectory module. It must output the time array of the trajectory sparse trajectory as the first output value in the tuple.

  • traj_args (list[float]) – List of arguments for the trajectory function. p is removed. Note: It must be a list, not a tuple because the new p values are inserted into the argument list.

  • bounds (list[float]) – Minimum and maximum values over which brentq will search for a root.

  • t_out (float) – The desired length of time for the waveform.

  • index_of_interest (int) – Index where to insert the new values in the traj_args list.

  • traj_kwargs (dict | None) – Keyword arguments for traj_module. Default is an empty dict.

  • xtol (float) – Absolute tolerance of the brentq root-finding - see :code: np.allclose() for details. Defaults to 2e-12 (scipy default).

  • rtol (float) – Relative tolerance of the brentq root-finding - see :code: np.allclose() for details. Defaults to ~8.8e-16 (scipy default).

Returns:

Root value.

Return type:

float

few.utils.utility.get_p_at_t(traj_module, t_out, traj_args, index_of_p=3, index_of_a=2, index_of_e=4, index_of_x=5, bounds=None, **kwargs)

Find the value of p that will give a specific length inspiral using Brent’s method.

If you want to generate an inspiral that is a specific length, you can adjust p accordingly. This function tells you what that value of p is based on the trajectory module and other input parameters at a desired time of observation.

This function uses scipy’s brentq routine to find the (presumed only) value of p that gives a trajectory of duration t_out.

Parameters:
  • traj_module (object) – Instantiated trajectory module. It must output the time array of the trajectory sparse trajectory as the first output value in the tuple.

  • t_out (float) – The desired length of time for the waveform.

  • traj_args (list[float]) – List of arguments for the trajectory function. p is removed. Note: It must be a list, not a tuple because the new p values are inserted into the argument list.

  • index_of_p (int) – Index where to insert the new p values in the traj_args list. Default is 3.

  • index_of_a (int) – Index of a in provided traj_module arguments. Default is 2.

  • index_of_e (int) – Index of e0 in provided traj_module arguments. Default is 4.

  • index_of_x (int) – Index of x0 in provided traj_module arguments. Default is 5.

  • bounds (list[float | None]) – Minimum and maximum values of p over which brentq will search for a root. If not given, will be set to [separatrix + 0.101, 50]. To supply only one of these two limits, set the other limit to None.

  • **kwargs – Keyword arguments for get_at_t().

Returns:

Value of p that creates the proper length trajectory.

Return type:

float

few.utils.utility.get_mu_at_t(traj_module, t_out, traj_args, index_of_mu=1, bounds=None, **kwargs)

Find the value of mu that will give a specific length inspiral using Brent’s method.

If you want to generate an inspiral that is a specific length, you can adjust mu accordingly. This function tells you what that value of mu is based on the trajectory module and other input parameters at a desired time of observation.

This function uses scipy’s brentq routine to find the (presumed only) value of mu that gives a trajectory of duration t_out.

Parameters:
  • traj_module (object) – Instantiated trajectory module. It must output the time array of the trajectory sparse trajectory as the first output value in the tuple.

  • t_out (float) – The desired length of time for the waveform.

  • traj_args (list[float]) – List of arguments for the trajectory function. p is removed. Note: It must be a list, not a tuple because the new p values are inserted into the argument list.

  • index_of_mu (int) – Index where to insert the new p values in the traj_args list. Default is 1.

  • bounds (list[float | None]) – Minimum and maximum values of p over which brentq will search for a root. If not given, will be set to [1e-1, 1e3]. To supply only one of these two limits, set the other limit to None.

  • **kwargs – Keyword arguments for get_at_t().

Returns:

Value of mu that creates the proper length trajectory.

Return type:

float

few.utils.utility.wrapper(*args, **kwargs)

Function to convert array and C/C++ class arguments to ptrs

This function checks the object type. If it is a cupy or numpy array, it will determine its pointer by calling the proper attributes. If you design a Cython class to be passed through python, it must have a ptr attribute.

If you use this function, you must convert input arrays to size_t data type in Cython and then properly cast the pointer as it enters the c++ function. See the Cython codes here for examples.

Parameters:
  • *args (list) – list of the arguments for a function.

  • **kwargs (dict) – dictionary of keyword arguments to be converted.

Returns:

(targs, tkwargs) where t indicates target (with pointer values

rather than python objects).

Return type:

Tuple

few.utils.utility.pointer_adjust(func)

Decorator function for cupy/numpy agnostic cython

This decorator applies few.utils.utility.wrapper() to functions via the decorator construction.

If you use this decorator, you must convert input arrays to size_t data type in Cython and then properly cast the pointer as it enters the c++ function. See the Cython codes here for examples.

Elliptic integrals (Legendre and Carlson symmetric forms)

Efficient implementations of complete elliptic integrals, obtained by applying the duplication theorem to compute Carlson’s symmetric forms.

Derived from Carlson’s algorithms in the SLATEC library.

few.utils.elliptic.RF(x, y, z, tol=0.0005)

Computes the Carlson symmetric form of the elliptic integral of the first kind, \(R_F\).

Parameters:
  • x (float) – First argument, \(x >= 0\).

  • y (float) – Second argument, \(y >= 0\).

  • z (float) – Third argument, \(z >= 0\).

  • tol (float) – Numerical tolerance parameter (defaults to 5e-4)

Returns:

The Carlson symmetric form \(R_F\).

Return type:

float

few.utils.elliptic.RC(x, y)

Computes the reduced Carlson symmetric form of the elliptic integral of the first kind, \(R_C\).

Unlike the other Carlson symmetric forms provided, this one can be computed analytically in terms of transcendental functions.

Parameters:
  • x (float) – First argument, \(x >= 0\).

  • y (float) – Second argument, \(y >= 0\).

Returns:

The Carlson symmetric form \(R_J\).

Return type:

float

few.utils.elliptic.RJ(x, y, z, p, tol=0.0005)

Computes the Carlson symmetric form of the elliptic integral of the second kind, \(R_J\).

No more than one of \((x, y, z)\) may be equal to zero.

Parameters:
  • x (float) – First argument, \(x >= 0\).

  • y (float) – Second argument, \(y >= 0\).

  • z (float) – Third argument, \(z >= 0\).

  • p (float) – Fourth argument, \(p > 0\).

  • tol (float) – Numerical tolerance parameter (defaults to 5e-4)

Returns:

The Carlson symmetric form \(R_J\).

Return type:

float

few.utils.elliptic.RD(x, y, z, tol=0.0005)

Computes the reduced Carlson symmetric form of the elliptic integral of the second kind, \(R_D\).

It is a requirement that \(x+y > 0\).

Parameters:
  • x (float) – First argument, \(x >= 0\).

  • y (float) – Second argument, \(y >= 0\).

  • z (float) – Third argument, \(z > 0\).

  • tol (float) – Numerical tolerance parameter (defaults to 5e-4)

Returns:

The Carlson symmetric form \(R_D\).

Return type:

float

few.utils.elliptic.EllipK(k, tol=0.0005)

Computes the complete elliptic integral of the first kind, \(K(k)\).

Switches to a polynomial approximation from Abramowitz & Stegun for \(k > 1 - 10^{-10}\).

Parameters:
  • k (float) – The elliptic modulus, where k in [0, 1].

  • tol (float) – Numerical tolerance parameter (defaults to 5e-4)

Returns:

The complete elliptic integral \(K(k)\).

Return type:

float

few.utils.elliptic.EllipE(k, tol=0.0005)

Computes the complete elliptic integral of the second kind, \(E(k)\).

Switches to a polynomial approximation from Abramowitz & Stegun for \(k > 1 - 10^{-10}\).

Parameters:
  • k (float) – The elliptic modulus, where k in [0, 1].

  • tol (float) – Numerical tolerance parameter (defaults to 5e-4)

Returns:

The complete elliptic integral \(E(k)\).

Return type:

float

few.utils.elliptic.EllipPi(n, k, tol=0.0005)

Computes the complete elliptic integral of the third kind, \(\Pi(n, k)\).

Parameters:
  • n (float) – The characteristic.

  • k (float) – The elliptic modulus, where k in [0, 1].

  • tol (float) – Numerical tolerance parameter (defaults to 5e-4)

Returns:

The complete elliptic integral \(\Pi(n, k)\).

Return type:

float

PN Parameter Mapping Utilities

few.utils.pn_map.Y_to_xI(a, p, e, Y)

Convert from \(Y=\cos{\iota}\) to \(x_I=\cos{I}\).

Converts between the two different inclination parameters. \(\cos{I}\equiv x_I\), where \(I\) describes the orbit’s inclination from the equatorial plane. \(\cos{\iota}\equiv Y\), where \(\cos{\iota}=L/\sqrt{L^2 + Q}\).

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • p (float | ndarray) – Values of separation, \(p\).

  • e (float | ndarray) – Values of eccentricity, \(e\).

  • Y (float | ndarray) – Values of PN cosine of inclination \(Y\).

Return type:

float | ndarray

returns: \(x=\cos{I}\) values with shape based on input shapes.

few.utils.pn_map.xI_to_Y(a, p, e, x)

Convert from \(x_I=\cos{I}\) to \(Y=\cos{\iota}\).

Converts between the two different inclination parameters. \(\cos{I}\equiv x_I\), where \(I\) describes the orbit’s inclination from the equatorial plane. \(\cos{\iota}\equiv Y\), where \(\cos{\iota}=L/\sqrt{L^2 + Q}\).

Parameters:
  • a (float | ndarray) – Dimensionless spin of massive black hole. If other parameters are arrays and the spin is scalar, it will be cast to a 1D array.

  • p (float | ndarray) – Values of separation, \(p\).

  • e (float | ndarray) – Values of eccentricity, \(e\).

  • x (float | ndarray) – Values of cosine of the inclination, \(x=\cos{I}\).

Returns:

\(Y\) values with shape based on input shapes.

Return type:

float | ndarray