pycbc.psd package
Submodules
pycbc.psd.analytical module
Provides reference PSDs from LALSimulation and pycbc.psd.analytical_space.
More information about how to use these ground-based detectors’ PSD can be found in the guide about Generating an Analytic PSD from lalsimulation. For space-borne ones, see pycbc.psd.analytical_space module.
- pycbc.psd.analytical.AdVBNSOptimizedSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVBNSOptimizedSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVDesignSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVDesignSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVEarlyHighSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVEarlyHighSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVEarlyLowSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVEarlyLowSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVLateHighSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVLateHighSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVLateLowSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVLateLowSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVMidHighSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVMidHighSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVMidLowSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVMidLowSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.AdVO3LowT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVO3LowT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.AdVO4IntermediateT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVO4IntermediateT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.AdVO4T1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdVO4T1800545 PSD from LALSimulation.
- pycbc.psd.analytical.AdvVirgo(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the AdvVirgo PSD from LALSimulation.
- pycbc.psd.analytical.CosmicExplorerP1600143(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the CosmicExplorerP1600143 PSD from LALSimulation.
- pycbc.psd.analytical.CosmicExplorerPessimisticP1600143(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the CosmicExplorerPessimisticP1600143 PSD from LALSimulation.
- pycbc.psd.analytical.CosmicExplorerWidebandP1600143(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the CosmicExplorerWidebandP1600143 PSD from LALSimulation.
- pycbc.psd.analytical.EinsteinTelescopeP1600143(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the EinsteinTelescopeP1600143 PSD from LALSimulation.
- pycbc.psd.analytical.GEO(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the GEO PSD from LALSimulation.
- pycbc.psd.analytical.GEOHF(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the GEOHF PSD from LALSimulation.
- pycbc.psd.analytical.KAGRA(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRA PSD from LALSimulation.
- pycbc.psd.analytical.KAGRA128MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRA128MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRA25MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRA25MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRA80MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRA80MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRADesignSensitivityT1600593(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRADesignSensitivityT1600593 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRAEarlySensitivityT1600593(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRAEarlySensitivityT1600593 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRALateSensitivityT1600593(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRALateSensitivityT1600593 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRAMidSensitivityT1600593(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRAMidSensitivityT1600593 PSD from LALSimulation.
- pycbc.psd.analytical.KAGRAOpeningSensitivityT1600593(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the KAGRAOpeningSensitivityT1600593 PSD from LALSimulation.
- pycbc.psd.analytical.TAMA(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the TAMA PSD from LALSimulation.
- pycbc.psd.analytical.Virgo(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the Virgo PSD from LALSimulation.
- pycbc.psd.analytical.aLIGO140MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGO140MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGO175MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGO175MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOAPlusDesignSensitivityT1800042(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOAPlusDesignSensitivityT1800042 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOAdVO3LowT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOAdVO3LowT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOAdVO4IntermediateT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOAdVO4IntermediateT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOAdVO4T1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOAdVO4T1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOBHBH20Deg(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOBHBH20Deg PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOBHBH20DegGWINC(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOBHBH20DegGWINC PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOBNSOptimizedSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOBNSOptimizedSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGODesignSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGODesignSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGODesignSensitivityT1800044(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGODesignSensitivityT1800044 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOEarlyHighSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOEarlyHighSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOEarlyLowSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOEarlyLowSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOHighFrequency(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOHighFrequency PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOHighFrequencyGWINC(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOHighFrequencyGWINC PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOKAGRA128MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOKAGRA128MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOKAGRA25MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOKAGRA25MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOKAGRA80MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOKAGRA80MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOLateHighSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOLateHighSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOLateLowSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOLateLowSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOMidHighSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOMidHighSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOMidLowSensitivityP1200087(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOMidLowSensitivityP1200087 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGONSNSOpt(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGONSNSOpt PSD from LALSimulation.
- pycbc.psd.analytical.aLIGONSNSOptGWINC(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGONSNSOptGWINC PSD from LALSimulation.
- pycbc.psd.analytical.aLIGONoSRMHighPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGONoSRMHighPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGONoSRMLowPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGONoSRMLowPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGONoSRMLowPowerGWINC(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGONoSRMLowPowerGWINC PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOO3LowT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOO3LowT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumBHBH20Deg(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumBHBH20Deg PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumHighFrequency(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumHighFrequency PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumNSNSOpt(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumNSNSOpt PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumNoSRMHighPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumNoSRMHighPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumNoSRMLowPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumNoSRMLowPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumZeroDetHighPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumZeroDetHighPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOQuantumZeroDetLowPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOQuantumZeroDetLowPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOThermal(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOThermal PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOZeroDetHighPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOZeroDetHighPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOZeroDetHighPowerGWINC(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOZeroDetHighPowerGWINC PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOZeroDetLowPower(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOZeroDetLowPower PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOZeroDetLowPowerGWINC(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOZeroDetLowPowerGWINC PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOaLIGO140MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOaLIGO140MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOaLIGO175MpcT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOaLIGO175MpcT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOaLIGODesignSensitivityT1800044(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOaLIGODesignSensitivityT1800044 PSD from LALSimulation.
- pycbc.psd.analytical.aLIGOaLIGOO3LowT1800545(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the aLIGOaLIGOO3LowT1800545 PSD from LALSimulation.
- pycbc.psd.analytical.eLIGOModel(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the eLIGOModel PSD from LALSimulation.
- pycbc.psd.analytical.eLIGOShot(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the eLIGOShot PSD from LALSimulation.
- pycbc.psd.analytical.flat_unity(length, delta_f, low_freq_cutoff)[source]
Returns a FrequencySeries of ones above the low_frequency_cutoff.
- Parameters:
- Returns:
Returns a FrequencySeries containing the unity PSD model.
- Return type:
- pycbc.psd.analytical.from_string(psd_name, length, delta_f, low_freq_cutoff, **kwargs)[source]
Generate a frequency series containing a LALSimulation or built-in space-borne detectors’ PSD specified by name.
- Parameters:
psd_name (string) – PSD name as found in LALSimulation (minus the SimNoisePSD prefix) or pycbc.psd.analytical_space.
length (int) – Length of the frequency series in samples.
delta_f (float) – Frequency resolution of the frequency series.
low_freq_cutoff (float) – Frequencies below this value are set to zero.
**kwargs – All other keyword arguments are passed to the PSD model.
- Returns:
psd – The generated frequency series.
- Return type:
- pycbc.psd.analytical.get_lalsim_psd_list()[source]
Return a list of available reference PSD functions from LALSimulation.
- pycbc.psd.analytical.get_psd_model_list()[source]
Returns a list of available reference PSD functions.
- Returns:
Returns a list of names of reference PSD functions.
- Return type:
- pycbc.psd.analytical.get_pycbc_psd_list()[source]
Return a list of available reference PSD functions coded in PyCBC.
- Returns:
Returns a list of names of all reference PSD functions coded in PyCBC.
- Return type:
- pycbc.psd.analytical.iLIGOModel(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the iLIGOModel PSD from LALSimulation.
- pycbc.psd.analytical.iLIGOSRD(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the iLIGOSRD PSD from LALSimulation.
- pycbc.psd.analytical.iLIGOSeismic(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the iLIGOSeismic PSD from LALSimulation.
- pycbc.psd.analytical.iLIGOShot(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the iLIGOShot PSD from LALSimulation.
- pycbc.psd.analytical.iLIGOThermal(length, delta_f, low_freq_cutoff)
Return a FrequencySeries containing the iLIGOThermal PSD from LALSimulation.
pycbc.psd.analytical_space module
This module provides (semi-)analytical PSDs and sensitivity curves for space borne detectors, such as LISA, Taiji, and TianQin. Based on LISA technical note <LISA-LCST-SGS-TN-001>, LDC manual <LISA-LCST-SGS-MAN-001>, paper <10.1088/1361-6382/ab1101>, <10.1088/0264-9381/33/3/035010>, <10.1103/PhysRevD.102.063021>, <10.1103/PhysRevD.100.043003>, and <10.1103/PhysRevD.107.064021>.
- pycbc.psd.analytical_space.analytical_csd_lisa_tdi_XY(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, tdi=None)[source]
The cross-spectrum density between LISA’s TDI channel X and Y.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The CSD between LISA’s TDI-1.5/2.0 channel X and Y.
- Return type:
Notes
Please see Eq.(56) in <LISA-LCST-SGS-MAN-001(Radler)> for more details.
- pycbc.psd.analytical_space.analytical_psd_lisa_tdi_AE(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, tdi=None)[source]
The PSD of LISA’s TDI-1.5/2.0 channel A and E.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The PSD of LISA’s TDI-1.5/2.0 channel A and E.
- Return type:
Notes
Please see Eq.(58) in <LISA-LCST-SGS-MAN-001(Radler)> for more details.
- pycbc.psd.analytical_space.analytical_psd_lisa_tdi_AE_confusion(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (A,E channel) for LISA with Galactic confusion noise.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
duration (float) – The duration of observation, between 0 and 10, in the unit of years.
tdi (string) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (A,E channel) for LISA with Galactic confusion noise.
- Return type:
- pycbc.psd.analytical_space.analytical_psd_lisa_tdi_T(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, tdi=None)[source]
The PSD of LISA’s TDI-1.5/2.0 channel T.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The PSD of LISA’s TDI-1.5/2.0 channel T.
- Return type:
Notes
Please see Eq.(59) in <LISA-LCST-SGS-MAN-001(Radler)> for more details.
- pycbc.psd.analytical_space.analytical_psd_lisa_tdi_XYZ(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, tdi=None)[source]
The TDI-1.5/2.0 analytical PSD (X,Y,Z channel) for LISA.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for LISA.
- Return type:
Notes
Please see Eq.(19-20) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.analytical_psd_taiji_confusion_noise(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (X,Y,Z channel) for Taiji Galactic confusion noise, no instrumental noise.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
duration (float) – The duration of observation, between 0 and 4, in the unit of years.
tdi (string) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for Taiji Galactic confusion noise, no instrumental noise.
- Return type:
- pycbc.psd.analytical_space.analytical_psd_taiji_tdi_AE(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, acc_noise_level=3e-15, oms_noise_level=8e-12, tdi=None)[source]
The PSD of Taiji’s TDI-1.5/2.0 channel A and E.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The PSD of Taiji’s TDI-1.5/2.0 channel A and E.
- Return type:
Notes
Please see <10.1103/PhysRevD.107.064021> for more details.
- pycbc.psd.analytical_space.analytical_psd_taiji_tdi_AE_confusion(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, acc_noise_level=3e-15, oms_noise_level=8e-12, duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (A,E channel) for Taiji with Galactic confusion noise.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
duration (float) – The duration of observation, between 0 and 4, in the unit of years.
tdi (string) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (A,E channel) for Taiji with Galactic confusion noise.
- Return type:
- pycbc.psd.analytical_space.analytical_psd_taiji_tdi_T(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, acc_noise_level=3e-15, oms_noise_level=8e-12, tdi=None)[source]
The PSD of Taiji’s TDI-1.5/2.0 channel T.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The PSD of Taiji’s TDI-1.5/2.0 channel T.
- Return type:
Notes
Please see <10.1103/PhysRevD.107.064021> for more details.
- pycbc.psd.analytical_space.analytical_psd_taiji_tdi_XYZ(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, acc_noise_level=3e-15, oms_noise_level=8e-12, tdi=None)[source]
The TDI-1.5/2.0 analytical PSD (X,Y,Z channel) for Taiji.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for Taiji.
- Return type:
Notes
Please see <10.1103/PhysRevD.107.064021> for more details.
- pycbc.psd.analytical_space.analytical_psd_tianqin_confusion_noise(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (X,Y,Z channel) for TianQin Galactic confusion noise, no instrumental noise.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
duration (float) – The duration of observation, between 0 and 5, in the unit of years.
tdi (string) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for TianQin Galactic confusion noise, no instrumental noise.
- Return type:
- pycbc.psd.analytical_space.analytical_psd_tianqin_tdi_AE(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, acc_noise_level=1e-15, oms_noise_level=1e-12, tdi=None)[source]
The PSD of TianQin’s TDI-1.5/2.0 channel A and E.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The PSD of TianQin’s TDI-1.5/2.0 channel A and E.
- Return type:
Notes
Please see Table(1) in <10.1088/0264-9381/33/3/035010> for more details.
- pycbc.psd.analytical_space.analytical_psd_tianqin_tdi_AE_confusion(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, acc_noise_level=1e-15, oms_noise_level=1e-12, duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (A,E channel) for TianQin with Galactic confusion noise.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
duration (float) – The duration of observation, between 0 and 5, in the unit of years.
tdi (string) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (A,E channel) for TianQin with Galactic confusion noise.
- Return type:
- pycbc.psd.analytical_space.analytical_psd_tianqin_tdi_T(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, acc_noise_level=1e-15, oms_noise_level=1e-12, tdi=None)[source]
The PSD of TianQin’s TDI-1.5/2.0 channel T.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The PSD of TianQin’s TDI-1.5/2.0 channel T.
- Return type:
Notes
Please see Table(1) in <10.1088/0264-9381/33/3/035010> for more details.
- pycbc.psd.analytical_space.analytical_psd_tianqin_tdi_XYZ(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, acc_noise_level=1e-15, oms_noise_level=1e-12, tdi=None)[source]
The TDI-1.5/2.0 analytical PSD (X,Y,Z channel) for TianQin.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
tdi (str) – The version of TDI. Choose from “1.5” or “2.0”.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for TianQin.
- Return type:
Notes
Please see Table(1) in <10.1088/0264-9381/33/3/035010> for more details.
- pycbc.psd.analytical_space.averaged_fplus_sq_approximated(f, len_arm=None)[source]
A simplified fit for TDI-based space-borne GW detectors’ squared antenna response function, averaged over sky and polarization angle.
\[<\left(F_{X}^{+}\right)^{2}>\approx \frac{3}{20} \frac{1}{ 1+0.6(\omega L)^{2}}\]- Parameters:
- Returns:
fp_sq_approx – The sky and polarization angle averaged squared antenna response.
- Return type:
float or numpy.array
Notes
Please see Eq.(9) in <10.1088/1361-6382/ab1101> for more details.
- pycbc.psd.analytical_space.averaged_lisa_fplus_sq_numerical(f, len_arm=2500000000.0)[source]
A numerical fit for LISA’s squared antenna response function, averaged over sky and polarization angle.
- Parameters:
- Returns:
fp_sq_numerical – The sky and polarization angle averaged squared antenna response.
- Return type:
float or numpy.array
Notes
Please see Eq.(36) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.averaged_response_lisa_tdi(f, len_arm=2500000000.0, tdi=None)[source]
LISA’s TDI-1.5/2.0 response function to GW, averaged over sky and polarization angle.
- Parameters:
- Returns:
response_tdi – The sky and polarization angle averaged TDI-1.5/2.0 response to GW.
- Return type:
float or numpy.array
Notes
Please see Eq.(39-40) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.averaged_response_taiji_tdi(f, len_arm=3000000000.0, tdi=None)[source]
Taiji’s TDI-1.5/2.0 response function to GW, averaged over sky and polarization angle.
- Parameters:
- Returns:
response_tdi – The sky and polarization angle averaged TDI-1.5/2.0 response to GW.
- Return type:
float or numpy.array
- pycbc.psd.analytical_space.averaged_response_tianqin_tdi(f, len_arm=173205080.75688773, tdi=None)[source]
TianQin’s TDI-1.5/2.0 response function to GW, averaged over sky and polarization angle.
- Parameters:
- Returns:
response_tdi – The sky and polarization angle averaged TDI-1.5/2.0 response to GW.
- Return type:
float or numpy.array
- pycbc.psd.analytical_space.averaged_tianqin_fplus_sq_numerical(f, len_arm=173205080.75688773)[source]
A numerical fit for TianQin’s squared antenna response function, averaged over sky and polarization angle.
- Parameters:
- Returns:
fp_sq_numerical – The sky and polarization angle averaged squared antenna response.
- Return type:
float or numpy.array
Notes
Please see Eq.(15-16) in <10.1103/PhysRevD.100.043003> for more details.
- pycbc.psd.analytical_space.confusion_fit_lisa(length, delta_f, low_freq_cutoff, duration=1.0)[source]
The LISA’s sensitivity curve for Galactic confusion noise, averaged over sky and polarization angle. No instrumental noise.
- Parameters:
- Returns:
fseries – The sky and polarization angle averaged LISA’s sensitivity curve for Galactic confusion noise. No instrumental noise.
- Return type:
Notes
Please see Eq.(85-86) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.confusion_fit_taiji(length, delta_f, low_freq_cutoff, duration=1.0)[source]
The Taiji’s sensitivity curve for Galactic confusion noise, averaged over sky and polarization angle. No instrumental noise. Only valid for 0.1 mHz < f < 10 mHz. Note that the results between 0.5, 1, 2, and 4 years are extrapolated, might be non-physical.
- Parameters:
- Returns:
fseries – The sky and polarization angle averaged Taiji’s sensitivity curve for Galactic confusion noise. No instrumental noise.
- Return type:
Notes
Please see Eq.(6) and Table(I) in <10.1103/PhysRevD.107.064021> for more details.
- pycbc.psd.analytical_space.confusion_fit_tianqin(length, delta_f, low_freq_cutoff, duration=1.0)[source]
The TianQin’s sensitivity curve for Galactic confusion noise, averaged over sky and polarization angle. No instrumental noise. Only valid for 0.5 mHz < f < 10 mHz. Note that the results between 0.5, 1, 2, 4, and 5 years are extrapolated, might be non-physical.
- Parameters:
- Returns:
fseries – The sky and polarization angle averaged TianQin’s sensitivity curve for Galactic confusion noise. No instrumental noise.
- Return type:
Notes
Please see Table(II) in <10.1103/PhysRevD.102.063021> for more details.
- pycbc.psd.analytical_space.lisa_psd_components(f, acc_noise_level=3e-15, oms_noise_level=1.5e-11)[source]
The PSD of LISA’s acceleration and OMS noise.
- Parameters:
- Returns:
The PSD value or array for acceleration and OMS noise.
- Return type:
low_freq_component, high_freq_component
- pycbc.psd.analytical_space.omega_length(f, len_arm=None)[source]
The function to calculate 2*pi*f*arm_length.
- pycbc.psd.analytical_space.psd_lisa_acc_noise(f, acc_noise_level=3e-15)[source]
The PSD of LISA’s acceleration noise.
- Parameters:
- Returns:
s_acc_nu – The PSD value or array for acceleration noise.
- Return type:
float or numpy.array
Notes
Please see Eq.(11-13) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.psd_lisa_oms_noise(f, oms_noise_level=1.5e-11)[source]
The PSD of LISA’s OMS noise.
- Parameters:
- Returns:
s_oms_nu – The PSD value or array for OMS noise.
- Return type:
float or numpy.array
Notes
Please see Eq.(9-10) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.psd_taiji_acc_noise(f, acc_noise_level=3e-15)[source]
The PSD of Taiji’s acceleration noise.
- Parameters:
- Returns:
s_acc_nu – The PSD value or array for acceleration noise.
- Return type:
float or numpy.array
Notes
Please see Eq.(2) in <10.1103/PhysRevD.107.064021> for more details.
- pycbc.psd.analytical_space.psd_taiji_oms_noise(f, oms_noise_level=8e-12)[source]
The PSD of Taiji’s OMS noise.
- Parameters:
- Returns:
s_oms_nu – The PSD value or array for OMS noise.
- Return type:
float or numpy.array
Notes
Please see Eq.(1) in <10.1103/PhysRevD.107.064021> for more details.
- pycbc.psd.analytical_space.psd_tianqin_acc_noise(f, acc_noise_level=1e-15)[source]
The PSD of TianQin’s acceleration noise.
- Parameters:
- Returns:
s_acc_nu – The PSD value or array for acceleration noise.
- Return type:
float or numpy.array
Notes
Please see Table(1) in <10.1088/0264-9381/33/3/035010> and that paper for more details.
- pycbc.psd.analytical_space.psd_tianqin_oms_noise(f, oms_noise_level=1e-12)[source]
The PSD of TianQin’s OMS noise.
- Parameters:
- Returns:
s_oms_nu – The PSD value or array for OMS noise.
- Return type:
float or numpy.array
Notes
Please see Table(1) in <10.1088/0264-9381/33/3/035010> and that paper for more details.
- pycbc.psd.analytical_space.semi_analytical_psd_lisa_confusion_noise(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (X,Y,Z channel) for LISA Galactic confusion noise, no instrumental noise.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
duration (float) – The duration of observation, between 0 and 10, in the unit of years.
tdi (string) – The version of TDI, currently only for 1.5 or 2.0.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for LISA Galactic confusion noise, no instrumental noise.
- Return type:
- pycbc.psd.analytical_space.sensitivity_curve_lisa_SciRD(length, delta_f, low_freq_cutoff)[source]
The analytical LISA’s sensitivity curve in SciRD, averaged over sky and polarization angle.
- Parameters:
- Returns:
fseries – The sky and polarization angle averaged analytical LISA’s sensitivity curve in SciRD.
- Return type:
Notes
Please see Eq.(114) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.sensitivity_curve_lisa_confusion(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, base_model='semi', duration=1.0)[source]
The LISA’s sensitivity curve with Galactic confusion noise, averaged over sky and polarization angle.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
base_model (string) – The base model of sensitivity curve, chosen from “semi” or “SciRD”.
duration (float) – The duration of observation, between 0 and 10, in the unit of years.
- Returns:
fseries – The sky and polarization angle averaged LISA’s sensitivity curve with Galactic confusion noise.
- Return type:
Notes
Please see Eq.(85-86) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.sensitivity_curve_lisa_semi_analytical(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11)[source]
The semi-analytical LISA’s sensitivity curve (6-links), averaged over sky and polarization angle.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
- Returns:
fseries – The sky and polarization angle averaged semi-analytical LISA’s sensitivity curve (6-links).
- Return type:
Notes
Please see Eq.(42-43) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.sensitivity_curve_taiji_analytical(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, acc_noise_level=3e-15, oms_noise_level=8e-12)[source]
The analytical Taiji’s sensitivity curve (6-links), averaged over sky and polarization angle.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
- Returns:
fseries – The sky and polarization angle averaged analytical Taiji’s sensitivity curve (6-links).
- Return type:
- pycbc.psd.analytical_space.sensitivity_curve_taiji_confusion(length, delta_f, low_freq_cutoff, len_arm=3000000000.0, acc_noise_level=3e-15, oms_noise_level=8e-12, duration=1.0)[source]
The Taiji’s sensitivity curve with Galactic confusion noise, averaged over sky and polarization angle.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of Taiji, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
duration (float) – The duration of observation, between 0 and 4, in the unit of years.
- Returns:
fseries – The sky and polarization angle averaged Taiji’s sensitivity curve with Galactic confusion noise.
- Return type:
- pycbc.psd.analytical_space.sensitivity_curve_tianqin_analytical(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, acc_noise_level=1e-15, oms_noise_level=1e-12)[source]
The analytical TianQin’s sensitivity curve (6-links), averaged over sky and polarization angle.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
- Returns:
fseries – The sky and polarization angle averaged analytical TianQin’s sensitivity curve (6-links).
- Return type:
- pycbc.psd.analytical_space.sensitivity_curve_tianqin_confusion(length, delta_f, low_freq_cutoff, len_arm=173205080.75688773, acc_noise_level=1e-15, oms_noise_level=1e-12, duration=1.0)[source]
The TianQin’s sensitivity curve with Galactic confusion noise, averaged over sky and polarization angle.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of TianQin, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
duration (float) – The duration of observation, between 0 and 5, in the unit of years.
- Returns:
fseries – The sky and polarization angle averaged TianQin’s sensitivity curve with Galactic confusion noise.
- Return type:
- pycbc.psd.analytical_space.sh_transformed_psd_lisa_tdi_XYZ(length, delta_f, low_freq_cutoff, len_arm=2500000000.0, acc_noise_level=3e-15, oms_noise_level=1.5e-11, base_model='semi', duration=1.0, tdi=None)[source]
The TDI-1.5/2.0 PSD (X,Y,Z channel) for LISA with Galactic confusion noise, transformed from LISA sensitivity curve.
- Parameters:
length (int) – Length of output Frequencyseries.
delta_f (float) – Frequency step for output FrequencySeries.
low_freq_cutoff (float) – Low-frequency cutoff for output FrequencySeries.
len_arm (float) – The arm length of LISA, in the unit of “m”.
acc_noise_level (float) – The level of acceleration noise.
oms_noise_level (float) – The level of OMS noise.
base_model (string) – The base model of sensitivity curve, chosen from “semi” or “SciRD”.
duration (float) – The duration of observation, between 0 and 10, in the unit of years.
tdi (string) – The version of TDI, currently only for 1.5 or 2.0.
- Returns:
fseries – The TDI-1.5/2.0 PSD (X,Y,Z channel) for LISA with Galactic confusion noise, transformed from LISA sensitivity curve.
- Return type:
Notes
Please see Eq.(7,41-43) in <LISA-LCST-SGS-TN-001> for more details.
- pycbc.psd.analytical_space.taiji_psd_components(f, acc_noise_level=3e-15, oms_noise_level=8e-12)[source]
The PSD of Taiji’s acceleration and OMS noise.
- Parameters:
- Returns:
The PSD value or array for acceleration and OMS noise.
- Return type:
low_freq_component, high_freq_component
pycbc.psd.estimate module
Utilites to estimate PSDs from data.
- pycbc.psd.estimate.bandlimited_interpolate(series, delta_f)[source]
Return a new PSD 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
- Returns:
interpolated series – A new FrequencySeries that has been interpolated.
- Return type:
- pycbc.psd.estimate.interpolate(series, delta_f, length=None)[source]
Return a new PSD 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
length (None or int) – The desired number of frequency samples. The default is None, so it will be calculated from the given series and delta_f. But this will cause an inconsistency issue of length sometimes, so if length is given, then just use it.
- Returns:
interpolated series – A new FrequencySeries that has been interpolated.
- Return type:
- pycbc.psd.estimate.inverse_spectrum_truncation(psd, max_filter_len, low_frequency_cutoff=None, trunc_method=None)[source]
Modify a PSD such that the impulse response associated with its inverse square root is no longer than max_filter_len time samples. In practice this corresponds to a coarse graining or smoothing of the PSD.
- Parameters:
psd (FrequencySeries) – PSD whose inverse spectrum is to be truncated.
max_filter_len (int) – Maximum length of the time-domain filter in samples.
low_frequency_cutoff ({None, int}) – Frequencies below low_frequency_cutoff are zeroed in the output.
trunc_method ({None, 'hann'}) – Function used for truncating the time-domain filter. None produces a hard truncation at max_filter_len.
- Returns:
psd – PSD whose inverse spectrum has been truncated.
- Return type:
- Raises:
ValueError – For invalid types or values of max_filter_len and low_frequency_cutoff.
Notes
See arXiv:gr-qc/0509116 for details.
- pycbc.psd.estimate.median_bias(n)[source]
Calculate the bias of the median average PSD computed from n segments.
- Parameters:
n (int) – Number of segments used in PSD estimation.
- Returns:
ans – Calculated bias.
- Return type:
- Raises:
ValueError – For non-integer or non-positive n.
Notes
See arXiv:gr-qc/0509116 appendix B for details.
- pycbc.psd.estimate.welch(timeseries, seg_len=4096, seg_stride=2048, window='hann', avg_method='median', num_segments=None, require_exact_data_fit=False)[source]
PSD estimator based on Welch’s method.
- Parameters:
timeseries (TimeSeries) – Time series for which the PSD is to be estimated.
seg_len (int) – Segment length in samples.
seg_stride (int) – Separation between consecutive segments, in samples.
window ({'hann', numpy.ndarray}) – Function used to window segments before Fourier transforming, or a numpy.ndarray that specifies the window.
avg_method ({'median', 'mean', 'median-mean'}) – Method used for averaging individual segment PSDs.
- Returns:
psd – Frequency series containing the estimated PSD.
- Return type:
- Raises:
ValueError – For invalid choices of seg_len, seg_stride window and avg_method and for inconsistent combinations of len(timeseries), seg_len and seg_stride.
Notes
See arXiv:gr-qc/0509116 for details.
pycbc.psd.read module
Utilities to read PSDs from files.
- pycbc.psd.read.from_numpy_arrays(freq_data, noise_data, length, delta_f, low_freq_cutoff)[source]
Interpolate n PSD (as two 1-dimensional arrays of frequency and data) to the desired length, delta_f and low frequency cutoff.
- Parameters:
freq_data (array) – Array of frequencies.
noise_data (array) – PSD values corresponding to frequencies in freq_arr.
length (int) – Length of the frequency series in samples.
delta_f (float) – Frequency resolution of the frequency series in Herz.
low_freq_cutoff (float) – Frequencies below this value are set to zero.
- Returns:
psd – The generated frequency series.
- Return type:
- pycbc.psd.read.from_txt(filename, length, delta_f, low_freq_cutoff, is_asd_file=True)[source]
Read an ASCII file containing one-sided ASD or PSD data and generate a frequency series with the corresponding PSD. The ASD or PSD data is interpolated in order to match the desired resolution of the generated frequency series.
- Parameters:
filename (string) – Path to a two-column ASCII file. The first column must contain the frequency (positive frequencies only) and the second column must contain the amplitude density OR power spectral density.
length (int) – Length of the frequency series in samples.
delta_f (float) – Frequency resolution of the frequency series in Herz.
low_freq_cutoff (float) – Frequencies below this value are set to zero.
is_asd_file (Boolean) – If false assume that the second column holds power spectral density. If true assume that the second column holds amplitude spectral density. Default: True
- Returns:
psd – The generated frequency series.
- Return type:
- Raises:
ValueError – If the ASCII file contains negative, infinite or NaN frequencies or amplitude densities.
- pycbc.psd.read.from_xml(filename, length, delta_f, low_freq_cutoff, ifo_string=None, root_name='psd')[source]
Read an ASCII file containing one-sided ASD or PSD data and generate a frequency series with the corresponding PSD. The ASD or PSD data is interpolated in order to match the desired resolution of the generated frequency series.
- Parameters:
filename (string) – Path to a two-column ASCII file. The first column must contain the frequency (positive frequencies only) and the second column must contain the amplitude density OR power spectral density.
length (int) – Length of the frequency series in samples.
delta_f (float) – Frequency resolution of the frequency series in Herz.
low_freq_cutoff (float) – Frequencies below this value are set to zero.
ifo_string (string) – Use the PSD in the file’s PSD dictionary with this ifo string. If not given and only one PSD present in the file return that, if not given and multiple (or zero) PSDs present an exception will be raised.
root_name (string (default='psd')) – If given use this as the root name for the PSD XML file. If this means nothing to you, then it is probably safe to ignore this option.
- Returns:
psd – The generated frequency series.
- Return type:
pycbc.psd.variation module
PSD Variation
- pycbc.psd.variation.calc_filt_psd_variation(strain, segment, short_segment, psd_long_segment, psd_duration, psd_stride, psd_avg_method, low_freq, high_freq)[source]
Calculates time series of PSD variability
This function first splits the segment up into 512 second chunks. It then calculates the PSD over this 512 second. The PSD is used to to create a filter that is the composition of three filters: 1. Bandpass filter between f_low and f_high. 2. Weighting filter which gives the rough response of a CBC template. 3. Whitening filter. Next it makes the convolution of this filter with the stretch of data. This new time series is given to the “mean_square” function, which computes the mean square of the timeseries within an 8 seconds window, once per second. The result, which is the variance of the S/N in that stride for the Parseval theorem, is then stored in a timeseries.
- Parameters:
strain (TimeSeries) – Input strain time series to estimate PSDs
segment ({float, 8}) – Duration of the segments for the mean square estimation in seconds.
short_segment ({float, 0.25}) – Duration of the short segments for the outliers removal.
psd_long_segment ({float, 512}) – Duration of the long segments for PSD estimation in seconds.
psd_duration ({float, 8}) – Duration of FFT segments for long term PSD estimation, in seconds.
psd_stride ({float, 4}) – Separation between FFT segments for long term PSD estimation, in seconds.
psd_avg_method ({string, 'median'}) – Method for averaging PSD estimation segments.
low_freq ({float, 20}) – Minimum frequency to consider the comparison between PSDs.
high_freq ({float, 480}) – Maximum frequency to consider the comparison between PSDs.
- Returns:
psd_var – Time series of the variability in the PSD estimation
- Return type:
- pycbc.psd.variation.create_full_filt(freqs, filt, plong, srate, psd_duration)[source]
Create a filter to convolve with strain data to find PSD variation.
- Parameters:
- Returns:
full_filt – The full filter used to calculate PSD variation.
- Return type:
numpy.ndarray
- pycbc.psd.variation.find_trigger_value(psd_var, idx, start, sample_rate)[source]
Find the PSD variation value at a particular time with the filter method. If the time is outside the timeseries bound, 1. is given.
- Parameters:
psd_var (TimeSeries) – Time series of the varaibility in the PSD estimation
idx (numpy.ndarray) – Time indices of the triggers
start (float) – GPS start time
sample_rate (float) – Sample rate defined in ini file
- Returns:
vals – PSD variation value at a particular time
- Return type:
- pycbc.psd.variation.live_calc_psd_variation(strain, full_filt, increment, data_trim=2.0, short_stride=0.25)[source]
Calculate the psd variation in the PyCBC Live search.
The Live strain data is convolved with the filter to produce a timeseries containing the PSD variation values for each sample. The mean square of the timeseries is calculated over the short_stride to find outliers caused by short duration glitches. Outliers are replaced with the average of adjacent elements in the array. This array is then further averaged every second to produce the PSD variation timeseries.
- Parameters:
strain (pycbc.timeseries) – Live data being searched through by the PyCBC Live search.
full_filt (numpy.ndarray) – A filter created by live_create_filter.
increment (float) – The number of seconds in each increment in the PyCBC Live search.
data_trim (float) – The number of seconds to be trimmed from either end of the convolved timeseries to prevent artefacts.
short_stride (float) – The number of seconds to average the PSD variation timeseries over to remove the effects of short duration glitches.
- Returns:
psd_var – A timeseries containing the PSD variation values.
- Return type:
pycbc.timeseries
- pycbc.psd.variation.live_create_filter(psd_estimated, psd_duration, sample_rate, low_freq=20, high_freq=480)[source]
Create a filter to be used in the calculation of the psd variation for the PyCBC Live search. This filter combines a bandpass between a lower and upper frequency and an estimated signal response so that the variance will be 1 when the filter is applied to white noise.
Within the PyCBC Live search this filter needs to be recreated every time the estimated psd is updated and needs to be unique for each detector.
- Parameters:
psd_estimated (pycbc.frequencyseries) – The current PyCBC Live PSD: variations are measured relative to this estimate.
psd_duration (float) – The duration of the estimation of the psd, in seconds.
sample_rate (int) – The sample rate of the strain data being search over.
low_freq (int (default = 20)) – The lower frequency to apply in the bandpass filter.
high_freq (int (default = 480)) – The upper frequency to apply in the bandpass filter.
- Returns:
full_filt – The complete filter to be convolved with the strain data to find the psd variation value.
- Return type:
numpy.ndarray
- pycbc.psd.variation.live_find_var_value(triggers, psd_var_timeseries)[source]
Extract the PSD variation values at trigger times by linear interpolation.
- Parameters:
triggers (dict) – Dictionary containing input trigger times.
psd_var_timeseries (pycbc.timeseries) – A timeseries containing the PSD variation value for each second of the latest increment in PyCBC Live.
- Returns:
psd_var_vals – Array of interpolated PSD variation values at trigger times.
- Return type:
numpy.ndarray
- pycbc.psd.variation.mean_square(data, delta_t, srate, short_stride, stride)[source]
Calculate mean square of given time series once per stride
First of all this function calculate the mean square of given time series once per short_stride. This is used to find and remove outliers due to short glitches. Here an outlier is defined as any element which is greater than two times the average of its closest neighbours. Every outlier is substituted with the average of the corresponding adjacent elements. Then, every second the function compute the mean square of the smoothed time series, within the stride.
- Parameters:
- Returns:
m_s – Mean square of given time series
- Return type:
List
Module contents
- pycbc.psd.associate_psds_to_multi_ifo_segments(opt, fd_segments, gwstrain, flen, delta_f, flow, ifos, dyn_range_factor=1.0, precision=None)[source]
Associate PSDs to segments for all ifos when using the multi-detector CLI
- pycbc.psd.associate_psds_to_segments(opt, fd_segments, gwstrain, flen, delta_f, flow, dyn_range_factor=1.0, precision=None)[source]
Generate a set of overlapping PSDs covering the data in GWstrain. Then associate these PSDs with the appropriate segment in strain_segments.
- Parameters:
opt (object) – Result of parsing the CLI with OptionParser, or any object with the required attributes (psd_model, psd_file, asd_file, psd_estimation, psd_segment_length, psd_segment_stride, psd_inverse_length, psd_output).
fd_segments (StrainSegments.fourier_segments() object) – The fourier transforms of the various analysis segments. The psd attribute of each segment is updated to point to the appropriate PSD.
gwstrain (Strain object) – The timeseries of raw data on which to estimate PSDs.
flen (int) – The length in samples of the output PSDs.
delta_f (float) – The frequency step of the output PSDs.
flow (float) – The low frequncy cutoff to use when calculating the PSD.
dyn_range_factor ({1, float}) – For PSDs taken from models or text files, if dyn_range_factor is not None, then the PSD is multiplied by dyn_range_factor ** 2.
precision (str, choices (None,'single','double')) – If not specified, or specified as None, the precision of the returned PSD will match the precision of the data, if measuring a PSD, or will match the default precision of the model if using an analytical PSD. If ‘single’ the PSD will be converted to float32, if not already in that precision. If ‘double’ the PSD will be converted to float64, if not already in that precision.
- pycbc.psd.associate_psds_to_single_ifo_segments(opt, fd_segments, gwstrain, flen, delta_f, flow, ifo, dyn_range_factor=1.0, precision=None)[source]
Associate PSDs to segments for a single ifo when using the multi-detector CLI
- pycbc.psd.from_cli(opt, length, delta_f, low_frequency_cutoff, strain=None, dyn_range_factor=1, precision=None)[source]
Parses the CLI options related to the noise PSD and returns a FrequencySeries with the corresponding PSD. If necessary, the PSD is linearly interpolated to achieve the resolution specified in the CLI.
- Parameters:
opt (object) – Result of parsing the CLI with OptionParser, or any object with the required attributes (psd_model, psd_file, asd_file, psd_estimation, psd_segment_length, psd_segment_stride, psd_inverse_length, psd_output).
length (int) – The length in samples of the output PSD.
delta_f (float) – The frequency step of the output PSD.
low_frequency_cutoff (float) – The low frequncy cutoff to use when calculating the PSD.
strain ({None, TimeSeries}) – Time series containing the data from which the PSD should be measured, when psd_estimation is in use.
dyn_range_factor ({1, float}) – For PSDs taken from models or text files, if dyn_range_factor is not None, then the PSD is multiplied by dyn_range_factor ** 2.
precision (str, choices (None,'single','double')) – If not specified, or specified as None, the precision of the returned PSD will match the precision of the data, if measuring a PSD, or will match the default precision of the model if using an analytical PSD. If ‘single’ the PSD will be converted to float32, if not already in that precision. If ‘double’ the PSD will be converted to float64, if not already in that precision.
- Returns:
psd – The frequency series containing the PSD.
- Return type:
- pycbc.psd.from_cli_multi_ifos(opt, length_dict, delta_f_dict, low_frequency_cutoff_dict, ifos, strain_dict=None, **kwargs)[source]
Get the PSD for all ifos when using the multi-detector CLI
- pycbc.psd.from_cli_single_ifo(opt, length, delta_f, low_frequency_cutoff, ifo, **kwargs)[source]
Get the PSD for a single ifo when using the multi-detector CLI
- pycbc.psd.generate_overlapping_psds(opt, gwstrain, flen, delta_f, flow, dyn_range_factor=1.0, precision=None)[source]
Generate a set of overlapping PSDs to cover a stretch of data. This allows one to analyse a long stretch of data with PSD measurements that change with time.
- Parameters:
opt (object) – Result of parsing the CLI with OptionParser, or any object with the required attributes (psd_model, psd_file, asd_file, psd_estimation, psd_segment_length, psd_segment_stride, psd_inverse_length, psd_output).
gwstrain (Strain object) – The timeseries of raw data on which to estimate PSDs.
flen (int) – The length in samples of the output PSDs.
delta_f (float) – The frequency step of the output PSDs.
flow (float) – The low frequncy cutoff to use when calculating the PSD.
dyn_range_factor ({1, float}) – For PSDs taken from models or text files, if dyn_range_factor is not None, then the PSD is multiplied by dyn_range_factor ** 2.
precision (str, choices (None,'single','double')) – If not specified, or specified as None, the precision of the returned PSD will match the precision of the data, if measuring a PSD, or will match the default precision of the model if using an analytical PSD. If ‘single’ the PSD will be converted to float32, if not already in that precision. If ‘double’ the PSD will be converted to float64, if not already in that precision.
- Returns:
psd_and_times – This is a list of tuples containing one entry for each PSD. The first and second entries (start, end) in each tuple represent the index range of the gwstrain data that was used to estimate that PSD. The third entry (psd) contains the PSD estimate between that interval.
- Return type:
list of (start, end, PSD) tuples
- pycbc.psd.insert_psd_option_group(parser, output=True, include_data_options=True)[source]
Adds the options used to call the pycbc.psd.from_cli function to an optparser as an OptionGroup. This should be used if you want to use these options in your code.
- Parameters:
parser (object) – OptionParser instance.
- pycbc.psd.insert_psd_option_group_multi_ifo(parser)[source]
Adds the options used to call the pycbc.psd.from_cli function to an optparser as an OptionGroup. This should be used if you want to use these options in your code.
- Parameters:
parser (object) – OptionParser instance.
- pycbc.psd.verify_psd_options(opt, parser)[source]
Parses the CLI options and verifies that they are consistent and reasonable.