scri

Module for operating on gravitational waveforms in various forms

Classes

WaveformBaseBase class

This is probably not needed directly; it is just used for inheritance by other objects.

WaveformModes: Complex spin-weighted spherical-harmonic modes

The modes must include all m values for a range of ell values. This is the “classic” version of a WaveformBase object we might normally think of.

WaveformGrid: Complex quantity evaluated along world lines of grid points on the sphere

To perform translations or boosts, we need to transform to physical space, along a series of selected world lines distributed evenly across the sphere. These values may need to be interpolated to new time values, and they will presumably need to be transformed back to WaveformModes.

WaveformInDetector: Real quantities as observed in an inertial detector

Detectors only measure one polarization, so they deal with real quantities. Also, data is measured in evenly spaced time steps. This object can be created from a WaveformModes object.

WaveformInDetectorFT: (Complex) Fourier transform of a WaveformInDetector

This contains only the positive-frequency values since the transformed data is real.

Functions

version_info()

Show version information about this module and various dependencies

class scri.WaveformGrid(*args, **kwargs)[source]
Attributes:
abs
arg
arg_unwrapped
conformal_weight
data_2d
data_ddot
data_dot
data_iint
data_int
data_type_latex
data_type_string
descriptor_string

Create a simple string describing the content of the waveform

frame_type_string
gamma_weight

Non-conformal effect of a boost.

is_valid
m_scaling
n_data_sets
n_phi
n_theta
n_times
num
r_scaling
spin_weight

Methods

SI_units(current_unit_mass_in_solar_masses)

Assuming current quantities are in geometric units, convert to SI units

compare(w_a[, min_time_step, min_time])

Return a waveform with differences between the two inputs

copy()

Return a (deep) copy of the object

copy_without_data()

Return a copy of the object, with empty t, frame, and data fields

deepcopy()

Return a deep copy of the object

ensure_validity([alter, assertions])

Try to ensure that the WaveformGrid object is valid

from_modes(w_modes, **kwargs)

Construct grid object from modes, with optional BMS transformation

interpolate(tprime)

Interpolate the frame and data onto the new set of time steps

max_norm_index([skip_fraction_of_data])

Index of time step with largest norm

max_norm_time([skip_fraction_of_data])

Return time at which largest norm occurs in data

norm([take_sqrt, indices])

L2 norm of the waveform

to_modes([ell_max, ell_min])

Transform to modes of a spin-weighted spherical harmonic expansion

transform(w_modes, **kwargs)

Transform modes by some BMS transformation

ensure_validity(alter=True, assertions=False)[source]

Try to ensure that the WaveformGrid object is valid

See WaveformBase.ensure_validity for the basic tests. This function also includes tests that data is complex, and consistent with the n_theta and n_phi values.

classmethod from_modes(w_modes, **kwargs)[source]

Construct grid object from modes, with optional BMS transformation

This “constructor” is designed with the goal of transforming the frame in which the modes are measured. If this is not desired, it can be called without those parameters.

It is important to note that the input transformation parameters are applied in the order listed in the parameter list below:

  1. (Super)Translations

  2. Rotation (about the origin of the translated system)

  3. Boost

All input parameters refer to the transformation required to take the mode’s inertial frame onto the inertial frame of the grid’s inertial observers. In what follows, the inertial frame of the modes will be unprimed, while the inertial frame of the grid will be primed. NOTE: These are passive transformations, e.g. supplying the option space_translation=[0, 0, 5] to a Schwarzschild spacetime will move the coordinates to z’=z+5 and so the center of mass will be shifted in the negative z direction by 5 in the new coordinates.

The translations (space, time, spacetime, or super) can be given in various ways, which may override each other. Ultimately, however, they are essentially combined into a single function alpha, representing the supertranslation, which transforms the asymptotic time variable u as

u’(theta, phi) = u - alpha(theta, phi)

A simple time translation would correspond to

alpha(theta, phi) = time_translation

A pure spatial translation would correspond to

alpha(theta, phi) = np.dot(space_translation, -nhat(theta, phi))

where np.dot is the usual dot product, and nhat is the unit vector in the given direction.

Parameters:
w_modesWaveformModes

The object storing the modes of the original waveform, which will be converted to values on a grid in this function. This is the only required argument to this function.

n_thetaint, optional
n_phiint, optional

