scri.asymptotic_bondi_data

Classes

AsymptoticBondiData(time, ell_max[, ...])

Class to store asymptotic Bondi data

class scri.asymptotic_bondi_data.AsymptoticBondiData(time, ell_max, multiplication_truncator=<built-in function sum>, frameType=1)[source]

Class to store asymptotic Bondi data

This class stores time data, along with the corresponding values of psi0 through psi4 and sigma. For simplicity, the data are stored as one contiguous array. That is, all values are stored at all times, even if they are zero, and all Modes objects are stored with ell_min=0, even when their spins are not zero.

The single contiguous array is then viewed as 6 separate ModesTimeSeries objects, which enables them to track their spin weights, and provides various convenient methods like eth and ethbar; dot and ddot for time-derivatives; int and iint for time-integrations; norm to take the norm of a function over the sphere; bar for conjugation of the functions (which is different from just conjugating the mode weights); etc. It also handles algebra correctly – particularly addition (which is disallowed when the spin weights differ) and multiplication (which can be delicate with regards to the resulting ell values).

This may lead to some headaches when the user tries to do things that are disabled by Modes objects. The goal is to create headaches if and only if the user is trying to do things that really should never be done (like conjugating mode weights, rather than the underlying function; adding modes with different spin weights; etc.). Please open issues for any situations that don’t meet this standard.

This class also provides various convenience methods for computing things like the mass aspect, the Bondi four-momentum, the Bianchi identities, etc.

Attributes:
LM
bondi_violation_norms

Compute norms of violations of Bondi-gauge conditions

bondi_violations

Compute violations of Bondi-gauge constraints

ell_max
ell_min
n_modes
n_times
psi0
psi1
psi2
psi3
psi4
sigma
t
time
u

Methods

CWWY_angular_momentum()

Compute the Chen/Wang/Wang/Yau angular momentum vector.

bianchi_0([lhs, rhs])

Return the left- and/or right-hand sides of the Psi0 component of the Bianchi identity

bianchi_1([lhs, rhs])

Return the left- and/or right-hand sides of the Psi1 component of the Bianchi identity

bianchi_2([lhs, rhs])

Return the left- and/or right-hand sides of the Psi2 component of the Bianchi identity

bondi_CoM_charge()

Compute the center-of-mass charge vector

bondi_angular_momentum()

Compute the total Bondi angular momentum vector

bondi_boost_charge()

Compute the Bondi boost charge vector

bondi_constraints([lhs, rhs])

Compute Bondi-gauge constraint equations

bondi_dimensionless_spin()

Compute the dimensionless Bondi spin vector

bondi_four_momentum()

Compute the Bondi four-momentum

bondi_rest_mass()

Compute the rest mass from the Bondi four-momentum

constraint_3([lhs, rhs])

Return the left- and/or right-hand sides of the Psi3 expression in Bondi gauge

constraint_4([lhs, rhs])

Return the left- and/or right-hand sides of the Psi4 expression in Bondi gauge

constraint_mass_aspect([lhs, rhs])

Return the left- and/or right-hand sides of the mass-aspect reality condition in Bondi gauge

from_initial_values(time[, ell_max, sigma0, ...])

Construct Bondi data from sigma as a function of time and optional initial values

map_to_abd_frame(target_abd[, t_0, ...])

Transform an abd object to a target abd object using data at t=0.

map_to_superrest_frame([t_0, ...])

Transform an abd object to the superrest frame.

mass_aspect([truncate_ell])

Compute the Bondi mass aspect of the AsymptoticBondiData.

supermomentum(supermomentum_def, **kwargs)

Compute the supermomentum

transform(**kwargs)

Apply BMS transformation to AsymptoticBondiData object

copy

interpolate

CWWY_angular_momentum()

Compute the Chen/Wang/Wang/Yau angular momentum vector.

See Eq. (5) of https://arxiv.org/abs/2102.03235

property LM
bianchi_0(lhs=True, rhs=True)

Return the left- and/or right-hand sides of the Psi0 component of the Bianchi identity

In Bondi coordinates, the Bianchi identities simplify, resulting in this expression (among others) for the time derivative of ψ₀:

ψ̇₀ = ðψ₁ + 3 σ ψ₂

Parameters:
lhs: bool [defaults to True]

If True, return the left-hand side of the equation above

rhs: bool [defaults to True]

If True, return the right-hand side of the equation above

