int main(int argc, char *argv[]) { LALFrStream *stream; REAL8TimeSeries *series; LIGOTimeGPS start; XLALSetErrorHandler(XLALAbortErrorHandler); parseargs(argc, argv); /* get the data */ stream = XLALFrStreamCacheOpen(cache); XLALGPSSetREAL8(&start, t0 - pad); series = XLALFrStreamInputREAL8TimeSeries(stream, channel, &start, dt + 2.0 * pad, 0); XLALFrStreamClose(stream); /* manipulate the data */ if (srate > 0) XLALResampleREAL8TimeSeries(series, 1.0 / srate); if (minfreq > 0) XLALHighPassREAL8TimeSeries(series, minfreq, 0.9, 8); if (maxfreq > 0) XLALLowPassREAL8TimeSeries(series, maxfreq, 0.9, 8); if (pad > 0) series = XLALResizeREAL8TimeSeries(series, pad / series->deltaT, dt / series->deltaT); if (df > 0) { /* we are computing a spectrum */ REAL8FrequencySeries *spectrum; REAL8FFTPlan *plan; REAL8Window *window; size_t seglen = 1.0 / (df * series->deltaT); /* make sure that the time series length is commensurate with seglen */ if (((2 * series->data->length) % seglen) != 0) { size_t newlen = ((2 * series->data->length) / seglen) * seglen; series = XLALResizeREAL8TimeSeries(series, 0, newlen); } spectrum = XLALCreateREAL8FrequencySeries(series->name, &series->epoch, 0.0, df, &lalDimensionlessUnit, seglen/2 + 1); plan = XLALCreateForwardREAL8FFTPlan(seglen, 0); window = XLALCreateHannREAL8Window(seglen); XLALREAL8AverageSpectrumWelch(spectrum, series, seglen, seglen/2, window, plan); if (minfreq > 0 || maxfreq > 0) { size_t first = minfreq / spectrum->deltaF; size_t last = maxfreq > 0 ? maxfreq / spectrum->deltaF : spectrum->data->length; spectrum = XLALResizeREAL8FrequencySeries(spectrum, first, last - first); } output_fs(outfile, spectrum); XLALDestroyREAL8Window(window); XLALDestroyREAL8FFTPlan(plan); XLALDestroyREAL8FrequencySeries(spectrum); } else { /* we are outputting a time series */ output_ts(outfile, series); } XLALDestroyREAL8TimeSeries(series); return 0; }
/** * Window and IFFT a FD waveform to TD, then window in TD. * Requires that the FD waveform be generated outside of f_min and f_max. * FD waveform is modified. */ static int FDToTD(REAL8TimeSeries **signalTD, const COMPLEX16FrequencySeries *signalFD, const REAL8 totalMass, const REAL8 deltaT, const REAL8 f_min, const REAL8 f_max, const REAL8 f_min_wide, const REAL8 f_max_wide) { const LIGOTimeGPS gpstime_zero = {0, 0}; const size_t nf = signalFD->data->length; const size_t nt = 2 * (nf - 1); /* Calculate the expected lengths of the FD and TD vectors from the * deltaF and deltaT passed in to this function */ //const size_t exp_nt = (size_t) 1. / (signalFD->deltaF * deltaT); //const size_t exp_nf = (exp_nt / 2) + 1; const REAL8 windowLength = 20. * totalMass * LAL_MTSUN_SI / deltaT; const REAL8 winFLo = 0.2*f_min_wide + 0.8*f_min; /* frequency used for tapering, slightly less than f_min to minimize FFT artifacts * equivalent to winFLo = f_min_wide + 0.8*(f_min - f_min_wide) */ REAL8 winFHi = f_max_wide; COMPLEX16 *FDdata = signalFD->data->data; REAL8FFTPlan *revPlan; REAL8 *TDdata; size_t k; /* check inputs */ if (f_min_wide >= f_min) XLAL_ERROR(XLAL_EDOM); /* apply the softening window function */ if (winFHi > 0.5 / deltaT) winFHi = 0.5 / deltaT; for (k = nf;k--;) { const REAL8 f = k / (deltaT * nt); REAL8 softWin = PlanckTaper(f, f_min_wide, winFLo) * (1.0 - PlanckTaper(f, f_max, winFHi)); FDdata[k] *= softWin; } /* allocate output */ *signalTD = XLALCreateREAL8TimeSeries("h", &gpstime_zero, 0.0, deltaT, &lalStrainUnit, nt); /* Inverse Fourier transform */ revPlan = XLALCreateReverseREAL8FFTPlan(nt, 1); if (!revPlan) { XLALDestroyREAL8TimeSeries(*signalTD); *signalTD = NULL; XLAL_ERROR(XLAL_EFUNC); } XLALREAL8FreqTimeFFT(*signalTD, signalFD, revPlan); XLALDestroyREAL8FFTPlan(revPlan); if (!(*signalTD)) XLAL_ERROR(XLAL_EFUNC); /* apply a linearly decreasing window at the end * of the waveform in order to avoid edge effects. */ if (windowLength > (*signalTD)->data->length) XLAL_ERROR(XLAL_ERANGE); TDdata = (*signalTD)->data->data; for (k = windowLength; k--;) TDdata[nt-k-1] *= k / windowLength; return XLAL_SUCCESS; }
REAL8 calculate_ligo_snr_from_strain_real8( REAL8TimeSeries *strain, const CHAR ifo[3]) { REAL8 ret = -1, snrSq, freq, psdValue; REAL8 deltaF; REAL8FFTPlan *pfwd; COMPLEX16FrequencySeries *fftData; UINT4 k; /* create the time series */ deltaF = strain->deltaT * strain->data->length; fftData = XLALCreateCOMPLEX16FrequencySeries( strain->name, &(strain->epoch), 0, deltaF, &lalDimensionlessUnit, strain->data->length/2 + 1 ); /* perform the fft */ pfwd = XLALCreateForwardREAL8FFTPlan( strain->data->length, 0 ); XLALREAL8TimeFreqFFT( fftData, strain, pfwd ); /* compute the SNR for initial LIGO at design */ for ( snrSq = 0, k = 0; k < fftData->data->length; k++ ) { freq = fftData->deltaF * k; if ( ifo[0] == 'V' ) { if (freq < 35) continue; LALVIRGOPsd( NULL, &psdValue, freq ); psdValue /= 9e-46; } else { if (freq < 40) continue; LALLIGOIPsd( NULL, &psdValue, freq ); } fftData->data->data[k] /= 3e-23; snrSq += creal(fftData->data->data[k]) * creal(fftData->data->data[k]) / psdValue; snrSq += cimag(fftData->data->data[k]) * cimag(fftData->data->data[k]) / psdValue; } snrSq *= 4*fftData->deltaF; XLALDestroyREAL8FFTPlan( pfwd ); XLALDestroyCOMPLEX16FrequencySeries( fftData ); ret = sqrt(snrSq); return ret; }
/* creates a waveform in the frequency domain; the waveform might be generated * in the time-domain and transformed */ int create_fd_waveform(COMPLEX16FrequencySeries ** htilde_plus, COMPLEX16FrequencySeries ** htilde_cross, struct params p) { clock_t timer_start = 0; double chirplen, deltaF; int chirplen_exp; /* length of the chirp in samples */ chirplen = imr_time_bound(p.f_min, p.m1, p.m2, p.s1z, p.s2z) * p.srate; /* make chirplen next power of two */ frexp(chirplen, &chirplen_exp); chirplen = ldexp(1.0, chirplen_exp); deltaF = p.srate / chirplen; if (p.verbose) fprintf(stderr, "using frequency resolution deltaF = %g Hz\n", deltaF); if (p.condition) { if (p.verbose) { fprintf(stderr, "generating waveform in frequency domain using XLALSimInspiralFD...\n"); timer_start = clock(); } XLALSimInspiralFD(htilde_plus, htilde_cross, p.m1, p.m2, p.s1x, p.s1y, p.s1z, p.s2x, p.s2y, p.s2z, p.distance, p.inclination, p.phiRef, p.longAscNodes, p.eccentricity, p.meanPerAno, deltaF, p.f_min, 0.5 * p.srate, p.fRef, p.params, p.approx); if (p.verbose) fprintf(stderr, "generation took %g seconds\n", (double)(clock() - timer_start) / CLOCKS_PER_SEC); } else if (p.domain == LAL_SIM_DOMAIN_FREQUENCY) { if (p.verbose) { fprintf(stderr, "generating waveform in frequency domain using XLALSimInspiralChooseFDWaveform...\n"); timer_start = clock(); } XLALSimInspiralChooseFDWaveform(htilde_plus, htilde_cross, p.m1, p.m2, p.s1x, p.s1y, p.s1z, p.s2x, p.s2y, p.s2z, p.distance, p.inclination, p.phiRef, p.longAscNodes, p.eccentricity, p.meanPerAno, deltaF, p.f_min, 0.5 * p.srate, p.fRef, p.params, p.approx); if (p.verbose) fprintf(stderr, "generation took %g seconds\n", (double)(clock() - timer_start) / CLOCKS_PER_SEC); } else { REAL8TimeSeries *h_plus = NULL; REAL8TimeSeries *h_cross = NULL; REAL8FFTPlan *plan; /* generate time domain waveform */ if (p.verbose) { fprintf(stderr, "generating waveform in time domain using XLALSimInspiralChooseTDWaveform...\n"); timer_start = clock(); } XLALSimInspiralChooseTDWaveform(&h_plus, &h_cross, p.m1, p.m2, p.s1x, p.s1y, p.s1z, p.s2x, p.s2y, p.s2z, p.distance, p.inclination, p.phiRef, p.longAscNodes, p.eccentricity, p.meanPerAno, 1.0 / p.srate, p.f_min, p.fRef, p.params, p.approx); if (p.verbose) fprintf(stderr, "generation took %g seconds\n", (double)(clock() - timer_start) / CLOCKS_PER_SEC); /* resize the waveforms to the required length */ XLALResizeREAL8TimeSeries(h_plus, h_plus->data->length - (size_t) chirplen, (size_t) chirplen); XLALResizeREAL8TimeSeries(h_cross, h_cross->data->length - (size_t) chirplen, (size_t) chirplen); /* put the waveform in the frequency domain */ /* (the units will correct themselves) */ if (p.verbose) { fprintf(stderr, "transforming waveform to frequency domain...\n"); timer_start = clock(); } *htilde_plus = XLALCreateCOMPLEX16FrequencySeries("htilde_plus", &h_plus->epoch, 0.0, deltaF, &lalDimensionlessUnit, (size_t) chirplen / 2 + 1); *htilde_cross = XLALCreateCOMPLEX16FrequencySeries("htilde_cross", &h_cross->epoch, 0.0, deltaF, &lalDimensionlessUnit, (size_t) chirplen / 2 + 1); plan = XLALCreateForwardREAL8FFTPlan((size_t) chirplen, 0); XLALREAL8TimeFreqFFT(*htilde_cross, h_cross, plan); XLALREAL8TimeFreqFFT(*htilde_plus, h_plus, plan); if (p.verbose) fprintf(stderr, "transformation took %g seconds\n", (double)(clock() - timer_start) / CLOCKS_PER_SEC); /* clean up */ XLALDestroyREAL8FFTPlan(plan); XLALDestroyREAL8TimeSeries(h_cross); XLALDestroyREAL8TimeSeries(h_plus); } return 0; }
REAL8 calculate_snr_from_strain_and_psd_real8( REAL8TimeSeries *strain, REAL8FrequencySeries *psd, REAL8 startFreq, const CHAR ifo[3]) { REAL8 ret = -1, snrSq, freq, psdValue; REAL8 deltaF; REAL8FFTPlan *pfwd; COMPLEX16FrequencySeries *fftData; UINT4 k; /* create the time series */ deltaF = strain->deltaT * strain->data->length; fftData = XLALCreateCOMPLEX16FrequencySeries( strain->name, &(strain->epoch), 0, deltaF, &lalDimensionlessUnit, strain->data->length/2 + 1 ); /* perform the fft */ pfwd = XLALCreateForwardREAL8FFTPlan( strain->data->length, 0 ); XLALREAL8TimeFreqFFT( fftData, strain, pfwd ); /* The PSD, if provided, comes in as it was in the original file */ /* since we don't know deltaF until we get here. Interpolate now */ if ( psd ) { psd = XLALInterpolatePSD(psd, 1.0 / deltaF); } /* compute the SNR for initial LIGO at design */ for ( snrSq = 0, k = 0; k < fftData->data->length; k++ ) { freq = fftData->deltaF * k; if ( psd ) { if ( freq < startFreq || k > psd->data->length ) continue; psdValue = psd->data->data[k]; psdValue /= 9e-46; } else if ( ifo[0] == 'V' ) { if (freq < 35) continue; LALVIRGOPsd( NULL, &psdValue, freq ); psdValue /= 9e-46; } else { if (freq < 40) continue; LALLIGOIPsd( NULL, &psdValue, freq ); } fftData->data->data[k] /= 3e-23; snrSq += creal(fftData->data->data[k]) * creal(fftData->data->data[k]) / psdValue; snrSq += cimag(fftData->data->data[k]) * cimag(fftData->data->data[k]) / psdValue; } snrSq *= 4*fftData->deltaF; XLALDestroyREAL8FFTPlan( pfwd ); XLALDestroyCOMPLEX16FrequencySeries( fftData ); if ( psd ) XLALDestroyREAL8FrequencySeries( psd ); ret = sqrt(snrSq); printf("Obtained snr=%f\n", ret); return ret; }