pycbc.filter package

Submodules

pycbc.filter.autocorrelation module

This modules provides functions for calculating the autocorrelation function and length of a data series.

pycbc.filter.autocorrelation.calculate_acf(data, delta_t=1.0, unbiased=False)[source]

Calculates the one-sided autocorrelation function.

Calculates the autocorrelation function (ACF) and returns the one-sided ACF. The ACF is defined as the autocovariance divided by the variance. The ACF can be estimated using

\[\hat{R}(k) = \frac{1}{n \sigma^{2}} \sum_{t=1}^{n-k} \left( X_{t} - \mu \right) \left( X_{t+k} - \mu \right)\]

Where \(\hat{R}(k)\) is the ACF, \(X_{t}\) is the data series at time t, \(\mu\) is the mean of \(X_{t}\), and \(\sigma^{2}\) is the variance of \(X_{t}\).

Parameters:
  • data (TimeSeries or numpy.array) – A TimeSeries or numpy.array of data.

  • delta_t (float) – The time step of the data series if it is not a TimeSeries instance.

  • unbiased (bool) – If True the normalization of the autocovariance function is n-k instead of n. This is called the unbiased estimation of the autocovariance. Note that this does not mean the ACF is unbiased.

Returns:

acf – If data is a TimeSeries then acf will be a TimeSeries of the one-sided ACF. Else acf is a numpy.array.

Return type:

numpy.array

pycbc.filter.autocorrelation.calculate_acl(data, m=5, dtype=<class 'int'>)[source]

Calculates the autocorrelation length (ACL).

Given a normalized autocorrelation function \(\rho[i]\) (by normalized, we mean that \(\rho[0] = 1\)), the ACL \(\tau\) is:

\[\tau = 1 + 2 \sum_{i=1}^{K} \rho[i].\]

The number of samples used \(K\) is found by using the first point such that:

\[m \tau[K] \leq K,\]

where \(m\) is a tuneable parameter (default = 5). If no such point exists, then the given data set it too short to estimate the ACL; in this case inf is returned.

This algorithm for computing the ACL is taken from:

  1. Madras and A.D. Sokal, J. Stat. Phys. 50, 109 (1988).

Parameters:
  • data (TimeSeries or array) – A TimeSeries of data.

  • m (int) – The number of autocorrelation lengths to use for determining the window size \(K\) (see above).

  • dtype (int or float) – The datatype of the output. If the dtype was set to int, then the ceiling is returned.

Returns:

acl – The autocorrelation length. If the ACL cannot be estimated, returns numpy.inf.

Return type:

int or float

pycbc.filter.matchedfilter module

This modules provides functions for matched filtering along with associated utilities.

class pycbc.filter.matchedfilter.LiveBatchMatchedFilter(templates, snr_threshold, chisq_bins, sg_chisq, maxelements=134217728, snr_abort_threshold=None, newsnr_threshold=None, max_triggers_in_batch=None)[source]

Bases: object

Calculate SNR and signal consistency tests in a batched progression

combine_results(results)[source]

Combine results from different batches of filtering

process_all()[source]

Process every batch group and return as single result

process_data(data_reader)[source]

Process the data for all of the templates

set_data(data)[source]

Set the data reader object to use

class pycbc.filter.matchedfilter.MatchedFilterControl(low_frequency_cutoff, high_frequency_cutoff, snr_threshold, tlen, delta_f, dtype, segment_list, template_output, use_cluster, downsample_factor=1, upsample_threshold=1, upsample_method='pruned_fft', gpu_callback_method='none', cluster_function='symmetric')[source]

Bases: object

full_matched_filter_and_cluster_fc(segnum, template_norm, window, epoch=None)[source]

Returns the complex snr timeseries, normalization of the complex snr, the correlation vector frequency series, the list of indices of the triggers, and the snr values at the trigger locations. Returns empty lists for these for points that are not above the threshold.

Calculated the matched filter, threshold, and cluster.

Parameters:
  • segnum (int) – Index into the list of segments at MatchedFilterControl construction against which to filter.

  • template_norm (float) – The htilde, template normalization factor.

  • window (int) – Size of the window over which to cluster triggers, in samples

