LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
util::LArFFTW Class Reference

#include "LArFFTW.h"

Public Types

using FloatVector = std::vector< float >
 
using DoubleVector = std::vector< double >
 
using ComplexVector = std::vector< std::complex< double >>
 

Public Member Functions

 LArFFTW (int transformSize, const void *fplan, const void *rplan, int fitbins)
 
 ~LArFFTW ()
 
template<class T >
void DoFFT (std::vector< T > &input)
 
template<class T >
void DoFFT (std::vector< T > &input, ComplexVector &output)
 
template<class T >
void DoInvFFT (std::vector< T > &output)
 
template<class T >
void DoInvFFT (ComplexVector &input, std::vector< T > &output)
 
template<class T >
void Convolute (std::vector< T > &func, const ComplexVector &kern)
 
template<class T >
void Convolute (std::vector< T > &func, std::vector< T > &resp)
 
template<class T >
void Deconvolute (std::vector< T > &func, const ComplexVector &kern)
 
template<class T >
void Deconvolute (std::vector< T > &func, std::vector< T > &resp)
 
template<class T >
void Correlate (std::vector< T > &func, const ComplexVector &kern)
 
template<class T >
void Correlate (std::vector< T > &func, std::vector< T > &resp)
 
void ShiftData (ComplexVector &input, double shift)
 
template<class T >
void ShiftData (std::vector< T > &input, double shift)
 
template<class T >
void AlignedSum (std::vector< T > &input, std::vector< T > &output, bool add=true)
 
template<class T >
PeakCorrelation (std::vector< T > &shape1, std::vector< T > &shape2)
 

Private Attributes

ComplexVector fKern
 
ComplexVector fCompTemp
 
std::vector< float > fConvHist
 
int fSize
 
int fFreqSize
 
void * fIn
 
void * fOut
 
const void * fPlan
 
void * rIn
 
void * rOut
 
const void * rPlan
 
int fFitBins
 
gshf::MarqFitAlgfMarqFitAlg
 

Detailed Description

Definition at line 18 of file LArFFTW.h.

Member Typedef Documentation

using util::LArFFTW::ComplexVector = std::vector<std::complex<double>>

Definition at line 23 of file LArFFTW.h.

using util::LArFFTW::DoubleVector = std::vector<double>

Definition at line 22 of file LArFFTW.h.

using util::LArFFTW::FloatVector = std::vector<float>

Definition at line 21 of file LArFFTW.h.

Constructor & Destructor Documentation

util::LArFFTW::LArFFTW ( int  transformSize,
const void *  fplan,
const void *  rplan,
int  fitbins 
)

Definition at line 5 of file LArFFTW.cxx.

References fCompTemp, fConvHist, fFitBins, fFreqSize, fIn, fKern, fOut, fSize, rIn, and rOut.

6  : fSize(transformSize), fPlan(fplan), rPlan(rplan), fFitBins(fitbins)
7 {
8 
9  fFreqSize = fSize / 2 + 1;
10 
11  // ... Real-Complex
12  fIn = fftw_malloc(sizeof(double) * fSize);
13  fOut = fftw_malloc(sizeof(fftw_complex) * fFreqSize);
14 
15  // ... Complex-Real
16  rIn = fftw_malloc(sizeof(fftw_complex) * fFreqSize);
17  rOut = fftw_malloc(sizeof(double) * fSize);
18 
19  // ... allocate other data vectors
20  fCompTemp.resize(fFreqSize);
21  fKern.resize(fFreqSize);
22  fConvHist.resize(fFitBins);
23 }
ComplexVector fKern
Definition: LArFFTW.h:65
void * fIn
Definition: LArFFTW.h:70
int fFreqSize
Definition: LArFFTW.h:69
const void * fPlan
Definition: LArFFTW.h:72
void * rOut
Definition: LArFFTW.h:74
std::vector< float > fConvHist
Definition: LArFFTW.h:67
int fFitBins
Definition: LArFFTW.h:76
void * rIn
Definition: LArFFTW.h:73
const void * rPlan
Definition: LArFFTW.h:75
void * fOut
Definition: LArFFTW.h:71
ComplexVector fCompTemp
Definition: LArFFTW.h:66
util::LArFFTW::~LArFFTW ( )

