LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
StandardHit3DBuilder_tool.cc
Go to the documentation of this file.
1 
8 // Framework Includes
15 #include "art_root_io/TFileDirectory.h"
16 #include "art_root_io/TFileService.h"
21 #include "cetlib/cpu_timer.h"
22 #include "fhiclcpp/ParameterSet.h"
24 
25 // LArSoft includes
30 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
31 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
33 
34 // Eigen
35 #include <Eigen/Core>
36 
37 // std includes
38 #include <iostream>
39 #include <memory>
40 #include <numeric> // std::accumulate
41 #include <string>
42 
43 // Ack!
44 #include "TH1F.h"
45 #include "TTree.h"
46 
47 //------------------------------------------------------------------------------------------------------------------------------------------
48 // implementation follows
49 
50 namespace lar_cluster3d {
51 
57  // forward declaration to define an ordering function for our hit set
58  struct Hit2DSetCompare {
59  bool operator()(const reco::ClusterHit2D*, const reco::ClusterHit2D*) const;
60  };
61 
62  using HitVector = std::vector<const reco::ClusterHit2D*>;
63  using PlaneToHitVectorMap = std::map<geo::PlaneID, HitVector>;
64  using TPCToPlaneToHitVectorMap = std::map<geo::TPCID, PlaneToHitVectorMap>;
65  using Hit2DList = std::list<reco::ClusterHit2D>;
66  using Hit2DSet = std::set<const reco::ClusterHit2D*, Hit2DSetCompare>;
67  using WireToHitSetMap = std::map<unsigned int, Hit2DSet>;
68  using PlaneToWireToHitSetMap = std::map<geo::PlaneID, WireToHitSetMap>;
69  using TPCToPlaneToWireToHitSetMap = std::map<geo::TPCID, PlaneToWireToHitSetMap>;
70  using HitVectorMap = std::map<size_t, HitVector>;
71 
72  //using HitPairVector = std::vector<std::unique_ptr<reco::ClusterHit3D>>;
73 
77  class StandardHit3DBuilder : virtual public IHit3DBuilder {
78  public:
84  explicit StandardHit3DBuilder(fhicl::ParameterSet const& pset);
85 
90  void produces(art::ProducesCollector&) override;
91 
92  void configure(const fhicl::ParameterSet&) override;
93 
100  void Hit3DBuilder(art::Event&, reco::HitPairList&, RecobHitToPtrMap&) override;
101 
105  float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
106  {
107  return m_timeVector[index];
108  }
109 
110  private:
116  void CollectArtHits(const art::Event& evt) const;
117 
121  void BuildHit3D(reco::HitPairList& hitPairList) const;
122 
126  void CreateNewRecobHitCollection(art::Event&,
128  std::vector<recob::Hit>&,
130 
134  void makeWireAssns(const art::Event&,
136  RecobHitToPtrMap&) const;
137 
141  void makeRawDigitAssns(const art::Event&,
143  RecobHitToPtrMap&) const;
144 
148  size_t BuildHitPairMap(PlaneToHitVectorMap& planeToHitVectorMap,
149  reco::HitPairList& hitPairList) const;
150 
155  std::vector<std::pair<HitVector::iterator, HitVector::iterator>>;
156 
157  size_t BuildHitPairMapByTPC(PlaneHitVectorItrPairVec& planeHitVectorItrPairVec,
158  reco::HitPairList& hitPairList) const;
159 
163  using HitMatchPair = std::pair<const reco::ClusterHit2D*, reco::ClusterHit3D>;
164  using HitMatchPairVec = std::vector<HitMatchPair>;
165  using HitMatchPairVecMap = std::map<geo::WireID, HitMatchPairVec>;
166 
167  int findGoodHitPairs(const reco::ClusterHit2D*,
170  HitMatchPairVecMap&) const;
171 
175  void findGoodTriplets(HitMatchPairVecMap&, HitMatchPairVecMap&, reco::HitPairList&) const;
176 
180  bool makeHitPair(reco::ClusterHit3D& pairOut,
181  const reco::ClusterHit2D* hit1,
182  const reco::ClusterHit2D* hit2,
183  float hitWidthSclFctr = 1.,
184  size_t hitPairCntr = 0) const;
185 
189  bool makeHitTriplet(reco::ClusterHit3D& pairOut,
190  const reco::ClusterHit3D& pairIn,
191  const reco::ClusterHit2D* hit2) const;
192 
196  bool makeDeadChannelPair(reco::ClusterHit3D& pairOut,
197  const reco::ClusterHit3D& pair,
198  size_t maxStatus,
199  size_t minStatus) const;
200 
204  const reco::ClusterHit2D* FindBestMatchingHit(const Hit2DSet& hit2DSet,
205  const reco::ClusterHit3D& pair,
206  float pairDeltaTimeLimits) const;
207 
211  int FindNumberInRange(const Hit2DSet& hit2DSet,
212  const reco::ClusterHit3D& pair,
213  float range) const;
214 
218  geo::WireID NearestWireID(const Eigen::Vector3f& position, const geo::WireID& wireID) const;
219 
223  float DistanceFromPointToHitWire(const Eigen::Vector3f& position,
224  const geo::WireID& wireID) const;
225 
229  void BuildChannelStatusVec() const;
230 
234  float chargeIntegral(float, float, float, int, int) const;
235 
239  using ChannelStatusVec = std::vector<size_t>;
240  using ChannelStatusByPlaneVec = std::vector<ChannelStatusVec>;
241 
245  void clear();
246 
250  std::vector<art::InputTag> m_hitFinderTagVec;
254  std::vector<int> m_invalidTPCVec;
255  float
259 
261  float m_wirePitch[3];
262  mutable std::vector<float> m_timeVector;
263 
265 
266  // Define some basic histograms
267  TTree* m_tupleTree;
268 
269  mutable std::vector<float> m_deltaTimeVec;
270  mutable std::vector<float> m_chiSquare3DVec;
271  mutable std::vector<float> m_maxPullVec;
272  mutable std::vector<float> m_overlapFractionVec;
273  mutable std::vector<float> m_overlapRangeVec;
274  mutable std::vector<float> m_maxDeltaPeakVec;
275  mutable std::vector<float> m_maxSideVecVec;
276  mutable std::vector<float> m_pairWireDistVec;
277  mutable std::vector<float> m_smallChargeDiffVec;
278  mutable std::vector<int> m_smallIndexVec;
279  mutable std::vector<float> m_qualityMetricVec;
280  mutable std::vector<float> m_spacePointChargeVec;
281  mutable std::vector<float> m_hitAsymmetryVec;
282 
283  // Get instances of the primary data structures needed
287 
289  mutable size_t m_numBadChannels;
290 
291  mutable bool m_weHaveAllBeenHereBefore = false;
292 
294  const lariov::ChannelStatusProvider* m_channelFilter;
295  };
296 
298  : m_geometry(art::ServiceHandle<geo::Geometry const>{}.get())
300  {
301  this->configure(pset);
302  }
303 
304  //------------------------------------------------------------------------------------------------------------------------------------------
305 
307  {
308  collector.produces<std::vector<recob::Hit>>();
311  }
312 
313  //------------------------------------------------------------------------------------------------------------------------------------------
314 
316  {
317  m_hitFinderTagVec = pset.get<std::vector<art::InputTag>>("HitFinderTagVec",
318  std::vector<art::InputTag>{"gaushit"});
319  m_enableMonitoring = pset.get<bool>("EnableMonitoring", true);
320  m_numSigmaPeakTime = pset.get<float>("NumSigmaPeakTime", 3.);
321  m_hitWidthSclFctr = pset.get<float>("HitWidthScaleFactor", 6.);
322  m_deltaPeakTimeSig = pset.get<float>("DeltaPeakTimeSig", 1.7);
323  m_zPosOffset = pset.get<float>("ZPosOffset", 0.0);
324  m_invalidTPCVec = pset.get<std::vector<int>>("InvalidTPCVec", std::vector<int>{});
325  m_wirePitchScaleFactor = pset.get<float>("WirePitchScaleFactor", 1.9);
326  m_maxHit3DChiSquare = pset.get<float>("MaxHitChiSquare", 6.0);
327  m_outputHistograms = pset.get<bool>("OutputHistograms", false);
328 
330 
331  // Returns the wire pitch per plane assuming they will be the same for all TPCs
332  constexpr geo::TPCID tpcid{0, 0};
333  m_wirePitch[0] = m_geometry->WirePitch(geo::PlaneID{tpcid, 0});
334  m_wirePitch[1] = m_geometry->WirePitch(geo::PlaneID{tpcid, 1});
335  m_wirePitch[2] = m_geometry->WirePitch(geo::PlaneID{tpcid, 2});
336 
337  // Access ART's TFileService, which will handle creating and writing
338  // histograms and n-tuples for us.
340 
341  if (m_outputHistograms) {
342  m_tupleTree = tfs->make<TTree>("Hit3DBuilderTree", "Tree by StandardHit3DBuilder");
343 
344  clear();
345 
346  m_tupleTree->Branch("DeltaTime2D", "std::vector<float>", &m_deltaTimeVec);
347  m_tupleTree->Branch("ChiSquare3D", "std::vector<float>", &m_chiSquare3DVec);
348  m_tupleTree->Branch("MaxPullValue", "std::vector<float>", &m_maxPullVec);
349  m_tupleTree->Branch("OverlapFraction", "std::vector<float>", &m_overlapFractionVec);
350  m_tupleTree->Branch("OverlapRange", "std::vector<float>", &m_overlapRangeVec);
351  m_tupleTree->Branch("MaxDeltaPeak", "std::vector<float>", &m_maxDeltaPeakVec);
352  m_tupleTree->Branch("MaxSideVec", "std::vector<float>", &m_maxSideVecVec);
353  m_tupleTree->Branch("PairWireDistVec", "std::vector<float>", &m_pairWireDistVec);
354  m_tupleTree->Branch("SmallChargeDiff", "std::vector<float>", &m_smallChargeDiffVec);
355  m_tupleTree->Branch("SmallChargeIdx", "std::vector<int>", &m_smallIndexVec);
356  m_tupleTree->Branch("QualityMetric", "std::vector<float>", &m_qualityMetricVec);
357  m_tupleTree->Branch("SPCharge", "std::vector<float>", &m_spacePointChargeVec);
358  m_tupleTree->Branch("HitAsymmetry", "std::vector<float>", &m_hitAsymmetryVec);
359  }
360 
361  return;
362  }
363 
365  {
366  m_deltaTimeVec.clear();
367  m_chiSquare3DVec.clear();
368  m_maxPullVec.clear();
369  m_overlapFractionVec.clear();
370  m_overlapRangeVec.clear();
371  m_maxDeltaPeakVec.clear();
372  m_maxSideVecVec.clear();
373  m_pairWireDistVec.clear();
374  m_smallChargeDiffVec.clear();
375  m_smallIndexVec.clear();
376  m_qualityMetricVec.clear();
377  m_spacePointChargeVec.clear();
378  m_hitAsymmetryVec.clear();
379 
380  return;
381  }
382 
384  {
385  // This is called each event, clear out the previous version and start over
386  if (!m_channelStatus.empty()) m_channelStatus.clear();
387 
388  m_numBadChannels = 0;
390 
391  // Loop through views/planes to set the wire length vectors
392  constexpr geo::TPCID tpcid{0, 0};
393  for (unsigned int idx = 0; idx < m_channelStatus.size(); idx++) {
395  }
396 
397  // Loop through the channels and mark those that are "bad"
398  for (size_t channel = 0; channel < m_geometry->Nchannels(); channel++) {
399  if (m_channelFilter->IsGood(channel)) continue;
400 
401  std::vector<geo::WireID> wireIDVec = m_geometry->ChannelToWire(channel);
402  geo::WireID wireID = wireIDVec[0];
403  lariov::ChannelStatusProvider::Status_t chanStat = m_channelFilter->Status(channel);
404 
405  m_channelStatus[wireID.Plane][wireID.Wire] = chanStat;
407  }
408  }
409 
412  {
413  return (*left).getAvePeakTime() < (*right).getAvePeakTime();
414  }
415 
416  struct HitPairClusterOrder {
419  {
420  // Watch out for the case where two clusters can have the same number of hits!
421  if (left->second.size() == right->second.size()) return left->first < right->first;
422 
423  return left->second.size() > right->second.size();
424  }
425  };
426 
428  reco::HitPairList& hitPairList,
429  RecobHitToPtrMap& clusterHitToArtPtrMap)
430  {
431  // Clear the internal data structures
432  m_clusterHit2DMasterList.clear();
433  m_planeToHitVectorMap.clear();
434  m_planeToWireToHitSetMap.clear();
435 
436  m_timeVector.resize(NUMTIMEVALUES, 0.);
437 
438  // Get a hit refiner
439  std::unique_ptr<std::vector<recob::Hit>> outputHitPtrVec(new std::vector<recob::Hit>);
440 
441  // Recover the 2D hits and then organize them into data structures which will be used in the
442  // DBscan algorithm for building the 3D clusters
443  this->CollectArtHits(evt);
444 
445  // If there are no hits in our view/wire data structure then do not proceed with the full analysis
446  if (!m_planeToWireToHitSetMap.empty()) {
447  // Call the algorithm that builds 3D hits
448  this->BuildHit3D(hitPairList);
449 
450  // If we built 3D points then attempt to output a new hit list as well
451  if (!hitPairList.empty())
452  CreateNewRecobHitCollection(evt, hitPairList, *outputHitPtrVec, clusterHitToArtPtrMap);
453  }
454 
455  // Set up to make the associations (if desired)
457  std::unique_ptr<art::Assns<recob::Wire, recob::Hit>> wireAssns(
459 
460  makeWireAssns(evt, *wireAssns, clusterHitToArtPtrMap);
461 
463  std::unique_ptr<art::Assns<raw::RawDigit, recob::Hit>> rawDigitAssns(
465 
466  makeRawDigitAssns(evt, *rawDigitAssns, clusterHitToArtPtrMap);
467 
468  // Move everything into the event
469  evt.put(std::move(outputHitPtrVec));
470  evt.put(std::move(wireAssns));
471  evt.put(std::move(rawDigitAssns));
472 
473  // Handle tree output too
474  if (m_outputHistograms) {
475  m_tupleTree->Fill();
476 
477  clear();
478  }
479  }
480 
482  {
487  cet::cpu_timer theClockMakeHits;
488 
489  if (m_enableMonitoring) theClockMakeHits.start();
490 
491  // The first task is to take the lists of input 2D hits (a map of view to sorted lists of 2D hits)
492  // and then to build a list of 3D hits to be used in downstream processing
494 
495  size_t numHitPairs = BuildHitPairMap(m_planeToHitVectorMap, hitPairList);
496 
497  if (m_enableMonitoring) {
498  theClockMakeHits.stop();
499 
500  m_timeVector[BUILDTHREEDHITS] = theClockMakeHits.accumulated_real_time();
501  }
502 
503  mf::LogDebug("Cluster3D") << ">>>>> 3D hit building done, found " << numHitPairs << " 3D Hits"
504  << std::endl;
505  }
506 
507  //------------------------------------------------------------------------------------------------------------------------------------------
508  namespace {
509 
510  class SetHitEarliestTimeOrder {
511  public:
512  SetHitEarliestTimeOrder() : m_numRMS(1.) {}
513  SetHitEarliestTimeOrder(float numRMS) : m_numRMS(numRMS) {}
514 
515  bool operator()(const reco::ClusterHit2D* left, const reco::ClusterHit2D* right) const
516  {
517  return left->getTimeTicks() - m_numRMS * left->getHit()->RMS() <
518  right->getTimeTicks() - m_numRMS * right->getHit()->RMS();
519  }
520 
521  public:
522  float m_numRMS;
523  };
524 
525  using HitVectorItrPair = std::pair<HitVector::iterator, HitVector::iterator>;
526 
527  class SetStartTimeOrder {
528  public:
529  SetStartTimeOrder() : m_numRMS(1.) {}
530  SetStartTimeOrder(float numRMS) : m_numRMS(numRMS) {}
531 
532  bool operator()(const HitVectorItrPair& left, const HitVectorItrPair& right) const
533  {
534  // Protect against possible issue?
535  if (left.first != left.second && right.first != right.second) {
536  // Sort by "modified start time" of pulse
537  return (*left.first)->getTimeTicks() - m_numRMS * (*left.first)->getHit()->RMS() <
538  (*right.first)->getTimeTicks() - m_numRMS * (*right.first)->getHit()->RMS();
539  }
540 
541  return left.first != left.second;
542  }
543 
544  private:
545  float m_numRMS;
546  };
547 
548  bool SetPairStartTimeOrder(const reco::ClusterHit3D& left, const reco::ClusterHit3D& right)
549  {
550  // Sort by "modified start time" of pulse
551  return left.getAvePeakTime() - left.getSigmaPeakTime() <
552  right.getAvePeakTime() - right.getSigmaPeakTime();
553  }
554  }
555 
556  //------------------------------------------------------------------------------------------------------------------------------------------
557 
559  reco::HitPairList& hitPairList) const
560  {
570  size_t totalNumHits(0);
571  size_t hitPairCntr(0);
572 
573  size_t nTriplets(0);
574  size_t nDeadChanHits(0);
575 
576  // Set up to loop over cryostats and tpcs...
577  for (size_t cryoIdx = 0; cryoIdx < m_geometry->Ncryostats(); cryoIdx++) {
578  for (size_t tpcIdx = 0; tpcIdx < m_geometry->NTPC(); tpcIdx++) {
580  planeToHitVectorMap.find(geo::PlaneID(cryoIdx, tpcIdx, 0));
582  planeToHitVectorMap.find(geo::PlaneID(cryoIdx, tpcIdx, 1));
584  planeToHitVectorMap.find(geo::PlaneID(cryoIdx, tpcIdx, 2));
585 
586  size_t nPlanesWithHits =
587  (mapItr0 != planeToHitVectorMap.end() && !mapItr0->second.empty() ? 1 : 0) +
588  (mapItr1 != planeToHitVectorMap.end() && !mapItr1->second.empty() ? 1 : 0) +
589  (mapItr2 != planeToHitVectorMap.end() && !mapItr2->second.empty() ? 1 : 0);
590 
591  if (nPlanesWithHits < 2) continue;
592 
593  HitVector& hitVector0 = mapItr0->second;
594  HitVector& hitVector1 = mapItr1->second;
595  HitVector& hitVector2 = mapItr2->second;
596 
597  // We are going to resort the hits into "start time" order...
598  std::sort(hitVector0.begin(),
599  hitVector0.end(),
600  SetHitEarliestTimeOrder(m_numSigmaPeakTime)); //SetHitStartTimeOrder);
601  std::sort(hitVector1.begin(),
602  hitVector1.end(),
603  SetHitEarliestTimeOrder(m_numSigmaPeakTime)); //SetHitStartTimeOrder);
604  std::sort(hitVector2.begin(),
605  hitVector2.end(),
606  SetHitEarliestTimeOrder(m_numSigmaPeakTime)); //SetHitStartTimeOrder);
607 
608  PlaneHitVectorItrPairVec hitItrVec = {
609  HitVectorItrPair(hitVector0.begin(), hitVector0.end()),
610  HitVectorItrPair(hitVector1.begin(), hitVector1.end()),
611  HitVectorItrPair(hitVector2.begin(), hitVector2.end())};
612 
613  totalNumHits += BuildHitPairMapByTPC(hitItrVec, hitPairList);
614  }
615  }
616 
617  // Return the hit pair list but sorted by z and y positions (faster traversal in next steps)
618  hitPairList.sort(SetPairStartTimeOrder);
619 
620  // Where are we?
621  mf::LogDebug("Cluster3D") << "Total number hits: " << totalNumHits << std::endl;
622  mf::LogDebug("Cluster3D") << "Created a total of " << hitPairList.size()
623  << " hit pairs, counted: " << hitPairCntr << std::endl;
624  mf::LogDebug("Cluster3D") << "-- Triplets: " << nTriplets
625  << ", dead channel pairs: " << nDeadChanHits << std::endl;
626 
627  return hitPairList.size();
628  }
629 
631  reco::HitPairList& hitPairList) const
632  {
643  // Define functions to set start/end iterators in the loop below
644  auto SetStartIterator =
645  [](HitVector::iterator startItr, HitVector::iterator endItr, float rms, float startTime) {
646  while (startItr != endItr) {
647  float numRMS(rms);
648  if ((*startItr)->getTimeTicks() + numRMS * (*startItr)->getHit()->RMS() < startTime)
649  startItr++;
650  else
651  break;
652  }
653  return startItr;
654  };
655 
656  auto SetEndIterator =
657  [](HitVector::iterator firstItr, HitVector::iterator endItr, float rms, float endTime) {
658  while (firstItr != endItr) {
659  float numRMS(rms);
660  if ((*firstItr)->getTimeTicks() - numRMS * (*firstItr)->getHit()->RMS() < endTime)
661  firstItr++;
662  else
663  break;
664  }
665  return firstItr;
666  };
667 
668  //*********************************************************************************
669  // Basically, we try to loop until done...
670  while (1) {
671  // Sort so that the earliest hit time will be the first element, etc.
672  std::sort(hitItrVec.begin(), hitItrVec.end(), SetStartTimeOrder(m_numSigmaPeakTime));
673 
674  // This loop iteration's golden hit
675  const reco::ClusterHit2D* goldenHit = *hitItrVec[0].first;
676 
677  // The range of history... (for this hit)
678  float goldenTimeStart = goldenHit->getTimeTicks() -
679  m_numSigmaPeakTime * goldenHit->getHit()->RMS() -
680  std::numeric_limits<float>::epsilon();
681  float goldenTimeEnd = goldenHit->getTimeTicks() +
682  m_numSigmaPeakTime * goldenHit->getHit()->RMS() +
683  std::numeric_limits<float>::epsilon();
684 
685  // Set iterators to insure we'll be in the overlap ranges
686  HitVector::iterator hitItr1Start = SetStartIterator(
687  hitItrVec[1].first, hitItrVec[1].second, m_numSigmaPeakTime, goldenTimeStart);
688  HitVector::iterator hitItr1End =
689  SetEndIterator(hitItr1Start, hitItrVec[1].second, m_numSigmaPeakTime, goldenTimeEnd);
690  HitVector::iterator hitItr2Start = SetStartIterator(
691  hitItrVec[2].first, hitItrVec[2].second, m_numSigmaPeakTime, goldenTimeStart);
692  HitVector::iterator hitItr2End =
693  SetEndIterator(hitItr2Start, hitItrVec[2].second, m_numSigmaPeakTime, goldenTimeEnd);
694 
695  // Since we'll use these many times in the internal loops, pre make the pairs for the second set of hits
696  HitMatchPairVecMap pair12Map;
697  HitMatchPairVecMap pair13Map;
698 
699  size_t n12Pairs = findGoodHitPairs(goldenHit, hitItr1Start, hitItr1End, pair12Map);
700  size_t n13Pairs = findGoodHitPairs(goldenHit, hitItr2Start, hitItr2End, pair13Map);
701 
702  if (n12Pairs > n13Pairs)
703  findGoodTriplets(pair12Map, pair13Map, hitPairList);
704  else
705  findGoodTriplets(pair13Map, pair12Map, hitPairList);
706 
707  hitItrVec[0].first++;
708 
709  int nPlanesWithHits(0);
710 
711  for (auto& pair : hitItrVec)
712  if (pair.first != pair.second) nPlanesWithHits++;
713 
714  if (nPlanesWithHits < 2) break;
715  }
716 
717  return hitPairList.size();
718  }
719 
721  HitVector::iterator& startItr,
722  HitVector::iterator& endItr,
723  HitMatchPairVecMap& hitMatchMap) const
724  {
725  int numPairs(0);
726 
727  // Loop through the input secon hits and make pairs
728  while (startItr != endItr) {
729  const reco::ClusterHit2D* hit = *startItr++;
730  reco::ClusterHit3D pair;
731 
732  // pair returned with a negative ave time is signal of failure
733  if (!makeHitPair(pair, goldenHit, hit, m_hitWidthSclFctr)) continue;
734 
735  geo::WireID wireID = hit->WireID();
736 
737  hitMatchMap[wireID].emplace_back(hit, pair);
738 
739  numPairs++;
740  }
741 
742  return numPairs;
743  }
744 
746  HitMatchPairVecMap& pair13Map,
747  reco::HitPairList& hitPairList) const
748  {
749  // Build triplets from the two lists of hit pairs
750  if (!pair12Map.empty()) {
751  // temporary container for dead channel hits
752  std::vector<reco::ClusterHit3D> tempDeadChanVec;
753  reco::ClusterHit3D deadChanPair;
754 
755  // Keep track of which third plane hits have been used
756  std::map<const reco::ClusterHit3D*, bool> usedPairMap;
757 
758  // Initial population of this map with the pair13Map hits
759  for (const auto& pair13 : pair13Map) {
760  for (const auto& hit2Dhit3DPair : pair13.second)
761  usedPairMap[&hit2Dhit3DPair.second] = false;
762  }
763 
764  // The outer loop is over all hit pairs made from the first two plane combinations
765  for (const auto& pair12 : pair12Map) {
766  if (pair12.second.empty()) continue;
767 
768  // This loop is over hit pairs that share the same first two plane wires but may have different
769  // hit times on those wires
770  for (const auto& hit2Dhit3DPair12 : pair12.second) {
771  const reco::ClusterHit3D& pair1 = hit2Dhit3DPair12.second;
772 
773  // populate the map with initial value
774  usedPairMap[&pair1] = false;
775 
776  // The simplest approach here is to loop over all possibilities and let the triplet builder weed out the weak candidates
777  for (const auto& pair13 : pair13Map) {
778  if (pair13.second.empty()) continue;
779 
780  for (const auto& hit2Dhit3DPair13 : pair13.second) {
781  const reco::ClusterHit2D* hit2 = hit2Dhit3DPair13.first;
782  const reco::ClusterHit3D& pair2 = hit2Dhit3DPair13.second;
783 
784  // If success try for the triplet
785  reco::ClusterHit3D triplet;
786 
787  if (makeHitTriplet(triplet, pair1, hit2)) {
788  triplet.setID(hitPairList.size());
789  hitPairList.emplace_back(triplet);
790  usedPairMap[&pair1] = true;
791  usedPairMap[&pair2] = true;
792  }
793  }
794  }
795  }
796  }
797 
798  // One more loop through the other pairs to check for sick channels
799  if (m_numBadChannels > 0) {
800  for (const auto& pairMapPair : usedPairMap) {
801  if (pairMapPair.second) continue;
802 
803  const reco::ClusterHit3D* pair = pairMapPair.first;
804 
805  // Here we look to see if we failed to make a triplet because the partner wire was dead/noisy/sick
806  if (makeDeadChannelPair(deadChanPair, *pair, 4, 0))
807  tempDeadChanVec.emplace_back(deadChanPair);
808  }
809 
810  // Handle the dead wire triplets
811  if (!tempDeadChanVec.empty()) {
812  // If we have many then see if we can trim down a bit by keeping those with time significance
813  if (tempDeadChanVec.size() > 1) {
814  // Sort by "significance" of agreement
815  std::sort(tempDeadChanVec.begin(),
816  tempDeadChanVec.end(),
817  [](const auto& left, const auto& right) {
818  return left.getDeltaPeakTime() / left.getSigmaPeakTime() <
819  right.getDeltaPeakTime() / right.getSigmaPeakTime();
820  });
821 
822  // What is the range of "significance" from first to last?
823  float firstSig = tempDeadChanVec.front().getDeltaPeakTime() /
824  tempDeadChanVec.front().getSigmaPeakTime();
825  float lastSig =
826  tempDeadChanVec.back().getDeltaPeakTime() / tempDeadChanVec.back().getSigmaPeakTime();
827  float sigRange = lastSig - firstSig;
828 
829  if (lastSig > 0.5 * m_deltaPeakTimeSig && sigRange > 0.5) {
830  // Declare a maximum of 1.5 * the average of the first and last pairs...
831  float maxSignificance = std::max(0.75 * (firstSig + lastSig), 1.0);
832 
833  std::vector<reco::ClusterHit3D>::iterator firstBadElem = std::find_if(
834  tempDeadChanVec.begin(),
835  tempDeadChanVec.end(),
836  [&maxSignificance](const auto& pair) {
837  return pair.getDeltaPeakTime() / pair.getSigmaPeakTime() > maxSignificance;
838  });
839 
840  // But only keep the best 10?
841  if (std::distance(tempDeadChanVec.begin(), firstBadElem) > 20)
842  firstBadElem = tempDeadChanVec.begin() + 20;
843  // Keep at least one hit...
844  else if (firstBadElem == tempDeadChanVec.begin())
845  firstBadElem++;
846 
847  tempDeadChanVec.resize(std::distance(tempDeadChanVec.begin(), firstBadElem));
848  }
849  }
850 
851  for (auto& pair : tempDeadChanVec) {
852  pair.setID(hitPairList.size());
853  hitPairList.emplace_back(pair);
854  }
855  }
856  }
857  }
858 
859  return;
860  }
861 
863  const reco::ClusterHit2D* hit1,
864  const reco::ClusterHit2D* hit2,
865  float hitWidthSclFctr,
866  size_t hitPairCntr) const
867  {
868  // Assume failure
869  bool result(false);
870 
871  // We assume in this routine that we are looking at hits in different views
872  // The first mission is to check that the wires intersect
873  const geo::WireID& hit1WireID = hit1->WireID();
874  const geo::WireID& hit2WireID = hit2->WireID();
875 
876  geo::WireIDIntersection widIntersect;
877 
878  if (m_geometry->WireIDsIntersect(hit1WireID, hit2WireID, widIntersect)) {
879  // Wires intersect so now we can check the timing
880  float hit1Peak = hit1->getTimeTicks();
881  float hit1Sigma = hit1->getHit()->RMS();
882 
883  float hit2Peak = hit2->getTimeTicks();
884  float hit2Sigma = hit2->getHit()->RMS();
885 
886  // "Long hits" are an issue... so we deal with these differently
887  int hit1NDF = hit1->getHit()->DegreesOfFreedom();
888  int hit2NDF = hit2->getHit()->DegreesOfFreedom();
889 
890  // Basically, allow the range to extend to the nearest end of the snippet
891  if (hit1NDF < 2)
892  hit1Sigma = std::min(hit1Peak - float(hit1->getHit()->StartTick()),
893  float(hit1->getHit()->EndTick()) - hit1Peak);
894  if (hit2NDF < 2)
895  hit2Sigma = std::min(hit2Peak - float(hit2->getHit()->StartTick()),
896  float(hit2->getHit()->EndTick()) - hit2Peak);
897 
898  // The "hit sigma" is the gaussian fit sigma of the hit, we need to expand a bit to allow hit overlap efficiency
899  float hit1Width = hitWidthSclFctr * hit1Sigma;
900  float hit2Width = hitWidthSclFctr * hit2Sigma;
901 
902  // Coarse check hit times are "in range"
903  if (fabs(hit1Peak - hit2Peak) <= (hit1Width + hit2Width)) {
904  // Check to see that hit peak times are consistent with each other
905  float hit1SigSq = hit1Sigma * hit1Sigma;
906  float hit2SigSq = hit2Sigma * hit2Sigma;
907  float deltaPeakTime = std::fabs(hit1Peak - hit2Peak);
908  float sigmaPeakTime = std::sqrt(hit1SigSq + hit2SigSq);
909 
910  // delta peak time consistency check here
911  if (deltaPeakTime < m_deltaPeakTimeSig *
912  sigmaPeakTime) // 2 sigma consistency? (do this way to avoid divide)
913  {
914  float oneOverWghts = hit1SigSq * hit2SigSq / (hit1SigSq + hit2SigSq);
915  float avePeakTime = (hit1Peak / hit1SigSq + hit2Peak / hit2SigSq) * oneOverWghts;
916  float totalCharge = hit1->getHit()->Integral() + hit2->getHit()->Integral();
917  float hitChiSquare = std::pow((hit1Peak - avePeakTime), 2) / hit1SigSq +
918  std::pow((hit2Peak - avePeakTime), 2) / hit2SigSq;
919 
920  float xPositionHit1(hit1->getXPosition());
921  float xPositionHit2(hit2->getXPosition());
922  float xPosition = (xPositionHit1 / hit1SigSq + xPositionHit2 / hit2SigSq) * hit1SigSq *
923  hit2SigSq / (hit1SigSq + hit2SigSq);
924 
925  Eigen::Vector3f position(
926  xPosition, float(widIntersect.y), float(widIntersect.z) - m_zPosOffset);
927 
928  // If to here then we need to sort out the hit pair code telling what views are used
929  unsigned statusBits = 1 << hit1->WireID().Plane | 1 << hit2->WireID().Plane;
930 
931  // handle status bits for the 2D hits
936 
939 
940  reco::ClusterHit2DVec hitVector;
941 
942  hitVector.resize(3, NULL);
943 
944  hitVector[hit1->WireID().Plane] = hit1;
945  hitVector[hit2->WireID().Plane] = hit2;
946 
947  unsigned int cryostatIdx = hit1->WireID().Cryostat;
948  unsigned int tpcIdx = hit1->WireID().TPC;
949 
950  // Initialize the wireIdVec
951  std::vector<geo::WireID> wireIDVec = {geo::WireID(cryostatIdx, tpcIdx, 0, 0),
952  geo::WireID(cryostatIdx, tpcIdx, 1, 0),
953  geo::WireID(cryostatIdx, tpcIdx, 2, 0)};
954 
955  wireIDVec[hit1->WireID().Plane] = hit1->WireID();
956  wireIDVec[hit2->WireID().Plane] = hit2->WireID();
957 
958  // For compiling at the moment
959  std::vector<float> hitDelTSigVec = {0., 0., 0.};
960 
961  hitDelTSigVec[hit1->WireID().Plane] = deltaPeakTime / sigmaPeakTime;
962  hitDelTSigVec[hit2->WireID().Plane] = deltaPeakTime / sigmaPeakTime;
963 
964  // Create the 3D cluster hit
965  hitPair.initialize(hitPairCntr,
966  statusBits,
967  position,
968  totalCharge,
969  avePeakTime,
970  deltaPeakTime,
971  sigmaPeakTime,
972  hitChiSquare,
973  0.,
974  0.,
975  0.,
976  0.,
977  hitVector,
978  hitDelTSigVec,
979  wireIDVec);
980 
981  result = true;
982  }
983  }
984  }
985 
986  // Send it back
987  return result;
988  }
989 
991  const reco::ClusterHit3D& pair,
992  const reco::ClusterHit2D* hit) const
993  {
994  // Assume failure
995  bool result(false);
996 
997  // We are going to force the wire pitch here, some time in the future we need to fix
998  static const double wirePitch =
999  0.5 * m_wirePitchScaleFactor * *std::max_element(m_wirePitch, m_wirePitch + 3);
1000 
1001  // Recover hit info
1002  float hitTimeTicks = hit->getTimeTicks();
1003  float hitSigma = hit->getHit()->RMS();
1004 
1005  // Special case check
1006  if (hitSigma > 2. * hit->getHit()->PeakAmplitude())
1007  hitSigma = 2. * hit->getHit()->PeakAmplitude();
1008 
1009  // Let's do a quick consistency check on the input hits to make sure we are in range...
1010  // Require the W hit to be "in range" with the UV Pair
1011  if (fabs(hitTimeTicks - pair.getAvePeakTime()) <
1012  m_hitWidthSclFctr * (pair.getSigmaPeakTime() + hitSigma)) {
1013  // Check the distance from the point to the wire the hit is on
1014  float hitWireDist = DistanceFromPointToHitWire(pair.getPosition(), hit->WireID());
1015 
1016  if (m_outputHistograms) m_maxSideVecVec.push_back(hitWireDist);
1017 
1018  // Reject hits that are not within range
1019  if (hitWireDist < wirePitch) {
1020  if (m_outputHistograms) m_pairWireDistVec.push_back(hitWireDist);
1021 
1022  // Use the existing code to see the U and W hits are willing to pair with the V hit
1023  reco::ClusterHit3D pair0h;
1024  reco::ClusterHit3D pair1h;
1025 
1026  // Recover all the hits involved
1027  const reco::ClusterHit2DVec& pairHitVec = pair.getHits();
1028  const reco::ClusterHit2D* hit0 = pairHitVec[0];
1029  const reco::ClusterHit2D* hit1 = pairHitVec[1];
1030 
1031  if (!hit0)
1032  hit0 = pairHitVec[2];
1033  else if (!hit1)
1034  hit1 = pairHitVec[2];
1035 
1036  // If good pairs made here then we can try to make a triplet
1037  if (makeHitPair(pair0h, hit0, hit, m_hitWidthSclFctr) &&
1038  makeHitPair(pair1h, hit1, hit, m_hitWidthSclFctr)) {
1039  // Get a copy of the input hit vector (note the order is by plane - by definition)
1040  reco::ClusterHit2DVec hitVector = pair.getHits();
1041 
1042  // include the new hit
1043  hitVector[hit->WireID().Plane] = hit;
1044 
1045  // Set up to get average peak time, hitChiSquare, etc.
1046  unsigned int statusBits(0x7);
1047  float avePeakTime(0.);
1048  float weightSum(0.);
1049  float xPosition(0.);
1050 
1051  // And get the wire IDs
1052  std::vector<geo::WireID> wireIDVec = {geo::WireID(), geo::WireID(), geo::WireID()};
1053 
1054  // First loop through the hits to get WireIDs and calculate the averages
1055  for (size_t planeIdx = 0; planeIdx < 3; planeIdx++) {
1056  const reco::ClusterHit2D* hit2D = hitVector[planeIdx];
1057 
1058  wireIDVec[planeIdx] = hit2D->WireID();
1059 
1062 
1064 
1065  float hitRMS = hit2D->getHit()->RMS();
1066  float peakTime = hit2D->getTimeTicks();
1067 
1068  // Basically, allow the range to extend to the nearest end of the snippet
1069  if (hit2D->getHit()->DegreesOfFreedom() < 2)
1070  hitRMS = std::min(hit2D->getTimeTicks() - float(hit2D->getHit()->StartTick()),
1071  float(hit2D->getHit()->EndTick()) - hit2D->getTimeTicks());
1072 
1073  float weight = 1. / (hitRMS * hitRMS);
1074 
1075  avePeakTime += peakTime * weight;
1076  xPosition += hit2D->getXPosition() * weight;
1077  weightSum += weight;
1078  }
1079 
1080  avePeakTime /= weightSum;
1081  xPosition /= weightSum;
1082 
1083  Eigen::Vector2f pair0hYZVec(pair0h.getPosition()[1], pair0h.getPosition()[2]);
1084  Eigen::Vector2f pair1hYZVec(pair1h.getPosition()[1], pair1h.getPosition()[2]);
1085  Eigen::Vector2f pairYZVec(pair.getPosition()[1], pair.getPosition()[2]);
1086  Eigen::Vector3f position(xPosition,
1087  float((pairYZVec[0] + pair0hYZVec[0] + pair1hYZVec[0]) / 3.),
1088  float((pairYZVec[1] + pair0hYZVec[1] + pair1hYZVec[1]) / 3.));
1089 
1090  // Armed with the average peak time, now get hitChiSquare and the sig vec
1091  float hitChiSquare(0.);
1092  float sigmaPeakTime(std::sqrt(1. / weightSum));
1093  std::vector<float> hitDelTSigVec;
1094 
1095  for (const auto& hit2D : hitVector) {
1096  float hitRMS = hit2D->getHit()->RMS();
1097 
1098  // Basically, allow the range to extend to the nearest end of the snippet
1099  if (hit2D->getHit()->DegreesOfFreedom() < 2)
1100  hitRMS = std::min(hit2D->getTimeTicks() - float(hit2D->getHit()->StartTick()),
1101  float(hit2D->getHit()->EndTick()) - hit2D->getTimeTicks());
1102 
1103  float combRMS = std::sqrt(hitRMS * hitRMS - sigmaPeakTime * sigmaPeakTime);
1104  float peakTime = hit2D->getTimeTicks();
1105  float deltaTime = peakTime - avePeakTime;
1106  float hitSig = deltaTime / combRMS;
1107 
1108  hitChiSquare += hitSig * hitSig;
1109 
1110  hitDelTSigVec.emplace_back(std::fabs(hitSig));
1111  }
1112 
1113  if (m_outputHistograms) m_chiSquare3DVec.push_back(hitChiSquare);
1114 
1115  int lowMinIndex(std::numeric_limits<int>::max());
1116  int lowMaxIndex(std::numeric_limits<int>::min());
1117  int hiMinIndex(std::numeric_limits<int>::max());
1118  int hiMaxIndex(std::numeric_limits<int>::min());
1119 
1120  // First task is to get the min/max values for the common overlap region
1121  for (const auto& hit2D : hitVector) {
1122  float range = 2. * hit2D->getHit()->RMS();
1123 
1124  // Basically, allow the range to extend to the nearest end of the snippet
1125  if (hit2D->getHit()->DegreesOfFreedom() < 2)
1126  range = std::min(hit2D->getTimeTicks() - float(hit2D->getHit()->StartTick()),
1127  float(hit2D->getHit()->EndTick()) - hit2D->getTimeTicks());
1128 
1129  int hitStart = hit2D->getHit()->PeakTime() - range - 0.5;
1130  int hitStop = hit2D->getHit()->PeakTime() + range + 0.5;
1131 
1132  lowMinIndex = std::min(hitStart, lowMinIndex);
1133  lowMaxIndex = std::max(hitStart, lowMaxIndex);
1134  hiMinIndex = std::min(hitStop + 1, hiMinIndex);
1135  hiMaxIndex = std::max(hitStop + 1, hiMaxIndex);
1136  }
1137 
1138  // Keep only "good" hits...
1139  if (hitChiSquare < m_maxHit3DChiSquare && hiMinIndex > lowMaxIndex) {
1140  // One more pass through hits to get charge
1141  std::vector<float> chargeVec;
1142 
1143  for (const auto& hit2D : hitVector)
1144  chargeVec.push_back(chargeIntegral(hit2D->getHit()->PeakTime(),
1145  hit2D->getHit()->PeakAmplitude(),
1146  hit2D->getHit()->RMS(),
1147  lowMaxIndex,
1148  hiMinIndex));
1149 
1150  float totalCharge =
1151  std::accumulate(chargeVec.begin(), chargeVec.end(), 0.) / float(chargeVec.size());
1152  float overlapRange = float(hiMinIndex - lowMaxIndex);
1153  float overlapFraction = overlapRange / float(hiMaxIndex - lowMinIndex);
1154 
1155  // Set up to compute the charge asymmetry
1156  std::vector<float> smallestChargeDiffVec;
1157  std::vector<float> chargeAveVec;
1158  float smallestDiff(std::numeric_limits<float>::max());
1159  float maxDeltaPeak(0.);
1160  size_t chargeIndex(0);
1161 
1162  for (size_t idx = 0; idx < 3; idx++) {
1163  size_t leftIdx = (idx + 2) % 3;
1164  size_t rightIdx = (idx + 1) % 3;
1165 
1166  smallestChargeDiffVec.push_back(std::abs(chargeVec[leftIdx] - chargeVec[rightIdx]));
1167  chargeAveVec.push_back(float(0.5 * (chargeVec[leftIdx] + chargeVec[rightIdx])));
1168 
1169  if (smallestChargeDiffVec.back() < smallestDiff) {
1170  smallestDiff = smallestChargeDiffVec.back();
1171  chargeIndex = idx;
1172  }
1173 
1174  // Take opportunity to look at peak time diff
1175  float deltaPeakTime =
1176  hitVector[leftIdx]->getTimeTicks() - hitVector[rightIdx]->getTimeTicks();
1177 
1178  if (std::abs(deltaPeakTime) > maxDeltaPeak) maxDeltaPeak = std::abs(deltaPeakTime);
1179 
1180  if (m_outputHistograms) m_deltaTimeVec.push_back(deltaPeakTime);
1181  }
1182 
1183  float chargeAsymmetry = (chargeAveVec[chargeIndex] - chargeVec[chargeIndex]) /
1184  (chargeAveVec[chargeIndex] + chargeVec[chargeIndex]);
1185 
1186  // If this is true there has to be a negative charge that snuck in somehow
1187  if (chargeAsymmetry < -1. || chargeAsymmetry > 1.) {
1188  const geo::WireID& hitWireID = hitVector[chargeIndex]->WireID();
1189 
1190  std::cout << "============> Charge asymmetry out of range: " << chargeAsymmetry
1191  << " <============" << std::endl;
1192  std::cout << " hit C: " << hitWireID.Cryostat << ", TPC: " << hitWireID.TPC
1193  << ", Plane: " << hitWireID.Plane << ", Wire: " << hitWireID.Wire
1194  << std::endl;
1195  std::cout << " charge: " << chargeVec[0] << ", " << chargeVec[1] << ", "
1196  << chargeVec[2] << std::endl;
1197  std::cout << " index: " << chargeIndex << ", smallest diff: " << smallestDiff
1198  << std::endl;
1199  return result;
1200  }
1201 
1202  // Usurping "deltaPeakTime" to be the maximum pull
1203  float deltaPeakTime = *std::max_element(hitDelTSigVec.begin(), hitDelTSigVec.end());
1204 
1205  if (m_outputHistograms) {
1206  m_smallChargeDiffVec.push_back(smallestDiff);
1207  m_smallIndexVec.push_back(chargeIndex);
1208  m_maxPullVec.push_back(deltaPeakTime);
1209  m_qualityMetricVec.push_back(hitChiSquare);
1210  m_spacePointChargeVec.push_back(totalCharge);
1211  m_overlapFractionVec.push_back(overlapFraction);
1212  m_overlapRangeVec.push_back(overlapRange);
1213  m_maxDeltaPeakVec.push_back(maxDeltaPeak);
1214  m_hitAsymmetryVec.push_back(chargeAsymmetry);
1215  }
1216 
1217  // Try to weed out cases where overlap doesn't match peak separation
1218  if (maxDeltaPeak > overlapRange) return result;
1219 
1220  // Create the 3D cluster hit
1221  hitTriplet.initialize(0,
1222  statusBits,
1223  position,
1224  totalCharge,
1225  avePeakTime,
1226  deltaPeakTime,
1227  sigmaPeakTime,
1228  hitChiSquare,
1229  overlapFraction,
1230  chargeAsymmetry,
1231  0.,
1232  0.,
1233  hitVector,
1234  hitDelTSigVec,
1235  wireIDVec);
1236 
1237  result = true;
1238  }
1239  }
1240  }
1241  }
1242 
1243  // return success/fail
1244  return result;
1245  }
1246 
1248  float peakAmp,
1249  float peakSigma,
1250  int low,
1251  int hi) const
1252  {
1253  float integral(0);
1254 
1255  for (int sigPos = low; sigPos < hi; sigPos++) {
1256  float arg = (float(sigPos) - peakMean + 0.5) / peakSigma;
1257  integral += peakAmp * std::exp(-0.5 * arg * arg);
1258  }
1259 
1260  return integral;
1261  }
1262 
1264  const reco::ClusterHit3D& pair,
1265  size_t maxChanStatus,
1266  size_t minChanStatus) const
1267  {
1268  // Assume failure (most common result)
1269  bool result(false);
1270 
1271  const reco::ClusterHit2D* hit0 = pair.getHits()[0];
1272  const reco::ClusterHit2D* hit1 = pair.getHits()[1];
1273 
1274  size_t missPlane(2);
1275 
1276  // u plane hit is missing
1277  if (!hit0) {
1278  hit0 = pair.getHits()[2];
1279  missPlane = 0;
1280  }
1281  // v plane hit is missing
1282  else if (!hit1) {
1283  hit1 = pair.getHits()[2];
1284  missPlane = 1;
1285  }
1286 
1287  // Which plane is missing?
1288  geo::WireID wireID0 = hit0->WireID();
1289  geo::WireID wireID1 = hit1->WireID();
1290 
1291  // Ok, recover the wireID expected in the third plane...
1292  geo::WireID wireIn(wireID0.Cryostat, wireID0.TPC, missPlane, 0);
1293  geo::WireID wireID = NearestWireID(pair.getPosition(), wireIn);
1294 
1295  // There can be a round off issue so check the next wire as well
1296  bool wireStatus = m_channelStatus[wireID.Plane][wireID.Wire] < maxChanStatus &&
1297  m_channelStatus[wireID.Plane][wireID.Wire] >= minChanStatus;
1298  bool wireOneStatus = m_channelStatus[wireID.Plane][wireID.Wire + 1] < maxChanStatus &&
1299  m_channelStatus[wireID.Plane][wireID.Wire + 1] >= minChanStatus;
1300 
1301  // Make sure they are of at least the minimum status
1302  if (wireStatus || wireOneStatus) {
1303  // Sort out which is the wire we're dealing with
1304  if (!wireStatus) wireID.Wire += 1;
1305 
1306  // Want to refine position since we "know" the missing wire
1307  geo::WireIDIntersection widIntersect0;
1308 
1309  if (m_geometry->WireIDsIntersect(wireID0, wireID, widIntersect0)) {
1310  geo::WireIDIntersection widIntersect1;
1311 
1312  if (m_geometry->WireIDsIntersect(wireID1, wireID, widIntersect1)) {
1313  Eigen::Vector3f newPosition(
1314  pair.getPosition()[0], pair.getPosition()[1], pair.getPosition()[2]);
1315 
1316  newPosition[1] = (newPosition[1] + widIntersect0.y + widIntersect1.y) / 3.;
1317  newPosition[2] =
1318  (newPosition[2] + widIntersect0.z + widIntersect1.z - 2. * m_zPosOffset) / 3.;
1319 
1320  pairOut = pair;
1321  pairOut.setWireID(wireID);
1322  pairOut.setPosition(newPosition);
1323 
1328 
1331 
1332  result = true;
1333  }
1334  }
1335  }
1336 
1337  return result;
1338  }
1339 
1341  const Hit2DSet& hit2DSet,
1342  const reco::ClusterHit3D& pair,
1343  float pairDeltaTimeLimits) const
1344  {
1345  static const float minCharge(0.);
1346 
1347  const reco::ClusterHit2D* bestVHit(0);
1348 
1349  float pairAvePeakTime(pair.getAvePeakTime());
1350 
1351  // Idea is to loop through the input set of hits and look for the best combination
1352  for (const auto& hit2D : hit2DSet) {
1353  if (hit2D->getHit()->Integral() < minCharge) continue;
1354 
1355  float hitVPeakTime(hit2D->getTimeTicks());
1356  float deltaPeakTime(pairAvePeakTime - hitVPeakTime);
1357 
1358  if (deltaPeakTime > pairDeltaTimeLimits) continue;
1359 
1360  if (deltaPeakTime < -pairDeltaTimeLimits) break;
1361 
1362  pairDeltaTimeLimits = fabs(deltaPeakTime);
1363  bestVHit = hit2D;
1364  }
1365 
1366  return bestVHit;
1367  }
1368 
1370  const reco::ClusterHit3D& pair,
1371  float range) const
1372  {
1373  static const float minCharge(0.);
1374 
1375  int numberInRange(0);
1376  float pairAvePeakTime(pair.getAvePeakTime());
1377 
1378  // Idea is to loop through the input set of hits and look for the best combination
1379  for (const auto& hit2D : hit2DSet) {
1380  if (hit2D->getHit()->Integral() < minCharge) continue;
1381 
1382  float hitVPeakTime(hit2D->getTimeTicks());
1383  float deltaPeakTime(pairAvePeakTime - hitVPeakTime);
1384 
1385  if (deltaPeakTime > range) continue;
1386 
1387  if (deltaPeakTime < -range) break;
1388 
1389  numberInRange++;
1390  }
1391 
1392  return numberInRange;
1393  }
1394 
1395  geo::WireID StandardHit3DBuilder::NearestWireID(const Eigen::Vector3f& position,
1396  const geo::WireID& wireIDIn) const
1397  {
1398  geo::WireID wireID = wireIDIn;
1399 
1400  // Embed the call to the geometry's services nearest wire id method in a try-catch block
1401  try {
1402  // Switch from NearestWireID to this method to avoid the roundoff error issues...
1403  double distanceToWire =
1404  m_geometry->Plane(wireIDIn).WireCoordinate(geo::vect::toPoint(position.data()));
1405 
1406  wireID.Wire = int(distanceToWire);
1407  }
1408  catch (std::exception& exc) {
1409  // This can happen, almost always because the coordinates are **just** out of range
1410  mf::LogWarning("Cluster3D") << "Exception caught finding nearest wire, position - "
1411  << exc.what() << std::endl;
1412 
1413  // Assume extremum for wire number depending on z coordinate
1414  if (position[2] < 0.5 * m_geometry->DetLength())
1415  wireID.Wire = 0;
1416  else
1417  wireID.Wire = m_geometry->Nwires(wireIDIn.asPlaneID()) - 1;
1418  }
1419 
1420  return wireID;
1421  }
1422 
1423  float StandardHit3DBuilder::DistanceFromPointToHitWire(const Eigen::Vector3f& position,
1424  const geo::WireID& wireIDIn) const
1425  {
1426  float distance;
1427 
1428  // Embed the call to the geometry's services nearest wire id method in a try-catch block
1429  try {
1430  // Get the wire endpoints
1431  Eigen::Vector3d wireStart;
1432  Eigen::Vector3d wireEnd;
1433 
1434  m_geometry->WireEndPoints(wireIDIn, &wireStart[0], &wireEnd[0]);
1435 
1436  // Want the hit position to have same x value as wire coordinates
1437  Eigen::Vector3d hitPosition(wireStart[0], position[1], position[2]);
1438 
1439  // Want the wire direction
1440  Eigen::Vector3d wireDir = wireEnd - wireStart;
1441 
1442  wireDir.normalize();
1443 
1444  // Get arc length to doca
1445  double arcLen = (hitPosition - wireStart).dot(wireDir);
1446 
1447  Eigen::Vector3d docaVec = hitPosition - (wireStart + arcLen * wireDir);
1448 
1449  distance = docaVec.norm();
1450  }
1451  catch (std::exception& exc) {
1452  // This can happen, almost always because the coordinates are **just** out of range
1453  mf::LogWarning("Cluster3D") << "Exception caught finding nearest wire, position - "
1454  << exc.what() << std::endl;
1455 
1456  // Assume extremum for wire number depending on z coordinate
1457  distance = 0.;
1458  }
1459 
1460  return distance;
1461  }
1462 
1463  //------------------------------------------------------------------------------------------------------------------------------------------
1464  bool SetHitTimeOrder(const reco::ClusterHit2D* left, const reco::ClusterHit2D* right)
1465  {
1466  // Sort by "modified start time" of pulse
1467  return left->getHit()->PeakTime() < right->getHit()->PeakTime();
1468  }
1469 
1471  const reco::ClusterHit2D* right) const
1472  {
1473  return left->getHit()->PeakTime() < right->getHit()->PeakTime();
1474  }
1475 
1476  //------------------------------------------------------------------------------------------------------------------------------------------
1478  {
1483  // Start by getting a vector of valid, non empty hit collections to make sure we really have something to do here...
1484  // Here is a container for the hits...
1485  std::vector<const recob::Hit*> recobHitVec;
1486 
1487  auto const clock_data =
1489  auto const det_prop =
1491 
1492  // Loop through the list of input sources
1493  for (const auto& inputTag : m_hitFinderTagVec) {
1494  art::Handle<std::vector<recob::Hit>> recobHitHandle;
1495  evt.getByLabel(inputTag, recobHitHandle);
1496 
1497  if (!recobHitHandle.isValid() || recobHitHandle->size() == 0) continue;
1498 
1499  recobHitVec.reserve(recobHitVec.size() + recobHitHandle->size());
1500 
1501  for (const auto& hit : *recobHitHandle)
1502  recobHitVec.push_back(&hit);
1503  }
1504 
1505  // If the vector is empty there is nothing to do
1506  if (recobHitVec.empty()) return;
1507 
1508  cet::cpu_timer theClockMakeHits;
1509 
1510  if (m_enableMonitoring) theClockMakeHits.start();
1511 
1512  // We'll want to correct the hit times for the plane offsets
1513  // (note this is already taken care of when converting to position)
1514  std::map<geo::PlaneID, double> planeOffsetMap;
1515 
1516  // Try to output a formatted string
1517  std::string debugMessage("");
1518 
1519  // Initialize the plane to hit vector map
1520  for (size_t cryoIdx = 0; cryoIdx < m_geometry->Ncryostats(); cryoIdx++) {
1521  for (size_t tpcIdx = 0; tpcIdx < m_geometry->NTPC(); tpcIdx++) {
1522  m_planeToHitVectorMap[geo::PlaneID(cryoIdx, tpcIdx, 0)] = HitVector();
1523  m_planeToHitVectorMap[geo::PlaneID(cryoIdx, tpcIdx, 1)] = HitVector();
1524  m_planeToHitVectorMap[geo::PlaneID(cryoIdx, tpcIdx, 2)] = HitVector();
1525 
1526  // What we want here are the relative offsets between the planes
1527  // Note that plane 0 is assumed the "first" plane and is the reference
1528  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 0)] = 0.;
1529  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 1)] =
1530  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 1)) -
1531  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0));
1532  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 2)] =
1533  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 2)) -
1534  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0));
1535 
1536  // Should we provide output?
1538  std::ostringstream outputString;
1539 
1540  outputString << "***> plane 0 offset: "
1541  << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 0)]
1542  << ", plane 1: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 1)]
1543  << ", plane 2: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 2)]
1544  << "\n";
1545  debugMessage += outputString.str();
1546  outputString << " Det prop plane 0: "
1547  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0))
1548  << ", plane 1: "
1549  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 1))
1550  << ", plane 2: "
1551  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 2))
1552  << ", Trig: " << trigger_offset(clock_data) << "\n";
1553  debugMessage += outputString.str();
1554  }
1555  }
1556  }
1557 
1559  mf::LogDebug("Cluster3D") << debugMessage << std::endl;
1560 
1562  }
1563 
1564  // Cycle through the recob hits to build ClusterHit2D objects and insert
1565  // them into the map
1566  for (const auto& recobHit : recobHitVec) {
1567  // Reject hits with negative charge, these are misreconstructed
1568  if (recobHit->Integral() < 0.) continue;
1569 
1570  // For some detectors we can have multiple wire ID's associated to a given channel.
1571  // So we recover the list of these wire IDs
1572  const std::vector<geo::WireID>& wireIDs = m_geometry->ChannelToWire(recobHit->Channel());
1573 
1574  // And then loop over all possible to build out our maps
1575  for (const auto& wireID : wireIDs) {
1576  // Check if this is an invalid TPC
1577  // (for example, in protoDUNE there are logical TPC's which see no signal)
1578  if (std::find(m_invalidTPCVec.begin(), m_invalidTPCVec.end(), wireID.TPC) !=
1579  m_invalidTPCVec.end())
1580  continue;
1581 
1582  // Note that a plane ID will define cryostat, TPC and plane
1583  const geo::PlaneID& planeID = wireID.planeID();
1584 
1585  double hitPeakTime(recobHit->PeakTime() - planeOffsetMap[planeID]);
1586  double xPosition(det_prop.ConvertTicksToX(
1587  recobHit->PeakTime(), planeID.Plane, planeID.TPC, planeID.Cryostat));
1588 
1589  m_clusterHit2DMasterList.emplace_back(0, 0., 0., xPosition, hitPeakTime, wireID, recobHit);
1590 
1591  m_planeToHitVectorMap[planeID].push_back(&m_clusterHit2DMasterList.back());
1592  m_planeToWireToHitSetMap[planeID][wireID.Wire].insert(&m_clusterHit2DMasterList.back());
1593  }
1594  }
1595 
1596  // Make a loop through to sort the recover hits in time order
1597  for (auto& hitVectorMap : m_planeToHitVectorMap)
1598  std::sort(hitVectorMap.second.begin(), hitVectorMap.second.end(), SetHitTimeOrder);
1599 
1600  if (m_enableMonitoring) {
1601  theClockMakeHits.stop();
1602 
1603  m_timeVector[COLLECTARTHITS] = theClockMakeHits.accumulated_real_time();
1604  }
1605 
1606  mf::LogDebug("Cluster3D") << ">>>>> Number of ART hits: " << m_clusterHit2DMasterList.size()
1607  << std::endl;
1608  }
1609 
1610  //------------------------------------------------------------------------------------------------------------------------------------------
1611 
1613  reco::HitPairList& hitPairList,
1614  std::vector<recob::Hit>& hitPtrVec,
1615  RecobHitToPtrMap& recobHitToPtrMap)
1616  {
1617  // Set up the timing
1618  cet::cpu_timer theClockBuildNewHits;
1619 
1620  if (m_enableMonitoring) theClockBuildNewHits.start();
1621 
1622  // We want to build a unique list of hits from the input 3D points which, we know, reuse 2D points frequently
1623  // At the same time we need to create a new 2D hit with the "correct" WireID and replace the old 2D hit in the
1624  // ClusterHit2D object with this new one... while keeping track of the use of the old ones. My head is spinning...
1625  // Declare a set which will allow us to keep track of those CusterHit2D objects we have seen already
1626  std::set<const reco::ClusterHit2D*> visitedHit2DSet;
1627 
1628  // Use this handy art utility to make art::Ptr objects to the new recob::Hits for use in the output phase
1629  art::PtrMaker<recob::Hit> ptrMaker(event);
1630 
1631  // Reserve enough memory to replace every recob::Hit we have considered (this is upper limit)
1632  hitPtrVec.reserve(m_clusterHit2DMasterList.size());
1633 
1634  // Scheme is to loop through all 3D hits, then through each associated ClusterHit2D object
1635  for (reco::ClusterHit3D& hit3D : hitPairList) {
1636  reco::ClusterHit2DVec& hit2DVec = hit3D.getHits();
1637 
1638  // The loop is over the index so we can recover the correct WireID to associate to the new hit when made
1639  for (size_t idx = 0; idx < hit3D.getHits().size(); idx++) {
1640  const reco::ClusterHit2D* hit2D = hit2DVec[idx];
1641 
1642  // Have we seen this 2D hit already?
1643  if (visitedHit2DSet.find(hit2D) == visitedHit2DSet.end()) {
1644  visitedHit2DSet.insert(hit2D);
1645 
1646  // Create and save the new recob::Hit with the correct WireID
1647  hitPtrVec.emplace_back(
1648  recob::HitCreator(*hit2D->getHit(), hit3D.getWireIDs()[idx]).copy());
1649 
1650  // Recover a pointer to it...
1651  recob::Hit* newHit = &hitPtrVec.back();
1652 
1653  // Create a mapping from this hit to an art Ptr representing it
1654  recobHitToPtrMap[newHit] = ptrMaker(hitPtrVec.size() - 1);
1655 
1656  // And set the pointer to this hit in the ClusterHit2D object
1657  const_cast<reco::ClusterHit2D*>(hit2D)->setHit(newHit);
1658  }
1659  }
1660  }
1661 
1662  size_t numNewHits = hitPtrVec.size();
1663 
1664  if (m_enableMonitoring) {
1665  theClockBuildNewHits.stop();
1666 
1667  m_timeVector[BUILDNEWHITS] = theClockBuildNewHits.accumulated_real_time();
1668  }
1669 
1670  mf::LogDebug("Cluster3D") << ">>>>> New output recob::Hit size: " << numNewHits << " (vs "
1671  << m_clusterHit2DMasterList.size() << " input)" << std::endl;
1672 
1673  return;
1674  }
1675 
1678  RecobHitToPtrMap& recobHitPtrMap) const
1679  {
1680  // Let's make sure the input associations container is empty
1682 
1683  // First task is to recover all of the previous wire <--> hit associations and map them by channel number
1684  // Create the temporary container
1685  std::unordered_map<raw::ChannelID_t, art::Ptr<recob::Wire>> channelToWireMap;
1686 
1687  // Go through the list of input sources and fill out the map
1688  for (const auto& inputTag : m_hitFinderTagVec) {
1690  evt.getValidHandle<std::vector<recob::Hit>>(inputTag);
1691 
1692  art::FindOneP<recob::Wire> hitToWireAssns(hitHandle, evt, inputTag);
1693 
1694  if (hitToWireAssns.isValid()) {
1695  for (size_t wireIdx = 0; wireIdx < hitToWireAssns.size(); wireIdx++) {
1696  art::Ptr<recob::Wire> wire = hitToWireAssns.at(wireIdx);
1697 
1698  channelToWireMap[wire->Channel()] = wire;
1699  }
1700  }
1701  }
1702 
1703  // Now fill the container
1704  for (const auto& hitPtrPair : recobHitPtrMap) {
1705  raw::ChannelID_t channel = hitPtrPair.first->Channel();
1706 
1707  std::unordered_map<raw::ChannelID_t, art::Ptr<recob::Wire>>::iterator chanWireItr =
1708  channelToWireMap.find(channel);
1709 
1710  if (!(chanWireItr != channelToWireMap.end())) {
1711  mf::LogDebug("Cluster3D") << "** Did not find channel to wire match! Skipping..."
1712  << std::endl;
1713  continue;
1714  }
1715 
1716  wireAssns.addSingle(chanWireItr->second, hitPtrPair.second);
1717  }
1718 
1719  return;
1720  }
1721 
1724  RecobHitToPtrMap& recobHitPtrMap) const
1725  {
1726  // Let's make sure the input associations container is empty
1727  rawDigitAssns = art::Assns<raw::RawDigit, recob::Hit>();
1728 
1729  // First task is to recover all of the previous wire <--> hit associations and map them by channel number
1730  // Create the temporary container
1731  std::unordered_map<raw::ChannelID_t, art::Ptr<raw::RawDigit>> channelToRawDigitMap;
1732 
1733  // Go through the list of input sources and fill out the map
1734  for (const auto& inputTag : m_hitFinderTagVec) {
1736  evt.getValidHandle<std::vector<recob::Hit>>(inputTag);
1737 
1738  art::FindOneP<raw::RawDigit> hitToRawDigitAssns(hitHandle, evt, inputTag);
1739 
1740  if (hitToRawDigitAssns.isValid()) {
1741  for (size_t rawDigitIdx = 0; rawDigitIdx < hitToRawDigitAssns.size(); rawDigitIdx++) {
1742  art::Ptr<raw::RawDigit> rawDigit = hitToRawDigitAssns.at(rawDigitIdx);
1743 
1744  channelToRawDigitMap[rawDigit->Channel()] = rawDigit;
1745  }
1746  }
1747  }
1748 
1749  // Now fill the container
1750  for (const auto& hitPtrPair : recobHitPtrMap) {
1751  raw::ChannelID_t channel = hitPtrPair.first->Channel();
1752 
1753  std::unordered_map<raw::ChannelID_t, art::Ptr<raw::RawDigit>>::iterator chanRawDigitItr =
1754  channelToRawDigitMap.find(channel);
1755 
1756  if (!(chanRawDigitItr != channelToRawDigitMap.end())) {
1757  mf::LogDebug("Cluster3D") << "** Did not find channel to wire match! Skipping..."
1758  << std::endl;
1759  continue;
1760  }
1761 
1762  rawDigitAssns.addSingle(chanRawDigitItr->second, hitPtrPair.second);
1763  }
1764 
1765  return;
1766  }
1767 
1768  //------------------------------------------------------------------------------------------------------------------------------------------
1769 
1771 } // namespace lar_cluster3d
void initialize(size_t id, unsigned int statusBits, const Eigen::Vector3f &position, float totalCharge, float avePeakTime, float deltaPeakTime, float sigmaPeakTime, float hitChiSquare, float overlapFraction, float chargeAsymmetry, float docaToAxis, float arclenToPoca, const ClusterHit2DVec &hitVec, const std::vector< float > &hitDelTSigVec, const std::vector< geo::WireID > &wireIDVec)
Definition: Cluster3D.cxx:156
intermediate_table::iterator iterator
bool makeHitTriplet(reco::ClusterHit3D &pairOut, const reco::ClusterHit3D &pairIn, const reco::ClusterHit2D *hit2) const
Make a 3D HitPair object by checking two hits.
void findGoodTriplets(HitMatchPairVecMap &, HitMatchPairVecMap &, reco::HitPairList &) const
This algorithm takes lists of hit pairs and finds good triplets.
float m_wirePitchScaleFactor
Scaling factor to determine max distance allowed between candidate pairs.
void setWireID(const geo::WireID &wid) const
Definition: Cluster3D.cxx:191
unsigned int NTPC(CryostatID const &cryoid=cryostat_zero) const
Returns the total number of TPCs in the specified cryostat.
Definition: GeometryCore.h:686
double z
z position of intersection
Definition: geo_types.h:792
std::list< reco::ClusterHit3D > HitPairList
Definition: Cluster3D.h:330
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:102
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:42
void CreateNewRecobHitCollection(art::Event &, reco::HitPairList &, std::vector< recob::Hit > &, RecobHitToPtrMap &)
Create a new 2D hit collection from hits associated to 3D space points.
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.
T * get() const
Definition: ServiceHandle.h:69
std::vector< size_t > ChannelStatusVec
define data structure for keeping track of channel status
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:155
bool SetHitTimeOrder(const reco::ClusterHit2D *left, const reco::ClusterHit2D *right)
float RMS() const
RMS of the hit shape, in tick units.
Definition: Hit.h:228
Declaration of signal hit object.
size_t BuildHitPairMapByTPC(PlaneHitVectorItrPairVec &planeHitVectorItrPairVec, reco::HitPairList &hitPairList) const
The data type to uniquely identify a Plane.
Definition: geo_types.h:463
bool WireIDsIntersect(WireID const &wid1, WireID const &wid2, Point_t &intersection) const
Computes the intersection between two wires.
const geo::WireID & WireID() const
Definition: Cluster3D.h:76
ChannelID_t Channel() const
DAQ channel this raw data was read from.
Definition: RawDigit.h:213
int findGoodHitPairs(const reco::ClusterHit2D *, HitVector::iterator &, HitVector::iterator &, HitMatchPairVecMap &) const
bool operator()(const reco::ClusterHit2D *, const reco::ClusterHit2D *) const
::geo::Point_t toPoint(Point const &p)
Convert the specified point into a geo::Point_t.
constexpr auto abs(T v)
Returns the absolute value of the argument.
int DegreesOfFreedom() const
Initial tdc tick for hit.
Definition: Hit.h:264
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
void BuildHit3D(reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
std::map< geo::WireID, HitMatchPairVec > HitMatchPairVecMap
bool makeHitPair(reco::ClusterHit3D &pairOut, const reco::ClusterHit2D *hit1, const reco::ClusterHit2D *hit2, float hitWidthSclFctr=1., size_t hitPairCntr=0) const
Make a HitPair object by checking two hits.
void produces(art::ProducesCollector &) override
Each algorithm may have different objects it wants "produced" so use this to let the top level produc...
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:244
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:563
double WireCoordinate(geo::Point_t const &point) const
Returns the coordinate of the point on the plane, in wire units.
Definition: PlaneGeo.h:797
float getSigmaPeakTime() const
Definition: Cluster3D.h:162
void Hit3DBuilder(art::Event &, reco::HitPairList &, RecobHitToPtrMap &) override
Given a set of recob hits, run DBscan to form 3D clusters.
std::map< geo::PlaneID, WireToHitSetMap > PlaneToWireToHitSetMap
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
Length_t DetLength(TPCID const &tpcid=tpc_zero) const
Returns the length of the active volume of the specified TPC.
bool SetPeakHitPairIteratorOrder(const reco::HitPairList::iterator &left, const reco::HitPairList::iterator &right)
const reco::ClusterHit2D * FindBestMatchingHit(const Hit2DSet &hit2DSet, const reco::ClusterHit3D &pair, float pairDeltaTimeLimits) const
A utility routine for finding a 2D hit closest in time to the given pair.
void BuildChannelStatusVec() const
Create the internal channel status vector (assume will eventually be event-by-event) ...
const recob::Hit * getHit() const
Definition: Cluster3D.h:77
unsigned int Nchannels() const
Returns the number of TPC readout channels in the detector.
float PeakAmplitude() const
The estimated amplitude of the hit at its peak, in ADC units.
Definition: Hit.h:232
IDparameter< geo::PlaneID > PlaneID
Member type of validated geo::PlaneID parameter.
float m_maxHit3DChiSquare
Provide ability to select hits based on "chi square".
bool isValid() const noexcept
Definition: Handle.h:203
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
float getAvePeakTime() const
Definition: Cluster3D.h:160
std::map< size_t, HitVector > HitVectorMap
Class managing the creation of a new recob::Hit object.
Definition: HitCreator.h:87
void makeWireAssns(const art::Event &, art::Assns< recob::Wire, recob::Hit > &, RecobHitToPtrMap &) const
Create recob::Wire to recob::Hit associations.
geo::WireID NearestWireID(const Eigen::Vector3f &position, const geo::WireID &wireID) const
Jacket the calls to finding the nearest wire in order to intercept the exceptions if out of range...
Helper functions to create a hit.
PlaneGeo const & Plane(PlaneID const &planeid) const
Returns the specified wire.
std::vector< art::InputTag > m_hitFinderTagVec
Data members to follow.
void produces(std::string const &instanceName={}, Persistable const persistable=Persistable::Yes)
std::list< reco::ClusterHit2D > Hit2DList
IDparameter< geo::WireID > WireID
Member type of validated geo::WireID parameter.
bool m_outputHistograms
Take the time to create and fill some histograms for diagnostics.
std::map< geo::PlaneID, HitVector > PlaneToHitVectorMap
raw::ChannelID_t Channel() const
Returns the ID of the channel (or InvalidChannelID)
Definition: Wire.h:223
float chargeIntegral(float, float, float, int, int) const
Perform charge integration between limits.
T get(std::string const &key) const
Definition: ParameterSet.h:314
StandardHit3DBuilder class definiton.
int FindNumberInRange(const Hit2DSet &hit2DSet, const reco::ClusterHit3D &pair, float range) const
A utility routine for returning the number of 2D hits from the list in a given range.
void CollectArtHits(const art::Event &evt) const
Extract the ART hits and the ART hit-particle relationships.
float getXPosition() const
Definition: Cluster3D.h:74
constexpr PlaneID const & asPlaneID() const
Conversion to PlaneID (for convenience of notation).
Definition: geo_types.h:520
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
The data type to uniquely identify a TPC.
Definition: geo_types.h:381
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:481
constexpr auto dot(Vector const &a, OtherVector const &b)
Return cross product of two vectors.
raw::TDCtick_t StartTick() const
Initial tdc tick for hit.
Definition: Hit.h:212
IHit3DBuilder interface class definiton.
Definition: IHit3DBuilder.h:37
Detector simulation of raw signals on wires.
bool makeDeadChannelPair(reco::ClusterHit3D &pairOut, const reco::ClusterHit3D &pair, size_t maxStatus, size_t minStatus) const
Make a 3D HitPair object from a valid pair and a dead channel in the missing plane.
float DistanceFromPointToHitWire(const Eigen::Vector3f &position, const geo::WireID &wireID) const
Jacket the calls to finding the nearest wire in order to intercept the exceptions if out of range...
void makeRawDigitAssns(const art::Event &, art::Assns< raw::RawDigit, recob::Hit > &, RecobHitToPtrMap &) const
Create raw::RawDigit to recob::Hit associations.
const_iterator end() const
Definition: Assns.h:514
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
std::set< const reco::ClusterHit2D *, Hit2DSetCompare > Hit2DSet
raw::TDCtick_t EndTick() const
Final tdc tick for hit.
Definition: Hit.h:216
const lariov::ChannelStatusProvider * m_channelFilter
std::map< geo::TPCID, PlaneToWireToHitSetMap > TPCToPlaneToWireToHitSetMap
std::vector< const reco::ClusterHit2D * > HitVector
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:94
double weight
Definition: plottest35.C:25
bool operator()(const reco::HitPairClusterMap::iterator &left, const reco::HitPairClusterMap::iterator &right)
std::map< geo::TPCID, PlaneToHitVectorMap > TPCToPlaneToHitVectorMap
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
std::pair< const reco::ClusterHit2D *, reco::ClusterHit3D > HitMatchPair
This builds a list of candidate hit pairs from lists of hits on two planes.
bool getByLabel(std::string const &label, std::string const &instance, Handle< PROD > &result) const
void setID(const size_t &id) const
Definition: Cluster3D.h:176
std::vector< ChannelStatusVec > ChannelStatusByPlaneVec
double y
y position of intersection
Definition: geo_types.h:791
std::vector< std::pair< HitVector::iterator, HitVector::iterator >> PlaneHitVectorItrPairVec
Given the ClusterHit2D objects, build the HitPairMap.
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
void addSingle(Ptr< left_t > const &left, Ptr< right_t > const &right, data_t const &data)
Definition: Assns.h:549
size_t BuildHitPairMap(PlaneToHitVectorMap &planeToHitVectorMap, reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
unsigned int Nwires(PlaneID const &planeid) const
Returns the total number of wires in the specified plane.
int trigger_offset(DetectorClocksData const &data)
Definition: MVAAlg.h:12
StandardHit3DBuilder(fhicl::ParameterSet const &pset)
Constructor.
constexpr WireID()=default
Default constructor: an invalid TPC ID.
unsigned int Nplanes(TPCID const &tpcid=tpc_zero) const
Returns the total number of planes in the specified TPC.
Definition: GeometryCore.h:977
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:46
TCEvent evt
Definition: DataStructs.cxx:8
unsigned int ChannelID_t
Type representing the ID of a readout channel.
Definition: RawTypes.h:28
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:399
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:82
Namespace collecting geometry-related classes utilities.
float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
If monitoring, recover the time to execute a particular function.
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:168
void clear()
clear the tuple vectors before processing next event
std::unordered_map< const recob::Hit *, art::Ptr< recob::Hit >> RecobHitToPtrMap
Defines a structure mapping art representation to internal.
Definition: IHit3DBuilder.h:60
void setPosition(const Eigen::Vector3f &pos) const
Definition: Cluster3D.h:183
Definition: fwd.h:26
bool SetPairStartTimeOrder(const reco::ClusterHit3D &left, const reco::ClusterHit3D &right)
void configure(const fhicl::ParameterSet &) override
Interface for configuring the particular algorithm tool.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
vec_iX clear()
Event finding and building.
void WireEndPoints(WireID const &wireid, double *xyzStart, double *xyzEnd) const
Fills two arrays with the coordinates of the wire end points.
std::map< unsigned int, Hit2DSet > WireToHitSetMap
unsigned getStatusBits() const
Definition: Cluster3D.h:71
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:79