pygwb.baseline.Baseline

class pygwb.baseline.Baseline(name, interferometer_1, interferometer_2, duration=None, frequencies=None, calibration_epsilon=0, notch_list_path='', coarse_grain_psd=False, coarse_grain_csd=True, overlap_factor_welch=0.5, overlap_factor=0.5, window_fftgram_dict={'window_fftgram': 'hann'}, window_fftgram_dict_welch={'window_fftgram': 'hann'}, N_average_segments_psd=2, sampling_frequency=None)[source]

Bases: object

__init__(name, interferometer_1, interferometer_2, duration=None, frequencies=None, calibration_epsilon=0, notch_list_path='', coarse_grain_psd=False, coarse_grain_csd=True, overlap_factor_welch=0.5, overlap_factor=0.5, window_fftgram_dict={'window_fftgram': 'hann'}, window_fftgram_dict_welch={'window_fftgram': 'hann'}, N_average_segments_psd=2, sampling_frequency=None)[source]

Instantiate a Baseline.

Parameters:
  • name (str) – Name for the baseline, e.g H1H2

  • interferometer_1/2 (bilby.Interferometer object) – The two detectors spanning the baseline.

  • duration (float, optional) – The duration in seconds of each data segment in the interferometers. None by default, in which case duration is inherited from the interferometers.

  • frequencies (array_like, optional) – The frequency array for the Baseline and interferometers.

  • calibration_epsilon (float, optional) – Calibration uncertainty for this baseline.

  • notch_list_path (str, optional) – File path of the baseline notch list.

  • coarse_grain_psd (bool, optional) – Whether to apply coarse graining to obtain PSD spectra. Default is False.

  • coarse_grain_csd (bool, optional) – Whether to apply coarse graining to obtain CSD spectra. Default is True.

  • overlap_factor_welch (float, optional) – Overlap factor to use when if using Welch’s method to estimate spectra (NOT coarsegraining). For “hann” window use 0.5 overlap_factor and for “boxcar”window use 0 overlap_factor. Default is 0.5 (50% overlap), which is optimal when using Welch’s method with a “hann” window.

  • overlap_factor (float, optional) – Factor by which to overlap the segments in the psd and csd estimation. Default is 1/2, if set to 0 no overlap is performed.

  • window_fftgram_dict (dict, optional) – Dictionary containing name and parameters describing which window to use when producing fftgrams for csds (and psds if these are coarse-grained). Default is “hann”.

  • window_fftgram_dict_welch (dict, optional) – Dictionary containing name and parameters describing which window to use when producing fftgrams with Welch’s method. Default is “hann”.

  • N_average_segments_psd (int, optional) – Number of segments used for PSD averaging (from both sides of the segment of interest) N_avg_segs should be even and >= 2.

__call__(*args, **kwargs)

Call self as a function.

Methods

__init__(name, interferometer_1, ...[, ...])

Instantiate a Baseline.

calc_baseline_orf([polarization, frequencies])

Calculate the overlap reduction function for this baseline.

calculate_delta_sigma_cut(delta_sigma_cut, ...)

Calculate the delta sigma cut using the naive and average psds, if set in the baseline.

check_sampling_frequencies_match_baseline_ifos(...)

Check that the sampling frequency of the two interferometers in this Baseline match the Baseline sampling frequency.

crop_frequencies_average_psd_csd(flow, fhigh)

Crop frequencies of average PSDs and CSDS.

from_interferometers(interferometers[, ...])

Load a Baseline from a list of interferometers.

from_parameters(interferometer_1, ...[, ...])

Load a Baseline from a Parameters object.

load_from_pickle(filename)

Load baseline object from pickle file.

reweight([new_alpha, new_fref])

Reweight all the frequency-weighted attributes of this Baseline, if these are set.

save_point_estimate_spectra(save_data_type, ...)

Save the overall point estimate Y, its error bar sigma, the frequency-dependent estimates and variances and the corresponding frequencies in the required save_data_type, which can be npz, pickle, json or hdf5.

save_psds_csds(save_data_type, filename)

Save the average and naive psds and csds and the corresponding frequencies in the required save_data_type, which can be npz, pickle, json or hdf5.

save_to_pickle(filename[, wipe])

Save baseline object to pickle file.

set_average_cross_spectral_density()

If csd has been calculated, sets the average csd for the baseline.

set_average_power_spectral_densities()

If psds have been calculated, sets the average psd in each ifo.

set_coherence_spectrum([flow, fhigh, ...])

Set the coherence spectrum between detectors, averaged over all data in the baseline.

set_cross_and_power_spectral_density(...)