Returns:

  • snr (TimeSeries) – A time series containing the complex snr.

  • norm (float) – The normalization of the complex snr.

  • correlation (FrequencySeries) – A frequency series containing the correlation vector.

  • idx (Array) – List of indices of the triggers.

  • snrv (Array) – The snr values at the trigger locations.

full_matched_filter_and_cluster_symm(segnum, template_norm, window, epoch=None)[source]

Returns the complex snr timeseries, normalization of the complex snr, the correlation vector frequency series, the list of indices of the triggers, and the snr values at the trigger locations. Returns empty lists for these for points that are not above the threshold.

Calculated the matched filter, threshold, and cluster.

Parameters:
  • segnum (int) – Index into the list of segments at MatchedFilterControl construction against which to filter.

  • template_norm (float) – The htilde, template normalization factor.

  • window (int) – Size of the window over which to cluster triggers, in samples

Returns:

  • snr (TimeSeries) – A time series containing the complex snr.

  • norm (float) – The normalization of the complex snr.

  • correlation (FrequencySeries) – A frequency series containing the correlation vector.

  • idx (Array) – List of indices of the triggers.

  • snrv (Array) – The snr values at the trigger locations.

full_matched_filter_thresh_only(segnum, template_norm, window=None, epoch=None)[source]

Returns the complex snr timeseries, normalization of the complex snr, the correlation vector frequency series, the list of indices of the triggers, and the snr values at the trigger locations. Returns empty lists for these for points that are not above the threshold.

Calculated the matched filter, threshold, and cluster.

Parameters:
  • segnum (int) – Index into the list of segments at MatchedFilterControl construction against which to filter.

  • template_norm (float) – The htilde, template normalization factor.

  • window (int) – Size of the window over which to cluster triggers, in samples. This is IGNORED by this function, and provided only for API compatibility.

Returns:

  • snr (TimeSeries) – A time series containing the complex snr.

  • norm (float) – The normalization of the complex snr.

  • correlation (FrequencySeries) – A frequency series containing the correlation vector.

  • idx (Array) – List of indices of the triggers.

  • snrv (Array) – The snr values at the trigger locations.

hierarchical_matched_filter_and_cluster(segnum, template_norm, window)[source]

Returns the complex snr timeseries, normalization of the complex snr, the correlation vector frequency series, the list of indices of the triggers, and the snr values at the trigger locations. Returns empty lists for these for points that are not above the threshold.

Calculated the matched filter, threshold, and cluster.

Parameters:
  • segnum (int) – Index into the list of segments at MatchedFilterControl construction

  • template_norm (float) – The htilde, template normalization factor.

  • window (int) – Size of the window over which to cluster triggers, in samples

Returns:

  • snr (TimeSeries) – A time series containing the complex snr at the reduced sample rate.

  • norm (float) – The normalization of the complex snr.

  • correlation (FrequencySeries) – A frequency series containing the correlation vector.

  • idx (Array) – List of indices of the triggers.

  • snrv (Array) – The snr values at the trigger locations.

class pycbc.filter.matchedfilter.MatchedFilterSkyMaxControl(low_frequency_cutoff, high_frequency_cutoff, snr_threshold, tlen, delta_f, dtype)[source]

Bases: object

full_matched_filter_and_cluster(hplus, hcross, hplus_norm, hcross_norm, psd, stilde, window)[source]

Return the complex snr and normalization.

Calculated the matched filter, threshold, and cluster.

Parameters:
  • h_quantities (Various) – FILL ME IN

  • stilde (FrequencySeries) – The strain data to be filtered.

  • window (int) – The size of the cluster window in samples.

Returns:

  • snr (TimeSeries) – A time series containing the complex snr.

  • norm (float) – The normalization of the complex snr.

  • correlation (FrequencySeries) – A frequency series containing the correlation vector.

  • idx (Array) – List of indices of the triggers.

  • snrv (Array) – The snr values at the trigger locations.

class pycbc.filter.matchedfilter.MatchedFilterSkyMaxControlNoPhase(low_frequency_cutoff, high_frequency_cutoff, snr_threshold, tlen, delta_f, dtype)[source]

Bases: MatchedFilterSkyMaxControl

pycbc.filter.matchedfilter.compute_followup_snr_series(data_reader, htilde, trig_time, duration=0.095, check_state=True, coinc_window=0.05)[source]

