LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
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 18 of file LArFFT.cc.

References fSize, InitializeFFT(), resetSizePerRun(), and art::ActivityRegistry::sPreBeginRun.

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

Definition at line 76 of file LArFFT.cc.

References fConvHist, fFFT, fInverseFFT, and fPeakFit.

77 {
78  delete fFFT;
79  delete fInverseFFT;
80  delete fPeakFit;
81  delete fConvHist;
82 }
TFFTRealComplex * fFFT
object to do FFT
Definition: LArFFT.h:82
TF1 * fPeakFit
Definition: LArFFT.h:77
TFFTComplexReal * fInverseFFT
object to do Inverse FF
Definition: LArFFT.h:83
TH1D * fConvHist
Definition: LArFFT.h:78

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 241 of file LArFFT.h.

References fSize, PeakCorrelation(), and ShiftData().

Referenced by caldata::CalWireAna::analyze().

242 {
243  double shift = PeakCorrelation(shape1, shape2);
244 
245  ShiftData(shape1, shift);
246 
247  if (add)
248  for (int i = 0; i < fSize; i++)
249  shape1[i] += shape2[i];
250 
251  return;
252 }
void ShiftData(std::vector< TComplex > &input, double shift)
Definition: LArFFT.cc:116
int fSize
Definition: LArFFT.h:73
T PeakCorrelation(std::vector< T > &shape1, std::vector< T > &shape2)
Definition: LArFFT.h:270
template<class T >
void util::LArFFT::Convolute ( std::vector< T > &  input,
std::vector< T > &  respFunc 
)
inline

Definition at line 170 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, fFreqSize, and fKern.

Referenced by util::SignalShaping::Convolute(), util::SignalShaping::Deconvolute(), caldata::CalWire::produce(), and caldata::CalWireT962::produce().

171 {
172  DoFFT(shape1, fKern);
173  DoFFT(shape2, fCompTemp);
174 
175  for (int i = 0; i < fFreqSize; i++)
176  fCompTemp[i] *= fKern[i];
177 
178  DoInvFFT(fCompTemp, shape1);
179 
180  return;
181 }
std::vector< TComplex > fKern
Definition: LArFFT.h:80
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::Convolute ( std::vector< T > &  input,
std::vector< TComplex > &  kern 
)
inline

Definition at line 189 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, and fFreqSize.

190 {
191  DoFFT(input, fCompTemp);
192 
193  for (int i = 0; i < fFreqSize; i++)
194  fCompTemp[i] *= kern[i];
195 
196  DoInvFFT(fCompTemp, input);
197 
198  return;
199 }
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::Correlate ( std::vector< T > &  input,
std::vector< T > &  respFunc 
)
inline

Definition at line 204 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, fFreqSize, and fKern.

Referenced by PeakCorrelation().

205 {
206  DoFFT(shape1, fKern);
207  DoFFT(shape2, fCompTemp);
208 
209  for (int i = 0; i < fFreqSize; i++)
210  fCompTemp[i] *= TComplex::Conjugate(fKern[i]);
211 
212  DoInvFFT(fCompTemp, shape1);
213 
214  return;
215 }
std::vector< TComplex > fKern
Definition: LArFFT.h:80
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::Correlate ( std::vector< T > &  input,
std::vector< TComplex > &  kern 
)
inline

Definition at line 223 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, and fFreqSize.

224 {
225  DoFFT(input, fCompTemp);
226 
227  for (int i = 0; i < fFreqSize; i++)
228  fCompTemp[i] *= TComplex::Conjugate(kern[i]);
229 
230  DoInvFFT(fCompTemp, input);
231 
232  return;
233 }
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::Deconvolute ( std::vector< T > &  input,
std::vector< T > &  respFunc 
)
inline

Definition at line 135 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, fFreqSize, and fKern.

