All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
opdet::OpDeconvolutionAlgWiener Class Reference
Inheritance diagram for opdet::OpDeconvolutionAlgWiener:
opdet::OpDeconvolutionAlg

Public Member Functions

 OpDeconvolutionAlgWiener (fhicl::ParameterSet const &p)
 
 ~OpDeconvolutionAlgWiener ()
 
std::vector< raw::OpDetWaveformRunDeconvolution (std::vector< raw::OpDetWaveform > const &wfHandle) override
 

Private Member Functions

void ApplyExpoAvSmoothing (std::vector< double > &wf)
 
void ApplyUnAvSmoothing (std::vector< double > &wf)
 
size_t WfSizeFFT (size_t n)
 
std::vector< double > ScintArrivalTimesShape (size_t n, detinfo::LArProperties const &lar_prop)
 
void SubtractBaseline (std::vector< double > &wf, double baseline)
 
void EstimateBaselineStdDev (std::vector< double > &wf, double &_mean, double &_stddev)
 
std::vector< TComplex > DeconvolutionKernel (size_t size, double baseline_stddev, double snr_scaling)
 
- Private Member Functions inherited from opdet::OpDeconvolutionAlg
virtual ~OpDeconvolutionAlg () noexcept=default
 

Private Attributes

bool fDebug
 
int fMaxFFTSizePow
 
std::vector< double > fSinglePEWave
 
bool fApplyExpoAvSmooth
 
bool fApplyUnAvSmooth
 
float fExpoAvSmoothPar
 
short unsigned int fUnAvNeighbours
 
double fHypoSignalTimeWindow
 
bool fHypoSignalCustom
 
std::vector< double > fHypoSignalTimeConsts
 
std::vector< double > fHypoSignalWeights
 
double fHypoSignalScale
 
double fPMTChargeToADC
 
double fDecoWaveformPrecision
 
short unsigned int fBaselineSample
 
std::string fOpDetDataFile
 
std::string fFilter
 
std::string fElectronics
 
bool fScaleHypoSignal
 
bool fUseParamFilter
 
std::vector< double > fFilterParams
 
double fNormUnAvSmooth
 
double fSamplingFreq
 
double fDaphne_Freq
 
size_t MaxBinsFFT
 
unsigned int NDecoWf
 
TF1 * fFilterTF1
 
std::vector< double > fSignalHypothesis
 
std::vector< double > fNoiseHypothesis
 
art::ServiceHandle
< art::TFileService > 
tfs
 
art::ServiceHandle< util::LArFFTfft_service
 

Detailed Description

Definition at line 36 of file OpDeconvolutionAlgWiener_tool.cc.

Constructor & Destructor Documentation

opdet::OpDeconvolutionAlgWiener::OpDeconvolutionAlgWiener ( fhicl::ParameterSet const &  p)
explicit

Definition at line 96 of file OpDeconvolutionAlgWiener_tool.cc.

