LArSoft  v07_13_02
Liquid Argon Software toolkit - http://larsoft.org/
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 
12 #ifndef OpHitFinder_H
13 #define OpHitFinder_H 1
14 
15 // LArSoft includes
17 #include "larcore/CoreUtils/ServiceUtil.h" // lar::providerFrom()
36 
37 // Framework includes
41 #include "fhiclcpp/ParameterSet.h"
45 
46 // ROOT includes
47 
48 // C++ Includes
49 #include <map>
50 #include <vector>
51 #include <string>
52 #include <memory>
53 
54 namespace opdet {
55 
56  class OpHitFinder : public art::EDProducer{
57  public:
58 
59  // Standard constructor and destructor for an ART module.
60  explicit OpHitFinder(const fhicl::ParameterSet&);
61  virtual ~OpHitFinder();
62 
63  void beginJob();
64  void endJob();
65  void reconfigure(fhicl::ParameterSet const& pset);
66 
67  // The producer routine, called once per event.
68  void produce(art::Event&);
69 
70  std::map< int, int > GetChannelMap();
71  std::vector< double > GetSPEScales();
72  std::vector< double > GetSPEShifts();
73 
74  private:
75 
76  // The parameters we'll read from the .fcl file.
77  std::string fInputModule; // Input tag for OpDetWaveform collection
78  std::string fGenModule;
79  std::vector< std::string > fInputLabels;
80  std::set< unsigned int > fChannelMasks;
81 
85 
86  Float_t fHitThreshold;
87  unsigned int fMaxOpChannel;
88 
89  calib::IPhotonCalibrator const* fCalib = nullptr;
90  };
91 
92 }
93 
94 namespace opdet {
96 }
97 
98 #endif
99 
100 namespace opdet {
101 
102  //----------------------------------------------------------------------------
103  // Constructor
105  fPulseRecoMgr()
106  {
107 
108  reconfigure(pset);
109 
110  // Initialize the hit finder algorithm
111  auto const hit_alg_pset = pset.get< fhicl::ParameterSet >("HitAlgoPset");
112  std::string threshAlgName = hit_alg_pset.get< std::string >("Name");
113  if (threshAlgName == "Threshold")
114  fThreshAlg = new pmtana::AlgoThreshold(hit_alg_pset);
115  else if (threshAlgName == "SiPM")
116  fThreshAlg = new pmtana::AlgoSiPM(hit_alg_pset);
117  else if (threshAlgName == "SlidingWindow")
118  fThreshAlg = new pmtana::AlgoSlidingWindow(hit_alg_pset);
119  else if (threshAlgName == "FixedWindow")
120  fThreshAlg = new pmtana::AlgoFixedWindow(hit_alg_pset);
121  else if (threshAlgName == "CFD" )
122  fThreshAlg = new pmtana::AlgoCFD(hit_alg_pset);
124  << "Cannot find implementation for "
125  << threshAlgName << " algorithm.\n";
126 
127  auto const ped_alg_pset = pset.get< fhicl::ParameterSet >("PedAlgoPset");
128  std::string pedAlgName = ped_alg_pset.get< std::string >("Name");
129  if (pedAlgName == "Edges")
130  fPedAlg = new pmtana::PedAlgoEdges(ped_alg_pset);
131  else if (pedAlgName == "RollingMean")
132  fPedAlg = new pmtana::PedAlgoRollingMean(ped_alg_pset);
133  else if (pedAlgName == "UB" )
134  fPedAlg = new pmtana::PedAlgoUB(ped_alg_pset);
136  << "Cannot find implementation for "
137  << pedAlgName << " algorithm.\n";
138 
139  produces< std::vector< recob::OpHit > >();
140 
143 
144  }
145 
146  //----------------------------------------------------------------------------
148  {
149 
150  // Indicate that the Input Module comes from .fcl
151  fInputModule = pset.get< std::string >("InputModule");
152  fGenModule = pset.get< std::string >("GenModule");
153  fInputLabels = pset.get< std::vector< std::string > >("InputLabels");
154 
155  for (auto const& ch : pset.get< std::vector< unsigned int > >
156  ("ChannelMasks", std::vector< unsigned int >()))
157  fChannelMasks.insert(ch);
158 
159  fHitThreshold = pset.get< float >("HitThreshold");
160  bool useCalibrator = pset.get< bool > ("UseCalibrator", false);
161 
162  auto const& geometry(*lar::providerFrom< geo::Geometry >());
163  fMaxOpChannel = geometry.MaxOpChannel();
164 
165  if (useCalibrator) {
166  // If useCalibrator, get it from ART
167  fCalib = lar::providerFrom<calib::IPhotonCalibratorService>();
168  }
169  else {
170  // If not useCalibrator, make an internal one based
171  // on fhicl settings to hit finder.
172  bool areaToPE = pset.get< bool > ("AreaToPE");
173  float SPEArea = pset.get< float >("SPEArea");
174  float SPEShift = pset.get< float >("SPEShift", 0.);
175 
176  // Reproduce behavior from GetSPEScales()
177  if (!areaToPE) SPEArea = 20;
178 
179  // Delete and replace if we are reconfiguring
180  if (fCalib) {
181  delete fCalib;
182  }
183 
184  fCalib = new calib::PhotonCalibratorStandard(SPEArea, SPEShift, areaToPE);
185  }
186  }
187 
188  //----------------------------------------------------------------------------
189  // Destructor
191  {
192 
193  delete fThreshAlg;
194  delete fPedAlg;
195 
196  }
197 
198  //----------------------------------------------------------------------------
200  {
201  }
202 
203  //----------------------------------------------------------------------------
205  {
206  }
207 
208  //----------------------------------------------------------------------------
210  {
211 
212  // These is the storage pointer we will put in the event
213  std::unique_ptr< std::vector< recob::OpHit > >
214  HitPtr(new std::vector< recob::OpHit >);
215 
216  std::vector< const sim::BeamGateInfo* > beamGateArray;
217  try
218  {
219  evt.getView(fGenModule, beamGateArray);
220  }
221  catch (art::Exception const& err)
222  {
223  if ( err.categoryCode() != art::errors::ProductNotFound ) throw;
224  }
225 
226  auto const& geometry(*lar::providerFrom< geo::Geometry >());
227  auto const& detectorClocks(*lar::providerFrom< detinfo::DetectorClocksService >());
228  auto const& calibrator(*fCalib);
229  //
230  // Get the pulses from the event
231  //
232 
233  // Reserve a large enough array
234  int totalsize = 0;
235  for (auto label : fInputLabels)
236  {
238  evt.getByLabel(fInputModule, label, wfHandle);
239  if (!wfHandle.isValid()) continue; // Skip non-existent collections
240  totalsize += wfHandle->size();
241  }
242 
243  // Load pulses into WaveformVector
244  std::vector< raw::OpDetWaveform > WaveformVector;
245  WaveformVector.reserve(totalsize);
246  for (auto label : fInputLabels)
247  {
249  evt.getByLabel(fInputModule, label, wfHandle);
250  if (!wfHandle.isValid()) continue; // Skip non-existent collections
251 
252  //WaveformVector.insert(WaveformVector.end(),
253  // wfHandle->begin(), wfHandle->end());
254  for(auto const& wf : *wfHandle)
255  {
256  if (fChannelMasks.find(wf.ChannelNumber())
257  != fChannelMasks.end()) continue;
258  WaveformVector.push_back(wf);
259  }
260  }
261 
262  RunHitFinder(WaveformVector,
263  *HitPtr,
265  *fThreshAlg,
266  geometry,
268  detectorClocks,
269  calibrator);
270 
271  // Store results into the event
272  evt.put(std::move(HitPtr));
273 
274  }
275 
276 } // namespace opdet
Utilities related to art service access.
std::set< unsigned int > fChannelMasks
void reconfigure(fhicl::ParameterSet const &pset)
unsigned int fMaxOpChannel
std::vector< double > GetSPEScales()
std::vector< std::string > fInputLabels
void AddRecoAlgo(pmtana::PMTPulseRecoBase *algo, PMTPedestalBase *ped_algo=nullptr)
A method to set pulse reconstruction algorithm.
OpHitFinder(const fhicl::ParameterSet &)
ProductID put(std::unique_ptr< PROD > &&product)
Definition: Event.h:102
bool isValid() const
Definition: Handle.h:190
std::size_t getView(std::string const &moduleLabel, std::string const &productInstanceName, std::vector< ELEMENT const * > &result) const
Definition: DataViewImpl.h:474
std::map< int, int > GetChannelMap()
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:42
Class definition file of PedAlgoRollingMean.
Class definition file of AlgoCFD.
T get(std::string const &key) const
Definition: ParameterSet.h:231
Class definition file of PedAlgoUB.
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::DetectorClocks const &detectorClocks, calib::IPhotonCalibrator const &calibrator)
Definition: OpHitAlg.cxx:18
pmtana::PMTPulseRecoBase * fThreshAlg
void produce(art::Event &)
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
Class definition file of AlgoFixedWindow.
calib::IPhotonCalibrator const * fCalib
std::vector< double > GetSPEShifts()
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
Class definition file of AlgoSlidingWindow.
Class definition file of AlgoThreshold.
pmtana::PulseRecoManager fPulseRecoMgr
void SetDefaultPedAlgo(pmtana::PMTPedestalBase *algo)
A method to set a choice of pedestal estimation method.
Class definition file of PMTPulseRecoBase.
Class definition file of PedAlgoEdges.
TCEvent evt
Definition: DataStructs.cxx:5
art framework interface to geometry description
pmtana::PMTPedestalBase * fPedAlg
Class definition file of PulseRecoManager.