Definition at line 25 of file LArFFTW.cxx.

References fIn, fOut, fPlan, rIn, rOut, and rPlan.

26 {
27  fPlan = 0;
28  fftw_free(fIn);
29  fIn = 0;
30  fftw_free((fftw_complex*)fOut);
31  fOut = 0;
32 
33  rPlan = 0;
34  fftw_free((fftw_complex*)rIn);
35  rIn = 0;
36  fftw_free(rOut);
37  rOut = 0;
38 }
void * fIn
Definition: LArFFTW.h:70
const void * fPlan
Definition: LArFFTW.h:72
void * rOut
Definition: LArFFTW.h:74
void * rIn
Definition: LArFFTW.h:73
const void * rPlan
Definition: LArFFTW.h:75
void * fOut
Definition: LArFFTW.h:71

Member Function Documentation

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

Definition at line 364 of file LArFFTW.h.

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

365 {
366  double shift = PeakCorrelation(shape1, shape2);
367 
368  ShiftData(shape1, shift);
369 
370  if (add)
371  for (int i = 0; i < fSize; i++)
372  shape1[i] += shape2[i];
373 
374  return;
375 }
T PeakCorrelation(std::vector< T > &shape1, std::vector< T > &shape2)
Definition: LArFFTW.h:382
void ShiftData(ComplexVector &input, double shift)
Definition: LArFFTW.cxx:43
template<class T >
void util::LArFFTW::Convolute ( std::vector< T > &  func,
const ComplexVector kern 
)
inline

Definition at line 170 of file LArFFTW.h.

References DoFFT(), DoInvFFT(), fFreqSize, fOut, fSize, n, and rIn.

171 {
172 
173  // ... Make sure that time series and kernel have the correct size.
174  int n = func.size();
175  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad time series size = " << n << "\n"; }
176  n = kern.size();
177  if (n != fFreqSize) { throw cet::exception("LArFFTW") << "Bad kernel size = " << n << "\n"; }
178 
179  DoFFT(func);
180 
181  // ..perform the convolution
182  for (int i = 0; i < fFreqSize; ++i) {
183  double re = ((fftw_complex*)fOut)[i][0];
184  double im = ((fftw_complex*)fOut)[i][1];
185  ((fftw_complex*)rIn)[i][0] = re * kern[i].real() - im * kern[i].imag();
186  ((fftw_complex*)rIn)[i][1] = re * kern[i].imag() + im * kern[i].real();
187  }
188 
189  DoInvFFT(func);
190 }
int fFreqSize
Definition: LArFFTW.h:69
void * rIn
Definition: LArFFTW.h:73
void * fOut
Definition: LArFFTW.h:71
Char_t n[5]
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<class T >
void util::LArFFTW::Convolute ( std::vector< T > &  func,
std::vector< T > &  resp 
)
inline

Definition at line 196 of file LArFFTW.h.

References DoFFT(), DoInvFFT(), fFreqSize, fKern, fOut, fSize, n, and rIn.

197 {
198 
199  // ... Make sure that time series has the correct size.
200  int n = func1.size();
201  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad 1st time series size = " << n << "\n"; }
202  n = func2.size();
203  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad 2nd time series size = " << n << "\n"; }
204 
205  DoFFT(func2);
206  for (int i = 0; i < fFreqSize; ++i) {
207  fKern[i].real(((fftw_complex*)fOut)[i][0]);
208  fKern[i].imag(((fftw_complex*)fOut)[i][1]);
209  }
210  DoFFT(func1);
211 
212  // ..perform the convolution
213  for (int i = 0; i < fFreqSize; ++i) {
214  double re = ((fftw_complex*)fOut)[i][0];
215  double im = ((fftw_complex*)fOut)[i][1];
216  ((fftw_complex*)rIn)[i][0] = re * fKern[i].real() - im * fKern[i].imag();
217  ((fftw_complex*)rIn)[i][1] = re * fKern[i].imag() + im * fKern[i].real();
218  }
219 
220  DoInvFFT(func1);
221 }
ComplexVector fKern
Definition: LArFFTW.h:65
int fFreqSize
Definition: LArFFTW.h:69
int func2(int i)
Definition: XFunc.cc:50
void * rIn
Definition: LArFFTW.h:73
int func1(int i)
Definition: XFunc.cc:39
void * fOut
Definition: LArFFTW.h:71
Char_t n[5]
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<class T >
void util::LArFFTW::Correlate ( std::vector< T > &  func,
const ComplexVector kern 
)
inline