If both parameters are True, a tuple with elements (lhs_value, rhs_value) is
returned; otherwise just the requested value is returned.
bianchi_1(lhs=True, rhs=True)

Return the left- and/or right-hand sides of the Psi1 component of the Bianchi identity

In Bondi coordinates, the Bianchi identities simplify, resulting in this expression (among others) for the time derivative of ψ₁:

ψ̇₁ = ðψ₂ + 2 σ ψ₃

Parameters:
lhs: bool [defaults to True]

If True, return the left-hand side of the equation above

rhs: bool [defaults to True]

If True, return the right-hand side of the equation above

If both parameters are True, a tuple with elements (lhs_value, rhs_value) is
returned; otherwise just the requested value is returned.
bianchi_2(lhs=True, rhs=True)

Return the left- and/or right-hand sides of the Psi2 component of the Bianchi identity

In Bondi coordinates, the Bianchi identities simplify, resulting in this expression (among others) for the time derivative of ψ₂:

ψ̇₂ = ðψ₃ + σ ψ₄

Parameters:
lhs: bool [defaults to True]

If True, return the left-hand side of the equation above

rhs: bool [defaults to True]

If True, return the right-hand side of the equation above

If both parameters are True, a tuple with elements (lhs_value, rhs_value) is
returned; otherwise just the requested value is returned.
bondi_CoM_charge()

Compute the center-of-mass charge vector

Gⁱ = Nⁱ + t Pⁱ = - [ψ₁ + σ ðσ̄ + ½ð(σ σ̄)]

where Nⁱ is the boost charge and Pⁱ is the momentum. See Eq. (3.4) of arXiv:1912.03164.

bondi_angular_momentum()

Compute the total Bondi angular momentum vector

i (ψ₁ + σ ðσ̄)

See Eq. (8) of Dray (1985) iopscience.iop.org/article/10.1088/0264-9381/2/1/002

bondi_boost_charge()

Compute the Bondi boost charge vector

  • [ψ₁ + σ ðσ̄ + ½ð(σ σ̄) - t ð ℜ{ψ₂ + σ ∂ₜσ̄}]

See Eq. (8) of Dray (1985) iopscience.iop.org/article/10.1088/0264-9381/2/1/002

bondi_constraints(lhs=True, rhs=True)

Compute Bondi-gauge constraint equations

Bondi gauge establishes some relations that the data must satisfy:

ψ̇₀ = ðψ₁ + 3 σ ψ₂

ψ̇₁ = ðψ₂ + 2 σ ψ₃

ψ̇₂ = ðψ₃ + 1 σ ψ₄

ψ₃ = -∂ðσ̄/∂u

ψ₄ = -∂²σ̄/∂u²

Im[ψ₂] = -Im[ð²σ̄ + σ ∂σ̄/∂u]

This function returns a 6-tuple of 2-tuples, corresponding to these 6 equations and their left- and right-hand sides.

bondi_dimensionless_spin()

Compute the dimensionless Bondi spin vector

bondi_four_momentum()

Compute the Bondi four-momentum

This is just the ell<2 component of the mass aspect, expressed as a four-vector.

bondi_rest_mass()

Compute the rest mass from the Bondi four-momentum

property bondi_violation_norms

Compute norms of violations of Bondi-gauge conditions

Bondi gauge establishes some relations that the data must satisfy:

ψ̇₀ = ðψ₁ + 3 σ ψ₂

ψ̇₁ = ðψ₂ + 2 σ ψ₃

ψ̇₂ = ðψ₃ + 1 σ ψ₄

ψ₃ = -∂ðσ̄/∂u

ψ₄ = -∂²σ̄/∂u²

Im[ψ₂] = -Im[ð²σ̄ + σ ∂σ̄/∂u]

This function returns a tuple of 6 arrays, corresponding to the norms of these 6 equations, in which the right-hand side is subtracted from the left-hand side, and then the squared magnitude of that result is integrated over the sphere. No integration is performed over time.

property bondi_violations

Compute violations of Bondi-gauge constraints

Bondi gauge establishes some relations that the data must satisfy:

ψ̇₀ = ðψ₁ + 3 σ ψ₂

ψ̇₁ = ðψ₂ + 2 σ ψ₃

ψ̇₂ = ðψ₃ + 1 σ ψ₄

ψ₃ = -∂ðσ̄/∂u

