Public Member Functions | Private Member Functions | Private Attributes | List of all members
util::LArFFT Class Reference

#include <LArFFT.h>

Public Member Functions

 LArFFT (fhicl::ParameterSet const &pset, art::ActivityRegistry &reg)
 
 ~LArFFT ()
 
template<class T >
void DoFFT (std::vector< T > &input, std::vector< TComplex > &output)
 
template<class T >
void DoInvFFT (std::vector< TComplex > &input, std::vector< T > &output)
 
template<class T >
void Deconvolute (std::vector< T > &input, std::vector< T > &respFunc)
 
template<class T >
void Deconvolute (std::vector< T > &input, std::vector< TComplex > &kern)
 
template<class T >
void Convolute (std::vector< T > &input, std::vector< T > &respFunc)
 
template<class T >
void Convolute (std::vector< T > &input, std::vector< TComplex > &kern)
 
template<class T >
void Correlate (std::vector< T > &input, std::vector< T > &respFunc)
 
template<class T >
void Correlate (std::vector< T > &input, std::vector< TComplex > &kern)
 
template<class T >
void AlignedSum (std::vector< T > &input, std::vector< T > &output, bool add=true)
 
void ShiftData (std::vector< TComplex > &input, double shift)
 
template<class T >
void ShiftData (std::vector< T > &input, double shift)
 
template<class T >
PeakCorrelation (std::vector< T > &shape1, std::vector< T > &shape2)
 
int FFTSize () const
 
std::string FFTOptions () const
 
int FFTFitBins () const
 
void ReinitializeFFT (int, std::string, int)
 

Private Member Functions

void InitializeFFT ()
 
void resetSizePerRun (art::Run const &)
 

Private Attributes

int fSize
 
int fFreqSize
 
std::string fOption
 
int fFitBins
 
TF1 * fPeakFit
 
TH1D * fConvHist
 
std::vector< TComplex > fCompTemp
 
std::vector< TComplex > fKern
 
TFFTRealComplex * fFFT
 object to do FFT More...
 
TFFTComplexReal * fInverseFFT
 object to do Inverse FF More...
 

Detailed Description

Definition at line 26 of file LArFFT.h.

Constructor & Destructor Documentation

util::LArFFT::LArFFT ( fhicl::ParameterSet const &  pset,
art::ActivityRegistry reg 
)

Definition at line 20 of file LArFFT_service.cc.

21  : fSize(pset.get<int>("FFTSize", 0))
22  , fOption(pset.get<std::string>("FFTOption"))
23  , fFitBins(pset.get<int>("FitBins"))
24 {
25  // Default to the readout window size if the user didn't input
26  // a specific size
27  if (fSize <= 0) {
28  // Creating a service handle to DetectorPropertiesService not only
29  // creates the service if it doesn't exist, it also guarantees
30  // that its callbacks are invoked before any of LArFFT's callbacks
31  // are invoked.
33  ->DataForJob()
34  .ReadOutWindowSize();
36  }
37  InitializeFFT();
38 }
void resetSizePerRun(art::Run const &)
std::string string
Definition: nybbler.cc:12
int fSize
Definition: LArFFT.h:77
GlobalSignal< detail::SignalResponseType::FIFO, void(Run const &)> sPreBeginRun
std::string fOption
Definition: LArFFT.h:79
int fFitBins
Definition: LArFFT.h:80
void InitializeFFT()
util::LArFFT::~LArFFT ( )

Definition at line 80 of file LArFFT_service.cc.

81 {
82  delete fFFT;
83  delete fInverseFFT;
84  delete fPeakFit;
85  delete fConvHist;
86 }
TFFTRealComplex * fFFT
object to do FFT
Definition: LArFFT.h:86
TF1 * fPeakFit
Definition: LArFFT.h:81
TFFTComplexReal * fInverseFFT
object to do Inverse FF
Definition: LArFFT.h:87
TH1D * fConvHist
Definition: LArFFT.h:82

Member Function Documentation

template<class T >
void util::LArFFT::AlignedSum ( std::vector< T > &  input,
std::vector< T > &  output,
bool  add = true 
)
inline

