proper phase correction for Re/Im Fourier transform (MUSR-206)
This commit is contained in:
parent
a41c518ecc
commit
a63ed33de7
@ -17,6 +17,7 @@ NEW the chi^2 calculation in single-histogram and asymmetry fits is parallelized
|
||||
if musrfit is built using a compiler supporting OpenMP (e.g. GCC >= 4.2)
|
||||
Using --disable-omp this feature can be disabled on the configure level.
|
||||
NEW any2many: force the user to define the exact NeXus ouput format (HDF4,HDF5,XML)
|
||||
FIXED proper phase correction for Re/Im Fourier transform (MUSR-206)
|
||||
FIXED when release or restore is called (command block commands), the
|
||||
corresponding parameter error is set to 2% of its parameter value (MUSR-188).
|
||||
FIXED the wrong normalization when using "SCALE_N0_BKG FALSE" together with view_packing (MUSR-184)
|
||||
|
@ -36,6 +36,7 @@ using namespace std;
|
||||
|
||||
#include "TH1F.h"
|
||||
#include "TF1.h"
|
||||
#include "TAxis.h"
|
||||
|
||||
//#include "TFile.h"
|
||||
|
||||
@ -179,9 +180,18 @@ void PFourier::Transform(UInt_t apodizationTag)
|
||||
fftw_execute(fFFTwPlan);
|
||||
|
||||
// correct the phase for tstart != 0.0
|
||||
// find the first bin >= fStartTime
|
||||
Double_t shiftTime = 0.0;
|
||||
for (Int_t i=1; i<fData->GetXaxis()->GetNbins(); i++) {
|
||||
if (fData->GetXaxis()->GetBinCenter(i) >= fStartTime) {
|
||||
shiftTime = fData->GetXaxis()->GetBinCenter(i);
|
||||
break;
|
||||
}
|
||||
}
|
||||
|
||||
Double_t phase, re, im;
|
||||
for (UInt_t i=0; i<fNoOfBins; i++) {
|
||||
phase = 2.0*PI/(fTimeResolution*fNoOfBins) * i * fStartTime;
|
||||
phase = 2.0*PI/(fTimeResolution*fNoOfBins) * i * shiftTime;
|
||||
re = fOut[i][0] * cos(phase) + fOut[i][1] * sin(phase);
|
||||
im = -fOut[i][0] * sin(phase) + fOut[i][1] * cos(phase);
|
||||
fOut[i][0] = re;
|
||||
|
Loading…
x
Reference in New Issue
Block a user