LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
SpacePointHit3DBuilder_tool.cc
Go to the documentation of this file.
1 
8 // Framework Includes
15 #include "art_root_io/TFileService.h"
19 #include "cetlib/cpu_timer.h"
20 #include "fhiclcpp/ParameterSet.h"
22 
23 // LArSoft includes
34 
35 // std includes
36 #include <cmath>
37 #include <iostream>
38 #include <map>
39 #include <memory>
40 #include <numeric> // std::accumulate
41 #include <unordered_map>
42 #include <utility>
43 #include <vector>
44 
45 // Ack!
46 #include "TTree.h"
47 
48 //------------------------------------------------------------------------------------------------------------------------------------------
49 // implementation follows
50 
51 namespace lar_cluster3d {
52 
56  class SpacePointHit3DBuilder : virtual public IHit3DBuilder {
57  public:
63  explicit SpacePointHit3DBuilder(fhicl::ParameterSet const& pset);
64 
69 
74  virtual void produces(art::ProducesCollector&) override;
75 
76  void configure(const fhicl::ParameterSet&) override;
77 
84  void Hit3DBuilder(art::Event& evt, reco::HitPairList& hitPairList, RecobHitToPtrMap&) override;
85 
89  float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
90  {
91  return fTimeVector.at(index);
92  }
93 
94  private:
98  void clear();
99 
103  float chargeIntegral(float, float, float, int, int) const;
104 
105  using Hit2DVector = std::vector<reco::ClusterHit2D>;
106 
116 
118  mutable std::vector<float> fTimeVector;
119 
120  // Define some basic histograms
121  TTree* m_tupleTree;
122 
123  mutable std::vector<float> m_deltaTimeVec;
124  mutable std::vector<float> m_chiSquare3DVec;
125  mutable std::vector<float> m_maxPullVec;
126  mutable std::vector<float> m_overlapFractionVec;
127  mutable std::vector<float> m_overlapRangeVec;
128  mutable std::vector<float> m_maxSideVecVec;
129  mutable std::vector<float> m_pairWireDistVec;
130  mutable std::vector<float> m_smallChargeDiffVec;
131  mutable std::vector<int> m_smallIndexVec;
132  mutable std::vector<float> m_qualityMetricVec;
133  mutable std::vector<float> m_spacePointChargeVec;
134  mutable std::vector<float> m_hitAsymmetryVec;
135 
136  // Get instances of the primary data structures needed
138 
140  };
141 
143  {
144  this->configure(pset);
145  }
146 
147  //------------------------------------------------------------------------------------------------------------------------------------------
148 
150 
152  {
153  collector.produces<std::vector<recob::Hit>>();
154 
157  }
158 
159  //------------------------------------------------------------------------------------------------------------------------------------------
160 
162  {
163  fSpacePointProducerLabel = pset.get<art::InputTag>("SpacePointProducerLabel");
164  fHitProducerLabel = pset.get<art::InputTag>("HitProducerLabel");
165  fDoWireAssns = pset.get<bool>("DoWireAssns", true);
166  fDoRawDigitAssns = pset.get<bool>("DoRawDigitAssns", true);
167  fEnableMonitoring = pset.get<bool>("EnableMonitoring", true);
168  m_maxHit3DChiSquare = pset.get<float>("MaxHitChiSquare", 6.0);
169  m_outputHistograms = pset.get<bool>("OutputHistograms", false);
170 
171  // Access ART's TFileService, which will handle creating and writing
172  // histograms and n-tuples for us.
174 
175  if (m_outputHistograms) {
176  m_tupleTree = tfs->make<TTree>("Hit3DBuilderTree", "Tree by StandardHit3DBuilder");
177 
178  clear();
179 
180  m_tupleTree->Branch("DeltaTime2D", "std::vector<float>", &m_deltaTimeVec);
181  m_tupleTree->Branch("ChiSquare3D", "std::vector<float>", &m_chiSquare3DVec);
182  m_tupleTree->Branch("MaxPullValue", "std::vector<float>", &m_maxPullVec);
183  m_tupleTree->Branch("OverlapFraction", "std::vector<float>", &m_overlapFractionVec);
184  m_tupleTree->Branch("OverlapRange", "std::vector<float>", &m_overlapRangeVec);
185  m_tupleTree->Branch("MaxSideVec", "std::vector<float>", &m_maxSideVecVec);
186  m_tupleTree->Branch("PairWireDistVec", "std::vector<float>", &m_pairWireDistVec);
187  m_tupleTree->Branch("SmallChargeDiff", "std::vector<float>", &m_smallChargeDiffVec);
188  m_tupleTree->Branch("SmallChargeIdx", "std::vector<int>", &m_smallIndexVec);
189  m_tupleTree->Branch("QualityMetric", "std::vector<float>", &m_qualityMetricVec);
190  m_tupleTree->Branch("SPCharge", "std::vector<float>", &m_spacePointChargeVec);
191  m_tupleTree->Branch("HitAsymmetry", "std::vector<float>", &m_hitAsymmetryVec);
192  }
193 
195  }
196 
198  {
199  m_deltaTimeVec.clear();
200  m_chiSquare3DVec.clear();
201  m_maxPullVec.clear();
202  m_overlapFractionVec.clear();
203  m_overlapRangeVec.clear();
204  m_maxSideVecVec.clear();
205  m_pairWireDistVec.clear();
206  m_smallChargeDiffVec.clear();
207  m_smallIndexVec.clear();
208  m_qualityMetricVec.clear();
209  m_spacePointChargeVec.clear();
210  m_hitAsymmetryVec.clear();
211 
212  return;
213  }
214 
216  reco::HitPairList& hitPairList,
217  RecobHitToPtrMap& recobHitToArtPtrMap)
218  {
223  fTimeVector.resize(NUMTIMEVALUES, 0.);
224 
225  cet::cpu_timer theClockMakeHits;
226 
227  if (fEnableMonitoring) theClockMakeHits.start();
228 
229  // Start by recovering the associations between space points and hits
231  evt.getByLabel(fSpacePointProducerLabel, hitSpacePointAssnsHandle);
232 
233  if (!hitSpacePointAssnsHandle.isValid()) return;
234 
235  // Get a hit refiner for the output hit collection
237 
238  // We need to spin through the associations first to build a map between the SpacePoints and
239  // the WireID associated to the collection plane... where for APA style TPCs this will be
240  // unambiguous (note that this is not true for ICARUS!)
241  using SpacePointToWireIDMap = std::unordered_map<const recob::SpacePoint*, geo::WireID>;
242 
243  SpacePointToWireIDMap spacePointToWireIDMap;
244 
245  for (const auto& assnPair : *hitSpacePointAssnsHandle) {
246  if (assnPair.first->SignalType() == geo::kCollection)
247  spacePointToWireIDMap[assnPair.second.get()] = assnPair.first->WireID();
248  }
249 
250  // First step is to loop through and get a mapping between space points and associated hits
251  // and, importantly, a list of unique hits (and mapping between art ptr and hit)
252  using OldHitToNewHitMap = std::map<const recob::Hit*, const recob::Hit*>;
253  using SpacePointHitVecMap = std::map<const recob::SpacePoint*, std::vector<const recob::Hit*>>;
254 
255  OldHitToNewHitMap oldHitToNewHitMap;
256  SpacePointHitVecMap spacePointHitVecMap;
257 
258  // We need a container for our new hits...
259  std::unique_ptr<std::vector<recob::Hit>> newHitVecPtr(new std::vector<recob::Hit>);
260 
261  // reserve a chunk of memory... cannot be more hits than 3 x # spacer points...
262  newHitVecPtr->reserve(3 * hitSpacePointAssnsHandle->size());
263 
264  // Use this handy art utility to make art::Ptr objects to the new recob::Hits for use in the output phase
265  art::PtrMaker<recob::Hit> ptrMaker(evt);
266 
267  for (auto& assnPair : *hitSpacePointAssnsHandle) {
268  const art::Ptr<recob::SpacePoint> spacePoint = assnPair.second;
269  const art::Ptr<recob::Hit>& recobHit = assnPair.first;
270 
271  // If we have seen this hit before then no need to create new hit
272  if (oldHitToNewHitMap.find(recobHit.get()) == oldHitToNewHitMap.end()) {
273  // Recover the reference WireID from our previous map
274  geo::WireID refWireID = spacePointToWireIDMap[spacePoint.get()];
275  geo::WireID thisWireID = recobHit.get()->WireID();
276 
277  // Recover the list of possible WireIDs from the geometry service
278  const std::vector<geo::WireID>& wireIDs =
279  fGeometry->ChannelToWire(recobHit.get()->Channel());
280 
281  // Loop to find match
282  for (const auto& wireID : wireIDs) {
283  if (wireID.TPC != refWireID.TPC || wireID.Cryostat != refWireID.Cryostat) continue;
284  thisWireID = wireID;
285  break;
286  }
287 
288  // Create and save the new recob::Hit with the correct WireID
289  newHitVecPtr->emplace_back(recob::HitCreator(*recobHit.get(), thisWireID).copy());
290 
291  // Recover a pointer to it...
292  recob::Hit* newHit = &newHitVecPtr->back();
293 
294  spacePointHitVecMap[spacePoint.get()].push_back(newHit);
295 
296  recobHitToArtPtrMap[newHit] = ptrMaker(newHitVecPtr->size() - 1);
297  oldHitToNewHitMap[recobHit.get()] = newHit;
298  }
299  else
300  spacePointHitVecMap[spacePoint.get()].push_back(oldHitToNewHitMap[recobHit.get()]);
301  }
302 
303  // We'll want to correct the hit times for the plane offsets
304  // (note this is already taken care of when converting to position)
305  std::map<geo::PlaneID, double> planeOffsetMap;
306 
307  auto const clock_data =
309  auto const det_prop =
311 
312  // Initialize the plane to hit vector map
313  for (size_t cryoIdx = 0; cryoIdx < fGeometry->Ncryostats(); cryoIdx++) {
314  for (size_t tpcIdx = 0; tpcIdx < fGeometry->NTPC(); tpcIdx++) {
315  // What we want here are the relative offsets between the planes
316  // Note that plane 0 is assumed the "first" plane and is the reference
317  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 0)] = 0.;
318  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 1)] =
319  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 1)) -
320  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0));
321  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 2)] =
322  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 2)) -
323  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0));
324 
325  std::cout << "***> plane 0 offset: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 0)]
326  << ", plane 1: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 1)]
327  << ", plane 2: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 2)] << std::endl;
328  std::cout << " Det prop plane 0: "
329  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0))
330  << ", plane 1: " << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 1))
331  << ", plane 2: " << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 2))
332  << ", Trig: " << trigger_offset(clock_data) << std::endl;
333  }
334  }
335 
336  // We need temporary mapping from recob::Hit's to our 2D hits
337  using RecobHitTo2DHitMap = std::map<const recob::Hit*, const reco::ClusterHit2D*>;
338 
339  RecobHitTo2DHitMap recobHitTo2DHitMap;
340 
341  // Set the size of the container for our hits
342  m_clusterHit2DMasterVec.clear();
343  m_clusterHit2DMasterVec.reserve(oldHitToNewHitMap.size());
344 
345  // Now go throught the list of unique hits and create the 2D hits we'll use
346  for (auto& hitPair : oldHitToNewHitMap) {
347  const recob::Hit* recobHit = hitPair.second;
348  const geo::WireID& hitWireID(recobHit->WireID());
349 
350  double hitPeakTime(
351  recobHit->PeakTime() -
352  planeOffsetMap.at(hitWireID.planeID())); //planeOffsetMap[hitWireID.planeID()]);
353  double xPosition(det_prop.ConvertTicksToX(
354  recobHit->PeakTime(), hitWireID.Plane, hitWireID.TPC, hitWireID.Cryostat));
355 
356  m_clusterHit2DMasterVec.emplace_back(0, 0., 0., xPosition, hitPeakTime, hitWireID, recobHit);
357 
358  recobHitTo2DHitMap[recobHit] = &m_clusterHit2DMasterVec.back();
359  }
360 
361  // Now we can go through the space points and build our 3D hits
362  for (auto& pointPair : spacePointHitVecMap) {
363  const recob::SpacePoint* spacePoint = pointPair.first;
364  const std::vector<const recob::Hit*>& recobHitVec = pointPair.second;
365 
366  if (recobHitVec.size() != 3) {
367  std::cout << "************>>>>>> do not have 3 hits associated to space point! "
368  << recobHitVec.size() << " ***************" << std::endl;
369  continue;
370  }
371 
372  reco::ClusterHit2DVec hitVector(recobHitVec.size());
373 
374  for (const auto& recobHit : recobHitVec) {
375  const reco::ClusterHit2D* hit2D = recobHitTo2DHitMap.at(recobHit);
376 
377  hitVector[hit2D->WireID().Plane] = hit2D;
378  }
379 
380  // Set up to get average peak time, hitChiSquare, etc.
381  unsigned int statusBits(0x7);
382  float avePeakTime(0.);
383  float weightSum(0.);
384 
385  // And get the wire IDs
386  std::vector<geo::WireID> wireIDVec = {geo::WireID(), geo::WireID(), geo::WireID()};
387 
388  // First loop through the hits to get WireIDs and calculate the averages
389  for (size_t planeIdx = 0; planeIdx < 3; planeIdx++) {
390  const reco::ClusterHit2D* hit2D = hitVector[planeIdx];
391 
392  wireIDVec[planeIdx] = hit2D->WireID();
393 
396 
398 
399  float hitRMS = hit2D->getHit()->RMS();
400  float weight = 1. / (hitRMS * hitRMS);
401  float peakTime = hit2D->getTimeTicks();
402 
403  avePeakTime += peakTime * weight;
404  weightSum += weight;
405  }
406 
407  avePeakTime /= weightSum;
408 
409  // Armed with the average peak time, now get hitChiSquare and the sig vec
410  float hitChiSquare(0.);
411  float sigmaPeakTime(std::sqrt(1. / weightSum));
412  std::vector<float> hitDelTSigVec;
413 
414  for (const auto& hit2D : hitVector) {
415  float hitRMS = hit2D->getHit()->RMS();
416  float combRMS = std::sqrt(hitRMS * hitRMS - sigmaPeakTime * sigmaPeakTime);
417  float peakTime = hit2D->getTimeTicks();
418  float deltaTime = peakTime - avePeakTime;
419  float hitSig = deltaTime / combRMS;
420 
421  hitChiSquare += hitSig * hitSig;
422 
423  hitDelTSigVec.emplace_back(std::fabs(hitSig));
424  }
425 
426  if (m_outputHistograms) m_chiSquare3DVec.push_back(hitChiSquare);
427 
428  // Need to determine the hit overlap ranges
429  int lowMinIndex(std::numeric_limits<int>::max());
430  int lowMaxIndex(std::numeric_limits<int>::min());
431  int hiMinIndex(std::numeric_limits<int>::max());
432  int hiMaxIndex(std::numeric_limits<int>::min());
433 
434  // This loop through hits to find min/max values for the common overlap region
435  for (const auto& hit2D : hitVector) {
436  int hitStart = hit2D->getHit()->PeakTime() - 2. * hit2D->getHit()->RMS() - 0.5;
437  int hitStop = hit2D->getHit()->PeakTime() + 2. * hit2D->getHit()->RMS() + 0.5;
438 
439  lowMinIndex = std::min(hitStart, lowMinIndex);
440  lowMaxIndex = std::max(hitStart, lowMaxIndex);
441  hiMinIndex = std::min(hitStop + 1, hiMinIndex);
442  hiMaxIndex = std::max(hitStop + 1, hiMaxIndex);
443  }
444 
445  // Keep only "good" hits...
446  if (hitChiSquare < m_maxHit3DChiSquare && hiMinIndex > lowMaxIndex) {
447  // One more pass through hits to get charge
448  std::vector<float> chargeVec;
449 
450  for (const auto& hit2D : hitVector)
451  chargeVec.push_back(chargeIntegral(hit2D->getHit()->PeakTime(),
452  hit2D->getHit()->PeakAmplitude(),
453  hit2D->getHit()->RMS(),
454  lowMaxIndex,
455  hiMinIndex));
456 
457  float totalCharge =
458  std::accumulate(chargeVec.begin(), chargeVec.end(), 0.) / float(chargeVec.size());
459  float overlapRange = float(hiMinIndex - lowMaxIndex);
460  float overlapFraction = overlapRange / float(hiMaxIndex - lowMinIndex);
461 
462  // Set up to compute the charge asymmetry
463  std::vector<float> smallestChargeDiffVec;
464  std::vector<float> chargeAveVec;
465  float smallestDiff(std::numeric_limits<float>::max());
466  size_t chargeIndex(0);
467 
468  for (size_t idx = 0; idx < 3; idx++) {
469  size_t leftIdx = (idx + 2) % 3;
470  size_t rightIdx = (idx + 1) % 3;
471 
472  smallestChargeDiffVec.push_back(std::abs(chargeVec[leftIdx] - chargeVec[rightIdx]));
473  chargeAveVec.push_back(float(0.5 * (chargeVec[leftIdx] + chargeVec[rightIdx])));
474 
475  if (smallestChargeDiffVec.back() < smallestDiff) {
476  smallestDiff = smallestChargeDiffVec.back();
477  chargeIndex = idx;
478  }
479 
480  // Take opportunity to look at peak time diff
481  if (m_outputHistograms) {
482  float deltaPeakTime =
483  hitVector[leftIdx]->getTimeTicks() - hitVector[rightIdx]->getTimeTicks();
484 
485  m_deltaTimeVec.push_back(deltaPeakTime);
486  }
487  }
488 
489  float chargeAsymmetry = (chargeAveVec[chargeIndex] - chargeVec[chargeIndex]) /
490  (chargeAveVec[chargeIndex] + chargeVec[chargeIndex]);
491 
492  // If this is true there has to be a negative charge that snuck in somehow
493  if (chargeAsymmetry < -1. || chargeAsymmetry > 1.) {
494  const geo::WireID& hitWireID = hitVector[chargeIndex]->WireID();
495 
496  std::cout << "============> Charge asymmetry out of range: " << chargeAsymmetry
497  << " <============" << std::endl;
498  std::cout << " hit C: " << hitWireID.Cryostat << ", TPC: " << hitWireID.TPC
499  << ", Plane: " << hitWireID.Plane << ", Wire: " << hitWireID.Wire << std::endl;
500  std::cout << " charge: " << chargeVec[0] << ", " << chargeVec[1] << ", "
501  << chargeVec[2] << std::endl;
502  std::cout << " index: " << chargeIndex << ", smallest diff: " << smallestDiff
503  << std::endl;
504  continue;
505  }
506 
507  // Usurping "deltaPeakTime" to be the maximum pull
508  float deltaPeakTime = *std::max_element(hitDelTSigVec.begin(), hitDelTSigVec.end());
509 
510  if (m_outputHistograms) {
511  m_smallChargeDiffVec.push_back(smallestDiff);
512  m_smallIndexVec.push_back(chargeIndex);
513  m_maxPullVec.push_back(deltaPeakTime);
514  m_qualityMetricVec.push_back(hitChiSquare);
515  m_spacePointChargeVec.push_back(totalCharge);
516  m_overlapFractionVec.push_back(overlapFraction);
517  m_overlapRangeVec.push_back(overlapRange);
518  m_hitAsymmetryVec.push_back(chargeAsymmetry);
519  }
520 
521  Eigen::Vector3f position(
522  float(spacePoint->XYZ()[0]), float(spacePoint->XYZ()[1]), float(spacePoint->XYZ()[2]));
523 
524  // Create the 3D cluster hit
525  hitPairList.emplace_back(0,
526  statusBits,
527  position,
528  totalCharge,
529  avePeakTime,
530  deltaPeakTime,
531  sigmaPeakTime,
532  hitChiSquare,
533  overlapFraction,
534  chargeAsymmetry,
535  0.,
536  0.,
537  hitVector,
538  hitDelTSigVec,
539  wireIDVec);
540  }
541  }
542 
543  // Now we give the new hits to the refinery
544  // Note that one advantage of using this utility is that it handles the
545  // Hit/Wire and Hit/RawDigit associations all behind the scenes for us
546  hitRefiner.use_hits(std::move(newHitVecPtr));
547 
548  // Output the new hit collection to the event
549  hitRefiner.put_into();
550 
551  // Handle tree output too
552  m_tupleTree->Fill();
553 
554  clear();
555 
556  if (fEnableMonitoring) {
557  theClockMakeHits.stop();
558 
559  fTimeVector[BUILDTHREEDHITS] = theClockMakeHits.accumulated_real_time();
560  }
561 
562  mf::LogDebug("Cluster3D") << ">>>>> 3D hit building done, found " << hitPairList.size()
563  << " 3D Hits" << std::endl;
564 
565  return;
566  }
567 
569  float peakAmp,
570  float peakSigma,
571  int low,
572  int hi) const
573  {
574  float integral(0);
575 
576  for (int sigPos = low; sigPos < hi; sigPos++) {
577  float arg = (float(sigPos) - peakMean + 0.5) / peakSigma;
578  integral += peakAmp * std::exp(-0.5 * arg * arg);
579  }
580 
581  return integral;
582  }
583 
584  //------------------------------------------------------------------------------------------------------------------------------------------
585  //------------------------------------------------------------------------------------------------------------------------------------------
586 
588 } // namespace lar_cluster3d
bool m_outputHistograms
Take the time to create and fill some histograms for diagnostics.
unsigned int NTPC(CryostatID const &cryoid=cryostat_zero) const
Returns the total number of TPCs in the specified cryostat.
Definition: GeometryCore.h:686
std::list< reco::ClusterHit3D > HitPairList
Definition: Cluster3D.h:330
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:42
art::InputTag fSpacePointProducerLabel
Data members to follow.
float getTimeTicks() const
Definition: Cluster3D.h:75
std::vector< WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
float RMS() const
RMS of the hit shape, in tick units.
Definition: Hit.h:228
float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
If monitoring, recover the time to execute a particular function.
Declaration of signal hit object.
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
float chargeIntegral(float, float, float, int, int) const
Perform charge integration between limits.
const geo::WireID & WireID() const
Definition: Cluster3D.h:76
std::vector< reco::ClusterHit2D > Hit2DVector
constexpr auto abs(T v)
Returns the absolute value of the argument.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
geo::WireID const & WireID() const
Initial tdc tick for hit.
Definition: Hit.h:280
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
void Hit3DBuilder(art::Event &evt, reco::HitPairList &hitPairList, RecobHitToPtrMap &) override
Given a set of recob hits, run DBscan to form 3D clusters.
const recob::Hit * getHit() const
Definition: Cluster3D.h:77
IDparameter< geo::PlaneID > PlaneID
Member type of validated geo::PlaneID parameter.
bool isValid() const noexcept
Definition: Handle.h:203
void clear()
clear the tuple vectors before processing next event
Class managing the creation of a new recob::Hit object.
Definition: HitCreator.h:87
Helper functions to create a hit.
void produces(std::string const &instanceName={}, Persistable const persistable=Persistable::Yes)
void use_hits(std::unique_ptr< std::vector< recob::Hit >> &&srchits)
Uses the specified collection as data product.
Definition: HitCreator.cxx:508
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
void put_into(art::Event &)
Moves the data into the event.
Definition: HitCreator.h:850
T get(std::string const &key) const
Definition: ParameterSet.h:314
A class handling a collection of hits and its associations.
Definition: HitCreator.h:787
SpacePointHit3DBuilder class definiton.
const Double32_t * XYZ() const
Definition: SpacePoint.h:78
TimeValues
enumerate the possible values for time checking if monitoring timing
Definition: IHit3DBuilder.h:73
The geometry of one entire detector, as served by art.
Definition: Geometry.h:181
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:481
Definition of data types for geometry description.
IHit3DBuilder interface class definiton.
Definition: IHit3DBuilder.h:37
This provides an art tool interface definition for tools which construct 3D hits used in 3D clusterin...
std::vector< const reco::ClusterHit2D * > ClusterHit2DVec
Definition: Cluster3D.h:90
virtual void produces(art::ProducesCollector &) override
Each algorithm may have different objects it wants "produced" so use this to let the top level produc...
double weight
Definition: plottest35.C:25
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
float m_maxHit3DChiSquare
Provide ability to select hits based on "chi square".
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
SpacePointHit3DBuilder(fhicl::ParameterSet const &pset)
Constructor.
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
int trigger_offset(DetectorClocksData const &data)
constexpr WireID()=default
Default constructor: an invalid TPC ID.
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:46
TCEvent evt
Definition: DataStructs.cxx:8
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
T const * get() const
Definition: Ptr.h:138
raw::ChannelID_t Channel() const
ID of the readout channel the hit was extracted from.
Definition: Hit.h:268
std::unordered_map< const recob::Hit *, art::Ptr< recob::Hit >> RecobHitToPtrMap
Defines a structure mapping art representation to internal.
Definition: IHit3DBuilder.h:60
art framework interface to geometry description
void configure(const fhicl::ParameterSet &) override
Interface for configuring the particular algorithm tool.
Signal from collection planes.
Definition: geo_types.h:152
unsigned getStatusBits() const
Definition: Cluster3D.h:71
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:79