Given a StrainBuffer, a template frequency series and a trigger time, compute a portion of the SNR time series centered on the trigger for its rapid sky localization and followup.

If the trigger time is too close to the boundary of the valid data segment the SNR series is calculated anyway and might be slightly contaminated by filter and wrap-around effects. For reasonable durations this will only affect a small fraction of the triggers and probably in a negligible way.

Parameters:
  • data_reader (StrainBuffer) – The StrainBuffer object to read strain data from.

  • htilde (FrequencySeries) – The frequency series containing the template waveform.

  • trig_time ({float, lal.LIGOTimeGPS}) – The trigger time.

  • duration (float (optional)) – Duration of the computed SNR series in seconds. If omitted, it defaults to twice the Earth light travel time plus 10 ms of timing uncertainty.

  • check_state (boolean) – If True, and the detector was offline or flagged for bad data quality at any point during the inspiral, then return (None, None) instead.

  • coinc_window (float (optional)) – Maximum possible time between coincident triggers at different detectors. This is needed to properly determine data padding.

Returns:

snr – The portion of SNR around the trigger. None if the detector is offline or has bad data quality, and check_state is True.

Return type:

TimeSeries

pycbc.filter.matchedfilter.compute_max_snr_over_sky_loc_stat(hplus, hcross, hphccorr, hpnorm=None, hcnorm=None, out=None, thresh=0, analyse_slice=None)[source]

Matched filter maximised over polarization and orbital phase.

This implements the statistic derived in 1603.02444. It is encouraged to read that work to understand the limitations and assumptions implicit in this statistic before using it.

Parameters:
  • hplus (TimeSeries) – This is the IFFTed complex SNR time series of (h+, data). If not normalized, supply the normalization factor so this can be done! It is recommended to normalize this before sending through this function

  • hcross (TimeSeries) – This is the IFFTed complex SNR time series of (hx, data). If not normalized, supply the normalization factor so this can be done!

  • hphccorr (float) – The real component of the overlap between the two polarizations Re[(h+, hx)]. Note that the imaginary component does not enter the detection statistic. This must be normalized and is sign-sensitive.

  • thresh (float) – Used for optimization. If we do not care about the value of SNR values below thresh we can calculate a quick statistic that will always overestimate SNR and then only calculate the proper, more expensive, statistic at points where the quick SNR is above thresh.

  • hpsigmasq (float) – The normalization factor (h+, h+). Default = None (=1, already normalized)

  • hcsigmasq (float) – The normalization factor (hx, hx). Default = None (=1, already normalized)

  • out (TimeSeries (optional, default=None)) – If given, use this array to store the output.

Returns:

det_stat – The SNR maximized over sky location

Return type:

TimeSeries

pycbc.filter.matchedfilter.compute_max_snr_over_sky_loc_stat_no_phase(hplus, hcross, hphccorr, hpnorm=None, hcnorm=None, out=None, thresh=0, analyse_slice=None)[source]

Matched filter maximised over polarization phase.

This implements the statistic derived in 1709.09181. It is encouraged to read that work to understand the limitations and assumptions implicit in this statistic before using it.

In contrast to compute_max_snr_over_sky_loc_stat this function performs no maximization over orbital phase, treating that as an intrinsic parameter. In the case of aligned-spin 2,2-mode only waveforms, this collapses to the normal statistic (at twice the computational cost!)

Parameters:
  • hplus (TimeSeries) – This is the IFFTed complex SNR time series of (h+, data). If not normalized, supply the normalization factor so this can be done! It is recommended to normalize this before sending through this function

  • hcross (TimeSeries) – This is the IFFTed complex SNR time series of (hx, data). If not normalized, supply the normalization factor so this can be done!

  • hphccorr (float) – The real component of the overlap between the two polarizations Re[(h+, hx)]. Note that the imaginary component does not enter the detection statistic. This must be normalized and is sign-sensitive.

  • thresh (float) – Used for optimization. If we do not care about the value of SNR values below thresh we can calculate a quick statistic that will always overestimate SNR and then only calculate the proper, more expensive, statistic at points where the quick SNR is above thresh.

  • hpsigmasq (float) – The normalization factor (h+, h+). Default = None (=1, already normalized)

  • hcsigmasq (float) – The normalization factor (hx, hx). Default = None (=1, already normalized)

  • out (TimeSeries (optional, default=None)) – If given, use this array to store the output.