ψ₄ = -∂²σ̄/∂u²

Im[ψ₂] = -Im[ð²σ̄ + σ ∂σ̄/∂u]

This function returns a tuple of 6 arrays, corresponding to these 6 equations, in which the right-hand side is subtracted from the left-hand side. No norms are taken.

constraint_3(lhs=True, rhs=True)

Return the left- and/or right-hand sides of the Psi3 expression in Bondi gauge

In Bondi coordinates, the value of ψ₃ is given by a time derivative and an angular derivative of the (conjugate) shear:

ψ₃ = -∂ðσ̄/∂u

Parameters:
lhs: bool [defaults to True]

If True, return the left-hand side of the equation above

rhs: bool [defaults to True]

If True, return the right-hand side of the equation above

If both parameters are True, a tuple with elements (lhs_value, rhs_value) is
returned; otherwise just the requested value is returned.
constraint_4(lhs=True, rhs=True)

Return the left- and/or right-hand sides of the Psi4 expression in Bondi gauge

In Bondi coordinates, the value of ψ₄ is given by the second time derivative of the (conjugate) shear:

ψ₄ = -∂²σ̄/∂u²

Parameters:
lhs: bool [defaults to True]

If True, return the left-hand side of the equation above

rhs: bool [defaults to True]

If True, return the right-hand side of the equation above

If both parameters are True, a tuple with elements (lhs_value, rhs_value) is
returned; otherwise just the requested value is returned.
constraint_mass_aspect(lhs=True, rhs=True)

Return the left- and/or right-hand sides of the mass-aspect reality condition in Bondi gauge

In Bondi coordinates, the Bondi mass aspect is always real, resulting in this relationship:

Im[ψ₂] = -Im[ð²σ̄ + σ ∂σ̄/∂u]

Parameters:
lhs: bool [defaults to True]

If True, return the left-hand side of the equation above

rhs: bool [defaults to True]

If True, return the right-hand side of the equation above

If both parameters are True, a tuple with elements (lhs_value, rhs_value) is
returned; otherwise just the requested value is returned.
copy()[source]
property ell_max
property ell_min
classmethod from_initial_values(time, ell_max=8, sigma0=0.0, sigmadot0=0.0, sigmaddot0=0.0, psi2=0.0, psi1=0.0, psi0=0.0)

Construct Bondi data from sigma as a function of time and optional initial values

The initial-value formulation for Bondi gauge is determined by these relations:

ψ̇₀ = ðψ₁ + 3 σ ψ₂ ψ̇₁ = ðψ₂ + 2 σ ψ₃ ψ̇₂ = ðψ₃ + 1 σ ψ₄ ψ₃ = -∂ðσ̄/∂u ψ₄ = -∂²σ̄/∂u²

We also have a constraint on the initial value of ψ₂:

Im[ψ₂] = -Im[ð²σ̄ + σ ∂σ̄/∂u]

Given these expressions, and the value of sigma as a function of time and direction, we can find the values of sigma and all the Weyl components for all time – the Bondi data. This function implements that algorithm.

Note that the interpretation of the input mode weights depends on the dimension of sigma0. If it has dimension 0 or 1, the problem is integrated exactly, and the rest of the inputs are assumed to represent the corresponding values at time u=0. If it has dimension 2, sigma0 is interpreted as a function of time, and the rest of the problem is integrated numerically via splines, with psi2, psi1, and psi0 representing the corresponding values at time u=`time[0]` (the first element of the input time array); sigmadot0 and sigmaddot0 are ignored in this case.

Parameters:
cls: class

Class to construct the AsymptoticBondiData object. This function will also be used as a classmethod in that class, in which case this parameter will be passed automatically by calling this function as a method.

time: array_like

Times at which to compute the Bondi data. Must be 1-dimensional.

ell_max: int

Maximum ell value to be stored in the data

sigma0: 0.0 or array_like [defaults to 0.0]

This represents the value of sigma as a function of time and direction, or simply its initial value (see discussion above). The input quantity must broadcast against an array of shape (n_times, LM_total_size(0, ell_max)).

sigmadot0: 0.0 or array_like [defaults to 0.0]

This represents the time-derivative of sigma as a function of time and direction. Note that this is ignored if sigma0 is a 2-dimensional array; instead, this quantity is computed by differentiating a spline of sigma0. This must be 0- or 1-dimensional, and is assumed to represent the derivative at time 0. (Or, if sigmaddot0 is 0, then the constant derivative at any time.)

sigmaddot0: 0.0 or array_like [defaults to 0.0]

Just like sigmadot0, except for the second derivative. Note that this represents the second derivative only at time u=0; for all other times we expand as a Taylor series if this is used.

psi2: 0.0 or array_like [defaults to 0.0]

This represents the initial value of the psi2 field. Its imaginary part is determined by the condition that the mass aspect must be real-valued, so any imaginary part in the input is discarded. If array_like, this parameter must be 1-dimensional and have size LM_total_size(0, ell_max).

psi1: 0.0 or array_like [defaults to 0.0]

This represents the initial value of the psi1 field. If array_like, this parameter must be 1-dimensional and have size LM_total_size(0, ell_max).

psi0: 0.0 or array_like [defaults to 0.0]

This represents the initial value of the psi0 field. If array_like, this parameter must be 1-dimensional and have size LM_total_size(0, ell_max).

interpolate(new_times)[source]
map_to_abd_frame(target_abd, t_0=0, padding_time=250, N_itr_maxes={'CoM_transformation': 10, 'abd': 2, 'rotation': 10, 'superrest': 2, 'supertranslation': 10}, rel_err_tols={'CoM_transformation': 1e-12, 'rotation': 1e-12, 'supertranslation': 1e-12}, order=['supertranslation', 'rotation', 'CoM_transformation'], ell_max=None, alpha_ell_max=None, fix_time_phase_freedom=True, nprocs=4, print_conv=False)

Transform an abd object to a target abd object using data at t=0.

This computes the transformations necessary to map an abd object to a target abd object. It uses the function map_to_bms_frame with the target charges computed from the target abd object.

Parameters:
target_abdAsymptoticBondiData

Target AsymptoticBondiData to map to.

t_0: float, optional

When to map to the target BMS frame. Default is 0.

padding_timefloat, optional

Amount by which to pad around t=0 to speed up computations, i.e., distance from t=0 in each direction to be included in self.interpolate(…). This also determines the range over which certain BMS charges will be computed. Default is 250.

N_itr_maxesdict, optional

Maximum number of iterations to perform for each transformation. For ‘abd’ , this is the number of iterations to use for the abd procedure. For ‘superrest’, this is the number of iterations to use for the superrest procedure. For the other options, these are the number of iterations to use for finding each individual transformation. Default is N_itr_maxes = {

‘abd’: 2, ‘superrest’: 2, ‘CoM_transformation’: 10, ‘rotation’: 10, ‘supertranslation’: 10,

}.

rel_err_tolsdict, optional

Relative error tolerances for each transformation. Default is rel_err_tols = {

‘CoM_transformation’: 1e-12, ‘rotation’: 1e-12, ‘supertranslation’: 1e-12

}.

orderlist, optional

Order in which to solve for the BMS transformations. Default is [“supertranslation”, “rotation”, “CoM_transformation”].

ell_maxint, optional

Maximum ell to use for SWSH/Grid transformations. Default is self.ell_max.

alpha_ell_maxint, optional

Maximum ell of the supertranslation to use. Default is self.ell_max.

fix_time_phase_freedombool, optional

Whether or not to fix the time and phase freedom using a 2d minimization scheme. Default is True.

nprocsint, optional

Number of cpus to use during parallelization for fixing the time and phase freedom. Default is 4. ‘None’ corresponds to the maximum number. ‘-1’ corresponds to no parallelization.

print_conv: bool, defaults to False

Whether or not to print the termination criterion. Default is False.

Returns:
abd_primeAsymptoticBondiData

Result of self.transform(…) where the input transformations are the transformations found in the BMSTransformations object.

transformationsBMSTransformation

BMS transformation to map to the target BMS frame.

map_to_superrest_frame(t_0=0, target_PsiM_input=None, target_strain_input=None, padding_time=250, N_itr_maxes={'CoM_transformation': 10, 'rotation': 10, 'superrest': 2, 'supertranslation': 10}, rel_err_tols={'CoM_transformation': 1e-12, 'rotation': 1e-12, 'supertranslation': 1e-12}, order=['supertranslation', 'rotation', 'CoM_transformation'], ell_max=None, alpha_ell_max=None, fix_time_phase_freedom=False, modes=None, print_conv=False)

Transform an abd object to the superrest frame.