97 {
98  //read fhicl paramters
99  fDebug = p.get< bool >("Debug");
100  fMaxFFTSizePow = p.get< int >("MaxFFTSizePow", 15);
101  fApplyExpoAvSmooth = p.get< bool >("ApplyExpoAvSmooth");
102  fApplyUnAvSmooth = p.get< bool >("ApplyUnAvSmooth");
103  fExpoAvSmoothPar = p.get< float >("ExpoAvSmoothPar");
104  fUnAvNeighbours = p.get< short unsigned int >("UnAvNeighbours");
105  fHypoSignalTimeWindow = p.get< double >("HypoSignalTimeWindow");
106  fHypoSignalCustom = p.get< bool >("HypoSignalCustom", false);
107  if(fHypoSignalCustom){
108  fHypoSignalTimeConsts = p.get< std::vector<double> >("HypoSignalTimeConsts");
109  fHypoSignalWeights = p.get< std::vector<double> >("HypoSignalWeights");
110  }
111  fHypoSignalScale = p.get< double >("HypoSignalScale");
112  fPMTChargeToADC = p.get< double >("PMTChargeToADC");
113  fDecoWaveformPrecision = p.get< double >("DecoWaveformPrecision");
114  fBaselineSample = p.get< short unsigned int >("BaselineSample");
115  fOpDetDataFile = p.get< std::string >("OpDetDataFile");
116  fFilter = p.get< std::string >("Filter");
117  fElectronics = p.get< std::string >("Electronics");
118  fDaphne_Freq = p.get< double >("DaphneFreq");
119  fScaleHypoSignal = p.get< bool >("ScaleHypoSignal");
120  fUseParamFilter = p.get< bool >("UseParamFilter");
121  fFilterParams = p.get< std::vector<double> >("FilterParams");
122 
123  fNormUnAvSmooth=1./(2*fUnAvNeighbours+1);
124  NDecoWf=0;
125  MaxBinsFFT=std::pow(2, fMaxFFTSizePow);
126 
127  auto const clockData = art::ServiceHandle<detinfo::DetectorClocksService const>()->DataForJob();
128  fSamplingFreq=clockData.OpticalClock().Frequency()/1000.;//in GHz
129  if (fElectronics=="Daphne") fSamplingFreq=fDaphne_Freq/1000.;//in GHz
130  auto const* lar_prop = lar::providerFrom<detinfo::LArPropertiesService>();
131 
132 
133  //Load SER
134  std::string fname;
135  cet::search_path sp("FW_SEARCH_PATH");
136  sp.find_file(fOpDetDataFile, fname);
137  TFile* file = TFile::Open(fname.c_str(), "READ");
138  std::vector<double>* SinglePEVec_p;
139  file->GetObject("SinglePEVec", SinglePEVec_p);
140  if (fElectronics=="Daphne") file->GetObject("SinglePEVec_40ftCable_Daphne", SinglePEVec_p);
141  fSinglePEWave = *SinglePEVec_p;
142  //for(size_t k=0; k<fSinglePEWave.size(); k++) std::cout<<" "<<k<<":"<<fSinglePEWave[k];std::cout<<std::endl;
143  mf::LogInfo("OpDeconvolutionAlg")<<"Loaded SER from "<<fOpDetDataFile<<"... size="<<fSinglePEWave.size()<<std::endl;
144  fSinglePEWave.resize(MaxBinsFFT, 0);
145  file->Close();
146 
147  if(fUseParamFilter){
148  fFilterTF1 = new TF1("FilterTemplate", fFilter.c_str());
149  for(size_t k=0; k<fFilterParams.size(); k++)
150  fFilterTF1->SetParameter(k, fFilterParams[k]);
151  mf::LogInfo("OpDeconvolutionAlg")<<"Creating parametrized filter... TF1:"<<fFilter<<std::endl;
152  }
153  else{
154  //Create light signal hypothesis for "on-fly" Wiener filter
155  fSignalHypothesis.resize(MaxBinsFFT, 0);
156  if(fFilter=="Wiener")
158  else if(fFilter=="Wiener1PE")
159  fSignalHypothesis[0]=1;
160  mf::LogInfo("OpDeconvolutionAlg")<<"Built light signal hypothesis... L="<<fFilter<<" size"<<fSignalHypothesis.size()<<std::endl;
161  }
162 }
string fname
Definition: demo.py:5
std::vector< double > ScintArrivalTimesShape(size_t n, detinfo::LArProperties const &lar_prop)
* file
Definition: file_to_url.sh:69
pdgs p
Definition: selectors.fcl:22
pdgs k
Definition: selectors.fcl:22
opdet::OpDeconvolutionAlgWiener::~OpDeconvolutionAlgWiener ( )
inline

Definition at line 40 of file OpDeconvolutionAlgWiener_tool.cc.

40 {}

Member Function Documentation

void opdet::OpDeconvolutionAlgWiener::ApplyExpoAvSmoothing ( std::vector< double > &  wf)
private

Definition at line 247 of file OpDeconvolutionAlgWiener_tool.cc.

247  {
248  std::transform (std::next(wf.begin(), 1), wf.end(), wf.begin(), std::next(wf.begin(), 1),
249  [&](double _x, double _y) { return fExpoAvSmoothPar*_x+ (1. - fExpoAvSmoothPar)*_y; } );
250 }
static constexpr Sample_t transform(Sample_t sample)
void opdet::OpDeconvolutionAlgWiener::ApplyUnAvSmoothing ( std::vector< double > &  wf)
private

