pycbc.waveform package

Submodules

pycbc.waveform.bank module

This module provides classes that describe banks of waveforms

class pycbc.waveform.bank.FilterBank(filename, filter_length, delta_f, dtype, out=None, max_template_length=None, approximant=None, parameters=None, enable_compressed_waveforms=True, low_frequency_cutoff=None, waveform_decompression_method=None, **kwds)[source]

Bases: TemplateBank

generate_with_delta_f_and_max_freq(t_num, max_freq, delta_f, low_frequency_cutoff=None, cached_mem=None)[source]

Generate the template with index t_num using custom length.

get_decompressed_waveform(tempout, index, f_lower=None, approximant=None, df=None)[source]

Returns a frequency domain decompressed waveform for the template in the bank corresponding to the index taken in as an argument. The decompressed waveform is obtained by interpolating in frequency space, the amplitude and phase points for the compressed template that are read in from the bank.

class pycbc.waveform.bank.FilterBankSkyMax(filename, filter_length, delta_f, dtype, out_plus=None, out_cross=None, max_template_length=None, parameters=None, low_frequency_cutoff=None, **kwds)[source]

Bases: TemplateBank

class pycbc.waveform.bank.LiveFilterBank(filename, sample_rate, minimum_buffer, approximant=None, increment=8, parameters=None, low_frequency_cutoff=None, **kwds)[source]

Bases: TemplateBank

get_template(index, min_buffer=None)[source]
getslice(sindex)[source]
id_from_param(param_tuple)[source]

Get the index of this template based on its param tuple

Parameters:

param_tuple (tuple) – Tuple of the parameters which uniquely identify this template

Returns:

index – The ordered index that this template has in the template bank.

Return type:

int

round_up(num)[source]

Determine the length to use for this waveform by rounding.

Parameters:

num (int) – Proposed size of waveform in samples.

Returns:

size – The rounded size to use for the waveform buffer in samples. This is calculated using an internal increment attribute, which determines the discreteness of the rounding.

Return type:

int

class pycbc.waveform.bank.TemplateBank(filename, approximant=None, parameters=None, **kwds)[source]

Bases: object

Class to provide some basic helper functions and information about elements of a template bank.

Parameters:
  • filename (string) – The name of the file to load. Must end in ‘.xml[.gz]’ or ‘.hdf’. If an hdf file, it should have a ‘parameters’ in its attrs which gives a list of the names of fields to load from the file. If no ‘parameters’ are found, all of the top-level groups in the file will assumed to be parameters (a warning will be printed to stdout in this case). If an xml file, it must have a SnglInspiral table.

  • approximant ({None, (list of) string(s)}) – Specify the approximant(s) for each template in the bank. If None provided, will try to load the approximant from the file. The approximant may either be a single string (in which case the same approximant will be used for all templates) or a list of strings and conditionals specifying where to use the approximant. See boolargs_from_apprxstr for syntax.

  • parameters ({None, (list of) sting(s)}) – Specify what parameters to load from the file. If None, all of the parameters in the file (if an xml file, this is all of the columns in the SnglInspiral table, if an hdf file, this is given by the parameters attribute in the file). The list may include parameters that are derived from the file’s parameters, or functions thereof. For a full list of possible parameters, see WaveformArray.default_fields. If a derived parameter is specified, only the parameters needed to compute that parameter will be loaded from the file. For example, if parameters=’mchirp’, then only mass1, mass2 will be loaded from the file. Note that derived parameters can only be used if the needed parameters are in the file; e.g., you cannot use chi_eff if spin1z, spin2z, mass1, and mass2 are in the input file.

  • **kwds – Any additional keyword arguments are stored to the extra_args attribute.

table

An instance of a WaveformArray containing all of the information about the parameters of the bank.

Type:

WaveformArray

has_compressed_waveforms

True if compressed waveforms are present in the the (hdf) file; False otherwise.

Type:

{False, bool}

indoc

If an xml file was provided, an in-memory representation of the xml. Otherwise, None.

Type:

{None, xmldoc}

filehandler

If an hdf file was provided, the file handler pointing to the hdf file (left open after initialization). Otherwise, None.

Type:

{None, h5py.File}

extra_args

Any extra keyword arguments that were provided on initialization.

Type:

{None, dict}

approximant(index)[source]

Return the name of the approximant ot use at the given index

end_frequency(index)[source]

Return the end frequency of the waveform at the given index value

ensure_hash()[source]

Ensure that there is a correctly populated template_hash.

Check for a correctly populated template_hash and create if it doesn’t already exist.

ensure_standard_filter_columns(low_frequency_cutoff=None)[source]

Initialize FilterBank common fields

Parameters:

low_frequency_cutoff ({float, None}, Optional) – A low frequency cutoff which overrides any given within the template bank file.

property parameters

The parameters loaded from the input file. Same as table.fieldnames.

Type:

tuple

parse_approximant(approximant)[source]

Parses the given approximant argument, returning the approximant to use for each template in self. This is done by calling parse_approximant_arg using self’s table as the array; see that function for more details.

template_thinning(inj_filter_rejector)[source]

Remove templates from bank that are far from all injections.

write_to_hdf(filename, start_index=None, stop_index=None, force=False, skip_fields=None, write_compressed_waveforms=True)[source]

Writes self to the given hdf file.

Parameters:
  • filename (str) – The name of the file to write to. Must be a recognised HDF5 file extension

  • start_index (If a specific slice of the template bank is to be) – written to the hdf file, this would specify the index of the first template in the slice

  • stop_index (If a specific slice of the template bank is to be) – written to the hdf file, this would specify the index of the last template in the slice

  • force ({False, bool}) – If the file already exists, it will be overwritten if True. Otherwise, an OSError is raised if the file exists.

  • skip_fields ({None, (list of) strings}) – Do not write the given fields to the hdf file. Default is None, in which case all fields in self.table.fieldnames are written.

  • write_compressed_waveforms ({True, bool}) – Write compressed waveforms to the output (hdf) file if this is True, which is the default setting. If False, do not write the compressed waveforms group, but only the template parameters to the output file.

Returns:

The file handler to the output hdf file (left open).

Return type:

h5py.File

pycbc.waveform.bank.add_approximant_arg(parser, default=None, help=None)[source]

Adds an approximant argument to the given parser.

Parameters:
  • parser (ArgumentParser) – The argument parser to add the argument to.

  • default ({None, str}) – Specify a default for the approximant argument. Defaults to None.

  • help ({None, str}) – Provide a custom help message. If None, will use a descriptive message on how to specify the approximant.

pycbc.waveform.bank.boolargs_from_apprxstr(approximant_strs)[source]

Parses a list of strings specifying an approximant and where that approximant should be used into a list that can be understood by FieldArray.parse_boolargs.

Parameters:

apprxstr ((list of) string(s)) – The strings to parse. Each string should be formatted APPRX:COND, where APPRX is the approximant and COND is a string specifying where it should be applied (see FieldArgs.parse_boolargs for examples of conditional strings). The last string in the list may exclude a conditional argument, which is the same as specifying ‘:else’.

Returns:

boolargs – A list of tuples giving the approximant and where to apply them. This can be passed directly to FieldArray.parse_boolargs.

Return type:

list

pycbc.waveform.bank.find_variable_start_frequency(approximant, parameters, f_start, max_length, delta_f=1)[source]

Find a frequency value above the starting frequency that results in a waveform shorter than max_length.

pycbc.waveform.bank.parse_approximant_arg(approximant_arg, warray)[source]

Given an approximant arg (see add_approximant_arg) and a field array, figures out what approximant to use for each template in the array.

Parameters:
  • approximant_arg (list) – The approximant argument to parse. Should be the thing returned by ArgumentParser when parsing the argument added by add_approximant_arg.

  • warray (FieldArray) – The array to parse. Must be an instance of a FieldArray, or a class that inherits from FieldArray.

Returns:

A numpy array listing the approximants to use for each element in the warray.

Return type:

array

pycbc.waveform.bank.sigma_cached(self, psd)[source]

Cache sigma calculate for use in tandem with the FilterBank class

pycbc.waveform.bank.tuple_to_hash(tuple_to_be_hashed)[source]

Return a hash for a numpy array, avoids native (unsafe) python3 hash function

Parameters:

tuple_to_be_hashed (tuple) – The tuple which is being hashed Must be convertible to a numpy array

Returns:

an integer representation of the hashed array

Return type:

int

pycbc.waveform.compress module

Utilities for handling frequency compressed an unequally spaced frequency domain waveforms.

class pycbc.waveform.compress.CompressedWaveform(sample_points, amplitude, phase, interpolation=None, tolerance=None, mismatch=None, precision='double', load_to_memory=True)[source]

Bases: object

Class that stores information about a compressed waveform.

Parameters:
  • sample_points ({array, h5py.Dataset}) – The frequency points at which the compressed waveform is sampled.

  • amplitude ({array, h5py.Dataset}) – The amplitude of the waveform at the given sample_points.

  • phase ({array, h5py.Dataset}) – The phase of the waveform at the given sample_points.

  • interpolation ({None, str}) – The interpolation that was used when compressing the waveform for computing tolerance. This is also the default interpolation used when decompressing; see decompress for details.

  • tolerance ({None, float}) – The tolerance that was used when compressing the waveform.

  • mismatch ({None, float}) – The actual mismatch between the decompressed waveform (using the given interpolation) and the full waveform.

  • precision ({'double', str}) – The precision used to generate the compressed waveform’s amplitude and phase points. Default is ‘double’.

  • load_to_memory ({True, bool}) – If sample_points, amplitude, and/or phase is an hdf dataset, they will be cached in memory the first time they are accessed. Default is True.

load_to_memory

Whether or not to load sample_points/amplitude/phase into memory the first time they are accessed, if they are hdf datasets. Can be set directly to toggle this behavior.

Type:

bool

interpolation

The interpolation that was used when compressing the waveform, for checking the mismatch. Also the default interpolation used when decompressing.

Type:

str

tolerance

The tolerance that was used when compressing the waveform.

Type:

{None, float}

mismatch

The mismatch between the decompressed waveform and the original waveform.

Type:

{None, float}

precision

The precision used to generate and store the compressed waveform points. Options are ‘double’ or ‘single’; default is ‘double

Type:

{‘double’, str}

property amplitude

The amplitude of the waveform at the sample_points.

This is always returned as an array; the same logic as for sample_points is used to determine whether or not to cache in memory.

Returns:

amplitude

Return type:

Array

clear_cache()[source]

Clear self’s cache of amplitude, phase, and sample_points.

decompress(out=None, df=None, f_lower=None, interpolation=None)[source]

Decompress self.

Parameters:
  • out ({None, FrequencySeries}) – Write the decompressed waveform to the given frequency series. The decompressed waveform will have the same delta_f as out. Either this or df must be provided.

  • df ({None, float}) – Decompress the waveform such that its delta_f has the given value. Either this or out must be provided.

  • f_lower ({None, float}) – The starting frequency at which to decompress the waveform. Cannot be less than the minimum frequency in sample_points. If None provided, will default to the minimum frequency in sample_points.

  • interpolation ({None, str}) – The interpolation to use for decompressing the waveform. If None provided, will default to self.interpolation.

Returns:

The decompressed waveform.

Return type:

FrequencySeries

classmethod from_hdf(fp, template_hash, root=None, load_to_memory=True, load_now=False)[source]

Load a compressed waveform from the given hdf file handler.

The waveform is retrieved from: fp[‘[{root}/]compressed_waveforms/{template_hash}/{param}’], where param is the sample_points, amplitude, and phase.

Parameters:
  • fp (h5py.File) – An open hdf file to write the compressed waveform to.

  • template_hash ({hash, int, str}) – The id of the waveform.

  • root ({None, str}) – Retrieve the compressed_waveforms group from the given string. If None, compressed_waveforms will be assumed to be in the top level.

  • load_to_memory ({True, bool}) – Set the load_to_memory attribute to the given value in the returned instance.

  • load_now ({False, bool}) – Immediately load the sample_points/amplitude/phase to memory.

Returns:

An instance of this class with parameters loaded from the hdf file.

Return type:

CompressedWaveform

property phase

The phase of the waveform as the sample_points.

This is always returned as an array; the same logic as for sample_points returned as an array; the same logic as for sample_points is used to determine whether or not to cache in memory.

Returns:

phase

Return type:

Array

property sample_points

The frequencies at which the compressed waveform is sampled.

This is always returned as an array, even if the stored sample_points is an hdf dataset. If load_to_memory is True and the stored points are an hdf dataset, the sample_points will cached in memory the first time this attribute is accessed.

Returns:

sample_points

Return type:

Array

write_to_hdf(fp, template_hash, root=None, precision=None)[source]

Write the compressed waveform to the given hdf file handler.

The waveform is written to: fp[‘[{root}/]compressed_waveforms/{template_hash}/{param}’], where param is the sample_points, amplitude, and phase. The interpolation, tolerance, mismatch and precision are saved to the group’s attributes.

Parameters:
  • fp (h5py.File) – An open hdf file to write the compressed waveform to.

  • template_hash ({hash, int, str}) – A hash, int, or string to map the template to the waveform.

  • root ({None, str}) – Put the compressed_waveforms group in the given directory in the hdf file. If None, compressed_waveforms will be the root directory.

  • precision ({None, str}) – Cast the saved parameters to the given precision before saving. If None provided, will use whatever their current precision is. This will raise an error if the parameters have single precision but the requested precision is double.

pycbc.waveform.compress.compress_waveform(htilde, sample_points, tolerance, interpolation, precision, decomp_scratch=None, psd=None)[source]

Retrieves the amplitude and phase at the desired sample points, and adds frequency points in order to ensure that the interpolated waveform has a mismatch with the full waveform that is <= the desired tolerance. The mismatch is computed by finding 1-overlap between htilde and the decompressed waveform; no maximimization over phase/time is done, a PSD may be used.

Note

The decompressed waveform is only garaunteed to have a true mismatch <= the tolerance for the given interpolation and for no PSD. However, since no maximization over time/phase is performed when adding points, the actual mismatch between the decompressed waveform and htilde is better than the tolerance, using no PSD. Using a PSD does increase the mismatch, and can lead to mismatches > than the desired tolerance, but typically by only a factor of a few worse.

