Public Member Functions | Private Member Functions | Private Attributes | List of all members
opdet::OpHitFinder Class Reference
Inheritance diagram for opdet::OpHitFinder:
art::EDProducer art::detail::Producer art::detail::LegacyModule art::Modifier art::ModuleBase art::ProductRegistryHelper

Public Member Functions

 OpHitFinder (const fhicl::ParameterSet &)
 
virtual ~OpHitFinder ()
 
void produce (art::Event &)
 
- Public Member Functions inherited from art::EDProducer
 EDProducer (fhicl::ParameterSet const &pset)
 
template<typename Config >
 EDProducer (Table< Config > const &config)
 
std::string workerType () const
 
- Public Member Functions inherited from art::detail::Producer
virtual ~Producer () noexcept
 
 Producer (fhicl::ParameterSet const &)
 
 Producer (Producer const &)=delete
 
 Producer (Producer &&)=delete
 
Produceroperator= (Producer const &)=delete
 
Produceroperator= (Producer &&)=delete
 
void doBeginJob (SharedResources const &resources)
 
void doEndJob ()
 
void doRespondToOpenInputFile (FileBlock const &fb)
 
void doRespondToCloseInputFile (FileBlock const &fb)
 
void doRespondToOpenOutputFiles (FileBlock const &fb)
 
void doRespondToCloseOutputFiles (FileBlock const &fb)
 
bool doBeginRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doEndRun (RunPrincipal &rp, ModuleContext const &mc)
 
bool doBeginSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEndSubRun (SubRunPrincipal &srp, ModuleContext const &mc)
 
bool doEvent (EventPrincipal &ep, ModuleContext const &mc, std::atomic< std::size_t > &counts_run, std::atomic< std::size_t > &counts_passed, std::atomic< std::size_t > &counts_failed)
 
- Public Member Functions inherited from art::Modifier
 ~Modifier () noexcept
 
 Modifier ()
 
 Modifier (Modifier const &)=delete
 
 Modifier (Modifier &&)=delete
 
Modifieroperator= (Modifier const &)=delete
 
Modifieroperator= (Modifier &&)=delete
 
- Public Member Functions inherited from art::ModuleBase
virtual ~ModuleBase () noexcept
 
 ModuleBase ()
 
ModuleDescription const & moduleDescription () const
 
void setModuleDescription (ModuleDescription const &)
 
std::array< std::vector< ProductInfo >, NumBranchTypes > const & getConsumables () const
 
void sortConsumables (std::string const &current_process_name)
 
template<typename T , BranchType BT>
ViewToken< T > consumesView (InputTag const &tag)
 
template<typename T , BranchType BT>
ViewToken< T > mayConsumeView (InputTag const &tag)
 

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::stringfInputLabels
 
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
 

Additional Inherited Members

- Public Types inherited from art::EDProducer
using ModuleType = EDProducer
 
using WorkerType = WorkerT< EDProducer >
 
- Public Types inherited from art::detail::Producer
template<typename UserConfig , typename KeysToIgnore = void>
using Table = Modifier::Table< UserConfig, KeysToIgnore >
 
- Public Types inherited from art::Modifier
template<typename UserConfig , typename UserKeysToIgnore = void>
using Table = ProducerTable< UserConfig, detail::ModuleConfig, UserKeysToIgnore >
 
- Static Public Member Functions inherited from art::EDProducer
static void commitEvent (EventPrincipal &ep, Event &e)
 