Definition at line 253 of file OpDeconvolutionAlgWiener_tool.cc.

253  {
254  std::vector<double> wf_aux(wf.begin(), wf.end());
255  for(size_t bin=fUnAvNeighbours; bin<wf.size()-fUnAvNeighbours; bin++){
256  double sum=0.;
257  for(size_t nbin=bin-fUnAvNeighbours; nbin<=bin+fUnAvNeighbours; nbin++)
258  sum+=wf_aux[nbin];
259  //std::cout<<bin<<" "<<sum<<" "<<sum*fNormUnAvSmooth<<std::endl;
260  wf[bin]=sum*fNormUnAvSmooth;
261  }
262 }
constexpr details::BinObj< T > bin(T value)
Returns a wrapper to print the specified data in binary format.
std::vector< TComplex > opdet::OpDeconvolutionAlgWiener::DeconvolutionKernel ( size_t  size,
double  baseline_stddev,
double  snr_scaling 
)
private

Definition at line 363 of file OpDeconvolutionAlgWiener_tool.cc.

363  {
364  //Initizalize kernel
365  size_t size=WfSizeFFT(wfsize);
366  TComplex kerinit(0,0,false);
367  std::vector<TComplex> kernel(size, kerinit);
368  fft_service->ReinitializeFFT (size, "", 20);
369 
370  //Prepare detector response FFT
371  std::vector<double> ser( fSinglePEWave.begin(), std::next(fSinglePEWave.begin(), size) );
372  std::vector<TComplex> serfft;
373  serfft.resize(size);
374  fft_service->DoFFT(ser, serfft);
375 
376  if(fUseParamFilter){
377  double freq_step=fSamplingFreq/size;
378  for(size_t k=0; k<size/2; k++){
379  kernel[k]= fFilterTF1->Eval(k*freq_step) / serfft[k] ;
380  //std::cout<<k<<" "<<freq<<" "<<filval<<" "<<fSamplingFreq<<"\n";
381  }
382  }
383  else{
384  //Build Wiener filter kernel: G = Conj(R) / ( |R|^2 + |N|^2/|L|^2)
385  //R=Detector resopnse FFT
386  //N=Noise mean spectral power
387  //L=True signal mean spectral power
388 
389  //Prepare L
390  std::vector<double> hypo( fSignalHypothesis.begin(), std::next(fSignalHypothesis.begin(), size) );
391  std::vector<TComplex> hypofft;
392  hypofft.resize(size);
393  fft_service->DoFFT(hypo, hypofft);
394 
395  //Prepare Noise Spectral Power
396  double noise_power=wfsize*baseline_stddev*baseline_stddev;
397  if(fScaleHypoSignal){
398  noise_power/=pow(snr_scaling, 2);
399  }
400 
401  for(size_t k=0; k<size/2; k++){
402  double den = pow(TComplex::Abs(serfft[k]), 2) + noise_power / pow(TComplex::Abs(hypofft[k]), 2) ;
403  kernel[k]= TComplex::Conjugate( serfft[k] ) / den;
404  //std::cout<<k<<":"<<serfft[k]<<":"<<hypofft[k]<<":"<<den << " ";
405  }
406  }
407 
408 
409  if(fDebug){
410  std::string name="h_wienerfilter_"+std::to_string(NDecoWf);
411  TH1F * hs_wiener = tfs->make< TH1F >
412  (name.c_str(),"Wiener Filter;Frequency Bin;Magnitude",size/2, 0, size/2);
413  for(size_t k=0; k<size/2; k++)
414  hs_wiener->SetBinContent(k, TComplex::Abs( kernel[k]*serfft[k] ) );
415  }
416 
417  return kernel;
418 }
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
art::ServiceHandle< util::LArFFT > fft_service
art::ServiceHandle< art::TFileService > tfs
std::string to_string(WindowPattern const &pattern)
then echo fcl name
pdgs k
Definition: selectors.fcl:22
void opdet::OpDeconvolutionAlgWiener::EstimateBaselineStdDev ( std::vector< double > &  wf,
double &  _mean,
double &  _stddev 
)
private

Definition at line 318 of file OpDeconvolutionAlgWiener_tool.cc.