Parameters:
  • htilde (FrequencySeries) – The waveform to compress.

  • sample_points (array) – The frequencies at which to store the amplitude and phase. More points may be added to this, depending on the desired tolerance.

  • tolerance (float) – The maximum mismatch to allow between a decompressed waveform and htilde.

  • interpolation (str) – The interpolation to use for decompressing the waveform when computing overlaps.

  • precision (str) – The precision being used to generate and store the compressed waveform points.

  • decomp_scratch ({None, FrequencySeries}) – Optionally provide scratch space for decompressing the waveform. The provided frequency series must have the same delta_f and length as htilde.

  • psd ({None, FrequencySeries}) – The psd to use for calculating the overlap between the decompressed waveform and the original full waveform.

Returns:

The compressed waveform data; see CompressedWaveform for details.

Return type:

CompressedWaveform

pycbc.waveform.compress.fd_decompress(amp, phase, sample_frequencies, out=None, df=None, f_lower=None, interpolation='inline_linear')[source]

Decompresses an FD waveform using the given amplitude, phase, and the frequencies at which they are sampled at.

Parameters:
  • amp (array) – The amplitude of the waveform at the sample frequencies.

  • phase (array) – The phase of the waveform at the sample frequencies.

  • sample_frequencies (array) – The frequency (in Hz) of the waveform at the sample frequencies.

  • out ({None, FrequencySeries}) – The output array to save the decompressed waveform to. If this contains slots for frequencies > the maximum frequency in sample_frequencies, the rest of the values are zeroed. If not provided, must provide a df.

  • df ({None, float}) – The frequency step to use for the decompressed waveform. Must be provided if out is None.

  • f_lower ({None, float}) – The frequency to start the decompression at. If None, will use whatever the lowest frequency is in sample_frequencies. All values at frequencies less than this will be 0 in the decompressed waveform.

  • interpolation ({'inline_linear', str}) – The interpolation to use for the amplitude and phase. Default is ‘inline_linear’. If ‘inline_linear’ a custom interpolater is used. Otherwise, scipy.interpolate.interp1d is used; for other options, see possible values for that function’s kind argument.

Returns:

out – If out was provided, writes to that array. Otherwise, a new FrequencySeries with the decompressed waveform.

Return type:

FrequencySeries

pycbc.waveform.compress.inline_linear_interp(amp, phase, sample_frequencies, output, df, f_lower, imin, start_index)[source]

Generate a frequency-domain waveform via linear interpolation from sampled amplitude and phase. The sample frequency locations for the amplitude and phase must be the same. This function may be less accurate than scipy’s linear interpolation, but should be much faster. Additionally, it is ‘schemed’ and so may run under either CPU or GPU schemes.

This function is not ordinarily called directly, but rather by giving the argument ‘interpolation’ the value ‘inline_linear’ when calling the function ‘fd_decompress’ below.

Parameters:
  • amp (array) – The amplitude of the waveform at the sample frequencies.

  • phase (array) – The phase of the waveform at the sample frequencies.

  • sample_frequencies (array) – The frequency (in Hz) of the waveform at the sample frequencies.

  • output ({None, FrequencySeries}) – The output array to save the decompressed waveform to. If this contains slots for frequencies > the maximum frequency in sample_frequencies, the rest of the values are zeroed. If not provided, must provide a df.

  • df ({None, float}) – The frequency step to use for the decompressed waveform. Must be provided if out is None.

  • f_lower (float) – The frequency to start the decompression at. All values at frequencies less than this will be 0 in the decompressed waveform.

  • imin (int) – The index at which to start in the sampled frequency series. Must therefore be 0 <= imin < len(sample_frequencies)

  • start_index (int) – The index at which to start in the output frequency; i.e., ceil(f_lower/df).

Returns:

output – If out was provided, writes to that array. Otherwise, a new FrequencySeries with the decompressed waveform.

Return type:

FrequencySeries

pycbc.waveform.compress.mchirp_compression(m1, m2, fmin, fmax, min_seglen=0.02, df_multiple=None)[source]

Return the frequencies needed to compress a waveform with the given chirp mass. This is based on the estimate in rough_time_estimate.

Parameters:
  • m1 (float) – mass of first component object in solar masses

  • m2 (float) – mass of second component object in solar masses

  • fmin (float) – The starting frequency of the compressed waveform.

  • fmax (float) – The ending frequency of the compressed waveform.

  • min_seglen (float) – The inverse of this gives the maximum frequency step that is used.

  • df_multiple ({None, float}) – Make the compressed sampling frequencies a multiple of the given value. If None provided, the returned sample points can have any floating point value.

Returns:

The frequencies at which to evaluate the compressed waveform.

Return type:

array

pycbc.waveform.compress.rough_time_estimate(m1, m2, flow, fudge_length=1.1, fudge_min=0.02)[source]

A very rough estimate of the duration of the waveform.

An estimate of the waveform duration starting from flow. This is intended to be fast but not necessarily accurate. It should be an overestimate of the length. It is derived from a simplification of the 0PN post-newtonian terms and includes a fudge factor for possible ringdown, etc.

Parameters:
  • m1 (float) – mass of first component object in solar masses

  • m2 (float) – mass of second component object in solar masses

  • flow (float) – starting frequency of the waveform

  • fudge_length (optional, {1.1, float}) – Factor to multiply length estimate by to ensure it is a convservative value

  • fudge_min (optional, {0.2, float}) – Minimum signal duration that can be returned. This should be long enough to encompass the ringdown and errors in the precise end time.

Returns:

time – Time from flow untill the end of the waveform

Return type:

float

pycbc.waveform.compress.spa_compression(htilde, fmin, fmax, min_seglen=0.02, sample_frequencies=None)[source]

Returns the frequencies needed to compress the given frequency domain waveform. This is done by estimating t(f) of the waveform using the stationary phase approximation.

Parameters:
  • htilde (FrequencySeries) – The waveform to compress.

  • fmin (float) – The starting frequency of the compressed waveform.

  • fmax (float) – The ending frequency of the compressed waveform.

  • min_seglen (float) – The inverse of this gives the maximum frequency step that is used.

  • sample_frequencies ({None, array}) – The frequencies that the waveform is evaluated at. If None, will retrieve the frequencies from the waveform’s sample_frequencies attribute.

Returns:

The frequencies at which to evaluate the compressed waveform.

Return type:

array

pycbc.waveform.compress.vecdiff(htilde, hinterp, sample_points, psd=None)[source]

Computes a statistic indicating between which sample points a waveform and the interpolated waveform differ the most.

pycbc.waveform.decompress_cpu module

Utilities for handling frequency compressed an unequally spaced frequency domain waveforms.

pycbc.waveform.decompress_cpu.inline_linear_interp(amp, phase, sample_frequencies, output, df, f_lower, imin, start_index)[source]

pycbc.waveform.decompress_cpu_cython module

pycbc.waveform.decompress_cpu_cython.decomp_ccode_double(h, delta_f, hlen, start_index, sample_frequencies, amp, phase, sflen, imin)
pycbc.waveform.decompress_cpu_cython.decomp_ccode_float(h, delta_f, hlen, start_index, sample_frequencies, amp, phase, sflen, imin)

pycbc.waveform.generator module

This modules provides classes for generating waveforms.

class pycbc.waveform.generator.BaseCBCGenerator(generator, variable_args=(), **frozen_params)[source]

Bases: BaseGenerator

Adds ability to convert from various derived parameters to parameters needed by the waveform generators.

possible_args = {'amplitude_order', 'approximant', 'coa_phase', 'dalpha1', 'dalpha2', 'dalpha3', 'dalpha4', 'dalpha5', 'dbeta1', 'dbeta2', 'dbeta3', 'dchi0', 'dchi1', 'dchi2', 'dchi3', 'dchi4', 'dchi5', 'dchi5l', 'dchi6', 'dchi6l', 'dchi7', 'delta_f', 'delta_t', 'distance', 'dquad_mon1', 'dquad_mon2', 'eccentricity', 'eccentricity_order', 'f_final', 'f_final_func', 'f_lower', 'f_ref', 'frame_axis', 'inclination', 'lambda1', 'lambda2', 'lambda_octu1', 'lambda_octu2', 'long_asc_nodes', 'mass1', 'mass2', 'mean_per_ano', 'mode_array', 'modes_choice', 'numrel_data', 'octufmode1', 'octufmode2', 'phase_order', 'quadfmode1', 'quadfmode2', 'side_bands', 'spin1x', 'spin1y', 'spin1z', 'spin2x', 'spin2y', 'spin2z', 'spin_order', 'taper', 'tidal_order'}

The set of names of arguments that may be used in the variable_args or frozen_params.

Type:

set

class pycbc.waveform.generator.BaseFDomainDetFrameGenerator(rFrameGeneratorClass, epoch, detectors=None, variable_args=(), recalib=None, gates=None, **frozen_params)[source]

Bases: object

Base generator for frquency-domain waveforms in a detector frame.