136 {
137  DoFFT(respFunction, fKern);
138  DoFFT(input, fCompTemp);
139 
140  for (int i = 0; i < fFreqSize; i++)
141  fCompTemp[i] /= fKern[i];
142 
143  DoInvFFT(fCompTemp, input);
144 
145  return;
146 }
std::vector< TComplex > fKern
Definition: LArFFT.h:80
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::Deconvolute ( std::vector< T > &  input,
std::vector< TComplex > &  kern 
)
inline

Definition at line 154 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, and fFreqSize.

155 {
156  DoFFT(input, fCompTemp);
157 
158  for (int i = 0; i < fFreqSize; i++)
159  fCompTemp[i] /= kern[i];
160 
161  DoInvFFT(fCompTemp, input);
162 
163  return;
164 }
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::DoFFT ( std::vector< T > &  input,
std::vector< TComplex > &  output 
)
inline

Definition at line 95 of file LArFFT.h.

References fFFT, and fFreqSize.

Referenced by util::SignalShaping::AddResponseFunction(), caldata::CalWireAna::analyze(), Convolute(), detsim::SimWire::ConvoluteResponseFunctions(), Correlate(), Deconvolute(), and ShiftData().

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

Definition at line 117 of file LArFFT.h.

References fFreqSize, fInverseFFT, and fSize.

Referenced by util::SignalShaping::AddResponseFunction(), util::SignalShaping::CalculateDeconvKernel(), Convolute(), detsim::SimWire::ConvoluteResponseFunctions(), Correlate(), Deconvolute(), detsim::SimWire::GenNoise(), ShiftData(), and util::SignalShaping::ShiftResponseTime().

118 {
119  for (int i = 0; i < fFreqSize; ++i)
120  fInverseFFT->SetPointComplex(i, input[i]);
121 
122  fInverseFFT->Transform();
123  double factor = 1.0 / (double)fSize;
124 
125  for (int i = 0; i < fSize; ++i)
126  output[i] = factor * fInverseFFT->GetPointReal(i, false);
127 
128  return;
129 }
int fSize
Definition: LArFFT.h:73
TFFTComplexReal * fInverseFFT
object to do Inverse FF
Definition: LArFFT.h:83
int fFreqSize
Definition: LArFFT.h:74
int util::LArFFT::FFTFitBins ( ) const
inline

Definition at line 68 of file LArFFT.h.

References fFitBins, and ReinitializeFFT().

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

Definition at line 67 of file LArFFT.h.

References fOption.

67 { return fOption; }
std::string fOption
Definition: LArFFT.h:75
void util::LArFFT::InitializeFFT ( )
private

Definition at line 48 of file LArFFT.cc.

References fCompTemp, fConvHist, fFFT, fFitBins, fFreqSize, fInverseFFT, fKern, fOption, fPeakFit, and fSize.

Referenced by LArFFT(), and ReinitializeFFT().

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

Definition at line 270 of file LArFFT.h.

References Correlate(), DECLARE_ART_SERVICE, fConvHist, fFitBins, fPeakFit, and fSize.

Referenced by AlignedSum(), detsim::SimWire::ConvoluteResponseFunctions(), and util::SignalShaping::SetPeakResponseTime().

