LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
RiseTimeThreshold_tool.cc
Go to the documentation of this file.
1 
13 #include "fhiclcpp/types/Atom.h"
14 
15 #include "RiseTimeCalculatorBase.h"
16 
17 namespace pmtana {
18 
20 
21  public:
22  //Configuration parameters
23  struct Config {
24 
26  };
27 
28  // Default constructor
29  explicit RiseTimeThreshold(art::ToolConfigTable<Config> const& config);
30 
31  // Method to calculate the OpFlash t0
32  double RiseTime(const pmtana::Waveform_t& wf_pulse,
33  const pmtana::PedestalMean_t& ped_pulse,
34  bool _positive) const override;
35 
36  private:
37  double fPeakRatio;
38  };
39 
41  : fPeakRatio{config().PeakRatio()}
42  {}
43 
45  const pmtana::PedestalMean_t& ped_pulse,
46  bool _positive) const
47  {
48 
49  // Pedestal-subtracted pulse
50  std::vector<double> wf_aux(ped_pulse);
51  if (_positive) {
52  for (size_t ix = 0; ix < wf_aux.size(); ix++) {
53  wf_aux[ix] = ((double)wf_pulse[ix]) - wf_aux[ix];
54  }
55  }
56  else {
57  for (size_t ix = 0; ix < wf_aux.size(); ix++) {
58  wf_aux[ix] = wf_aux[ix] - ((double)wf_pulse[ix]);
59  }
60  }
61 
62  auto it_max = max_element(wf_aux.begin(), wf_aux.end());
63  size_t rise = std::lower_bound(wf_aux.begin(), it_max, fPeakRatio * (*it_max)) - wf_aux.begin();
64 
65  return rise;
66  }
67 
68 }
69 
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:42
RiseTimeThreshold(art::ToolConfigTable< Config > const &config)
Interfacce class for a tool to calculate the pulse rise time.
double RiseTime(const pmtana::Waveform_t &wf_pulse, const pmtana::PedestalMean_t &ped_pulse, bool _positive) const override
std::vector< short > Waveform_t
std::vector< double > PedestalMean_t