- Protected Member Functions inherited from art::ModuleBase
ConsumesCollectorconsumesCollector ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > consumes (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > consumesView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void consumesMany ()
 
template<typename T , BranchType = InEvent>
ProductToken< T > mayConsume (InputTag const &)
 
template<typename Element , BranchType = InEvent>
ViewToken< Element > mayConsumeView (InputTag const &)
 
template<typename T , BranchType = InEvent>
void mayConsumeMany ()
 

Detailed Description

Definition at line 50 of file OpHitFinder_module.cc.

Constructor & Destructor Documentation

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

Definition at line 91 of file OpHitFinder_module.cc.

91  : EDProducer{pset}, fPulseRecoMgr()
92  {
93  // Indicate that the Input Module comes from .fcl
94  fInputModule = pset.get<std::string>("InputModule");
95  fGenModule = pset.get<std::string>("GenModule");
96  fInputLabels = pset.get<std::vector<std::string>>("InputLabels");
97  fUseStartTime = pset.get<bool>("UseStartTime", false);
98 
99  for (auto const& ch :
100  pset.get<std::vector<unsigned int>>("ChannelMasks", std::vector<unsigned int>()))
101  fChannelMasks.insert(ch);
102 
103  fHitThreshold = pset.get<float>("HitThreshold");
104  bool useCalibrator = pset.get<bool>("UseCalibrator", false);
105 
106  auto const& geometry(*lar::providerFrom<geo::Geometry>());
107  fMaxOpChannel = geometry.MaxOpChannel();
108 
109  if (useCalibrator) {
110  // If useCalibrator, get it from ART
111  fCalib = lar::providerFrom<calib::IPhotonCalibratorService>();
112  }
113  else {
114  // If not useCalibrator, make an internal one based
115  // on fhicl settings to hit finder.
116  bool areaToPE = pset.get<bool>("AreaToPE");
117  float SPEArea = pset.get<float>("SPEArea");
118  float SPEShift = pset.get<float>("SPEShift", 0.);
119 
120  // Reproduce behavior from GetSPEScales()
121  if (!areaToPE) SPEArea = 20;
122 
123  // Delete and replace if we are reconfiguring
124  if (fCalib) { delete fCalib; }
125 
126  fCalib = new calib::PhotonCalibratorStandard(SPEArea, SPEShift, areaToPE);
127  }
128 
129  // Initialize the hit finder algorithm
130  auto const hit_alg_pset = pset.get<fhicl::ParameterSet>("HitAlgoPset");
131  std::string threshAlgName = hit_alg_pset.get<std::string>("Name");
132  if (threshAlgName == "Threshold")
133  fThreshAlg = new pmtana::AlgoThreshold(hit_alg_pset);
134  else if (threshAlgName == "SiPM")
135  fThreshAlg = new pmtana::AlgoSiPM(hit_alg_pset);
136  else if (threshAlgName == "SlidingWindow")
137  fThreshAlg = new pmtana::AlgoSlidingWindow(hit_alg_pset);
138  else if (threshAlgName == "FixedWindow")
139  fThreshAlg = new pmtana::AlgoFixedWindow(hit_alg_pset);
140  else if (threshAlgName == "CFD")
141  fThreshAlg = new pmtana::AlgoCFD(hit_alg_pset);
142  else
144  << "Cannot find implementation for " << threshAlgName << " algorithm.\n";
145 
146  auto const ped_alg_pset = pset.get<fhicl::ParameterSet>("PedAlgoPset");
147  std::string pedAlgName = ped_alg_pset.get<std::string>("Name");
148  if (pedAlgName == "Edges")
149  fPedAlg = new pmtana::PedAlgoEdges(ped_alg_pset);
150  else if (pedAlgName == "RollingMean")
151  fPedAlg = new pmtana::PedAlgoRollingMean(ped_alg_pset);
152  else if (pedAlgName == "UB")
153  fPedAlg = new pmtana::PedAlgoUB(ped_alg_pset);
154  else
156  << "Cannot find implementation for " << pedAlgName << " algorithm.\n";
157 
158  produces<std::vector<recob::OpHit>>();
159 
162  }
unsigned int fMaxOpChannel
std::string string
Definition: nybbler.cc:12
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
void AddRecoAlgo(pmtana::PMTPulseRecoBase *algo, PMTPedestalBase *ped_algo=nullptr)
A method to set pulse reconstruction algorithm.
std::set< unsigned int > fChannelMasks
T get(std::string const &key) const
Definition: ParameterSet.h:271
pmtana::PMTPulseRecoBase * fThreshAlg
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
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 166 of file OpHitFinder_module.cc.

167  {
168 
169  delete fThreshAlg;
170  delete fPedAlg;
171  }
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)
virtual

Implements art::EDProducer.