Definition at line 292 of file LArFFTW.h.

References DoFFT(), DoInvFFT(), fFreqSize, fOut, fSize, n, and rIn.

Referenced by PeakCorrelation().

293 {
294 
295  // ... Make sure that time series and kernel have the correct size.
296  int n = func.size();
297  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad time series size = " << n << "\n"; }
298  n = kern.size();
299  if (n != fFreqSize) { throw cet::exception("LArFFTW") << "Bad kernel size = " << n << "\n"; }
300 
301  DoFFT(func);
302 
303  // ..perform the correlation
304  for (int i = 0; i < fFreqSize; ++i) {
305  double re = ((fftw_complex*)fOut)[i][0];
306  double im = ((fftw_complex*)fOut)[i][1];
307  ((fftw_complex*)rIn)[i][0] = re * kern[i].real() + im * kern[i].imag();
308  ((fftw_complex*)rIn)[i][1] = -re * kern[i].imag() + im * kern[i].real();
309  }
310 
311  DoInvFFT(func);
312 }
int fFreqSize
Definition: LArFFTW.h:69
void * rIn
Definition: LArFFTW.h:73
void * fOut
Definition: LArFFTW.h:71
Char_t n[5]
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<class T >
void util::LArFFTW::Correlate ( std::vector< T > &  func,
std::vector< T > &  resp 
)
inline

Definition at line 318 of file LArFFTW.h.

References DoFFT(), DoInvFFT(), fFreqSize, fKern, fOut, fSize, n, and rIn.

319 {
320 
321  // ... Make sure that time series has the correct size.
322  int n = func1.size();
323  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad 1st time series size = " << n << "\n"; }
324  n = func2.size();
325  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad 2nd time series size = " << n << "\n"; }
326 
327  DoFFT(func2);
328  for (int i = 0; i < fFreqSize; ++i) {
329  fKern[i].real(((fftw_complex*)fOut)[i][0]);
330  fKern[i].imag(((fftw_complex*)fOut)[i][1]);
331  }
332  DoFFT(func1);
333 
334  // ..perform the correlation
335  for (int i = 0; i < fFreqSize; ++i) {
336  double re = ((fftw_complex*)fOut)[i][0];
337  double im = ((fftw_complex*)fOut)[i][1];
338  ((fftw_complex*)rIn)[i][0] = re * fKern[i].real() + im * fKern[i].imag();
339  ((fftw_complex*)rIn)[i][1] = -re * fKern[i].imag() + im * fKern[i].real();
340  }
341 
342  DoInvFFT(func1);
343 }
ComplexVector fKern
Definition: LArFFTW.h:65
int fFreqSize
Definition: LArFFTW.h:69
int func2(int i)
Definition: XFunc.cc:50
void * rIn
Definition: LArFFTW.h:73
int func1(int i)
Definition: XFunc.cc:39
void * fOut
Definition: LArFFTW.h:71
Char_t n[5]
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<class T >
void util::LArFFTW::Deconvolute ( std::vector< T > &  func,
const ComplexVector kern 
)
inline

Definition at line 227 of file LArFFTW.h.

References d, DoFFT(), DoInvFFT(), e, fFreqSize, fOut, fSize, n, and rIn.

228 {
229 
230  // ... Make sure that time series and kernel have the correct size.
231  int n = func.size();
232  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad time series size = " << n << "\n"; }
233  n = kern.size();
234  if (n != fFreqSize) { throw cet::exception("LArFFTW") << "Bad kernel size = " << n << "\n"; }
235 
236  DoFFT(func);
237 
238  // ..perform the deconvolution
239  double a, b, c, d, e;
240  for (int i = 0; i < fFreqSize; ++i) {
241  a = ((fftw_complex*)fOut)[i][0];
242  b = ((fftw_complex*)fOut)[i][1];
243  c = kern[i].real();
244  d = kern[i].imag();
245  e = 1. / (c * c + d * d);
246  ((fftw_complex*)rIn)[i][0] = (a * c + b * d) * e;
247  ((fftw_complex*)rIn)[i][1] = (b * c - a * d) * e;
248  }
249 
250  DoInvFFT(func);
251 }
int fFreqSize
Definition: LArFFTW.h:69
void * rIn
Definition: LArFFTW.h:73
Float_t d
Definition: plot.C:235
void * fOut
Definition: LArFFTW.h:71
Char_t n[5]
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
Float_t e
Definition: plot.C:35
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<class T >
void util::LArFFTW::Deconvolute ( std::vector< T > &  func,
std::vector< T > &  resp 
)
inline

Definition at line 257 of file LArFFTW.h.

References d, DoFFT(), DoInvFFT(), e, fFreqSize, fKern, fOut, fSize, n, and rIn.

258 {
259 
260  // ... Make sure that time series has the correct size.
261  int n = func.size();
262  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad 1st time series size = " << n << "\n"; }
263  n = resp.size();
264  if (n != fSize) { throw cet::exception("LArFFTW") << "Bad 2nd time series size = " << n << "\n"; }
265 
266  DoFFT(resp);
267  for (int i = 0; i < fFreqSize; ++i) {
268  fKern[i].real(((fftw_complex*)fOut)[i][0]);
269  fKern[i].imag(((fftw_complex*)fOut)[i][1]);
270  }
271  DoFFT(func);
272 
273  // ..perform the deconvolution
274  double a, b, c, d, e;
275  for (int i = 0; i < fFreqSize; ++i) {
276  a = ((fftw_complex*)fOut)[i][0];
277  b = ((fftw_complex*)fOut)[i][1];
278  c = fKern[i].real();
279  d = fKern[i].imag();
280  e = 1. / (c * c + d * d);
281  ((fftw_complex*)rIn)[i][0] = (a * c + b * d) * e;
282  ((fftw_complex*)rIn)[i][1] = (b * c - a * d) * e;
283  }
284 
285  DoInvFFT(func);
286 }
ComplexVector fKern
Definition: LArFFTW.h:65
int fFreqSize
Definition: LArFFTW.h:69
void * rIn
Definition: LArFFTW.h:73
Float_t d
Definition: plot.C:235
void * fOut
Definition: LArFFTW.h:71
Char_t n[5]
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
Float_t e
Definition: plot.C:35
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
template<class T >
void util::LArFFTW::DoFFT ( std::vector< T > &  input)
inline

Definition at line 87 of file LArFFTW.h.

References fIn, fOut, and fPlan.

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

88 {
89  // ..set point
90  for (size_t p = 0; p < input.size(); ++p) {
91  ((double*)fIn)[p] = input[p];
92  }
93 
94  // ..transform (using the New-array Execute Functions)
95  fftw_execute_dft_r2c((fftw_plan)fPlan, (double*)fIn, (fftw_complex*)fOut);
96 
97  return;
98 }
void * fIn
Definition: LArFFTW.h:70
const void * fPlan
Definition: LArFFTW.h:72
void * fOut
Definition: LArFFTW.h:71
template<class T >
void util::LArFFTW::DoFFT ( std::vector< T > &  input,
ComplexVector output 
)
inline

Definition at line 104 of file LArFFTW.h.

References fFreqSize, fIn, fOut, and fPlan.

105 {
106  // ..set point
107  for (size_t p = 0; p < input.size(); ++p) {
108  ((double*)fIn)[p] = input[p];
109  }
110 
111  // ..transform (using the New-array Execute Functions)
112  fftw_execute_dft_r2c((fftw_plan)fPlan, (double*)fIn, (fftw_complex*)fOut);
113 
114  for (int i = 0; i < fFreqSize; ++i) {
115  output[i].real(((fftw_complex*)fOut)[i][0]);
116  output[i].imag(((fftw_complex*)fOut)[i][1]);
117  }
118 
119  return;
120 }
void * fIn
Definition: LArFFTW.h:70
int fFreqSize
Definition: LArFFTW.h:69
const void * fPlan
Definition: LArFFTW.h:72
void * fOut
Definition: LArFFTW.h:71
template<class T >
void util::LArFFTW::DoInvFFT ( std::vector< T > &  output)
inline