Number of points in the equi-angular grid in the colatitude (theta) and azimuth (phi) directions. Each defaults to 2*ell_max+1, which is optimal for accuracy and speed. However, this ell_max will typically be greater than the input waveform’s ell_max by at least one, or the ell_max of the input supertranslation (whichever is greater). This is to minimally account for the power at higher orders that such a supertranslation introduces. You may wish to increase this further if the spatial size of your supertranslation is large compared to the smallest wavelength you wish to capture in your data [e.g., ell_max*Omega_orbital_max/speed_of_light], or if your boost speed is close to the speed of light.

time_translationfloat, optional

Defaults to zero. Nonzero overrides spacetime_translation and supertranslation.

space_translationfloat array of length 3, optional

Defaults to empty (no translation). Non-empty overrides spacetime_translation and supertranslation.

spacetime_translationfloat array of length 4, optional

Defaults to empty (no translation). Non-empty overrides supertranslation.

supertranslationcomplex array, optional

This gives the complex components of the spherical-harmonic expansion of the supertranslation in standard form, starting from ell=0 up to some ell_max, which may be different from the ell_max of the input WaveformModes object. Supertranslations must be real, so these values must obey the condition

alpha^{ell,m} = (-1)^m ar{alpha}^{ell,-m}

Defaults to empty, which causes no supertranslation.

frame_rotationquaternion, optional

Transformation applied to (x,y,z) basis of the mode’s inertial frame. For example, the basis z vector of the new grid frame may be written as

z’ = frame_rotation * z * frame_rotation.inverse()

Defaults to 1, corresponding to the identity transformation (no frame_rotation).

boost_velocityfloat array of length 3, optional

This is the three-velocity vector of the grid frame relative to the mode frame. The norm of this vector is checked to ensure that it is smaller than 1. Defaults to [], corresponding to no boost.

psi4_modesWaveformModes, required only if w_modes is type psi3, psi2, psi1, or psi0
psi3_modesWaveformModes, required only if w_modes is type psi2, psi1, or psi0
psi2_modesWaveformModes, required only if w_modes is type psi1 or psi0
psi1_modesWaveformModes, required only if w_modes is type psi0

The objects storing the modes of the original waveforms of the same spacetime that w_modes belongs to. A BMS transformation of an asymptotic Weyl scalar requires mode information from all higher index Weyl scalars. E.g. if w_modes is of type scri.psi2, then psi4_modes and psi3_modes will be required. Note: the WaveformModes objects supplied to these arguments will themselves NOT be transformed. Please use the AsymptoticBondiData class to efficiently transform all asymptotic data at once.

Returns:
WaveformGrid
property n_phi
property n_theta
to_modes(ell_max=None, ell_min=None)[source]

Transform to modes of a spin-weighted spherical harmonic expansion

Parameters:
selfWaveformGrid object

This is the object to be transformed to SWSH modes

ell_maxint, optional

The largest ell value to include in the output data. Default value is deduced from n_theta and n_phi.

ell_minint, optional

The smallest ell value to include in the output data. Default value is abs(spin_weight).

classmethod transform(w_modes, **kwargs)[source]

Transform modes by some BMS transformation

This simply applies the WaveformGrid.from_modes function, followed by the WaveformGrid.to_modes function. See their respective docstrings for more details. However, note that the ell_max parameter used in the second function call defaults here to the ell_max value in the input waveform. This is slightly different from the usual default, because WaveformGrid.from_modes usually increases the effective ell value by 1.

class scri.WaveformModes(*args, **kwargs)[source]

Object containing time, frame, and data, along with related information

This object collects all the data needed to manipulate time-series of spin-weighted spherical-harmonic modes of gravitational waves, as well as a few useful informational members. Various methods exist to manipulate the data and extract important information.

As much as possible, completeness of modes is enforced – meaning that if any (ell,m) mode with a given ell is included in the data, all modes with abs(m)<=ell must also be included. The data are stored as complex numbers, and assumed to be in standard order corresponding to

[f(ell,m) for ell in range(ell_min, ell_max+1) for m in range(-ell,ell+1)]

There is no automated check that can be done to assure that the order is correct. (Precession breaks all symmetries, for example.) The only check that is done ensures that the input data array has the correct dimensions. It is up to the user to ensure that the order is correct.

Attributes:
tfloat array

Time steps corresponding to other data

framequaternion array

Rotors taking static basis onto decomposition basis

data2-d complex array

Mode values of the spin-weighted spherical-harmonic decomposition. The array is assumed to have a first dimension equal to the length of t, and a second dimension equal to the number of modes as deduced from the values of ell_min and ell_max below. As noted above, a very particular order is assumed for the data. Input is treated as a reference, so in-place transformations may modify your original array. To avoid this behavior, pass a copy of your array.

