scri.extrapolation

Functions

ErroredExtrapolations(TopLevelOutputDir, ...)

Find directories with errors.

FindPossibleExtrapolationsToRun(TopLevelInputDir)

Find all possible extrapolations.

NewerDataThanExtrapolation(TopLevelInputDir, ...)

Find newer data than extrapolation.

RunExtrapolation(TopLevelInputDir, ...)

StartedButUnfinishedExtrapolations(...)

Find directories with extrapolations that started but didn't finish.

UnstartedExtrapolations(TopLevelOutputDir, ...)

Find unstarted extrapolation directories.

datatype_from_filename(filename)

extrapolate(**kwargs)

Perform extrapolations from finite-radius data.

intersection(t1, t2[, min_step, min_time, ...])

Return the intersection of two time sequences.

monotonic_indices(T[, MinTimeStep])

Given an array of times, return the indices that make the array strictly monotonic.

pick_Ch_mass([filename])

Deduce the best Christodoulou mass by finding the statistical "mode" (after binning).

read_finite_radius_data([ChMass, filename, ...])

Read data at various radii, and offset by tortoise coordinate.

read_finite_radius_waveform(filename, ...)

read_finite_radius_waveform_nrar(filename, ...)

read_finite_radius_waveform_rpxmb_or_rpdmb(...)

This is just a worker function defined for read_finite_radius_data, below, reading a single waveform from an h5 file of many waveforms.

set_common_time(Ws, Radii[, MinTimeStep, ...])

Interpolate Waveforms and radius data to a common set of times.

validate_group_of_waveforms(h5file, ...)

validate_single_waveform(h5file, filename, ...)

scri.extrapolation.ErroredExtrapolations(TopLevelOutputDir, SubdirectoriesAndDataFiles)[source]

Find directories with errors.

scri.extrapolation.FindPossibleExtrapolationsToRun(TopLevelInputDir)[source]

Find all possible extrapolations.

scri.extrapolation.NewerDataThanExtrapolation(TopLevelInputDir, TopLevelOutputDir, SubdirectoriesAndDataFiles)[source]

Find newer data than extrapolation.

scri.extrapolation.RunExtrapolation(TopLevelInputDir, TopLevelOutputDir, Subdirectory, DataFile, Template)[source]
scri.extrapolation.StartedButUnfinishedExtrapolations(TopLevelOutputDir, SubdirectoriesAndDataFiles)[source]

Find directories with extrapolations that started but didn’t finish.

scri.extrapolation.UnstartedExtrapolations(TopLevelOutputDir, SubdirectoriesAndDataFiles)[source]

Find unstarted extrapolation directories.

scri.extrapolation.datatype_from_filename(filename)[source]
scri.extrapolation.extrapolate(**kwargs)[source]

Perform extrapolations from finite-radius data.

See arXiv:2010.15200 for specific details.

Parameters:
InputDirectorystr, (Default: ‘./’)

Where to find the input data. Can be relative or absolute.

OutputDirectorystr, (Default: ‘./’)

This directory will be made if it does not exist. If this is passed but is the empty string, no files of any kind will be output.

DataFilestr, (Default: ‘rh_FiniteRadii_CodeUnits.h5’)

Input file holding the data from all the radii.

ChMassfloat, (Default: 0.0)

Christodoulou mass in the same units as the rest of the data. All the data will be rescaled into units such that this is one. If this is zero, the Christodoulou mass will be extracted automatically from the horizons file below.

HorizonsFilestr, (Default: ‘Horizons.h5’)

File name to read for horizon data (if ChMass is 0.0).

CoordRadiilist of str or list of int, (Default: [])

List of strings containing the radii to use, or of (integer) indices of the list of waveform names. If this is a list of indices, the order is just the order output by the command list(h5py.File(DataFile)) which should be the same as h5ls. If the list is empty, all radii that can be found are used.

ExtrapolationOrderslist of int, (Default: [-1, 2, 3, 4, 5, 6])

Negative numbers correspond to extracted data, counting down from the outermost extraction radius (which is -1).

UseOmegabool, (Default: False)

