All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
OpHitFinder_module.cc
Go to the documentation of this file.
1 // -*- mode: c++; c-basic-offset: 2; -*-
2 // Gleb Sinev, Duke, 2016
3 //
4 // This module finds periods of time-localized activity
5 // on each channel of the optical system and creates OpHits.
6 //
7 // Split from OpFlashFinder_module.cc
8 // by Ben Jones, MIT, 2013
9 //
10 
11 // LArSoft includes
23 #include "larcore/CoreUtils/ServiceUtil.h" // lar::providerFrom()
32 
33 // Framework includes
34 #include "art/Framework/Core/EDProducer.h"
35 #include "art/Framework/Core/ModuleMacros.h"
36 #include "art/Framework/Principal/Event.h"
37 #include "art/Framework/Principal/Handle.h"
38 #include "art/Utilities/make_tool.h"
39 #include "canvas/Utilities/Exception.h"
40 #include "fhiclcpp/ParameterSet.h"
41 
42 // ROOT includes
43 
44 // C++ Includes
45 #include <map>
46 #include <memory>
47 #include <string>
48 
49 namespace {
50  template <typename T>
51  pmtana::PMTPulseRecoBase* thresholdAlgorithm(fhicl::ParameterSet const& hit_alg_pset,
52  std::optional<fhicl::ParameterSet> const& rise_alg_pset)
53  {
54  if (rise_alg_pset)
55  return new T(hit_alg_pset, art::make_tool<pmtana::RiseTimeCalculatorBase>(*rise_alg_pset) );
56  else
57  return new T(hit_alg_pset, nullptr);
58  }
59 }
60 
61 namespace opdet {
62 
63  class OpHitFinder : public art::EDProducer {
64  public:
65  // Standard constructor and destructor for an ART module.
66  explicit OpHitFinder(const fhicl::ParameterSet&);
67  virtual ~OpHitFinder();
68 
69  // The producer routine, called once per event.
70  void produce(art::Event&);
71 
72  private:
73  std::map<int, int> GetChannelMap();
74  std::vector<double> GetSPEScales();
75  std::vector<double> GetSPEShifts();
76 
77  // The parameters we'll read from the .fcl file.
78  std::string fInputModule; // Input tag for OpDetWaveform collection
79  std::string fGenModule;
80  std::vector<std::string> fInputLabels;
81  std::set<unsigned int> fChannelMasks;
82 
86 
87  Float_t fHitThreshold;
88  unsigned int fMaxOpChannel;
90 
91  calib::IPhotonCalibrator const* fCalib = nullptr;
92  };
93 
94 }
95 
96 namespace opdet {
97  DEFINE_ART_MODULE(OpHitFinder)
98 }
99 
100 namespace opdet {
101 
102  //----------------------------------------------------------------------------
103  // Constructor
104  OpHitFinder::OpHitFinder(const fhicl::ParameterSet& pset) : EDProducer{pset}, fPulseRecoMgr()
105  {
106  // Indicate that the Input Module comes from .fcl
107  fInputModule = pset.get<std::string>("InputModule");
108  fGenModule = pset.get<std::string>("GenModule");
109  fInputLabels = pset.get<std::vector<std::string>>("InputLabels");
110  fUseStartTime = pset.get<bool>("UseStartTime", false);
111 
112  for (auto const& ch :
113  pset.get<std::vector<unsigned int>>("ChannelMasks", std::vector<unsigned int>()))
114  fChannelMasks.insert(ch);
115 
116  fHitThreshold = pset.get<float>("HitThreshold");
117  bool useCalibrator = pset.get<bool>("UseCalibrator", false);
118 
119  auto const& geometry(*lar::providerFrom<geo::Geometry>());
120  fMaxOpChannel = geometry.MaxOpChannel();
121 
122  if (useCalibrator) {
123  // If useCalibrator, get it from ART
124  fCalib = lar::providerFrom<calib::IPhotonCalibratorService>();
125  }
126  else {
127  // If not useCalibrator, make an internal one based
128  // on fhicl settings to hit finder.
129  bool areaToPE = pset.get<bool>("AreaToPE");
130  float SPEArea = pset.get<float>("SPEArea");
131  float SPEShift = pset.get<float>("SPEShift", 0.);
132 
133  // Reproduce behavior from GetSPEScales()
134  if (!areaToPE) SPEArea = 20;
135 
136  // Delete and replace if we are reconfiguring
137  if (fCalib) { delete fCalib; }
138 
139  fCalib = new calib::PhotonCalibratorStandard(SPEArea, SPEShift, areaToPE);
140  }
141 
142  // Initialize the rise time calculator tool
143  auto const rise_alg_pset = pset.get_if_present<fhicl::ParameterSet>("RiseTimeCalculator");
144 
145  // Initialize the hit finder algorithm
146  auto const hit_alg_pset = pset.get<fhicl::ParameterSet>("HitAlgoPset");
147  std::string threshAlgName = hit_alg_pset.get<std::string>("Name");
148  if (threshAlgName == "Threshold")
149  fThreshAlg = thresholdAlgorithm<pmtana::AlgoThreshold>(hit_alg_pset, rise_alg_pset);
150  else if (threshAlgName == "SiPM")
151  fThreshAlg = thresholdAlgorithm<pmtana::AlgoSiPM>(hit_alg_pset, rise_alg_pset);
152  else if (threshAlgName == "SlidingWindow")
153  fThreshAlg = thresholdAlgorithm<pmtana::AlgoSlidingWindow>(hit_alg_pset, rise_alg_pset);
154  else if (threshAlgName == "FixedWindow")
155  fThreshAlg = thresholdAlgorithm<pmtana::AlgoFixedWindow>(hit_alg_pset, rise_alg_pset);
156  else if (threshAlgName == "CFD")
157  fThreshAlg = thresholdAlgorithm<pmtana::AlgoCFD>(hit_alg_pset, rise_alg_pset);
158  else
159  throw art::Exception(art::errors::UnimplementedFeature)
160  << "Cannot find implementation for " << threshAlgName << " algorithm.\n";
161 
162  // Initialize the pedestal estimation algorithm
163  auto const ped_alg_pset = pset.get<fhicl::ParameterSet>("PedAlgoPset");
164  std::string pedAlgName = ped_alg_pset.get<std::string>("Name");
165  if (pedAlgName == "Edges")
166  fPedAlg = new pmtana::PedAlgoEdges(ped_alg_pset);
167  else if (pedAlgName == "RollingMean")
168  fPedAlg = new pmtana::PedAlgoRollingMean(ped_alg_pset);
169  else if (pedAlgName == "UB")
170  fPedAlg = new pmtana::PedAlgoUB(ped_alg_pset);
171  else
172  throw art::Exception(art::errors::UnimplementedFeature)
173  << "Cannot find implementation for " << pedAlgName << " algorithm.\n";
174 
175  produces<std::vector<recob::OpHit>>();
176 
177  fPulseRecoMgr.AddRecoAlgo(fThreshAlg);
178  fPulseRecoMgr.SetDefaultPedAlgo(fPedAlg);
179  }
180 
181  //----------------------------------------------------------------------------
182  // Destructor
184  {
185 
186  delete fThreshAlg;
187  delete fPedAlg;
188  }
189 
190  //----------------------------------------------------------------------------
191  void
193  {
194 
195  // These is the storage pointer we will put in the event
196  std::unique_ptr<std::vector<recob::OpHit>> HitPtr(new std::vector<recob::OpHit>);
197 
198  std::vector<const sim::BeamGateInfo*> beamGateArray;
199  try {
200  evt.getView(fGenModule, beamGateArray);
201  }
202  catch (art::Exception const& err) {
203  if (err.categoryCode() != art::errors::ProductNotFound) throw;
204  }
205 
206  auto const& geometry(*lar::providerFrom<geo::Geometry>());
207  auto const clock_data =
208  art::ServiceHandle<detinfo::DetectorClocksService const>()->DataFor(evt);
209  auto const& calibrator(*fCalib);
210  //
211  // Get the pulses from the event
212  //
213 
214  // Load pulses into WaveformVector
215  if (fChannelMasks.empty() && fInputLabels.size() < 2) {
216  art::Handle<std::vector<raw::OpDetWaveform>> wfHandle;
217  if (fInputLabels.empty())
218  evt.getByLabel(fInputModule, wfHandle);
219  else
220  evt.getByLabel(fInputModule, fInputLabels.front(), wfHandle);
221  assert(wfHandle.isValid());
222  RunHitFinder(*wfHandle,
223  *HitPtr,
225  *fThreshAlg,
226  geometry,
228  clock_data,
229  calibrator,
230  fUseStartTime);
231  }
232  else {
233 
234  // Reserve a large enough array
235  int totalsize = 0;
236  for (auto label : fInputLabels) {
237  art::Handle<std::vector<raw::OpDetWaveform>> wfHandle;
238  evt.getByLabel(fInputModule, label, wfHandle);
239  if (!wfHandle.isValid()) continue; // Skip non-existent collections
240  totalsize += wfHandle->size();
241  }
242 
243  std::vector<raw::OpDetWaveform> WaveformVector;
244  WaveformVector.reserve(totalsize);
245 
246  for (auto label : fInputLabels) {
247  art::Handle<std::vector<raw::OpDetWaveform>> wfHandle;
248  evt.getByLabel(fInputModule, label, wfHandle);
249  if (!wfHandle.isValid()) continue; // Skip non-existent collections
250 
251  //WaveformVector.insert(WaveformVector.end(),
252  // wfHandle->begin(), wfHandle->end());
253  for (auto const& wf : *wfHandle) {
254  if (fChannelMasks.find(wf.ChannelNumber()) != fChannelMasks.end()) continue;
255  WaveformVector.push_back(wf);
256  }
257  }
258 
259  RunHitFinder(WaveformVector,
260  *HitPtr,
262  *fThreshAlg,
263  geometry,
265  clock_data,
266  calibrator,
267  fUseStartTime);
268  }
269  // Store results into the event
270  evt.put(std::move(HitPtr));
271  }
272 
273 } // namespace opdet
Utilities related to art service access.
const geo::GeometryCore * geometry
unsigned int fMaxOpChannel
std::map< int, int > GetChannelMap()
EResult err(const char *call)
OpHitFinder(const fhicl::ParameterSet &)
std::set< unsigned int > fChannelMasks
BEGIN_PROLOG xarapuca_vuv SBNDDecoOpHitFinderXArapuca SPEArea
Class definition file of PedAlgoRollingMean.
Class definition file of AlgoCFD.
std::vector< double > GetSPEShifts()
Class definition file of PedAlgoUB.
unsigned int MaxOpChannel() const
Largest optical channel number.
pmtana::PMTPulseRecoBase * fThreshAlg
void produce(art::Event &)
Class definition file of AlgoFixedWindow.
calib::IPhotonCalibrator const * fCalib
void RunHitFinder(std::vector< raw::OpDetWaveform > const &opDetWaveformVector, std::vector< recob::OpHit > &hitVector, pmtana::PulseRecoManager const &pulseRecoMgr, pmtana::PMTPulseRecoBase const &threshAlg, geo::GeometryCore const &geometry, float hitThreshold, detinfo::DetectorClocksData const &clocksData, calib::IPhotonCalibrator const &calibrator, bool use_start_time)
Definition: OpHitAlg.cxx:30
Class definition file of AlgoSlidingWindow.
Class definition file of AlgoThreshold.
pmtana::PulseRecoManager fPulseRecoMgr
std::vector< double > GetSPEScales()
Class definition file of PMTPulseRecoBase.
Class definition file of PedAlgoEdges.
TCEvent evt
Definition: DataStructs.cxx:8
std::vector< std::string > fInputLabels
art framework interface to geometry description
pmtana::PMTPedestalBase * fPedAlg
Class definition file of PulseRecoManager.