/** * Unit-Test for function XLALSFTVectorToLFT(). * Generates random data (timeseries + corresponding SFTs), * then feeds the SFTs into XLALSFTVectorToLFT() * and checks correctness of output Fourier transform by * comparing to FT of original real-valued timeseries. * * \note This indirectly also checks XLALSFTVectorToCOMPLEX8TimeSeries() * which is used by XLALSFTVectorToLFT(). * * returns XLAL_SUCCESS on success, XLAL-error otherwise */ int test_XLALSFTVectorToLFT ( void ) { // ----- generate real-valued random timeseries and corresponding SFTs REAL4TimeSeries *tsR4 = NULL; SFTVector *sfts0 = NULL; XLAL_CHECK ( XLALgenerateRandomData ( &tsR4, &sfts0 ) == XLAL_SUCCESS, XLAL_EFUNC ); UINT4 numSamplesR4 = tsR4->data->length; REAL8 dt_R4 = tsR4->deltaT; REAL8 TspanR4 = numSamplesR4 * dt_R4; // ----- consider only the frequency band [3Hz, 4Hz] REAL8 out_fMin = 3; REAL8 out_Band = 1; SFTVector *sftsBand; XLAL_CHECK ( (sftsBand = XLALExtractBandFromSFTVector ( sfts0, out_fMin, out_Band )) != NULL, XLAL_EFUNC ); XLALDestroySFTVector ( sfts0 ); // ----- 1) compute FFT on original REAL4 timeseries ----- REAL4FFTPlan *planR4; XLAL_CHECK ( (planR4 = XLALCreateForwardREAL4FFTPlan ( numSamplesR4, 0 )) != NULL, XLAL_EFUNC ); SFTtype *lft0; XLAL_CHECK ( (lft0 = XLALCreateSFT ( numSamplesR4/2+1 )) != NULL, XLAL_EFUNC ); XLAL_CHECK ( XLALREAL4ForwardFFT ( lft0->data, tsR4->data, planR4 ) == XLAL_SUCCESS, XLAL_EFUNC ); strcpy ( lft0->name, tsR4->name ); lft0->f0 = 0; lft0->epoch = tsR4->epoch; REAL8 dfLFT = 1.0 / TspanR4; lft0->deltaF = dfLFT; // ----- extract frequency band of interest SFTtype *lftR4 = NULL; XLAL_CHECK ( XLALExtractBandFromSFT ( &lftR4, lft0, out_fMin, out_Band ) == XLAL_SUCCESS, XLAL_EFUNC ); for ( UINT4 k=0; k < lftR4->data->length; k ++ ) { lftR4->data->data[k] *= dt_R4; } // ----- 2) compute LFT directly from SFTs ---------- SFTtype *lftSFTs0; XLAL_CHECK ( (lftSFTs0 = XLALSFTVectorToLFT ( sftsBand, 1 )) != NULL, XLAL_EFUNC ); XLALDestroySFTVector ( sftsBand ); // ----- re-extract frequency band of interest SFTtype *lftSFTs = NULL; XLAL_CHECK ( XLALExtractBandFromSFT ( &lftSFTs, lftSFTs0, out_fMin, out_Band ) == XLAL_SUCCESS, XLAL_EFUNC ); if ( lalDebugLevel & LALINFO ) { // ----- write debug output XLAL_CHECK ( XLALdumpREAL4TimeSeries ( "TS_R4.dat", tsR4 ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( write_SFTdata ("LFT_R4T.dat", lftR4 ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( write_SFTdata ("LFT_SFTs.dat", lftSFTs ) == XLAL_SUCCESS, XLAL_EFUNC ); } // end: debug output // ========== compare resulting LFTs ========== VectorComparison XLAL_INIT_DECL(tol0); XLALPrintInfo ("Comparing LFT with itself: should give 0 for all measures\n"); XLAL_CHECK ( XLALCompareSFTs ( lftR4, lftR4, &tol0 ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( XLALCompareSFTs ( lftSFTs, lftSFTs, &tol0 ) == XLAL_SUCCESS, XLAL_EFUNC ); VectorComparison XLAL_INIT_DECL(tol); tol.relErr_L1 = 4e-2; tol.relErr_L2 = 5e-2; tol.angleV = 5e-2; // rad tol.relErr_atMaxAbsx = 1.2e-2; tol.relErr_atMaxAbsy = 1.2e-2; XLALPrintInfo ("Comparing LFT from REAL4-timeseries, to LFT from heterodyned COMPLEX8-timeseries:\n"); XLAL_CHECK ( XLALCompareSFTs ( lftR4, lftSFTs, &tol ) == XLAL_SUCCESS, XLAL_EFUNC ); // ---------- free memory ---------- XLALDestroyREAL4TimeSeries ( tsR4 ); XLALDestroyREAL4FFTPlan ( planR4 ); XLALDestroySFT ( lft0 ); XLALDestroySFT ( lftR4 ); XLALDestroySFT ( lftSFTs ); XLALDestroySFT ( lftSFTs0 ); return XLAL_SUCCESS; } // test_XLALSFTVectorToLFT()
int test_XLALSincInterpolateSFT ( void ) { REAL8 f0 = 0; // heterodyning frequency REAL8 sigmaN = 0.001; REAL8 dt = 0.1; // sampling frequency = 10Hz LIGOTimeGPS epoch = { 100, 0 }; REAL8 tStart = XLALGPSGetREAL8 ( &epoch ); UINT4 numSamples = 1000; REAL8 Tspan = numSamples * dt; REAL8 df = 1.0 / Tspan; UINT4 numSamples0padded = 3 * numSamples; REAL8 Tspan0padded = numSamples0padded * dt; REAL8 df0padded = 1.0 / Tspan0padded; UINT4 numBins = NhalfPosDC ( numSamples ); UINT4 numBins0padded = NhalfPosDC ( numSamples0padded ); // original timeseries REAL4TimeSeries* ts; XLAL_CHECK ( (ts = XLALCreateREAL4TimeSeries ( "test TS_in", &epoch, f0, dt, &emptyLALUnit, numSamples )) != NULL, XLAL_EFUNC ); for ( UINT4 j = 0; j < numSamples; j ++ ) { ts->data->data[j] = crealf ( testSignal ( tStart + j * dt, sigmaN ) ); } // for j < numSamples // zero-padded to double length REAL4TimeSeries* ts0padded; XLAL_CHECK ( (ts0padded = XLALCreateREAL4TimeSeries ( "test TS_padded", &epoch, f0, dt, &emptyLALUnit, numSamples0padded )) != NULL, XLAL_EFUNC ); memcpy ( ts0padded->data->data, ts->data->data, numSamples * sizeof(ts0padded->data->data[0]) ); memset ( ts0padded->data->data + numSamples, 0, (numSamples0padded - numSamples) * sizeof(ts0padded->data->data[0]) ); // compute FFT on ts and ts0padded REAL4FFTPlan *plan, *plan0padded; XLAL_CHECK ( (plan = XLALCreateForwardREAL4FFTPlan ( numSamples, 0 )) != NULL, XLAL_EFUNC ); XLAL_CHECK ( (plan0padded = XLALCreateForwardREAL4FFTPlan ( numSamples0padded, 0 )) != NULL, XLAL_EFUNC ); COMPLEX8Vector *fft, *fft0padded; XLAL_CHECK ( (fft = XLALCreateCOMPLEX8Vector ( numBins )) != NULL, XLAL_ENOMEM ); XLAL_CHECK ( (fft0padded = XLALCreateCOMPLEX8Vector ( numBins0padded )) != NULL, XLAL_ENOMEM ); XLAL_CHECK ( XLALREAL4ForwardFFT ( fft, ts->data, plan ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( XLALREAL4ForwardFFT ( fft0padded, ts0padded->data, plan0padded ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyREAL4TimeSeries ( ts ); XLALDestroyREAL4TimeSeries ( ts0padded ); XLALDestroyREAL4FFTPlan( plan ); XLALDestroyREAL4FFTPlan( plan0padded ); SFTtype XLAL_INIT_DECL(tmp); tmp.f0 = f0; tmp.deltaF = df; tmp.data = fft; REAL8 Band = 0.5/dt - f0; SFTtype *sft = NULL; XLAL_CHECK ( XLALExtractBandFromSFT ( &sft, &tmp, f0, Band ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyCOMPLEX8Vector ( fft ); tmp.f0 = f0; tmp.deltaF = df0padded; tmp.data = fft0padded; SFTtype *sft0padded = NULL; XLAL_CHECK ( XLALExtractBandFromSFT ( &sft0padded, &tmp, f0, Band ) == XLAL_SUCCESS, XLAL_EFUNC ); XLALDestroyCOMPLEX8Vector ( fft0padded ); // ---------- interpolate input SFT onto frequency bins of padded-ts FFT for comparison UINT4 Dterms = 16; REAL8 safetyBins = (Dterms + 1.0); // avoid truncated interpolation to minimize errors, set to 0 for seeing boundary-effects [they're not so bad...] REAL8 fMin = f0 + safetyBins * df; fMin = round(fMin / df0padded) * df0padded; UINT4 numBinsOut = numBins0padded - 3 * safetyBins; REAL8 BandOut = (numBinsOut-1) * df0padded; SFTtype *sftUpsampled = NULL; XLAL_CHECK ( XLALExtractBandFromSFT ( &sftUpsampled, sft0padded, fMin + 0.5*df0padded, BandOut - df0padded ) == XLAL_SUCCESS, XLAL_EFUNC ); SFTtype *sftInterpolated; XLAL_CHECK ( (sftInterpolated = XLALSincInterpolateSFT ( sft, fMin, df0padded, numBinsOut, Dterms )) != NULL, XLAL_EFUNC ); // ----- out debug info if ( lalDebugLevel & LALINFO ) { XLAL_CHECK ( write_SFTdata ( "SFT_in.dat", sft ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( write_SFTdata ( "SFT_in0padded.dat", sft0padded ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( write_SFTdata ( "SFT_upsampled.dat", sftUpsampled ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( write_SFTdata ( "SFT_interpolated.dat", sftInterpolated ) == XLAL_SUCCESS, XLAL_EFUNC ); } // if LALINFO // ---------- check accuracy of interpolation VectorComparison XLAL_INIT_DECL(tol); tol.relErr_L1 = 8e-2; tol.relErr_L2 = 8e-2; tol.angleV = 8e-2; tol.relErr_atMaxAbsx = 4e-3; tol.relErr_atMaxAbsy = 4e-3; XLALPrintInfo ("Comparing Dirichlet SFT interpolation with upsampled SFT:\n"); XLAL_CHECK ( XLALCompareSFTs ( sftUpsampled, sftInterpolated, &tol ) == XLAL_SUCCESS, XLAL_EFUNC ); // ---------- free memory XLALDestroySFT ( sft ); XLALDestroySFT ( sft0padded ); XLALDestroySFT ( sftUpsampled ); XLALDestroySFT ( sftInterpolated ); return XLAL_SUCCESS; } // test_XLALSincInterpolateSFT()
int main ( void ) { const char *fn = __func__; //LALStatus status = empty_status; SFTtype *mySFT; LIGOTimeGPS epoch = { 731210229, 0 }; REAL8 dFreq = 1.0 / 1800.0; REAL8 f0 = 150.0 - 2.0 * dFreq; /* init data array */ COMPLEX8 vals[] = { crectf( -1.249241e-21, 1.194085e-21 ), crectf( 2.207420e-21, 2.472366e-22 ), crectf( 1.497939e-21, 6.593609e-22 ), crectf( 3.544089e-20, -9.365807e-21 ), crectf( 1.292773e-21, -1.402466e-21 ) }; UINT4 numBins = sizeof ( vals ) / sizeof(vals[0] ); if ( (mySFT = XLALCreateSFT ( numBins )) == NULL ) { XLALPrintError ("%s: Failed to create test-SFT using XLALCreateSFT(), xlalErrno = %d\n", fn, xlalErrno ); return XLAL_EFAILED; } /* init header */ strcpy ( mySFT->name, "H1;testSFTRngmed" ); mySFT->epoch = epoch; mySFT->f0 = f0; mySFT->deltaF = dFreq; /* we simply copy over these data-values into the SFT */ UINT4 iBin; for ( iBin = 0; iBin < numBins; iBin ++ ) mySFT->data->data[iBin] = vals[iBin]; /* get memory for running-median vector */ REAL8FrequencySeries rngmed; INIT_MEM ( rngmed ); XLAL_CHECK ( (rngmed.data = XLALCreateREAL8Vector ( numBins )) != NULL, XLAL_EFUNC, "Failed XLALCreateREAL8Vector ( %d )", numBins ); // ---------- Test running-median PSD estimation in simple blocksize cases // ------------------------------------------------------------ // TEST 1: odd blocksize = 3 // ------------------------------------------------------------ UINT4 blockSize3 = 3; /* reference result for 3-bin block running-median computed in octave: octave> sft = [ \ -1.249241e-21 + 1.194085e-21i, \ 2.207420e-21 + 2.472366e-22i, \ 1.497939e-21 + 6.593609e-22i, \ 3.544089e-20 - 9.365807e-21i, \ 1.292773e-21 - 1.402466e-21i \ ]; octave> periodo = abs(sft).^2; octave> m1 = median ( periodo(1:3) ); m2 = median ( periodo(2:4) ); m3 = median ( periodo (3:5 ) ); octave> rngmed = [ m1, m1, m2, m3, m3 ]; octave> printf ("rngmedREF3 = { %.16g, %.16g, %.16g, %.16g, %.16g };\n", rngmed ); rngmedREF3[] = { 2.986442063306e-42, 2.986442063306e-42, 4.933828992779561e-42, 3.638172910684999e-42, 3.638172910684999e-42 }; */ REAL8 rngmedREF3[] = { 2.986442063306e-42, 2.986442063306e-42, 4.933828992779561e-42, 3.638172910684999e-42, 3.638172910684999e-42 }; /* compute running median */ XLAL_CHECK ( XLALSFTtoRngmed ( &rngmed, mySFT, blockSize3 ) == XLAL_SUCCESS, XLAL_EFUNC, "XLALSFTtoRngmed() failed."); /* get median->mean bias correction, needed for octave-reference results, to make * them comparable to the bias-corrected results from LALSFTtoRngmed() */ REAL8 medianBias3 = XLALRngMedBias ( blockSize3 ); XLAL_CHECK ( xlalErrno == 0, XLAL_EFUNC, "XLALRngMedBias() failed."); BOOLEAN pass = 1; const CHAR *passStr; printf ("%4s %22s %22s %8s <%g\n", "Bin", "rngmed(LAL)", "rngmed(Octave)", "relError", tol); for (iBin=0; iBin < numBins; iBin ++ ) { REAL8 rngmedVAL = rngmed.data->data[iBin]; REAL8 rngmedREF = rngmedREF3[iBin] / medianBias3; // apply median-bias correction REAL8 relErr = REL_ERR ( rngmedREF, rngmedVAL ); if ( relErr > tol ) { pass = 0; passStr = "fail"; } else { passStr = "OK."; } printf ("%4d %22.16g %22.16g %8.1g %s\n", iBin, rngmedVAL, rngmedREF, relErr, passStr ); } /* for iBin < numBins */ // ------------------------------------------------------------ // TEST 2: even blocksize = 4 // ------------------------------------------------------------ UINT4 blockSize4 = 4; /* reference result for 4-bin block running-median computed in octave: octave> m1 = median ( periodo(1:4) ); m2 = median ( periodo(2:5) ); octave> rngmed = [ m1, m1, m1, m2, m2 ]; octave> printf ("rngmedREF4[] = { %.16g, %.16g, %.16g, %.16g, %.16g };\n", rngmed ); rngmedREF4[] = { 3.96013552804278e-42, 3.96013552804278e-42, 3.96013552804278e-42, 4.28600095173228e-42, 4.28600095173228e-42 }; */ REAL8 rngmedREF4[] = { 3.96013552804278e-42, 3.96013552804278e-42, 3.96013552804278e-42, 4.28600095173228e-42, 4.28600095173228e-42 }; /* compute running median */ XLAL_CHECK ( XLALSFTtoRngmed ( &rngmed, mySFT, blockSize4 ) == XLAL_SUCCESS, XLAL_EFUNC, "XLALSFTtoRngmed() failed."); /* get median->mean bias correction, needed for octave-reference results, to make * them comparable to the bias-corrected results from LALSFTtoRngmed() */ REAL8 medianBias4 = XLALRngMedBias ( blockSize4 ); XLAL_CHECK ( xlalErrno == 0, XLAL_EFUNC, "XLALRngMedBias() failed."); printf ("%4s %22s %22s %8s <%g\n", "Bin", "rngmed(LAL)", "rngmed(Octave)", "relError", tol); for (iBin=0; iBin < numBins; iBin ++ ) { REAL8 rngmedVAL = rngmed.data->data[iBin]; REAL8 rngmedREF = rngmedREF4[iBin] / medianBias4; // apply median-bias correction REAL8 relErr = REL_ERR ( rngmedREF, rngmedVAL ); if ( relErr > tol ) { pass = 0; passStr = "fail"; } else { passStr = "OK."; } printf ("%4d %22.16g %22.16g %8.1g %s\n", iBin, rngmedVAL, rngmedREF, relErr, passStr ); } /* for iBin < numBins */ /* free memory */ XLALDestroyREAL8Vector ( rngmed.data ); XLALDestroySFT ( mySFT ); LALCheckMemoryLeaks(); if ( !pass ) { printf ("Test failed! Difference exceeded tolerance.\n"); return XLAL_EFAILED; } else { printf ("Test passed.\n"); return XLAL_SUCCESS; } } /* main() */