All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Private Member Functions | Private Attributes | List of all members
opdet::OpHitFinder Class Reference
Inheritance diagram for opdet::OpHitFinder:

Public Member Functions

 OpHitFinder (const fhicl::ParameterSet &)
 
virtual ~OpHitFinder ()
 
void produce (art::Event &)
 

Private Member Functions

std::map< int, int > GetChannelMap ()
 
std::vector< double > GetSPEScales ()
 
std::vector< double > GetSPEShifts ()
 

Private Attributes

std::string fInputModule
 
std::string fGenModule
 
std::vector< std::string > fInputLabels
 
std::set< unsigned int > fChannelMasks
 
pmtana::PulseRecoManager fPulseRecoMgr
 
pmtana::PMTPulseRecoBasefThreshAlg
 
pmtana::PMTPedestalBasefPedAlg
 
Float_t fHitThreshold
 
unsigned int fMaxOpChannel
 
bool fUseStartTime
 
calib::IPhotonCalibrator const * fCalib = nullptr
 

Detailed Description

Definition at line 63 of file OpHitFinder_module.cc.

Constructor & Destructor Documentation

opdet::OpHitFinder::OpHitFinder ( const fhicl::ParameterSet &  pset)
explicit

Definition at line 104 of file OpHitFinder_module.cc.

104  : 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>());
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 
179  }
const geo::GeometryCore * geometry
unsigned int fMaxOpChannel
void AddRecoAlgo(pmtana::PMTPulseRecoBase *algo, PMTPedestalBase *ped_algo=nullptr)
A method to set pulse reconstruction algorithm.
std::set< unsigned int > fChannelMasks
BEGIN_PROLOG xarapuca_vuv SBNDDecoOpHitFinderXArapuca SPEArea
unsigned int MaxOpChannel() const
Largest optical channel number.
pmtana::PMTPulseRecoBase * fThreshAlg
calib::IPhotonCalibrator const * fCalib
pmtana::PulseRecoManager fPulseRecoMgr
void SetDefaultPedAlgo(pmtana::PMTPedestalBase *algo)
A method to set a choice of pedestal estimation method.
std::vector< std::string > fInputLabels
pmtana::PMTPedestalBase * fPedAlg
opdet::OpHitFinder::~OpHitFinder ( )
virtual

Definition at line 183 of file OpHitFinder_module.cc.

184  {
185 
186  delete fThreshAlg;
187  delete fPedAlg;
188  }
pmtana::PMTPulseRecoBase * fThreshAlg
pmtana::PMTPedestalBase * fPedAlg

Member Function Documentation

std::map<int, int> opdet::OpHitFinder::GetChannelMap ( )
private
std::vector<double> opdet::OpHitFinder::GetSPEScales ( )
private
std::vector<double> opdet::OpHitFinder::GetSPEShifts ( )
private
void opdet::OpHitFinder::produce ( art::Event &  evt)

Definition at line 192 of file OpHitFinder_module.cc.

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  }
const geo::GeometryCore * geometry
EResult err(const char *call)
std::set< unsigned int > fChannelMasks
pmtana::PMTPulseRecoBase * fThreshAlg
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
pmtana::PulseRecoManager fPulseRecoMgr
TCEvent evt
Definition: DataStructs.cxx:8
std::vector< std::string > fInputLabels

Member Data Documentation

calib::IPhotonCalibrator const* opdet::OpHitFinder::fCalib = nullptr
private

Definition at line 91 of file OpHitFinder_module.cc.

std::set<unsigned int> opdet::OpHitFinder::fChannelMasks
private

Definition at line 81 of file OpHitFinder_module.cc.

std::string opdet::OpHitFinder::fGenModule
private

Definition at line 79 of file OpHitFinder_module.cc.

Float_t opdet::OpHitFinder::fHitThreshold
private

Definition at line 87 of file OpHitFinder_module.cc.

std::vector<std::string> opdet::OpHitFinder::fInputLabels
private

Definition at line 80 of file OpHitFinder_module.cc.

std::string opdet::OpHitFinder::fInputModule
private

Definition at line 78 of file OpHitFinder_module.cc.

unsigned int opdet::OpHitFinder::fMaxOpChannel
private

Definition at line 88 of file OpHitFinder_module.cc.

pmtana::PMTPedestalBase* opdet::OpHitFinder::fPedAlg
private

Definition at line 85 of file OpHitFinder_module.cc.

pmtana::PulseRecoManager opdet::OpHitFinder::fPulseRecoMgr
private

Definition at line 83 of file OpHitFinder_module.cc.

pmtana::PMTPulseRecoBase* opdet::OpHitFinder::fThreshAlg
private

Definition at line 84 of file OpHitFinder_module.cc.

bool opdet::OpHitFinder::fUseStartTime
private

Definition at line 89 of file OpHitFinder_module.cc.


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