Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALSimulation 6.2.0.1-8a6b96f
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

Detailed Description

Routines to compute a stochastic gravitational-wave background power spectrum and to produce a continuous stream of simulated gravitational-wave detector strain.

Prototypes

REAL8FrequencySeriesXLALSimSGWBOmegaGWFlatSpectrum (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...
 
REAL8FrequencySeriesXLALSimSGWBOmegaGWPowerLawSpectrum (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...
 
REAL8FrequencySeriesXLALSimSGWBOmegaGWNumericalSpectrumFromFile (const char *fname, size_t length)
 

Function Documentation

◆ XLALSimSGWBOmegaGWFlatSpectrum()

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.

Parameters
[in]Omega0sgwb spectrum power (dimensionless)
[in]flowlow frequncy cutoff of SGWB spectrum (Hz)
[in]deltaFfrequency bin width (Hz)
[in]lengthnumber of frequency bins

Definition at line 147 of file LALSimSGWB.c.

◆ XLALSimSGWBOmegaGWPowerLawSpectrum()

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.

Parameters
[in]Omegarefsgwb spectrum power at reference frequency (dimensionless)
[in]alphasgwb spectrum power law power
[in]frefreference frequency (Hz)
[in]flowlow frequncy cutoff of SGWB spectrum (Hz)
[in]deltaFfrequency bin width (Hz)
[in]lengthnumber of frequency bins

Definition at line 178 of file LALSimSGWB.c.

◆ XLALSimSGWB()

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.

#include <stdio.h>
#include <gsl/gsl_rng.h>
#include <lal/LALStdlib.h>
#include <lal/LALDetectors.h>
#include <lal/FrequencySeries.h>
#include <lal/TimeSeries.h>
#include <lal/Units.h>
#include <lal/LALSimSGWB.h>
int mksgwbdata(void)
{
const double flow = 40.0; // 40 Hz low frequency cutoff
const double duration = 16.0; // 16 second segments
const double srate = 16384.0; // sampling rate in Hertz
const double Omega0 = 1e-6; // fraction of critical energy in GWs
const double H0 = 0.72 * LAL_H0FAC_SI; // Hubble's constant in seconds
LALDetector detectors[3] = {H1, L1, V1}; // the network of detectors
size_t length = duration * srate; // segment length
size_t stride = length / 2; // stride between segments
LIGOTimeGPS epoch = { 0, 0 };
REAL8FrequencySeries *OmegaGW; // the spectrum of the SGWB
REAL8TimeSeries *h[3]; // the strain induced in the network of detectors
gsl_rng *rng;
gsl_rng_env_setup();
rng = gsl_rng_alloc(gsl_rng_default);
h[0] = XLALCreateREAL8TimeSeries("H1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
h[1] = XLALCreateREAL8TimeSeries("L1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
h[2] = XLALCreateREAL8TimeSeries("V1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
OmegaGW = XLALSimSGWBOmegaGWFlatSpectrum(Omega0, flow, deltaF, seglen/2 + 1);
XLALSimSGWB(h, detectors, 3, 0, Omega0, flow, H0, rng); // first time to initialize
while (1) { // infinite loop
size_t j;
for (j = 0; j < stride; ++j) // output first stride points
printf("%.9f\t%e\t%e\t%e\n", XLALGPSGetREAL8(&h[0]->epoch) + j / srate, h[0]->data->data[j], h[1]->data->data[j], h[2]->data->data[j]);
XLALSimSGWB(h, detectors, 3, stride, Omega0, flow, H0, rng); // make more data
}
double e
Definition: bh_ringdown.c:117
double srate
double duration
double flow
sigmaKerr data[0]
const LALDetector lalCachedDetectors[LAL_NUM_DETECTORS]
#define LAL_H0FAC_SI
LAL_LLO_4K_DETECTOR
LAL_VIRGO_DETECTOR
LAL_LHO_4K_DETECTOR
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 s...
Definition: LALSimSGWB.c:147
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 ...
Definition: LALSimSGWB.c:276
REAL8TimeSeries * XLALCreateREAL8TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
const LALUnit lalStrainUnit
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
LALDetector detectors[LAL_NUM_DETECTORS]
Definition: sgwb.c:112
double Omega0
Definition: sgwb.c:105
LIGOTimeGPS epoch
Definition: unicorn.c:20

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).

Note
  • If stride = 0, initialize h by generating one (periodic) realization of noise; subsequent calls should have non-zero stride.
  • If stride = h->data->length then generate one segment of non-periodic noise by generating two different realizations and feathering them together.
Warning
Only the first stride points are valid.
Parameters
h[in/out] array of sgwb timeseries for detector network
[in]detectorsarray of detectors in network
[in]numDetectorsnumber of detectors in network
[in]stridestride (samples)
[in]OmegaGWsgwb spectrum frequeny series
[in]H0Hubble's constant (s)
[in]rngGSL random number generator

Definition at line 276 of file LALSimSGWB.c.

◆ XLALSimSGWBFlatSpectrum()

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.

#include <stdio.h>
#include <gsl/gsl_rng.h>
#include <lal/LALStdlib.h>
#include <lal/LALDetectors.h>
#include <lal/FrequencySeries.h>
#include <lal/TimeSeries.h>
#include <lal/Units.h>
#include <lal/LALSimSGWB.h>
int mkgwbdata_flat(void)
{
const double flow = 40.0; // 40 Hz low frequency cutoff
const double duration = 16.0; // 16 second segments
const double srate = 16384.0; // sampling rate in Hertz
const double Omega0 = 1e-6; // fraction of critical energy in GWs
const double H0 = 0.72 * LAL_H0FAC_SI; // Hubble's constant in seconds
LALDetector detectors[3] = {H1, L1, V1}; // the network of detectors
size_t length = duration * srate; // segment length
size_t stride = length / 2; // stride between segments
LIGOTimeGPS epoch = { 0, 0 };
REAL8TimeSeries *h[3]; // the strain induced in the network of detectors
gsl_rng *rng;
gsl_rng_env_setup();
rng = gsl_rng_alloc(gsl_rng_default);
h[0] = XLALCreateREAL8TimeSeries("H1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
h[1] = XLALCreateREAL8TimeSeries("L1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
h[2] = XLALCreateREAL8TimeSeries("V1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
XLALSimSGWBFlatSpectrum(h, detectors, 3, 0, Omega0, flow, H0, rng); // first time to initialize
while (1) { // infinite loop
size_t j;
for (j = 0; j < stride; ++j) // output first stride points
printf("%.9f\t%e\t%e\t%e\n", XLALGPSGetREAL8(&h[0]->epoch) + j / srate, h[0]->data->data[j], h[1]->data->data[j], h[2]->data->data[j]);
XLALSimSGWBFlatSpectrum(h, detectors, 3, stride, Omega0, flow, H0, rng); // make more data
}
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 ...
Definition: LALSimSGWB.c:429

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).

Note
  • If stride = 0, initialize h by generating one (periodic) realization of noise; subsequent calls should have non-zero stride.
  • If stride = h->data->length then generate one segment of non-periodic noise by generating two different realizations and feathering them together.
Warning
Only the first stride points are valid.
Parameters
h[in/out] array of sgwb timeseries for detector network
[in]detectorsarray of detectors in network
[in]numDetectorsnumber of detectors in network
[in]stridestride (samples)
[in]Omega0flat sgwb spectrum power (dimensionless)
[in]flowlow frequency cutoff (Hz)
[in]H0Hubble's constant (s)
[in]rngGSL random number generator

Definition at line 429 of file LALSimSGWB.c.

◆ XLALSimSGWBPowerLawSpectrum()

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.

#include <stdio.h>
#include <gsl/gsl_rng.h>
#include <lal/LALStdlib.h>
#include <lal/LALDetectors.h>
#include <lal/FrequencySeries.h>
#include <lal/TimeSeries.h>
#include <lal/Units.h>
#include <lal/LALSimSGWB.h>
int mksgwbdata_powerlaw(void)
{
const double flow = 40.0; // 40 Hz low frequency cutoff
const double duration = 16.0; // 16 second segments
const double srate = 16384.0; // sampling rate in Hertz
const double Omegaref = 1e-6; // fraction of critical energy in GWs
const double fref = 100; // reference frequency in Hertz
const double alpha = 3.0; // sgwb spectrum power law power
const double H0 = 0.72 * LAL_H0FAC_SI; // Hubble's constant in seconds
LALDetector detectors[3] = {H1, L1, V1}; // the network of detectors
size_t length = duration * srate; // segment length
size_t stride = length / 2; // stride between segments
LIGOTimeGPS epoch = { 0, 0 };
REAL8TimeSeries *h[3]; // the strain induced in the network of detectors
gsl_rng *rng;
gsl_rng_env_setup();
rng = gsl_rng_alloc(gsl_rng_default);
h[0] = XLALCreateREAL8TimeSeries("H1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
h[1] = XLALCreateREAL8TimeSeries("L1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
h[2] = XLALCreateREAL8TimeSeries("V1:STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
XLALSimSGWBPowerLawSpectrum(h, detectors, 3, 0, Omegaref, alpha, fref, flow, H0, rng); // first time to initialize
while (1) { // infinite loop
size_t j;
for (j = 0; j < stride; ++j) // output first stride points
printf("%.9f\t%e\t%e\t%e\n", XLALGPSGetREAL8(&h[0]->epoch) + j / srate, h[0]->data->data[j], h[1]->data->data[j], h[2]->data->data[j]);
XLALSimSGWBPowerLawSpectrum(h, detectors, 3, stride, Omegaref, alpha, fref, flow, H0, rng); // make more data
}
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 ...
Definition: LALSimSGWB.c:529
double alpha
Definition: sgwb.c:106
double fref
Definition: sgwb.c:107

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).

Note
  • If stride = 0, initialize h by generating one (periodic) realization of noise; subsequent calls should have non-zero stride.
  • If stride = h->data->length then generate one segment of non-periodic noise by generating two different realizations and feathering them together.
Warning
Only the first stride points are valid.
Parameters
h[in/out] array of sgwb timeseries for detector network
[in]detectorsarray of detectors in network
[in]numDetectorsnumber of detectors in network
[in]stridestride (samples)
[in]Omegarefsgwb spectrum power at reference frequency (dimensionless)
[in]alphasgwb spectrum power power law
[in]frefsgwb spectrum reference frequency (Hz)
[in]flowlow frequency cutoff (Hz)
[in]H0Hubble's constant (s)
[in]rngGSL random number generator

Definition at line 529 of file LALSimSGWB.c.

◆ XLALSimSGWBOmegaGWNumericalSpectrumFromFile()

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.