#include <LArFFTW.h>
|
| 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 > |
T | PeakCorrelation (std::vector< T > &shape1, std::vector< T > &shape2) |
|
Definition at line 18 of file LArFFTW.h.
util::LArFFTW::LArFFTW |
( |
int |
transformSize, |
|
|
const void * |
fplan, |
|
|
const void * |
rplan, |
|
|
int |
fitbins |
|
) |
| |
Definition at line 5 of file LArFFTW.cxx.
6 :
fSize (transformSize)
15 fIn = fftw_malloc(
sizeof(
double)*
fSize);
std::vector< float > fConvHist
util::LArFFTW::~LArFFTW |
( |
| ) |
|
Definition at line 28 of file LArFFTW.cxx.
33 fftw_free((fftw_complex*)
fOut);
37 fftw_free((fftw_complex*)
rIn);
template<class T >
void util::LArFFTW::AlignedSum |
( |
std::vector< T > & |
input, |
|
|
std::vector< T > & |
output, |
|
|
bool |
add = true |
|
) |
| |
|
inline |
Definition at line 376 of file LArFFTW.h.
384 if(add)
for(
int i = 0; i <
fSize; i++) shape1[i]+=shape2[i];
T PeakCorrelation(std::vector< T > &shape1, std::vector< T > &shape2)
void ShiftData(ComplexVector &input, double shift)
template<class T >
void util::LArFFTW::Convolute |
( |
std::vector< T > & |
func, |
|
|
const ComplexVector & |
kern |
|
) |
| |
|
inline |
Definition at line 156 of file LArFFTW.h.
162 throw cet::exception(
"LArFFTW") <<
"Bad time series size = " << n <<
"\n";
166 throw cet::exception(
"LArFFTW") <<
"Bad kernel size = " << n <<
"\n";
173 double re = ((fftw_complex*)
fOut)[i][0];
174 double im = ((fftw_complex*)
fOut)[i][1];
175 ((fftw_complex*)
rIn)[i][0] = re*kern[i].real()-im*kern[i].imag();
176 ((fftw_complex*)
rIn)[i][1] = re*kern[i].imag()+im*kern[i].real();
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
template<class T >
void util::LArFFTW::Convolute |
( |
std::vector< T > & |
func, |
|
|
std::vector< T > & |
resp |
|
) |
| |
|
inline |
Definition at line 186 of file LArFFTW.h.
190 int n = func1.size();
192 throw cet::exception(
"LArFFTW") <<
"Bad 1st time series size = " << n <<
"\n";
196 throw cet::exception(
"LArFFTW") <<
"Bad 2nd time series size = " << n <<
"\n";
201 fKern[i].real(((fftw_complex*)
fOut)[i][0]);
202 fKern[i].imag(((fftw_complex*)
fOut)[i][1]);
208 double re = ((fftw_complex*)
fOut)[i][0];
209 double im = ((fftw_complex*)
fOut)[i][1];
210 ((fftw_complex*)
rIn)[i][0] = re*
fKern[i].real()-im*
fKern[i].imag();
211 ((fftw_complex*)
rIn)[i][1] = re*
fKern[i].imag()+im*
fKern[i].real();
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
template<class T >
void util::LArFFTW::Correlate |
( |
std::vector< T > & |
func, |
|
|
const ComplexVector & |
kern |
|
) |
| |
|
inline |
Definition at line 295 of file LArFFTW.h.
301 throw cet::exception(
"LArFFTW") <<
"Bad time series size = " << n <<
"\n";
305 throw cet::exception(
"LArFFTW") <<
"Bad kernel size = " << n <<
"\n";
312 double re = ((fftw_complex*)
fOut)[i][0];
313 double im = ((fftw_complex*)
fOut)[i][1];
314 ((fftw_complex*)
rIn)[i][0] = re*kern[i].real()+im*kern[i].imag();
315 ((fftw_complex*)
rIn)[i][1] = -re*kern[i].imag()+im*kern[i].real();
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
template<class T >
void util::LArFFTW::Correlate |
( |
std::vector< T > & |
func, |
|
|
std::vector< T > & |
resp |
|
) |
| |
|
inline |
Definition at line 326 of file LArFFTW.h.
330 int n = func1.size();
332 throw cet::exception(
"LArFFTW") <<
"Bad 1st time series size = " << n <<
"\n";
336 throw cet::exception(
"LArFFTW") <<
"Bad 2nd time series size = " << n <<
"\n";
341 fKern[i].real(((fftw_complex*)
fOut)[i][0]);
342 fKern[i].imag(((fftw_complex*)
fOut)[i][1]);
348 double re = ((fftw_complex*)
fOut)[i][0];
349 double im = ((fftw_complex*)
fOut)[i][1];
350 ((fftw_complex*)
rIn)[i][0] = re*
fKern[i].real()+im*
fKern[i].imag();
351 ((fftw_complex*)
rIn)[i][1] = -re*
fKern[i].imag()+im*
fKern[i].real();
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
template<class T >
void util::LArFFTW::Deconvolute |
( |
std::vector< T > & |
func, |
|
|
const ComplexVector & |
kern |
|
) |
| |
|
inline |
Definition at line 221 of file LArFFTW.h.
227 throw cet::exception(
"LArFFTW") <<
"Bad time series size = " << n <<
"\n";
231 throw cet::exception(
"LArFFTW") <<
"Bad kernel size = " << n <<
"\n";
239 a = ((fftw_complex*)
fOut)[i][0];
240 b = ((fftw_complex*)
fOut)[i][1];
244 ((fftw_complex*)
rIn)[i][0] = (a*c+b*d)*e;
245 ((fftw_complex*)
rIn)[i][1] = (b*c-a*d)*e;
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
template<class T >
void util::LArFFTW::Deconvolute |
( |
std::vector< T > & |
func, |
|
|
std::vector< T > & |
resp |
|
) |
| |
|
inline |
Definition at line 255 of file LArFFTW.h.
261 throw cet::exception(
"LArFFTW") <<
"Bad 1st time series size = " << n <<
"\n";
265 throw cet::exception(
"LArFFTW") <<
"Bad 2nd time series size = " << n <<
"\n";
270 fKern[i].real(((fftw_complex*)
fOut)[i][0]);
271 fKern[i].imag(((fftw_complex*)
fOut)[i][1]);
278 a = ((fftw_complex*)
fOut)[i][0];
279 b = ((fftw_complex*)
fOut)[i][1];
283 ((fftw_complex*)
rIn)[i][0] = (a*c+b*d)*e;
284 ((fftw_complex*)
rIn)[i][1] = (b*c-a*d)*e;
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
template<class T >
void util::LArFFTW::DoFFT |
( |
std::vector< T > & |
input | ) |
|
|
inline |
Definition at line 76 of file LArFFTW.h.
79 for(
size_t p = 0;
p < input.size(); ++
p){
80 ((
double *)
fIn)[
p] = input[
p];
84 fftw_execute_dft_r2c((fftw_plan)
fPlan,(
double*)
fIn,(fftw_complex*)
fOut);
template<class T >
void util::LArFFTW::DoFFT |
( |
std::vector< T > & |
input, |
|
|
ComplexVector & |
output |
|
) |
| |
|
inline |
Definition at line 92 of file LArFFTW.h.
95 for(
size_t p = 0;
p < input.size(); ++
p){
96 ((
double *)
fIn)[
p] = input[
p];
100 fftw_execute_dft_r2c((fftw_plan)
fPlan,(
double*)
fIn,(fftw_complex*)
fOut);
BEGIN_PROLOG sequence::SlidingWindowTriggerPatternsOppositeWindows END_PROLOG simSlidingORM6O6 effSlidingORW output
template<class T >
void util::LArFFTW::DoInvFFT |
( |
std::vector< T > & |
output | ) |
|
|
inline |
Definition at line 113 of file LArFFTW.h.
116 fftw_execute_dft_c2r((fftw_plan)
rPlan,(fftw_complex*)
rIn,(
double*)
rOut);
119 double factor = 1.0/(double)
fSize;
120 const double *
array = (
const double*)(
rOut);
121 for(
int i = 0; i <
fSize; ++i){
122 output[i] = factor*array[i];
then echo ***************************************echo array
template<class T >
void util::LArFFTW::DoInvFFT |
( |
ComplexVector & |
input, |
|
|
std::vector< T > & |
output |
|
) |
| |
|
inline |
Definition at line 131 of file LArFFTW.h.
135 ((fftw_complex*)
rIn)[i][0] = input[i].real();
136 ((fftw_complex*)
rIn)[i][1] = input[i].imag();
140 fftw_execute_dft_c2r((fftw_plan)
rPlan,(fftw_complex*)
rIn,(
double*)
rOut);
143 double factor = 1.0/(double)
fSize;
144 const double *
array = (
const double*)(
rOut);
145 for(
int i = 0; i <
fSize; ++i){
146 output[i] = factor*array[i];
then echo ***************************************echo array
template<class T >
T util::LArFFTW::PeakCorrelation |
( |
std::vector< T > & |
shape1, |
|
|
std::vector< T > & |
shape2 |
|
) |
| |
|
inline |
Definition at line 393 of file LArFFTW.h.
396 float chiSqr = std::numeric_limits<float>::max();
397 float dchiSqr = std::numeric_limits<float>::max();
398 const float chiCut = 1
e-3;
400 std::vector<float>
p;
402 std::vector<T> holder = shape1;
405 int maxT = max_element(holder.begin(), holder.end())-holder.begin();
410 if(startT+i < 0) offset=
fSize;
413 if(holder[i+startT+offset]<=0.) {
432 mf::LogWarning(
"LArFFTW") <<
"Peak Correlation Fitting failed";
434 }
else if (trial>100){
438 while (fabs(dchiSqr) >= chiCut);
439 if (!fitResult)p1=p[1];
441 return p1 + 0.5 + startT;
BEGIN_PROLOG TPC Trig offset(g4 rise time) ProjectToHeight
double lambda(double a, double b, double c)
std::vector< float > fConvHist
void Correlate(std::vector< T > &func, const ComplexVector &kern)
gshf::MarqFitAlg * fMarqFitAlg
physics associatedGroupsWithLeft p1
void util::LArFFTW::ShiftData |
( |
ComplexVector & |
input, |
|
|
double |
shift |
|
) |
| |
Definition at line 46 of file LArFFTW.cxx.
48 double factor = -2.0*std::acos(-1)*
shift/(double)
fSize;
51 input[i] *=
std::exp(std::complex<double>(0,factor*(
double)i));
template<class T >
void util::LArFFTW::ShiftData |
( |
std::vector< T > & |
input, |
|
|
double |
shift |
|
) |
| |
|
inline |
Definition at line 362 of file LArFFTW.h.
void ShiftData(ComplexVector &input, double shift)
void DoInvFFT(std::vector< T > &output)
void DoFFT(std::vector< T > &input)
std::vector<float> util::LArFFTW::fConvHist |
|
private |
int util::LArFFTW::fFitBins |
|
private |
int util::LArFFTW::fFreqSize |
|
private |
gshf::MarqFitAlg* util::LArFFTW::fMarqFitAlg |
|
private |
void* util::LArFFTW::fOut |
|
private |
const void* util::LArFFTW::fPlan |
|
private |
void* util::LArFFTW::rOut |
|
private |
const void* util::LArFFTW::rPlan |
|
private |
The documentation for this class was generated from the following files: