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

Detailed Description

Routines to produce a continuous stream of simulated gravitational-wave detector noise.

Prototypes

int XLALSimNoise (REAL8TimeSeries *s, size_t stride, REAL8FrequencySeries *psd, gsl_rng *rng)
 Routine that may be used to generate sequential segments of data with a specified stride from one segment to the next. More...
 

Function Documentation

◆ XLALSimNoise()

int XLALSimNoise ( REAL8TimeSeries s,
size_t  stride,
REAL8FrequencySeries psd,
gsl_rng *  rng 
)

Routine that may be used to generate sequential segments of data with a specified stride from one segment to the next.

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 detector noise with an Initial LIGO spectrum above 40 Hz:

#include <stdio.h>
#include <gsl/gsl_rng.h>
#include <lal/LALStdlib.h>
#include <lal/FrequencySeries.h>
#include <lal/TimeSeries.h>
#include <lal/Units.h>
#include <lal/LALSimNoise.h>
void mkligodata(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
size_t length = duration * srate; // segment length
size_t stride = length / 2; // stride between segments
LIGOTimeGPS epoch = { 0, 0 };
gsl_rng *rng;
gsl_rng_env_setup();
rng = gsl_rng_alloc(gsl_rng_default);
seg = XLALCreateREAL8TimeSeries("STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
psd = XLALCreateREAL8FrequencySeries("LIGO SRD", &epoch, 0.0, 1.0/duration, &lalSecondUnit, length/2 + 1);
XLALSimNoise(seg, 0, psd, rng); // first time to initialize
while (1) { // infinite loop
double t0 = XLALGPSGetREAL8(&seg->epoch);
size_t j;
for (j = 0; j < stride; ++j) // output first stride points
printf("%.9f\t%e\n", t0 + j*seg->deltaT, seg->data->data[j]);
XLALSimNoise(seg, stride, psd, rng); // make more data
}
}
double srate
double duration
double flow
REAL8FrequencySeries * XLALCreateREAL8FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
int XLALSimNoise(REAL8TimeSeries *s, size_t stride, REAL8FrequencySeries *psd, gsl_rng *rng)
Routine that may be used to generate sequential segments of data with a specified stride from one seg...
Definition: LALSimNoise.c:144
int XLALSimNoisePSD(REAL8FrequencySeries *psd, double flow, double(*psdfunc)(double))
Evaluates a power spectral density function, psdfunc, at the frequencies required to populate the fre...
double XLALSimNoisePSDiLIGOSRD(double f)
Provides the noise power spectrum based on a phenomenological fit to the SRD curve for iLIGO.
REAL8TimeSeries * XLALCreateREAL8TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
const LALUnit lalStrainUnit
const LALUnit lalSecondUnit
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
REAL8Sequence * data
LIGOTimeGPS epoch
REAL8 * data
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
s[in/out] noise time series
[in]stridestride (samples)
[in]psdpower spectrum frequency series
[in]rngGSL random number generator

Definition at line 144 of file LALSimNoise.c.