/*-------------------------------------------------- * calculate the difference of two SFT-vectors *--------------------------------------------------*/ SFTVector * subtractSFTVectors ( const SFTVector *sftvect1, const SFTVector *sftvect2 ) { XLAL_CHECK_NULL ( (sftvect1 != NULL) && (sftvect2 != NULL), XLAL_EINVAL ); XLAL_CHECK_NULL ( sftvect1->length == sftvect2->length, XLAL_EINVAL ); UINT4 numSFTs = sftvect1->length; UINT4 numBins = sftvect1->data[0].data->length; SFTVector *vect; XLAL_CHECK_NULL ( (vect = XLALCreateSFTVector ( numSFTs, numBins )) != NULL, XLAL_EFUNC ); for ( UINT4 alpha = 0; alpha < numSFTs; alpha ++ ) { for ( UINT4 j=0; j < numBins; j++ ) { vect->data[alpha].data->data[j] = sftvect1->data[alpha].data->data[j] - sftvect2->data[alpha].data->data[j]; } } /* for alpha < numSFTs */ return vect; } /* subtractSFTVectors() */
/// /// Make SFTs from given REAL8TimeSeries at given timestamps, potentially applying a time-domain window on each timestretch first /// SFTVector * XLALMakeSFTsFromREAL8TimeSeries ( const REAL8TimeSeries *timeseries, //!< input time-series const LIGOTimeGPSVector *timestamps, //!< timestamps to produce SFTs for (can be NULL), if given must all lies within timeseries' time-span const char *windowType, //!< optional time-domain window function to apply before FFTing REAL8 windowBeta //!< window parameter, if any ) { XLAL_CHECK_NULL ( timeseries != NULL, XLAL_EINVAL, "Invalid NULL input 'timeseries'\n"); XLAL_CHECK_NULL ( timestamps != NULL, XLAL_EINVAL, "Invalid NULL input 'timestamps'\n"); REAL8 dt = timeseries->deltaT; // timeseries timestep */ REAL8 Tsft = timestamps->deltaT; REAL8 df = 1.0 / Tsft; // SFT frequency spacing // make sure that number of timesamples/SFT is an integer (up to possible rounding error 'eps') REAL8 timestepsSFT0 = Tsft / dt; UINT4 timestepsSFT = lround ( timestepsSFT0 ); XLAL_CHECK_NULL ( fabs ( timestepsSFT0 - timestepsSFT ) / timestepsSFT0 < eps, XLAL_ETOL, "Inconsistent sampling-step (dt=%g) and Tsft=%g: must be integer multiple Tsft/dt = %g >= %g\n", dt, Tsft, timestepsSFT0, eps ); // prepare window function if requested REAL8Window *window = NULL; if ( windowType != NULL ) { XLAL_CHECK_NULL ( (window = XLALCreateNamedREAL8Window ( windowType, windowBeta, timestepsSFT )) != NULL, XLAL_EFUNC ); } // ---------- Prepare FFT ---------- REAL8Vector *timeStretchCopy; // input array of length N XLAL_CHECK_NULL ( (timeStretchCopy = XLALCreateREAL8Vector ( timestepsSFT )) != NULL, XLAL_EFUNC, "XLALCreateREAL4Vector(%d) failed.\n", timestepsSFT ); UINT4 numSFTBins = timestepsSFT / 2 + 1; // number of positive frequency-bins + 'DC' to be stored in SFT fftw_complex *fftOut; // output array of length N/2 + 1 XLAL_CHECK_NULL ( (fftOut = fftw_malloc ( numSFTBins * sizeof(fftOut[0]) )) != NULL, XLAL_ENOMEM, "fftw_malloc(%d*sizeof(complex)) failed\n", numSFTBins ); fftw_plan fftplan; // FFTW plan LAL_FFTW_WISDOM_LOCK; XLAL_CHECK_NULL ( (fftplan = fftw_plan_dft_r2c_1d ( timestepsSFT, timeStretchCopy->data, fftOut, FFTW_ESTIMATE)) != NULL, XLAL_EFUNC ); // FIXME: or try FFTW_MEASURE LAL_FFTW_WISDOM_UNLOCK; LIGOTimeGPS tStart = timeseries->epoch; // get last possible start-time for an SFT REAL8 duration = round ( timeseries->data->length * dt ); // rounded to seconds LIGOTimeGPS tLast = tStart; XLALGPSAdd( &tLast, duration - Tsft ); // check that all timestamps lie within [tStart, tLast] for ( UINT4 i = 0; i < timestamps->length; i ++ ) { char buf1[256], buf2[256]; XLAL_CHECK_NULL ( XLALGPSDiff ( &tStart, &(timestamps->data[i]) ) <= 0, XLAL_EDOM, "Timestamp i=%d: %s before start-time %s\n", i, XLALGPSToStr ( buf1, &(timestamps->data[i]) ), XLALGPSToStr ( buf2, &tStart ) ); XLAL_CHECK_NULL ( XLALGPSDiff ( &tLast, &(timestamps->data[i]) ) >=0, XLAL_EDOM, "Timestamp i=%d: %s after last start-time %s\n", i, XLALGPSToStr ( buf1, &(timestamps->data[i]) ), XLALGPSToStr ( buf2, &tLast ) ); } UINT4 numSFTs = timestamps->length; // prepare output SFT-vector SFTVector *sftvect; XLAL_CHECK_NULL ( (sftvect = XLALCreateSFTVector ( numSFTs, numSFTBins )) != NULL, XLAL_EFUNC, "XLALCreateSFTVector(numSFTs=%d, numBins=%d) failed.\n", numSFTs, numSFTBins ); // main loop: apply FFT to the requested time-stretches and store in output SFTs for ( UINT4 iSFT = 0; iSFT < numSFTs; iSFT++ ) { SFTtype *thisSFT = &(sftvect->data[iSFT]); // point to current SFT-slot to store output in // find the start-bin for this SFT in the time-series REAL8 offset = XLALGPSDiff ( &(timestamps->data[iSFT]), &tStart ); INT4 offsetBins = lround ( offset / dt ); // copy timeseries-data for that SFT into local buffer memcpy ( timeStretchCopy->data, timeseries->data->data + offsetBins, timeStretchCopy->length * sizeof(timeStretchCopy->data[0]) ); // window the current time series stretch if required REAL8 sigma_window = 1; if ( window != NULL ) { sigma_window = sqrt ( window->sumofsquares / window->data->length ); for( UINT4 iBin = 0; iBin < timeStretchCopy->length; iBin++ ) { timeStretchCopy->data[iBin] *= window->data->data[iBin]; } } // if window // FFT this time-stretch fftw_execute ( fftplan ); // fill the header of the i'th output SFT */ strcpy ( thisSFT->name, timeseries->name ); thisSFT->epoch = timestamps->data[iSFT]; thisSFT->f0 = timeseries->f0; // SFT starts at heterodyning frequency thisSFT->deltaF = df; // normalize DFT-data to conform to v2 specification ==> multiply DFT by (dt/sigma{window}) // the SFT normalization in case of windowing follows the conventions detailed in the SFTv2 specification, // namely LIGO-T040164, and in particular Eqs.(3),(4) and (6) in T010095-00.pdf // https://dcc.ligo.org/cgi-bin/private/DocDB/ShowDocument?.submit=Number&docid=T010095 // https://dcc.ligo.org/DocDB/0026/T010095/000/T010095-00.pdf REAL8 norm = dt / sigma_window; for ( UINT4 k = 0; k < numSFTBins ; k ++ ) { thisSFT->data->data[k] = (COMPLEX8) ( norm * fftOut[k] ); } // correct heterodyning-phase, IF NECESSARY: ie if (fHet * tStart) is not an integer, such that phase-corr = multiple of 2pi if ( ( (INT4)timeseries->f0 != timeseries->f0 ) || (timeseries->epoch.gpsNanoSeconds != 0) || (thisSFT->epoch.gpsNanoSeconds != 0) ) { XLAL_CHECK_NULL ( XLALcorrect_phase ( thisSFT, timeseries->epoch) == XLAL_SUCCESS, XLAL_EFUNC ); } } // for iSFT < numSFTs // free memory fftw_free ( fftOut ); LAL_FFTW_WISDOM_LOCK; fftw_destroy_plan ( fftplan ); LAL_FFTW_WISDOM_UNLOCK; XLALDestroyREAL8Vector ( timeStretchCopy ); XLALDestroyREAL8Window ( window ); return sftvect; } // XLALMakeSFTsFromREAL8TimeSeries()