Set the power spectral density in each interferometer and the cross spectral density for the baseline object when data are available.

set_frequency_mask([notch_list_path, ...])

Set frequency mask to frequencies attribute.

set_point_estimate_sigma([badtimes, alpha, ...])

Set point estimate sigma based on a set of parameters.

set_point_estimate_sigma_spectrogram([...])

Set point estimate and sigma spectrogram.

set_point_estimate_sigma_spectrum([...])

Set time-integrated point estimate spectrum and variance in each frequency bin.

Attributes

badGPStimes

GPS times flagged by delta sigma cut.

coherence_dict

Coherence dictionary based on data in this baseline.

coherence_spectrum

Coherence spectrum calculated using data in this baseline.

csd_segment_offset

CSD segment offset to use for this baseline.

delta_sigmas

Values of delta sigmas for data segments in the baseline.

duration

Duration in seconds of a unit segment of data stored in the baseline detectors.

frequencies

Frequency array associated to this baseline.

gamma_v

Overlap reduction function for asymmetrically polarised backgrounds, as described in https://arxiv.org/pdf/0707.0535.pdf

orf_polarization

Overlap reduction function polarization

overlap_reduction_function

Overlap reduction function associated to this baseline, calculated for the requested polarization.

point_estimate

Point estimate (in Omega units) calculated using data in this baseline.

point_estimate_spectrogram

Point estimate spectrogram (in Omega units) calculated using data in this baseline.

point_estimate_spectrum

Point estimate spectrum (in Omega units) calculated using data in this baseline.

sampling_frequency

Sampling frequency of the data stored in this baseline.

scalar_overlap_reduction_function

Overlap reduction function calculated for scalar polarization.

sigma

Sigma (in Omega units) calculated using data in this baseline.

sigma_spectrogram

Sigma spectrogram (in Omega units) calculated using data in this baseline.

sigma_spectrum

Sigma spectrum (in Omega units) calculated using data in this baseline.

tensor_overlap_reduction_function

Overlap reduction function calculated for tensor polarization.

vector_overlap_reduction_function

Overlap reduction function calculated for vector polarization.

property badGPStimes

GPS times flagged by delta sigma cut.

calc_baseline_orf(polarization='tensor', frequencies=None)[source]

Calculate the overlap reduction function for this baseline. Wraps the ORF module.

Parameters:
  • polarization (str, optional) – Polarization of the signal to consider (scalar, vector, tensor) for the ORF calculation. Default is tensor.

  • frequencies (array_like, optional) – Frequency array to use in the calculation of the ORF. By default, self.frequencies is used.

Returns:
orf: array_like

Overlap reduction function for the required polarization.

See also

pygwb.orfs.calc_orf

Method to compute the overlap reduction function.

calculate_delta_sigma_cut(delta_sigma_cut, alphas, fref, flow=20, fhigh=1726, notch_list_path='', polarization='tensor', return_naive_and_averaged_sigmas: bool = False)[source]

Calculate the delta sigma cut using the naive and average psds, if set in the baseline.

Parameters:
  • delta_sigma_cut (float) – The cutoff to implement in the delta sigma cut.

  • alphas (list) – Set of spectral indices to use in the delta sigma cut calculation.

  • flow (float, optional) – Low frequency. Default is 20 Hz.

  • fhigh (float, optional) – High frequency. Default is 1726 Hz.

  • notch_list_path (str, optional) – File path of the baseline notch list.

  • fref (int) – Reference frequency (Hz).

  • return_naive_and_averaged_sigmas (bool, optional) – Option to return naive and sliding sigmas. Default is False.

  • polarization (str, optional) – Polarization of the signal to consider (scalar, vector, tensor) for the orf calculation. Default is tensor.

See also

pygwb.delta_sigma_cut.run_dsc

Function used to run the delta sigma cut.

check_sampling_frequencies_match_baseline_ifos(sampling_frequency)[source]

Check that the sampling frequency of the two interferometers in this Baseline match the Baseline sampling frequency.

Parameters:

sampling_frequency (float) – The sampling frequency that is being set for the Baseline.

Notes

If the sampling frequency passed is None, the Baseline sampling frequency will be set to that of the interferometers, if these match. If these don’t match, an error will be raised. If the sampling frequency of the interferometers is also None, then no sampling frequency will be set, and the user can set it at a later time.

property coherence_dict

Coherence dictionary based on data in this baseline.

property coherence_spectrum

Coherence spectrum calculated using data in this baseline.

crop_frequencies_average_psd_csd(flow, fhigh)[source]

