int test_XLALSincInterpolateCOMPLEX8TimeSeries ( void ) { COMPLEX8TimeSeries* tsIn; REAL8 f0 = 100; // heterodyning frequency REAL8 dt = 0.1; // sampling frequency = 10Hz LIGOTimeGPS epoch = { 100, 0 }; REAL8 tStart = XLALGPSGetREAL8 ( &epoch ); UINT4 numSamples = 1000; REAL8 Tspan = numSamples * dt; XLAL_CHECK ( (tsIn = XLALCreateCOMPLEX8TimeSeries ( "test TS_in", &epoch, f0, dt, &emptyLALUnit, numSamples )) != NULL, XLAL_EFUNC ); for ( UINT4 j = 0; j < numSamples; j ++ ) { tsIn->data->data[j] = testSignal ( tStart + j * dt, 0 ); } // for j < numSamples // ---------- interpolate this onto new time-samples UINT4 Dterms = 16; REAL8 safety = (Dterms+1.0) * dt; // avoid truncated interpolation to minimize errors, set to 0 for seeing boundary-effects [they're not so bad...] LIGOTimeGPS epochOut = epoch; XLALGPSAdd ( &epochOut, safety ); REAL8 TspanOut = Tspan - 2 * safety; REAL8 dtOut = dt / 10; UINT4 numSamplesOut = lround ( TspanOut / dtOut ); COMPLEX8TimeSeries *tsOut; XLAL_CHECK ( (tsOut = XLALCreateCOMPLEX8TimeSeries ( "test TS_out", &epochOut, f0, dtOut, &emptyLALUnit, numSamplesOut )) != NULL, XLAL_EFUNC ); REAL8 tStartOut = XLALGPSGetREAL8 ( &epochOut ); REAL8Vector *times_out; XLAL_CHECK ( (times_out = XLALCreateREAL8Vector ( numSamplesOut )) != NULL, XLAL_EFUNC ); for ( UINT4 j = 0; j < numSamplesOut; j ++ ) { REAL8 t_j = tStartOut + j * dtOut; times_out->data[j] = t_j; } // for j < numSamplesOut XLAL_CHECK ( XLALSincInterpolateCOMPLEX8TimeSeries ( tsOut->data, times_out, tsIn, Dterms ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyREAL8Vector ( times_out ); // ---------- check accuracy of interpolation COMPLEX8TimeSeries *tsFull; XLAL_CHECK ( (tsFull = XLALCreateCOMPLEX8TimeSeries ( "test TS_full", &epochOut, f0, dtOut, &emptyLALUnit, numSamplesOut )) != NULL, XLAL_EFUNC ); for ( UINT4 j = 0; j < numSamplesOut; j ++ ) { tsFull->data->data[j] = testSignal ( tStartOut + j * dtOut, 0 ); } // for j < numSamplesOut // ----- out debug info if ( lalDebugLevel & LALINFO ) { XLAL_CHECK ( XLALdumpCOMPLEX8TimeSeries ( "TS_in.dat", tsIn ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( XLALdumpCOMPLEX8TimeSeries ( "TS_out.dat", tsOut ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( XLALdumpCOMPLEX8TimeSeries ( "TS_full.dat", tsFull ) == XLAL_SUCCESS, XLAL_EFUNC ); } // if LALINFO VectorComparison XLAL_INIT_DECL(tol); tol.relErr_L1 = 2e-2; tol.relErr_L2 = 2e-2; tol.angleV = 2e-2; tol.relErr_atMaxAbsx = 2e-2; tol.relErr_atMaxAbsy = 2e-2; XLALPrintInfo ("Comparing sinc-interpolated timeseries to exact signal timeseries:\n"); VectorComparison XLAL_INIT_DECL(cmp); XLAL_CHECK ( XLALCompareCOMPLEX8Vectors ( &cmp, tsOut->data, tsFull->data, &tol ) == XLAL_SUCCESS, XLAL_EFUNC ); // ---------- free memory XLALDestroyCOMPLEX8TimeSeries ( tsIn ); XLALDestroyCOMPLEX8TimeSeries ( tsOut ); XLALDestroyCOMPLEX8TimeSeries ( tsFull ); return XLAL_SUCCESS; } // test_XLALSincInterpolateCOMPLEX8TimeSeries()
/** * Turn the given SFTvector into one long time-series, properly dealing with gaps. */ COMPLEX8TimeSeries * XLALSFTVectorToCOMPLEX8TimeSeries ( const SFTVector *sftsIn /**< [in] SFT vector */ ) { // check input sanity XLAL_CHECK_NULL ( (sftsIn !=NULL) && (sftsIn->length > 0), XLAL_EINVAL ); // create a local copy of the input SFTs, as they will be locally modified! SFTVector *sfts; XLAL_CHECK_NULL ( (sfts = XLALDuplicateSFTVector ( sftsIn )) != NULL, XLAL_EFUNC ); /* define some useful shorthands */ UINT4 numSFTs = sfts->length; SFTtype *firstSFT = &(sfts->data[0]); SFTtype *lastSFT = &(sfts->data[numSFTs-1]); UINT4 numFreqBinsSFT = firstSFT->data->length; REAL8 dfSFT = firstSFT->deltaF; REAL8 Tsft = 1.0 / dfSFT; REAL8 deltaT = Tsft / numFreqBinsSFT; // complex FFT: numSamplesSFT = numFreqBinsSFT REAL8 f0SFT = firstSFT->f0; /* if the start and end input pointers are NOT NULL then determine start and time-span of the final long time-series */ LIGOTimeGPS start = firstSFT->epoch; LIGOTimeGPS end = lastSFT->epoch; XLALGPSAdd ( &end, Tsft ); /* determine output time span */ REAL8 Tspan; XLAL_CHECK_NULL ( (Tspan = XLALGPSDiff ( &end, &start ) ) > 0, XLAL_EINVAL ); UINT4 numSamples = lround ( Tspan / deltaT ); /* determine the heterodyning frequency */ /* fHet = DC of our internal DFTs */ UINT4 NnegSFT = NhalfNeg ( numFreqBinsSFT ); REAL8 fHet = f0SFT + 1.0 * NnegSFT * dfSFT; /* ----- Prepare invFFT of SFTs: compute plan for FFTW */ COMPLEX8FFTPlan *SFTplan; XLAL_CHECK_NULL ( (SFTplan = XLALCreateReverseCOMPLEX8FFTPlan( numFreqBinsSFT, 0 )) != NULL, XLAL_EFUNC ); /* ----- Prepare short time-series holding ONE invFFT of a single SFT */ LIGOTimeGPS XLAL_INIT_DECL(epoch); COMPLEX8TimeSeries *sTS; XLAL_CHECK_NULL ( (sTS = XLALCreateCOMPLEX8TimeSeries ( "short timeseries", &epoch, 0, deltaT, &emptyLALUnit, numFreqBinsSFT )) != NULL, XLAL_EFUNC ); /* ----- prepare long TimeSeries container ---------- */ COMPLEX8TimeSeries *lTS; XLAL_CHECK_NULL ( (lTS = XLALCreateCOMPLEX8TimeSeries ( firstSFT->name, &start, fHet, deltaT, &emptyLALUnit, numSamples )) != NULL, XLAL_EFUNC ); memset ( lTS->data->data, 0, numSamples * sizeof(*lTS->data->data)); /* set all time-samples to zero (in case there are gaps) */ /* ---------- loop over all SFTs and inverse-FFT them ---------- */ for ( UINT4 n = 0; n < numSFTs; n ++ ) { SFTtype *thisSFT = &(sfts->data[n]); /* find bin in long timeseries corresponding to starttime of *this* SFT */ REAL8 offset_n = XLALGPSDiff ( &(thisSFT->epoch), &start ); UINT4 bin0_n = lround ( offset_n / deltaT ); /* round to closest bin */ REAL8 nudge_n = bin0_n * deltaT - offset_n; /* rounding error */ nudge_n = 1e-9 * round ( nudge_n * 1e9 ); /* round to closest nanosecond */ /* nudge SFT into integer timestep bin if necessary */ XLAL_CHECK_NULL ( XLALTimeShiftSFT ( thisSFT, nudge_n ) == XLAL_SUCCESS, XLAL_EFUNC ); /* determine heterodyning phase-correction for this SFT */ REAL8 offset = XLALGPSDiff ( &thisSFT->epoch, &start ); // updated value after time-shift // fHet * Tsft is an integer by construction, because fHet was chosen as a frequency-bin of the input SFTs // therefore we only need the remainder (offset % Tsft) REAL8 offsetEff = fmod ( offset, Tsft ); REAL8 hetCycles = fmod ( fHet * offsetEff, 1); // heterodyning phase-correction for this SFT if ( nudge_n != 0 ){ XLALPrintInfo("n = %d, offset_n = %g, nudge_n = %g, offset = %g, offsetEff = %g, hetCycles = %g\n", n, offset_n, nudge_n, offset, offsetEff, hetCycles ); } REAL4 hetCorrection_re, hetCorrection_im; XLAL_CHECK_NULL ( XLALSinCos2PiLUT ( &hetCorrection_im, &hetCorrection_re, -hetCycles ) == XLAL_SUCCESS, XLAL_EFUNC ); COMPLEX8 hetCorrection = crectf( hetCorrection_re, hetCorrection_im ); /* Note: we also bundle the overall normalization of 'df' into the het-correction. * This ensures that the resulting timeseries will have the correct normalization, according to * x_l = invFT[sft]_l = df * sum_{k=0}^{N-1} xt_k * e^(i 2pi k l / N ) * where x_l is the l-th timestamp, and xt_k is the k-th frequency bin of the SFT. * See the LAL-conventions on FFTs: http://www.ligo.caltech.edu/docs/T/T010095-00.pdf * (the FFTw convention does not contain the factor of 'df', which is why we need to * apply it ourselves) * */ hetCorrection *= dfSFT; XLAL_CHECK_NULL ( XLALReorderSFTtoFFTW (thisSFT->data) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK_NULL ( XLALCOMPLEX8VectorFFT( sTS->data, thisSFT->data, SFTplan ) == XLAL_SUCCESS, XLAL_EFUNC ); for ( UINT4 j=0; j < sTS->data->length; j++) { sTS->data->data[j] *= hetCorrection; } // for j < numFreqBinsSFT // copy the short (shifted) heterodyned timeseries into correct location within long timeseries UINT4 binsLeft = numSamples - bin0_n; UINT4 copyLen = MYMIN ( numFreqBinsSFT, binsLeft ); /* make sure not to write past the end of the long TS */ memcpy ( &lTS->data->data[bin0_n], sTS->data->data, copyLen * sizeof(lTS->data->data[0]) ); } /* for n < numSFTs */ // cleanup memory XLALDestroySFTVector ( sfts ); XLALDestroyCOMPLEX8TimeSeries ( sTS ); XLALDestroyCOMPLEX8FFTPlan ( SFTplan ); return lTS; } // XLALSFTVectorToCOMPLEX8TimeSeries()