Definition at line 175 of file OpHitFinder_module.cc.

176  {
177 
178  // These is the storage pointer we will put in the event
179  std::unique_ptr<std::vector<recob::OpHit>> HitPtr(new std::vector<recob::OpHit>);
180 
181  std::vector<const sim::BeamGateInfo*> beamGateArray;
182  try {
183  evt.getView(fGenModule, beamGateArray);
184  }
185  catch (art::Exception const& err) {
186  if (err.categoryCode() != art::errors::ProductNotFound) throw;
187  }
188 
189  auto const& geometry(*lar::providerFrom<geo::Geometry>());
190  auto const clock_data =
192  auto const& calibrator(*fCalib);
193  //
194  // Get the pulses from the event
195  //
196 
197  // Load pulses into WaveformVector
198  if (fChannelMasks.empty() && fInputLabels.size() < 2) {
200  if (fInputLabels.empty())
201  evt.getByLabel(fInputModule, wfHandle);
202  else
203  evt.getByLabel(fInputModule, fInputLabels.front(), wfHandle);
204  assert(wfHandle.isValid());
205  RunHitFinder(*wfHandle,
206  *HitPtr,
208  *fThreshAlg,
209  geometry,
211  clock_data,
212  calibrator,
213  fUseStartTime);
214  }
215  else {
216 
217  // Reserve a large enough array
218  int totalsize = 0;
219  for (auto label : fInputLabels) {
221  evt.getByLabel(fInputModule, label, wfHandle);
222  if (!wfHandle.isValid()) continue; // Skip non-existent collections
223  totalsize += wfHandle->size();
224  }
225 
226  std::vector<raw::OpDetWaveform> WaveformVector;
227  WaveformVector.reserve(totalsize);
228 
229  for (auto label : fInputLabels) {
231  evt.getByLabel(fInputModule, label, wfHandle);
232  if (!wfHandle.isValid()) continue; // Skip non-existent collections
233 
234  //WaveformVector.insert(WaveformVector.end(),
235  // wfHandle->begin(), wfHandle->end());
236  for (auto const& wf : *wfHandle) {
237  if (fChannelMasks.find(wf.ChannelNumber()) != fChannelMasks.end()) continue;
238  WaveformVector.push_back(wf);
239  }
240  }
241 
242  RunHitFinder(WaveformVector,
243  *HitPtr,
245  *fThreshAlg,
246  geometry,
248  clock_data,
249  calibrator,
250  fUseStartTime);
251  }
252  // Store results into the event
253  evt.put(std::move(HitPtr));
254  }
std::set< unsigned int > fChannelMasks
bool isValid() const noexcept
Definition: Handle.h:191
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
Definition: DataViewImpl.h:633
def move(depos, offset)
Definition: depos.py:107
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
void err(const char *fmt,...)
Definition: message.cpp:226
pmtana::PMTPulseRecoBase * fThreshAlg
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
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:29
pmtana::PulseRecoManager fPulseRecoMgr
std::size_t getView(std::string const &moduleLabel, std::string const &productInstanceName, std::string const &processName, std::vector< ELEMENT const * > &result) const
Definition: DataViewImpl.h:500
std::vector< std::string > fInputLabels

Member Data Documentation

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

Definition at line 78 of file OpHitFinder_module.cc.

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

Definition at line 68 of file OpHitFinder_module.cc.

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

Definition at line 66 of file OpHitFinder_module.cc.

Float_t opdet::OpHitFinder::fHitThreshold
private

Definition at line 74 of file OpHitFinder_module.cc.

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

Definition at line 67 of file OpHitFinder_module.cc.

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

Definition at line 65 of file OpHitFinder_module.cc.

unsigned int opdet::OpHitFinder::fMaxOpChannel
private

Definition at line 75 of file OpHitFinder_module.cc.

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

Definition at line 72 of file OpHitFinder_module.cc.

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

Definition at line 70 of file OpHitFinder_module.cc.

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

Definition at line 71 of file OpHitFinder_module.cc.

bool opdet::OpHitFinder::fUseStartTime
private

Definition at line 76 of file OpHitFinder_module.cc.


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