318  {
319  double minADC=*min_element(wf.begin(), wf.end());
320  double maxADC=*max_element(wf.begin(), wf.end());
321  unsigned nbins=25*ceil(maxADC-minADC);
322  TH1F h_std = TH1F("",";;", nbins, 0, (maxADC-minADC)/2);
323  TH1F h_mean = TH1F("",";;", nbins, minADC, maxADC);
324 
325  for(size_t ix=0; ix<wf.size()-fBaselineSample; ix++){
326  double sum2=0, sum=0;
327  for(size_t jx=ix; jx<ix+fBaselineSample; jx++){
328  sum = sum + wf.at(jx);
329  }
330  sum/=fBaselineSample;
331  for(size_t jx=ix; jx<ix+fBaselineSample; jx++){
332  sum2 = sum2 + pow( wf.at(jx)-sum, 2 );
333  }
334  h_std.Fill( std::sqrt(sum2/fBaselineSample) );
335  h_mean.Fill( sum );
336  //std::cout<<ix<<":"<<wf[ix]<<":"<<sum<<" ";
337  }
338 
339  _stddev=h_std.GetXaxis()->GetBinCenter(h_std.GetMaximumBin());
340  _mean=h_mean.GetXaxis()->GetBinCenter(h_mean.GetMaximumBin());
341 
342  if(fDebug){
343  std::cout<<" -- Estimating baseline...StdDev: "<<_stddev<<" Bias="<<_mean<<std::endl;
344  std::cout<<" .. "<<minADC<<" "<<maxADC<<" "<<maxADC-minADC<<" "<<ceil(log10(maxADC-minADC))<<" "<<nbins<<"\n";
345 
346  std::string name="h_baselinestddev_"+std::to_string(NDecoWf)+std::to_string(_mean);
347  TH1F * hs_std = tfs->make< TH1F > (name.c_str(),"Baseline StdDev;ADC;# entries",
348  h_std.GetNbinsX(), h_std.GetXaxis()->GetXmin(), h_std.GetXaxis()->GetXmax());
349  for(int k=1; k<=h_std.GetNbinsX(); k++)
350  hs_std->SetBinContent(k, h_std.GetBinContent(k));
351 
352  name="h_baselinemean_"+std::to_string(NDecoWf)+std::to_string(_mean);
353  TH1F * hs_mean = tfs->make< TH1F >(name.c_str(),"Baseline Mean;ADC;# entries",
354  h_mean.GetNbinsX(), h_mean.GetXaxis()->GetXmin(), h_mean.GetXaxis()->GetXmax());
355  for(int k=1; k<=h_mean.GetNbinsX(); k++)
356  hs_mean->SetBinContent(k, h_mean.GetBinContent(k));
357  }
358 
359  return;
360 }
art::ServiceHandle< art::TFileService > tfs
std::string to_string(WindowPattern const &pattern)
then echo fcl name
pdgs k
Definition: selectors.fcl:22
BEGIN_PROLOG could also be cout
std::vector< raw::OpDetWaveform > opdet::OpDeconvolutionAlgWiener::RunDeconvolution ( std::vector< raw::OpDetWaveform > const &  wfHandle)
overridevirtual

Implements opdet::OpDeconvolutionAlg.

Definition at line 165 of file OpDeconvolutionAlgWiener_tool.cc.