Definition at line 126 of file LArFFTW.h.

References lar::dump::array(), fSize, rIn, rOut, and rPlan.

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

127 {
128  // ..transform (using the New-array Execute Functions)
129  fftw_execute_dft_c2r((fftw_plan)rPlan, (fftw_complex*)rIn, (double*)rOut);
130 
131  // ..get point real
132  double factor = 1.0 / (double)fSize;
133  const double* array = (const double*)(rOut);
134  for (int i = 0; i < fSize; ++i) {
135  output[i] = factor * array[i];
136  }
137 
138  return;
139 }
void * rOut
Definition: LArFFTW.h:74
auto array(Array const &a)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:250
void * rIn
Definition: LArFFTW.h:73
const void * rPlan
Definition: LArFFTW.h:75
template<class T >
void util::LArFFTW::DoInvFFT ( ComplexVector input,
std::vector< T > &  output 
)
inline

Definition at line 145 of file LArFFTW.h.

References lar::dump::array(), fFreqSize, fSize, rIn, rOut, and rPlan.

146 {
147  // ..set point complex
148  for (int i = 0; i < fFreqSize; ++i) {
149  ((fftw_complex*)rIn)[i][0] = input[i].real();
150  ((fftw_complex*)rIn)[i][1] = input[i].imag();
151  }
152 
153  // ..transform (using the New-array Execute Functions)
154  fftw_execute_dft_c2r((fftw_plan)rPlan, (fftw_complex*)rIn, (double*)rOut);
155 
156  // ..get point real
157  double factor = 1.0 / (double)fSize;
158  const double* array = (const double*)(rOut);
159  for (int i = 0; i < fSize; ++i) {
160  output[i] = factor * array[i];
161  }
162 
163  return;
164 }
int fFreqSize
Definition: LArFFTW.h:69
void * rOut
Definition: LArFFTW.h:74
auto array(Array const &a)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:250
void * rIn
Definition: LArFFTW.h:73
const void * rPlan
Definition: LArFFTW.h:75
template<class T >
T util::LArFFTW::PeakCorrelation ( std::vector< T > &  shape1,
std::vector< T > &  shape2 
)
inline

Definition at line 382 of file LArFFTW.h.

References Correlate(), e, fConvHist, fFitBins, fMarqFitAlg, fSize, and gshf::MarqFitAlg::mrqdtfit().

Referenced by AlignedSum().

383 {
384  float chiSqr = std::numeric_limits<float>::max();
385  float dchiSqr = std::numeric_limits<float>::max();
386  const float chiCut = 1e-3;
387  float lambda = 0.001; // Marquardt damping parameter
388  std::vector<float> p;
389 
390  std::vector<T> holder = shape1;
391  Correlate(holder, shape2);
392 
393  int maxT = max_element(holder.begin(), holder.end()) - holder.begin();
394  float startT = maxT - fFitBins / 2;
395  int offset = 0;
396 
397  for (int i = 0; i < fFitBins; i++) {
398  if (startT + i < 0)
399  offset = fSize;
400  else if (startT + i > fSize)
401  offset = -fSize;
402  else
403  offset = 0;
404  if (holder[i + startT + offset] <= 0.) { fConvHist[i] = 0.; }
405  else {
406  fConvHist[i] = holder[i + startT + offset];
407  }
408  }
409 
410  p[0] = *max_element(fConvHist.begin(), fConvHist.end());
411  p[1] = fFitBins / 2;
412  p[2] = fFitBins / 2;
413  float p1 = p[1]; // save initial p[1] guess
414 
415  int fitResult{-1};
416  int trial = 0;
417  lambda = -1.; // initialize lambda on first call
418  do {
419  fitResult = fMarqFitAlg->mrqdtfit(lambda, &p[0], &fConvHist[0], 3, fFitBins, chiSqr, dchiSqr);
420  trial++;
421  if (fitResult) {
422  mf::LogWarning("LArFFTW") << "Peak Correlation Fitting failed";
423  break;
424  }
425  else if (trial > 100) {
426  break;
427  }
428  } while (fabs(dchiSqr) >= chiCut);
429  if (!fitResult) p1 = p[1]; // if fit succeeded, use fit result
430 
431  return p1 + 0.5 + startT;
432 }
int mrqdtfit(float &lambda, float p[], float y[], const int nParam, const int nData, float &chiSqr, float &dchiSqr)
Definition: MarqFitAlg.cxx:294
std::vector< float > fConvHist
Definition: LArFFTW.h:67
int fFitBins
Definition: LArFFTW.h:76
void Correlate(std::vector< T > &func, const ComplexVector &kern)
Definition: LArFFTW.h:292
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
gshf::MarqFitAlg * fMarqFitAlg
Definition: LArFFTW.h:78
Float_t e
Definition: plot.C:35
void util::LArFFTW::ShiftData ( ComplexVector input,
double  shift 
)

