LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
pmtana::RiseTimeGaussFit Class Reference
Inheritance diagram for pmtana::RiseTimeGaussFit:
pmtana::RiseTimeCalculatorBase

Classes

struct  Config
 

Public Member Functions

 RiseTimeGaussFit (art::ToolConfigTable< Config > const &config)
 
double RiseTime (const pmtana::Waveform_t &wf_pulse, const pmtana::PedestalMean_t &ped_pulse, bool _positive) const override
 
std::size_t findFirstMax (const std::vector< double > &arr, double threshold) const
 

Private Attributes

double fMinAmp
 
double fInitSigma
 
double fTolerance
 
int fNbins
 

Detailed Description

Definition at line 24 of file RiseTimeGaussFit_tool.cc.

Constructor & Destructor Documentation

pmtana::RiseTimeGaussFit::RiseTimeGaussFit ( art::ToolConfigTable< Config > const &  config)
explicit

Definition at line 53 of file RiseTimeGaussFit_tool.cc.

References fInitSigma, fNbins, and fTolerance.

54  : fMinAmp{config().MinAmp()}
55  , fInitSigma{config().InitSigma()}
56  , fTolerance{config().Tolerance()}
57  , fNbins{config().Nbins()}
58  {}

Member Function Documentation

std::size_t pmtana::RiseTimeGaussFit::findFirstMax ( const std::vector< double > &  arr,
double  threshold 
) const

Linear search, O(N), 1st peak should be close to the start of the vector for scintillation LAr Signals. Returns the position of the first local max

Definition at line 106 of file RiseTimeGaussFit_tool.cc.

References DEFINE_ART_CLASS_TOOL, and n.

Referenced by RiseTime().

107  {
111  double max = arr[0];
112  for (std::size_t i = 1, n = arr.size(); i < n; ++i) {
113  if (arr[i] >= max)
114  max = arr[i];
115 
116  else if (max < threshold) //didn't pass the threshold, keep searching
117  max = arr[i];
118 
119  else
120  return i - 1;
121  }
122 
123  //No local maxima found.
124  mf::LogInfo("RiseTimeGaussFit") << "No local max found above fixed threshold: " << threshold;
125  return 0;
126  }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
Char_t n[5]
double pmtana::RiseTimeGaussFit::RiseTime ( const pmtana::Waveform_t wf_pulse,
const pmtana::PedestalMean_t ped_pulse,
bool  _positive 
) const
overridevirtual

Implements pmtana::RiseTimeCalculatorBase.

Definition at line 60 of file RiseTimeGaussFit_tool.cc.

References util::abs(), f, findFirstMax(), fInitSigma, fMinAmp, fNbins, and fTolerance.

63  {
64 
65  // Pedestal-subtracted pulse
66  std::vector<double> wf_aux(ped_pulse);
67  if (_positive) {
68  for (size_t ix = 0; ix < wf_aux.size(); ix++) {
69  wf_aux[ix] = ((double)wf_pulse[ix]) - wf_aux[ix];
70  }
71  }
72  else {
73  for (size_t ix = 0; ix < wf_aux.size(); ix++) {
74  wf_aux[ix] = wf_aux[ix] - ((double)wf_pulse[ix]);
75  }
76  }
77 
78  // Find first local maximum
79  size_t first_max = findFirstMax(wf_aux, fMinAmp);
80 
81  // Create & fill th1
82  TH1F* h_aux = new TH1F("aux", "aux", wf_aux.size(), -0.5, wf_aux.size() - 0.5);
83  for (long unsigned int j = 0; j < wf_aux.size(); j++)
84  h_aux->SetBinContent(j + 1, wf_aux[j]);
85 
86  // Function & initial values for the fit, ROOT Gaussian: [p0]*e**(-0.5 (x-[p1])**2 / [p2]**2)
87  TF1* f = new TF1("f", "gaus", double(first_max) - fNbins, double(first_max) + fNbins);
88  f->SetParameters(wf_aux[first_max], first_max, fInitSigma);
89  // Fit
90  h_aux->Fit(f, "q", "SAME", first_max - fNbins, first_max + fNbins);
91  double t_fit = f->GetParameter(1);
92  double peak_time;
93  if (
94  std::abs(t_fit - first_max) <
95  fTolerance) { //check fit is close in distance to the original max, use max peak as time otherwise
96  peak_time = t_fit;
97  }
98  else {
99  peak_time = first_max;
100  mf::LogInfo("RiseTimeGaussFit") << "No good fit found, keeping 1st max bin instead";
101  }
102 
103  return peak_time;
104  }
MaybeLogger_< ELseverityLevel::ELsev_info, false > LogInfo
constexpr auto abs(T v)
Returns the absolute value of the argument.
std::size_t findFirstMax(const std::vector< double > &arr, double threshold) const
TFile f
Definition: plotHisto.C:6

Member Data Documentation

double pmtana::RiseTimeGaussFit::fInitSigma
private

Definition at line 48 of file RiseTimeGaussFit_tool.cc.

Referenced by RiseTime(), and RiseTimeGaussFit().

double pmtana::RiseTimeGaussFit::fMinAmp
private

Definition at line 47 of file RiseTimeGaussFit_tool.cc.

Referenced by RiseTime().

int pmtana::RiseTimeGaussFit::fNbins
private

Definition at line 50 of file RiseTimeGaussFit_tool.cc.

Referenced by RiseTime(), and RiseTimeGaussFit().

double pmtana::RiseTimeGaussFit::fTolerance
private

Definition at line 49 of file RiseTimeGaussFit_tool.cc.

Referenced by RiseTime(), and RiseTimeGaussFit().


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