271 {
272  fConvHist->Reset("ICE");
273  std::vector<T> holder = shape1;
274  Correlate(holder, shape2);
275 
276  int maxT = max_element(holder.begin(), holder.end()) - holder.begin();
277  float startT = maxT - fFitBins / 2;
278  int offset = 0;
279 
280  for (int i = 0; i < fFitBins; i++) {
281  if (startT + i < 0)
282  offset = fSize;
283  else if (startT + i > fSize)
284  offset = -fSize;
285  else
286  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:73
void Correlate(std::vector< T > &input, std::vector< T > &respFunc)
Definition: LArFFT.h:204
TF1 * fPeakFit
Definition: LArFFT.h:77
TH1D * fConvHist
Definition: LArFFT.h:78
int fFitBins
Definition: LArFFT.h:76
void util::LArFFT::ReinitializeFFT ( int  size,
std::string  option,
int  fitbins 
)

Definition at line 85 of file LArFFT.cc.

References fConvHist, fFFT, fFitBins, fInverseFFT, fOption, fPeakFit, fSize, InitializeFFT(), and util::size().

Referenced by FFTFitBins(), and resetSizePerRun().

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

Definition at line 39 of file LArFFT.cc.

References fFitBins, fOption, fSize, and ReinitializeFFT().

Referenced by LArFFT().

40 {
42  ->DataForJob()
43  .ReadOutWindowSize();
45 }
int fSize
Definition: LArFFT.h:73
std::string fOption
Definition: LArFFT.h:75
int fFitBins
Definition: LArFFT.h:76
void ReinitializeFFT(int, std::string, int)
Definition: LArFFT.cc:85
void util::LArFFT::ShiftData ( std::vector< TComplex > &  input,
double  shift 
)

Definition at line 116 of file LArFFT.cc.

References fFreqSize, and fSize.

Referenced by AlignedSum(), detsim::SimWire::ConvoluteResponseFunctions(), ShiftData(), and util::SignalShaping::ShiftResponseTime().

117 {
118  double factor = -2.0 * TMath::Pi() * shift / (double)fSize;
119 
120  for (int i = 0; i < fFreqSize; i++)
121  input[i] *= TComplex::Exp(TComplex(0, factor * (double)i));
122 
123  return;
124 }
int fSize
Definition: LArFFT.h:73
int fFreqSize
Definition: LArFFT.h:74
template<class T >
void util::LArFFT::ShiftData ( std::vector< T > &  input,
double  shift 
)
inline

Definition at line 257 of file LArFFT.h.

References DoFFT(), DoInvFFT(), fCompTemp, and ShiftData().

258 {
259  DoFFT(input, fCompTemp);
260  ShiftData(fCompTemp, shift);
261  DoInvFFT(fCompTemp, input);
262 
263  return;
264 }
void ShiftData(std::vector< TComplex > &input, double shift)
Definition: LArFFT.cc:116
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
Definition: LArFFT.h:95
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
Definition: LArFFT.h:117
std::vector< TComplex > fCompTemp
Definition: LArFFT.h:79

Member Data Documentation

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

Definition at line 79 of file LArFFT.h.

Referenced by Convolute(), Correlate(), Deconvolute(), InitializeFFT(), and ShiftData().

TH1D* util::LArFFT::fConvHist
private

Definition at line 78 of file LArFFT.h.

Referenced by InitializeFFT(), PeakCorrelation(), ReinitializeFFT(), and ~LArFFT().

TFFTRealComplex* util::LArFFT::fFFT
private

object to do FFT

Definition at line 82 of file LArFFT.h.

Referenced by DoFFT(), InitializeFFT(), ReinitializeFFT(), and ~LArFFT().

int util::LArFFT::fFitBins
private

Definition at line 76 of file LArFFT.h.

Referenced by FFTFitBins(), InitializeFFT(), PeakCorrelation(), ReinitializeFFT(), and resetSizePerRun().

int util::LArFFT::fFreqSize
private

Definition at line 74 of file LArFFT.h.

Referenced by Convolute(), Correlate(), Deconvolute(), DoFFT(), DoInvFFT(), InitializeFFT(), and ShiftData().

TFFTComplexReal* util::LArFFT::fInverseFFT
private

object to do Inverse FF

Definition at line 83 of file LArFFT.h.

Referenced by DoInvFFT(), InitializeFFT(), ReinitializeFFT(), and ~LArFFT().

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

Definition at line 80 of file LArFFT.h.

Referenced by Convolute(), Correlate(), Deconvolute(), and InitializeFFT().

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

Definition at line 75 of file LArFFT.h.

Referenced by FFTOptions(), InitializeFFT(), ReinitializeFFT(), and resetSizePerRun().

TF1* util::LArFFT::fPeakFit
private

Definition at line 77 of file LArFFT.h.

Referenced by InitializeFFT(), PeakCorrelation(), ReinitializeFFT(), and ~LArFFT().

int util::LArFFT::fSize
private

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