Utilities

Mode Filtering

class few.utils.modeselector.ModeSelector(amplitude_generator, ylm_generator=None, mode_selection=None, include_minus_mkn=True, mode_selection_threshold=1e-05, sensitivity_fn=None, modeinds_map=None, force_backend=None, **kwargs)

Bases: ParallelModuleBase

Filter teukolsky amplitudes based on power contribution.

This module generates teukolsky modes and their associated ylms given an input trajectory, then filters these modes according to either a list of requested modes or the the power contribution from each mode. Mode filtering is performed to (roughly) achieve the requested mismatch with the minimum number of modes. If a sensitivity curve is provided, mode filtering is weighted according to this sensitivity curve.

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

Parameters:
  • amplitude_generator (object) – Object that generates the teukolsky amplitudes for the modes. This should be an instance of a class that has a __call__ method that takes the parameters \((a, p, e, xI)\) and returns the teukolsky amplitudes for the modes. It should also have mode-related attributes obtained from subclassing few.utils.baseclasses.SphericalHarmonic.

  • ylm_generator (object | None) – Object that generates the Ylm values for the modes. This should be an instance of a class that has a __call__ method that takes the parameters \((l, m, \theta, \phi)\) and returns the Ylm values for the modes. Default is few.utils.ylm.GetYlms, which generates Ylm values for the modes based on the l and m indices.

  • mode_selection (str | list | ndarray | None) – Determines the type of mode filtering to perform. If None, use default mode filtering provided by mode_selector. If ‘all’, it will run all modes without filtering. If ‘threshold’ it will override other options to filter by the threshold value set by mode_selection_threshold. If a list of tuples (or lists) of mode indices (e.g. [(\(l_1,m_1,k_1,n_1\)), (\(l_2,m_2,k_2,n_2\))]) is provided, it will return those modes combined into a single waveform. Default is None.

  • include_minus_mkn (bool | None) – If True, then include \((-m, -k, -n)\) mode when computing a \((m, k, n)\) mode. This only affects modes if mode_selection is a list of specific modes when the class is called. If True, this list of modes provided at call time must only contain \(m\geq 0\). Default is True.

  • mode_selection_threshold (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.

  • 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.

  • modeinds_map (ndarray | None) – Map of mode indices to Teukolsky amplitude data. This is a 4D array of shape (l, m, k, n) that maps the mode indices to one-dimensional indices from the amplitude module output. This is used to efficiently select the modes from the amplitude module output. By default this is obtained from the supplied amplitude module.

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

  • force_backend (str | Backend | None)

l_arr

l-mode indices for each mode index, given by amplitude module.

Type:

array

m_arr

m-mode indices for each mode index, given by amplitude module.

Type:

array

k_arr

k-mode indices for each mode index, given by amplitude module.

Type:

array

n_arr

n-mode indices for each mode index, given by amplitude module.

Type:

array

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

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

classmethod all_citations()

Return all the citations from the registry as printable BibTeX string

Return type:

str

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:

Iterable[REFERENCE | str]

property xp: ModuleType

Return the module providing ndarray capabilities

include_minus_mkn

Whether to include modes with m < 0 in the output.

Type:

bool

negative_m_flag

Specifies whether there are negative m-modes in mode_selection.

Type:

bool

mode_selection_threshold

Target mismatch threshold for mode selection. Modes will be selected such that the mismatch is approximately equal to this value.

Type:

float

classmethod supported_backends()

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

__call__(t, a, p, e, xI, theta, phi, online_mode_selection_args=None, mode_selection=None, include_minus_mkn=None, mode_selection_threshold=None, return_sort_inds=False)

Call to sort and filer teukolsky modes.

This is the call function that takes an inspiral trajectory and returns the teukolsky modes and ylms required to produce a waveform of either a given mismatch or the requested mode content.

Parameters:
  • t – Time array for the input trajectory.

  • a – Dimensionless spin parameter of the primary black hole.

  • p – Semi-latus rectum of the trajectory.

  • e – Eccentricity of the trajectory.

  • xI – Initial cosine(inclination) for the trajectory.

  • theta – Polar source-frame viewing angle.

  • phi – Azimuthal source-frame viewing angle.

  • online_mode_selection_args (dict | None) – Dictionary of arguments necessary to determine the mode frequencies along the trajectory. This should contain the keys ‘f_phi’, ‘f_theta’, and ‘f_r’, which are the fundamental frequencies of the modes along the trajectory.

  • mode_selection (str | list | ndarray | None) – Determines the type of mode filtering to perform. If None, use default mode filtering provided by mode_selector. If ‘all’, it will run all modes without filtering. If ‘threshold’ it will override other options to filter by the threshold value set by mode_selection_threshold. If a list of tuples (or lists) of mode indices (e.g. [(\(l_1,m_1,k_1,n_1\)), (\(l_2,m_2,k_2,n_2\))]) is provided, it will return those modes combined into a single waveform. If include_minus_mkn = True, we require that \(m \geq 0\) for this list. Default is None.

  • modeinds_map – Map of mode indices to Teukolsky amplitude data. Only required if mode_selection is a list of specific mode. Default is None.

  • include_minus_mkn (bool | None) – If True, then include \((-m, -k, -n)\) mode when computing a \((m, k, n)\) mode. This only affects modes if mode_selection is a list of specific modes. Default is True.

  • mode_selection_threshold (float) – Target waveform mismatch used to determine the contributing modes. Lowering this value will calculate more modes, slowing the waveform generation 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_sort_inds (bool) – If True, also return the indices sorting the modes according to their contribution. Only used when filtering in this mode. Default is False.

Return type:

tuple[ndarray]

few.utils.modeselector.get_selected_modes_from_initial_conditions(mode_selector_module, traj_module, m1, m2, a, p0, e0, xI0, theta, phi, traj_args=None, traj_kwargs=None, mode_selector_kwargs=None)

Get the selected modes from the initial conditions of an EMRI trajectory. This function uses the provided trajectory module to generate the trajectory and the mode selector module to select the modes based on the trajectory.

Parameters:
  • mode_selector_module (ModeSelector) – Instance of the ModeSelector class.

  • traj_module (EMRIInspiral) – Instance of the EMRIInspiral class.

  • m1 (float) – Mass of the primary black hole in solar masses.

  • m2 (float) – Mass of the secondary black hole in solar masses.

  • a (float) – Dimensionless spin parameter of the primary black hole.

  • p0 (float) – Initial semi-latus rectum of the trajectory.

  • e0 (float) – Initial eccentricity of the trajectory.

  • xI0 (float) – Initial cosine(inclination) for the trajectory.

  • theta (float) – Polar source-frame viewing angle.

  • phi (float) – Azimuthal source-frame viewing angle.

  • traj_args (list | None) – Additional arguments to pass to the trajectory module.

  • traj_kwargs (dict | None) – Additional keyword arguments to pass to the trajectory module.

  • mode_selector_kwargs (dict | None) – Additional keyword arguments to pass to the mode selector module.

Returns:

A dict containing the selected teuk_modes, ylms, their indices (ls, ms, ks, ns), the trajectory traj and the online_mode_selection_args.

Return type:

dict

(-2) Spin-Weighted Spherical Harmonics

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

Bases: ParallelModuleBase

(-2) Spin-weighted Spherical Harmonics

The class generates (-2) spin-weighted spherical harmonics, \(Y_{lm}(\Theta,\phi)\).

Parameters:
  • include_minus_m (bool) – Set True if only providing \(m\geq0\), it will return twice the number of requested modes with the second half as modes with \(m<0\) for array inputs of \(l,m\). 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, include_minus_m=None)

Call method for Ylms.

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

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

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

  • theta (float) – Polar viewing angle.

  • phi (float) – Azimuthal viewing angle.

  • include_minus_m (bool | None)

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

classmethod all_citations()

Return all the citations from the registry as printable BibTeX string

Return type:

str

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:

Iterable[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.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, lower_bound_buffer=1e-06, upper_bound_maximum=1000000.0, **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 the minimum and maximum values of p for the trajectory modules.

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

  • lower_bound_buffer (float) – A float that offsets the lower bound by a small number to prevent starting on the separatrix If not provided, it will default to 1e-6

  • upper_bound_maximum (float) – A float that sets the maximum value of the upper bound to a large finite number if max_p returns inf If not provided, it will default to 1e6

Returns:

Value of p that creates the proper length trajectory.

Return type:

float

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

Find the value of m2 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 m2 accordingly. This function tells you what that value of m2 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 m2 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_m2 (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 m2 over which brentq will search for a root. If not given, will be set to [1e-1, m1]. To supply only one of these two limits, set the other limit to None.

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

Returns:

Value of m2 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.mappings.pn.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.mappings.pn.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