166 {
167  std::vector<raw::OpDetWaveform> wfDeco;
168  wfDeco.reserve(wfVector.size());
169 
170  for(auto const& wf : wfVector)
171  {
172  //Read waveform
173  size_t wfsize=wf.Waveform().size();
174  if(wfsize>MaxBinsFFT){
175  mf::LogWarning("OpDeconvolutionAlg")<<"Skipping waveform...waveform size is"<<wfsize<<"...maximum allowed FFT size is="<<MaxBinsFFT<<std::endl;
176  continue;
177  }
178  mf::LogInfo("OpDeconvolutionAlg")<<"Deconvolving waveform:"<<NDecoWf<<"...size="<<wfsize<<std::endl;
179  size_t wfsizefft=WfSizeFFT(wfsize);
180 
181  std::vector<double> wave;
182  wave.reserve(wfsizefft);
183  wave.assign(wf.Waveform().begin(), wf.Waveform().end());
184 
185  //Reserve minimum ADCC value for Wiener filter
186  double minADC=*min_element(wave.begin(), wave.end());
187 
188  //Apply waveform smoothing
190  ApplyExpoAvSmoothing(wave);
191  if(fApplyUnAvSmooth)
192  ApplyUnAvSmoothing(wave);
193 
194  //Estimate baseline standrd deviation
195  double baseline_mean=0., baseline_stddev=1.;
196  EstimateBaselineStdDev(wave, baseline_mean, baseline_stddev);
197  double wfPeakPE=fHypoSignalScale*(baseline_mean-minADC)/fPMTChargeToADC;
198  SubtractBaseline(wave, baseline_mean);
199 
200  //Create deconvolution kernel
201  wave.resize(wfsizefft, 0);
202  std::vector<TComplex> fDeconvolutionKernel=DeconvolutionKernel(wfsize, baseline_stddev, wfPeakPE);
203 
204  //Deconvolve raw signal (covolve with kernel)
205  fft_service->ReinitializeFFT(wfsizefft, "", 20);
206  fft_service->Convolute(wave, fDeconvolutionKernel);
207  wave.resize(wfsize);
208 
209  //Set deconvlved waveform precision and restore baseline before saving
210  EstimateBaselineStdDev(wave, baseline_mean, baseline_stddev);
211  SubtractBaseline(wave, baseline_mean);
212  double fDecoWfScaleFactor=1./fDecoWaveformPrecision;
213  std::transform(wave.begin(), wave.end(), wave.begin(), [fDecoWfScaleFactor](double &dec){ return fDecoWfScaleFactor*dec; } );
214 
215  //Debbuging and save wf in hist file
216  if(fDebug){
217  std::cout<<".....Debbuging.....\n";
218  auto minADC_ix=min_element(wave.begin(), wave.end());
219  std::cout<<"Stamp="<<wf.TimeStamp()<<" OpCh"<<wf.ChannelNumber()<<" MinADC="<<minADC<<" (";
220  std::cout<<minADC_ix-wave.begin()<<") Size="<<wf.Waveform().size()<<" ScFactor="<<wfPeakPE<<"\n\n";
221 
222  std::string name="h_deco"+std::to_string(NDecoWf)+"_"+std::to_string(wf.ChannelNumber())+"_"+std::to_string(wf.TimeStamp());
223  TH1F * h_deco = tfs->make< TH1F >(name.c_str(),";Bin;#PE", MaxBinsFFT, 0, MaxBinsFFT);
224  for(size_t k=0; k<wave.size(); k++){
225  //if(fDebug) std::cout<<k<<":"<<wave[k]<<":"<<rawsignal[k]<<" ";
226  h_deco->Fill(k, wave[k]);
227  }
228 
229  name="h_raw"+std::to_string(NDecoWf)+"_"+std::to_string(wf.ChannelNumber())+"_"+std::to_string(wf.TimeStamp());
230  TH1F * h_raw = tfs->make< TH1F >(name.c_str(),";Bin;ADC", MaxBinsFFT, 0, MaxBinsFFT);
231  for(size_t k=0; k<wf.Waveform().size(); k++){
232  //if(fDebug) std::cout<<k<<":"<<wave[k]<<":"<<rawsignal[k]<<" ";
233  h_raw->Fill(k, wf.Waveform()[k]);
234  }
235  }
236 
237  //raw::OpDetWaveform decowf(wf.TimeStamp(), wf.ChannelNumber(), std::vector<short unsigned int> (wave.begin(), std::next(wave.begin(), wf.Waveform().size()) ) );
238  raw::OpDetWaveform decowf( wf.TimeStamp(), wf.ChannelNumber(), std::vector<short unsigned int> (wave.begin(), wave.end()) );
239  wfDeco.push_back(decowf);
240  NDecoWf++;
241  }
242 
243  return wfDeco;
244 }
static constexpr Sample_t transform(Sample_t sample)
void SubtractBaseline(std::vector< double > &wf, double baseline)
std::vector< TComplex > DeconvolutionKernel(size_t size, double baseline_stddev, double snr_scaling)
void EstimateBaselineStdDev(std::vector< double > &wf, double &_mean, double &_stddev)
void ApplyExpoAvSmoothing(std::vector< double > &wf)
art::ServiceHandle< util::LArFFT > fft_service
art::ServiceHandle< art::TFileService > tfs
std::string to_string(WindowPattern const &pattern)
then echo fcl name
pdgs k
Definition: selectors.fcl:22
BEGIN_PROLOG could also be cout
void ApplyUnAvSmoothing(std::vector< double > &wf)
std::vector< double > opdet::OpDeconvolutionAlgWiener::ScintArrivalTimesShape ( size_t  n,
detinfo::LArProperties const &  lar_prop 
)
private

