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
LALSimNoise.c
Go to the documentation of this file.
1/*
2* Copyright (C) 2007 Jolien Creighton
3*
4* This program is free software; you can redistribute it and/or modify
5* it under the terms of the GNU General Public License as published by
6* the Free Software Foundation; either version 2 of the License, or
7* (at your option) any later version.
8*
9* This program is distributed in the hope that it will be useful,
10* but WITHOUT ANY WARRANTY; without even the implied warranty of
11* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12* GNU General Public License for more details.
13*
14* You should have received a copy of the GNU General Public License
15* along with with program; see the file COPYING. If not, write to the
16* Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17* MA 02110-1301 USA
18*/
19
20#include <complex.h>
21#include <math.h>
22#include <string.h>
23#include <gsl/gsl_rng.h>
24#include <gsl/gsl_randist.h>
25
26#include <lal/Date.h>
27#include <lal/LALConstants.h>
28#include <lal/LALStdlib.h>
29#include <lal/FrequencySeries.h>
30#include <lal/Sequence.h>
31#include <lal/TimeSeries.h>
32#include <lal/TimeFreqFFT.h>
33#include <lal/Units.h>
34#include <lal/LALSimNoise.h>
35
36
37/*
38 * This routine generates a single segment of data. Note that this segment is
39 * generated in the frequency domain and is inverse Fourier transformed into
40 * the time domain; consequently the data is periodic in the time domain.
41 */
43{
44 size_t k;
45 REAL8FFTPlan *plan;
47
48 plan = XLALCreateReverseREAL8FFTPlan(s->data->length, 0);
49 if (! plan)
51
52 stilde = XLALCreateCOMPLEX16FrequencySeries("STILDE", &s->epoch, 0.0, 1.0/(s->data->length * s->deltaT), &lalDimensionlessUnit, s->data->length/2 + 1);
53 if (! stilde) {
56 }
57
58 /* correct units: [stilde] = sqrt([psd] * seconds) */
60 XLALUnitSqrt(&stilde->sampleUnits, &stilde->sampleUnits);
61
62 stilde->data->data[0] = 0.0;
63 for (k = 0; k < s->data->length/2 + 1; ++k) {
64 double sigma = 0.5 * sqrt(psd->data->data[k] / psd->deltaF);
65 stilde->data->data[k] = gsl_ran_gaussian_ziggurat(rng, sigma);
66 stilde->data->data[k] += I * gsl_ran_gaussian_ziggurat(rng, sigma);
67 }
68
69 XLALREAL8FreqTimeFFT(s, stilde, plan);
70
73 return 0;
74}
75
76/**
77 * @addtogroup LALSimNoise_c
78 * @brief Routines to produce a continuous stream of simulated
79 * gravitational-wave detector noise.
80 * @{
81 */
82
83/**
84 * @brief Routine that may be used to generate sequential segments of data with
85 * a specified stride from one segment to the next.
86 *
87 * Calling instructions: for the first call, set stride = 0; subsequent calls
88 * should pass the same time series and have non-zero stride. This routine
89 * will advance the time series by an amount given by the stride and will
90 * generate new data so that the data is continuous from one segment to the
91 * next. For example: the following routine will output a continuous stream of
92 * detector noise with an Initial LIGO spectrum above 40 Hz:
93 *
94 * @code
95 * #include <stdio.h>
96 * #include <gsl/gsl_rng.h>
97 * #include <lal/LALStdlib.h>
98 * #include <lal/FrequencySeries.h>
99 * #include <lal/TimeSeries.h>
100 * #include <lal/Units.h>
101 * #include <lal/LALSimNoise.h>
102 * void mkligodata(void)
103 * {
104 * const double flow = 40.0; // 40 Hz low frequency cutoff
105 * const double duration = 16.0; // 16 second segments
106 * const double srate = 16384.0; // sampling rate in Hertz
107 * size_t length = duration * srate; // segment length
108 * size_t stride = length / 2; // stride between segments
109 * LIGOTimeGPS epoch = { 0, 0 };
110 * REAL8FrequencySeries *psd;
111 * REAL8TimeSeries *seg;
112 * gsl_rng *rng;
113 * gsl_rng_env_setup();
114 * rng = gsl_rng_alloc(gsl_rng_default);
115 * seg = XLALCreateREAL8TimeSeries("STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
116 * psd = XLALCreateREAL8FrequencySeries("LIGO SRD", &epoch, 0.0, 1.0/duration, &lalSecondUnit, length/2 + 1);
117 * XLALSimNoisePSD(psd, flow, XLALSimNoisePSDiLIGOSRD);
118 * XLALSimNoise(seg, 0, psd, rng); // first time to initialize
119 * while (1) { // infinite loop
120 * double t0 = XLALGPSGetREAL8(&seg->epoch);
121 * size_t j;
122 * for (j = 0; j < stride; ++j) // output first stride points
123 * printf("%.9f\t%e\n", t0 + j*seg->deltaT, seg->data->data[j]);
124 * XLALSimNoise(seg, stride, psd, rng); // make more data
125 * }
126 * }
127 * @endcode
128 *
129 * If only one single segment of data is required, set stride to be the length
130 * of the timeseries data vector. This will make a single segment of data
131 * that is *not* periodic (also, in this case it will not advance the epoch of
132 * the timeseries).
133 *
134 * @note
135 * - If stride = 0, initialize h by generating one (periodic) realization of
136 * noise; subsequent calls should have non-zero stride.
137 *
138 * - If stride = h->data->length then generate one segment of non-periodic
139 * noise by generating two different realizations and feathering them
140 * together.
141 *
142 * @warning Only the first stride points are valid.
143 */
145 REAL8TimeSeries *s, /**< [in/out] noise time series */
146 size_t stride, /**< [in] stride (samples) */
147 REAL8FrequencySeries *psd, /**< [in] power spectrum frequency series */
148 gsl_rng *rng /**< [in] GSL random number generator */
149)
150{
151 REAL8Vector *overlap;
152 size_t j;
153
154 /* Use a default RNG if a NULL pointer was passed in */
155 if (!rng)
156 rng = gsl_rng_alloc(gsl_rng_default);
157
158 /* make sure that the resolution of the frequency series is
159 * commensurate with the requested time series */
160 if (s->data->length/2 + 1 != psd->data->length
161 || (size_t)floor(0.5 + 1.0/(s->deltaT * psd->deltaF)) != s->data->length)
163
164 /* stride cannot be longer than data length */
165 if (stride > s->data->length)
167
168 if (stride == 0) { /* generate segment with no feathering */
169 XLALSimNoiseSegment(s, psd, rng);
170 return 0;
171 } else if (stride == s->data->length) {
172 /* will generate two independent noise realizations
173 * and feather them together with full overlap */
174 XLALSimNoiseSegment(s, psd, rng);
175 stride = 0;
176 }
177
178 overlap = XLALCreateREAL8Sequence(s->data->length - stride);
179
180 /* copy overlap region between the old and the new data to temporary storage */
181 memcpy(overlap->data, s->data->data + stride, overlap->length*sizeof(*overlap->data));
182
183 /* generate the new data */
184 XLALSimNoiseSegment(s, psd, rng);
185
186 /* feather old data in overlap region with new data */
187 for (j = 0; j < overlap->length; ++j) {
188 double x = cos(LAL_PI*j/(2.0 * overlap->length));
189 double y = sin(LAL_PI*j/(2.0 * overlap->length));
190 s->data->data[j] = x*overlap->data[j] + y*s->data->data[j];
191 }
192
194
195 /* advance time */
196 XLALGPSAdd(&s->epoch, stride * s->deltaT);
197 return 0;
198}
199
200/** @} */
201
202/*
203 *
204 * TEST CODE
205 *
206 */
207
208#if 0
209
210/* Example routine listed in documentation. */
211void mkligodata(void)
212{
213 const double flow = 40.0; // 40 Hz low frequency cutoff
214 const double duration = 16.0; // 16 second segments
215 const double srate = 16384.0; // sampling rate in Hertz
216 size_t length = duration * srate; // segment length
217 size_t stride = length / 2; // stride between segments
218 LIGOTimeGPS epoch = { 0, 0 };
220 REAL8TimeSeries *seg;
221 gsl_rng *rng;
222 gsl_rng_env_setup();
223 rng = gsl_rng_alloc(gsl_rng_default);
224 seg = XLALCreateREAL8TimeSeries("STRAIN", &epoch, 0.0, 1.0/srate, &lalStrainUnit, length);
225 psd = XLALCreateREAL8FrequencySeries("LIGO SRD", &epoch, 0.0, 1.0/duration, &lalSecondUnit, length/2 + 1);
227 XLALSimNoise(seg, 0, psd, rng); // first time to initialize
228 while (1) { // infinite loop
229 double t0 = XLALGPSGetREAL8(&seg->epoch);
230 size_t j;
231 for (j = 0; j < stride; ++j) // output first stride points
232 printf("%.9f\t%e\n", t0 + j*seg->deltaT, seg->data->data[j]);
233 XLALSimNoise(seg, stride, psd, rng); // make more data
234 }
235}
236
237
238/*
239 * Test routine that generates 1024 seconds of data in 8 second segments with
240 * 50% overlap. Also produced is the PSD for this data for comparison with the
241 * PSD used to generate the data.
242 */
243int test_noise(void)
244{
245 const double recdur = 1024.0; // duration of the entire data record in seconds
246 const double segdur = 8.0; // duration of a segment in seconds
247 const double srate = 4096.0; // sample rate in Hertz
248 const double flow = 9.0; // low frequency cutoff in Hertz
249 LIGOTimeGPS epoch = {0, 0};
250 REAL8TimeSeries *seg;
251 REAL8TimeSeries *rec;
253 REAL8Window *window;
254 REAL8FFTPlan *plan;
255 gsl_rng *rng;
256 double deltaT = 1.0 / srate;
257 double deltaF = 1.0 / segdur;
258 size_t reclen = recdur * srate;
259 size_t seglen = segdur * srate;
260 size_t stride = seglen / 2;
261 size_t numseg = 1 + (reclen - seglen)/stride;
262 size_t klow = flow / deltaF;
263 size_t i, j, k;
264 FILE *fp;
265
266 if (reclen != (numseg - 1)*stride + seglen) {
267 fprintf(stderr, "warning: numseg, seglen, reclen, and stride not commensurate\n");
268 reclen = (numseg - 1)*stride + seglen;
269 }
270
271 rng = gsl_rng_alloc(gsl_rng_default);
272 seg = XLALCreateREAL8TimeSeries("STRAIN", &epoch, 0.0, deltaT, &lalStrainUnit, seglen);
273 rec = XLALCreateREAL8TimeSeries("STRAIN", &epoch, 0.0, deltaT, &lalStrainUnit, reclen);
274 psd = XLALCreateREAL8FrequencySeries("PSD", &epoch, 0.0, deltaF, &lalSecondUnit, seglen/2 + 1);
275
277
278 fp = fopen("psd1.dat", "w");
279 for (k = klow; k < psd->data->length - 1; ++k)
280 fprintf(fp, "%e\t%e\n", k * psd->deltaF, sqrt(psd->data->data[k]));
281 fclose(fp);
282
283 for (i = 0; i < numseg; ++i) {
284 /* first time, initialize; subsequent times, stride */
285 XLALSimNoise(seg, i ? stride : 0, psd, rng);
286 /* copy stride amount of data or seglen on last time */
287 memcpy(rec->data->data + i*stride, seg->data->data, (i == numseg - 1 ? seglen : stride) * sizeof(*rec->data->data));
288 }
289
290 fp = fopen("noise.dat", "w");
291 for (j = 0; j < rec->data->length; ++j) {
292 epoch = rec->epoch;
293 XLALGPSAdd(&epoch, j*rec->deltaT);
294 fprintf(fp, "%.9f\t%e\n", XLALGPSGetREAL8(&epoch), rec->data->data[j]);
295 }
296 fclose(fp);
297
298 plan = XLALCreateForwardREAL8FFTPlan(seglen, 0);
299 window = XLALCreateHannREAL8Window(seglen);
300 XLALREAL8AverageSpectrumWelch(psd, rec, seglen, stride, window, plan);
301
302 fp = fopen("psd2.dat", "w");
303 for (k = klow; k < psd->data->length - 1; ++k)
304 fprintf(fp, "%e\t%e\n", k * psd->deltaF, sqrt(psd->data->data[k]));
305 fclose(fp);
306
312 gsl_rng_free(rng);
313 return 0;
314}
315
316int main(void)
317{
319 gsl_rng_env_setup();
320 // mkligodata();
321 test_noise();
323 return 0;
324}
325
326#endif
void LALCheckMemoryLeaks(void)
static double sigma(const double a, const double b, const sysq *system)
Internal function that computes the spin-spin couplings.
static int XLALSimNoiseSegment(REAL8TimeSeries *s, REAL8FrequencySeries *psd, gsl_rng *rng)
Definition: LALSimNoise.c:42
#define fprintf
int main(int argc, char *argv[])
Definition: bh_qnmode.c:226
int s
Definition: bh_qnmode.c:137
double i
Definition: bh_ringdown.c:118
double segdur
double srate
double duration
double flow
COMPLEX16FrequencySeries * XLALCreateCOMPLEX16FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
REAL8FrequencySeries * XLALCreateREAL8FrequencySeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaF, const LALUnit *sampleUnits, size_t length)
void XLALDestroyREAL8FrequencySeries(REAL8FrequencySeries *series)
void XLALDestroyCOMPLEX16FrequencySeries(COMPLEX16FrequencySeries *series)
#define LAL_PI
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 XLALSimNoisePSDaLIGOHighFrequency(double f)
Provides the noise power spectrum for aLIGO under the configuration tuned to narrow-band high-frequen...
double XLALSimNoisePSDiLIGOSRD(double f)
Provides the noise power spectrum based on a phenomenological fit to the SRD curve for iLIGO.
void XLALDestroyREAL8FFTPlan(REAL8FFTPlan *plan)
REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan(UINT4 size, int measurelvl)
REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan(UINT4 size, int measurelvl)
void XLALDestroyREAL8Sequence(REAL8Sequence *sequence)
REAL8Sequence * XLALCreateREAL8Sequence(size_t length)
int XLALREAL8AverageSpectrumWelch(REAL8FrequencySeries *spectrum, const REAL8TimeSeries *tseries, UINT4 seglen, UINT4 stride, const REAL8Window *window, const REAL8FFTPlan *plan)
int XLALREAL8FreqTimeFFT(REAL8TimeSeries *tser, const COMPLEX16FrequencySeries *freq, const REAL8FFTPlan *plan)
void XLALDestroyREAL8TimeSeries(REAL8TimeSeries *series)
REAL8TimeSeries * XLALCreateREAL8TimeSeries(const CHAR *name, const LIGOTimeGPS *epoch, REAL8 f0, REAL8 deltaT, const LALUnit *sampleUnits, size_t length)
const LALUnit lalStrainUnit
const LALUnit lalSecondUnit
const LALUnit lalDimensionlessUnit
LALUnit * XLALUnitMultiply(LALUnit *output, const LALUnit *unit1, const LALUnit *unit2)
LALUnit * XLALUnitSqrt(LALUnit *output, const LALUnit *input)
void XLALDestroyREAL8Window(REAL8Window *window)
REAL8Window * XLALCreateHannREAL8Window(UINT4 length)
void XLALAbortErrorHandler(const char *func, const char *file, int line, int errnum)
XLALErrorHandlerType * XLALSetErrorHandler(XLALErrorHandlerType *newHandler)
#define XLAL_ERROR(...)
XLAL_EFUNC
XLAL_EINVAL
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
REAL8 XLALGPSGetREAL8(const LIGOTimeGPS *epoch)
list y
x
COMPLEX16Sequence * data
COMPLEX16 * data
REAL8Sequence * data
REAL8Sequence * data
LIGOTimeGPS epoch
REAL8 * data
FILE * fp
LIGOTimeGPS epoch
Definition: unicorn.c:20
double deltaT
Definition: unicorn.c:24