ell_minint

Smallest ell value present in data.

ell_maxint

Largest ell value present in data.

LMint array (read only)

Array of (ell,m) values in Waveform data

historylist of strings

As far as possible, all functions applied to the object are recorded in the history variable. In fact, the object should almost be able to be recreated using the commands in the history list. Commands taking large arrays, however, are shortened – so the data will not be entirely reconstructable.

frameTypeint

Index corresponding to scri.FrameType appropriate for data.

dataTypeint

Index corresponding to scri.DataType appropriate for data.

r_is_scaled_outbool

True if the data have been multiplied by the appropriate power of radius so that the asymptotic value can be finite and nonzero.

m_is_scaled_outbool

True if the data have been scaled by the appropriate value of the total mass so that they are dimensionless.

numint (read only)

Automatically assigned number of this object. The constructor of this type keeps count of the number of objects it has created, to assign each object a more-or-less unique ID for use in the history strings. This counter is reset at the beginning of each python session. Subclasses should automatically have a different counter.

Methods

LLComparisonMatrix(W2)

Calculate the <LL> quantity with respect to the modes of two Waveforms

LLDominantEigenvector([RoughDirection, ...])

Calculate the principal axis of the LL matrix

LLMatrix()

Calculate the <LL> quantity with respect to the modes

LVector(W2)

Calculate the <L> quantity with respect to the modes

LdtVector()

Calculate the <Ldt> quantity with respect to the modes

SI_units(current_unit_mass_in_solar_masses)

Assuming current quantities are in geometric units, convert to SI units

align_decomposition_frame_to_modes(t_fid[, ...])

Fix the attitude of the corotating frame.

angular_momentum_flux([hdot])

Compute angular momentum flux from waveform

angular_velocity([include_frame_velocity])

Angular velocity of Waveform

apply_eth(operations[, eth_convention])

Apply spin raising/lowering operators to waveform mode data in a specified order.

boost_flux([hdot])

Computes the boost flux from the waveform.

compare(w_a[, min_time_step, min_time])

Return a waveform with differences between the two inputs

convert_from_conjugate_pairs()

Convert modes from conjugate-pair format in place

convert_to_conjugate_pairs()

Convert modes to conjugate-pair format in place

copy()

Return a (deep) copy of the object

copy_without_data()

Return a copy of the object, with empty t, frame, and data fields

deepcopy()

Return a deep copy of the object

energy_flux()

Compute energy flux from waveform

ensure_validity([alter, assertions])

Try to ensure that the WaveformModes object is valid

from_grid(w_grid, ell_max)

Transform to modes of a spin-weighted spherical harmonic expansion

from_sxs(w_sxs[, ...])

Construct this object from an sxs.WaveformModes object

index(ell, m)

Index of given (ell,m) mode in the data

indices(ell_m)

Indices of given (ell,m) modes in the data

inner_product(b[, t1, t2, allow_LM_differ, ...])

Compute the all-angles inner product <self, b>.

interpolate(tprime)

Interpolate the frame and data onto the new set of time steps

ladder_factor(operations, s, ell[, ...])

Compute the 'ladder factor' for applying a sequence of spin raising/lower operations on a SWSH of given s, ell.

max_norm_index([skip_fraction_of_data])

Index of time step with largest norm

max_norm_time([skip_fraction_of_data])

Return time at which largest norm occurs in data

momentum_flux()

Compute momentum flux from waveform

norm([take_sqrt, indices])

L2 norm of the waveform

poincare_fluxes([hdot])

Compute fluxes of energy, momentum, angular momemntum, and boost

rotate_decomposition_basis(R_basis)

Rotate a Waveform in place

rotate_physical_system(R_phys)

Rotate a Waveform in place

to_coprecessing_frame([RoughDirection, ...])

Transform waveform (in place) to a coprecessing frame

to_corotating_frame([R0, tolerance, ...])

Transform waveform (in place) to a corotating frame

to_grid(**kwargs)

Construct grid object from modes, with optional BMS transformation

transform(**kwargs)

Transform modes by some BMS transformation

truncate([tol])

Truncate the precision of this object's data in place

to_inertial_frame

LLComparisonMatrix(W2)

Calculate the <LL> quantity with respect to the modes of two Waveforms

The matrix is given in the (possibly rotating) mode frame (X,Y,Z), rather than the inertial frame (x,y,z).