Definition at line 244 of file LArFFT.h.

247 {
248  double shift = PeakCorrelation(shape1,shape2);
249 
250  ShiftData(shape1,shift);
251 
252  if(add)for(int i = 0; i < fSize; i++) shape1[i]+=shape2[i];
253 
254  return;
255 }
void ShiftData(std::vector< TComplex > &input, double shift)
int fSize
Definition: LArFFT.h:77
Coord add(Coord c1, Coord c2)
Definition: restypedef.cpp:23
T PeakCorrelation(std::vector< T > &shape1, std::vector< T > &shape2)
Definition: LArFFT.h:272
template<class T >
void util::LArFFT::Convolute ( std::vector< T > &  input,
std::vector< T > &  respFunc 
)
inline

Definition at line 173 of file LArFFT.h.

175 {
176  DoFFT(shape1, fKern);
177  DoFFT(shape2, fCompTemp);
178 
179  for(int i = 0; i < fFreqSize; i++)
180  fCompTemp[i]*=fKern[i];
181 
182  DoInvFFT(fCompTemp, shape1);
183 
184  return;
185 }
std::vector< TComplex > fKern
Definition: LArFFT.h:84
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::Convolute ( std::vector< T > &  input,
std::vector< TComplex > &  kern 
)
inline

Definition at line 192 of file LArFFT.h.

194 {
195  DoFFT(input, fCompTemp);
196 
197  for(int i = 0; i < fFreqSize; i++)
198  fCompTemp[i]*=kern[i];
199 
200  DoInvFFT(fCompTemp, input);
201 
202  return;
203 }
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::Correlate ( std::vector< T > &  input,
std::vector< T > &  respFunc 
)
inline

Definition at line 207 of file LArFFT.h.

209 {
210  DoFFT(shape1, fKern);
211  DoFFT(shape2, fCompTemp);
212 
213  for(int i = 0; i < fFreqSize; i++)
214  fCompTemp[i]*=TComplex::Conjugate(fKern[i]);
215 
216  DoInvFFT(fCompTemp, shape1);
217 
218  return;
219 }
std::vector< TComplex > fKern
Definition: LArFFT.h:84
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::Correlate ( std::vector< T > &  input,
std::vector< TComplex > &  kern 
)
inline

Definition at line 226 of file LArFFT.h.

228 {
229  DoFFT(input, fCompTemp);
230 
231  for(int i = 0; i < fFreqSize; i++)
232  fCompTemp[i]*=TComplex::Conjugate(kern[i]);
233 
234  DoInvFFT(fCompTemp, input);
235 
236  return;
237 }
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::Deconvolute ( std::vector< T > &  input,
std::vector< T > &  respFunc 
)
inline

Definition at line 138 of file LArFFT.h.

140 {
141  DoFFT(respFunction, fKern);
142  DoFFT(input, fCompTemp);
143 
144  for(int i = 0; i < fFreqSize; i++)
145  fCompTemp[i]/=fKern[i];
146 
147  DoInvFFT(fCompTemp, input);
148 
149  return;
150 }
std::vector< TComplex > fKern
Definition: LArFFT.h:84
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::Deconvolute ( std::vector< T > &  input,
std::vector< TComplex > &  kern 
)
inline

Definition at line 157 of file LArFFT.h.

159 {
160  DoFFT(input, fCompTemp);
161 
162  for(int i = 0; i < fFreqSize; i++)
163  fCompTemp[i]/=kern[i];
164 
165  DoInvFFT(fCompTemp, input);
166 
167  return;
168 }
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::DoFFT ( std::vector< T > &  input,
std::vector< TComplex > &  output 
)
inline

Definition at line 98 of file LArFFT.h.

100 {
101  double real = 0.; //real value holder
102  double imaginary = 0.; //imaginary value hold
103 
104  // set the points
105  for(size_t p = 0; p < input.size(); ++p)
106  fFFT->SetPoint(p, input[p]);
107 
108  fFFT->Transform();
109 
110  for(int i = 0; i < fFreqSize; ++i){
111  fFFT->GetPointComplex(i, real, imaginary);
112  output[i]=TComplex(real, imaginary);
113  }
114 
115  return;
116 }
TFFTRealComplex * fFFT
object to do FFT
Definition: LArFFT.h:86
p
Definition: test.py:223
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::DoInvFFT ( std::vector< TComplex > &  input,
std::vector< T > &  output 
)
inline

