LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
SnippetHit3DBuilder_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
31 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
32 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
34 
35 // Eigen
36 #include <Eigen/Core>
37 
38 // std includes
39 #include <iostream>
40 #include <memory>
41 #include <numeric> // std::accumulate
42 #include <string>
43 
44 // Ack!
45 #include "TH1F.h"
46 #include "TTree.h"
47 
48 //------------------------------------------------------------------------------------------------------------------------------------------
49 // implementation follows
50 
51 namespace lar_cluster3d {
52 
58  // forward declaration to define an ordering function for our hit set
59  struct Hit2DSetCompare {
60  bool operator()(const reco::ClusterHit2D*, const reco::ClusterHit2D*) const;
61  };
62 
63  using HitVector = std::vector<const reco::ClusterHit2D*>;
64  using HitStartEndPair = std::pair<raw::TDCtick_t, raw::TDCtick_t>;
65  using SnippetHitMap = std::map<HitStartEndPair, HitVector>;
66  using PlaneToSnippetHitMap = std::map<geo::PlaneID, SnippetHitMap>;
67  using TPCToPlaneToSnippetHitMap = std::map<geo::TPCID, PlaneToSnippetHitMap>;
68  using Hit2DList = std::list<reco::ClusterHit2D>;
69  using Hit2DSet = std::set<const reco::ClusterHit2D*, Hit2DSetCompare>;
70  using WireToHitSetMap = std::map<unsigned int, Hit2DSet>;
71  using PlaneToWireToHitSetMap = std::map<geo::PlaneID, WireToHitSetMap>;
72  using TPCToPlaneToWireToHitSetMap = std::map<geo::TPCID, PlaneToWireToHitSetMap>;
73  using HitVectorMap = std::map<size_t, HitVector>;
74  using SnippetHitMapItrPair = std::pair<SnippetHitMap::iterator, SnippetHitMap::iterator>;
75  using PlaneSnippetHitMapItrPairVec = std::vector<SnippetHitMapItrPair>;
76 
80  class SnippetHit3DBuilder : virtual public IHit3DBuilder {
81  public:
87  explicit SnippetHit3DBuilder(fhicl::ParameterSet const& pset);
88 
93  virtual void produces(art::ProducesCollector&) override;
94 
95  virtual void configure(const fhicl::ParameterSet&) override;
96 
103  virtual void Hit3DBuilder(art::Event&, reco::HitPairList&, RecobHitToPtrMap&) override;
104 
108  virtual float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
109  {
110  return m_timeVector[index];
111  }
112 
113  private:
119  void CollectArtHits(const art::Event& evt) const;
120 
124  void BuildHit3D(reco::HitPairList& hitPairList) const;
125 
129  void CreateNewRecobHitCollection(art::Event&,
131  std::vector<recob::Hit>&,
133 
137  void makeWireAssns(const art::Event&,
139  RecobHitToPtrMap&) const;
140 
144  void makeRawDigitAssns(const art::Event&,
146  RecobHitToPtrMap&) const;
147 
151  size_t BuildHitPairMap(PlaneToSnippetHitMap& planeToHitVectorMap,
152  reco::HitPairList& hitPairList) const;
153 
157  size_t BuildHitPairMapByTPC(PlaneSnippetHitMapItrPairVec& planeSnippetHitMapItrPairVec,
158  reco::HitPairList& hitPairList) const;
159 
163  using HitMatchTriplet =
164  std::tuple<const reco::ClusterHit2D*, const reco::ClusterHit2D*, const reco::ClusterHit3D>;
165  using HitMatchTripletVec = std::vector<HitMatchTriplet>;
166  using HitMatchTripletVecMap = std::map<geo::WireID, HitMatchTripletVec>;
167 
168  int findGoodHitPairs(SnippetHitMap::iterator&,
171  HitMatchTripletVecMap&) const;
172 
176  void findGoodTriplets(HitMatchTripletVecMap&, HitMatchTripletVecMap&, reco::HitPairList&) const;
177 
181  bool makeHitPair(reco::ClusterHit3D& pairOut,
182  const reco::ClusterHit2D* hit1,
183  const reco::ClusterHit2D* hit2,
184  float hitWidthSclFctr = 1.,
185  size_t hitPairCntr = 0) const;
186 
190  bool makeHitTriplet(reco::ClusterHit3D& pairOut,
191  const reco::ClusterHit3D& pairIn,
192  const reco::ClusterHit2D* hit2) const;
193 
197  bool makeDeadChannelPair(reco::ClusterHit3D& pairOut,
198  const reco::ClusterHit3D& pair,
199  size_t maxStatus,
200  size_t minStatus) const;
201 
206  bool WireIDsIntersect(const geo::WireID&, const geo::WireID&, geo::WireIDIntersection&) const;
207 
211  float closestApproach(const Eigen::Vector3f&,
212  const Eigen::Vector3f&,
213  const Eigen::Vector3f&,
214  const Eigen::Vector3f&,
215  float&,
216  float&) const;
217 
221  const reco::ClusterHit2D* FindBestMatchingHit(const Hit2DSet& hit2DSet,
222  const reco::ClusterHit3D& pair,
223  float pairDeltaTimeLimits) const;
224 
228  int FindNumberInRange(const Hit2DSet& hit2DSet,
229  const reco::ClusterHit3D& pair,
230  float range) const;
231 
235  geo::WireID NearestWireID(const Eigen::Vector3f& position, const geo::WireID& wireID) const;
236 
240  float DistanceFromPointToHitWire(const Eigen::Vector3f& position,
241  const geo::WireID& wireID) const;
242 
246  void BuildChannelStatusVec() const;
247 
251  float chargeIntegral(float, float, float, int, int) const;
252 
256  using ChannelStatusVec = std::vector<size_t>;
257  using ChannelStatusByPlaneVec = std::vector<ChannelStatusVec>;
258 
262  void clear();
263 
267  std::vector<art::InputTag> m_hitFinderTagVec;
274  std::vector<int> m_invalidTPCVec;
275  float
279 
281  float m_wirePitch[3];
282  mutable std::vector<float> m_timeVector;
283 
285 
286  // Define some basic histograms
287  TTree* m_tupleTree;
288 
289  mutable std::vector<float> m_deltaTimeVec;
290  mutable std::vector<float> m_chiSquare3DVec;
291  mutable std::vector<float> m_maxPullVec;
292  mutable std::vector<float> m_overlapFractionVec;
293  mutable std::vector<float> m_overlapRangeVec;
294  mutable std::vector<float> m_maxDeltaPeakVec;
295  mutable std::vector<float> m_maxSideVecVec;
296  mutable std::vector<float> m_pairWireDistVec;
297  mutable std::vector<float> m_smallChargeDiffVec;
298  mutable std::vector<int> m_smallIndexVec;
299  mutable std::vector<float> m_qualityMetricVec;
300  mutable std::vector<float> m_spacePointChargeVec;
301  mutable std::vector<float> m_hitAsymmetryVec;
302 
303  // Get instances of the primary data structures needed
307 
309  mutable size_t m_numBadChannels;
310 
311  mutable bool m_weHaveAllBeenHereBefore = false;
312 
313  const geo::Geometry* m_geometry; //< pointer to the Geometry service
314  const lariov::ChannelStatusProvider* m_channelFilter;
315  };
316 
318  : m_channelFilter(&art::ServiceHandle<lariov::ChannelStatusService const>()->GetProvider())
319 
320  {
321  this->configure(pset);
322  }
323 
324  //------------------------------------------------------------------------------------------------------------------------------------------
325 
327  {
328  collector.produces<std::vector<recob::Hit>>();
331  }
332 
333  //------------------------------------------------------------------------------------------------------------------------------------------
334 
336  {
337  m_hitFinderTagVec = pset.get<std::vector<art::InputTag>>(
338  "HitFinderTagVec", std::vector<art::InputTag>() = {"gaushit"});
339  m_enableMonitoring = pset.get<bool>("EnableMonitoring", true);
340  m_hitWidthSclFctr = pset.get<float>("HitWidthScaleFactor", 6.);
341  m_rangeNumSig = pset.get<float>("RangeNumSigma", 3.);
342  m_LongHitStretchFctr = pset.get<float>("LongHitsStretchFactor", 1.5);
343  m_pulseHeightFrac = pset.get<float>("PulseHeightFraction", 0.5);
344  m_PHLowSelection = pset.get<float>("PHLowSelection", 20.);
345  m_deltaPeakTimeSig = pset.get<float>("DeltaPeakTimeSig", 1.7);
346  m_zPosOffset = pset.get<float>("ZPosOffset", 0.0);
347  m_invalidTPCVec = pset.get<std::vector<int>>("InvalidTPCVec", std::vector<int>());
348  m_wirePitchScaleFactor = pset.get<float>("WirePitchScaleFactor", 1.9);
349  m_maxHit3DChiSquare = pset.get<float>("MaxHitChiSquare", 6.0);
350  m_outputHistograms = pset.get<bool>("OutputHistograms", false);
351 
353 
354  // Returns the wire pitch per plane assuming they will be the same for all TPCs
355  constexpr geo::TPCID tpcid{0, 0};
356  m_wirePitch[0] = m_geometry->WirePitch(geo::PlaneID{tpcid, 0});
357  m_wirePitch[1] = m_geometry->WirePitch(geo::PlaneID{tpcid, 1});
358  m_wirePitch[2] = m_geometry->WirePitch(geo::PlaneID{tpcid, 2});
359 
360  // Access ART's TFileService, which will handle creating and writing
361  // histograms and n-tuples for us.
362  if (m_outputHistograms) {
364 
365  m_tupleTree = tfs->make<TTree>("Hit3DBuilderTree", "Tree by SnippetHit3DBuilder");
366 
367  clear();
368 
369  m_tupleTree->Branch("DeltaTime2D", "std::vector<float>", &m_deltaTimeVec);
370  m_tupleTree->Branch("ChiSquare3D", "std::vector<float>", &m_chiSquare3DVec);
371  m_tupleTree->Branch("MaxPullValue", "std::vector<float>", &m_maxPullVec);
372  m_tupleTree->Branch("OverlapFraction", "std::vector<float>", &m_overlapFractionVec);
373  m_tupleTree->Branch("OverlapRange", "std::vector<float>", &m_overlapRangeVec);
374  m_tupleTree->Branch("MaxDeltaPeak", "std::vector<float>", &m_maxDeltaPeakVec);
375  m_tupleTree->Branch("MaxSideVec", "std::vector<float>", &m_maxSideVecVec);
376  m_tupleTree->Branch("PairWireDistVec", "std::vector<float>", &m_pairWireDistVec);
377  m_tupleTree->Branch("SmallChargeDiff", "std::vector<float>", &m_smallChargeDiffVec);
378  m_tupleTree->Branch("SmallChargeIdx", "std::vector<int>", &m_smallIndexVec);
379  m_tupleTree->Branch("QualityMetric", "std::vector<float>", &m_qualityMetricVec);
380  m_tupleTree->Branch("SPCharge", "std::vector<float>", &m_spacePointChargeVec);
381  m_tupleTree->Branch("HitAsymmetry", "std::vector<float>", &m_hitAsymmetryVec);
382  }
383  }
384 
386  {
387  m_deltaTimeVec.clear();
388  m_chiSquare3DVec.clear();
389  m_maxPullVec.clear();
390  m_overlapFractionVec.clear();
391  m_overlapRangeVec.clear();
392  m_maxDeltaPeakVec.clear();
393  m_maxSideVecVec.clear();
394  m_pairWireDistVec.clear();
395  m_smallChargeDiffVec.clear();
396  m_smallIndexVec.clear();
397  m_qualityMetricVec.clear();
398  m_spacePointChargeVec.clear();
399  m_hitAsymmetryVec.clear();
400  }
401 
403  {
404  // This is called each event, clear out the previous version and start over
405  m_channelStatus.clear();
406 
407  m_numBadChannels = 0;
409 
410  // Loop through views/planes to set the wire length vectors
411  constexpr geo::TPCID tpcid{0, 0};
412  for (unsigned int idx = 0; idx < m_channelStatus.size(); idx++) {
414  }
415 
416  // Loop through the channels and mark those that are "bad"
417  for (size_t channel = 0; channel < m_geometry->Nchannels(); channel++) {
418  if (m_channelFilter->IsGood(channel)) continue;
419 
420  std::vector<geo::WireID> wireIDVec = m_geometry->ChannelToWire(channel);
421  geo::WireID wireID = wireIDVec[0];
422  lariov::ChannelStatusProvider::Status_t chanStat = m_channelFilter->Status(channel);
423 
424  m_channelStatus[wireID.Plane][wireID.Wire] = chanStat;
426  }
427  }
428 
431  {
432  return (*left).getAvePeakTime() < (*right).getAvePeakTime();
433  }
434 
435  struct HitPairClusterOrder {
438  {
439  // Watch out for the case where two clusters can have the same number of hits!
440  if (left->second.size() == right->second.size()) return left->first < right->first;
441 
442  return left->second.size() > right->second.size();
443  }
444  };
445 
447  reco::HitPairList& hitPairList,
448  RecobHitToPtrMap& clusterHitToArtPtrMap)
449  {
450  // Clear the internal data structures
451  m_clusterHit2DMasterList.clear();
452  m_planeToSnippetHitMap.clear();
453  m_planeToWireToHitSetMap.clear();
454 
455  m_timeVector.resize(NUMTIMEVALUES, 0.);
456 
457  // Get a hit refiner
458  std::unique_ptr<std::vector<recob::Hit>> outputHitPtrVec(new std::vector<recob::Hit>);
459 
460  // Recover the 2D hits and then organize them into data structures which will be used in the
461  // DBscan algorithm for building the 3D clusters
462  this->CollectArtHits(evt);
463 
464  // If there are no hits in our view/wire data structure then do not proceed with the full analysis
465  if (!m_planeToWireToHitSetMap.empty()) {
466  // Call the algorithm that builds 3D hits
467  this->BuildHit3D(hitPairList);
468 
469  // If we built 3D points then attempt to output a new hit list as well
470  if (!hitPairList.empty())
471  CreateNewRecobHitCollection(evt, hitPairList, *outputHitPtrVec, clusterHitToArtPtrMap);
472  }
473 
474  // Set up to make the associations (if desired)
476  std::unique_ptr<art::Assns<recob::Wire, recob::Hit>> wireAssns(
478 
479  makeWireAssns(evt, *wireAssns, clusterHitToArtPtrMap);
480 
482  std::unique_ptr<art::Assns<raw::RawDigit, recob::Hit>> rawDigitAssns(
484 
485  makeRawDigitAssns(evt, *rawDigitAssns, clusterHitToArtPtrMap);
486 
487  // Move everything into the event
488  evt.put(std::move(outputHitPtrVec));
489  evt.put(std::move(wireAssns));
490  evt.put(std::move(rawDigitAssns));
491 
492  // Handle tree output too
493  if (m_outputHistograms) {
494  m_tupleTree->Fill();
495 
496  clear();
497  }
498  }
499 
501  {
506  cet::cpu_timer theClockMakeHits;
507 
508  if (m_enableMonitoring) theClockMakeHits.start();
509 
510  // The first task is to take the lists of input 2D hits (a map of view to sorted lists of 2D hits)
511  // and then to build a list of 3D hits to be used in downstream processing
513 
514  size_t numHitPairs = BuildHitPairMap(m_planeToSnippetHitMap, hitPairList);
515 
516  if (m_enableMonitoring) {
517  theClockMakeHits.stop();
518 
519  m_timeVector[BUILDTHREEDHITS] = theClockMakeHits.accumulated_real_time();
520  }
521 
522  mf::LogDebug("Cluster3D") << ">>>>> 3D hit building done, found " << numHitPairs << " 3D Hits"
523  << std::endl;
524  }
525 
526  //------------------------------------------------------------------------------------------------------------------------------------------
529  {
530  // Special case handling, there is nothing to compare for the left or right
531  if (left.first == left.second) return false;
532  if (right.first == right.second) return true;
533 
534  // de-referencing a bunch of pairs here...
535  return left.first->first.first < right.first->first.first;
536  }
537  };
538 
540  {
541  // Sort by "modified start time" of pulse
542  return left.getAvePeakTime() - left.getSigmaPeakTime() <
543  right.getAvePeakTime() - right.getSigmaPeakTime();
544  }
545 
546  //------------------------------------------------------------------------------------------------------------------------------------------
547 
549  reco::HitPairList& hitPairList) const
550  {
560  size_t totalNumHits(0);
561  size_t hitPairCntr(0);
562 
563  size_t nTriplets(0);
564  size_t nDeadChanHits(0);
565 
566  // Set up to loop over cryostats and tpcs...
567  for (size_t cryoIdx = 0; cryoIdx < m_geometry->Ncryostats(); cryoIdx++) {
568  for (size_t tpcIdx = 0; tpcIdx < m_geometry->NTPC(); tpcIdx++) {
570  planeToSnippetHitMap.find(geo::PlaneID(cryoIdx, tpcIdx, 0));
572  planeToSnippetHitMap.find(geo::PlaneID(cryoIdx, tpcIdx, 1));
574  planeToSnippetHitMap.find(geo::PlaneID(cryoIdx, tpcIdx, 2));
575 
576  size_t nPlanesWithHits =
577  (mapItr0 != planeToSnippetHitMap.end() && !mapItr0->second.empty() ? 1 : 0) +
578  (mapItr1 != planeToSnippetHitMap.end() && !mapItr1->second.empty() ? 1 : 0) +
579  (mapItr2 != planeToSnippetHitMap.end() && !mapItr2->second.empty() ? 1 : 0);
580 
581  if (nPlanesWithHits < 2) continue;
582 
583  SnippetHitMap& snippetHitMap0 = mapItr0->second;
584  SnippetHitMap& snippetHitMap1 = mapItr1->second;
585  SnippetHitMap& snippetHitMap2 = mapItr2->second;
586 
587  PlaneSnippetHitMapItrPairVec hitItrVec = {
588  SnippetHitMapItrPair(snippetHitMap0.begin(), snippetHitMap0.end()),
589  SnippetHitMapItrPair(snippetHitMap1.begin(), snippetHitMap1.end()),
590  SnippetHitMapItrPair(snippetHitMap2.begin(), snippetHitMap2.end())};
591 
592  totalNumHits += BuildHitPairMapByTPC(hitItrVec, hitPairList);
593  }
594  }
595 
596  // Return the hit pair list but sorted by z and y positions (faster traversal in next steps)
597  hitPairList.sort(SetPairStartTimeOrder);
598 
599  // Where are we?
600  mf::LogDebug("Cluster3D") << "Total number hits: " << totalNumHits << std::endl;
601  mf::LogDebug("Cluster3D") << "Created a total of " << hitPairList.size()
602  << " hit pairs, counted: " << hitPairCntr << std::endl;
603  mf::LogDebug("Cluster3D") << "-- Triplets: " << nTriplets
604  << ", dead channel pairs: " << nDeadChanHits << std::endl;
605 
606  return hitPairList.size();
607  }
608 
610  PlaneSnippetHitMapItrPairVec& snippetHitMapItrVec,
611  reco::HitPairList& hitPairList) const
612  {
623  // Define functions to set start/end iterators in the loop below
624  auto SetStartIterator =
625  [](SnippetHitMap::iterator startItr, SnippetHitMap::iterator endItr, float startTime) {
626  while (startItr != endItr) {
627  if (startItr->first.second < startTime)
628  startItr++;
629  else
630  break;
631  }
632  return startItr;
633  };
634 
635  auto SetEndIterator =
636  [](SnippetHitMap::iterator lastItr, SnippetHitMap::iterator endItr, float endTime) {
637  while (lastItr != endItr) {
638  if (lastItr->first.first < endTime)
639  lastItr++;
640  else
641  break;
642  }
643  return lastItr;
644  };
645 
646  //*********************************************************************************
647  // Basically, we try to loop until done...
648  while (1) {
649  // Sort so that the earliest hit time will be the first element, etc.
650  std::sort(snippetHitMapItrVec.begin(), snippetHitMapItrVec.end(), SetStartTimeOrder());
651 
652  // Make sure there are still hits on at least
653  int nPlanesWithHits(0);
654 
655  for (auto& pair : snippetHitMapItrVec)
656  if (pair.first != pair.second) nPlanesWithHits++;
657 
658  if (nPlanesWithHits < 2) break;
659 
660  // End condition: no more hit snippets
661  // if (snippetHitMapItrVec.front().first == snippetHitMapItrVec.front().second) break;
662 
663  // This loop iteration's snippet iterator
664  SnippetHitMap::iterator firstSnippetItr = snippetHitMapItrVec.front().first;
665 
666  // Set iterators to insure we'll be in the overlap ranges
667  SnippetHitMap::iterator snippetHitMapItr1Start = SetStartIterator(
668  snippetHitMapItrVec[1].first, snippetHitMapItrVec[1].second, firstSnippetItr->first.first);
669  SnippetHitMap::iterator snippetHitMapItr1End = SetEndIterator(
670  snippetHitMapItr1Start, snippetHitMapItrVec[1].second, firstSnippetItr->first.second);
671  SnippetHitMap::iterator snippetHitMapItr2Start = SetStartIterator(
672  snippetHitMapItrVec[2].first, snippetHitMapItrVec[2].second, firstSnippetItr->first.first);
673  SnippetHitMap::iterator snippetHitMapItr2End = SetEndIterator(
674  snippetHitMapItr2Start, snippetHitMapItrVec[2].second, firstSnippetItr->first.second);
675 
676  // Since we'll use these many times in the internal loops, pre make the pairs for the second set of hits
677  HitMatchTripletVecMap pair12Map;
678  HitMatchTripletVecMap pair13Map;
679 
680  size_t n12Pairs =
681  findGoodHitPairs(firstSnippetItr, snippetHitMapItr1Start, snippetHitMapItr1End, pair12Map);
682  size_t n13Pairs =
683  findGoodHitPairs(firstSnippetItr, snippetHitMapItr2Start, snippetHitMapItr2End, pair13Map);
684 
685  if (n12Pairs > n13Pairs)
686  findGoodTriplets(pair12Map, pair13Map, hitPairList);
687  else
688  findGoodTriplets(pair13Map, pair12Map, hitPairList);
689 
690  snippetHitMapItrVec.front().first++;
691  }
692 
693  return hitPairList.size();
694  }
695 
697  SnippetHitMap::iterator& startItr,
698  SnippetHitMap::iterator& endItr,
699  HitMatchTripletVecMap& hitMatchMap) const
700  {
701  int numPairs(0);
702 
703  HitVector::iterator firstMaxItr =
704  std::max_element(firstSnippetItr->second.begin(),
705  firstSnippetItr->second.end(),
706  [](const auto& left, const auto& right) {
707  return left->getHit()->PeakAmplitude() < right->getHit()->PeakAmplitude();
708  });
709  float firstPHCut = firstMaxItr != firstSnippetItr->second.end() ?
710  m_pulseHeightFrac * (*firstMaxItr)->getHit()->PeakAmplitude() :
711  4096.;
712 
713  // Loop through the hits on the first snippet
714  for (const auto& hit1 : firstSnippetItr->second) {
715  // Let's focus on the largest hit in the chain
716  if (hit1->getHit()->DegreesOfFreedom() > 1 && hit1->getHit()->PeakAmplitude() < firstPHCut &&
717  hit1->getHit()->PeakAmplitude() < m_PHLowSelection)
718  continue;
719 
720  // Inside loop iterator
721  SnippetHitMap::iterator secondHitItr = startItr;
722 
723  // Loop through the input secon hits and make pairs
724  while (secondHitItr != endItr) {
725  HitVector::iterator secondMaxItr = std::max_element(
726  secondHitItr->second.begin(),
727  secondHitItr->second.end(),
728  [](const auto& left, const auto& right) {
729  return left->getHit()->PeakAmplitude() < right->getHit()->PeakAmplitude();
730  });
731  float secondPHCut = secondMaxItr != secondHitItr->second.end() ?
732  m_pulseHeightFrac * (*secondMaxItr)->getHit()->PeakAmplitude() :
733  0.;
734 
735  for (const auto& hit2 : secondHitItr->second) {
736  // Again, focus on the large hits
737  if (hit2->getHit()->DegreesOfFreedom() > 1 &&
738  hit2->getHit()->PeakAmplitude() < secondPHCut &&
739  hit2->getHit()->PeakAmplitude() < m_PHLowSelection)
740  continue;
741 
742  reco::ClusterHit3D pair;
743 
744  // pair returned with a negative ave time is signal of failure
745  if (!makeHitPair(pair, hit1, hit2, m_hitWidthSclFctr)) continue;
746 
747  std::vector<const recob::Hit*> recobHitVec = {nullptr, nullptr, nullptr};
748 
749  recobHitVec[hit1->WireID().Plane] = hit1->getHit();
750  recobHitVec[hit2->WireID().Plane] = hit2->getHit();
751 
752  geo::WireID wireID = hit2->WireID();
753 
754  hitMatchMap[wireID].emplace_back(hit1, hit2, pair);
755 
756  numPairs++;
757  }
758 
759  secondHitItr++;
760  }
761  }
762 
763  return numPairs;
764  }
765 
767  HitMatchTripletVecMap& pair13Map,
768  reco::HitPairList& hitPairList) const
769  {
770  // Build triplets from the two lists of hit pairs
771  if (!pair12Map.empty()) {
772  // temporary container for dead channel hits
773  std::vector<reco::ClusterHit3D> tempDeadChanVec;
774  reco::ClusterHit3D deadChanPair;
775 
776  // Keep track of which third plane hits have been used
777  std::map<const reco::ClusterHit3D*, bool> usedPairMap;
778 
779  // Initial population of this map with the pair13Map hits
780  for (const auto& pair13 : pair13Map) {
781  for (const auto& hit2Dhit3DPair : pair13.second)
782  usedPairMap[&std::get<2>(hit2Dhit3DPair)] = false;
783  }
784 
785  // The outer loop is over all hit pairs made from the first two plane combinations
786  for (const auto& pair12 : pair12Map) {
787  if (pair12.second.empty()) continue;
788 
789  // This loop is over hit pairs that share the same first two plane wires but may have different
790  // hit times on those wires
791  for (const auto& hit2Dhit3DPair12 : pair12.second) {
792  const reco::ClusterHit3D& pair1 = std::get<2>(hit2Dhit3DPair12);
793 
794  // populate the map with initial value
795  usedPairMap[&pair1] = false;
796 
797  // The simplest approach here is to loop over all possibilities and let the triplet builder weed out the weak candidates
798  for (const auto& pair13 : pair13Map) {
799  if (pair13.second.empty()) continue;
800 
801  for (const auto& hit2Dhit3DPair13 : pair13.second) {
802  // Protect against double counting
803  if (std::get<0>(hit2Dhit3DPair12) != std::get<0>(hit2Dhit3DPair13)) continue;
804 
805  const reco::ClusterHit2D* hit2 = std::get<1>(hit2Dhit3DPair13);
806  const reco::ClusterHit3D& pair2 = std::get<2>(hit2Dhit3DPair13);
807 
808  // If success try for the triplet
809  reco::ClusterHit3D triplet;
810 
811  if (makeHitTriplet(triplet, pair1, hit2)) {
812  triplet.setID(hitPairList.size());
813  hitPairList.emplace_back(triplet);
814  usedPairMap[&pair1] = true;
815  usedPairMap[&pair2] = true;
816  }
817  }
818  }
819  }
820  }
821 
822  // One more loop through the other pairs to check for sick channels
823  if (m_numBadChannels > 0) {
824  for (const auto& pairMapPair : usedPairMap) {
825  if (pairMapPair.second) continue;
826 
827  const reco::ClusterHit3D* pair = pairMapPair.first;
828 
829  // Here we look to see if we failed to make a triplet because the partner wire was dead/noisy/sick
830  if (makeDeadChannelPair(deadChanPair, *pair, 4, 0))
831  tempDeadChanVec.emplace_back(deadChanPair);
832  }
833 
834  // Handle the dead wire triplets
835  if (!tempDeadChanVec.empty()) {
836  // If we have many then see if we can trim down a bit by keeping those with time significance
837  if (tempDeadChanVec.size() > 1) {
838  // Sort by "significance" of agreement
839  std::sort(tempDeadChanVec.begin(),
840  tempDeadChanVec.end(),
841  [](const auto& left, const auto& right) {
842  return left.getDeltaPeakTime() / left.getSigmaPeakTime() <
843  right.getDeltaPeakTime() / right.getSigmaPeakTime();
844  });
845 
846  // What is the range of "significance" from first to last?
847  float firstSig = tempDeadChanVec.front().getDeltaPeakTime() /
848  tempDeadChanVec.front().getSigmaPeakTime();
849  float lastSig =
850  tempDeadChanVec.back().getDeltaPeakTime() / tempDeadChanVec.back().getSigmaPeakTime();
851  float sigRange = lastSig - firstSig;
852 
853  if (lastSig > 0.5 * m_deltaPeakTimeSig && sigRange > 0.5) {
854  // Declare a maximum of 1.5 * the average of the first and last pairs...
855  float maxSignificance = std::max(0.75 * (firstSig + lastSig), 1.0);
856 
857  std::vector<reco::ClusterHit3D>::iterator firstBadElem = std::find_if(
858  tempDeadChanVec.begin(),
859  tempDeadChanVec.end(),
860  [&maxSignificance](const auto& pair) {
861  return pair.getDeltaPeakTime() / pair.getSigmaPeakTime() > maxSignificance;
862  });
863 
864  // But only keep the best 10?
865  if (std::distance(tempDeadChanVec.begin(), firstBadElem) > 20)
866  firstBadElem = tempDeadChanVec.begin() + 20;
867  // Keep at least one hit...
868  else if (firstBadElem == tempDeadChanVec.begin())
869  firstBadElem++;
870 
871  tempDeadChanVec.resize(std::distance(tempDeadChanVec.begin(), firstBadElem));
872  }
873  }
874 
875  for (auto& pair : tempDeadChanVec) {
876  pair.setID(hitPairList.size());
877  hitPairList.emplace_back(pair);
878  }
879  }
880  }
881  }
882 
883  return;
884  }
885 
887  const reco::ClusterHit2D* hit1,
888  const reco::ClusterHit2D* hit2,
889  float hitWidthSclFctr,
890  size_t hitPairCntr) const
891  {
892  // Assume failure
893  bool result(false);
894 
895  // Start by checking time consistency since this is fastest
896  // Wires intersect so now we can check the timing
897  float hit1Peak = hit1->getTimeTicks();
898  float hit1Sigma = hit1->getHit()->RMS();
899 
900  float hit2Peak = hit2->getTimeTicks();
901  float hit2Sigma = hit2->getHit()->RMS();
902 
903  // "Long hits" are an issue... so we deal with these differently
904  int hit1NDF = hit1->getHit()->DegreesOfFreedom();
905  int hit2NDF = hit2->getHit()->DegreesOfFreedom();
906 
907  // Basically, allow the range to extend to the nearest end of the snippet
908  if (hit1NDF < 2)
909  hit1Sigma *= m_LongHitStretchFctr; // This sets the range to the width of the pulse
910  if (hit2NDF < 2) hit2Sigma *= m_LongHitStretchFctr;
911 
912  // The "hit sigma" is the gaussian fit sigma of the hit, we need to expand a bit to allow hit overlap efficiency
913  float hit1Width = hitWidthSclFctr * hit1Sigma;
914  float hit2Width = hitWidthSclFctr * hit2Sigma;
915 
916  // Coarse check hit times are "in range"
917  if (fabs(hit1Peak - hit2Peak) <= (hit1Width + hit2Width)) {
918  // Check to see that hit peak times are consistent with each other
919  float hit1SigSq = hit1Sigma * hit1Sigma;
920  float hit2SigSq = hit2Sigma * hit2Sigma;
921  float deltaPeakTime = std::fabs(hit1Peak - hit2Peak);
922  float sigmaPeakTime = std::sqrt(hit1SigSq + hit2SigSq);
923 
924  // delta peak time consistency check here
925  if (deltaPeakTime <
926  m_deltaPeakTimeSig * sigmaPeakTime) // 2 sigma consistency? (do this way to avoid divide)
927  {
928  // We assume in this routine that we are looking at hits in different views
929  // The first mission is to check that the wires intersect
930  const geo::WireID& hit1WireID = hit1->WireID();
931  const geo::WireID& hit2WireID = hit2->WireID();
932 
933  geo::WireIDIntersection widIntersect;
934 
935  if (WireIDsIntersect(hit1WireID, hit2WireID, widIntersect)) {
936  float oneOverWghts = hit1SigSq * hit2SigSq / (hit1SigSq + hit2SigSq);
937  float avePeakTime = (hit1Peak / hit1SigSq + hit2Peak / hit2SigSq) * oneOverWghts;
938  float totalCharge = hit1->getHit()->Integral() + hit2->getHit()->Integral();
939  float hitChiSquare = std::pow((hit1Peak - avePeakTime), 2) / hit1SigSq +
940  std::pow((hit2Peak - avePeakTime), 2) / hit2SigSq;
941 
942  float xPositionHit1(hit1->getXPosition());
943  float xPositionHit2(hit2->getXPosition());
944  float xPosition = (xPositionHit1 / hit1SigSq + xPositionHit2 / hit2SigSq) * hit1SigSq *
945  hit2SigSq / (hit1SigSq + hit2SigSq);
946 
947  Eigen::Vector3f position(
948  xPosition, float(widIntersect.y), float(widIntersect.z) - m_zPosOffset);
949 
950  // If to here then we need to sort out the hit pair code telling what views are used
951  unsigned statusBits = 1 << hit1->WireID().Plane | 1 << hit2->WireID().Plane;
952 
953  // handle status bits for the 2D hits
958 
961 
962  reco::ClusterHit2DVec hitVector;
963 
964  hitVector.resize(3, NULL);
965 
966  hitVector[hit1->WireID().Plane] = hit1;
967  hitVector[hit2->WireID().Plane] = hit2;
968 
969  unsigned int cryostatIdx = hit1->WireID().Cryostat;
970  unsigned int tpcIdx = hit1->WireID().TPC;
971 
972  // Initialize the wireIdVec
973  std::vector<geo::WireID> wireIDVec = {geo::WireID(cryostatIdx, tpcIdx, 0, 0),
974  geo::WireID(cryostatIdx, tpcIdx, 1, 0),
975  geo::WireID(cryostatIdx, tpcIdx, 2, 0)};
976 
977  wireIDVec[hit1->WireID().Plane] = hit1->WireID();
978  wireIDVec[hit2->WireID().Plane] = hit2->WireID();
979 
980  // For compiling at the moment
981  std::vector<float> hitDelTSigVec = {0., 0., 0.};
982 
983  hitDelTSigVec[hit1->WireID().Plane] = deltaPeakTime / sigmaPeakTime;
984  hitDelTSigVec[hit2->WireID().Plane] = deltaPeakTime / sigmaPeakTime;
985 
986  // Create the 3D cluster hit
987  hitPair.initialize(hitPairCntr,
988  statusBits,
989  position,
990  totalCharge,
991  avePeakTime,
992  deltaPeakTime,
993  sigmaPeakTime,
994  hitChiSquare,
995  0.,
996  0.,
997  0.,
998  0.,
999  hitVector,
1000  hitDelTSigVec,
1001  wireIDVec);
1002 
1003  result = true;
1004  }
1005  }
1006  // else std::cout << "-MakeHitPair, deltaPeakTime: " << deltaPeakTime << ", scl fctr: " << m_deltaPeakTimeSig << ", sigmaPeakTime: " << sigmaPeakTime << std::endl;
1007  }
1008  // else std::cout << "-MakeHitPair, delta peak: " << hit1Peak - hit2Peak << ", hit1Width: " << hit1Width << ", hit2Width: " << hit2Width << std::endl;
1009 
1010  // Send it back
1011  return result;
1012  }
1013 
1015  const reco::ClusterHit3D& pair,
1016  const reco::ClusterHit2D* hit) const
1017  {
1018  // Assume failure
1019  bool result(false);
1020 
1021  // We are going to force the wire pitch here, some time in the future we need to fix
1022  static const double wirePitch =
1023  0.5 * m_wirePitchScaleFactor * *std::max_element(m_wirePitch, m_wirePitch + 3);
1024 
1025  // Recover hit info
1026  float hitTimeTicks = hit->getTimeTicks();
1027  float hitSigma = hit->getHit()->RMS();
1028 
1029  // Special case check
1030  if (hitSigma > 2. * hit->getHit()->PeakAmplitude())
1031  hitSigma = 2. * hit->getHit()->PeakAmplitude();
1032 
1033  // Let's do a quick consistency check on the input hits to make sure we are in range...
1034  // Require the W hit to be "in range" with the UV Pair
1035  if (fabs(hitTimeTicks - pair.getAvePeakTime()) <
1036  m_hitWidthSclFctr * (pair.getSigmaPeakTime() + hitSigma)) {
1037  // Check the distance from the point to the wire the hit is on
1038  float hitWireDist = DistanceFromPointToHitWire(pair.getPosition(), hit->WireID());
1039 
1040  if (m_outputHistograms) m_maxSideVecVec.push_back(hitWireDist);
1041 
1042  // Reject hits that are not within range
1043  if (hitWireDist < wirePitch) {
1044  if (m_outputHistograms) m_pairWireDistVec.push_back(hitWireDist);
1045 
1046  // Use the existing code to see the U and W hits are willing to pair with the V hit
1047  reco::ClusterHit3D pair0h;
1048  reco::ClusterHit3D pair1h;
1049 
1050  // Recover all the hits involved
1051  const reco::ClusterHit2DVec& pairHitVec = pair.getHits();
1052  const reco::ClusterHit2D* hit0 = pairHitVec[0];
1053  const reco::ClusterHit2D* hit1 = pairHitVec[1];
1054 
1055  if (!hit0)
1056  hit0 = pairHitVec[2];
1057  else if (!hit1)
1058  hit1 = pairHitVec[2];
1059 
1060  // If good pairs made here then we can try to make a triplet
1061  if (makeHitPair(pair0h, hit0, hit, m_hitWidthSclFctr) &&
1062  makeHitPair(pair1h, hit1, hit, m_hitWidthSclFctr)) {
1063  // Get a copy of the input hit vector (note the order is by plane - by definition)
1064  reco::ClusterHit2DVec hitVector = pair.getHits();
1065 
1066  // include the new hit
1067  hitVector[hit->WireID().Plane] = hit;
1068 
1069  // Set up to get average peak time, hitChiSquare, etc.
1070  unsigned int statusBits(0x7);
1071  float avePeakTime(0.);
1072  float weightSum(0.);
1073  float xPosition(0.);
1074 
1075  // And get the wire IDs
1076  std::vector<geo::WireID> wireIDVec = {geo::WireID(), geo::WireID(), geo::WireID()};
1077 
1078  // First loop through the hits to get WireIDs and calculate the averages
1079  for (size_t planeIdx = 0; planeIdx < 3; planeIdx++) {
1080  const reco::ClusterHit2D* hit2D = hitVector[planeIdx];
1081 
1082  wireIDVec[planeIdx] = hit2D->WireID();
1083 
1086 
1088 
1089  float hitRMS = hit2D->getHit()->RMS();
1090  float peakTime = hit2D->getTimeTicks();
1091 
1092  // Basically, allow the range to extend to the nearest end of the snippet
1093  if (hit2D->getHit()->DegreesOfFreedom() < 2) hitRMS *= m_LongHitStretchFctr;
1094  //hitRMS = std::min(hit2D->getTimeTicks() - float(hit2D->getHit()->StartTick()),float(hit2D->getHit()->EndTick())-hit2D->getTimeTicks());
1095 
1096  float weight = 1. / (hitRMS * hitRMS);
1097 
1098  avePeakTime += peakTime * weight;
1099  xPosition += hit2D->getXPosition() * weight;
1100  weightSum += weight;
1101  }
1102 
1103  avePeakTime /= weightSum;
1104  xPosition /= weightSum;
1105 
1106  Eigen::Vector2f pair0hYZVec(pair0h.getPosition()[1], pair0h.getPosition()[2]);
1107  Eigen::Vector2f pair1hYZVec(pair1h.getPosition()[1], pair1h.getPosition()[2]);
1108  Eigen::Vector2f pairYZVec(pair.getPosition()[1], pair.getPosition()[2]);
1109  Eigen::Vector3f position(xPosition,
1110  float((pairYZVec[0] + pair0hYZVec[0] + pair1hYZVec[0]) / 3.),
1111  float((pairYZVec[1] + pair0hYZVec[1] + pair1hYZVec[1]) / 3.));
1112 
1113  // Armed with the average peak time, now get hitChiSquare and the sig vec
1114  float hitChiSquare(0.);
1115  float sigmaPeakTime(std::sqrt(1. / weightSum));
1116  std::vector<float> hitDelTSigVec;
1117 
1118  for (const auto& hit2D : hitVector) {
1119  float hitRMS = hit2D->getHit()->RMS();
1120 
1121  // Basically, allow the range to extend to the nearest end of the snippet
1122  if (hit2D->getHit()->DegreesOfFreedom() < 2) hitRMS *= m_LongHitStretchFctr;
1123  //hitRMS = std::min(hit2D->getTimeTicks() - float(hit2D->getHit()->StartTick()),float(hit2D->getHit()->EndTick())-hit2D->getTimeTicks());
1124 
1125  float combRMS = std::sqrt(hitRMS * hitRMS - sigmaPeakTime * sigmaPeakTime);
1126  float peakTime = hit2D->getTimeTicks();
1127  float deltaTime = peakTime - avePeakTime;
1128  float hitSig = deltaTime / combRMS;
1129 
1130  hitChiSquare += hitSig * hitSig;
1131 
1132  hitDelTSigVec.emplace_back(std::fabs(hitSig));
1133  }
1134 
1135  if (m_outputHistograms) m_chiSquare3DVec.push_back(hitChiSquare);
1136 
1137  int lowMinIndex(std::numeric_limits<int>::max());
1138  int lowMaxIndex(std::numeric_limits<int>::min());
1139  int hiMinIndex(std::numeric_limits<int>::max());
1140  int hiMaxIndex(std::numeric_limits<int>::min());
1141 
1142  // First task is to get the min/max values for the common overlap region
1143  for (const auto& hit2D : hitVector) {
1144  float range = m_rangeNumSig * hit2D->getHit()->RMS();
1145 
1146  // Basically, allow the range to extend to the nearest end of the snippet
1147  if (hit2D->getHit()->DegreesOfFreedom() < 2) range *= m_LongHitStretchFctr;
1148  //range = std::min(hit2D->getTimeTicks() - float(hit2D->getHit()->StartTick()),float(hit2D->getHit()->EndTick())-hit2D->getTimeTicks());
1149 
1150  int hitStart = hit2D->getHit()->PeakTime() - range - 0.5;
1151  int hitStop = hit2D->getHit()->PeakTime() + range + 0.5;
1152 
1153  lowMinIndex = std::min(hitStart, lowMinIndex);
1154  lowMaxIndex = std::max(hitStart, lowMaxIndex);
1155  hiMinIndex = std::min(hitStop + 1, hiMinIndex);
1156  hiMaxIndex = std::max(hitStop + 1, hiMaxIndex);
1157  }
1158 
1159  // Keep only "good" hits...
1160  if (hitChiSquare < m_maxHit3DChiSquare && hiMinIndex > lowMaxIndex) {
1161  // One more pass through hits to get charge
1162  std::vector<float> chargeVec;
1163 
1164  for (const auto& hit2D : hitVector)
1165  chargeVec.push_back(chargeIntegral(hit2D->getHit()->PeakTime(),
1166  hit2D->getHit()->PeakAmplitude(),
1167  hit2D->getHit()->RMS(),
1168  lowMaxIndex,
1169  hiMinIndex));
1170 
1171  float totalCharge =
1172  std::accumulate(chargeVec.begin(), chargeVec.end(), 0.) / float(chargeVec.size());
1173  float overlapRange = float(hiMinIndex - lowMaxIndex);
1174  float overlapFraction = overlapRange / float(hiMaxIndex - lowMinIndex);
1175 
1176  // Set up to compute the charge asymmetry
1177  std::vector<float> smallestChargeDiffVec;
1178  std::vector<float> chargeAveVec;
1179  float smallestDiff(std::numeric_limits<float>::max());
1180  float maxDeltaPeak(0.);
1181  size_t chargeIndex(0);
1182 
1183  for (size_t idx = 0; idx < 3; idx++) {
1184  size_t leftIdx = (idx + 2) % 3;
1185  size_t rightIdx = (idx + 1) % 3;
1186 
1187  smallestChargeDiffVec.push_back(std::abs(chargeVec[leftIdx] - chargeVec[rightIdx]));
1188  chargeAveVec.push_back(float(0.5 * (chargeVec[leftIdx] + chargeVec[rightIdx])));
1189 
1190  if (smallestChargeDiffVec.back() < smallestDiff) {
1191  smallestDiff = smallestChargeDiffVec.back();
1192  chargeIndex = idx;
1193  }
1194 
1195  // Take opportunity to look at peak time diff
1196  float deltaPeakTime =
1197  hitVector[leftIdx]->getTimeTicks() - hitVector[rightIdx]->getTimeTicks();
1198 
1199  if (std::abs(deltaPeakTime) > maxDeltaPeak) maxDeltaPeak = std::abs(deltaPeakTime);
1200 
1201  if (m_outputHistograms) m_deltaTimeVec.push_back(deltaPeakTime);
1202  }
1203 
1204  float chargeAsymmetry = (chargeAveVec[chargeIndex] - chargeVec[chargeIndex]) /
1205  (chargeAveVec[chargeIndex] + chargeVec[chargeIndex]);
1206 
1207  // If this is true there has to be a negative charge that snuck in somehow
1208  if (chargeAsymmetry < -1. || chargeAsymmetry > 1.) {
1209  const geo::WireID& hitWireID = hitVector[chargeIndex]->WireID();
1210 
1211  std::cout << "============> Charge asymmetry out of range: " << chargeAsymmetry
1212  << " <============" << std::endl;
1213  std::cout << " hit C: " << hitWireID.Cryostat << ", TPC: " << hitWireID.TPC
1214  << ", Plane: " << hitWireID.Plane << ", Wire: " << hitWireID.Wire
1215  << std::endl;
1216  std::cout << " charge: " << chargeVec[0] << ", " << chargeVec[1] << ", "
1217  << chargeVec[2] << std::endl;
1218  std::cout << " index: " << chargeIndex << ", smallest diff: " << smallestDiff
1219  << std::endl;
1220  return result;
1221  }
1222 
1223  // Usurping "deltaPeakTime" to be the maximum pull
1224  float deltaPeakTime = *std::max_element(hitDelTSigVec.begin(), hitDelTSigVec.end());
1225 
1226  if (m_outputHistograms) {
1227  m_smallChargeDiffVec.push_back(smallestDiff);
1228  m_smallIndexVec.push_back(chargeIndex);
1229  m_maxPullVec.push_back(deltaPeakTime);
1230  m_qualityMetricVec.push_back(hitChiSquare);
1231  m_spacePointChargeVec.push_back(totalCharge);
1232  m_overlapFractionVec.push_back(overlapFraction);
1233  m_overlapRangeVec.push_back(overlapRange);
1234  m_maxDeltaPeakVec.push_back(maxDeltaPeak);
1235  m_hitAsymmetryVec.push_back(chargeAsymmetry);
1236  }
1237 
1238  // Try to weed out cases where overlap doesn't match peak separation
1239  if (maxDeltaPeak > 3. * overlapRange) return result;
1240 
1241  // Create the 3D cluster hit
1242  hitTriplet.initialize(0,
1243  statusBits,
1244  position,
1245  totalCharge,
1246  avePeakTime,
1247  deltaPeakTime,
1248  sigmaPeakTime,
1249  hitChiSquare,
1250  overlapFraction,
1251  chargeAsymmetry,
1252  0.,
1253  0.,
1254  hitVector,
1255  hitDelTSigVec,
1256  wireIDVec);
1257 
1258  result = true;
1259  }
1260  // else std::cout << "-Rejecting triple with chiSquare: " << hitChiSquare << " and hiMinIndex: " << hiMinIndex << ", loMaxIndex: " << lowMaxIndex << std::endl;
1261  }
1262  }
1263  }
1264  // else std::cout << "-MakeTriplet hit cut, delta: " << hitTimeTicks - pair.getAvePeakTime() << ", min scale fctr: " <<m_hitWidthSclFctr << ", pair sig: " << pair.getSigmaPeakTime() << ", hitSigma: " << hitSigma << std::endl;
1265 
1266  // return success/fail
1267  return result;
1268  }
1269 
1271  const geo::WireID& wireID1,
1272  geo::WireIDIntersection& widIntersection) const
1273  {
1274  bool success(false);
1275 
1276  // Do quick check that things are in the same logical TPC
1277  if (wireID0.Cryostat != wireID1.Cryostat || wireID0.TPC != wireID1.TPC ||
1278  wireID0.Plane == wireID1.Plane)
1279  return success;
1280 
1281  // Recover wire geometry information for each wire
1282  const geo::WireGeo& wireGeo0 = m_geometry->WireIDToWireGeo(wireID0);
1283  const geo::WireGeo& wireGeo1 = m_geometry->WireIDToWireGeo(wireID1);
1284 
1285  // Get wire position and direction for first wire
1286  auto wirePosArr = wireGeo0.GetCenter();
1287 
1288  Eigen::Vector3f wirePos0(wirePosArr.X(), wirePosArr.Y(), wirePosArr.Z());
1289  Eigen::Vector3f wireDir0(
1290  wireGeo0.Direction().X(), wireGeo0.Direction().Y(), wireGeo0.Direction().Z());
1291 
1292  //*********************************
1293  // Kludge
1294  // if (wireID0.Plane > 0) wireDir0[2] = -wireDir0[2];
1295 
1296  // And now the second one
1297  wirePosArr = wireGeo1.GetCenter();
1298 
1299  Eigen::Vector3f wirePos1(wirePosArr.X(), wirePosArr.Y(), wirePosArr.Z());
1300  Eigen::Vector3f wireDir1(
1301  wireGeo1.Direction().X(), wireGeo1.Direction().Y(), wireGeo1.Direction().Z());
1302 
1303  //**********************************
1304  // Kludge
1305  // if (wireID1.Plane > 0) wireDir1[2] = -wireDir1[2];
1306 
1307  // Get the distance of closest approach
1308  float arcLen0;
1309  float arcLen1;
1310 
1311  if (closestApproach(wirePos0, wireDir0, wirePos1, wireDir1, arcLen0, arcLen1)) {
1312  // Now check that arc lengths are within range
1313  if (std::abs(arcLen0) < wireGeo0.HalfL() && std::abs(arcLen1) < wireGeo1.HalfL()) {
1314  Eigen::Vector3f poca0 = wirePos0 + arcLen0 * wireDir0;
1315 
1316  widIntersection.y = poca0[1];
1317  widIntersection.z = poca0[2];
1318 
1319  success = true;
1320  }
1321  }
1322 
1323  return success;
1324  }
1325 
1326  float SnippetHit3DBuilder::closestApproach(const Eigen::Vector3f& P0,
1327  const Eigen::Vector3f& u0,
1328  const Eigen::Vector3f& P1,
1329  const Eigen::Vector3f& u1,
1330  float& arcLen0,
1331  float& arcLen1) const
1332  {
1333  // Technique is to compute the arclength to each point of closest approach
1334  Eigen::Vector3f w0 = P0 - P1;
1335  float a(1.);
1336  float b(u0.dot(u1));
1337  float c(1.);
1338  float d(u0.dot(w0));
1339  float e(u1.dot(w0));
1340  float den(a * c - b * b);
1341 
1342  arcLen0 = (b * e - c * d) / den;
1343  arcLen1 = (a * e - b * d) / den;
1344 
1345  Eigen::Vector3f poca0 = P0 + arcLen0 * u0;
1346  Eigen::Vector3f poca1 = P1 + arcLen1 * u1;
1347 
1348  return (poca0 - poca1).norm();
1349  }
1350 
1352  float peakAmp,
1353  float peakSigma,
1354  int low,
1355  int hi) const
1356  {
1357  float integral(0);
1358 
1359  for (int sigPos = low; sigPos < hi; sigPos++) {
1360  float arg = (float(sigPos) - peakMean + 0.5) / peakSigma;
1361  integral += peakAmp * std::exp(-0.5 * arg * arg);
1362  }
1363 
1364  return integral;
1365  }
1366 
1368  const reco::ClusterHit3D& pair,
1369  size_t maxChanStatus,
1370  size_t minChanStatus) const
1371  {
1372  // Assume failure (most common result)
1373  bool result(false);
1374 
1375  const reco::ClusterHit2D* hit0 = pair.getHits()[0];
1376  const reco::ClusterHit2D* hit1 = pair.getHits()[1];
1377 
1378  size_t missPlane(2);
1379 
1380  // u plane hit is missing
1381  if (!hit0) {
1382  hit0 = pair.getHits()[2];
1383  missPlane = 0;
1384  }
1385  // v plane hit is missing
1386  else if (!hit1) {
1387  hit1 = pair.getHits()[2];
1388  missPlane = 1;
1389  }
1390 
1391  // Which plane is missing?
1392  geo::WireID wireID0 = hit0->WireID();
1393  geo::WireID wireID1 = hit1->WireID();
1394 
1395  // Ok, recover the wireID expected in the third plane...
1396  geo::WireID wireIn(wireID0.Cryostat, wireID0.TPC, missPlane, 0);
1397  geo::WireID wireID = NearestWireID(pair.getPosition(), wireIn);
1398 
1399  // There can be a round off issue so check the next wire as well
1400  bool wireStatus = m_channelStatus[wireID.Plane][wireID.Wire] < maxChanStatus &&
1401  m_channelStatus[wireID.Plane][wireID.Wire] >= minChanStatus;
1402  bool wireOneStatus = m_channelStatus[wireID.Plane][wireID.Wire + 1] < maxChanStatus &&
1403  m_channelStatus[wireID.Plane][wireID.Wire + 1] >= minChanStatus;
1404 
1405  // Make sure they are of at least the minimum status
1406  if (wireStatus || wireOneStatus) {
1407  // Sort out which is the wire we're dealing with
1408  if (!wireStatus) wireID.Wire += 1;
1409 
1410  // Want to refine position since we "know" the missing wire
1411  geo::WireIDIntersection widIntersect0;
1412 
1413  if (m_geometry->WireIDsIntersect(wireID0, wireID, widIntersect0)) {
1414  geo::WireIDIntersection widIntersect1;
1415 
1416  if (m_geometry->WireIDsIntersect(wireID1, wireID, widIntersect1)) {
1417  Eigen::Vector3f newPosition(
1418  pair.getPosition()[0], pair.getPosition()[1], pair.getPosition()[2]);
1419 
1420  newPosition[1] = (newPosition[1] + widIntersect0.y + widIntersect1.y) / 3.;
1421  newPosition[2] =
1422  (newPosition[2] + widIntersect0.z + widIntersect1.z - 2. * m_zPosOffset) / 3.;
1423 
1424  pairOut = pair;
1425  pairOut.setWireID(wireID);
1426  pairOut.setPosition(newPosition);
1427 
1432 
1435 
1436  result = true;
1437  }
1438  }
1439  }
1440 
1441  return result;
1442  }
1443 
1445  const Hit2DSet& hit2DSet,
1446  const reco::ClusterHit3D& pair,
1447  float pairDeltaTimeLimits) const
1448  {
1449  static const float minCharge(0.);
1450 
1451  const reco::ClusterHit2D* bestVHit(0);
1452 
1453  float pairAvePeakTime(pair.getAvePeakTime());
1454 
1455  // Idea is to loop through the input set of hits and look for the best combination
1456  for (const auto& hit2D : hit2DSet) {
1457  if (hit2D->getHit()->Integral() < minCharge) continue;
1458 
1459  float hitVPeakTime(hit2D->getTimeTicks());
1460  float deltaPeakTime(pairAvePeakTime - hitVPeakTime);
1461 
1462  if (deltaPeakTime > pairDeltaTimeLimits) continue;
1463 
1464  if (deltaPeakTime < -pairDeltaTimeLimits) break;
1465 
1466  pairDeltaTimeLimits = fabs(deltaPeakTime);
1467  bestVHit = hit2D;
1468  }
1469 
1470  return bestVHit;
1471  }
1472 
1474  const reco::ClusterHit3D& pair,
1475  float range) const
1476  {
1477  static const float minCharge(0.);
1478 
1479  int numberInRange(0);
1480  float pairAvePeakTime(pair.getAvePeakTime());
1481 
1482  // Idea is to loop through the input set of hits and look for the best combination
1483  for (const auto& hit2D : hit2DSet) {
1484  if (hit2D->getHit()->Integral() < minCharge) continue;
1485 
1486  float hitVPeakTime(hit2D->getTimeTicks());
1487  float deltaPeakTime(pairAvePeakTime - hitVPeakTime);
1488 
1489  if (deltaPeakTime > range) continue;
1490 
1491  if (deltaPeakTime < -range) break;
1492 
1493  numberInRange++;
1494  }
1495 
1496  return numberInRange;
1497  }
1498 
1499  geo::WireID SnippetHit3DBuilder::NearestWireID(const Eigen::Vector3f& position,
1500  const geo::WireID& wireIDIn) const
1501  {
1502  geo::WireID wireID = wireIDIn;
1503 
1504  // Embed the call to the geometry's services nearest wire id method in a try-catch block
1505  try {
1506  // Switch from NearestWireID to this method to avoid the roundoff error issues...
1507  double distanceToWire =
1508  m_geometry->Plane(wireIDIn).WireCoordinate(geo::vect::toPoint(position.data()));
1509 
1510  wireID.Wire = int(distanceToWire);
1511  }
1512  catch (std::exception& exc) {
1513  // This can happen, almost always because the coordinates are **just** out of range
1514  mf::LogWarning("Cluster3D") << "Exception caught finding nearest wire, position - "
1515  << exc.what() << std::endl;
1516 
1517  // Assume extremum for wire number depending on z coordinate
1518  if (position[2] < 0.5 * m_geometry->DetLength())
1519  wireID.Wire = 0;
1520  else
1521  wireID.Wire = m_geometry->Nwires(wireIDIn.asPlaneID()) - 1;
1522  }
1523 
1524  return wireID;
1525  }
1526 
1527  float SnippetHit3DBuilder::DistanceFromPointToHitWire(const Eigen::Vector3f& position,
1528  const geo::WireID& wireIDIn) const
1529  {
1530  float distance = std::numeric_limits<float>::max();
1531 
1532  // Embed the call to the geometry's services nearest wire id method in a try-catch block
1533  try {
1534  // Recover wire geometry information for each wire
1535  const geo::WireGeo& wireGeo = m_geometry->WireIDToWireGeo(wireIDIn);
1536 
1537  // Get wire position and direction for first wire
1538  auto const wirePosArr = wireGeo.GetCenter();
1539 
1540  Eigen::Vector3f wirePos(wirePosArr.X(), wirePosArr.Y(), wirePosArr.Z());
1541  Eigen::Vector3f wireDir(
1542  wireGeo.Direction().X(), wireGeo.Direction().Y(), wireGeo.Direction().Z());
1543 
1544  //*********************************
1545  // Kludge
1546  // if (wireIDIn.Plane > 0) wireDir[2] = -wireDir[2];
1547 
1548  // Want the hit position to have same x value as wire coordinates
1549  Eigen::Vector3f hitPosition(wirePos[0], position[1], position[2]);
1550 
1551  // Get arc length to doca
1552  double arcLen = (hitPosition - wirePos).dot(wireDir);
1553 
1554  // Make sure arclen is in range
1555  if (abs(arcLen) < wireGeo.HalfL()) {
1556  Eigen::Vector3f docaVec = hitPosition - (wirePos + arcLen * wireDir);
1557 
1558  distance = docaVec.norm();
1559  }
1560  }
1561  catch (std::exception& exc) {
1562  // This can happen, almost always because the coordinates are **just** out of range
1563  mf::LogWarning("Cluster3D") << "Exception caught finding nearest wire, position - "
1564  << exc.what() << std::endl;
1565 
1566  // Assume extremum for wire number depending on z coordinate
1567  distance = 0.;
1568  }
1569 
1570  return distance;
1571  }
1572 
1573  //------------------------------------------------------------------------------------------------------------------------------------------
1575  {
1576  // Sort by "modified start time" of pulse
1577  return left->getHit()->PeakTime() < right->getHit()->PeakTime();
1578  }
1579 
1581  const reco::ClusterHit2D* right) const
1582  {
1583  return left->getHit()->PeakTime() < right->getHit()->PeakTime();
1584  }
1585 
1586  //------------------------------------------------------------------------------------------------------------------------------------------
1588  {
1593  // Start by getting a vector of valid, non empty hit collections to make sure we really have something to do here...
1594  // Here is a container for the hits...
1595  std::vector<const recob::Hit*> recobHitVec;
1596 
1597  // Loop through the list of input sources
1598  for (const auto& inputTag : m_hitFinderTagVec) {
1599  art::Handle<std::vector<recob::Hit>> recobHitHandle;
1600  evt.getByLabel(inputTag, recobHitHandle);
1601 
1602  if (!recobHitHandle.isValid() || recobHitHandle->size() == 0) continue;
1603 
1604  recobHitVec.reserve(recobHitVec.size() + recobHitHandle->size());
1605 
1606  for (const auto& hit : *recobHitHandle)
1607  recobHitVec.push_back(&hit);
1608  }
1609 
1610  // If the vector is empty there is nothing to do
1611  if (recobHitVec.empty()) return;
1612 
1613  cet::cpu_timer theClockMakeHits;
1614 
1615  if (m_enableMonitoring) theClockMakeHits.start();
1616 
1617  // We'll want to correct the hit times for the plane offsets
1618  // (note this is already taken care of when converting to position)
1619  std::map<geo::PlaneID, double> planeOffsetMap;
1620 
1621  // Need the detector properties which needs the clocks
1622  auto const clock_data =
1624  auto const det_prop =
1626 
1627  // Try to output a formatted string
1628  std::string debugMessage("");
1629 
1630  // Initialize the plane to hit vector map
1631  for (size_t cryoIdx = 0; cryoIdx < m_geometry->Ncryostats(); cryoIdx++) {
1632  for (size_t tpcIdx = 0; tpcIdx < m_geometry->NTPC(); tpcIdx++) {
1633  m_planeToSnippetHitMap[geo::PlaneID(cryoIdx, tpcIdx, 0)] = SnippetHitMap();
1634  m_planeToSnippetHitMap[geo::PlaneID(cryoIdx, tpcIdx, 1)] = SnippetHitMap();
1635  m_planeToSnippetHitMap[geo::PlaneID(cryoIdx, tpcIdx, 2)] = SnippetHitMap();
1636 
1637  // What we want here are the relative offsets between the planes
1638  // Note that plane 0 is assumed the "first" plane and is the reference
1639  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 0)] = 0.;
1640  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 1)] =
1641  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 1)) -
1642  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0));
1643  planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 2)] =
1644  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 2)) -
1645  det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0));
1646 
1647  // Should we provide output?
1649  std::ostringstream outputString;
1650 
1651  outputString << "***> plane 0 offset: "
1652  << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 0)]
1653  << ", plane 1: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 1)]
1654  << ", plane 2: " << planeOffsetMap[geo::PlaneID(cryoIdx, tpcIdx, 2)]
1655  << "\n";
1656  debugMessage += outputString.str();
1657  outputString << " Det prop plane 0: "
1658  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 0))
1659  << ", plane 1: "
1660  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 1))
1661  << ", plane 2: "
1662  << det_prop.GetXTicksOffset(geo::PlaneID(cryoIdx, tpcIdx, 2))
1663  << ", Trig: " << trigger_offset(clock_data) << "\n";
1664  debugMessage += outputString.str();
1665  }
1666  }
1667  }
1668 
1670  mf::LogDebug("Cluster3D") << debugMessage << std::endl;
1671 
1672  std::cout << debugMessage << std::endl;
1673 
1675  }
1676 
1677  // Cycle through the recob hits to build ClusterHit2D objects and insert
1678  // them into the map
1679  for (const auto& recobHit : recobHitVec) {
1680  // Reject hits with negative charge, these are misreconstructed
1681  if (recobHit->Integral() < 0.) continue;
1682 
1683  // For some detectors we can have multiple wire ID's associated to a given channel.
1684  // So we recover the list of these wire IDs
1685  const std::vector<geo::WireID>& wireIDs = m_geometry->ChannelToWire(recobHit->Channel());
1686 
1687  // Start/End ticks to identify the snippet
1688  HitStartEndPair hitStartEndPair(recobHit->StartTick(), recobHit->EndTick());
1689 
1690  // Can this really happen?
1691  if (hitStartEndPair.second <= hitStartEndPair.first) {
1692  std::cout << "Yes, found a hit with end time less than start time: "
1693  << hitStartEndPair.first << "/" << hitStartEndPair.second
1694  << ", mult: " << recobHit->Multiplicity() << std::endl;
1695  continue;
1696  }
1697 
1698  // And then loop over all possible to build out our maps
1699  //for(const auto& wireID : wireIDs)
1700  for (auto wireID : wireIDs) {
1701  // Check if this is an invalid TPC
1702  // (for example, in protoDUNE there are logical TPC's which see no signal)
1703  if (std::find(m_invalidTPCVec.begin(), m_invalidTPCVec.end(), wireID.TPC) !=
1704  m_invalidTPCVec.end())
1705  continue;
1706 
1707  // Note that a plane ID will define cryostat, TPC and plane
1708  const geo::PlaneID& planeID = wireID.planeID();
1709 
1710  double hitPeakTime(recobHit->PeakTime() - planeOffsetMap[planeID]);
1711  double xPosition(det_prop.ConvertTicksToX(
1712  recobHit->PeakTime(), planeID.Plane, planeID.TPC, planeID.Cryostat));
1713 
1714  m_clusterHit2DMasterList.emplace_back(0, 0., 0., xPosition, hitPeakTime, wireID, recobHit);
1715 
1716  m_planeToSnippetHitMap[planeID][hitStartEndPair].emplace_back(
1717  &m_clusterHit2DMasterList.back());
1718  m_planeToWireToHitSetMap[planeID][wireID.Wire].insert(&m_clusterHit2DMasterList.back());
1719  }
1720  }
1721 
1722  // Make a loop through to sort the recover hits in time order
1723  // for(auto& hitVectorMap : m_planeToSnippetHitMap)
1724  // std::sort(hitVectorMap.second.begin(), hitVectorMap.second.end(), SetHitTimeOrder);
1725 
1726  if (m_enableMonitoring) {
1727  theClockMakeHits.stop();
1728 
1729  m_timeVector[COLLECTARTHITS] = theClockMakeHits.accumulated_real_time();
1730  }
1731 
1732  mf::LogDebug("Cluster3D") << ">>>>> Number of ART hits: " << m_clusterHit2DMasterList.size()
1733  << std::endl;
1734  }
1735 
1736  //------------------------------------------------------------------------------------------------------------------------------------------
1737 
1739  reco::HitPairList& hitPairList,
1740  std::vector<recob::Hit>& hitPtrVec,
1741  RecobHitToPtrMap& recobHitToPtrMap)
1742  {
1743  // Set up the timing
1744  cet::cpu_timer theClockBuildNewHits;
1745 
1746  if (m_enableMonitoring) theClockBuildNewHits.start();
1747 
1748  // We want to build a unique list of hits from the input 3D points which, we know, reuse 2D points frequently
1749  // At the same time we need to create a new 2D hit with the "correct" WireID and replace the old 2D hit in the
1750  // ClusterHit2D object with this new one... while keeping track of the use of the old ones. My head is spinning...
1751  // Declare a set which will allow us to keep track of those CusterHit2D objects we have seen already
1752  std::set<const reco::ClusterHit2D*> visitedHit2DSet;
1753 
1754  // Use this handy art utility to make art::Ptr objects to the new recob::Hits for use in the output phase
1755  art::PtrMaker<recob::Hit> ptrMaker(event);
1756 
1757  // Reserve enough memory to replace every recob::Hit we have considered (this is upper limit)
1758  hitPtrVec.reserve(m_clusterHit2DMasterList.size());
1759 
1760  // Scheme is to loop through all 3D hits, then through each associated ClusterHit2D object
1761  for (reco::ClusterHit3D& hit3D : hitPairList) {
1762  reco::ClusterHit2DVec& hit2DVec = hit3D.getHits();
1763 
1764  // The loop is over the index so we can recover the correct WireID to associate to the new hit when made
1765  for (size_t idx = 0; idx < hit3D.getHits().size(); idx++) {
1766  const reco::ClusterHit2D* hit2D = hit2DVec[idx];
1767 
1768  // Have we seen this 2D hit already?
1769  if (visitedHit2DSet.find(hit2D) == visitedHit2DSet.end()) {
1770  visitedHit2DSet.insert(hit2D);
1771 
1772  // Create and save the new recob::Hit with the correct WireID
1773  hitPtrVec.emplace_back(
1774  recob::HitCreator(*hit2D->getHit(), hit3D.getWireIDs()[idx]).copy());
1775 
1776  // Recover a pointer to it...
1777  recob::Hit* newHit = &hitPtrVec.back();
1778 
1779  // Create a mapping from this hit to an art Ptr representing it
1780  recobHitToPtrMap[newHit] = ptrMaker(hitPtrVec.size() - 1);
1781 
1782  // And set the pointer to this hit in the ClusterHit2D object
1783  const_cast<reco::ClusterHit2D*>(hit2D)->setHit(newHit);
1784  }
1785  }
1786  }
1787 
1788  size_t numNewHits = hitPtrVec.size();
1789 
1790  if (m_enableMonitoring) {
1791  theClockBuildNewHits.stop();
1792 
1793  m_timeVector[BUILDNEWHITS] = theClockBuildNewHits.accumulated_real_time();
1794  }
1795 
1796  mf::LogDebug("Cluster3D") << ">>>>> New output recob::Hit size: " << numNewHits << " (vs "
1797  << m_clusterHit2DMasterList.size() << " input)" << std::endl;
1798 
1799  return;
1800  }
1801 
1804  RecobHitToPtrMap& recobHitPtrMap) const
1805  {
1806  // Let's make sure the input associations container is empty
1808 
1809  // First task is to recover all of the previous wire <--> hit associations and map them by channel number
1810  // Create the temporary container
1811  std::unordered_map<raw::ChannelID_t, art::Ptr<recob::Wire>> channelToWireMap;
1812 
1813  // Go through the list of input sources and fill out the map
1814  for (const auto& inputTag : m_hitFinderTagVec) {
1816  evt.getValidHandle<std::vector<recob::Hit>>(inputTag);
1817 
1818  art::FindOneP<recob::Wire> hitToWireAssns(hitHandle, evt, inputTag);
1819 
1820  if (hitToWireAssns.isValid()) {
1821  for (size_t wireIdx = 0; wireIdx < hitToWireAssns.size(); wireIdx++) {
1822  art::Ptr<recob::Wire> wire = hitToWireAssns.at(wireIdx);
1823 
1824  channelToWireMap[wire->Channel()] = wire;
1825  }
1826  }
1827  }
1828 
1829  // Now fill the container
1830  for (const auto& hitPtrPair : recobHitPtrMap) {
1831  raw::ChannelID_t channel = hitPtrPair.first->Channel();
1832 
1833  std::unordered_map<raw::ChannelID_t, art::Ptr<recob::Wire>>::iterator chanWireItr =
1834  channelToWireMap.find(channel);
1835 
1836  if (!(chanWireItr != channelToWireMap.end())) {
1837  //mf::LogDebug("Cluster3D") << "** Did not find channel to wire match! Skipping..." << std::endl;
1838  continue;
1839  }
1840 
1841  wireAssns.addSingle(chanWireItr->second, hitPtrPair.second);
1842  }
1843 
1844  return;
1845  }
1846 
1849  RecobHitToPtrMap& recobHitPtrMap) const
1850  {
1851  // Let's make sure the input associations container is empty
1852  rawDigitAssns = art::Assns<raw::RawDigit, recob::Hit>();
1853 
1854  // First task is to recover all of the previous wire <--> hit associations and map them by channel number
1855  // Create the temporary container
1856  std::unordered_map<raw::ChannelID_t, art::Ptr<raw::RawDigit>> channelToRawDigitMap;
1857 
1858  // Go through the list of input sources and fill out the map
1859  for (const auto& inputTag : m_hitFinderTagVec) {
1861  evt.getValidHandle<std::vector<recob::Hit>>(inputTag);
1862 
1863  art::FindOneP<raw::RawDigit> hitToRawDigitAssns(hitHandle, evt, inputTag);
1864 
1865  if (hitToRawDigitAssns.isValid()) {
1866  for (size_t rawDigitIdx = 0; rawDigitIdx < hitToRawDigitAssns.size(); rawDigitIdx++) {
1867  art::Ptr<raw::RawDigit> rawDigit = hitToRawDigitAssns.at(rawDigitIdx);
1868 
1869  channelToRawDigitMap[rawDigit->Channel()] = rawDigit;
1870  }
1871  }
1872  }
1873 
1874  // Now fill the container
1875  for (const auto& hitPtrPair : recobHitPtrMap) {
1876  raw::ChannelID_t channel = hitPtrPair.first->Channel();
1877 
1878  std::unordered_map<raw::ChannelID_t, art::Ptr<raw::RawDigit>>::iterator chanRawDigitItr =
1879  channelToRawDigitMap.find(channel);
1880 
1881  if (chanRawDigitItr == channelToRawDigitMap.end()) {
1882  //mf::LogDebug("Cluster3D") << "** Did not find channel to wire match! Skipping..." << std::endl;
1883  continue;
1884  }
1885 
1886  rawDigitAssns.addSingle(chanRawDigitItr->second, hitPtrPair.second);
1887  }
1888 
1889  return;
1890  }
1891 
1892  //------------------------------------------------------------------------------------------------------------------------------------------
1893 
1895 } // namespace lar_cluster3d
TTree * m_tupleTree
output analysis tree
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
Geometry description of a TPC wireThe wire is a single straight segment on a wire plane...
Definition: WireGeo.h:114
void clear()
clear the tuple vectors before processing next event
std::map< geo::TPCID, PlaneToSnippetHitMap > TPCToPlaneToSnippetHitMap
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
std::pair< raw::TDCtick_t, raw::TDCtick_t > HitStartEndPair
Point_t const & GetCenter() const
Returns the world coordinate of the center of the wire [cm].
Definition: WireGeo.h:221
double z
z position of intersection
Definition: geo_types.h:792
std::list< reco::ClusterHit3D > HitPairList
Definition: Cluster3D.h:330
WireGeo const & WireIDToWireGeo(WireID const &wireid) const
Returns the specified wire.
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
float getTimeTicks() const
Definition: Cluster3D.h:75
SnippetHit3DBuilder class definiton.
std::vector< WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
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.
std::vector< ChannelStatusVec > ChannelStatusByPlaneVec
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 BuildChannelStatusVec() const
Create the internal channel status vector (assume will eventually be event-by-event) ...
T * get() const
Definition: ServiceHandle.h:69
void CreateNewRecobHitCollection(art::Event &, reco::HitPairList &, std::vector< recob::Hit > &, RecobHitToPtrMap &)
Create a new 2D hit collection from hits associated to 3D space points.
const Eigen::Vector3f getPosition() const
Definition: Cluster3D.h:155
bool SetHitTimeOrder(const reco::ClusterHit2D *left, const reco::ClusterHit2D *right)
const lariov::ChannelStatusProvider * m_channelFilter
Float_t den
Definition: plot.C:35
float RMS() const
RMS of the hit shape, in tick units.
Definition: Hit.h:228
What follows are several highly useful typedefs which we want to expose to the outside world...
Declaration of signal hit object.
std::vector< size_t > ChannelStatusVec
define data structure for keeping track of channel status
float chargeIntegral(float, float, float, int, int) const
Perform charge integration between limits.
std::map< geo::PlaneID, SnippetHitMap > PlaneToSnippetHitMap
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
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.
std::map< geo::WireID, HitMatchTripletVec > HitMatchTripletVecMap
bool makeHitTriplet(reco::ClusterHit3D &pairOut, const reco::ClusterHit3D &pairIn, const reco::ClusterHit2D *hit2) const
Make a 3D HitPair object by checking two hits.
constexpr auto abs(T v)
Returns the absolute value of the argument.
void BuildHit3D(reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
size_t BuildHitPairMapByTPC(PlaneSnippetHitMapItrPairVec &planeSnippetHitMapItrPairVec, reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
int DegreesOfFreedom() const
Initial tdc tick for hit.
Definition: Hit.h:264
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:211
SnippetHit3DBuilder(fhicl::ParameterSet const &pset)
Constructor.
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
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...
std::map< geo::PlaneID, WireToHitSetMap > PlaneToWireToHitSetMap
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
Definition: GeometryCore.h:430
float m_wirePitchScaleFactor
Scaling factor to determine max distance allowed between candidate pairs.
Length_t DetLength(TPCID const &tpcid=tpc_zero) const
Returns the length of the active volume of the specified TPC.
float m_maxHit3DChiSquare
Provide ability to select hits based on "chi square".
bool SetPeakHitPairIteratorOrder(const reco::HitPairList::iterator &left, const reco::HitPairList::iterator &right)
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.
bool isValid() const noexcept
Definition: Handle.h:203
PutHandle< PROD > put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: Event.h:77
void findGoodTriplets(HitMatchTripletVecMap &, HitMatchTripletVecMap &, reco::HitPairList &) const
This algorithm takes lists of hit pairs and finds good triplets.
std::tuple< const reco::ClusterHit2D *, const reco::ClusterHit2D *, const reco::ClusterHit3D > HitMatchTriplet
This builds a list of candidate hit pairs from lists of hits on two planes.
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
int findGoodHitPairs(SnippetHitMap::iterator &, SnippetHitMap::iterator &, SnippetHitMap::iterator &, HitMatchTripletVecMap &) const
Helper functions to create a hit.
PlaneGeo const & Plane(PlaneID const &planeid) const
Returns the specified wire.
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.
Vector_t Direction() const
Definition: WireGeo.h:289
raw::ChannelID_t Channel() const
Returns the ID of the channel (or InvalidChannelID)
Definition: Wire.h:223
T get(std::string const &key) const
Definition: ParameterSet.h:314
virtual void configure(const fhicl::ParameterSet &) override
Interface for configuring the particular algorithm tool.
std::vector< HitMatchTriplet > HitMatchTripletVec
void makeWireAssns(const art::Event &, art::Assns< recob::Wire, recob::Hit > &, RecobHitToPtrMap &) const
Create recob::Wire to recob::Hit associations.
bool operator()(const SnippetHitMapItrPair &left, const SnippetHitMapItrPair &right) const
Float_t d
Definition: plot.C:235
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
void makeRawDigitAssns(const art::Event &, art::Assns< raw::RawDigit, recob::Hit > &, RecobHitToPtrMap &) const
Create raw::RawDigit to recob::Hit associations.
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.
IHit3DBuilder interface class definiton.
Definition: IHit3DBuilder.h:37
std::map< HitStartEndPair, HitVector > SnippetHitMap
Detector simulation of raw signals on wires.
std::vector< art::InputTag > m_hitFinderTagVec
Data members to follow.
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
std::map< geo::TPCID, PlaneToWireToHitSetMap > TPCToPlaneToWireToHitSetMap
std::vector< const reco::ClusterHit2D * > HitVector
double HalfL() const
Returns half the length of the wire [cm].
Definition: WireGeo.h:172
Filters for channels, events, etc.
size_t BuildHitPairMap(PlaneToSnippetHitMap &planeToHitVectorMap, reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
std::vector< SnippetHitMapItrPair > PlaneSnippetHitMapItrPairVec
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...
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)
float PeakTime() const
Time of the signal peak, in tick units.
Definition: Hit.h:220
ValidHandle< PROD > getValidHandle(InputTag const &tag) const
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.
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
virtual void produces(art::ProducesCollector &) override
Each algorithm may have different objects it wants "produced" so use this to let the top level produc...
bool WireIDsIntersect(const geo::WireID &, const geo::WireID &, geo::WireIDIntersection &) const
function to detemine if two wires "intersect" (in the 2D sense)
void CollectArtHits(const art::Event &evt) const
Extract the ART hits and the ART hit-particle relationships.
double y
y position of intersection
Definition: geo_types.h:791
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
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
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
virtual void Hit3DBuilder(art::Event &, reco::HitPairList &, RecobHitToPtrMap &) override
Given a set of recob hits, run DBscan to form 3D clusters.
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
Float_t e
Definition: plot.C:35
second_as<> second
Type of time stored in seconds, in double precision.
Definition: spacetime.h:82
float closestApproach(const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, const Eigen::Vector3f &, float &, float &) const
function to compute the distance of closest approach and the arc length to the points of closest appr...
bool m_outputHistograms
Take the time to create and fill some histograms for diagnostics.
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:168
std::pair< SnippetHitMap::iterator, SnippetHitMap::iterator > SnippetHitMapItrPair
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
virtual float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
If monitoring, recover the time to execute a particular function.
bool SetPairStartTimeOrder(const reco::ClusterHit3D &left, const reco::ClusterHit3D &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.
art framework interface to geometry description
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
vec_iX clear()
Event finding and building.
std::map< unsigned int, Hit2DSet > WireToHitSetMap
unsigned getStatusBits() const
Definition: Cluster3D.h:71
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:79