<LL>^{ab} = sum_{ell,m,m’} bar{f}^{ell,m’} < ell,m’ | L_a L_b | ell,m > g^{ell,m}

LLDominantEigenvector(RoughDirection=array([0., 0., 1.]), RoughDirectionIndex=0)

Calculate the principal axis of the LL matrix

param Lmodes L modes to evaluate (optional) param RoughDirection Vague guess about the preferred initial (optional)

If Lmodes is empty (default), all L modes are used. Setting Lmodes to [2] or [2,3,4], for example, restricts the range of the sum.

Ell is the direction of the angular velocity for a PN system, so some rough guess about that direction allows us to choose the direction of the eigenvectors output by this function to be more parallel than anti-parallel to that direction. The default is to simply choose the z axis, since this is most often the correct choice anyway.

The vector is given in the (possibly rotating) mode frame (X,Y,Z), rather than the inertial frame (x,y,z).

LLMatrix()

Calculate the <LL> quantity with respect to the modes

The matrix is given in the (possibly rotating) mode frame (X,Y,Z), rather than the inertial frame (x,y,z).

This quantity is as prescribed by O’Shaughnessy et al. <http://arxiv.org/abs/1109.5224>, except that no normalization is imposed, and this operates on whatever type of data is input.

<LL>^{ab} = Re { sum_{ell,m,m’} bar{f}^{ell,m’} < ell,m’ | L_a L_b | ell,m > f^{ell,m} }

property LM

Array of (ell,m) values in Waveform data

This array is just a flat array of [ell,m] pairs. It is automatically recomputed each time ell_min or ell_max is changed. Specifically, it is

np.array([[ell,m]

for ell in range(self.ell_min,self.ell_max+1) for m in range(-ell,ell+1)])

LVector(W2)

Calculate the <L> quantity with respect to the modes

The vector is given in the (possibly rotating) mode frame (X,Y,Z), rather than the inertial frame (x,y,z).

<L>^{a} = sum_{ell,m,m’} bar{f}^{ell,m’} < ell,m’ | L_a | ell,m > g^{ell,m}

LdtVector()

Calculate the <Ldt> quantity with respect to the modes

The vector is given in the (possibly rotating) mode frame (X,Y,Z), rather than the inertial frame (x,y,z).

<Ldt>^{a} = sum_{ell,m,m’} bar{f}^{ell,m’} < ell,m’ | L_a | ell,m > (df/dt)^{ell,m}

align_decomposition_frame_to_modes(t_fid, nHat_t_fid=quaternion(0, 1, 0, 0), ell_max=None)

Fix the attitude of the corotating frame.

The corotating frame is only defined up to some constant rotor R_eps; if R_corot is corotating, then so is R_corot*R_eps. This function uses that freedom to ensure that the frame is aligned with the Waveform modes at the fiducial time. In particular, it ensures that the Z axis of the frame in which the decomposition is done is along the dominant eigenvector of the <LL> matrix (suggested by O’Shaughnessy et al.), and the phase of the (2,2) mode is zero.

If ell_max is None (default), all ell modes are used.

Parameters:
w: WaveformModes

Object to be aligned

t_fid: float

Fiducial time at which the alignment should happen

nHat_t_fid: quaternion (optional)

The approximate direction of nHat at t_fid; defaults to x

ell_max: int

Highest ell mode to use in computing the <LL> matrix

angular_momentum_flux(hdot=None)

Compute angular momentum flux from waveform

This implements Eq. (2.24) from Ruiz et al. (2008) [0707.4654] by using matrix_expectation_value with (internal) helper functions j_z, j_plus, and j_minus.

angular_velocity(include_frame_velocity=False)

Angular velocity of Waveform

This function calculates the angular velocity of a Waveform object from its modes. This was introduced in Sec. II of “Angular velocity of gravitational radiation and the corotating frame” <http://arxiv.org/abs/1302.2919>. Essentially, this is the angular velocity of the rotating frame in which the time dependence of the modes is minimized.

The vector is given in the (possibly rotating) mode frame (X,Y,Z), which is not necessarily equal to the inertial frame (x,y,z).

apply_eth(operations, eth_convention='NP')[source]

Apply spin raising/lowering operators to waveform mode data in a specified order. This does not modify the original waveform object.

Parameters:
operations: list or str composed of +1, -1, ‘+’, ‘-’, ‘ð’, or ‘ð̅’

The order of eth (+1, ‘+’, or ‘ð’) and ethbar (-1, ‘-’, or ‘ð̅’) operations to perform, applied from right to left. Example, operations=’–+’ will perform on WaveformModes data f the operation ethbar(ethbar(eth(f))). Note that the order of operations is right-to-left.