Definition at line 120 of file LArFFT.h.

122 {
123  for(int i = 0; i < fFreqSize; ++i)
124  fInverseFFT->SetPointComplex(i, input[i]);
125 
126  fInverseFFT->Transform();
127  double factor = 1.0/(double) fSize;
128 
129  for(int i = 0; i < fSize; ++i)
130  output[i] = factor*fInverseFFT->GetPointReal(i,false);
131 
132  return;
133 }
int fSize
Definition: LArFFT.h:77
TFFTComplexReal * fInverseFFT
object to do Inverse FF
Definition: LArFFT.h:87
int fFreqSize
Definition: LArFFT.h:78
int util::LArFFT::FFTFitBins ( ) const
inline

Definition at line 71 of file LArFFT.h.

71 { return fFitBins; }
int fFitBins
Definition: LArFFT.h:80
std::string util::LArFFT::FFTOptions ( ) const
inline

Definition at line 70 of file LArFFT.h.

70 { return fOption; }
std::string fOption
Definition: LArFFT.h:79
int util::LArFFT::FFTSize ( ) const
inline

Definition at line 69 of file LArFFT.h.

69 { return fSize; }
int fSize
Definition: LArFFT.h:77
void util::LArFFT::InitializeFFT ( )
private

Definition at line 52 of file LArFFT_service.cc.

53 {
54  int i;
55  for (i = 1; i < fSize; i *= 2) {}
56  fSize = i;
57  fFreqSize = fSize / 2 + 1;
58 
59  // allocate and setup Transform objects
60  fFFT = new TFFTRealComplex(fSize, false);
61  fInverseFFT = new TFFTComplexReal(fSize, false);
62 
63  int dummy[1] = {0};
64  // appears to be dummy argument from root page
65  fFFT->Init(fOption.c_str(), -1, dummy);
66  fInverseFFT->Init(fOption.c_str(), 1, dummy);
67 
68  fPeakFit = new TF1("fPeakFit", "gaus"); //allocate function used for peak fitting
69  fConvHist = new TH1D("fConvHist",
70  "Convolution Peak Data",
71  fFitBins,
72  0,
73  fFitBins); //allocate histogram for peak fitting
74  //allocate other data vectors
75  fCompTemp.resize(fFreqSize);
76  fKern.resize(fFreqSize);
77 }
int fSize
Definition: LArFFT.h:77
std::vector< TComplex > fKern
Definition: LArFFT.h:84
TFFTRealComplex * fFFT
object to do FFT
Definition: LArFFT.h:86
TF1 * fPeakFit
Definition: LArFFT.h:81
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83
cet::LibraryManager dummy("noplugin")
std::string fOption
Definition: LArFFT.h:79
TFFTComplexReal * fInverseFFT
object to do Inverse FF
Definition: LArFFT.h:87
TH1D * fConvHist
Definition: LArFFT.h:82
int fFitBins
Definition: LArFFT.h:80
int fFreqSize
Definition: LArFFT.h:78
template<class T >
T util::LArFFT::PeakCorrelation ( std::vector< T > &  shape1,
std::vector< T > &  shape2 
)
inline

Definition at line 272 of file LArFFT.h.

