LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
AlgoFixedWindow.cxx
Go to the documentation of this file.
1 //
3 // AlgoFixedWindow source
4 //
6 
8 
9 #include "AlgoFixedWindow.h"
10 
11 namespace pmtana {
12 
13  //*******************************************************************************
14  AlgoFixedWindow::AlgoFixedWindow(const std::string name) : PMTPulseRecoBase(name)
15  //*******************************************************************************
16  {
17  Reset();
18 
19  _index_start = 0;
20 
21  _index_end = 0;
22  }
23 
24  //****************************************************************************************
26  const fhicl::ParameterSet& pset,
27  std::unique_ptr<pmtana::RiseTimeCalculatorBase> risetimecalculator,
28  //AlgoFixedWindow::AlgoFixedWindow(const ::fcllite::PSet& pset,
29  const std::string name)
30  : PMTPulseRecoBase(name)
31  //****************************************************************************************
32  {
33  _risetime_calc_ptr = std::move(risetimecalculator);
34 
35  Reset();
36 
37  _index_start = pset.get<size_t>("StartIndex");
38 
39  _index_end = pset.get<size_t>("EndIndex");
40  }
41 
42  //***************************************************************
44  //***************************************************************
45  {
46  if (!(_pulse_v.size())) _pulse_v.push_back(_pulse);
47 
48  _pulse_v[0].reset_param();
49  }
50 
51  //***************************************************************
53  const PedestalMean_t& mean_v,
54  const PedestalSigma_t& sigma_v)
55  //***************************************************************
56  {
57  this->Reset();
58 
59  if (_index_start >= wf.size()) return true;
60 
61  _pulse_v[0].t_start = (double)(_index_start);
62 
63  _pulse_v[0].ped_mean = mean_v.front();
64 
65  _pulse_v[0].ped_sigma = sigma_v.front();
66 
67  if (!_index_end)
68 
69  _pulse_v[0].t_end = (double)(wf.size() - 1);
70 
71  else if (_index_end < wf.size())
72 
73  _pulse_v[0].t_end = (double)_index_end;
74 
75  else
76 
77  _pulse_v[0].t_end = wf.size() - 1;
78 
79  _pulse_v[0].t_max =
80  PMTPulseRecoBase::Max(wf, _pulse_v[0].peak, _index_start, _pulse_v[0].t_end);
81 
82  _pulse_v[0].peak -= mean_v.front();
83 
85 
86  _pulse_v[0].area =
87  _pulse_v[0].area - (_pulse_v[0].t_end - _pulse_v[0].t_start + 1) * mean_v.front();
88 
90  _pulse_v[0].t_rise = _risetime_calc_ptr->RiseTime(
91  {wf.begin() + _pulse.t_start, wf.begin() + _pulse.t_end},
92  {mean_v.begin() + _pulse.t_start, mean_v.begin() + _pulse.t_end},
93  true);
94 
95  return true;
96  }
97 
98 }
bool Integral(const std::vector< short > &wf, double &result, size_t begin=0, size_t end=0) const
std::vector< double > PedestalSigma_t
size_t _index_end
index marker for the end of pulse time window
size_t _index_start
index marker for the beginning of the pulse time window
bool RecoPulse(const pmtana::Waveform_t &, const pmtana::PedestalMean_t &, const pmtana::PedestalSigma_t &)
Implementation of AlgoFixedWindow::reco() method.
pulse_param _pulse
A subject pulse_param object to be filled with the last reconstructed pulse parameters.
size_t Max(const std::vector< short > &wf, double &result, size_t begin=0, size_t end=0) const
T get(std::string const &key) const
Definition: ParameterSet.h:314
std::vector< short > Waveform_t
std::unique_ptr< pmtana::RiseTimeCalculatorBase > _risetime_calc_ptr
Tool for rise time calculation.
Class definition file of AlgoFixedWindow.
AlgoFixedWindow(const std::string name="FixedWindow")
Default ctor.
std::vector< double > PedestalMean_t
pulse_param_array _pulse_v
A container array of pulse_param struct objects to store (possibly multiple) reconstructed pulse(s)...
void Reset()
Implementation of AlgoFixedWindow::reset() method.