eth_convention: either ‘NP’ or ‘GHP’ [default: ‘NP’]

Choose between Newman-Penrose or Geroch-Held-Penrose convention

Returns:
mode_data: array of complex

Note that the returned data has the same shape as this object’s data attribute, and the modes correspond to the same (ell, m) values. In particular, if the spin weight changes, the output data will no longer satisfy the expectation that ell_min == abs(s).

boost_flux(hdot=None)

Computes the boost flux from the waveform.

This implements Eq. (C.1) from Flanagan & Nichols (2016) [1510.03386] for the boost vector field.

Notes

Boost flux is analytically calculated using the Wald & Zoupas formalism, following the calculation by Flanagan & Nichols [1510.03386]. Then, the NP formalism is implemented to make it computationally feasible. The conventions of Ruiz et al. (2008) [0707.4654] are followed.

The expression for boost flux is

boost_flux(h) = (-1/32π) (

(1/8) [⟨ð̄N|χ|ð̄h⟩ - ⟨ðN|χ|ðh⟩ + ⟨ð̄h|χ|ð̄N⟩ - ⟨ðh|χ|ðN⟩ + 6⟨N|χ|h⟩ + 6⟨h|χ|N⟩] - (1/4) [⟨ðN|ðχ|h⟩ + ⟨h|ð̄χ|ðN⟩] - (u/2) ⟨N|χ|N⟩

)

h and N have spin weight s=-2. χ is an l=1, s=0 function characterizing the direction of the boost. ð and ð̄ are spin raising and lowering operators.

The matrix elements for most of the terms are identical to pᶻ, p⁺, and p⁻ matrix elements, except that input for s is different. For example, ⟨ð̄N|χ|ð̄h⟩ will have s = -3. There are two terms ⟨ðN|ðχ|h⟩ and ⟨h|ð̄χ|ðN⟩ for which matrix elements have been defined below.

convert_from_conjugate_pairs()[source]

Convert modes from conjugate-pair format in place

This function reverses the effects of convert_to_conjugate_pairs. See that function’s docstring for details.

convert_to_conjugate_pairs()[source]

Convert modes to conjugate-pair format in place

This function alters this object’s modes to store the sum and difference of pairs with opposite m values. If we denote the modes f[l, m], then we define

s[l, m] = (f[l, m] + f̄[l, -m]) / √2 d[l, m] = (f[l, m] - f̄[l, -m]) / √2

For m<0 we replace the mode data with d[l, -m], for m=0 we do nothing, and for m>0 we replace the mode data with s[l, m]. That is, the mode data on output look like this:

[…, d[2, 2], d[2, 1], f[2, 0], s[2, 1], s[2, 2], d[3, 3], d[3, 2], …]

The factor of √2 is chosen so that the sum of the magnitudes squared at each time for this data is the same as it is for the original data.

copy_without_data()[source]

Return a copy of the object, with empty t, frame, and data fields

Note that subclasses may override this to set some related data members. For example, WaveformModes.copy_without_data sets the ell_min and ell_max fields appropriately. If you wish to only skip t, frame, and data, you can simply use WaveformBase.copy_without_data(W). The latter is useful if, for example, you will only be making changes to those three fields, and want everything else preserved.

Also note that some slicing operations can achieve similar – but different – goals. For example, w = w[:, :0] will simply empty data and ells, without affecting the time and frame.

property ell_max
property ell_min
property ells

Return self.ell_min,self.ell_max

energy_flux()

Compute energy flux from waveform

This implements Eq. (2.8) from Ruiz et al. (2008) [0707.4654].

ensure_validity(alter=True, assertions=False)[source]

Try to ensure that the WaveformModes object is valid

See WaveformBase.ensure_validity for the basic tests. This function also includes tests that data is complex, and consistent with the ell_min and ell_max values.

property eth

Returns the spin-raised waveform mode data.

property ethbar

Returns the spin-lowered waveform mode data.

classmethod from_grid(w_grid, ell_max)

Transform to modes of a spin-weighted spherical harmonic expansion

Parameters:
selfWaveformGrid object

This is the object to be transformed to SWSH modes

ell_maxint, optional

The largest ell value to include in the output data. Default value is deduced from n_theta and n_phi.

ell_minint, optional

The smallest ell value to include in the output data. Default value is abs(spin_weight).

classmethod from_sxs(w_sxs, override_exception_from_invalidity=False)[source]

Construct this object from an sxs.WaveformModes object

Note that the resulting object will likely contain references to the same underlying data contained in the original object; modifying one will modify the other. You can make a copy of the result — using code like WaveformModes.from_sxs(w_sxs).copy() — to obtain separate data.

index(ell, m)[source]

Index of given (ell,m) mode in the data

Parameters:
ell: int
m: int
Returns:
idx: int

Index such that self.LM[idx] is [ell, m]

indices(ell_m)[source]

Indices of given (ell,m) modes in the data

Parameters:
ell_m: Nx2 array-like, int
Returns:
idcs: array of int

Index such that self.LM[idcs] is ell_m

inner_product(b, t1=None, t2=None, allow_LM_differ=False, allow_times_differ=False)[source]

Compute the all-angles inner product <self, b>.

self and b should have the same spin-weight, set of (ell,m) modes, and times. Differing sets of modes and times can optionally be turned on.

Parameters:
bWaveformModes object

The other set of modes to inner product with.

t1float, optional [default: None]
t2float, optional [default: None]

Lower and higher bounds of time integration

allow_LM_differbool, optional [default: False]

If True and if the set of (ell,m) modes between self and b differ, then the inner product will be computed using the intersection of the set of modes.

allow_times_differ: bool, optional [default: False]

If True and if the set of times between self and b differ, then both WaveformModes will be interpolated to the intersection of the set of times.

Returns:
inner_productcomplex

<self, b>

ladder_factor(operations, s, ell, eth_convention='NP')[source]

Compute the ‘ladder factor’ for applying a sequence of spin raising/lower operations on a SWSH of given s, ell.

Parameters:
operations: list or str composed of +1, -1, ‘+’, ‘-’, ‘ð’, or ‘ð̅’

The order of eth (+1, ‘+’, or ‘ð’) and ethbar (-1, ‘-’, or ‘ð̅’) operations to perform, applied from right to left. Example, operations=’–+’ will perform on WaveformModes data f the operation ethbar(ethbar(eth(f))). Note that the order of operations is right-to-left.

eth_convention: either ‘NP’ or ‘GHP’ [default: ‘NP’]

Choose between Newman-Penrose or Geroch-Held-Penrose convention

Returns:
float
momentum_flux()

Compute momentum flux from waveform

This implements Eq. (2.11) from Ruiz et al. (2008) [0707.4654] by using matrix_expectation_value with p_z, p_plus, and p_minus.

property n_modes
property parity_antisymmetric_part

Component of waveform that changes sign under parity_conjugate

property parity_conjugate

Reflect modes along all axes

See “Gravitational-wave modes from precessing black-hole binaries” by Boyle et al. (2014) for more details.

property parity_symmetric_part

Component of waveform invariant under parity_conjugate

property parity_violation_normalized

Norm of parity-antisymmetric component divided by norm of waveform

property parity_violation_squared

(Squared) norm of parity-antisymmetric component of waveform

poincare_fluxes(hdot=None)

Compute fluxes of energy, momentum, angular momemntum, and boost

This function will compute the time derivative 1 or 0 times (if an optional argument is passed) so is more efficient than separate calls to energy_flux, momentum_flux, angular_momentum_flux, and boost_flux.

Parameters:
hWaveformModes object

Must have type h

hdotWaveformModes object, optional [default: None]

The time derivative of h. If None, computed from h. Must have type hdot.

Returns:
edot, pdot, jdot, boost_fluxndarray, ndarray, ndarray, ndarray
rotate_decomposition_basis(R_basis)

Rotate a Waveform in place

This function takes a Waveform object W and either a quaternion or array of quaternions R_basis. It applies that rotation to the decomposition basis of the modes in the Waveform. The change in basis is also recorded in the Waveform’s frame data.

For more information on the analytical details, see http://moble.github.io/spherical_functions/SWSHs.html#rotating-swshs

rotate_physical_system(R_phys)

Rotate a Waveform in place

This just rotates the decomposition basis by the inverse of the input rotor(s). See rotate_decomposition_basis.

For more information on the analytical details, see http://moble.github.io/spherical_functions/SWSHs.html#rotating-swshs

to_coprecessing_frame(RoughDirection=array([0., 0., 1.]), RoughDirectionIndex=None, transition_times=None)

Transform waveform (in place) to a coprecessing frame

Parameters:
W: waveform

Waveform object to be transformed in place.

RoughDirection: 3-array [defaults to np.array([0.0, 0.0, 1.0])]

Vague guess about the preferred initial axis, to choose the sign of the eigenvectors.

RoughDirectionIndex: int or None [defaults to None]

Time index at which to apply RoughDirection guess.

transition_times2-tuple or None [defaults to None]

If a 2-tuple of floats is given, these will be interpreted as the beginning and ending times (respectively) to transition from using the coprecessing frame to stopping rotation altogether. This can be helpful for ensuring that the frame doesn’t fluctuate wildly during ringdown.

to_corotating_frame(R0=quaternion(1, 0, 0, 0), tolerance=1e-12, z_alignment_region=None, return_omega=False, truncate_log_frame=False)

Transform waveform (in place) to a corotating frame

Parameters:
W: waveform

Waveform object to be transformed in place

R0: quaternion [defaults to 1]

Initial value of frame when integrating angular velocity

tolerance: float [defaults to 1e-12]

Absolute tolerance used in integration of angular velocity

z_alignment_region: None or 2-tuple of floats [defaults to None]

If not None, the dominant eigenvector of the <LL> matrix is aligned with the z axis, averaging over this portion of the data. The first and second elements of the input are considered fractions of the inspiral at which to begin and end the average. For example, (0.1, 0.9) would lead to starting 10% of the time from the first time step to the max norm time, and ending at 90% of that time.

return_omega: bool [defaults to False]

If True, return a 2-tuple consisting of the waveform in the corotating frame (the usual returned object) and the angular-velocity data. That is frequently also needed, so this is just a more efficient way of getting the data.

truncate_log_frame: bool [defaults to False]

If True, set bits of log(frame) with lower significance than tolerance to zero, and use exp(truncated(log(frame))) to rotate the waveform. Also returns log_frame along with the waveform (and optionally omega)

to_grid(**kwargs)

Construct grid object from modes, with optional BMS transformation

This “constructor” is designed with the goal of transforming the frame in which the modes are measured. If this is not desired, it can be called without those parameters.

It is important to note that the input transformation parameters are applied in the order listed in the parameter list below:

  1. (Super)Translations

  2. Rotation (about the origin of the translated system)

  3. Boost

All input parameters refer to the transformation required to take the mode’s inertial frame onto the inertial frame of the grid’s inertial observers. In what follows, the inertial frame of the modes will be unprimed, while the inertial frame of the grid will be primed. NOTE: These are passive transformations, e.g. supplying the option space_translation=[0, 0, 5] to a Schwarzschild spacetime will move the coordinates to z’=z+5 and so the center of mass will be shifted in the negative z direction by 5 in the new coordinates.

The translations (space, time, spacetime, or super) can be given in various ways, which may override each other. Ultimately, however, they are essentially combined into a single function alpha, representing the supertranslation, which transforms the asymptotic time variable u as

u’(theta, phi) = u - alpha(theta, phi)

A simple time translation would correspond to

alpha(theta, phi) = time_translation

A pure spatial translation would correspond to

alpha(theta, phi) = np.dot(space_translation, -nhat(theta, phi))

where np.dot is the usual dot product, and nhat is the unit vector in the given direction.

Parameters:
w_modesWaveformModes

The object storing the modes of the original waveform, which will be converted to values on a grid in this function. This is the only required argument to this function.

n_thetaint, optional
n_phiint, optional

Number of points in the equi-angular grid in the colatitude (theta) and azimuth (phi) directions. Each defaults to 2*ell_max+1, which is optimal for accuracy and speed. However, this ell_max will typically be greater than the input waveform’s ell_max by at least one, or the ell_max of the input supertranslation (whichever is greater). This is to minimally account for the power at higher orders that such a supertranslation introduces. You may wish to increase this further if the spatial size of your supertranslation is large compared to the smallest wavelength you wish to capture in your data [e.g., ell_max*Omega_orbital_max/speed_of_light], or if your boost speed is close to the speed of light.

time_translationfloat, optional

Defaults to zero. Nonzero overrides spacetime_translation and supertranslation.

space_translationfloat array of length 3, optional

Defaults to empty (no translation). Non-empty overrides spacetime_translation and supertranslation.

spacetime_translationfloat array of length 4, optional

Defaults to empty (no translation). Non-empty overrides supertranslation.

supertranslationcomplex array, optional

This gives the complex components of the spherical-harmonic expansion of the supertranslation in standard form, starting from ell=0 up to some ell_max, which may be different from the ell_max of the input WaveformModes object. Supertranslations must be real, so these values must obey the condition

alpha^{ell,m} = (-1)^m ar{alpha}^{ell,-m}