Parameters:
  • rFrameGeneratorClass (class) – The class to use for generating the waveform in the radiation frame, e.g., FDomainCBCGenerator. This should be the class, not an instance of the class (the class will be initialized with the appropriate arguments internally).

  • detectors ({None, list of strings}) – The names of the detectors to use. If provided, all location parameters must be included in either the variable args or the frozen params. If None, the generate function will just return the plus polarization returned by the rFrameGeneratorClass shifted by any desired time shift.

  • epoch ({float, lal.LIGOTimeGPS) – The epoch start time to set the waveform to. A time shift = tc - epoch is applied to waveforms before returning.

  • variable_args ({(), list or tuple}) – A list or tuple of strings giving the names and order of parameters that will be passed to the generate function.

  • **frozen_params – Keyword arguments setting the parameters that will not be changed from call-to-call of the generate function.

detectors

The dictionary of detectors that antenna patterns are calculated for on each call of generate. If no detectors were provided, will be {'RF': None}, where “RF” means “radiation frame”.

Type:

dict

detector_names

The list of detector names. If no detectors were provided, then this will be [‘RF’] for “radiation frame”.

Type:

list

current_params

A dictionary of name, value pairs of the arguments that were last used by the generate function.

Type:

dict

rframe_generator

The instance of the radiation-frame generator that is used for waveform generation. All parameters in current_params except for the location params are passed to this class’s generate function.

Type:

instance of rFrameGeneratorClass

frozen_location_args

Any location parameters that were included in the frozen_params.

Type:

dict

variable_args

The list of names of arguments that are passed to the generate function.

Type:

tuple

property epoch

The GPS start time of the frequency series returned by the generate function. A time shift is applied to the waveform equal to tc-epoch. Update by using set_epoch

abstract generate(**kwargs)[source]

The function that generates the waveforms.

location_args = {}

Should be overriden by children classes with a set of parameters that set the binary’s location.

Type:

Set

abstract select_rframe_generator(approximant)[source]

Method to select waveform generator based on an approximant.

set_epoch(epoch)[source]

Sets the epoch; epoch should be a float or a LIGOTimeGPS.

property static_args

Returns a dictionary of the static arguments.

class pycbc.waveform.generator.BaseGenerator(generator, variable_args=(), record_failures=False, **frozen_params)[source]

Bases: object

A wrapper class to call a waveform generator with a set of frozen parameters and a set of variable parameters. The frozen parameters and values, along with a list of variable parameter names, are set at initialization. This way, repeated calls can be made to the underlying generator by simply passing a list of values for the variable parameters to this class’s generate function.

Parameters:
  • generator (function) – The function that is called for waveform generation.

  • variable_args ({(), list}) – A tuple or list of strings giving the names and order of variable parameters that will be passed to the waveform generator when the generate function is called.

  • record_failures (boolean) – Store output files containing the parameters of failed waveform generation. Default is False.

  • **frozen_params – These keyword arguments are the ones that will be frozen in the waveform generator. For a list of possible parameters, see pycbc.waveform.cbc_parameters.

generator

The function that is called for waveform generation.

Type:

function

variable_args

The list of names of variable arguments.

Type:

tuple

frozen_params

A dictionary of the frozen keyword arguments that are always passed to the waveform generator function.

Type:

dict

current_params

A dictionary of the frozen keyword arguments and variable arguments that were last passed to the waveform generator.

Type:

dict

generate(**kwargs)[source]

Generates a waveform from the keyword args. The current params are updated with the given kwargs, then the generator is called.

property static_args

Returns a dictionary of the static arguments.

class pycbc.waveform.generator.FDomainCBCGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseCBCGenerator

Generates frequency-domain CBC waveforms in the radiation frame.

Uses waveform.get_fd_waveform as a generator function to create frequency- domain CBC waveforms in the radiation frame; i.e., with no detector response function applied. For more details, see BaseGenerator.

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import FDomainCBCGenerator
>>> generator = FDomainCBCGenerator(variable_args=['mass1', 'mass2'], delta_f=1./32, f_lower=30., approximant='TaylorF2')

Create a waveform with the variable arguments (in this case, mass1, mass2):

>>> generator.generate(mass1=1.4, mass2=1.4)
    (<pycbc.types.frequencyseries.FrequencySeries at 0x1110c1450>,
     <pycbc.types.frequencyseries.FrequencySeries at 0x1110c1510>)
class pycbc.waveform.generator.FDomainCBCModesGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseCBCGenerator

Generates frequency-domain CBC waveform modes.

Uses waveform_modes.get_fd_waveform_modes() as a generator function to create frequency-domain CBC waveforms mode-by-mode, without applying spherical harmonics.

For details, on methods and arguments, see BaseGenerator.

class pycbc.waveform.generator.FDomainDetFrameGenerator(rFrameGeneratorClass, epoch, detectors=None, variable_args=(), recalib=None, gates=None, **frozen_params)[source]

Bases: BaseFDomainDetFrameGenerator

Generates frequency-domain waveform in a specific frame.

Generates a waveform using the given radiation frame generator class, and applies the detector response function and appropriate time offset.

Parameters:
  • rFrameGeneratorClass (class) – The class to use for generating the waveform in the radiation frame, e.g., FDomainCBCGenerator. This should be the class, not an instance of the class (the class will be initialized with the appropriate arguments internally).

  • detectors ({None, list of strings}) – The names of the detectors to use. If provided, all location parameters must be included in either the variable args or the frozen params. If None, the generate function will just return the plus polarization returned by the rFrameGeneratorClass shifted by any desired time shift.

  • epoch ({float, lal.LIGOTimeGPS) – The epoch start time to set the waveform to. A time shift = tc - epoch is applied to waveforms before returning.

  • variable_args ({(), list or tuple}) – A list or tuple of strings giving the names and order of parameters that will be passed to the generate function.

  • **frozen_params – Keyword arguments setting the parameters that will not be changed from call-to-call of the generate function.

detectors

The dictionary of detectors that antenna patterns are calculated for on each call of generate. If no detectors were provided, will be {'RF': None}, where “RF” means “radiation frame”.

Type:

dict

detector_names

The list of detector names. If no detectors were provided, then this will be [‘RF’] for “radiation frame”.

Type:

list

epoch

The GPS start time of the frequency series returned by the generate function. A time shift is applied to the waveform equal to tc-epoch. Update by using set_epoch.

Type:

lal.LIGOTimeGPS

current_params

A dictionary of name, value pairs of the arguments that were last used by the generate function.

Type:

dict

rframe_generator

The instance of the radiation-frame generator that is used for waveform generation. All parameters in current_params except for the location params are passed to this class’s generate function.

Type:

instance of rFrameGeneratorClass

frozen_location_args

Any location parameters that were included in the frozen_params.

Type:

dict

variable_args

The list of names of arguments that are passed to the generate function.

Type:

tuple

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import FDomainDetFrameGenerator
>>> generator = FDomainDetFrameGenerator(waveform.FDomainCBCGenerator, 0., variable_args=['mass1', 'mass2', 'spin1z', 'spin2z', 'tc', 'ra', 'dec', 'polarization'], detectors=['H1', 'L1'], delta_f=1./64, f_lower=20., approximant='SEOBNRv2_ROM_DoubleSpin')

Generate a waveform:

>>> generator.generate(mass1=38.6, mass2=29.3, spin1z=0.33, spin2z=-0.94, tc=2.43, ra=1.37, dec=-1.26, polarization=2.76)
{'H1': <pycbc.types.frequencyseries.FrequencySeries at 0x116637350>,
 'L1': <pycbc.types.frequencyseries.FrequencySeries at 0x116637a50>}
generate(**kwargs)[source]

Generates a waveform, applies a time shift and the detector response function from the given kwargs.

location_args = {'dec', 'polarization', 'ra', 'tc'}

set([‘tc’, ‘ra’, ‘dec’, ‘polarization’]): The set of location parameters. These are not passed to the rFrame generator class; instead, they are used to apply the detector response function and/or shift the waveform in time. The parameters are:

  • tc: The GPS time of coalescence (should be geocentric time).

  • ra: Right ascension.

  • dec: declination

  • polarization: polarization.

All of these must be provided in either the variable args or the frozen params if detectors is not None. If detectors is None, tc may optionally be provided.

static select_rframe_generator(approximant)[source]

Returns a radiation frame generator class based on the approximant string.

class pycbc.waveform.generator.FDomainDetFrameModesGenerator(rFrameGeneratorClass, epoch, detectors=None, variable_args=(), recalib=None, gates=None, **frozen_params)[source]

Bases: BaseFDomainDetFrameGenerator

Generates frequency-domain waveform modes in a specific frame.

Generates both polarizations of every waveform mode using the given radiation frame generator class, and applies the time shift. Detector response functions are not applied.

Parameters:
  • rFrameGeneratorClass (class) – The class to use for generating the waveform modes in the radiation frame, e.g., FDomainCBCModesGenerator. This should be the class, not an instance of the class (the class will be initialized with the appropriate arguments internally). The class should have a generate function that returns a dictionary of waveforms keyed by the modes.

  • detectors ({None, list of strings}) – The names of the detectors to use. If provided, all location parameters must be included in either the variable args or the frozen params. If None, the generate function will just return the plus polarization returned by the rFrameGeneratorClass shifted by any desired time shift.

  • epoch ({float, lal.LIGOTimeGPS) – The epoch start time to set the waveform to. A time shift = tc - epoch is applied to waveforms before returning.

  • variable_args ({(), list or tuple}) – A list or tuple of strings giving the names and order of parameters that will be passed to the generate function.

  • **frozen_params – Keyword arguments setting the parameters that will not be changed from call-to-call of the generate function.

detectors

The dictionary of detectors that antenna patterns are calculated for on each call of generate. If no detectors were provided, will be {'RF': None}, where “RF” means “radiation frame”.

Type:

dict

detector_names

The list of detector names. If no detectors were provided, then this will be [‘RF’] for “radiation frame”.

Type:

list

epoch

The GPS start time of the frequency series returned by the generate function. A time shift is applied to the waveform equal to tc-epoch. Update by using set_epoch.

Type:

lal.LIGOTimeGPS

current_params

A dictionary of name, value pairs of the arguments that were last used by the generate function.

Type:

dict

rframe_generator

The instance of the radiation-frame generator that is used for waveform generation. All parameters in current_params except for the location params are passed to this class’s generate function.

Type:

instance of rFrameGeneratorClass

frozen_location_args

Any location parameters that were included in the frozen_params.

Type:

dict

variable_args

The list of names of arguments that are passed to the generate function.

Type:

tuple

generate(**kwargs)[source]

Generates and returns a waveform decompsed into separate modes.

Returns:

Dictionary of detector names -> modes -> (ulm, vlm), where ulm, vlm are the frequency-domain representations of the real and imaginary parts, respectively, of the complex time series representation of the hlm.

Return type:

dict

location_args = {'dec', 'ra', 'tc'}

set([‘tc’, ‘ra’, ‘dec’]): The set of location parameters. These are not passed to the rFrame generator class; instead, they are used to apply the detector response function and/or shift the waveform in time. The parameters are:

  • tc: The GPS time of coalescence (should be geocentric time).

  • ra: Right ascension.

  • dec: declination

All of these must be provided in either the variable args or the frozen params if detectors is not None. If detectors is None, tc may optionally be provided.

static select_rframe_generator(approximant)[source]

Returns a radiation frame generator class based on the approximant string.

class pycbc.waveform.generator.FDomainDetFrameTwoPolGenerator(rFrameGeneratorClass, epoch, detectors=None, variable_args=(), recalib=None, gates=None, **frozen_params)[source]

Bases: BaseFDomainDetFrameGenerator

Generates frequency-domain waveform in a specific frame.

Generates both polarizations of a waveform using the given radiation frame generator class, and applies the time shift. Detector response functions are not applied.

Parameters:
  • rFrameGeneratorClass (class) – The class to use for generating the waveform in the radiation frame, e.g., FDomainCBCGenerator. This should be the class, not an instance of the class (the class will be initialized with the appropriate arguments internally).

  • detectors ({None, list of strings}) – The names of the detectors to use. If provided, all location parameters must be included in either the variable args or the frozen params. If None, the generate function will just return the plus polarization returned by the rFrameGeneratorClass shifted by any desired time shift.

  • epoch ({float, lal.LIGOTimeGPS) – The epoch start time to set the waveform to. A time shift = tc - epoch is applied to waveforms before returning.

  • variable_args ({(), list or tuple}) – A list or tuple of strings giving the names and order of parameters that will be passed to the generate function.

  • **frozen_params – Keyword arguments setting the parameters that will not be changed from call-to-call of the generate function.

detectors

The dictionary of detectors that antenna patterns are calculated for on each call of generate. If no detectors were provided, will be {'RF': None}, where “RF” means “radiation frame”.

Type:

dict

detector_names

The list of detector names. If no detectors were provided, then this will be [‘RF’] for “radiation frame”.

Type:

list

epoch

The GPS start time of the frequency series returned by the generate function. A time shift is applied to the waveform equal to tc-epoch. Update by using set_epoch.

Type:

lal.LIGOTimeGPS

current_params

A dictionary of name, value pairs of the arguments that were last used by the generate function.

Type:

dict

rframe_generator

The instance of the radiation-frame generator that is used for waveform generation. All parameters in current_params except for the location params are passed to this class’s generate function.

Type:

instance of rFrameGeneratorClass

frozen_location_args

Any location parameters that were included in the frozen_params.

Type:

dict

variable_args

The list of names of arguments that are passed to the generate function.

Type:

tuple

generate(**kwargs)[source]

Generates a waveform polarizations and applies a time shift.

Returns:

Dictionary of detector names -> (hp, hc), where hp, hc are the plus and cross polarization, respectively.

Return type:

dict

location_args = {'dec', 'ra', 'tc'}

set([‘tc’, ‘ra’, ‘dec’]): The set of location parameters. These are not passed to the rFrame generator class; instead, they are used to apply the detector response function and/or shift the waveform in time. The parameters are:

  • tc: The GPS time of coalescence (should be geocentric time).

  • ra: Right ascension.

  • dec: declination

All of these must be provided in either the variable args or the frozen params if detectors is not None. If detectors is None, tc may optionally be provided.

static select_rframe_generator(approximant)[source]

Returns a radiation frame generator class based on the approximant string.

class pycbc.waveform.generator.FDomainDetFrameTwoPolNoRespGenerator(rFrameGeneratorClass, epoch, detectors=None, variable_args=(), recalib=None, gates=None, **frozen_params)[source]

Bases: BaseFDomainDetFrameGenerator

Generates frequency-domain waveform in a specific frame.

Generates both polarizations of a waveform using the given radiation frame generator class, and applies the time shift. Detector response functions are not applied.

Parameters:
  • rFrameGeneratorClass (class) – The class to use for generating the waveform in the radiation frame, e.g., FDomainCBCGenerator. This should be the class, not an instance of the class (the class will be initialized with the appropriate arguments internally).

  • detectors ({None, list of strings}) – The names of the detectors to use. If provided, all location parameters must be included in either the variable args or the frozen params. If None, the generate function will just return the plus polarization returned by the rFrameGeneratorClass shifted by any desired time shift.

  • epoch ({float, lal.LIGOTimeGPS) – The epoch start time to set the waveform to. A time shift = tc - epoch is applied to waveforms before returning.

  • variable_args ({(), list or tuple}) – A list or tuple of strings giving the names and order of parameters that will be passed to the generate function.

  • **frozen_params – Keyword arguments setting the parameters that will not be changed from call-to-call of the generate function.

detectors

The dictionary of detectors that antenna patterns are calculated for on each call of generate. If no detectors were provided, will be {'RF': None}, where “RF” means “radiation frame”.

Type:

dict

detector_names

The list of detector names. If no detectors were provided, then this will be [‘RF’] for “radiation frame”.

Type:

list

epoch

The GPS start time of the frequency series returned by the generate function. A time shift is applied to the waveform equal to tc-epoch. Update by using set_epoch.

Type:

lal.LIGOTimeGPS

current_params

A dictionary of name, value pairs of the arguments that were last used by the generate function.

Type:

dict

rframe_generator

The instance of the radiation-frame generator that is used for waveform generation. All parameters in current_params except for the location params are passed to this class’s generate function.

Type:

instance of rFrameGeneratorClass

frozen_location_args

Any location parameters that were included in the frozen_params.

Type:

dict

variable_args

The list of names of arguments that are passed to the generate function.

Type:

tuple

generate(**kwargs)[source]

Generates a waveform polarizations

Returns:

Dictionary of detector names -> (hp, hc), where hp, hc are the plus and cross polarization, respectively.

Return type:

dict

static select_rframe_generator(approximant)[source]

Returns a radiation frame generator class based on the approximant string.

class pycbc.waveform.generator.FDomainFreqTauRingdownGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseGenerator

Uses ringdown.get_fd_from_freqtau as a generator function to create frequency-domain ringdown waveforms with higher modes in the radiation frame; i.e., with no detector response function applied. For more details, see BaseGenerator.

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import FDomainFreqTauRingdownGenerator
>>> generator = FDomainFreqTauRingdownGenerator(variable_args=['f_220',
                'tau_220','f_210','tau_210','amp220','amp210','phi220','phi210'],
                lmns=['221','211'], delta_f=1./32, f_lower=30., f_final=500)

Create a ringdown with the variable arguments:

>>> generator.generate(f_220=317., tau_220=0.003, f_210=274., tau_210=0.003,
                       amp220=1e-21, amp210=1./10, phi220=0., phi210=0.)
    (<pycbc.types.frequencyseries.FrequencySeries at 0x51614d0>,
     <pycbc.types.frequencyseries.FrequencySeries at 0x5161550>)
class pycbc.waveform.generator.FDomainMassSpinRingdownGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseGenerator

Uses ringdown.get_fd_from_final_mass_spin as a generator function to create frequency-domain ringdown waveforms with higher modes in the radiation frame; i.e., with no detector response function applied. For more details, see BaseGenerator.

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import FDomainMassSpinRingdownGenerator
>>> generator = FDomainMassSpinRingdownGenerator(variable_args=['final_mass',
                'final_spin','amp220','amp210','phi220','phi210'], lmns=['221','211'],
                delta_f=1./32, f_lower=30., f_final=500)

Create a ringdown with the variable arguments:

>>> generator.generate(final_mass=65., final_spin=0.7,
                       amp220=1e-21, amp210=1./10, phi220=0., phi210=0.)
    (<pycbc.types.frequencyseries.FrequencySeries at 0x51614d0>,
     <pycbc.types.frequencyseries.FrequencySeries at 0x5161550>)
class pycbc.waveform.generator.TDomainCBCGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseCBCGenerator

Create time domain CBC waveforms in the radiation frame.

Uses waveform.get_td_waveform as a generator function to create time- domain CBC waveforms in the radiation frame; i.e., with no detector response function applied. For more details, see BaseGenerator.

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import TDomainCBCGenerator
>>> generator = TDomainCBCGenerator(variable_args=['mass1', 'mass2'], delta_t=1./4096, f_lower=30., approximant='TaylorT4')

Create a waveform with the variable arguments (in this case, mass1, mass2):

>>> generator.generate(mass1=2., mass2=1.3)
    (<pycbc.types.timeseries.TimeSeries at 0x10e546710>,
     <pycbc.types.timeseries.TimeSeries at 0x115f37690>)
class pycbc.waveform.generator.TDomainCBCModesGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseCBCGenerator

Generates time domain CBC waveform modes.

Uses waveform_modes.get_td_waveform_modes() as a generator function to create time-domain CBC waveforms mode-by-mode, without applying spherical harmonics. The generate function returns a dictionary of modes -> (real, imag) part of the complex time series.

For details, on methods and arguments, see BaseGenerator.

class pycbc.waveform.generator.TDomainFreqTauRingdownGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseGenerator

Uses ringdown.get_td_from_freqtau as a generator function to create time-domain ringdown waveforms with higher modes in the radiation frame; i.e., with no detector response function applied. For more details, see BaseGenerator.

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import FDomainFreqTauRingdownGenerator
>>> generator = TDomainFreqTauRingdownGenerator(variable_args=['f_220',
                'tau_220','f_210','tau_210','amp220','amp210','phi220','phi210'],
                lmns=['221','211'], delta_t=1./2048)

Create a ringdown with the variable arguments:

>>> generator.generate(f_220=317., tau_220=0.003, f_210=274., tau_210=0.003,
                       amp220=1e-21, amp210=1./10, phi220=0., phi210=0.)
    (<pycbc.types.frequencyseries.FrequencySeries at 0x51614d0>,
     <pycbc.types.frequencyseries.FrequencySeries at 0x5161550>)
class pycbc.waveform.generator.TDomainMassSpinRingdownGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseGenerator

Uses ringdown.get_td_from_final_mass_spin as a generator function to create time-domain ringdown waveforms with higher modes in the radiation frame; i.e., with no detector response function applied. For more details, see BaseGenerator.

Examples

Initialize a generator:

>>> from pycbc.waveform.generator import TDomainMassSpinRingdownGenerator
>>> generator = TDomainMassSpinRingdownGenerator(variable_args=['final_mass',
                'final_spin','amp220','amp210','phi220','phi210'], lmns=['221','211'],
                delta_t=1./2048)

Create a ringdown with the variable arguments:

>>> generator.generate(final_mass=65., final_spin=0.7,
                       amp220=1e-21, amp210=1./10, phi220=0., phi210=0.)
    (<pycbc.types.frequencyseries.FrequencySeries at 0x51614d0>,
     <pycbc.types.frequencyseries.FrequencySeries at 0x5161550>)
class pycbc.waveform.generator.TDomainSupernovaeGenerator(variable_args=(), **frozen_params)[source]

Bases: BaseGenerator

Uses supernovae.py to create time domain core-collapse supernovae waveforms using a set of Principal Components provided in a .hdf file.

pycbc.waveform.generator.select_waveform_generator(approximant)[source]

Returns the single-IFO generator for the approximant.

Parameters:

approximant (str) – Name of waveform approximant. Valid names can be found using pycbc.waveform methods.

Returns:

generator – A waveform generator object.

Return type:

(PyCBC generator instance)

Examples

Get a list of available approximants: >>> from pycbc import waveform >>> waveform.fd_approximants() >>> waveform.td_approximants() >>> from pycbc.waveform import ringdown >>> ringdown.ringdown_fd_approximants.keys()

Get generator object: >>> from pycbc.waveform.generator import select_waveform_generator >>> select_waveform_generator(waveform.fd_approximants()[0])

pycbc.waveform.generator.select_waveform_modes_generator(approximant)[source]

Returns the single-IFO modes generator for the approximant.

Parameters:

approximant (str) – Name of waveform approximant. Valid names can be found using pycbc.waveform methods.

Returns:

generator – A waveform generator object.

Return type:

(PyCBC generator instance)

pycbc.waveform.multiband module

Tools and functions to calculate interpolate waveforms using multi-banding

pycbc.waveform.multiband.multiband_fd_waveform(bands=None, lengths=None, overlap=0, **p)[source]

Generate a fourier domain waveform using multibanding

Speed up generation of a fouerier domain waveform using multibanding. This allows for multi-rate sampling of the frequeny space. Each band is smoothed and stitched together to produce the final waveform. The base approximant must support ‘f_ref’ and ‘f_final’. The other parameters must be chosen carefully by the user.

Parameters:
  • bands (list or str) – The frequencies to split the waveform by. These should be chosen so that the corresponding length include all the waveform’s frequencies within this band.

  • lengths (list or str) – The corresponding length for each frequency band. This sets the resolution of the subband and should be chosen carefully so that it is sufficiently long to include all of the bands frequency content.

  • overlap (float) – The frequency width to apply tapering between bands.

  • params (dict) – The remaining keyworkd arguments passed to the base approximant waveform generation.

Returns:

  • hp (pycbc.types.FrequencySeries) – Plus polarization

  • hc (pycbc.type.FrequencySeries) – Cross polarization

pycbc.waveform.nltides module

Utilities for introducing nonlinear tidal effects into waveform approximants

pycbc.waveform.nltides.nltides_fourier_phase_difference(f, delta_f, f0, amplitude, n, m1, m2)[source]

Calculate the change to the Fourier phase change due to non-linear tides. Note that the Fourier phase Psi(f) is not the same as the gravitational-wave phase phi(f) and is computed by Delta Psi(f) = 2 pi f Delta t(f) - Delta phi(f)

Parameters:
  • f (numpy.array) – Array of frequency values to calculate the fourier phase difference

  • delta_f (float) – Frequency resolution of f array

  • f0 (float) – Frequency that NL effects switch on

  • amplitude (float) – Amplitude of effect

  • n (float) – Growth dependence of effect

  • m1 (float) – Mass of component 1

  • m2 (float) – Mass of component 2

Returns:

delta_psi – Fourier phase as a function of frequency

Return type:

numpy.array

pycbc.waveform.nltides.nonlinear_tidal_spa(**kwds)[source]

Generates a frequency-domain waveform that implements the TaylorF2+NL tide model described in https://arxiv.org/abs/1808.07013

pycbc.waveform.parameters module

Classes to define common parameters used for waveform generation.

class pycbc.waveform.parameters.Parameter(name, dtype=None, default=None, label=None, description='No description.')[source]

Bases: str

A class that stores information about a parameter. This is done by sub-classing string, adding additional attributes.

docstr(prefix='', include_label=True)[source]

Returns a string summarizing the parameter. Format is: <prefix>``name`` : {default, dtype} <prefix> description Label: label.

class pycbc.waveform.parameters.ParameterList(initlist=None)[source]

Bases: UserList

A list of parameters. Each element in the list is expected to be a Parameter instance.

property asdict

Returns a dictionary of the parameters keyed by the parameters.

property aslist

Cast to basic list.

default_dict()[source]

Returns a dictionary of the name and default value of each parameter.

defaults()[source]

Returns a list of the name and default value of each parameter, as tuples.

property description_dict

Return a dictionary of the name and description of each parameter.

property descriptions

Returns a list of the name and description of each parameter, as tuples.

docstr(prefix='', include_label=True)[source]

Returns the docstr of each parameter joined together.

property dtype_dict

Returns a dictionary of the name and dtype of each parameter.

property dtypes

Returns a list of the name and dtype of each parameter, as tuples.

property label_dict

Return a dictionary of the name and label of each parameter.

property labels

Returns a list of each parameter and its label, as tuples.

property names

Returns a list of the names of each parameter.

property nodefaults

Returns a ParameterList of the parameters that have None for defaults.

pycbc.waveform.plugin module

Utilities for handling waveform plugins

pycbc.waveform.plugin.add_custom_waveform(approximant, function, domain, sequence=False, has_det_response=False, force=False)[source]

Make custom waveform available to pycbc

Parameters:
  • approximant (str) – The name of the waveform

  • function (function) – The function to generate the waveform

  • domain (str) – Either ‘frequency’ or ‘time’ to indicate the domain of the waveform.

  • sequence (bool, False) – Function evaluates waveform at only chosen points (instead of a equal-spaced grid).

  • has_det_response (bool, False) – Check if waveform generator has built-in detector response.

pycbc.waveform.plugin.add_length_estimator(approximant, function)[source]

Add length estimator for an approximant

Parameters:
  • approximant (str) – Name of approximant

  • function (function) – A function which takes kwargs and returns the waveform length

pycbc.waveform.plugin.retrieve_waveform_plugins()[source]

Process external waveform plugins

pycbc.waveform.premerger module

Waveform approximants for the pre-merger detection of gravitational waves

pycbc.waveform.premerger.premerger_taylorf2(**p)[source]

Generate time-shifted TaylorF2

pycbc.waveform.pycbc_phenomC_tmplt module

pycbc.waveform.ringdown module

Generate ringdown templates in the time and frequency domain.

pycbc.waveform.ringdown.Kerr_factor(final_mass, distance)[source]

Return the factor final_mass/distance (in dimensionless units) for Kerr ringdowns

pycbc.waveform.ringdown.fd_damped_sinusoid(f_0, tau, amp, phi, freqs, t_0=0.0, l=2, m=2, n=0, inclination=0.0, azimuthal=0.0, harmonics='spherical', final_spin=None, pol=None, polnm=None)[source]

Return the frequency domain version of a damped sinusoid.

This is the frequency domain version td_damped_sinusoid() without a taper if an infinite sample rate were used to resolve the step function that turns on the damped sinusoid. See td_damped_sinusoid() for details.

Note

This function currently does not support using a different amplitude and phase for the -m modes (equivalent to setting dphi = dbeta = 0 in td_damped_sinusoid().

Parameters:
  • f_0 (float) – The central frequency of the damped sinusoid, in Hz.

  • tau (float) – The damping time, in seconds.

  • amp (float) – The intrinsic amplitude of the QNM (\(A^0_{lmn}\) in the above).

  • phi (float) – The reference phase of the QNM (\(\phi_{lmn}\)) in the above.

  • freqs (array) – Array of frequencies to evaluate the damped sinusoid over.

  • t_0 (float, optional) – The start time of ringdown. Default (0.) corresponds to the ringdown starting at the beginning of the equivalent segment in the time domain. A non-zero value for t_0 will shift the ringdown in time to the corresponding number of seconds from the start of the segment.

  • l (int, optional) – The l index; default is 2.

  • m (int, optional) – The m index; default is 2.

  • n (int, optional) – The overtone index; default is 0.

  • inclination (float, optional) – The inclination angle (\(\theta\) in the above). Ignored if harmonics='arbitrary'. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle (\(\varphi\) in the above). Ignored if harmonics='arbitrary'. Default is 0.

  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – Which harmonics to use. See spher_harms() for details. Default is spherical.

  • final_spin (float, optional) – The dimensionless spin of the black hole. Only needed if harmonics='spheroidal'.

  • pol (float, optional) – Angle to use for +m arbitrary harmonics. Only needed if harmonics='arbitrary'. See spher_harms() for details.

  • polnm (float, optional) – Angle to use for -m arbitrary harmonics. Only needed if harmonics='arbitrary'. See spher_harms() for details.

Returns:

  • hptilde (numpy.ndarray) – The plus polarization.

  • hctilde (numpy.ndarray) – The cross polarization.

pycbc.waveform.ringdown.fd_output_vector(freqs, damping_times, delta_f=None, f_final=None)[source]

Return an empty FrequencySeries with the appropriate size to fit all the quasi-normal modes present in freqs, damping_times

pycbc.waveform.ringdown.format_lmns(lmns)[source]

Checks if the format of the parameter lmns is correct, returning the appropriate format if not, and raise an error if nmodes=0.

The required format for the ringdown approximants is a list of lmn modes as a single whitespace-separated string, with n the number of overtones desired. Alternatively, a list object may be used, containing individual strings as elements for the lmn modes. For instance, lmns = ‘223 331’ are the modes 220, 221, 222, and 330. Giving lmns = [‘223’, ‘331’] is equivalent. The ConfigParser of a workflow might convert that to a single string (case 1 below) or a list with a single string (case 2), and this function will return the appropriate list of strings. If a different format is given, raise an error.

pycbc.waveform.ringdown.get_fd_from_final_mass_spin(template=None, **kwargs)[source]

Return frequency domain ringdown with all the modes specified.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • final_mass (float) – Mass of the final black hole in solar masses.

  • final_spin (float) – Dimensionless spin of the final black hole.

  • distance ({None, float}, optional) – Luminosity distance of the system. If specified, the returned ringdown will include the Kerr factor (final_mass/distance).

  • lmns (list) – Desired lmn modes as strings. All modes up to l = m = 7 are available. The n specifies the number of overtones desired for the corresponding lm pair, not the overtone number; maximum n=8. Example: lmns = [‘223’,’331’] are the modes 220, 221, 222, and 330

  • ref_amp (str, optional) – Which mode to use as the reference for computing amplitudes. Must be ‘amp220’ if distance is given. Default is ‘amp220’. The amplitude of the reference mode should be specified directly, while all other amplitudes are specified as ratios with respect to that mode. For example, if ref_amp = 'amp220', lmns = ['221', '331'], and no distance is provided, then amp220 = 1e-22, amp330 = 0.1 would result in the 220 mode having a strain amplitude of 1e-22 and the 330 mode having a strain amplitude of 1e-23. If distance is given, the amplitude of the reference mode will have a completely different order of magnitude. See table II in https://arxiv.org/abs/1107.0854 for an estimate. An amplitude for the reference mode must always be provided, even if that mode is not being generated. For example, if ref_amp = 'amp220' and lmns = ['331'] then both a amp220 and amp330 must be provided even though only the 330 mode will be created.

  • amplmn (float) – The amplitude of each mode, required for all modes specifed plus the reference mode. As described above, amplitudes should be specified relative to the reference mode.

  • philmn (float) – Phase of the lmn overtone, as many as the number of modes.

  • inclination (float) – Inclination of the system in radians. Ignored if harmonics='arbitrary'. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle in radians. Ignored if harmonics='arbitrary'. Usually this is not necessary to specify since it is degenerate with initial phase philmn; i.e., this is only useful if you have an expectation for what the phase of each mode is. Default is 0.

  • pollmn (float, optional) – Angle to use for +m arbitrary harmonics of the lmn mode in radians (example: pol220 = 0.1. Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • polnmlmn (float, optional) – Angle to use for -m arbitrary harmonics of the lmn mode in radians (example: polnm220 = 0.1). Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – Which harmonics to use. See spher_harms() for details. Default is spherical.

  • delta_flmn ({None, float}, optional) – GR deviation for the frequency of the lmn mode. If given, the lmn frequency will be converted to new_flmn = flmn + delta_flmn * flmn, with flmn the GR predicted value for the corresponding mass and spin.

  • delta_taulmn ({None, float}, optional) – GR deviation for the damping time of the lmn mode. If given, the lmn tau will be converted to new_taulmn = taulmn + delta_taulmn * taulmn, with taulmn the GR predicted value for the corresponding mass and spin.

  • delta_f ({None, float}, optional) – The frequency step used to generate the ringdown (not to be confused with the delta_flmn parameter that simulates GR violations). If None, it will be set to the inverse of the time at which the amplitude is 1/1000 of the peak amplitude (the minimum of all modes).

  • f_lower ({None, float}, optional) – The starting frequency of the output frequency series. If None, it will be set to delta_f.

  • f_final ({None, float}, optional) – The ending frequency of the output frequency series. If None, it will be set to the frequency at which the amplitude is 1/1000 of the peak amplitude (the maximum of all modes).

Returns:

  • hplustilde (FrequencySeries) – The plus phase of a ringdown with the lm modes specified and n overtones in frequency domain.

  • hcrosstilde (FrequencySeries) – The cross phase of a ringdown with the lm modes specified and n overtones in frequency domain.

pycbc.waveform.ringdown.get_fd_from_freqtau(template=None, **kwargs)[source]

Return frequency domain ringdown with all the modes specified.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • lmns (list) – Desired lmn modes as strings. All modes up to l = m = 7 are available. The n specifies the number of overtones desired for the corresponding lm pair, not the overtone number; maximum n=8. Example: lmns = [‘223’,’331’] are the modes 220, 221, 222, and 330

  • f_lmn (float) – Central frequency of the lmn overtone, as many as number of modes.

  • tau_lmn (float) – Damping time of the lmn overtone, as many as number of modes.

  • ref_amp (str, optional) – Which mode to use as the reference for computing amplitudes. Must be ‘amp220’ if distance is given. Default is ‘amp220’. The amplitude of the reference mode should be specified directly, while all other amplitudes are specified as ratios with respect to that mode. For example, if ref_amp = 'amp220', lmns = ['221', '331'], and no distance is provided, then amp220 = 1e-22, amp330 = 0.1 would result in the 220 mode having a strain amplitude of 1e-22 and the 330 mode having a strain amplitude of 1e-23. If distance is given, the amplitude of the reference mode will have a completely different order of magnitude. See table II in https://arxiv.org/abs/1107.0854 for an estimate. An amplitude for the reference mode must always be provided, even if that mode is not being generated. For example, if ref_amp = 'amp220' and lmns = ['331'] then both a amp220 and amp330 must be provided even though only the 330 mode will be created.

  • amplmn (float) – The amplitude of each mode, required for all modes specifed plus the reference mode. As described above, amplitudes should be specified relative to the reference mode.

  • philmn (float) – Phase of the lmn overtone, as many as the number of modes.

  • inclination (float) – Inclination of the system in radians. Ignored if harmonics='arbitrary'. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle in radians. Ignored if harmonics='arbitrary'. Usually this is not necessary to specify since it is degenerate with initial phase philmn; i.e., this is only useful if you have an expectation for what the phase of each mode is. Default is 0.

  • dphi[lmn] (float, optional) – The difference in phase between the +m and -m mode. See the documentation for dphi in td_damped_sinusoid() for details. You may specify a dphi{lmn} (ex. dphi220) separately for each mode, and/or a single dphi (without any lmn) for all modes that do not have dphi specified. Default is to use 0 for all modes.

  • dbeta[lmn] (float, optional) – The angular difference in the amplitudes of the +m and -m mode. See the documentation for dbeta in td_damped_sinusoid() for details. You may specify a dbeta{lmn} (ex. dbeta220) separately for each mode, and/or a single dbeta (without any lmn) for all modes that do not have dbeta specified. Default is to use 0 for all modes.

  • pollmn (float, optional) – Angle to use for +m arbitrary harmonics of the lmn mode in radians (example: pol220 = 0.1. Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • polnmlmn (float, optional) – Angle to use for -m arbitrary harmonics of the lmn mode in radians (example: polnm220 = 0.1). Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – Which harmonics to use. See spher_harms() for details. Default is spherical.

  • final_spin (float, optional) – Dimensionless spin of the final black hole. This is required if harmonics='spheroidal'. Ignored otherwise.

  • delta_f ({None, float}, optional) – The frequency step used to generate the ringdown. If None, it will be set to the inverse of the time at which the amplitude is 1/1000 of the peak amplitude (the minimum of all modes).

  • f_lower ({None, float}, optional) – The starting frequency of the output frequency series. If None, it will be set to delta_f.

  • f_final ({None, float}, optional) – The ending frequency of the output frequency series. If None, it will be set to the frequency at which the amplitude is 1/1000 of the peak amplitude (the maximum of all modes).

Returns:

  • hplustilde (FrequencySeries) – The plus phase of a ringdown with the lm modes specified and n overtones in frequency domain.

  • hcrosstilde (FrequencySeries) – The cross phase of a ringdown with the lm modes specified and n overtones in frequency domain.

pycbc.waveform.ringdown.get_td_from_final_mass_spin(template=None, **kwargs)[source]

Return time domain ringdown with all the modes specified.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • final_mass (float) – Mass of the final black hole in solar masses.

  • final_spin (float) – Dimensionless spin of the final black hole.

  • distance ({None, float}, optional) – Luminosity distance of the system. If specified, the returned ringdown will include the Kerr factor (final_mass/distance).

  • lmns (list) – Desired lmn modes as strings. All modes up to l = m = 7 are available. The n specifies the number of overtones desired for the corresponding lm pair, not the overtone number; maximum n=8. Example: lmns = [‘223’,’331’] are the modes 220, 221, 222, and 330

  • ref_amp (str, optional) – Which mode to use as the reference for computing amplitudes. Must be ‘amp220’ if distance is given. Default is ‘amp220’. The amplitude of the reference mode should be specified directly, while all other amplitudes are specified as ratios with respect to that mode. For example, if ref_amp = 'amp220', lmns = ['221', '331'], and no distance is provided, then amp220 = 1e-22, amp330 = 0.1 would result in the 220 mode having a strain amplitude of 1e-22 and the 330 mode having a strain amplitude of 1e-23. If distance is given, the amplitude of the reference mode will have a completely different order of magnitude. See table II in https://arxiv.org/abs/1107.0854 for an estimate. An amplitude for the reference mode must always be provided, even if that mode is not being generated. For example, if ref_amp = 'amp220' and lmns = ['331'] then both a amp220 and amp330 must be provided even though only the 330 mode will be created.

  • amplmn (float) – The amplitude of each mode, required for all modes specifed plus the reference mode. As described above, amplitudes should be specified relative to the reference mode.

  • philmn (float) – Phase of the lmn overtone, as many as the number of modes.

  • inclination (float) – Inclination of the system in radians. Ignored if harmonics='arbitrary'. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle in radians. Ignored if harmonics='arbitrary'. Usually this is not necessary to specify since it is degenerate with initial phase philmn; i.e., this is only useful if you have an expectation for what the phase of each mode is. Default is 0.

  • dphi[lmn] (float, optional) – The difference in phase between the +m and -m mode. See the documentation for dphi in td_damped_sinusoid() for details. You may specify a dphi{lmn} (ex. dphi220) separately for each mode, and/or a single dphi (without any lmn) for all modes that do not have dphi specified. Default is to use 0 for all modes.

  • dbeta[lmn] (float, optional) – The angular difference in the amplitudes of the +m and -m mode. See the documentation for dbeta in td_damped_sinusoid() for details. You may specify a dbeta{lmn} (ex. dbeta220) separately for each mode, and/or a single dbeta (without any lmn) for all modes that do not have dbeta specified. Default is to use 0 for all modes.

  • pollmn (float, optional) – Angle to use for +m arbitrary harmonics of the lmn mode in radians (example: pol220 = 0.1. Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • polnmlmn (float, optional) – Angle to use for -m arbitrary harmonics of the lmn mode in radians (example: polnm220 = 0.1). Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – Which harmonics to use. See spher_harms() for details. Default is spherical.

  • delta_flmn ({None, float}, optional) – GR deviation for the frequency of the lmn mode. If given, the lmn frequency will be converted to new_flmn = flmn + delta_flmn * flmn, with flmn the GR predicted value for the corresponding mass and spin.

  • delta_taulmn ({None, float}, optional) – GR deviation for the damping time of the lmn mode. If given, the lmn tau will be converted to new_taulmn = taulmn + delta_taulmn * taulmn, with taulmn the GR predicted value for the corresponding mass and spin.

  • delta_t ({None, float}, optional) – The time step used to generate the ringdown. If None, it will be set to the inverse of the frequency at which the amplitude is 1/1000 of the peak amplitude (the minimum of all modes).

  • t_final ({None, float}, optional) – The ending time of the output frequency series. If None, it will be set to the time at which the amplitude is 1/1000 of the peak amplitude (the maximum of all modes).

  • taper (bool, optional) – Add a rapid ringup with timescale tau/10 at the beginning of the waveform to avoid the abrupt turn on of the ringdown. Each mode and overtone will have a different taper depending on its tau, the final taper being the superposition of all the tapers. Default is False.

Returns:

  • hplus (TimeSeries) – The plus phase of a ringdown with the lm modes specified and n overtones in time domain.

  • hcross (TimeSeries) – The cross phase of a ringdown with the lm modes specified and n overtones in time domain.

pycbc.waveform.ringdown.get_td_from_freqtau(template=None, **kwargs)[source]

Return time domain ringdown with all the modes specified.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • lmns (list) – Desired lmn modes as strings. All modes up to l = m = 7 are available. The n specifies the number of overtones desired for the corresponding lm pair, not the overtone number; maximum n=8. Example: lmns = [‘223’,’331’] are the modes 220, 221, 222, and 330

  • f_lmn (float) – Central frequency of the lmn overtone, as many as number of modes.

  • tau_lmn (float) – Damping time of the lmn overtone, as many as number of modes.

  • ref_amp (str, optional) – Which mode to use as the reference for computing amplitudes. Must be ‘amp220’ if distance is given. Default is ‘amp220’. The amplitude of the reference mode should be specified directly, while all other amplitudes are specified as ratios with respect to that mode. For example, if ref_amp = 'amp220', lmns = ['221', '331'], and no distance is provided, then amp220 = 1e-22, amp330 = 0.1 would result in the 220 mode having a strain amplitude of 1e-22 and the 330 mode having a strain amplitude of 1e-23. If distance is given, the amplitude of the reference mode will have a completely different order of magnitude. See table II in https://arxiv.org/abs/1107.0854 for an estimate. An amplitude for the reference mode must always be provided, even if that mode is not being generated. For example, if ref_amp = 'amp220' and lmns = ['331'] then both a amp220 and amp330 must be provided even though only the 330 mode will be created.

  • amplmn (float) – The amplitude of each mode, required for all modes specifed plus the reference mode. As described above, amplitudes should be specified relative to the reference mode.

  • philmn (float) – Phase of the lmn overtone, as many as the number of modes.

  • inclination (float) – Inclination of the system in radians. Ignored if harmonics='arbitrary'. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle in radians. Ignored if harmonics='arbitrary'. Usually this is not necessary to specify since it is degenerate with initial phase philmn; i.e., this is only useful if you have an expectation for what the phase of each mode is. Default is 0.

  • dphi[lmn] (float, optional) – The difference in phase between the +m and -m mode. See the documentation for dphi in td_damped_sinusoid() for details. You may specify a dphi{lmn} (ex. dphi220) separately for each mode, and/or a single dphi (without any lmn) for all modes that do not have dphi specified. Default is to use 0 for all modes.

  • dbeta[lmn] (float, optional) – The angular difference in the amplitudes of the +m and -m mode. See the documentation for dbeta in td_damped_sinusoid() for details. You may specify a dbeta{lmn} (ex. dbeta220) separately for each mode, and/or a single dbeta (without any lmn) for all modes that do not have dbeta specified. Default is to use 0 for all modes.

  • pollmn (float, optional) – Angle to use for +m arbitrary harmonics of the lmn mode in radians (example: pol220 = 0.1. Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • polnmlmn (float, optional) – Angle to use for -m arbitrary harmonics of the lmn mode in radians (example: polnm220 = 0.1). Only needed if harmonics='arbitrary', ignored otherwise. See spher_harms() for details.

  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – Which harmonics to use. See spher_harms() for details. Default is spherical.

  • final_spin (float, optional) – Dimensionless spin of the final black hole. This is required if harmonics='spheroidal'. Ignored otherwise.

  • delta_t ({None, float}, optional) – The time step used to generate the ringdown. If None, it will be set to the inverse of the frequency at which the amplitude is 1/1000 of the peak amplitude (the minimum of all modes).

  • t_final ({None, float}, optional) – The ending time of the output frequency series. If None, it will be set to the time at which the amplitude is 1/1000 of the peak amplitude (the maximum of all modes).

  • taper (bool, optional) – Add a rapid ringup with timescale tau/10 at the beginning of the waveform to avoid the abrupt turn on of the ringdown. Each mode and overtone will have a different taper depending on its tau, the final taper being the superposition of all the tapers. Default is False.

Returns:

  • hplus (TimeSeries) – The plus phase of a ringdown with the lm modes specified and n overtones in time domain.

  • hcross (TimeSeries) – The cross phase of a ringdown with the lm modes specified and n overtones in time domain.

pycbc.waveform.ringdown.lm_amps_phases(**kwargs)[source]

Takes input_params and return dictionaries with amplitudes and phases of each overtone of a specific lm mode, checking that all of them are given. Will also look for dbetas and dphis. If (dphi|dbeta) (i.e., without a mode suffix) are provided, they will be used for all modes that don’t explicitly set a (dphi|dbeta){lmn}.

pycbc.waveform.ringdown.lm_arbitrary_harmonics(**kwargs)[source]

Take input_params and return dictionaries with arbitrary harmonics for each mode.

pycbc.waveform.ringdown.lm_deltaf(damping_times)[source]

Return the minimum delta_f of all the modes given, with delta_f given by the inverse of the time at which the amplitude of the ringdown falls to 1/1000 of the peak amplitude.

pycbc.waveform.ringdown.lm_deltat(freqs, damping_times)[source]

Return the minimum delta_t of all the modes given, with delta_t given by the inverse of the frequency at which the amplitude of the ringdown falls to 1/1000 of the peak amplitude.

pycbc.waveform.ringdown.lm_ffinal(freqs, damping_times)[source]

Return the maximum f_final of the modes given, with f_final the frequency at which the amplitude falls to 1/1000 of the peak amplitude

pycbc.waveform.ringdown.lm_freqs_taus(**kwargs)[source]

Take input_params and return dictionaries with frequencies and damping times of each overtone of a specific lm mode, checking that all of them are given.

pycbc.waveform.ringdown.lm_tfinal(damping_times)[source]

Return the maximum t_final of the modes given, with t_final the time at which the amplitude falls to 1/1000 of the peak amplitude

pycbc.waveform.ringdown.multimode_base(input_params, domain, freq_tau_approximant=False)[source]

Return a superposition of damped sinusoids in either time or frequency domains with parameters set by input_params.

Parameters:
  • input_params (dict) – Dictionary of parameters to generate the ringdowns with. See td_damped_sinusoid() and fd_damped_sinusoid() for supported parameters.

  • domain (string) – Choose domain of the waveform, either ‘td’ for time domain or ‘fd’ for frequency domain. If ‘td’ (‘fd’), the damped sinusoids will be generated with td_damped_sinusoid() (fd_damped_sinusoid()).

  • freq_tau_approximant ({False, bool}, optional) – Choose choose the waveform approximant to use. Either based on mass/spin (set to False, default), or on frequencies/damping times of the modes (set to True).

Returns:

  • hplus (TimeSeries) – The plus phase of a ringdown with the lm modes specified and n overtones in the chosen domain (time or frequency).

  • hcross (TimeSeries) – The cross phase of a ringdown with the lm modes specified and n overtones in the chosen domain (time or frequency).

pycbc.waveform.ringdown.parse_mode(lmn)[source]

Extracts overtones from an lmn.

pycbc.waveform.ringdown.props(obj, required, domain_args, **kwargs)[source]

Return a dictionary built from the combination of defaults, kwargs, and the attributes of the given object.

pycbc.waveform.ringdown.qnm_freq_decay(f_0, tau, decay)[source]

Return the frequency at which the amplitude of the ringdown falls to decay of the peak amplitude.

Parameters:
  • f_0 (float) – The ringdown-frequency, which gives the peak amplitude.

  • tau (float) – The damping time of the sinusoid.

  • decay (float) – The fraction of the peak amplitude.

Returns:

f_decay – The frequency at which the amplitude of the frequency-domain ringdown falls to decay of the peak amplitude.

Return type:

float

pycbc.waveform.ringdown.qnm_time_decay(tau, decay)[source]

Return the time at which the amplitude of the ringdown falls to decay of the peak amplitude.

Parameters:
  • tau (float) – The damping time of the sinusoid.

  • decay (float) – The fraction of the peak amplitude.

Returns:

t_decay – The time at which the amplitude of the time-domain ringdown falls to decay of the peak amplitude.

Return type:

float

pycbc.waveform.ringdown.spher_harms(harmonics='spherical', l=None, m=None, n=0, inclination=0.0, azimuthal=0.0, spin=None, pol=None, polnm=None)[source]

Return the +/-m harmonic polarizations.

This will return either spherical, spheroidal, or an arbitrary complex number depending on what harmonics is set to. If harmonics is set to arbitrary, then the “(+/-)m” harmonic will be \(e^{i \psi_(+/-)}\), where \(\psi_{+/-}\) are arbitrary angles provided by the user (using the pol(nm) arguments).

Parameters:
  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – The type of harmonic to generate. Default is spherical.

  • l (int, optional) – The l index. Must be provided if harmonics is ‘spherical’ or ‘spheroidal’.

  • m (int, optional) – The m index. Must be provided if harmonics is ‘spherical’ or ‘spheroidal’.

  • n (int, optional) – The overtone number. Only used if harmonics is ‘spheroidal’. Default is 0.

  • inclination (float, optional) – The inclination angle. Only used if harmonics is ‘spherical’ or ‘spheroidal’. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle. Only used if harmonics is ‘spherical’ or ‘spheroidal’. Default is 0.

  • spin (float, optional) – The dimensionless spin of the black hole. Must be provided if harmonics is ‘spheroidal’. Ignored otherwise.

  • pol (float, optional) – Angle (in radians) to use for arbitrary “+m” harmonic. Must be provided if harmonics is ‘arbitrary’. Ignored otherwise.

  • polnm (float, optional) – Angle (in radians) to use for arbitrary “-m” harmonic. Must be provided if harmonics is ‘arbitrary’. Ignored otherwise.

Returns:

  • xlm (complex) – The harmonic of the +m mode.

  • xlnm (complex) – The harmonic of the -m mode.

pycbc.waveform.ringdown.td_damped_sinusoid(f_0, tau, amp, phi, times, l=2, m=2, n=0, inclination=0.0, azimuthal=0.0, dphi=0.0, dbeta=0.0, harmonics='spherical', final_spin=None, pol=None, polnm=None)[source]

Return a time domain damped sinusoid (plus and cross polarizations) with central frequency f_0, damping time tau, amplitude amp and phase phi.

This returns the plus and cross polarization of the QNM, defined as \(h^{+,\times}_{l|m|n} = (\Re, \Im) \{ h_{l|m|n}(t)\right}\), where

\[\begin{split}h_{l|m|n}(t) &:= A_{lmn} X_{lmn}(\theta, \varphi) e^{-t/\tau_{lmn} + i(2\pi f_{lmn}t + \phi_{lmn})} \\ & + A_{l-mn} X_{l-mn}(\theta, \varphi) e^{-t/\tau_{lmn} - i(2\pi f_{lmn}t + \phi_{l-mn})}\end{split}\]

Here, the \(X_{lmn}(\theta, \varphi)\) are either the spherical or spheroidal harmonics, or an arbitrary complex number, depending on the input arguments. This uses the convention that the +/-m modes are related to each other via \(f_{l-mn} = -f_{lmn}\) and \(\tau_{l-mn} = \tau_{lmn}\). The amplitudes \(A_{l(-)mn}\) and phases \(\phi_{l(-)mn}\) of the +/-m modes are related to each other by:

\[\phi_{l-mn} = l\pi + \Delta \phi_{lmn} - \phi_{lmn}\]

and

\[\begin{split}A_{lmn} &= A^0_{lmn} \sqrt{2} \cos(\pi/4 + \Delta \beta_{lmn})\\ A_{lmn} &= A^0_{lmn} \sqrt{2} \sin(\pi/4 + \Delta \beta_{lmn}).\end{split}\]

Here, \(A^0_{lmn}\) is an overall fiducial amplitude (set by the amp) parameter, and

\[\begin{split}\Delta \beta_{lmn} &\in [-pi/4, pi/4], \\ \Delta \phi_{lmn} &\in (-pi, pi)\end{split}\]

are parameters that define the deviation from circular polarization. Circular polarization occurs when both \(\Delta \beta_{lmn}\) and \(\Delta \phi_{lmn}\) are zero (this is equivalent to assuming that \(h_{l-mn} = (-1)^l h_{lmn}^*\)).

Parameters:
  • f_0 (float) – The central frequency of the damped sinusoid, in Hz.

  • tau (float) – The damping time, in seconds.

  • amp (float) – The intrinsic amplitude of the QNM (\(A^0_{lmn}\) in the above).

  • phi (float) – The reference phase of the QNM (\(\phi_{lmn}\)) in the above.

  • times (array) – Array of times to use, where t=0 is considered the start of the ringdown. Times are assumed to be monotonically increasing. A taper of 10x the damping time will be used for any negative times.

  • l (int, optional) – The l index; default is 2.

  • m (int, optional) – The m index; default is 2.

  • n (int, optional) – The overtone index; default is 0.

  • inclination (float, optional) – The inclination angle (\(\theta\) in the above). Ignored if harmonics='arbitrary'. Default is 0.

  • azimuthal (float, optional) – The azimuthal angle (\(\varphi\) in the above). Ignored if harmonics='arbitrary'. Default is 0.

  • dphi (float, optional) – The difference in phase between the +m and -m mode (\(\Delta \phi_{lmn}\) in the above). Default is 0.

  • dbeta (float, optional) – The angular difference in the amplitudes of the +m and -m mode (\(\Delta \beta_{lmn}\) in the above). Default is 0. If this and dphi are both 0, will have circularly polarized waves.

  • harmonics ({'spherical', 'spheroidal', 'arbitrary'}, optional) – Which harmonics to use. See spher_harms() for details. Default is spherical.

  • final_spin (float, optional) – The dimensionless spin of the black hole. Only needed if harmonics='spheroidal'.

  • pol (float, optional) – Angle to use for +m arbitrary harmonics. Only needed if harmonics='arbitrary'. See spher_harms() for details.

  • polnm (float, optional) – Angle to use for -m arbitrary harmonics. Only needed if harmonics='arbitrary'. See spher_harms() for details.

Returns:

  • hplus (numpy.ndarray) – The plus polarization.

  • hcross (numpy.ndarray) – The cross polarization.

pycbc.waveform.ringdown.td_output_vector(freqs, damping_times, taper=False, delta_t=None, t_final=None)[source]

Return an empty TimeSeries with the appropriate size to fit all the quasi-normal modes present in freqs, damping_times

pycbc.waveform.sinegauss module

Generation of sine-Gaussian bursty type things

pycbc.waveform.sinegauss.fd_sine_gaussian(amp, quality, central_frequency, fmin, fmax, delta_f)[source]

Generate a Fourier domain sine-Gaussian

Parameters:
  • amp (float) – Amplitude of the sine-Gaussian

  • quality (float) – The quality factor

  • central_frequency (float) – The central frequency of the sine-Gaussian

  • fmin (float) – The minimum frequency to generate the sine-Gaussian. This determines the length of the output vector.

  • fmax (float) – The maximum frequency to generate the sine-Gaussian

  • delta_f (float) – The size of the frequency step

Returns:

sg – A Fourier domain sine-Gaussian

Return type:

pycbc.types.Frequencyseries

pycbc.waveform.spa_tmplt module

This module contains functions for generating common SPA template precalculated vectors.

pycbc.waveform.spa_tmplt.findchirp_chirptime(m1, m2, fLower, porder)[source]
pycbc.waveform.spa_tmplt.spa_amplitude_factor(**kwds)[source]
pycbc.waveform.spa_tmplt.spa_distance(psd, mass1, mass2, lower_frequency_cutoff, snr=8)[source]

Return the distance at a given snr (default=8) of the SPA TaylorF2 template.

pycbc.waveform.spa_tmplt.spa_length_in_time(**kwds)[source]

Returns the length in time of the template, based on the masses, PN order, and low-frequency cut-off.

pycbc.waveform.spa_tmplt.spa_tmplt(**kwds)[source]

Generate a minimal TaylorF2 approximant with optimizations for the sin/cos

pycbc.waveform.spa_tmplt.spa_tmplt_end(**kwds)[source]
pycbc.waveform.spa_tmplt.spa_tmplt_engine(htilde, kmin, phase_order, delta_f, piM, pfaN, pfa2, pfa3, pfa4, pfa5, pfl5, pfa6, pfl6, pfa7, amp_factor)[source]

Calculate the spa tmplt phase

pycbc.waveform.spa_tmplt.spa_tmplt_norm(psd, length, delta_f, f_lower)[source]
pycbc.waveform.spa_tmplt.spa_tmplt_precondition(length, delta_f, kmin=0)[source]

Return the amplitude portion of the TaylorF2 approximant, used to precondition the strain data. The result is cached, and so should not be modified, only read.

pycbc.waveform.spa_tmplt_cpu module

pycbc.waveform.spa_tmplt_cpu.cbrt_lookup(vmax, delta)
pycbc.waveform.spa_tmplt_cpu.get_cbrt(vmax, delta)
pycbc.waveform.spa_tmplt_cpu.get_log(vmax, delta)
pycbc.waveform.spa_tmplt_cpu.logv_lookup(vmax, delta)
pycbc.waveform.spa_tmplt_cpu.spa_tmplt_engine(htilde, kmin, phase_order, delta_f, piM, pfaN, pfa2, pfa3, pfa4, pfa5, pfl5, pfa6, pfl6, pfa7, amp_factor)

Calculate the spa tmplt phase

pycbc.waveform.spa_tmplt_cpu.spa_tmplt_inline_sequence(piM, pfaN, pfa2, pfa3, pfa4, pfa5, pfl5, pfa6, pfl6, pfa7, ampc, _fvals, _htilde)

pycbc.waveform.supernovae module

Generate core-collapse supernovae waveform for core bounce and subsequent postbounce oscillations.

pycbc.waveform.supernovae.get_corecollapse_bounce(**kwargs)[source]

Generates core bounce and postbounce waveform by using principal component basis vectors from a .hdf file. The waveform parameters are the coefficients of the principal components and the distance. The number of principal components used can also be varied.

pycbc.waveform.utils module

This module contains convenience utilities for manipulating waveforms

pycbc.waveform.utils.amplitude_from_frequencyseries(htilde)[source]

Returns the amplitude of the given frequency-domain waveform as a FrequencySeries.

Parameters:

htilde (FrequencySeries) – The waveform to get the amplitude of.

Returns:

The amplitude of the waveform as a function of frequency.

Return type:

FrequencySeries

pycbc.waveform.utils.amplitude_from_polarizations(h_plus, h_cross)[source]

Return gravitational wave amplitude

Return the gravitation-wave amplitude from the h_plus and h_cross polarizations of the waveform.

Parameters:
  • h_plus (TimeSeries) – An PyCBC TmeSeries vector that contains the plus polarization of the gravitational waveform.

  • h_cross (TimeSeries) – A PyCBC TmeSeries vector that contains the cross polarization of the gravitational waveform.

Returns:

GWAmplitude – A TimeSeries containing the gravitational wave amplitude.

Return type:

TimeSeries

Examples

>>> from pycbc.waveform import get_td_waveform, phase_from_polarizations
>>> hp, hc = get_td_waveform(approximant="TaylorT4", mass1=10, mass2=10,
                     f_lower=30, delta_t=1.0/4096)
>>> amp = amplitude_from_polarizations(hp, hc)
pycbc.waveform.utils.apply_fd_time_shift(htilde, shifttime, kmin=0, fseries=None, copy=True)[source]

Shifts a frequency domain waveform in time. The shift applied is shiftime - htilde.epoch.

Parameters:
  • htilde (FrequencySeries) – The waveform frequency series.

  • shifttime (float) – The time to shift the frequency series to.

  • kmin ({0, int}) – The starting index of htilde to apply the time shift. Default is 0.

  • fseries ({None, numpy array}) – The frequencies of each element in htilde. This is only needed if htilde is not sampled at equal frequency steps.

  • copy ({True, bool}) – Make a copy of htilde before applying the time shift. If False, the time shift will be applied to htilde’s data.

Returns:

A frequency series with the waveform shifted to the new time. If makecopy is True, will be a new frequency series; if makecopy is False, will be the same as htilde.

Return type:

FrequencySeries

pycbc.waveform.utils.apply_fseries_time_shift(htilde, dt, kmin=0, copy=True)[source]

Shifts a frequency domain waveform in time. The waveform is assumed to be sampled at equal frequency intervals.

pycbc.waveform.utils.ceilpow2(n)[source]

convenience function to determine a power-of-2 upper frequency limit

pycbc.waveform.utils.coalign_waveforms(h1, h2, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, resize=True)[source]

Return two time series which are aligned in time and phase.

The alignment is only to the nearest sample point and all changes to the phase are made to the first input waveform. Waveforms should not be split accross the vector boundary. If it is, please use roll or cyclic time shift to ensure that the entire signal is contiguous in the time series.

Parameters:
  • h1 (pycbc.types.TimeSeries) – The first waveform to align.

  • h2 (pycbc.types.TimeSeries) – The second waveform to align.

  • psd ({None, pycbc.types.FrequencySeries}) – A psd to weight the alignment

  • low_frequency_cutoff ({None, float}) – The low frequency cutoff to weight the matching in Hz.

  • high_frequency_cutoff ({None, float}) – The high frequency cutoff to weight the matching in Hz.

  • resize (Optional, {True, boolean}) – If true, the vectors will be resized to match each other. If false, they must be the same length and even in length

Returns:

  • h1 (pycbc.types.TimeSeries) – The shifted waveform to align with h2

  • h2 (pycbc.type.TimeSeries) – The resized (if necessary) waveform to align with h1.

pycbc.waveform.utils.fd_taper(out, start, end, beta=8, side='left')[source]

Applies a taper to the given FrequencySeries.

A half-kaiser window is used for the roll-off.

Parameters:
  • out (FrequencySeries) – The FrequencySeries to taper.

  • start (float) – The frequency (in Hz) to start the taper window.

  • end (float) – The frequency (in Hz) to end the taper window.

  • beta (int, optional) – The beta parameter to use for the Kaiser window. See scipy.signal.kaiser for details. Default is 8.

  • side ({'left', 'right'}) – The side to apply the taper to. If 'left' ('right'), the taper will roll up (down) between start and end, with all values before start (after end) set to zero. Default is 'left'.

Returns:

The tapered frequency series.

Return type:

FrequencySeries

pycbc.waveform.utils.fd_to_td(htilde, delta_t=None, left_window=None, right_window=None, left_beta=8, right_beta=8)[source]

Converts a FD waveform to TD.

A window can optionally be applied using fd_taper to the left or right side of the waveform before being converted to the time domain.

Parameters:
  • htilde (FrequencySeries) – The waveform to convert.

  • delta_t (float, optional) – Make the returned time series have the given delta_t.

  • left_window (tuple of float, optional) – A tuple giving the start and end frequency of the FD taper to apply on the left side. If None, no taper will be applied on the left.

  • right_window (tuple of float, optional) – A tuple giving the start and end frequency of the FD taper to apply on the right side. If None, no taper will be applied on the right.

  • left_beta (int, optional) – The beta parameter to use for the left taper. See fd_taper for details. Default is 8.

  • right_beta (int, optional) – The beta parameter to use for the right taper. Default is 8.

Returns:

The time-series representation of htilde.

Return type:

TimeSeries

pycbc.waveform.utils.frequency_from_polarizations(h_plus, h_cross)[source]

Return gravitational wave frequency

Return the gravitation-wave frequency as a function of time from the h_plus and h_cross polarizations of the waveform. It is 1 bin shorter than the input vectors and the sample times are advanced half a bin.

Parameters:
  • h_plus (TimeSeries) – A PyCBC TimeSeries vector that contains the plus polarization of the gravitational waveform.

  • h_cross (TimeSeries) – A PyCBC TimeSeries vector that contains the cross polarization of the gravitational waveform.

Returns:

GWFrequency – A TimeSeries containing the gravitational wave frequency as a function of time.

Return type:

TimeSeries

Examples

>>> from pycbc.waveform import get_td_waveform, phase_from_polarizations
>>> hp, hc = get_td_waveform(approximant="TaylorT4", mass1=10, mass2=10,
                     f_lower=30, delta_t=1.0/4096)
>>> freq = frequency_from_polarizations(hp, hc)
pycbc.waveform.utils.phase_from_frequencyseries(htilde, remove_start_phase=True)[source]

Returns the phase from the given frequency-domain waveform. This assumes that the waveform has been sampled finely enough that the phase cannot change by more than pi radians between each step.

Parameters:
  • htilde (FrequencySeries) – The waveform to get the phase for; must be a complex frequency series.

  • remove_start_phase ({True, bool}) – Subtract the initial phase before returning.

Returns:

The phase of the waveform as a function of frequency.

Return type:

FrequencySeries

pycbc.waveform.utils.phase_from_polarizations(h_plus, h_cross, remove_start_phase=True)[source]

Return gravitational wave phase

Return the gravitation-wave phase from the h_plus and h_cross polarizations of the waveform. The returned phase is always positive and increasing with an initial phase of 0.

Parameters:
  • h_plus (TimeSeries) – An PyCBC TmeSeries vector that contains the plus polarization of the gravitational waveform.

  • h_cross (TimeSeries) – A PyCBC TmeSeries vector that contains the cross polarization of the gravitational waveform.

Returns:

  • GWPhase (TimeSeries) – A TimeSeries containing the gravitational wave phase.

  • Examples

  • ——–s

  • >>> from pycbc.waveform import get_td_waveform, phase_from_polarizations

  • >>> hp, hc = get_td_waveform(approximant=”TaylorT4”, mass1=10, mass2=10, – f_lower=30, delta_t=1.0/4096)

  • >>> phase = phase_from_polarizations(hp, hc)

pycbc.waveform.utils.taper_timeseries(tsdata, tapermethod=None, return_lal=False)[source]

Taper either or both ends of a time series using wrapped LALSimulation functions

Parameters:
  • tsdata (TimeSeries) – Series to be tapered, dtype must be either float32 or float64

  • tapermethod (string) – Should be one of (‘TAPER_NONE’, ‘TAPER_START’, ‘TAPER_END’, ‘TAPER_STARTEND’, ‘start’, ‘end’, ‘startend’) - NB ‘TAPER_NONE’ will not change the series!

  • return_lal (Boolean) – If True, return a wrapped LAL time series object, else return a PyCBC time series.

pycbc.waveform.utils.td_taper(out, start, end, beta=8, side='left')[source]

Applies a taper to the given TimeSeries.

A half-kaiser window is used for the roll-off.

Parameters:
  • out (TimeSeries) – The TimeSeries to taper.

  • start (float) – The time (in s) to start the taper window.

  • end (float) – The time (in s) to end the taper window.

  • beta (int, optional) – The beta parameter to use for the Kaiser window. See scipy.signal.kaiser for details. Default is 8.

  • side ({'left', 'right'}) – The side to apply the taper to. If 'left' ('right'), the taper will roll up (down) between start and end, with all values before start (after end) set to zero. Default is 'left'.

Returns:

The tapered time series.

Return type:

TimeSeries

pycbc.waveform.utils.time_from_frequencyseries(htilde, sample_frequencies=None, discont_threshold=3.1101767270538954)[source]

Computes time as a function of frequency from the given frequency-domain waveform. This assumes the stationary phase approximation. Any frequencies lower than the first non-zero value in htilde are assigned the time at the first non-zero value. Times for any frequencies above the next-to-last non-zero value in htilde will be assigned the time of the next-to-last non-zero value.

Note

Some waveform models (e.g., SEOBNRv2_ROM_DoubleSpin) can have discontinuities in the phase towards the end of the waveform due to numerical error. We therefore exclude any points that occur after a discontinuity in the phase, as the time estimate becomes untrustworthy beyond that point. What determines a discontinuity in the phase is set by the discont_threshold. To turn this feature off, just set discont_threshold to a value larger than pi (due to the unwrapping of the phase, no two points can have a difference > pi).

Parameters:
  • htilde (FrequencySeries) – The waveform to get the time evolution of; must be complex.

  • sample_frequencies ({None, array}) – The frequencies at which the waveform is sampled. If None, will retrieve from htilde.sample_frequencies.

  • discont_threshold ({0.99*pi, float}) – If the difference in the phase changes by more than this threshold, it is considered to be a discontinuity. Default is 0.99*pi.

Returns:

The time evolution of the waveform as a function of frequency.

Return type:

FrequencySeries

pycbc.waveform.utils_cpu module

This module contains the CPU-specific code for convenience utilities for manipulating waveforms

pycbc.waveform.utils_cpu.apply_fseries_time_shift(htilde, dt, kmin=0, copy=True)

Shifts a frequency domain waveform in time. The waveform is assumed to be sampled at equal frequency intervals.

pycbc.waveform.utils_cpu.fstimeshift(freqseries, phi, kmin, kmax)
pycbc.waveform.utils_cpu.fstimeshift32(freqseries, phi, kmin, kmax)

pycbc.waveform.waveform module

Convenience functions to genenerate gravitational wave templates and waveforms.

exception pycbc.waveform.waveform.FailedWaveformError[source]

Bases: Exception

This should be raised if a waveform fails to generate.

exception pycbc.waveform.waveform.NoWaveformError[source]

Bases: Exception

This should be raised if generating a waveform would just result in all zeros being returned, e.g., if a requested f_final is <= f_lower.

pycbc.waveform.waveform.fd_approximants(scheme=<pycbc.scheme.DefaultScheme object>)[source]

Return a list containing the available fourier domain approximants for the given processing scheme.

pycbc.waveform.waveform.filter_approximants(scheme=<pycbc.scheme.DefaultScheme object>)[source]

Return a list of fourier domain approximants including those written specifically as templates.

pycbc.waveform.waveform.get_fd_det_waveform(template=None, **kwargs)[source]

Return a frequency domain gravitational waveform. The waveform generator includes detector response.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. An example would be a row in an xml table.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • coa_phase ({0.0, float}) – Coalesence phase of the binary (in rad).

  • inclination ({0.0, float}) – Inclination (rad), defined as the angle between the orbital angular momentum L and the line-of-sight at the reference frequency.

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • delta_f ({None, float}) – The frequency step used to generate the waveform (in Hz).

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

  • f_final ({0, float}) – The ending frequency of the waveform. The default (0) indicates that the choice is made by the respective approximant.

  • f_final_func ({, str}) – Use the given frequency function to compute f_final based on the parameters of the waveform.

Returns:

The detector-frame waveform (with detector response) in frequency domain. Keys are requested data channels, values are FrequencySeries.

Return type:

dict

pycbc.waveform.waveform.get_fd_det_waveform_sequence(template=None, **kwds)[source]

Return values of the waveform evaluated at the sequence of frequency points. The waveform generator includes detector response.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • coa_phase ({0.0, float}) – Coalesence phase of the binary (in rad).

  • inclination ({0.0, float}) – Inclination (rad), defined as the angle between the orbital angular momentum L and the line-of-sight at the reference frequency.

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • sample_points ({None, Array}) – An array of the frequencies (in Hz) at which to generate the waveform.

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

Returns:

The detector-frame waveform (with detector response) in frequency domain evaluated at the frequency points. Keys are requested data channels, values are FrequencySeries.

Return type:

dict

pycbc.waveform.waveform.get_fd_waveform(template=None, **kwargs)[source]

Return a frequency domain gravitational waveform.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • coa_phase ({0.0, float}) – Coalesence phase of the binary (in rad).

  • inclination ({0.0, float}) – Inclination (rad), defined as the angle between the orbital angular momentum L and the line-of-sight at the reference frequency.

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • delta_f ({None, float}) – The frequency step used to generate the waveform (in Hz).

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

  • f_final ({0, float}) – The ending frequency of the waveform. The default (0) indicates that the choice is made by the respective approximant.

  • f_final_func ({, str}) – Use the given frequency function to compute f_final based on the parameters of the waveform.

Returns:

  • hplustilde (FrequencySeries) – The plus phase of the waveform in frequency domain.

  • hcrosstilde (FrequencySeries) – The cross phase of the waveform in frequency domain.

pycbc.waveform.waveform.get_fd_waveform_from_td(**params)[source]

Return time domain version of fourier domain approximant.

This returns a frequency domain version of a fourier domain approximant, with padding and tapering at the start of the waveform.

Parameters:

params (dict) – The parameters defining the waveform to generator. See get_td_waveform.

Returns:

  • hp (pycbc.types.FrequencySeries) – Plus polarization time series

  • hc (pycbc.types.FrequencySeries) – Cross polarization time series

pycbc.waveform.waveform.get_fd_waveform_sequence(template=None, **kwds)[source]

Return values of the waveform evaluated at the sequence of frequency points. The waveform generator doesn’t include detector response.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • coa_phase ({0.0, float}) – Coalesence phase of the binary (in rad).

  • inclination ({0.0, float}) – Inclination (rad), defined as the angle between the orbital angular momentum L and the line-of-sight at the reference frequency.

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • sample_points ({None, Array}) – An array of the frequencies (in Hz) at which to generate the waveform.

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

Returns:

  • hplustilde (Array) – The plus phase of the waveform in frequency domain evaluated at the

  • frequency points.

  • hcrosstilde (Array) – The cross phase of the waveform in frequency domain evaluated at the

  • frequency points.

pycbc.waveform.waveform.get_sgburst_waveform(template=None, **kwargs)[source]

Return the plus and cross polarizations of a time domain sine-Gaussian burst waveform.

Parameters:
  • template (object) – An object that has attached properties. This can be used to subsitute for keyword arguments. A common example would be a row in an xml table.

  • approximant (string) – A string that indicates the chosen approximant. See td_approximants for available options.

  • q (float) – The quality factor of a sine-Gaussian burst

  • frequency (float) – The centre-frequency of a sine-Gaussian burst

  • delta_t (float) – The time step used to generate the waveform

  • hrss (float) – The strain rss

  • amplitude (float) – The strain amplitude

Returns:

  • hplus (TimeSeries) – The plus polarization of the waveform.

  • hcross (TimeSeries) – The cross polarization of the waveform.

pycbc.waveform.waveform.get_td_det_waveform_from_fd_det(template=None, rwrap=None, **params)[source]

Return time domain version of fourier domain approximant which includes detector response, with padding and tapering at the start of the waveform.

Parameters:
  • rwrap (float) – Cyclic time shift parameter in seconds. A fudge factor to ensure that the entire time series is contiguous in the array and not wrapped around the end.

  • params (dict) – The parameters defining the waveform to generator. See get_fd_det_waveform.

Returns:

The detector-frame waveform (with detector response) in time domain. Keys are requested data channels.

Return type:

dict

pycbc.waveform.waveform.get_td_waveform(template=None, **kwargs)[source]

Return the plus and cross polarizations of a time domain waveform.

Parameters:
  • template (object) – An object that has attached properties. This can be used to subsitute for keyword arguments. A common example would be a row in an xml table.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • coa_phase ({0.0, float}) – Coalesence phase of the binary (in rad).

  • inclination ({0.0, float}) – Inclination (rad), defined as the angle between the orbital angular momentum L and the line-of-sight at the reference frequency.

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • delta_t ({None, float}) – The time step used to generate the waveform (in s).

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

  • numrel_data ({, str}) – Sets the NR flags; only needed for NR waveforms.

  • frame_axis – Allow to choose among orbital_l, view and total_j

  • modes_choice – Allow to turn on among orbital_l, view and total_j

  • side_bands – Flag for generating sidebands

  • mode_array – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

Returns:

  • hplus (TimeSeries) – The plus polarization of the waveform.

  • hcross (TimeSeries) – The cross polarization of the waveform.

pycbc.waveform.waveform.get_td_waveform_from_fd(rwrap=None, **params)[source]

Return time domain version of fourier domain approximant.

This returns a time domain version of a fourier domain approximant, with padding and tapering at the start of the waveform.

Parameters:
  • rwrap (float) – Cyclic time shift parameter in seconds. A fudge factor to ensure that the entire time series is contiguous in the array and not wrapped around the end.

  • params (dict) – The parameters defining the waveform to generator. See get_fd_waveform.

Returns:

  • hp (pycbc.types.TimeSeries) – Plus polarization time series

  • hc (pycbc.types.TimeSeries) – Cross polarization time series

pycbc.waveform.waveform.get_template_amplitude_norm(template=None, **kwargs)[source]

Return additional constant template normalization. This only affects the effective distance calculation. Returns None for all templates with a physically meaningful amplitude.

pycbc.waveform.waveform.get_two_pol_waveform_filter(outplus, outcross, template, **kwargs)[source]

Return a frequency domain waveform filter for the specified approximant. Unlike get_waveform_filter this function returns both h_plus and h_cross components of the waveform, which are needed for searches where h_plus and h_cross are not related by a simple phase shift.

pycbc.waveform.waveform.get_waveform_end_frequency(template=None, **kwargs)[source]

Return the stop frequency of a template

pycbc.waveform.waveform.get_waveform_filter(out, template=None, **kwargs)[source]

Return a frequency domain waveform filter for the specified approximant

pycbc.waveform.waveform.get_waveform_filter_length_in_time(approximant, template=None, **kwargs)[source]

For filter templates, return the length in time of the template.

pycbc.waveform.waveform.get_waveform_filter_norm(approximant, psd, length, delta_f, f_lower)[source]

Return the normalization vector for the approximant

pycbc.waveform.waveform.print_fd_approximants()[source]
pycbc.waveform.waveform.print_sgburst_approximants()[source]
pycbc.waveform.waveform.print_td_approximants()[source]
pycbc.waveform.waveform.sgburst_approximants(scheme=<pycbc.scheme.DefaultScheme object>)[source]

Return a list containing the available time domain sgbursts for the given processing scheme.

pycbc.waveform.waveform.td_approximants(scheme=<pycbc.scheme.DefaultScheme object>)[source]

Return a list containing the available time domain approximants for the given processing scheme.

pycbc.waveform.waveform.td_waveform_to_fd_waveform(waveform, out=None, length=None, buffer_length=100)[source]

Convert a time domain into a frequency domain waveform by FFT. As a waveform is assumed to “wrap” in the time domain one must be careful to ensure the waveform goes to 0 at both “boundaries”. To ensure this is done correctly the waveform must have the epoch set such the merger time is at t=0 and the length of the waveform should be shorter than the desired length of the FrequencySeries (times 2 - 1) so that zeroes can be suitably pre- and post-pended before FFTing. If given, out is a memory array to be used as the output of the FFT. If not given memory is allocated internally. If present the length of the returned FrequencySeries is determined from the length out. If out is not given the length can be provided expicitly, or it will be chosen as the nearest power of 2. If choosing length explicitly the waveform length + buffer_length is used when choosing the nearest binary number so that some zero padding is always added.

pycbc.waveform.waveform.waveform_norm_exists(approximant)[source]

pycbc.waveform.waveform_modes module

Provides functions and utilities for generating waveforms mode-by-mode.

pycbc.waveform.waveform_modes.default_modes(approximant)[source]

Returns the default modes for the given approximant.

pycbc.waveform.waveform_modes.fd_waveform_mode_approximants()[source]

Frequency domain approximants that will return separate modes.

pycbc.waveform.waveform_modes.get_fd_waveform_modes(template=None, **kwargs)[source]

Generates frequency domain waveforms, but does not sum over the modes.

The returned values are the frequency-domain equivalents of the real and imaginary parts of the complex \(\mathfrak{h}_{\ell m}(t)\) time series. In other words, the returned values are equivalent to the Fourier Transform of the two time series returned by get_td_waveform_modes(); see that function for more details.

Parameters:
  • template (object) – An object that has attached properties. This can be used to subsitute for keyword arguments.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • delta_f ({None, float}) – The frequency step used to generate the waveform (in Hz).

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

  • f_final ({0, float}) – The ending frequency of the waveform. The default (0) indicates that the choice is made by the respective approximant.

  • f_final_func ({, str}) – Use the given frequency function to compute f_final based on the parameters of the waveform.

Returns:

  • ulm (dict) – Dictionary of mode tuples -> fourier transform of the real part of the hlm time series, as a pycbc.types.FrequencySeries.

  • vlm (dict) – Dictionary of mode tuples -> fourier transform of the imaginary part of the hlm time series, as a pycbc.types.FrequencySeries.

pycbc.waveform.waveform_modes.get_glm(l, m, theta)[source]

The maginitude of the \({}_{-2}Y_{\ell m}\).

The spin-weighted spherical harmonics can be written as \({}_{-2}Y_{\ell m}(\theta, \phi) = g_{\ell m}(\theta)e^{i m \phi}\). This returns the g_{ell m}(theta) part. Note that this is real.

Parameters:
  • l (int) – The \(\ell\) index of the spherical harmonic.

  • m (int) – The \(m\) index of the spherical harmonic.

  • theta (float) – The polar angle (in radians).

Returns:

The amplitude of the harmonic at the given polar angle.

Return type:

float

pycbc.waveform.waveform_modes.get_imrphenomxh_modes(**params)[source]

Generates IMRPhenomXHM waveforms mode-by-mode.

pycbc.waveform.waveform_modes.get_nrhybsur_modes(**params)[source]

Generates NRHybSur3dq8 waveform mode-by-mode.

All waveform parameters should be provided as keyword arguments. Recognized parameters are listed below. Unrecognized arguments are ignored.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • approximant (str) – The approximant to generate. Must be one of the NRHyb* models.

  • delta_t ({None, float}) – The time step used to generate the waveform (in s).

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • f_ref ({0, float}) – The reference frequency.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

Returns:

Dictionary of (l, m) -> (h_+, -h_x) TimeSeries.

Return type:

dict

pycbc.waveform.waveform_modes.get_nrsur_modes(**params)[source]

Generates NRSurrogate waveform mode-by-mode.

All waveform parameters should be provided as keyword arguments. Recognized parameters are listed below. Unrecognized arguments are ignored.

Parameters:
  • template (object) – An object that has attached properties. This can be used to substitute for keyword arguments. A common example would be a row in an xml table.

  • approximant (str) – The approximant to generate. Must be one of the NRSur* models.

  • delta_t ({None, float}) – The time step used to generate the waveform (in s).

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • f_ref ({0, float}) – The reference frequency.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

Returns:

Dictionary of (l, m) -> (h_+, -h_x) TimeSeries.

Return type:

dict

pycbc.waveform.waveform_modes.get_td_waveform_modes(template=None, **kwargs)[source]

Generates time domain waveforms, but does not sum over the modes.

The returned values are the real and imaginary parts of the complex \(\mathfrak{h}_{\ell m}(t)\). These are defined such that the plus and cross polarizations \(h_{+,\times}\) are:

\[h_{+,\times}(\theta, \phi; t) = (\Re, -\Im) \sum_{\ell m} {}_{-2}Y_{\ell m}(\theta, \phi) \mathfrak{h}_{\ell m}(t).\]
Parameters:
  • template (object) – An object that has attached properties. This can be used to subsitute for keyword arguments.

  • mass1 ({None, float}) – The mass of the first component object in the binary (in solar masses).

  • mass2 ({None, float}) – The mass of the second component object in the binary (in solar masses).

  • spin1x ({0.0, float}) – The x component of the first binary component’s dimensionless spin.

  • spin1y ({0.0, float}) – The y component of the first binary component’s dimensionless spin.

  • spin1z ({0.0, float}) – The z component of the first binary component’s dimensionless spin.

  • spin2x ({0.0, float}) – The x component of the second binary component’s dimensionless spin.

  • spin2y ({0.0, float}) – The y component of the second binary component’s dimensionless spin.

  • spin2z ({0.0, float}) – The z component of the second binary component’s dimensionless spin.

  • eccentricity ({0.0, float}) – Eccentricity.

  • lambda1 ({None, float}) – The dimensionless tidal deformability parameter of object 1.

  • lambda2 ({None, float}) – The dimensionless tidal deformability parameter of object 2.

  • dquad_mon1 ({None, float}) – Quadrupole-monopole parameter / m_1^5 -1.

  • dquad_mon2 ({None, float}) – Quadrupole-monopole parameter / m_2^5 -1.

  • lambda_octu1 ({None, float}) – The octupolar tidal deformability parameter of object 1.

  • lambda_octu2 ({None, float}) – The octupolar tidal deformability parameter of object 2.

  • quadfmode1 ({None, float}) – The quadrupolar f-mode angular frequency of object 1.

  • quadfmode2 ({None, float}) – The quadrupolar f-mode angular frequency of object 2.

  • octufmode1 ({None, float}) – The octupolar f-mode angular frequency of object 1.

  • octufmode2 ({None, float}) – The octupolar f-mode angular frequency of object 2.

  • dchi0 ({0.0, float}) – 0PN testingGR parameter.

  • dchi1 ({0.0, float}) – 0.5PN testingGR parameter.

  • dchi2 ({0.0, float}) – 1PN testingGR parameter.

  • dchi3 ({0.0, float}) – 1.5PN testingGR parameter.

  • dchi4 ({0.0, float}) – 2PN testingGR parameter.

  • dchi5 ({0.0, float}) – 2.5PN testingGR parameter.

  • dchi5l ({0.0, float}) – 2.5PN logrithm testingGR parameter.

  • dchi6 ({0.0, float}) – 3PN testingGR parameter.

  • dchi6l ({0.0, float}) – 3PN logrithm testingGR parameter.

  • dchi7 ({0.0, float}) – 3.5PN testingGR parameter.

  • dalpha1 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha2 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha3 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha4 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dalpha5 ({0.0, float}) – Merger-ringdown testingGR parameter.

  • dbeta1 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta2 ({0.0, float}) – Intermediate testingGR parameter.

  • dbeta3 ({0.0, float}) – Intermediate testingGR parameter.

  • distance ({1.0, float}) – Luminosity distance to the binary (in Mpc).

  • long_asc_nodes ({0.0, float}) – Longitude of ascending nodes axis (rad).

  • mean_per_ano ({0.0, float}) – Mean anomaly of the periastron (rad).

  • delta_t ({None, float}) – The time step used to generate the waveform (in s).

  • f_lower ({None, float}) – The starting frequency of the waveform (in Hz).

  • approximant ({None, str}) – A string that indicates the chosen approximant.

  • f_ref ({0, float}) – The reference frequency.

  • phase_order ({-1, int}) – The pN order of the orbital phase. The default of -1 indicates that all implemented orders are used.

  • spin_order ({-1, int}) – The pN order of the spin corrections. The default of -1 indicates that all implemented orders are used.

  • tidal_order ({-1, int}) – The pN order of the tidal corrections. The default of -1 indicates that all implemented orders are used.

  • amplitude_order ({-1, int}) – The pN order of the amplitude. The default of -1 indicates that all implemented orders are used.

  • eccentricity_order ({-1, int}) – The pN order of the eccentricity corrections.The default of -1 indicates that all implemented orders are used.

  • frame_axis ({0, int}) – Allow to choose among orbital_l, view and total_j

  • modes_choice ({0, int}) – Allow to turn on among orbital_l, view and total_j

  • side_bands ({0, int}) – Flag for generating sidebands

  • mode_array ({None, list}) – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

  • numrel_data ({, str}) – Sets the NR flags; only needed for NR waveforms.

  • frame_axis – Allow to choose among orbital_l, view and total_j

  • modes_choice – Allow to turn on among orbital_l, view and total_j

  • side_bands – Flag for generating sidebands

  • mode_array – Choose which (l,m) modes to include when generating a waveform. Only if approximant supports this feature.By default pass None and let lalsimulation use it’s default behaviour.Example: mode_array = [ [2,2], [2,-2] ]

Returns:

  • ulm (dict) – Dictionary of mode tuples -> real part of the hlm, as a pycbc.types.TimeSeries.

  • vlm (dict) – Dictionary of mode tuples -> imaginary part of the hlm, as a pycbc.types.TimeSeries.

pycbc.waveform.waveform_modes.sum_modes(hlms, inclination, phi)[source]

Applies spherical harmonics and sums modes to produce a plus and cross polarization.

Parameters:
  • hlms (dict) – Dictionary of (l, m) -> complex hlm. The hlm may be a complex number or array, or complex TimeSeries. All modes in the dictionary will be summed.

  • inclination (float) – The inclination to use.

  • phi (float) – The phase to use.

Returns:

The plus and cross polarization as a complex number. The real part gives the plus, the negative imaginary part the cross.

Return type:

complex float or array

pycbc.waveform.waveform_modes.td_waveform_mode_approximants()[source]

Time domain approximants that will return separate modes.

Module contents