Returns:

det_stat – The SNR maximized over sky location

Return type:

TimeSeries

pycbc.filter.matchedfilter.compute_u_val_for_sky_loc_stat(hplus, hcross, hphccorr, hpnorm=None, hcnorm=None, indices=None)[source]

The max-over-sky location detection statistic maximizes over a phase, an amplitude and the ratio of F+ and Fx, encoded in a variable called u. Here we return the value of u for the given indices.

pycbc.filter.matchedfilter.compute_u_val_for_sky_loc_stat_no_phase(hplus, hcross, hphccorr, hpnorm=None, hcnorm=None, indices=None)[source]

The max-over-sky location (no phase) detection statistic maximizes over an amplitude and the ratio of F+ and Fx, encoded in a variable called u. Here we return the value of u for the given indices.

pycbc.filter.matchedfilter.correlate(x, y, z)[source]
pycbc.filter.matchedfilter.followup_event_significance(ifo, data_reader, bank, template_id, coinc_times, coinc_threshold=0.005, lookback=150, duration=0.095)[source]

Given a detector, a template waveform and a set of candidate event times in different detectors, perform an on-source/off-source analysis to determine if the SNR in the first detector has a significant peak in the on-source window. The significance is given in terms of a p-value. See Dal Canton et al. 2021 (https://arxiv.org/abs/2008.07494) for details. A portion of the SNR time series around the on-source window is also returned for use in BAYESTAR.

If the calculation cannot be carried out, for example because ifo is not in observing mode at the requested time, then None is returned. Otherwise, the dict contains the following keys. snr_series is a TimeSeries object with the SNR time series for BAYESTAR. peak_time is the time of maximum SNR in the on-source window. pvalue is the p-value for the maximum on-source SNR compared to the off-source realizations. pvalue_saturated is a bool indicating whether the p-value is limited by the number of off-source realizations, i.e. whether the maximum on-source SNR is larger than all the off-source ones. sigma2 is the SNR normalization (squared) for the given template and detector.

Parameters:
  • ifo (str) – Which detector is being used for the calculation.

  • data_reader (StrainBuffer) – StrainBuffer object providing the data for the given detector.

  • bank (LiveFilterBank) – Template bank object providing the template related quantities.

  • template_id (int) – Index of the template in the bank.

  • coinc_times (dict) – Dictionary keyed by detector names reporting the coalescence times of a candidate measured at the different detectors. Used to define the on-source window of the candidate in ifo.

  • coinc_threshold (float) – Nominal statistical uncertainty in coinc_times; expands the on-source window by twice the given amount.

  • lookback (float) – Nominal amount of time to use for the calculation of the onsource and offsource SNR time series. The actual time may be reduced depending on the duration of the template and the strain buffer in the data reader (if so, a warning is logged).

  • duration (float) – Duration of the SNR time series to be reported to BAYESTAR.

Returns:

followup_info – Results of the followup calculation (see above) or None if ifo did not have usable data.

Return type:

dict or None

pycbc.filter.matchedfilter.get_cutoff_indices(flow, fhigh, df, N)[source]

Gets the indices of a frequency series at which to stop an overlap calculation.

Parameters:
  • flow (float) – The frequency (in Hz) of the lower index.

  • fhigh (float) – The frequency (in Hz) of the upper index.

  • df (float) – The frequency step (in Hz) of the frequency series.

  • N (int) – The number of points in the time series. Can be odd or even.

Returns:

  • kmin (int)

  • kmax (int)

pycbc.filter.matchedfilter.make_frequency_series(vec)[source]

Return a frequency series of the input vector.

If the input is a frequency series it is returned, else if the input vector is a real time series it is fourier transformed and returned as a frequency series.

Parameters:

vector (TimeSeries or FrequencySeries) –

Returns:

Frequency Series – A frequency domain version of the input vector.

Return type:

FrequencySeries

pycbc.filter.matchedfilter.match(vec1, vec2, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, v1_norm=None, v2_norm=None, subsample_interpolation=False, return_phase=False)[source]

Return the match between the two TimeSeries or FrequencySeries.

Return the match between two waveforms. This is equivalent to the overlap maximized over time and phase.

The maximization is only performed with discrete time-shifts, or a quadratic interpolation of them if the subsample_interpolation option is turned on; for a more precise computation of the match between two waveforms, use the optimized_match function. The accuracy of this function is guaranteed up to the fourth decimal place.

Parameters:
  • vec1 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • vec2 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • psd (Frequency Series) – A power spectral density to weight the overlap.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin the match.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop the match.

  • v1_norm ({None, float}, optional) – The normalization of the first waveform. This is equivalent to its sigmasq value. If None, it is internally calculated.

  • v2_norm ({None, float}, optional) – The normalization of the second waveform. This is equivalent to its sigmasq value. If None, it is internally calculated.

  • subsample_interpolation ({False, bool}, optional) – If True the peak will be interpolated between samples using a simple quadratic fit. This can be important if measuring matches very close to 1 and can cause discontinuities if you don’t use it as matches move between discrete samples. If True the index returned will be a float.

  • return_phase ({False, bool}, optional) – If True, also return the phase shift that gives the match.

Returns:

  • match (float)

  • index (int) – The number of samples to shift to get the match.

  • phi (float) – Phase to rotate complex waveform to get the match, if desired.

pycbc.filter.matchedfilter.matched_filter(template, data, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, sigmasq=None)[source]

Return the complex snr.

Return the complex snr, along with its associated normalization of the template, matched filtered against the data.

Parameters:
  • template (TimeSeries or FrequencySeries) – The template waveform

  • data (TimeSeries or FrequencySeries) – The strain data to be filtered.

  • psd (FrequencySeries) – The noise weighting of the filter.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin the filter calculation. If None, begin at the first frequency after DC.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop the filter calculation. If None, continue to the the nyquist frequency.

  • sigmasq ({None, float}, optional) – The template normalization. If none, this value is calculated internally.

Returns:

snr – A time series containing the complex snr.

Return type:

TimeSeries

pycbc.filter.matchedfilter.matched_filter_core(template, data, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, h_norm=None, out=None, corr_out=None)[source]

Return the complex snr and normalization.

Return the complex snr, along with its associated normalization of the template, matched filtered against the data.

Parameters:
  • template (TimeSeries or FrequencySeries) – The template waveform

  • data (TimeSeries or FrequencySeries) – The strain data to be filtered.

  • psd ({FrequencySeries}, optional) – The noise weighting of the filter.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin the filter calculation. If None, begin at the first frequency after DC.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop the filter calculation. If None, continue to the the nyquist frequency.

  • h_norm ({None, float}, optional) – The template normalization. If none, this value is calculated internally.

  • out ({None, Array}, optional) – An array to use as memory for snr storage. If None, memory is allocated internally.

  • corr_out ({None, Array}, optional) – An array to use as memory for correlation storage. If None, memory is allocated internally. If provided, management of the vector is handled externally by the caller. No zero’ing is done internally.

Returns:

  • snr (TimeSeries) – A time series containing the complex snr.

  • correlation (FrequencySeries) – A frequency series containing the correlation vector.

  • norm (float) – The normalization of the complex snr.

pycbc.filter.matchedfilter.optimized_match(vec1, vec2, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, v1_norm=None, v2_norm=None, return_phase=False)[source]

Given two waveforms (as numpy arrays), compute the optimized match between them, making use of scipy.minimize_scalar.

This function computes the same quantities as “match”; it is more accurate and slower.

Parameters:
  • vec1 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • vec2 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • psd (FrequencySeries) – A power spectral density to weight the overlap.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin the match.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop the match.

  • v1_norm ({None, float}, optional) – The normalization of the first waveform. This is equivalent to its sigmasq value. If None, it is internally calculated.

  • v2_norm ({None, float}, optional) – The normalization of the second waveform. This is equivalent to its sigmasq value. If None, it is internally calculated.

  • return_phase ({False, bool}, optional) – If True, also return the phase shift that gives the match.

Returns:

  • match (float)

  • index (int) – The number of samples to shift to get the match.

  • phi (float) – Phase to rotate complex waveform to get the match, if desired.

pycbc.filter.matchedfilter.overlap(vec1, vec2, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, normalized=True)[source]

Return the overlap between the two TimeSeries or FrequencySeries.

Parameters:
  • vec1 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • vec2 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • psd (Frequency Series) – A power spectral density to weight the overlap.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin the overlap.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop the overlap.

  • normalized ({True, boolean}, optional) – Set if the overlap is normalized. If true, it will range from 0 to 1.

Returns:

overlap

Return type:

float

pycbc.filter.matchedfilter.overlap_cplx(vec1, vec2, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None, normalized=True)[source]

Return the complex overlap between the two TimeSeries or FrequencySeries.

Parameters:
  • vec1 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • vec2 (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • psd (Frequency Series) – A power spectral density to weight the overlap.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin the overlap.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop the overlap.

  • normalized ({True, boolean}, optional) – Set if the overlap is normalized. If true, it will range from 0 to 1.

Returns:

overlap

Return type:

complex

pycbc.filter.matchedfilter.sigma(htilde, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None)[source]

Return the sigma of the waveform. See sigmasq for more details.

Parameters:
  • htilde (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • psd ({None, FrequencySeries}, optional) – The psd used to weight the accumulated power.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin considering waveform power.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop considering waveform power.

Returns:

sigmasq

Return type:

float

pycbc.filter.matchedfilter.sigmasq(htilde, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None)[source]

Return the loudness of the waveform. This is defined (see Duncan Brown’s thesis) as the unnormalized matched-filter of the input waveform, htilde, with itself. This quantity is usually referred to as (sigma)^2 and is then used to normalize matched-filters with the data.

Parameters:
  • htilde (TimeSeries or FrequencySeries) – The input vector containing a waveform.

  • psd ({None, FrequencySeries}, optional) – The psd used to weight the accumulated power.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin considering waveform power.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop considering waveform power.

Returns:

sigmasq

Return type:

float

pycbc.filter.matchedfilter.sigmasq_series(htilde, psd=None, low_frequency_cutoff=None, high_frequency_cutoff=None)[source]

Return a cumulative sigmasq frequency series.

Return a frequency series containing the accumulated power in the input up to that frequency.

Parameters:
  • htilde (TimeSeries or FrequencySeries) – The input vector

  • psd ({None, FrequencySeries}, optional) – The psd used to weight the accumulated power.

  • low_frequency_cutoff ({None, float}, optional) – The frequency to begin accumulating power. If None, start at the beginning of the vector.

  • high_frequency_cutoff ({None, float}, optional) – The frequency to stop considering accumulated power. If None, continue until the end of the input vector.

Returns:

Frequency Series – A frequency series containing the cumulative sigmasq.

Return type:

FrequencySeries

pycbc.filter.matchedfilter_cpu module

class pycbc.filter.matchedfilter_cpu.CPUCorrelator(x, y, z)

Bases: _BaseCorrelator

correlate(self)
pycbc.filter.matchedfilter_cpu.batch_correlate_execute(self, y)
pycbc.filter.matchedfilter_cpu.correlate(x, y, z)
pycbc.filter.matchedfilter_cpu.correlate_numpy(x, y, z)

pycbc.filter.matchedfilter_numpy module

pycbc.filter.matchedfilter_numpy.correlate(x, y, z)[source]

pycbc.filter.qtransform module

This module retrives a timeseries and then calculates the q-transform of that time series

pycbc.filter.qtransform.deltam_f(mismatch)[source]

Fractional mismatch between neighbouring tiles

Parameters:

mismatch ('float') – percentage of desired fractional mismatch

Return type:

type: ‘float’

pycbc.filter.qtransform.qplane(qplane_tile_dict, fseries, return_complex=False)[source]
Performs q-transform on each tile for each q-plane and selects

tile with the maximum energy. Q-transform can then be interpolated to a desired frequency and time resolution.

Parameters:
  • qplane_tile_dict – Dictionary containing a list of q-tile tupples for each q-plane

  • fseries ('pycbc FrequencySeries') – frequency-series data set

  • return_complex ({False, bool}) – Return the raw complex series instead of the normalized power.

Returns:

  • q (float) – The q of the maximum q plane

  • times (numpy.ndarray) – The time that the qtransform is sampled.

  • freqs (numpy.ndarray) – The frequencies that the qtransform is samled.

  • qplane (numpy.ndarray (2d)) – The two dimensional interpolated qtransform of this time series.

pycbc.filter.qtransform.qseries(fseries, Q, f0, return_complex=False)[source]

Calculate the energy ‘TimeSeries’ for the given fseries

Parameters:
  • fseries ('pycbc FrequencySeries') – frequency-series data set

  • Q – q value

  • f0 – central frequency

  • return_complex ({False, bool}) – Return the raw complex series instead of the normalized power.

Returns:

energy – A ‘TimeSeries’ of the normalized energy from the Q-transform of this tile against the data.

Return type:

‘~pycbc.types.TimeSeries’

pycbc.filter.qtransform.qtiling(fseries, qrange, frange, mismatch=0.2)[source]

Iterable constructor of QTile tuples

Parameters:
  • fseries ('pycbc FrequencySeries') – frequency-series data set

  • qrange – upper and lower bounds of q range

  • frange – upper and lower bounds of frequency range

  • mismatch – percentage of desired fractional mismatch

Returns:

qplane_tile_dict – dictionary containing Q-tile tuples for a set of Q-planes

Return type:

‘dict’

pycbc.filter.resample module

pycbc.filter.resample.fir_zero_filter(coeff, timeseries)[source]

Filter the timeseries with a set of FIR coefficients

Parameters:
  • coeff (numpy.ndarray) – FIR coefficients. Should be and odd length and symmetric.

  • timeseries (pycbc.types.TimeSeries) – Time series to be filtered.

Returns:

  • filtered_series (pycbc.types.TimeSeries) – Return the filtered timeseries, which has been properly shifted to account

  • for the FIR filter delay and the corrupted regions zeroed out.

pycbc.filter.resample.highpass(timeseries, frequency, filter_order=8, attenuation=0.1)[source]

Return a new timeseries that is highpassed.

Return a new time series that is highpassed above the frequency.

Parameters:
  • Series (Time) – The time series to be high-passed.

  • frequency (float) – The frequency below which is suppressed.

  • filter_order ({8, int}, optional) – The order of the filter to use when high-passing the time series.

  • attenuation ({0.1, float}, optional) – The attenuation of the filter.

Returns:

Time Series – A new TimeSeries that has been high-passed.

Return type:

TimeSeries

Raises:
  • TypeError: – time_series is not an instance of TimeSeries.

  • TypeError: – time_series is not real valued

pycbc.filter.resample.highpass_fir(timeseries, frequency, order, beta=5.0)[source]

Highpass filter the time series using an FIR filtered generated from the ideal response passed through a kaiser window (beta = 5.0)

Parameters:
  • Series (Time) – The time series to be high-passed.

  • frequency (float) – The frequency below which is suppressed.

  • order (int) – Number of corrupted samples on each side of the time series

  • beta (float) – Beta parameter of the kaiser window that sets the side lobe attenuation.

pycbc.filter.resample.interpolate_complex_frequency(series, delta_f, zeros_offset=0, side='right')[source]

Interpolate complex frequency series to desired delta_f.

Return a new complex frequency series that has been interpolated to the desired delta_f.

Parameters:
  • series (FrequencySeries) – Frequency series to be interpolated.

  • delta_f (float) – The desired delta_f of the output

  • zeros_offset (optional, {0, int}) – Number of sample to delay the start of the zero padding

  • side (optional, {'right', str}) – The side of the vector to zero pad

Returns:

interpolated series – A new FrequencySeries that has been interpolated.

Return type:

FrequencySeries

pycbc.filter.resample.lowpass(timeseries, frequency, filter_order=8, attenuation=0.1)[source]

Return a new timeseries that is lowpassed.

Return a new time series that is lowpassed below the frequency.

Parameters:
  • Series (Time) – The time series to be low-passed.

  • frequency (float) – The frequency above which is suppressed.

  • filter_order ({8, int}, optional) – The order of the filter to use when low-passing the time series.

  • attenuation ({0.1, float}, optional) – The attenuation of the filter.

Returns:

Time Series – A new TimeSeries that has been low-passed.

Return type:

TimeSeries

Raises:
  • TypeError: – time_series is not an instance of TimeSeries.

  • TypeError: – time_series is not real valued

pycbc.filter.resample.lowpass_fir(timeseries, frequency, order, beta=5.0)[source]

Lowpass filter the time series using an FIR filtered generated from the ideal response passed through a kaiser window (beta = 5.0)

Parameters:
  • Series (Time) – The time series to be low-passed.

  • frequency (float) – The frequency below which is suppressed.

  • order (int) – Number of corrupted samples on each side of the time series

  • beta (float) – Beta parameter of the kaiser window that sets the side lobe attenuation.

pycbc.filter.resample.notch_fir(timeseries, f1, f2, order, beta=5.0)[source]

notch filter the time series using an FIR filtered generated from the ideal response passed through a time-domain kaiser window (beta = 5.0)

The suppression of the notch filter is related to the bandwidth and the number of samples in the filter length. For a few Hz bandwidth, a length corresponding to a few seconds is typically required to create significant suppression in the notched band. To achieve frequency resolution df at sampling frequency fs, order should be at least fs/df.

Parameters:
  • Series (Time) – The time series to be notched.

  • f1 (float) – The start of the frequency suppression.

  • f2 (float) – The end of the frequency suppression.

  • order (int) – Number of corrupted samples on each side of the time series (Extent of the filter on either side of zero)

  • beta (float) – Beta parameter of the kaiser window that sets the side lobe attenuation.

pycbc.filter.resample.resample_to_delta_t(timeseries, delta_t, method='butterworth')[source]

Resmple the time_series to delta_t

Resamples the TimeSeries instance time_series to the given time step, delta_t. Only powers of two and real valued time series are supported at this time. Additional restrictions may apply to particular filter methods.

Parameters:
  • time_series (TimeSeries) – The time series to be resampled

  • delta_t (float) – The desired time step

Returns:

Time Series – A TimeSeries that has been resampled to delta_t.

Return type:

TimeSeries

Raises:
  • TypeError: – time_series is not an instance of TimeSeries.

  • TypeError: – time_series is not real valued

Examples

>>> h_plus_sampled = resample_to_delta_t(h_plus, 1.0/2048)

pycbc.filter.simd_correlate module

pycbc.filter.simd_correlate.correlate_parallel(ht, st, qt)[source]
pycbc.filter.simd_correlate.correlate_simd(ht, st, qt)[source]

pycbc.filter.simd_correlate_cython module

pycbc.filter.simd_correlate_cython.ccorrf_parallel(inconj, innoconj, out, arrlen, segsize)
pycbc.filter.simd_correlate_cython.ccorrf_simd(inconj, innoconj, out, len)

pycbc.filter.zpk module

pycbc.filter.zpk.filter_zpk(timeseries, z, p, k)[source]

Return a new timeseries that was filtered with a zero-pole-gain filter. The transfer function in the s-domain looks like: .. math:: frac{H(s) = (s - s_1) * (s - s_3) * … * (s - s_n)}{(s - s_2) * (s - s_4) * … * (s - s_m)}, m >= n

The zeroes, and poles entered in Hz are converted to angular frequency, along the imaginary axis in the s-domain s=i*omega. Then the zeroes, and poles are bilinearly transformed via: .. math:: z(s) = frac{(1 + s*T/2)}{(1 - s*T/2)}

Where z is the z-domain value, s is the s-domain value, and T is the sampling period. After the poles and zeroes have been bilinearly transformed, then the second-order sections are found and filter the data using scipy.

Parameters:
  • timeseries (TimeSeries) – The TimeSeries instance to be filtered.

  • z (array) – Array of zeros to include in zero-pole-gain filter design. In units of Hz.

  • p (array) – Array of poles to include in zero-pole-gain filter design. In units of Hz.

  • k (float) – Gain to include in zero-pole-gain filter design. This gain is a constant multiplied to the transfer function.

Returns:

Time Series – A new TimeSeries that has been filtered.

Return type:

TimeSeries

Examples

To apply a 5 zeroes at 100Hz, 5 poles at 1Hz, and a gain of 1e-10 filter to a TimeSeries instance, do: >>> filtered_data = zpk_filter(timeseries, [100]*5, [1]*5, 1e-10)

Module contents