Definition at line 277 of file OpDeconvolutionAlgWiener_tool.cc.

277  {
279  mf::LogWarning("OpDeconvolutionAlg")<<"HypoSignalWeights and HypoSignalTimeConsts must have the same size"<<
280  "...Switching to 2-exponential (default values in LArProperties)"<<std::endl;
281  }
283  fHypoSignalTimeConsts={lar_prop.ScintFastTimeConst(), lar_prop.ScintSlowTimeConst()};
284  fHypoSignalWeights={lar_prop.ScintYieldRatio(), 1.-lar_prop.ScintYieldRatio()};
285  }
286 
287  mf::LogInfo("OpDeconvolutionAlg")<<"Creating hypothesis signal for Wiener filter..."<<std::endl;
288  mf::LogInfo("OpDeconvolutionAlg")<<"Using a "<<fHypoSignalTimeConsts.size()<<"-exponential"<<std::endl;
289  mf::LogInfo("OpDeconvolutionAlg")<<"Time constants [ns]:";
290  for(auto &tau:fHypoSignalTimeConsts) mf::LogInfo("OpDeconvolutionAlg")<<tau<<" ";
291  mf::LogInfo("OpDeconvolutionAlg")<<"\nWeights:";
292  for(auto &w:fHypoSignalWeights) mf::LogInfo("OpDeconvolutionAlg")<<w<<" ";
293  mf::LogInfo("OpDeconvolutionAlg")<<"\n";
294 
295  std::vector<double> v(n, 0);
296  size_t maxbin=fSamplingFreq*fHypoSignalTimeWindow;
297  if(n<maxbin)
298  maxbin=n;
299 
300  for(size_t k=0; k<maxbin; k++){
301  double t = (double)(k) / fSamplingFreq; //in ns
302  double value=0;
303  for(size_t ix=0; ix<fHypoSignalTimeConsts.size(); ix++){
304  value+=fHypoSignalWeights[ix]*std::exp(-1.*t/fHypoSignalTimeConsts[ix]);
305  }
306  v[k]=value;
307  }
308  return v;
309 }
temporary value
pdgs k
Definition: selectors.fcl:22
void opdet::OpDeconvolutionAlgWiener::SubtractBaseline ( std::vector< double > &  wf,
double  baseline 
)
private

Definition at line 312 of file OpDeconvolutionAlgWiener_tool.cc.

312  {
313  std::transform (wf.begin(), wf.end(), wf.begin(), [&](double _x) { return _x-baseline; } );
314  return;
315 }
static constexpr Sample_t transform(Sample_t sample)
BEGIN_PROLOG baseline
size_t opdet::OpDeconvolutionAlgWiener::WfSizeFFT ( size_t  n)
private

Definition at line 265 of file OpDeconvolutionAlgWiener_tool.cc.

265  {
266  if (n && !(n & (n - 1)))
267  return n;
268  size_t cont=0;
269  while(n>0){
270  cont++;
271  n=(n>>1);
272  }
273  return 1<<cont;
274 }

Member Data Documentation

bool opdet::OpDeconvolutionAlgWiener::fApplyExpoAvSmooth
private

Definition at line 49 of file OpDeconvolutionAlgWiener_tool.cc.

bool opdet::OpDeconvolutionAlgWiener::fApplyUnAvSmooth
private

Definition at line 50 of file OpDeconvolutionAlgWiener_tool.cc.

short unsigned int opdet::OpDeconvolutionAlgWiener::fBaselineSample
private

Definition at line 60 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fDaphne_Freq
private

