mayawaves.radiation module
- class mayawaves.radiation.RadiationBundle(radiation_spheres: dict)
Class for interacting with all radiative information from the simulation.
- create_extrapolated_sphere(order: int = 1)
Create a RadiationSphere object extrapolated from the RadiationSphere at the provided radius for extrapolation.
Extrapolate the \(\Psi_4\) of all RadiationModes in the RadiationSphere at the set radius for extrapolation to infinite radius and create and store a new RadiationSphere with those extrapolated RadiationModes. Uses the extrapolation method described in https://arxiv.org/abs/1008.4360 and https://arxiv.org/abs/1108.4421.
- Parameters:
order (
int
, optional) – the extrapolation order. Defaults to 1.
- static create_radiation_bundle(radiation_group: Group)
Create a RadiationBundle
- Parameters:
radiation_group (h5py.Group) – group containing all the data pertaining to \(\Psi_4\) for all modes and all extraction radii.
- Returns:
a RadiationBundle consisting of RadiationSpheres for each extraction radius
- Return type:
- property extrapolated_sphere
The RadiationSphere with extrapolated radius. All \(\Psi_4\) data has been extrapolated to infinite radius using the method described in https://arxiv.org/abs/1008.4360 and https://arxiv.org/abs/1108.4421.
- property frame: Frame
The frame for the radiation modes
This can be either the raw frame (Frame.RAW) or the frame which corrects for the drift of the center of mass (FRAME.COM_CORRECTED).
- get_dEnergy_dt_radiated(extraction_radius: float | None = None, **kwargs) tuple
Rate at which energy is radiated, \(dE/dt\)
Uses the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
extraction_radius (
float
, optional) – radius for gravitational wave extraction. If not provided, extrapolates to infinite radius.**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minimum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value
- Returns:
time, \(dE/dt\) for the given modes
- Return type:
tuple
- get_dP_dt_radiated(extraction_radius: float | None = None, **kwargs) tuple
Rate at which linear momentum is radiated
Rate at which linear momentum is radiated through gravitational waves, computed using the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
extraction_radius (
float
, optional) – radius for gravitational wave extraction. If not provided, extrapolates to infinite radius.**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value
- Returns:
time, dP/dt for the given modes
- Return type:
tuple
- get_energy_radiated(extraction_radius: float | None = None, **kwargs) tuple
Cumulative radiated energy \(E\)
Uses the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
extraction_radius (
float
, optional) – radius for gravitational wave extraction. If not provided, extrapolates to infinite radius.**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value
- Returns:
time, cumulative energy radiated for the given modes
- Return type:
tuple
- get_linear_momentum_radiated(extraction_radius: float | None = None, **kwargs) tuple
Linear momentum radiated
Cummulative linear momentum radiated through gravitational waves, computed using the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
extraction_radius (
float
, optional) – radius for gravitational wave extraction. If not provided, extrapolates to infinite radius.**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value
- Returns:
time, cummulative linear momentum radiated for the given modes
- Return type:
tuple
- get_psi4_amplitude_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Amplitude of \(\Psi_4\) for a given mode and extraction radius.
Returns the amplitude of \(\Psi_4\) for a given mode and extraction radius. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for \(\Psi_4\) data. If 0, the data at infinity is provided.
- Returns:
\(\Psi_4\) amplitude for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_psi4_imaginary_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Imaginary component of \(\Psi_4\) for a given mode and extraction radius.
Returns the imaginary part of \(\Psi_4\) for a given mode and extraction radius. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for \(\Psi_4\) data. If 0, the data at infinity is provided.
- Returns:
\(\Psi_4\) imaginary component for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_psi4_max_time_for_mode(l: int, m: int, extraction_radius: float = 0) int
Time of maximum \(\Psi_4\) amplitude for a given mode and extraction radius.
The time at which the amplitude of \(\Psi_4\) reaches its peak. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – radius at which \(\Psi_4\) was extracted
- Returns:
time of \(\Psi_4\) max for given mode and extraction radius
- Return type:
float
- get_psi4_phase_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Phase of \(\Psi_4\) for a given mode and extraction radius.
Returns the phase of \(\Psi_4\) for a given mode and extraction radius. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for \(\Psi_4\) data. If 0, the data at infinity is provided.
- Returns:
\(\Psi_4\) phase for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_psi4_real_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Real component of \(\Psi_4\) for a given mode and extraction radius.
Returns the real part of \(\Psi_4\) for a given mode and extraction radius. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for \(\Psi_4\) data. If 0, the data at infinity is provided.
- Returns:
\(\Psi_4\) real component for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_strain_amplitude_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Amplitude of \(rh\) for a given mode and extraction radius.
Returns the amplitude of strain for a given mode and extraction radius. The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for strain data. If 0, the data at infinity is provided.
- Returns:
amplitude of \(rh\) for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_strain_cross_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Cross component of \(rh\) for a given mode and extraction radius.
Returns the cross component of strain for a given mode and extraction radius. The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for strain data. If 0, the data at infinity is provided.
- Returns:
\(rh_{\times}\) for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_strain_phase_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Phase of \(rh\) for a given mode and extraction radius.
Returns the phase of strain for a given mode and extraction radius. The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for strain data. If 0, the data at infinity is provided.
- Returns:
phase of \(rh\) for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_strain_plus_for_mode(l: int, m: int, extraction_radius: float = 0) ndarray
Plus component of \(rh\) for a given mode and extraction radius.
Returns the plus component of strain for a given mode and extraction radius. The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration. If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
extraction_radius (
float
, optional) – extraction radius for strain data. If 0, the data at infinity is provided.
- Returns:
\(rh_+\) for a given mode and extraction radius.
- Return type:
numpy.ndarray
- get_strain_recomposed_at_sky_location(theta: float, phi: float, extraction_radius: float = 0) tuple
Plus and cross components of strain recomposed at a given sky location
The strain is recomposed by summing up the modes using spin weighted spherical harmonics as \(h(t,\theta,\phi) = \sum_{\ell,m} {}_{-2}Y_{\ell,m}(\theta, \phi) h_{ \ell,m}(t)\). If the extraction radius is 0, the data is extrapolated to infinite radius.
- Parameters:
theta (float) – \(0 \leq \theta \lt \pi\)
phi (float) – \(0 \leq \phi \lt 2\pi\)
extraction_radius (
float
, optional) – extraction radius for strain data. If 0, the data at infinity is provided.
- Returns:
\(rh_{+}\) and \(rh_{\times}\) recomposed at a given sky location
- Return type:
tuple
- get_time(extraction_radius: float = 0) ndarray
Time array associated with all radiation timeseries at the given radius.
- Parameters:
extraction_radius (
float
, optional) – Extraction radius for which the time data is requested.- Returns:
array containing the time data for all radiation information at the given radius
- Return type:
numpy.ndarray
- property included_modes: list
\(\Psi_4\) is decomposed using spherical harmonics labeled by (l, m). This provides a list of all (l,m) modes included.
- property included_radii: list
List of all extraction radii included.
- property l_max: int
Maximum l mode included.
- property radiation_spheres: dict
Dictionary containing all radius -> RadiationSphere pairs
- property radius_for_extrapolation: float
The radius from which extrapolation to infinite radius will be computed.
- set_frame(new_frame: Frame, time: ndarray | None = None, center_of_mass: ndarray | None = None)
Set the frame to use when decomposing the modes.
Options given by the frame enum and are the raw, original frame or the frame corrected for center of mass drift.
- Parameters:
new_frame (Frame) – Frame to transform to
time (
numpy.ndarray
, optional) – Time stamps for center of mass. Only necessary if moving to center of mass corrected frame.center_of_mass (
numpy.ndarray
, optional) – Time series of center of mass. Only necessary if moving to center of mass corrected frame.
- class mayawaves.radiation.RadiationSphere(mode_dict: dict, time: ndarray, radius: float, extrapolated: bool = False)
Class for interacting with the radiative data associated with a single extraction sphere.
- static create_radiation_sphere(radius_group: dict, radius: float, extrapolated: bool = False)
Create a RadiationSphere
- Parameters:
radius_group (h5py.Group) – group containing all the data pertaining to \(\Psi_4\) for all modes at the given extraction radius.
radius (float) – the radius at which all the \(\Psi_4\) data has been extracted.
extrapolated (
bool
, optional) – Whether the data has been extrapolated to infinite radius. Defaults to false.
- Returns:
a RadiationSphere consisting of all RadiationModes at the given extraction radius
- Return type:
- property extrapolated: bool
Whether the \(\Psi_4\) data has been extrapolated to infinite radius.
- property frame: Frame
Frame the modes are decomposed in.
Options are either the raw inertial frame or the frame corrected for center of mass drift.
- get_dEnergy_dt_radiated(**kwargs) tuple
Rate at which energy is radiated, \(dE/dt\)
Uses the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value- Returns:
time, \(dE/dt\) for given modes.
- Return type:
tuple
- get_dP_dt_radiated(**kwargs) tuple
Rate at which linear momentum is radiated
Rate at which linear momentum is radiated through gravitational waves, computed using the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value- Returns:
time, \(dP/dt\) for given modes.
- Return type:
tuple
- get_energy_radiated(**kwargs) tuple
Cumulative radiated energy \(E\)
Uses the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value- Returns:
time, cumulative energy radiated for given modes.
- Return type:
tuple
- get_extrapolated_sphere(order: int = 1)
RadiationSphere object extrapolated from this RadiationSphere to infinite radius.
Extrapolate the \(\Psi_4\) of all RadiationModes in this RadiationSphere to infinite radius and create and return a new RadiationSphere with those extrapolated RadiationModes. Uses the method described in https://arxiv.org/abs/1008.4360 and https://arxiv.org/abs/1108.4421.
- Parameters:
order (
int
, optional) – the extrapolation order. Defaults to 1.- Returns:
A RadiationSphere with the same properties as this one but with \(\Psi_4\) extrapolated to infinite radius.
- Return type:
- get_linear_momentum_radiated(**kwargs) tuple
Linear momentum radiated
Linear momentum radiated through gravitational waves computed using the method described in https://arxiv.org/abs/0707.4654. If no lmin, lmax, l, or m are provided, it computes the total sum of all modes.
- Parameters:
**kwargs –
parameters to specify the modes you would like to sum over
Valid keyword arguments are:
lmin (
int
, optional): minumum value of l rangelmax (
int
, optional): maximum value of l rangel (
int
, optional): specific l valuem (
int
, optional): specific m value- Returns:
time, cumulative linear momentum radiated for given modes
- Return type:
tuple
- get_mode(l: int, m: int)
RadiationMode with a given (l, m).
Returns the RadiationMode extracted on this RadiationSphere with the given spherical harmonic decomposition mode.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
RadiationMode associated with the given (l, m)
- Return type:
- get_psi4_amplitude_for_mode(l: int, m: int) ndarray
Amplitude of \(\Psi_4\) for a given mode.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
\(\Psi_4\) amplitude for a given mode
- Return type:
numpy.ndarray
- get_psi4_imaginary_for_mode(l: int, m: int) ndarray
Imaginary component of \(\Psi_4\) for a given mode.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
\(\Psi_4\) imaginary component for a given mode
- Return type:
numpy.ndarray
- get_psi4_max_time_for_mode(l: int, m: int) ndarray
Time of maximum \(\Psi_4\) amplitude for a given mode.
The time at which the amplitude of \(\Psi_4\) reaches its peak.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
time of \(\Psi_4\) max for given mode
- Return type:
float
- get_psi4_phase_for_mode(l: int, m: int) ndarray
Phase of \(\Psi_4\) for a given mode.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
\(\Psi_4\) phase for a given mode
- Return type:
numpy.ndarray
- get_psi4_real_for_mode(l: int, m: int) ndarray
Real component of \(\Psi_4\) for a given mode.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
\(\Psi_4\) real component for a given mode
- Return type:
numpy.ndarray
- get_strain_amplitude_for_mode(l: int, m: int) ndarray
Amplitude of \(rh\) for a given mode.
The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
amplitude of \(rh\) for a given mode
- Return type:
numpy.ndarray
- get_strain_cross_for_mode(l: int, m: int) ndarray
Cross component of \(rh\) for a given mode.
The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
\(rh_{\times}\) for a given mode
- Return type:
numpy.ndarray
- get_strain_phase_for_mode(l: int, m: int) ndarray
Phase of \(rh\) for a given mode.
The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
phase of \(rh\) for a given mode
- Return type:
numpy.ndarray
- get_strain_plus_for_mode(l: int, m: int) ndarray
Plus component of \(rh\) for a given mode.
The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- Parameters:
l (int) – l value of mode
m (int) – m value of mode
- Returns:
\(rh_+\) for a given mode
- Return type:
numpy.ndarray
- get_strain_recomposed_at_sky_location(theta: float, phi: float) tuple
Plus and cross components of strain recomposed at a given sky location
The strain is recomposed by summing up the modes using spin weighted spherical harmonics as \(h(t,\theta,\phi) = \sum_{\ell,m} {}_{-2}Y_{\ell,m}(\theta, \phi) h_{ \ell,m}(t)\)
- Parameters:
theta (float) – \(0 \leq \theta \lt \pi\)
phi (float) – \(0 \leq \phi \lt 2\pi\)
- Returns:
\(rh_+\) and \(rh_{\times}\) recomposed at a given sky location
- Return type:
tuple
- property included_modes: list
\(\Psi_4\) is decomposed into spherical harmonics labeled by (l, m). This provides a list of all (l, m) modes included.
- property l_max: int
Maximum l mode included.
- property modes: dict
Dictionary containing all the (l, m) -> RadiationMode objects in the current frame.
- property radius: float
Radius at which the \(\Psi_4\) data was extracted.
- property raw_modes: dict
Dictionary containing all the (l, m) -> RadiationMode objects in their original frame.
- set_frame(new_frame: Frame, time: ndarray | None = None, center_of_mass: ndarray | None = None, alpha: ndarray | None = None, beta: ndarray | None = None)
Set the frame to use when decomposing the modes.
Options given by the frame enum and are the raw, original frame or corrected for center of mass drift.
- Parameters:
new_frame (Frame) – Frame to transform to
time (
numpy.ndarray
, optional) – Time stamps for center of mass. Only necessary if moving to center of mass corrected frame and not providing alpha and beta.center_of_mass (
numpy.ndarray
, optional) – Time series of center of mass. Only necessary if moving to center of mass corrected frame and not providing alpha and beta.alpha (
numpy.ndarray
, optional) – Offset for center of mass correction. Only necessary if moving to center of mass corrected frame and not providing the center of mass timeseries.beta (
numpy.ndarray
, optional) – Boost for center of mass correction. Only necessary if moving to center of mass corrected frame and not providing the center of mass timeseries.
- property time: ndarray
Time array associated with all timeseries provided by this RadiationSphere.
- class mayawaves.radiation.RadiationMode(l: int, m: int, rad: float, time: ndarray, psi4_group: Group | None = None, extrapolated: bool = False, radiation_sphere: RadiationSphere | None = None, psi4_real: ndarray | None = None, psi4_imaginary: ndarray | None = None, strain_plus: ndarray | None = None, strain_cross: ndarray | None = None, center_of_mass_corrected: bool = False)
Class for interacting with a single mode of radiation information when expressed in terms of spin weighted spherical harmonics.
- compute_and_store_strain()
Computes and stores the strain data from the \(\Psi_4\) data.
Computes the strain using the fixed-frequency integration method to compute the double time integral of \(\Psi_4\) using the fourier transform.
- static create_radiation_mode(psi4_group: Group | None = None, l: int = 0, m: int = 0, rad: float = 0, time: ndarray | None = None, extrapolated: bool = False, radiation_sphere: RadiationSphere | None = None)
Create a RadiationMode from psi4 group.
- Parameters:
psi4_group (
h5py.Group
, optional) – group containing all the data pertaining to \(\Psi_4\) for the given model (
int
, optional) – value of l for the spherical harmonic mode. Defaults to 0.m (
int
, optional) – value of m for the spherical harmonic mode. Defaults to 0.rad (
float
, optional) – radius at which the given mode was extracted. Defaults to 0.time (
numpy.ndarray
, optional) – array containing the time value associated with each data point of \(\Psi_4\). Defaults to empty array.extrapolated (
bool
, optional) – whether the data has been extrapolated to infinite radius. Defaults to false.radiation_sphere (
RadiationSphere
, optional) – the RadiationSphere this mode is associated with. Defaults to None.
- Returns:
The constructed RadiationMode
- Return type:
- property dEnergy_dt_radiated: ndarray
Time derivative of the energy radiated.
https://arxiv.org/pdf/0707.4654.pdf eq 3.8
- Returns:
\(dE/dt\) radiated
- Return type:
numpy.ndarray
- property dP_dt_radiated: ndarray
Rate at which linear momentum is radiated as gravitational waves
https://arxiv.org/pdf/0707.4654.pdf eqs 3.14, 3.15
- Returns:
\(dP/dt\) radiated
- Return type:
numpy.ndarray
- property energy_radiated: ndarray
Energy radiated, \(E\).
https://arxiv.org/pdf/0707.4654.pdf eq 3.8
- Returns:
cumulative \(E\) radiated
- Return type:
numpy.ndarray
- property extrapolated: bool
Whether the data within this mode has been extrapolated to infinite radius.
- get_mode_with_extrapolated_radius(order: int = 1)
RadiationMode object extrapolated from this RadiationMode to infinite radius.
Extrapolate the \(\Psi_4\) of this RadiationMode to infinite radius and create and return a new RadiationMode with that \(\Psi_4\) data. Uses the method described in https://arxiv.org/abs/1008.4360 and https://arxiv.org/abs/1108.4421.
- Parameters:
order (
int
, optional) – the extrapolation order. Defaults to 1.- Returns:
A RadiationMode with the same properties as this one but with \(\Psi_4\) extrapolated to infinite radius.
- Return type:
- property h_cross_dot
The imaginary part of the first time derivative of the strain.
- property h_plus_dot
The real part of the first time derivative of the strain.
- property l_value: int
l value in the (l, m) pair that denotes a mode of the spherical harmonic decomposition.
- property linear_momentum_radiated: ndarray
Linear momentum radiated as gravitational waves
https://arxiv.org/pdf/0707.4654.pdf eqs 3.14, 3.15
- Returns:
linear momentum radiated
- Return type:
numpy.ndarray
- property m_value: int
m value in the (l, m) pair that denotes a mode of the spherical harmonic decomposition.
- property omega_start
The starting frequency of \(\Psi_4\) where the frequency is computed as the time derivative of the phase.
- property psi4_amplitude: ndarray
The amplitude of the Weyl Scalar (\(\Psi_4\)) as a timeseries, computed from the real and imaginary parts as \(\sqrt{\mathcal{Re}(\Psi_4)^2 + \mathcal{Im}(\Psi_4)^2}\).
- property psi4_imaginary: ndarray
The imaginary part of the Weyl Scalar (\(\Psi_4\)) as a timeseries.
- property psi4_max_time
The time at which the amplitude of \(\Psi_4\) is at a maximum
- property psi4_omega
The frequency of \(\Psi_4\) computed as the time derivative of the phase.
- property psi4_phase: ndarray
The phase of the Weyl Scalar (\(\Psi_4\)) as a timeseries, computed from the real and imaginary parts as \(\textrm{tan}^{-1}( \mathcal{Im}(\Psi_4) / \mathcal{Re}(\Psi_4) )\).
- property psi4_real: ndarray
The real part of the Weyl Scalar (\(\Psi_4\)) as a timeseries.
- property radiation_sphere: RadiationSphere
The RadiationSphere this mode is associated with.
- property radius: float
Radius at which \(\Psi_4\) for this mode was extracted.
- property strain_amplitude: ndarray
The amplitude of \(rh\) as a timeseries, computed from the real and imaginary parts as \(\sqrt{\mathcal{Re}(rh)^2 + \mathcal{Im}(rh)^2}\). The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- property strain_cross: ndarray
The cross component of \(rh\) as a timeseries. The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- property strain_phase: ndarray
The phase of \(rh\) as a timeseries, computed from the real and imaginary parts as \(\textrm{tan}^{-1}( \mathcal{Im}(rh) / \mathcal{Re}(rh) )\). The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- property strain_plus: ndarray
The plus component of \(rh\) as a timeseries. The strain is the second time integral of \(\Psi_4\) computed using fixed-frequency integration.
- property time: ndarray
Array of time stamps associated with the \(\Psi_4\) and strain data.
- static ylm(l, m, theta, phi) float
Spherical harmonic for l, m at given angles theta and phi
- Parameters:
l – l value of mode
m – m value of mode
theta – angle defined from the +z-axis
phi – angle defined in the x-y plane moving counterclockwise from the +x-axis
- Returns:
The spherical harmonic value at the given l, m, theta, phi
- Return type:
float