Crop frequencies of average PSDs and CSDS. Done in place. This is not completely implemented yet.

Parameters:
  • flow (float) – Low frequency to crop.

  • fhigh (float) – High frequency to crop.

property csd_segment_offset

CSD segment offset to use for this baseline.

property delta_sigmas

Values of delta sigmas for data segments in the baseline.

property duration

Duration in seconds of a unit segment of data stored in the baseline detectors.

property frequencies

Frequency array associated to this baseline.

classmethod from_interferometers(interferometers, duration=None, calibration_epsilon=0)[source]

Load a Baseline from a list of interferometers.

Parameters:
  • interferometers (list) – List of two bilby Interferometer objects.

  • duration (float, optional) – Segment duration in seconds. Default is None.

  • calibration_epsilon (float, optional) – Calibration uncertainty for this baseline. Default is 0.

Returns:
Baseline: cls

Baseline class

classmethod from_parameters(interferometer_1, interferometer_2, parameters, frequencies=None)[source]

Load a Baseline from a Parameters object.

Parameters:
  • interferometer_1/2 (bilby.Interferometer object) – The two detectors spanning this baseline.

  • parameters (pygwb.parameters object) – Parameters object containing necessary parameters for the instantiation of the baseline, and subsequent analyses.

  • frequencies (array_like, optional) – Frequency array to use in the instantiation of this baseline. Default is None.

Returns:
Baseline: cls

Baseline class

property gamma_v

Overlap reduction function for asymmetrically polarised backgrounds, as described in https://arxiv.org/pdf/0707.0535.pdf

classmethod load_from_pickle(filename)[source]

Load baseline object from pickle file.

Parameters:

filename (str) – Filename (inclusive of path) to load the pickled baseline from.

Returns:
Baseline: cls

Baseline class.

property orf_polarization

Overlap reduction function polarization

property overlap_reduction_function

Overlap reduction function associated to this baseline, calculated for the requested polarization.

property point_estimate

Point estimate (in Omega units) calculated using data in this baseline.

property point_estimate_spectrogram

Point estimate spectrogram (in Omega units) calculated using data in this baseline.

property point_estimate_spectrum

Point estimate spectrum (in Omega units) calculated using data in this baseline.

reweight(new_alpha=None, new_fref=None)[source]

Reweight all the frequency-weighted attributes of this Baseline, if these are set.

Parameters:
  • new_alpha (float, optional) – New alpha to weight the spectra to. Default is None.

  • new_fref (float, optional) – New reference frequency to refer the spectra to. Default is None.

property sampling_frequency

Sampling frequency of the data stored in this baseline. This must match the sampling frequency stored in this baseline’s interferometers.

save_point_estimate_spectra(save_data_type, filename)[source]

Save the overall point estimate Y, its error bar sigma, the frequency-dependent estimates and variances and the corresponding frequencies in the required save_data_type, which can be npz, pickle, json or hdf5. You can call upon this data afterwards when loaoding in using the [‘key’] dictionary format.

Parameters:
  • save_data_type (str) – The required type of data file where the information will be stored.

  • filename (str) – The path/name of the file in which you want to save.

save_psds_csds(save_data_type, filename)[source]

Save the average and naive psds and csds and the corresponding frequencies in the required save_data_type, which can be npz, pickle, json or hdf5. One can call upon the data afterwards when loading in using the [‘key’] dictionary format.

Parameters:
  • save_data_type (str) – The required type of data file where the information will be stored.

  • filename (str) – The path/name of the file in which you want to save.

save_to_pickle(filename, wipe=True)[source]

Save baseline object to pickle file.

Parameters:

filename (str) – Filename (inclusive of path) to save the pickled baseline to.

property scalar_overlap_reduction_function

Overlap reduction function calculated for scalar polarization.

set_average_cross_spectral_density()[source]

If csd has been calculated, sets the average csd for the baseline.

set_average_power_spectral_densities()[source]

If psds have been calculated, sets the average psd in each ifo.

set_coherence_spectrum(flow=20, fhigh=1726, badtimes=None, apply_dsc=True)[source]

Set the coherence spectrum between detectors, averaged over all data in the baseline.

Parameters:
  • flow (float, optional) – Low frequency. Default is 20 Hz.

  • fhigh (float, optional) – High frequency. Default is 1726 Hz.

  • badtimes (array_like, optional) – Array of times to exclude from the coherence calculation. Default is None.

  • apply_dsc (bool, optional) – Apply delta sigma cut flag; if True, removes the badGPStimes from the spectra calculations. Default is True.

Notes

The coherence calculation uses averaged naive PSD estimates as the coherence is calculated using CSD and PSD estimates of each individual segment, calculated “on shell”.

set_cross_and_power_spectral_density(frequency_resolution)[source]

Set the power spectral density in each interferometer and the cross spectral density for the baseline object when data are available.

Parameters:

frequency_resolution (float) – The frequency resolution at which the cross and power spectral densities are calculated.

set_frequency_mask(notch_list_path='', apply_notches=True)[source]

Set frequency mask to frequencies attribute.

Parameters:
  • notch_list_path (str, optional) – Path to notch list to apply to frequency array. If not provided, no notching will be applied at this point.

  • apply_notches (bool, optional) – If True, apply frequency notches. Default is True.

See also

pygwb.notch.StochNotchList

Used to read in the frequency notches.

set_point_estimate_sigma(badtimes=None, alpha=0.0, fref=25, flow=20, fhigh=1726, notch_list_path='', polarization='tensor', apply_dsc=True, apply_notches=True)[source]

Set point estimate sigma based on a set of parameters. This is estimate of omega_gw in each frequency bin.

Parameters:
  • badtimes (array_like, optional) – Array of times to exclude from point estimate/sigma calculation. If no times are passed, none will be excluded.

  • alpha (float, optional) – Spectral index to use in the re-weighting. Default is 0.

  • fref (float, optional) – Reference frequency to use in the re-weighting. Default is 25.

  • flow (float, optional) – Low frequency. Default is 20 Hz.

  • fhigh (float, optional) – High frequency. Default is 1726 Hz.

  • notch_list_path (str, optional) – Path to the notch list to use in the spectrum; if the notch_list isn’t set in the baseline, user can pass it directly here. If it is not set and if none is passed no notches will be applied.

  • polarization (str, optional) – Polarization of the signal to consider (scalar, vector, tensor) for the orf calculation. Default is Tensor.

  • apply_dsc (bool, optional) – Apply delta sigma cut flag; if True, removes the badGPStimes from the spectra calculations. Default is True.

  • apply_notches (bool, optional) – Apply spectral notches flag; if True, remove the notches specified in the notch_list from the spectra calculations. Default is True.

set_point_estimate_sigma_spectrogram(alpha=0.0, fref=25, flow=20, fhigh=1726, polarization='tensor')[source]

Set point estimate and sigma spectrogram. Resulting spectrogram does not include frequency weighting for alpha.

Parameters:
  • alpha (float, optional) – Spectral index to use in the weighting. Default is 0.

  • fref (float, optional) – Reference frequency to use in the weighting calculation. Final result refers to this frequency. Default is 25 Hz.

  • flow (float, optional) – Lowest frequency to consider. Default is 20 Hz.

  • fhigh (float, optional) – Highest frequency to consider. Default is 1726 Hz.

  • polarization (str, optional) – Polarization of the signal to consider (scalar, vector, tensor) for the ORF calculation. Default is tensor.

set_point_estimate_sigma_spectrum(badtimes=None, alpha=0.0, fref=25, flow=20, fhigh=1726, notch_list_path='', polarization='tensor', apply_dsc=True, apply_notches=True)[source]

Set time-integrated point estimate spectrum and variance in each frequency bin.

Parameters:
  • badtimes (array_like, optional) – Array of times to exclude from point estimate/sigma calculation. If no times are passed, none will be excluded.

  • alpha (float, optional) – Spectral index to use in the re-weighting. Default is 0.

  • fref (float, optional) – Reference frequency to use in the re-weighting. Default is 25 Hz.

  • flow (float, optional) – Low frequency. Default is 20 Hz.

  • fhigh (float, optional) – High frequency. Default is 1726 Hz.

  • notch_list_path (str, optional) – Path to the notch list to use in the spectrum. Default is empty string.

  • polarization (str, optional) – Polarization of the signal to consider (scalar, vector, tensor) for the orf calculation. Default is tensor.

  • apply_dsc (bool, optional) – Apply delta sigma cut flag; if True, removes the badGPStimes from the spectra calculations. Default is True.

  • apply_notches (bool, optional) – Apply spectral notches flag; if True, remove the notches specified in the notch_list from the spectra calculations. Default is True.

property sigma

Sigma (in Omega units) calculated using data in this baseline.

property sigma_spectrogram

Sigma spectrogram (in Omega units) calculated using data in this baseline.

property sigma_spectrum

Sigma spectrum (in Omega units) calculated using data in this baseline.

property tensor_overlap_reduction_function

Overlap reduction function calculated for tensor polarization.

property vector_overlap_reduction_function

Overlap reduction function calculated for vector polarization.