Routines to compute a stochastic gravitational-wave background power spectrum and to produce a continuous stream of simulated gravitational-wave detector strain.
Prototypes | |
REAL8FrequencySeries * | XLALSimSGWBOmegaGWFlatSpectrum (double Omega0, double flow, double deltaF, size_t length) |
Creates a frequency series that contains a flat SGWB spectrum with the specified power Omega0 above some low frequency cutoff flow. More... | |
REAL8FrequencySeries * | XLALSimSGWBOmegaGWPowerLawSpectrum (double Omegaref, double alpha, double fref, double flow, double deltaF, size_t length) |
Creates a frequency series that contains a power law SGWB spectrum with the specified power Omegaref at reference frequency fref and specified power law power alpha above some low frequency cutoff flow. More... | |
int | XLALSimSGWB (REAL8TimeSeries **h, const LALDetector *detectors, size_t numDetectors, size_t stride, const REAL8FrequencySeries *OmegaGW, double H0, gsl_rng *rng) |
Routine that may be used to generate sequential segments of stochastic background gravitational wave signals for a network of detectors with a specified stride from one segment to the next. More... | |
int | XLALSimSGWBFlatSpectrum (REAL8TimeSeries **h, const LALDetector *detectors, size_t numDetectors, size_t stride, double Omega0, double flow, double H0, gsl_rng *rng) |
Routine that may be used to generate sequential segments of stochastic background gravitational wave signals for a network of detectors with a specified stride from one segment to the next. More... | |
int | XLALSimSGWBPowerLawSpectrum (REAL8TimeSeries **h, const LALDetector *detectors, size_t numDetectors, size_t stride, double Omegaref, double alpha, double fref, double flow, double H0, gsl_rng *rng) |
Routine that may be used to generate sequential segments of stochastic background gravitational wave signals for a network of detectors with a specified stride from one segment to the next. More... | |
REAL8FrequencySeries * | XLALSimSGWBOmegaGWNumericalSpectrumFromFile (const char *fname, size_t length) |
REAL8FrequencySeries * XLALSimSGWBOmegaGWFlatSpectrum | ( | double | Omega0, |
double | flow, | ||
double | deltaF, | ||
size_t | length | ||
) |
Creates a frequency series that contains a flat SGWB spectrum with the specified power Omega0 above some low frequency cutoff flow.
[in] | Omega0 | sgwb spectrum power (dimensionless) |
[in] | flow | low frequncy cutoff of SGWB spectrum (Hz) |
[in] | deltaF | frequency bin width (Hz) |
[in] | length | number of frequency bins |
Definition at line 147 of file LALSimSGWB.c.
REAL8FrequencySeries * XLALSimSGWBOmegaGWPowerLawSpectrum | ( | double | Omegaref, |
double | alpha, | ||
double | fref, | ||
double | flow, | ||
double | deltaF, | ||
size_t | length | ||
) |
Creates a frequency series that contains a power law SGWB spectrum with the specified power Omegaref at reference frequency fref and specified power law power alpha above some low frequency cutoff flow.
[in] | Omegaref | sgwb spectrum power at reference frequency (dimensionless) |
[in] | alpha | sgwb spectrum power law power |
[in] | fref | reference frequency (Hz) |
[in] | flow | low frequncy cutoff of SGWB spectrum (Hz) |
[in] | deltaF | frequency bin width (Hz) |
[in] | length | number of frequency bins |
Definition at line 178 of file LALSimSGWB.c.
int XLALSimSGWB | ( | REAL8TimeSeries ** | h, |
const LALDetector * | detectors, | ||
size_t | numDetectors, | ||
size_t | stride, | ||
const REAL8FrequencySeries * | OmegaGW, | ||
double | H0, | ||
gsl_rng * | rng | ||
) |
Routine that may be used to generate sequential segments of stochastic background gravitational wave signals for a network of detectors with a specified stride from one segment to the next.
The spectrum is specified by the frequency series OmegaGW.
Calling instructions: for the first call, set stride = 0; subsequent calls should pass the same time series and have non-zero stride. This routine will advance the time series by an amount given by the stride and will generate new data so that the data is continuous from one segment to the next. For example: the following routine will output a continuous stream of stochastic background signals with a "flat" (OmegaGW = const) spectrum for the HLV network.
If only one single segment of data is required, set stride to be the length of the timeseries data vector. This will make a single segment of data that is not periodic (also, in this case it will not advance the epoch of the timeseries).
h | [in/out] array of sgwb timeseries for detector network | |
[in] | detectors | array of detectors in network |
[in] | numDetectors | number of detectors in network |
[in] | stride | stride (samples) |
[in] | OmegaGW | sgwb spectrum frequeny series |
[in] | H0 | Hubble's constant (s) |
[in] | rng | GSL random number generator |
Definition at line 276 of file LALSimSGWB.c.
int XLALSimSGWBFlatSpectrum | ( | REAL8TimeSeries ** | h, |
const LALDetector * | detectors, | ||
size_t | numDetectors, | ||
size_t | stride, | ||
double | Omega0, | ||
double | flow, | ||
double | H0, | ||
gsl_rng * | rng | ||
) |
Routine that may be used to generate sequential segments of stochastic background gravitational wave signals for a network of detectors with a specified stride from one segment to the next.
The spectrum is flat for frequencies above flow with power given by Omega0, and zero for frequencies below the low frequency cutoff flow.
Calling instructions: for the first call, set stride = 0; subsequent calls should pass the same time series and have non-zero stride. This routine will advance the time series by an amount given by the stride and will generate new data so that the data is continuous from one segment to the next. For example: the following routine will output a continuous stream of stochastic background signals with a "flat" (OmegaGW = const) spectrum for the HLV network.
If only one single segment of data is required, set stride to be the length of the timeseries data vector. This will make a single segment of data that is not periodic (also, in this case it will not advance the epoch of the timeseries).
h | [in/out] array of sgwb timeseries for detector network | |
[in] | detectors | array of detectors in network |
[in] | numDetectors | number of detectors in network |
[in] | stride | stride (samples) |
[in] | Omega0 | flat sgwb spectrum power (dimensionless) |
[in] | flow | low frequency cutoff (Hz) |
[in] | H0 | Hubble's constant (s) |
[in] | rng | GSL random number generator |
Definition at line 429 of file LALSimSGWB.c.
int XLALSimSGWBPowerLawSpectrum | ( | REAL8TimeSeries ** | h, |
const LALDetector * | detectors, | ||
size_t | numDetectors, | ||
size_t | stride, | ||
double | Omegaref, | ||
double | alpha, | ||
double | fref, | ||
double | flow, | ||
double | H0, | ||
gsl_rng * | rng | ||
) |
Routine that may be used to generate sequential segments of stochastic background gravitational wave signals for a network of detectors with a specified stride from one segment to the next.
The spectrum is a power law with power alpha for frequencies above flow with power given by Omegaref at the reference frequency fref, and zero for frequencies below the low frequency cutoff flow.
Calling instructions: for the first call, set stride = 0; subsequent calls should pass the same time series and have non-zero stride. This routine will advance the time series by an amount given by the stride and will generate new data so that the data is continuous from one segment to the next. For example: the following routine will output a continuous stream of stochastic background signals with a "flat" (OmegaGW = const) spectrum for the HLV network.
If only one single segment of data is required, set stride to be the length of the timeseries data vector. This will make a single segment of data that is not periodic (also, in this case it will not advance the epoch of the timeseries).
h | [in/out] array of sgwb timeseries for detector network | |
[in] | detectors | array of detectors in network |
[in] | numDetectors | number of detectors in network |
[in] | stride | stride (samples) |
[in] | Omegaref | sgwb spectrum power at reference frequency (dimensionless) |
[in] | alpha | sgwb spectrum power power law |
[in] | fref | sgwb spectrum reference frequency (Hz) |
[in] | flow | low frequency cutoff (Hz) |
[in] | H0 | Hubble's constant (s) |
[in] | rng | GSL random number generator |
Definition at line 529 of file LALSimSGWB.c.
REAL8FrequencySeries * XLALSimSGWBOmegaGWNumericalSpectrumFromFile | ( | const char * | fname, |
size_t | length | ||
) |
< [in] sgwb numerical power spectrum
< [in] low frequncy cutoff of SGWB spectrum (Hz)
< [in] frequency bin width (Hz)
Definition at line 558 of file LALSimSGWB.c.