#include <LArFFT.h>
|
| LArFFT (fhicl::ParameterSet const &pset, art::ActivityRegistry ®) |
|
| ~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 > |
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) |
|
Definition at line 26 of file LArFFT.h.
util::LArFFT::LArFFT |
( |
fhicl::ParameterSet const & |
pset, |
|
|
art::ActivityRegistry & |
reg |
|
) |
| |
Definition at line 18 of file LArFFT.cc.
19 :
fSize(pset.get<
int>(
"FFTSize", 0))
20 ,
fOption(pset.get<std::string>(
"FFTOption"))
30 fSize = art::ServiceHandle<detinfo::DetectorPropertiesService const>()
void resetSizePerRun(art::Run const &)
util::LArFFT::~LArFFT |
( |
| ) |
|
Definition at line 78 of file LArFFT.cc.
TFFTRealComplex * fFFT
object to do FFT
TFFTComplexReal * fInverseFFT
object to do Inverse FF
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.
252 if(add)
for(
int i = 0; i <
fSize; i++) shape1[i]+=shape2[i];
void ShiftData(std::vector< TComplex > &input, double shift)
T PeakCorrelation(std::vector< T > &shape1, std::vector< T > &shape2)
template<class T >
void util::LArFFT::Convolute |
( |
std::vector< T > & |
input, |
|
|
std::vector< T > & |
respFunc |
|
) |
| |
|
inline |
Definition at line 173 of file LArFFT.h.
std::vector< TComplex > fKern
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
template<class T >
void util::LArFFT::Convolute |
( |
std::vector< T > & |
input, |
|
|
std::vector< TComplex > & |
kern |
|
) |
| |
|
inline |
Definition at line 192 of file LArFFT.h.
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
template<class T >
void util::LArFFT::Correlate |
( |
std::vector< T > & |
input, |
|
|
std::vector< T > & |
respFunc |
|
) |
| |
|
inline |
Definition at line 207 of file LArFFT.h.
std::vector< TComplex > fKern
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
template<class T >
void util::LArFFT::Correlate |
( |
std::vector< T > & |
input, |
|
|
std::vector< TComplex > & |
kern |
|
) |
| |
|
inline |
Definition at line 226 of file LArFFT.h.
232 fCompTemp[i]*=TComplex::Conjugate(kern[i]);
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
template<class T >
void util::LArFFT::Deconvolute |
( |
std::vector< T > & |
input, |
|
|
std::vector< T > & |
respFunc |
|
) |
| |
|
inline |
Definition at line 138 of file LArFFT.h.
std::vector< TComplex > fKern
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
template<class T >
void util::LArFFT::Deconvolute |
( |
std::vector< T > & |
input, |
|
|
std::vector< TComplex > & |
kern |
|
) |
| |
|
inline |
Definition at line 157 of file LArFFT.h.
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
template<class T >
void util::LArFFT::DoFFT |
( |
std::vector< T > & |
input, |
|
|
std::vector< TComplex > & |
output |
|
) |
| |
|
inline |
Definition at line 98 of file LArFFT.h.
102 double imaginary = 0.;
105 for(
size_t p = 0;
p < input.size(); ++
p)
106 fFFT->SetPoint(
p, input[
p]);
111 fFFT->GetPointComplex(i, real, imaginary);
112 output[i]=TComplex(real, imaginary);
TFFTRealComplex * fFFT
object to do FFT
template<class T >
void util::LArFFT::DoInvFFT |
( |
std::vector< TComplex > & |
input, |
|
|
std::vector< T > & |
output |
|
) |
| |
|
inline |
Definition at line 120 of file LArFFT.h.
127 double factor = 1.0/(double)
fSize;
129 for(
int i = 0; i <
fSize; ++i)
130 output[i] = factor*
fInverseFFT->GetPointReal(i,
false);
TFFTComplexReal * fInverseFFT
object to do Inverse FF
int util::LArFFT::FFTFitBins |
( |
| ) |
const |
|
inline |
std::string util::LArFFT::FFTOptions |
( |
| ) |
const |
|
inline |
int util::LArFFT::FFTSize |
( |
| ) |
const |
|
inline |
void util::LArFFT::InitializeFFT |
( |
| ) |
|
|
private |
Definition at line 50 of file LArFFT.cc.
53 for (i = 1; i <
fSize; i *= 2) {}
58 fFFT =
new TFFTRealComplex(fSize,
false);
66 fPeakFit =
new TF1(
"fPeakFit",
"gaus");
68 "Convolution Peak Data",
std::vector< TComplex > fKern
TFFTRealComplex * fFFT
object to do FFT
std::vector< TComplex > fCompTemp
TFFTComplexReal * fInverseFFT
object to do Inverse FF
template<class T >
T util::LArFFT::PeakCorrelation |
( |
std::vector< T > & |
shape1, |
|
|
std::vector< T > & |
shape2 |
|
) |
| |
|
inline |
Definition at line 272 of file LArFFT.h.
276 std::vector<T> holder = shape1;
279 int maxT = max_element(holder.begin(), holder.end())-holder.begin();
284 if(startT+i < 0) offset=
fSize;
287 fConvHist->Fill(i,holder[i+startT+offset]);
292 return fPeakFit->GetParameter(1)+startT;
BEGIN_PROLOG TPC Trig offset(g4 rise time) ProjectToHeight
void Correlate(std::vector< T > &input, std::vector< T > &respFunc)
void util::LArFFT::ReinitializeFFT |
( |
int |
size, |
|
|
std::string |
option, |
|
|
int |
fitbins |
|
) |
| |
Definition at line 88 of file LArFFT.cc.
TFFTRealComplex * fFFT
object to do FFT
decltype(auto) constexpr size(T &&obj)
ADL-aware version of std::size.
TFFTComplexReal * fInverseFFT
object to do Inverse FF
void util::LArFFT::resetSizePerRun |
( |
art::Run const & |
| ) |
|
|
private |
Definition at line 40 of file LArFFT.cc.
42 fSize = art::ServiceHandle<detinfo::DetectorPropertiesService const>()
void ReinitializeFFT(int, std::string, int)
void util::LArFFT::ShiftData |
( |
std::vector< TComplex > & |
input, |
|
|
double |
shift |
|
) |
| |
Definition at line 120 of file LArFFT.cc.
122 double factor = -2.0 * TMath::Pi() *
shift / (double)
fSize;
125 input[i] *= TComplex::Exp(TComplex(0, factor * (
double)i));
template<class T >
void util::LArFFT::ShiftData |
( |
std::vector< T > & |
input, |
|
|
double |
shift |
|
) |
| |
|
inline |
Definition at line 259 of file LArFFT.h.
void ShiftData(std::vector< TComplex > &input, double shift)
void DoFFT(std::vector< T > &input, std::vector< TComplex > &output)
void DoInvFFT(std::vector< TComplex > &input, std::vector< T > &output)
std::vector< TComplex > fCompTemp
std::vector<TComplex> util::LArFFT::fCompTemp |
|
private |
TH1D* util::LArFFT::fConvHist |
|
private |
TFFTRealComplex* util::LArFFT::fFFT |
|
private |
object to do FFT
Definition at line 86 of file LArFFT.h.
int util::LArFFT::fFitBins |
|
private |
int util::LArFFT::fFreqSize |
|
private |
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 |
std::string util::LArFFT::fOption |
|
private |
TF1* util::LArFFT::fPeakFit |
|
private |
The documentation for this class was generated from the following files: