/**
 * Change frequency-bin order from 'SFT' to fftw-convention
 * ie. from f[-(N-1)/2], ... f[-1], f[0], f[1], .... f[N/2]
 * to FFTW: f[0], f[1],...f[N/2], f[-(N-1)/2], ... f[-2], f[-1]
 */
int
XLALReorderSFTtoFFTW (COMPLEX8Vector *X)
{
  XLAL_CHECK ( (X != NULL) && (X->length > 0), XLAL_EINVAL );

  UINT4 N = X->length;
  UINT4 Npos_and_DC = NhalfPosDC ( N );
  UINT4 Nneg = NhalfNeg ( N );

  /* allocate temporary storage for swap */
  COMPLEX8 *tmp;
  XLAL_CHECK ( (tmp = XLALMalloc ( N * sizeof(*tmp) )) != NULL, XLAL_ENOMEM );

  memcpy ( tmp, X->data, N * sizeof(*tmp) );

  /* Copy second half of FFTW: 'negative' frequencies */
  memcpy ( X->data + Npos_and_DC, tmp, Nneg * sizeof(*tmp) );

  /* Copy first half of FFTW: 'DC + positive' frequencies */
  memcpy ( X->data, tmp + Nneg, Npos_and_DC * sizeof(*tmp) );

  XLALFree(tmp);

  return XLAL_SUCCESS;

}  // XLALReorderSFTtoFFTW()
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()