Defaults to empty, which causes no supertranslation.

frame_rotationquaternion, optional

Transformation applied to (x,y,z) basis of the mode’s inertial frame. For example, the basis z vector of the new grid frame may be written as

z’ = frame_rotation * z * frame_rotation.inverse()

Defaults to 1, corresponding to the identity transformation (no frame_rotation).

boost_velocityfloat array of length 3, optional

This is the three-velocity vector of the grid frame relative to the mode frame. The norm of this vector is checked to ensure that it is smaller than 1. Defaults to [], corresponding to no boost.

psi4_modesWaveformModes, required only if w_modes is type psi3, psi2, psi1, or psi0
psi3_modesWaveformModes, required only if w_modes is type psi2, psi1, or psi0
psi2_modesWaveformModes, required only if w_modes is type psi1 or psi0
psi1_modesWaveformModes, required only if w_modes is type psi0

The objects storing the modes of the original waveforms of the same spacetime that w_modes belongs to. A BMS transformation of an asymptotic Weyl scalar requires mode information from all higher index Weyl scalars. E.g. if w_modes is of type scri.psi2, then psi4_modes and psi3_modes will be required. Note: the WaveformModes objects supplied to these arguments will themselves NOT be transformed. Please use the AsymptoticBondiData class to efficiently transform all asymptotic data at once.

Returns:
WaveformGrid
to_inertial_frame()
property to_sxs

Convert this object to an sxs.WaveformModes object

Note that the resulting object will likely contain references to the same underlying data contained in the original object; modifying one will modify the other. You can make a copy of this object before calling this function — using code like w.copy().to_sxs — to obtain separate data.

transform(**kwargs)[source]

Transform modes by some BMS transformation

This simply applies the WaveformGrid.from_modes function, followed by the WaveformGrid.to_modes function. See their respective docstrings for more details.

However, note that the ell_max parameter used in the second function call defaults here to the ell_max value in the input waveform. This is slightly different from the usual default, because WaveformGrid.from_modes usually increases the effective ell value by 1.

truncate(tol=1e-10)[source]

Truncate the precision of this object’s data in place

This function sets bits in self.data to 0 when they will typically contribute less than tol times the norm of the Waveform at that instant in time. Here, “typically” means that we divide tol by the square-root of the number of modes in the waveform, so that if each mode contributes a random amount of error at that level, we can expect a total error of roughly tol times the norm.

property x_parity_antisymmetric_part

Return component of waveform that changes sign under x_parity_conjugate

property x_parity_conjugate

Reflect modes across y-z plane (along x axis)

See “Gravitational-wave modes from precessing black-hole binaries” by Boyle et al. (2014) for more details.

property x_parity_symmetric_part

Return component of waveform invariant under x_parity_conjugate

property x_parity_violation_normalized

Norm of x-parity-antisymmetric component divided by norm of waveform

property x_parity_violation_squared

(Squared) norm of x-parity-antisymmetric component of waveform

property y_parity_antisymmetric_part

Component of waveform that changes sign under y_parity_conjugate

property y_parity_conjugate

Reflect modes across x-z plane (along y axis)

See “Gravitational-wave modes from precessing black-hole binaries” by Boyle et al. (2014) for more details.

property y_parity_symmetric_part

Component of waveform invariant under y_parity_conjugate

property y_parity_violation_normalized

Norm of y-parity-antisymmetric component divided by norm of waveform

property y_parity_violation_squared

(Squared) norm of y-parity-antisymmetric component of waveform

property z_parity_antisymmetric_part

Component of waveform that changes sign under z_parity_conjugate

property z_parity_conjugate

Reflect modes across x-y plane (along z axis)

See “Gravitational-wave modes from precessing black-hole binaries” by Boyle et al. (2014) for more details.

property z_parity_symmetric_part

Component of waveform invariant under z_parity_conjugate

property z_parity_violation_normalized

Norm of z-parity-antisymmetric component divided by norm of waveform

property z_parity_violation_squared

(Squared) norm of z-parity-antisymmetric component of waveform

Modules

scri.LVC

Submodule for operating on LVC-format waveform files

scri.SpEC

Submodule for operating on SpEC waveform files

scri.asymptotic_bondi_data

scri.bms_transformations

scri.extrapolation

scri.flux

scri.mode_calculations

scri.modes_time_series

scri.plotting

scri.pn

scri.rotations

scri.sample_waveforms

scri.utilities

scri.waveform_base

scri.waveform_grid

scri.waveform_modes