This computes the transformations necessary to map an abd object to the superrest frame by iteratively minimizing various BMS charges at a certain time. This can be performed either without a target (in which case it maps to the superrest frame at the input time t_0) or with a target (in which case it maps to the frame of the target). For example, one could supply a t_0 which is toward the end of the simulation in which case this code would map to the superrest frame of the remnant BH. Or, one could supply a target PN strain and supermomentum and a t_0 during the early inspiral in which case this code would map to the PN BMS frame.

This wholly fixes the BMS frame of the abd object up to a time translation and a phase rotation.

Parameters:
t_0float, optional

When to map to the superrest frame. Default is 0.

target_PsiM_inputWaveformModes, optional

Target Moreschi supermomentum to map to. Default is 0.

target_strain_inputWaveformModes, optional

Target strain used to constrain the rotation freedom via the angular momentum flux. Default is aligning to the z-axis.

padding_timefloat, optional

Amount by which to pad around t_0 to speed up computations, i.e., distance from t_0 in each direction to be included in self.interpolate(…). This also determines the range over which certain BMS charges will be computed. Default is 250.

N_itr_maxesdict, optional

Maximum number of iterations to perform for each transformation. For ‘superrest’, this is the number of iterations to use for the superrest procedure. For the other options, these are the number of iterations to use for finding each individual transformation. Default is N_itr_maxes = {

‘superrest’: 2, ‘CoM_transformation’: 10, ‘rotation’: 10, ‘supertranslation’: 10,

}.

rel_err_tolsdict, optional

Relative error tolerances for each transformation. Default is rel_err_tols = {

‘CoM_transformation’: 1e-12, ‘rotation’: 1e-12, ‘supertranslation’: 1e-12

}.

orderlist, optional

Order in which to solve for the BMS transformations. Default is [“rotation”, “CoM_transformation”, “supertranslation”].

ell_maxint, optional

Maximum ell to use for SWSH/Grid transformations. Default is self.ell_max.

alpha_ell_maxint, optional

Maximum ell of the supertranslation to use. Default is self.ell_max.

fix_time_phase_freedombool, optional

Whether or not to fix the time and phase freedom using a 2d minimization scheme. Default is True.

modeslist, optional

List of modes to include when performing the 2d alignment. Default is every mode.

print_conv: bool, defaults to False

Whether or not to print the termination criterion. Default is False.

Returns:
abd_primeAsymptoticBondiData

Result of self.transform(…) where the input transformations are the transformations found in the BMSTransformations object.

transformationsBMSTransformation

BMS transformation to map to the target BMS frame.

mass_aspect(truncate_ell=<built-in function max>)

Compute the Bondi mass aspect of the AsymptoticBondiData.

The Bondi mass aspect is given by

M = -ℜ{ψ₂ + σ ∂ₜσ̄}

Note that the last term is a product between two fields. If, for example, these both have ell_max=8, then their full product would have ell_max=16, meaning that we would go from tracking 81 modes to 289. This shows that deciding how to truncate the output ell is important, which is why this function has the extra argument that it does.

Parameters:
truncate_ell: int, or callable [defaults to `max`]

Determines how the ell_max value of the output is determined. If an integer is passed, each term in the output is truncated to have at most that ell_max. (In particular, terms that will not be used in the output are simply not computed, without incurring any errors due to aliasing.) If a callable is passed, it is passed on to the spherical_functions.Modes.multiply method. See that function’s docstring for details. The default behavior will result in the output having ell_max equal to the largest of any of the individual Modes objects in the equation for M above – but not the product.

property n_modes
property n_times
property psi0
property psi1
property psi2
property psi3
property psi4
property sigma
supermomentum(supermomentum_def, **kwargs)

Compute the supermomentum

This function allows for several different definitions of the supermomentum. These differences only apply to ell > 1 modes, so they do not affect the Bondi four-momentum. See Eqs. (7-9) in arXiv:1404.2475 for the different supermomentum definitions and links to further references.

In the literature, there is an ambiguity of vocabulary. When it comes to other BMS charges, we clearly distinuish between the “charge” and the “aspect”. However, the term “supermomentum” is used for both. Accordingly, this function provides two ways to compute the supermomentum.

  1. By default, the supermomentum will be computed as

    Ψ = ψ₂ + σ ∂ₜσ̄ + f(θ, ϕ)

    (See below for the definitions of f.)

  2. By passing the option integrated=True, the supermomentum will instead be computed as

    Pₗₘ = - (1/4π) ∫ Ψ(θ, ϕ) Yₗₘ(θ, ϕ) dΩ

Parameters:
supermomentum_defstr

The definition of the supermomentum to be computed. One of the following (case-insensitive) options can be specified:

  • ‘Bondi-Sachs’ or ‘BS’ for f = 0

  • ‘Moreschi’ or ‘M’ for f = ð²σ̄

  • ‘Geroch’ or ‘G’ for f = ½ (ð²σ̄ - ð̄²σ)

  • ‘Geroch-Winicour’ or ‘GW’ for f = - ð̄²σ

integratedbool, optional

If True, then return the integrated form of the supermomentum — see Eq. (6) in arXiv:1404.2475. Default is False

working_ell_max: int, optional

The value of ell_max to be used to define the computation grid. The number of theta points and the number of phi points are set to 2*working_ell_max+1. Defaults to 2*self.ell_max.

Returns:
ModesTimeSeries
property t
property time
transform(**kwargs)

Apply BMS transformation to AsymptoticBondiData object

It is important to note that the input transformation parameters are applied in this order:

  1. (Super)Translations

  2. Rotation (about the origin)

  3. Boost (about the origin)

All input parameters refer to the transformation required to take the input data’s inertial frame onto the inertial frame of the output data’s inertial observers. In what follows, the coordinates of and functions in the input inertial frame will be unprimed, while corresponding values of the output inertial frame will be primed.

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 α, representing the supertranslation, which transforms the asymptotic time variable u as

u’(u, θ, ϕ) = u(u, θ, ϕ) - α(θ, ϕ)

A simple time translation by δt would correspond to

α(θ, ϕ) = δt # Independent of (θ, ϕ)

A pure spatial translation δx would correspond to

α(θ, ϕ) = -δx · n̂(θ, ϕ)

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

Parameters:
abd: AsymptoticBondiData

The object storing the modes of the original data, which will be transformed in this function. This is the only required argument to this function.

time_translation: float, optional

Defaults to zero. Nonzero overrides corresponding components of spacetime_translation and supertranslation parameters. Note that this is the actual change in the coordinate value, rather than the corresponding mode weight (which is what supertranslation represents).

space_translationfloat array of length 3, optional

Defaults to empty (no translation). Non-empty overrides corresponding components of spacetime_translation and supertranslation parameters. Note that this is the actual change in the coordinate value, rather than the corresponding mode weight (which is what supertranslation represents).

spacetime_translationfloat array of length 4, optional

Defaults to empty (no translation). Non-empty overrides corresponding components of supertranslation. Note that this is the actual change in the coordinate value, rather than the corresponding mode weight (which is what supertranslation represents).

supertranslationcomplex array [defaults to 0]

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 abd object. Supertranslations must be real, so these values should obey the condition

α^{ℓ,m} = (-1)^m ᾱ^{ℓ,-m}

This condition is actually imposed on the input data, so imaginary parts of α(θ, ϕ) will essentially be discarded. Defaults to empty, which causes no supertranslation. Note that some components may be overridden by the parameters above.

frame_rotationquaternion [defaults to 1]

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

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

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

boost_velocityfloat array of length 3 [defaults to (0, 0, 0)]

This is the three-velocity vector of the new frame relative to the input frame. The norm of this vector is required to be smaller than 1.

output_ell_max: int [defaults to abd.ell_max]

Maximum ell value in the output data.

working_ell_max: int [defaults to 2 * abd.ell_max]

Maximum ell value to use during the intermediate calculations. Rotations and time translations do not require this to be any larger than abd.ell_max, but other transformations will require more values of ell for accurate results. In particular, boosts are multiplied by time, meaning that a large boost of data with large values of time will lead to very large power in higher modes. Similarly, large (super)translations will couple power through a lot of modes. To avoid aliasing, this value should be large, to accomodate power in higher modes.

Returns:
abdprime: AsymptoticBondiData

Object representing the transformed data.

property u

Modules

scri.asymptotic_bondi_data.bms_charges

scri.asymptotic_bondi_data.constraints

scri.asymptotic_bondi_data.from_initial_values

scri.asymptotic_bondi_data.map_to_abd_frame

scri.asymptotic_bondi_data.map_to_superrest_frame

scri.asymptotic_bondi_data.transformations