Definition at line 70 of file OpDeconvolutionAlgWiener_tool.cc.

bool opdet::OpDeconvolutionAlgWiener::fDebug
private

Definition at line 46 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fDecoWaveformPrecision
private

Definition at line 59 of file OpDeconvolutionAlgWiener_tool.cc.

std::string opdet::OpDeconvolutionAlgWiener::fElectronics
private

Definition at line 63 of file OpDeconvolutionAlgWiener_tool.cc.

float opdet::OpDeconvolutionAlgWiener::fExpoAvSmoothPar
private

Definition at line 51 of file OpDeconvolutionAlgWiener_tool.cc.

std::string opdet::OpDeconvolutionAlgWiener::fFilter
private

Definition at line 62 of file OpDeconvolutionAlgWiener_tool.cc.

std::vector<double> opdet::OpDeconvolutionAlgWiener::fFilterParams
private

Definition at line 66 of file OpDeconvolutionAlgWiener_tool.cc.

TF1* opdet::OpDeconvolutionAlgWiener::fFilterTF1
private

Definition at line 74 of file OpDeconvolutionAlgWiener_tool.cc.

art::ServiceHandle<util::LArFFT> opdet::OpDeconvolutionAlgWiener::fft_service
private

Definition at line 92 of file OpDeconvolutionAlgWiener_tool.cc.

bool opdet::OpDeconvolutionAlgWiener::fHypoSignalCustom
private

Definition at line 54 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fHypoSignalScale
private

Definition at line 57 of file OpDeconvolutionAlgWiener_tool.cc.

std::vector<double> opdet::OpDeconvolutionAlgWiener::fHypoSignalTimeConsts
private

Definition at line 55 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fHypoSignalTimeWindow
private

Definition at line 53 of file OpDeconvolutionAlgWiener_tool.cc.

std::vector<double> opdet::OpDeconvolutionAlgWiener::fHypoSignalWeights
private

Definition at line 56 of file OpDeconvolutionAlgWiener_tool.cc.

int opdet::OpDeconvolutionAlgWiener::fMaxFFTSizePow
private

Definition at line 47 of file OpDeconvolutionAlgWiener_tool.cc.

std::vector<double> opdet::OpDeconvolutionAlgWiener::fNoiseHypothesis
private

Definition at line 76 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fNormUnAvSmooth
private

Definition at line 68 of file OpDeconvolutionAlgWiener_tool.cc.

std::string opdet::OpDeconvolutionAlgWiener::fOpDetDataFile
private

Definition at line 61 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fPMTChargeToADC
private

Definition at line 58 of file OpDeconvolutionAlgWiener_tool.cc.

double opdet::OpDeconvolutionAlgWiener::fSamplingFreq
private

Definition at line 69 of file OpDeconvolutionAlgWiener_tool.cc.

bool opdet::OpDeconvolutionAlgWiener::fScaleHypoSignal
private

Definition at line 64 of file OpDeconvolutionAlgWiener_tool.cc.

std::vector<double> opdet::OpDeconvolutionAlgWiener::fSignalHypothesis
private

Definition at line 75 of file OpDeconvolutionAlgWiener_tool.cc.

std::vector<double> opdet::OpDeconvolutionAlgWiener::fSinglePEWave
private

Definition at line 48 of file OpDeconvolutionAlgWiener_tool.cc.

short unsigned int opdet::OpDeconvolutionAlgWiener::fUnAvNeighbours
private

Definition at line 52 of file OpDeconvolutionAlgWiener_tool.cc.

bool opdet::OpDeconvolutionAlgWiener::fUseParamFilter
private

Definition at line 65 of file OpDeconvolutionAlgWiener_tool.cc.

size_t opdet::OpDeconvolutionAlgWiener::MaxBinsFFT
private

Definition at line 71 of file OpDeconvolutionAlgWiener_tool.cc.

unsigned int opdet::OpDeconvolutionAlgWiener::NDecoWf
private

Definition at line 72 of file OpDeconvolutionAlgWiener_tool.cc.

art::ServiceHandle<art::TFileService> opdet::OpDeconvolutionAlgWiener::tfs
private

Definition at line 90 of file OpDeconvolutionAlgWiener_tool.cc.


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