Whether or not to extrapolate as a function of lambda/r = 1/(r*m*omega), where omega is the instantaneous angular frequency of rotation. If this is True, the extrapolation will usually not converge for high N; if this is False, SVD will generally cause the convergence to appear to fall to roundoff, though the accuracy presumably is not so great.

OutputFrame{scri.Inertial, scri.Corotating}

Transform to this frame before comparison and output.

ExtrapolatedFilesstr, (Default: ‘Extrapolated_N{N}.h5’)
DifferenceFilesstr, (Default: ‘ExtrapConvergence_N{N}-N{Nm1}.h5’)

These are python-formatted output file names, where the extrapolation order N is substituted for ‘{N}’, and the previous extrapolation order is substituted for ‘{Nm1}’. The data-type inferred from the DataFile name is prepended. If DifferenceFiles is empty, the corresponding file is not output.

UseStupidNRARFormatbool, (Default: False)

If True (and ExtrapolatedFiles does not end in ‘.dat’), then the h5 output format will be that stupid, old NRAR/NINJA format that doesn’t convey enough information, is slow, and uses 33% more space than it needs to. But you know, if you’re into that kind of thing, whatever. Who am I to judge?

PlotFormatstr, (Default: ‘’)

The format of output plots. This can be the empty string, in which case no plotting is done. Or, these can be any of the formats supported by your installation of matplotlib.

MinTimeStepfloat, (Default: 0.005)

The smallest allowed time step in the output data.

EarliestTimefloat, (Default: -3.0e300)

The earliest time in the output data. For values less than 0, some of the data corresponds to times when only junk radiation is present.

LatestTimefloat, (Default: 3.0e300)

The latest time in the output data.

AlignmentTimefloat, (Default: None)

The time at which to align the Waveform with the dominant eigenvector of <LL>. If the input value is None or is outside of the input data, it will be reset to the midpoint of the waveform: (W_outer.T(0)+W_outer.T(-1))/2

NoiseFloorfloat, (Default: None)

ONLY USED FOR Psi1 AND Psi0 WAVEFORMS. The effective noise floor of the SpEC simulation. If Psi1 or Psi0 (NOT scaled by radius) falls beneath this value for certain extraction radii, then those radii are not included in the extrapolation. The value 1.0e-9 seems to work well for a few BBH systems.

scri.extrapolation.intersection(t1, t2, min_step=None, min_time=None, max_time=None)[source]

Return the intersection of two time sequences.

The time step at each point is the minimum of the time steps in t1 and t2 at that instant, or min_step, whichever is greater. The output starts at the earliest moment common to t1 and t2, or min_time, whichever is greater.

Parameters:
t1: 1-d float array
t2: 1-d float array
min_step: float
min_time: float
max_time: float
Returns:
1-d float array
scri.extrapolation.monotonic_indices(T, MinTimeStep=0.001)[source]

Given an array of times, return the indices that make the array strictly monotonic.

scri.extrapolation.pick_Ch_mass(filename='Horizons.h5')[source]

Deduce the best Christodoulou mass by finding the statistical “mode” (after binning).

scri.extrapolation.read_finite_radius_data(ChMass=0.0, filename='rh_FiniteRadii_CodeUnits.h5', CoordRadii=[])[source]

Read data at various radii, and offset by tortoise coordinate.

scri.extrapolation.read_finite_radius_waveform(filename, groupname, WaveformName, ChMass)[source]
scri.extrapolation.read_finite_radius_waveform_nrar(filename, WaveformName)[source]
scri.extrapolation.read_finite_radius_waveform_rpxmb_or_rpdmb(filename, groupname, WaveformName)[source]

This is just a worker function defined for read_finite_radius_data, below, reading a single waveform from an h5 file of many waveforms. You probably don’t need to call this directly.

scri.extrapolation.set_common_time(Ws, Radii, MinTimeStep=0.005, EarliestTime=-3e+300, LatestTime=3e+300)[source]

Interpolate Waveforms and radius data to a common set of times.

scri.extrapolation.validate_group_of_waveforms(h5file, filename, WaveformNames)[source]
scri.extrapolation.validate_single_waveform(h5file, filename, WaveformName, ExpectedNModes, ExpectedNTimes, LModes)[source]