Definition at line 43 of file LArFFTW.cxx.

References fFreqSize, and fSize.

Referenced by AlignedSum(), and ShiftData().

44 {
45  double factor = -2.0 * std::acos(-1) * shift / (double)fSize;
46 
47  for (int i = 0; i < fFreqSize; i++) {
48  input[i] *= std::exp(std::complex<double>(0, factor * (double)i));
49  }
50 
51  return;
52 }
int fFreqSize
Definition: LArFFTW.h:69
template<class T >
void util::LArFFTW::ShiftData ( std::vector< T > &  input,
double  shift 
)
inline

Definition at line 349 of file LArFFTW.h.

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

350 {
351  DoFFT(input, fCompTemp);
352  ShiftData(fCompTemp, shift);
353  DoInvFFT(fCompTemp, input);
354 
355  return;
356 }
void ShiftData(ComplexVector &input, double shift)
Definition: LArFFTW.cxx:43
ComplexVector fCompTemp
Definition: LArFFTW.h:66
void DoInvFFT(std::vector< T > &output)
Definition: LArFFTW.h:126
void DoFFT(std::vector< T > &input)
Definition: LArFFTW.h:87

Member Data Documentation

ComplexVector util::LArFFTW::fCompTemp
private

Definition at line 66 of file LArFFTW.h.

Referenced by LArFFTW(), and ShiftData().

std::vector<float> util::LArFFTW::fConvHist
private

Definition at line 67 of file LArFFTW.h.

Referenced by LArFFTW(), and PeakCorrelation().

int util::LArFFTW::fFitBins
private

Definition at line 76 of file LArFFTW.h.

Referenced by LArFFTW(), and PeakCorrelation().

int util::LArFFTW::fFreqSize
private

Definition at line 69 of file LArFFTW.h.

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

void* util::LArFFTW::fIn
private

Definition at line 70 of file LArFFTW.h.

Referenced by DoFFT(), LArFFTW(), and ~LArFFTW().

ComplexVector util::LArFFTW::fKern
private

Definition at line 65 of file LArFFTW.h.

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

gshf::MarqFitAlg* util::LArFFTW::fMarqFitAlg
private

Definition at line 78 of file LArFFTW.h.

Referenced by PeakCorrelation().

void* util::LArFFTW::fOut
private

Definition at line 71 of file LArFFTW.h.

Referenced by Convolute(), Correlate(), Deconvolute(), DoFFT(), LArFFTW(), and ~LArFFTW().

const void* util::LArFFTW::fPlan
private

Definition at line 72 of file LArFFTW.h.

Referenced by DoFFT(), and ~LArFFTW().

int util::LArFFTW::fSize
private
void* util::LArFFTW::rIn
private

Definition at line 73 of file LArFFTW.h.

Referenced by Convolute(), Correlate(), Deconvolute(), DoInvFFT(), LArFFTW(), and ~LArFFTW().

void* util::LArFFTW::rOut
private

Definition at line 74 of file LArFFTW.h.

Referenced by DoInvFFT(), LArFFTW(), and ~LArFFTW().

const void* util::LArFFTW::rPlan
private

Definition at line 75 of file LArFFTW.h.

Referenced by DoInvFFT(), and ~LArFFTW().


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