274 {
275  fConvHist->Reset("ICE");
276  std::vector<T> holder = shape1;
277  Correlate(holder,shape2);
278 
279  int maxT = max_element(holder.begin(), holder.end())-holder.begin();
280  float startT = maxT-fFitBins/2;
281  int offset = 0;
282 
283  for(int i = 0; i < fFitBins; i++) {
284  if(startT+i < 0) offset=fSize;
285  else if(startT+i > fSize) offset=-fSize;
286  else offset = 0;
287  fConvHist->Fill(i,holder[i+startT+offset]);
288  }
289 
290  fPeakFit->SetParameters(fConvHist->GetMaximum(),fFitBins/2,fFitBins/2);
291  fConvHist->Fit(fPeakFit,"QWNR","",0,fFitBins);
292  return fPeakFit->GetParameter(1)+startT;
293 }
int fSize
Definition: LArFFT.h:77
void Correlate(std::vector< T > &input, std::vector< T > &respFunc)
Definition: LArFFT.h:207
TF1 * fPeakFit
Definition: LArFFT.h:81
TH1D * fConvHist
Definition: LArFFT.h:82
int fFitBins
Definition: LArFFT.h:80
void util::LArFFT::ReinitializeFFT ( int  size,
std::string  option,
int  fitbins 
)

Definition at line 90 of file LArFFT_service.cc.

91 {
92  //delete these, which will be remade
93  delete fFFT;
94  delete fInverseFFT;
95  delete fPeakFit;
96  delete fConvHist;
97 
98  //set members
99  fSize = size;
100  fOption = option;
101  fFitBins = fitbins;
102 
103  //now initialize
104  InitializeFFT();
105 }
int fSize
Definition: LArFFT.h:77
TFFTRealComplex * fFFT
object to do FFT
Definition: LArFFT.h:86
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
Definition: StdUtils.h:92
TF1 * fPeakFit
Definition: LArFFT.h:81
std::string fOption
Definition: LArFFT.h:79
TFFTComplexReal * fInverseFFT
object to do Inverse FF
Definition: LArFFT.h:87
TH1D * fConvHist
Definition: LArFFT.h:82
int fFitBins
Definition: LArFFT.h:80
void InitializeFFT()
void util::LArFFT::resetSizePerRun ( art::Run const &  )
private

Definition at line 42 of file LArFFT_service.cc.

43 {
45  ->DataForJob()
46  .ReadOutWindowSize();
48 }
int fSize
Definition: LArFFT.h:77
std::string fOption
Definition: LArFFT.h:79
int fFitBins
Definition: LArFFT.h:80
void ReinitializeFFT(int, std::string, int)
void util::LArFFT::ShiftData ( std::vector< TComplex > &  input,
double  shift 
)

Definition at line 122 of file LArFFT_service.cc.

123 {
124  double factor = -2.0 * TMath::Pi() * shift / (double)fSize;
125 
126  for (int i = 0; i < fFreqSize; i++)
127  input[i] *= TComplex::Exp(TComplex(0, factor * (double)i));
128 
129  return;
130 }
int fSize
Definition: LArFFT.h:77
int fFreqSize
Definition: LArFFT.h:78
template<class T >
void util::LArFFT::ShiftData ( std::vector< T > &  input,
double  shift 
)
inline

Definition at line 259 of file LArFFT.h.

261 {
262  DoFFT(input,fCompTemp);
263  ShiftData(fCompTemp,shift);
264  DoInvFFT(fCompTemp,input);
265 
266  return;
267 }
void ShiftData(std::vector< TComplex > &input, double shift)
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:98
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:120
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:83

Member Data Documentation

std::vector<TComplex> util::LArFFT::fCompTemp
private

Definition at line 83 of file LArFFT.h.

TH1D* util::LArFFT::fConvHist
private

Definition at line 82 of file LArFFT.h.

TFFTRealComplex* util::LArFFT::fFFT
private

object to do FFT

Definition at line 86 of file LArFFT.h.

int util::LArFFT::fFitBins
private

Definition at line 80 of file LArFFT.h.

int util::LArFFT::fFreqSize
private

Definition at line 78 of file LArFFT.h.

TFFTComplexReal* util::LArFFT::fInverseFFT
private

object to do Inverse FF

Definition at line 87 of file LArFFT.h.

std::vector<TComplex> util::LArFFT::fKern
private

Definition at line 84 of file LArFFT.h.

std::string util::LArFFT::fOption
private

Definition at line 79 of file LArFFT.h.

TF1* util::LArFFT::fPeakFit
private

Definition at line 81 of file LArFFT.h.

int util::LArFFT::fSize
private

Definition at line 77 of file LArFFT.h.


The documentation for this class was generated from the following files: