scri.waveform_grid

Functions

process_transformation_kwargs(ell_max, **kwargs)

Classes

WaveformGrid(*args, **kwargs)

Attributes:

class scri.waveform_grid.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.

scri.waveform_grid.process_transformation_kwargs(ell_max, **kwargs)[source]