LArSoft  v06_85_00
Liquid Argon Software toolkit - http://larsoft.org/
StandardHit3DBuilder_tool.cc
Go to the documentation of this file.
1 
8 // Framework Includes
10 #include "cetlib/search_path.h"
11 #include "cetlib/cpu_timer.h"
13 
15 
16 // LArSoft includes
21 #include "larevt/CalibrationDBI/Interface/ChannelStatusService.h"
22 #include "larevt/CalibrationDBI/Interface/ChannelStatusProvider.h"
23 
24 // std includes
25 #include <string>
26 #include <functional>
27 #include <iostream>
28 #include <memory>
29 
30 //------------------------------------------------------------------------------------------------------------------------------------------
31 // implementation follows
32 
33 namespace lar_cluster3d {
34 
40 // forward declaration to define an ordering function for our hit set
42 {
43  bool operator() (const reco::ClusterHit2D*, const reco::ClusterHit2D*) const;
44 };
45 
46 using HitVector = std::vector<reco::ClusterHit2D*>;
47 using PlaneToHitVectorMap = std::map<geo::PlaneID, HitVector>;
48 using TPCToPlaneToHitVectorMap = std::map<geo::TPCID, PlaneToHitVectorMap>;
49 using Hit2DVector = std::vector<reco::ClusterHit2D>;
50 using Hit2DSet = std::set<const reco::ClusterHit2D*, Hit2DSetCompare>;
51 using WireToHitSetMap = std::map<unsigned int, Hit2DSet>;
52 using PlaneToWireToHitSetMap = std::map<geo::PlaneID, WireToHitSetMap>;
53 using TPCToPlaneToWireToHitSetMap = std::map<geo::TPCID, PlaneToWireToHitSetMap>;
54 using HitVectorMap = std::map<size_t, HitVector>;
55 
56 using HitPairVector = std::vector<std::unique_ptr<reco::ClusterHit3D>>;
57 
61 class StandardHit3DBuilder : virtual public IHit3DBuilder
62 {
63 public:
69  explicit StandardHit3DBuilder(fhicl::ParameterSet const &pset);
70 
75 
76  void configure(const fhicl::ParameterSet&) override;
77 
84  void Hit3DBuilder(const art::Event &evt, reco::HitPairList& hitPairList, RecobHitToPtrMap&) const override;
85 
89  float getTimeToExecute(IHit3DBuilder::TimeValues index) const override {return m_timeVector.at(index);}
90 
91 private:
92 
102  void CollectArtHits(const art::Event& evt,
103  RecobHitToPtrMap& hitToPtrMap) const;
104 
108  void BuildHit3D(reco::HitPairList& hitPairList) const;
109 
113  size_t BuildHitPairMap(PlaneToHitVectorMap& planeToHitVectorMap, reco::HitPairList& hitPairList) const;
114 
118  using PlaneHitVectorItrPairVec = std::vector<std::pair<HitVector::iterator,HitVector::iterator>>;
119 
120  size_t BuildHitPairMapByTPC(PlaneHitVectorItrPairVec& planeHitVectorItrPairVec, reco::HitPairList& hitPairList) const;
121 
125  using HitMatchPair = std::pair<const reco::ClusterHit2D*,reco::ClusterHit3D>;
126  using HitMatchPairVec = std::vector<HitMatchPair>;
127  using HitMatchPairVecMap = std::map<geo::WireID,HitMatchPairVec>;
128 
129  int findGoodHitPairs(const reco::ClusterHit2D*, HitVector::iterator&, HitVector::iterator&, HitMatchPairVecMap&) const;
130 
134  void findGoodTriplets(HitMatchPairVecMap&, HitMatchPairVecMap&, reco::HitPairList&, bool = false) const;
135 
139  bool makeHitPair(reco::ClusterHit3D& pairOut,
140  const reco::ClusterHit2D* hit1,
141  const reco::ClusterHit2D* hit2,
142  float hitWidthSclFctr = 1.,
143  size_t hitPairCntr = 0) const;
144 
148  bool makeHitTriplet(reco::ClusterHit3D& pairOut,
149  const reco::ClusterHit3D& pairIn,
150  const reco::ClusterHit2D* hit2) const;
151 
155  bool makeDeadChannelPair(reco::ClusterHit3D& pairOut, const reco::ClusterHit3D& pair, size_t maxStatus = 4, size_t minStatus = 0, float minOverlap=0.2) const;
156 
160  const reco::ClusterHit2D* FindBestMatchingHit(const Hit2DSet& hit2DSet, const reco::ClusterHit3D& pair, float pairDeltaTimeLimits) const;
161 
165  int FindNumberInRange(const Hit2DSet& hit2DSet, const reco::ClusterHit3D& pair, float range) const;
166 
170  geo::WireID NearestWireID(const float* position, const geo::WireID& wireID) const;
171 
175  void BuildChannelStatusVec(PlaneToWireToHitSetMap& planeToWiretoHitSetMap) const;
176 
180  using ChannelStatusVec = std::vector<size_t>;
181  using ChannelStatusByPlaneVec = std::vector<ChannelStatusVec>;
182 
190 
192  float m_wirePitch[3];
193  mutable std::vector<float> m_timeVector;
194 
196 
197  // Get instances of the primary data structures needed
201 
202 
204  mutable size_t m_numBadChannels;
205 
206  geo::Geometry* m_geometry; //< pointer to the Geometry service
207  const detinfo::DetectorProperties* m_detector; //< Pointer to the detector properties
208  const lariov::ChannelStatusProvider* m_channelFilter;
209 };
210 
212  m_channelFilter(&art::ServiceHandle<lariov::ChannelStatusService>()->GetProvider())
213 
214 {
215  this->configure(pset);
216 }
217 
218 //------------------------------------------------------------------------------------------------------------------------------------------
219 
221 {
222 }
223 
224 //------------------------------------------------------------------------------------------------------------------------------------------
225 
227 {
228  m_hitFinderTag = pset.get<art::InputTag>("HitFinderTag");
229  m_enableMonitoring = pset.get<bool> ("EnableMonitoring", true);
230  m_numSigmaPeakTime = pset.get<float> ("NumSigmaPeakTime", 3. );
231  m_hitWidthSclFctr = pset.get<float> ("HitWidthScaleFactor", 6. );
232  m_deltaPeakTimeSig = pset.get<float> ("DeltaPeakTimeSig", 1.7 );
233  m_zPosOffset = pset.get<float> ("ZPosOffset", 0.0 );
234 
236 
237  m_geometry = &*geometry;
238  m_detector = lar::providerFrom<detinfo::DetectorPropertiesService>();
239 
241  m_wirePitch[1] = m_geometry->WirePitch(1);
242  m_wirePitch[2] = m_geometry->WirePitch(2);
243 }
244 
246 {
247  // This is called each event, clear out the previous version and start over
248  if (!m_channelStatus.empty()) m_channelStatus.clear();
249 
250  m_numBadChannels = 0;
252 
253  // Loop through views/planes to set the wire length vectors
254  for(size_t idx = 0; idx < m_channelStatus.size(); idx++)
255  {
257  }
258 
259  // Loop through the channels and mark those that are "bad"
260  for(size_t channel = 0; channel < m_geometry->Nchannels(); channel++)
261  {
262  if( !m_channelFilter->IsGood(channel))
263  {
264  std::vector<geo::WireID> wireIDVec = m_geometry->ChannelToWire(channel);
265  geo::WireID wireID = wireIDVec[0];
266  lariov::ChannelStatusProvider::Status_t chanStat = m_channelFilter->Status(channel);
267 
268  m_channelStatus[wireID.Plane][wireID.Wire] = chanStat;
270  }
271  }
272 
273  // add quiet wires in U plane for microboone (this will done "correctly" in near term)
274  // PlaneToWireToHitSetMap::iterator plane0HitItr = planeToWireToHitSetMap.find(geo::PlaneID(0,0,0));
275 
276  // if (plane0HitItr != planeToWireToHitSetMap.end())
277  // {
279 
280  // for(size_t idx = 2016; idx < 2096; idx++) m_channelStatus[0][idx] = 3;
281  // for(size_t idx = 2192; idx < 2304; idx++) m_channelStatus[0][idx] = 3;
282  // for(size_t idx = 2352; idx < 2384; idx++) m_channelStatus[0][idx] = 3;
283  // //for(size_t idx = 2016; idx < 2096; idx++) if (wireToHitSetMap.find(idx) == wireToHitSetMap.end()) m_channelStatus[0][idx] = 3;
284  // //for(size_t idx = 2192; idx < 2304; idx++) if (wireToHitSetMap.find(idx) == wireToHitSetMap.end()) m_channelStatus[0][idx] = 3;
285  // //for(size_t idx = 2352; idx < 2384; idx++) if (wireToHitSetMap.find(idx) == wireToHitSetMap.end()) m_channelStatus[0][idx] = 3;
287  // }
288 
289  return;
290 }
291 
292 
294 {
295  return (*left)->getAvePeakTime() < (*right)->getAvePeakTime();
296 }
297 
298 struct HitPairClusterOrder
299 {
301  {
302  // Watch out for the case where two clusters can have the same number of hits!
303  if (left->second.size() == right->second.size())
304  return left->first < right->first;
305 
306  return left->second.size() > right->second.size();
307  }
308 };
309 
310 void StandardHit3DBuilder::Hit3DBuilder(const art::Event& evt, reco::HitPairList& hitPairList, RecobHitToPtrMap& clusterHitToArtPtrMap) const
311 {
312  // Clear the internal data structures
313  m_clusterHit2DMasterVec.clear();
314  m_planeToHitVectorMap.clear();
315  m_planeToWireToHitSetMap.clear();
316 
317  m_timeVector.resize(NUMTIMEVALUES, 0.);
318 
319  // Recover the 2D hits and then organize them into data structures which will be used in the
320  // DBscan algorithm for building the 3D clusters
321  this->CollectArtHits(evt, clusterHitToArtPtrMap);
322 
323 // if (m_enableMonitoring) theClockArtHits.stop();
324 
325  // If there are no hits in our view/wire data structure then do not proceed with the full analysis
326  if (!m_planeToWireToHitSetMap.empty())
327  {
328  // Call the algorithm that builds 3D hits
329  this->BuildHit3D(hitPairList);
330  }
331 
332  return;
333 }
334 
336 {
341  cet::cpu_timer theClockMakeHits;
342 
343  if (m_enableMonitoring) theClockMakeHits.start();
344 
345  // The first task is to take the lists of input 2D hits (a map of view to sorted lists of 2D hits)
346  // and then to build a list of 3D hits to be used in downstream processing
348 
349  size_t numHitPairs = BuildHitPairMap(m_planeToHitVectorMap, hitPairList);
350 
351  if (m_enableMonitoring)
352  {
353  theClockMakeHits.stop();
354 
355  m_timeVector[BUILDTHREEDHITS] = theClockMakeHits.accumulated_real_time();
356  }
357 
358  mf::LogDebug("Cluster3D") << ">>>>> 3D hit building done, found " << numHitPairs << " 3D Hits" << std::endl;
359 
360  return;
361 }
362 
363 //------------------------------------------------------------------------------------------------------------------------------------------
364 namespace {
365 bool SetHitStartTimeOrder(const reco::ClusterHit2D* left, const reco::ClusterHit2D* right)
366 {
367  // Sort by "modified start time" of pulse
368  return left->getHit().PeakTime() - left->getHit().RMS() < right->getHit().PeakTime() - right->getHit().RMS();
369 }
370 
371 using HitVectorItrPair = std::pair<HitVector::iterator,HitVector::iterator>;
372 
373 class SetStartTimeOrder
374 {
375 public:
376  SetStartTimeOrder() : m_numRMS(1.) {}
377  SetStartTimeOrder(float numRMS) : m_numRMS(numRMS) {}
378 
379  bool operator()(const HitVectorItrPair& left, const HitVectorItrPair& right) const
380  {
381  // Protect against possible issue?
382  if (left.first != left.second && right.first != right.second)
383  {
384  // Sort by "modified start time" of pulse
385  return (*left.first)->getTimeTicks() - m_numRMS*(*left.first)->getHit().RMS() < (*right.first)->getTimeTicks() - m_numRMS*(*right.first)->getHit().RMS();
386  }
387 
388  return left.first != left.second;
389  }
390 
391 private:
392  float m_numRMS;
393 };
394 
395 bool SetPairStartTimeOrder(const std::unique_ptr<reco::ClusterHit3D>& left, const std::unique_ptr<reco::ClusterHit3D>& right)
396 {
397  // Sort by "modified start time" of pulse
398  return left->getAvePeakTime() - left->getSigmaPeakTime() < right->getAvePeakTime() - right->getSigmaPeakTime();
399 }
400 }
401 
402 //------------------------------------------------------------------------------------------------------------------------------------------
403 
404 size_t StandardHit3DBuilder::BuildHitPairMap(PlaneToHitVectorMap& planeToHitVectorMap, reco::HitPairList& hitPairList) const
405 {
415  size_t totalNumHits(0);
416  size_t hitPairCntr(0);
417 
418  size_t nTriplets(0);
419  size_t nDeadChanHits(0);
420 
421  // Set up to loop over cryostats and tpcs...
422  for(size_t cryoIdx = 0; cryoIdx < m_geometry->Ncryostats(); cryoIdx++)
423  {
424  for(size_t tpcIdx = 0; tpcIdx < m_geometry->NTPC(); tpcIdx++)
425  {
426  PlaneToHitVectorMap::iterator mapItr0 = planeToHitVectorMap.find(geo::PlaneID(cryoIdx,tpcIdx,0));
427  PlaneToHitVectorMap::iterator mapItr1 = planeToHitVectorMap.find(geo::PlaneID(cryoIdx,tpcIdx,1));
428  PlaneToHitVectorMap::iterator mapItr2 = planeToHitVectorMap.find(geo::PlaneID(cryoIdx,tpcIdx,2));
429 
430  size_t nPlanesWithHits = (mapItr0 != planeToHitVectorMap.end() && !mapItr0->second.empty() ? 1 : 0)
431  + (mapItr1 != planeToHitVectorMap.end() && !mapItr1->second.empty() ? 1 : 0)
432  + (mapItr2 != planeToHitVectorMap.end() && !mapItr2->second.empty() ? 1 : 0);
433 
434  if (nPlanesWithHits < 2) continue;
435 
436  HitVector& hitVector0 = mapItr0->second;
437  HitVector& hitVector1 = mapItr1->second;
438  HitVector& hitVector2 = mapItr2->second;
439 
440  // We are going to resort the hits into "start time" order...
441  std::sort(hitVector0.begin(), hitVector0.end(), SetHitStartTimeOrder);
442  std::sort(hitVector1.begin(), hitVector1.end(), SetHitStartTimeOrder);
443  std::sort(hitVector2.begin(), hitVector2.end(), SetHitStartTimeOrder);
444 
445  PlaneHitVectorItrPairVec hitItrVec = {HitVectorItrPair(hitVector0.begin(),hitVector0.end()),
446  HitVectorItrPair(hitVector1.begin(),hitVector1.end()),
447  HitVectorItrPair(hitVector2.begin(),hitVector2.end())};
448 
449  totalNumHits += BuildHitPairMapByTPC(hitItrVec, hitPairList);
450  }
451  }
452 
453  // Return the hit pair list but sorted by z and y positions (faster traversal in next steps)
454  hitPairList.sort(SetPairStartTimeOrder);
455 
456  // Where are we?
457  mf::LogDebug("Cluster3D") << "Total number hits: " << totalNumHits << std::endl;
458  mf::LogDebug("Cluster3D") << "Created a total of " << hitPairList.size() << " hit pairs, counted: " << hitPairCntr << std::endl;
459  mf::LogDebug("Cluster3D") << "-- Triplets: " << nTriplets << ", dead channel pairs: " << nDeadChanHits << std::endl;
460 
461  return hitPairList.size();
462 }
463 
465 {
476  // Define functions to set start/end iterators in the loop below
477  auto SetStartIterator = [](HitVector::iterator startItr, HitVector::iterator endItr, float rms, float startTime)
478  {
479  while(startItr != endItr)
480  {
481  float numRMS(rms);
482  if ((*startItr)->getTimeTicks() + numRMS * (*startItr)->getHit().RMS() < startTime) startItr++;
483  else break;
484  }
485  return startItr;
486  };
487 
488  auto SetEndIterator = [](HitVector::iterator firstItr, HitVector::iterator endItr, float rms, float endTime)
489  {
490  while(firstItr != endItr)
491  {
492  float numRMS(rms);
493  if ((*firstItr)->getTimeTicks() - numRMS * (*firstItr)->getHit().RMS() < endTime) firstItr++;
494  else break;
495  }
496  return firstItr;
497  };
498 
499  size_t nTriplets(0);
500  size_t nDeadChanHits(0);
501 
502  //*********************************************************************************
503  // Basically, we try to loop until done...
504  while(1)
505  {
506  // Sort so that the earliest hit time will be the first element, etc.
507  std::sort(hitItrVec.begin(),hitItrVec.end(),SetStartTimeOrder(m_numSigmaPeakTime));
508 
509  // This loop iteration's golden hit
510  const reco::ClusterHit2D* goldenHit = *hitItrVec[0].first;
511 
512  // The range of history... (for this hit)
513  float goldenTimeStart = goldenHit->getTimeTicks() - m_numSigmaPeakTime * goldenHit->getHit().RMS() - 0.1;
514  float goldenTimeEnd = goldenHit->getTimeTicks() + m_numSigmaPeakTime * goldenHit->getHit().RMS() + 0.1;
515 
516  // Set iterators to insure we'll be in the overlap ranges
517  HitVector::iterator hitItr1Start = SetStartIterator(hitItrVec[1].first, hitItrVec[1].second, m_numSigmaPeakTime, goldenTimeStart);
518  HitVector::iterator hitItr1End = SetEndIterator( hitItr1Start, hitItrVec[1].second, m_numSigmaPeakTime, goldenTimeEnd);
519  HitVector::iterator hitItr2Start = SetStartIterator(hitItrVec[2].first, hitItrVec[2].second, m_numSigmaPeakTime, goldenTimeStart);
520  HitVector::iterator hitItr2End = SetEndIterator( hitItr2Start, hitItrVec[2].second, m_numSigmaPeakTime, goldenTimeEnd);
521 
522  // Since we'll use these many times in the internal loops, pre make the pairs for the second set of hits
523  size_t curHitListSize(hitPairList.size());
524  HitMatchPairVecMap pair12Map;
525  HitMatchPairVecMap pair13Map;
526 
527  size_t n12Pairs = findGoodHitPairs(goldenHit, hitItr1Start, hitItr1End, pair12Map);
528  size_t n13Pairs = findGoodHitPairs(goldenHit, hitItr2Start, hitItr2End, pair13Map);
529 
530  nDeadChanHits += hitPairList.size() - curHitListSize;
531  curHitListSize = hitPairList.size();
532 
533  if (n12Pairs > n13Pairs) findGoodTriplets(pair12Map, pair13Map, hitPairList);
534  else findGoodTriplets(pair13Map, pair12Map, hitPairList);
535 
536  nTriplets += hitPairList.size() - curHitListSize;
537 
538  hitItrVec[0].first++;
539 
540  int nPlanesWithHits(0);
541 
542  for(auto& pair : hitItrVec)
543  if (pair.first != pair.second) nPlanesWithHits++;
544 
545  if (nPlanesWithHits < 2) break;
546  }
547 
548  return hitPairList.size();
549 }
550 
552  HitVector::iterator& startItr,
553  HitVector::iterator& endItr,
554  HitMatchPairVecMap& hitMatchMap) const
555 {
556  int numPairs(0);
557 
558  // Loop through the input secon hits and make pairs
559  while(startItr != endItr)
560  {
561  reco::ClusterHit2D* hit = *startItr++;
562  reco::ClusterHit3D pair;
563 
564  // pair returned with a negative ave time is signal of failure
565  if (!makeHitPair(pair, goldenHit, hit, m_hitWidthSclFctr)) continue;
566 
567  geo::WireID wireID = hit->getHit().WireID();
568 
569  hitMatchMap[wireID].emplace_back(hit,pair);
570 
571  numPairs++;
572  }
573 
574  return numPairs;
575 }
576 
577 void StandardHit3DBuilder::findGoodTriplets(HitMatchPairVecMap& pair12Map, HitMatchPairVecMap& pair13Map, reco::HitPairList& hitPairList, bool tagged) const
578 {
579  // Build triplets from the two lists of hit pairs
580  if (!pair12Map.empty())
581  {
582  // temporary container for dead channel hits
583  std::vector<reco::ClusterHit3D> tempDeadChanVec;
584  reco::ClusterHit3D deadChanPair;
585 
586  // Keep track of which third plane hits have been used
587  std::map<const reco::ClusterHit3D*,bool> usedPairMap;
588 
589  // Initial population of this map with the pair13Map hits
590  for(const auto& pair13 : pair13Map)
591  {
592  for(const auto& hit2Dhit3DPair : pair13.second) usedPairMap[&hit2Dhit3DPair.second] = false;
593  }
594 
595  // The outer loop is over all hit pairs made from the first two plane combinations
596  for(const auto& pair12 : pair12Map)
597  {
598  if (pair12.second.empty()) continue;
599 
600  // Use the planeID for the first hit
601  geo::WireID missingPlaneID = pair12.first;
602 
603  // "Discover" the missing view (and we can't rely on assuming there are hits in the pair13Map at this point)
604  size_t missPlane = 0;
605 
606  if (!pair12.second.front().second.getHits()[1]) missPlane = 1;
607  else if (!pair12.second.front().second.getHits()[2]) missPlane = 2;
608 
609  missingPlaneID.Plane = missPlane;
610 
611  // This loop is over hit pairs that share the same first two plane wires but may have different
612  // hit times on those wires
613  for(const auto& hit2Dhit3DPair : pair12.second)
614  {
615  const reco::ClusterHit3D& pair1 = hit2Dhit3DPair.second;
616 
617  // Get the wire ID for the nearest wire to the position of this hit
618  geo::WireID wireID = NearestWireID(pair1.getPosition(), missingPlaneID);
619 
620  // populate the map with initial value
621  usedPairMap[&pair1] = false;
622 
623  // Now look up the hit pairs on the wire which matches the current hit pair
624  HitMatchPairVecMap::iterator thirdPlaneHitMapItr = pair13Map.find(wireID);
625 
626  // It may be that the "nearest" wire is actually the next wire (which can happen
627  // because we are right between two wires in this plane and may have some round off
628  // error pointing us a the wrong wire)
629  if (thirdPlaneHitMapItr == pair13Map.end())
630  {
631  wireID.Wire += 1;
632 
633  thirdPlaneHitMapItr = pair13Map.find(wireID);
634  }
635 
636  // Loop over third plane hits and try to form a triplet
637  if (thirdPlaneHitMapItr != pair13Map.end())
638  {
639  for(const auto& thirdPlaneHitItr : thirdPlaneHitMapItr->second)
640  {
641  const reco::ClusterHit2D* hit2 = thirdPlaneHitItr.first;
642  const reco::ClusterHit3D& pair2 = thirdPlaneHitItr.second;
643 
644  // If success try for the triplet
645  reco::ClusterHit3D triplet;
646 
647  if (makeHitTriplet(triplet, pair1, hit2))
648  {
649  triplet.setID(hitPairList.size());
650  hitPairList.emplace_back(new reco::ClusterHit3D(triplet));
651  usedPairMap[&pair1] = true;
652  usedPairMap[&pair2] = true;
653  }
654  }
655  }
656  }
657  }
658 
659  // One more loop through the other pairs to check for sick channels
660  if (m_numBadChannels > 0)
661  {
662  for(const auto& pairMapPair : usedPairMap)
663  {
664  if (pairMapPair.second) continue;
665 
666  const reco::ClusterHit3D* pair = pairMapPair.first;
667 
668  // Here we look to see if we failed to make a triplet because the partner wire was dead/noisy/sick
669  if (makeDeadChannelPair(deadChanPair, *pair, 4, 0, 0.)) tempDeadChanVec.emplace_back(deadChanPair);
670  }
671 
672  // Handle the dead wire triplets
673  if(!tempDeadChanVec.empty())
674  {
675  // If we have many then see if we can trim down a bit by keeping those with time significance
676  if (tempDeadChanVec.size() > 1)
677  {
678  // Sort by "significance" of agreement
679  std::sort(tempDeadChanVec.begin(),tempDeadChanVec.end(),[](const auto& left, const auto& right){return left.getDeltaPeakTime()/left.getSigmaPeakTime() < right.getDeltaPeakTime()/right.getSigmaPeakTime();});
680 
681  // What is the range of "significance" from first to last?
682  float firstSig = tempDeadChanVec.front().getDeltaPeakTime() / tempDeadChanVec.front().getSigmaPeakTime();
683  float lastSig = tempDeadChanVec.back().getDeltaPeakTime() / tempDeadChanVec.back().getSigmaPeakTime();
684  float sigRange = lastSig - firstSig;
685 
686  if (lastSig > 0.5 * m_deltaPeakTimeSig && sigRange > 0.5)
687  {
688  // Declare a maximum of 1.5 * the average of the first and last pairs...
689  float maxSignificance = std::max(0.75 * (firstSig + lastSig),1.0);
690 
691  std::vector<reco::ClusterHit3D>::iterator firstBadElem = std::find_if(tempDeadChanVec.begin(),tempDeadChanVec.end(),[&maxSignificance](const auto& pair){return pair.getDeltaPeakTime()/pair.getSigmaPeakTime() > maxSignificance;});
692 
693  // But only keep the best 10?
694  if (std::distance(tempDeadChanVec.begin(),firstBadElem) > 20) firstBadElem = tempDeadChanVec.begin() + 20;
695  // Keep at least one hit...
696  else if (firstBadElem == tempDeadChanVec.begin()) firstBadElem++;
697 
698  tempDeadChanVec.resize(std::distance(tempDeadChanVec.begin(),firstBadElem));
699  }
700  }
701 
702  for(auto& pair : tempDeadChanVec)
703  {
704  pair.setID(hitPairList.size());
705  hitPairList.emplace_back(new reco::ClusterHit3D(pair));
706  }
707  }
708  }
709  }
710 
711  return;
712 }
713 
715  const reco::ClusterHit2D* hit1,
716  const reco::ClusterHit2D* hit2,
717  float hitWidthSclFctr,
718  size_t hitPairCntr) const
719 {
720  // Assume failure
721  bool result(false);
722 
723  // We assume in this routine that we are looking at hits in different views
724  // The first mission is to check that the wires intersect
725  const geo::WireID& hit1WireID = hit1->getHit().WireID();
726  const geo::WireID& hit2WireID = hit2->getHit().WireID();
727 
728  geo::WireIDIntersection widIntersect;
729 
730  if (m_geometry->WireIDsIntersect(hit1WireID, hit2WireID, widIntersect))
731  {
732  // Wires intersect so now we can check the timing
733  float hit1Peak = hit1->getTimeTicks();
734  float hit1Sigma = hit1->getHit().RMS();
735 
736  float hit2Peak = hit2->getTimeTicks();
737  float hit2Sigma = hit2->getHit().RMS();
738 
739  // ad hoc correction for most bad fits...
740  if (hit1Sigma > 2. * hit1->getHit().PeakAmplitude()) hit1Sigma = 2. * hit1->getHit().PeakAmplitude();
741  if (hit2Sigma > 2. * hit2->getHit().PeakAmplitude()) hit2Sigma = 2. * hit2->getHit().PeakAmplitude();
742 
743  float hit1Width = hitWidthSclFctr * hit1Sigma;
744  float hit2Width = hitWidthSclFctr * hit2Sigma;
745 
746  // Coarse check hit times are "in range"
747  if (fabs(hit1Peak - hit2Peak) <= (hit1Width + hit2Width))
748  {
749  // Check to see that hit peak times are consistent with each other
750  float hit1SigSq = hit1Sigma * hit1Sigma;
751  float hit2SigSq = hit2Sigma * hit2Sigma;
752  float avePeakTime = (hit1Peak / hit1SigSq + hit2Peak / hit2SigSq) * hit1SigSq * hit2SigSq / (hit1SigSq + hit2SigSq);
753  float deltaPeakTime = std::fabs(hit1Peak - hit2Peak);
754  float sigmaPeakTime = std::sqrt(hit1SigSq + hit2SigSq);
755 
756  // delta peak time consistency check here
757  if (deltaPeakTime < m_deltaPeakTimeSig * sigmaPeakTime) // 2 sigma consistency? (do this way to avoid divide)
758  {
759  float totalCharge = hit1->getHit().Integral() + hit2->getHit().Integral();
760  float hitChiSquare = std::pow((hit1Peak - avePeakTime),2) / hit1SigSq
761  + std::pow((hit2Peak - avePeakTime),2) / hit2SigSq;
762 
763  float xPositionHit1(hit1->getXPosition());
764  float xPositionHit2(hit2->getXPosition());
765  float xPosition = (xPositionHit1 / hit1SigSq + xPositionHit2 / hit2SigSq) * hit1SigSq * hit2SigSq / (hit1SigSq + hit2SigSq);
766 
767  float position[] = {xPosition, float(widIntersect.y), float(widIntersect.z)-m_zPosOffset};
768 
769  // If to here then we need to sort out the hit pair code telling what views are used
770  unsigned statusBits = 1 << hit1->getHit().WireID().Plane | 1 << hit2->getHit().WireID().Plane;
771 
772  // handle status bits for the 2D hits
775 
778 
779  reco::ClusterHit2DVec hitVector;
780 
781  hitVector.resize(3, NULL);
782 
783  hitVector.at(hit1->getHit().WireID().Plane) = hit1;
784  hitVector.at(hit2->getHit().WireID().Plane) = hit2;
785 
786  unsigned int cryostatIdx = hit1->getHit().WireID().Cryostat;
787  unsigned int tpcIdx = hit1->getHit().WireID().TPC;
788 
789  // Initialize the wireIdVec
790  std::vector<geo::WireID> wireIDVec = {geo::WireID(cryostatIdx,tpcIdx,0,0),
791  geo::WireID(cryostatIdx,tpcIdx,1,0),
792  geo::WireID(cryostatIdx,tpcIdx,2,0)};
793 
794  wireIDVec[hit1->getHit().WireID().Plane] = hit1->getHit().WireID();
795  wireIDVec[hit2->getHit().WireID().Plane] = hit2->getHit().WireID();
796 
797  // For compiling at the moment
798  std::vector<float> hitDelTSigVec = {0.,0.,0.};
799 
800  hitDelTSigVec.at(hit1->getHit().WireID().Plane) = deltaPeakTime / sigmaPeakTime;
801  hitDelTSigVec.at(hit2->getHit().WireID().Plane) = deltaPeakTime / sigmaPeakTime;
802 
803  // Create the 3D cluster hit
804  hitPair.initialize(hitPairCntr,
805  statusBits,
806  position,
807  totalCharge,
808  avePeakTime,
809  deltaPeakTime,
810  sigmaPeakTime,
811  hitChiSquare,
812  0.,
813  0.,
814  hitVector,
815  hitDelTSigVec,
816  wireIDVec);
817 
818  result = true;
819  }
820  }
821  }
822 
823  // Send it back
824  return result;
825 }
826 
827 
829  const reco::ClusterHit3D& pair,
830  const reco::ClusterHit2D* hit) const
831 {
832  // Assume failure
833  bool result(false);
834 
835  static const float rmsToSig(1.0); //0.75); //0.57735027);
836 
837  // Recover hit info
838  float hitTimeTicks = hit->getTimeTicks();
839  float hitSigma = hit->getHit().RMS();
840 
841  // Special case check
842  if (hitSigma > 2. * hit->getHit().PeakAmplitude()) hitSigma = 2. * hit->getHit().PeakAmplitude();
843 
844  // Let's do a quick consistency check on the input hits to make sure we are in range...
845  // Require the W hit to be "in range" with the UV Pair
846  if (fabs(hitTimeTicks - pair.getAvePeakTime()) < m_hitWidthSclFctr * (pair.getSigmaPeakTime() + hitSigma))
847  {
848  // Timing in range, now check that the input hit wire "intersects" with the input pair's wires
849  geo::WireID wireID = NearestWireID(pair.getPosition(), hit->getHit().WireID());
850 
851  // There is an interesting round off issue that we need to watch for...
852  if (wireID.Wire == hit->getHit().WireID().Wire || wireID.Wire + 1 == hit->getHit().WireID().Wire)
853  {
854  // Use the existing code to see the U and W hits are willing to pair with the V hit
855  reco::ClusterHit3D pair0h;
856  reco::ClusterHit3D pair1h;
857 
858  // Recover all the hits involved
859  const reco::ClusterHit2DVec& pairHitVec = pair.getHits();
860  const reco::ClusterHit2D* hit0 = pairHitVec.at(0);
861  const reco::ClusterHit2D* hit1 = pairHitVec.at(1);
862 
863  if (!hit0) hit0 = pairHitVec.at(2);
864  else if (!hit1) hit1 = pairHitVec.at(2);
865 
866  // If good pairs made here then we can try to make a triplet
867  if (makeHitPair(pair0h, hit0, hit, m_hitWidthSclFctr) && makeHitPair(pair1h, hit1, hit, m_hitWidthSclFctr))
868  {
869  // We want to make sure the 3 sets of pair are really consistent - ie the difference in position in the Y-Z
870  // plane (the wire plane) is less than the wire pitch. Note that we are forming an equalateral triangle so we
871  // really only need compute one side of it to get the answer we're looking for...
872  float deltaPairY = pair1h.getPosition()[1] - pair0h.getPosition()[1];
873  float deltaPairZ = pair1h.getPosition()[2] - pair0h.getPosition()[2];
874  float yzdistance = std::sqrt(deltaPairY * deltaPairY + deltaPairZ * deltaPairZ);
875 
876  // The intersection of wires on 3 planes is actually an equilateral triangle... Each pair will have its position at one of the
877  // corners, the difference in distance along the z axis will be 1/2 wire spacing, the difference along the y axis is
878  // 1/2 wire space / cos(pitch)
879 // if (std::fabs(std::fabs(deltaZ_w) - 0.5 * m_wirePitch[2]) < .05 && std::fabs(std::fabs(deltaY_uv) - 0.5774 * m_wirePitch[2]) < 0.05)
880  if (yzdistance < 0.3)
881  {
882  // Get a copy of the input hit vector (note the order is by plane - by definition)
883  reco::ClusterHit2DVec hitVector = pair.getHits();
884 
885  // include the new hit
886  hitVector.at(hit->getHit().WireID().Plane) = hit;
887 
888  // Set up to get average peak time, hitChiSquare, etc.
889  unsigned int statusBits(0x7);
890  float totalCharge(0.);
891  float avePeakTime(0.);
892  float weightSum(0.);
893  float xPosition(0.);
894 
895  // And get the wire IDs
896  std::vector<geo::WireID> wireIDVec = {geo::WireID(0,0,geo::kU,0), geo::WireID(0,0,geo::kV,0), geo::WireID(0,0,geo::kW,0)};
897 
898  // First loop through the hits to get WireIDs and calculate the averages
899  for(size_t planeIdx = 0; planeIdx < 3; planeIdx++)
900  {
901  const reco::ClusterHit2D* hit2D = hitVector.at(planeIdx);
902 
903  wireIDVec.at(planeIdx) = hit2D->getHit().WireID();
904 
906 
908 
909  float hitRMS = rmsToSig * hit2D->getHit().RMS();
910  float weight = 1. / (hitRMS * hitRMS);
911  float peakTime = hit2D->getTimeTicks();
912 
913  avePeakTime += peakTime * weight;
914  xPosition += hit2D->getXPosition() * weight;
915  weightSum += weight;
916  totalCharge += hit2D->getHit().Integral();
917  }
918 
919  avePeakTime /= weightSum;
920  xPosition /= weightSum;
921 
922  float position[] = { xPosition,
923  float((pair.getPosition()[1] + pair0h.getPosition()[1] + pair1h.getPosition()[1]) / 3.),
924  float((pair.getPosition()[2] + pair0h.getPosition()[2] + pair1h.getPosition()[2]) / 3.)};
925 
926  // Armed with the average peak time, now get hitChiSquare and the sig vec
927  float hitChiSquare(0.);
928  float sigmaPeakTime(std::sqrt(1./weightSum));
929  std::vector<float> hitDelTSigVec;
930 
931  for(const auto& hit2D : hitVector)
932  {
933  float hitRMS = rmsToSig * hit2D->getHit().RMS();
934  float combRMS = std::sqrt(hitRMS*hitRMS - sigmaPeakTime*sigmaPeakTime);
935  float peakTime = hit2D->getTimeTicks();
936  float deltaTime = peakTime - avePeakTime;
937  float hitSig = deltaTime / combRMS; //hitRMS;
938 
939  hitChiSquare += hitSig * hitSig;
940 
941  hitDelTSigVec.emplace_back(std::fabs(hitSig));
942  }
943 
944  // Usurping "deltaPeakTime" to be the maximum pull
945  float deltaPeakTime = *std::max_element(hitDelTSigVec.begin(),hitDelTSigVec.end());
946 
947  // Create the 3D cluster hit
948  hitTriplet.initialize(0,
949  statusBits,
950  position,
951  totalCharge,
952  avePeakTime,
953  deltaPeakTime,
954  sigmaPeakTime,
955  hitChiSquare,
956  0.,
957  0.,
958  hitVector,
959  hitDelTSigVec,
960  wireIDVec);
961 
962  result = true;
963  }
964  }
965  }
966  }
967 
968  // return success/fail
969  return result;
970 }
971 
973  const reco::ClusterHit3D& pair,
974  size_t maxChanStatus,
975  size_t minChanStatus,
976  float minOverlap) const
977 {
978  // Assume failure (most common result)
979  bool result(false);
980 
981  const reco::ClusterHit2D* hit0 = pair.getHits().at(0);
982  const reco::ClusterHit2D* hit1 = pair.getHits().at(1);
983 
984  size_t missPlane(2);
985 
986  // u plane hit is missing
987  if (!hit0)
988  {
989  hit0 = pair.getHits().at(2);
990  missPlane = 0;
991  }
992  // v plane hit is missing
993  else if (!hit1)
994  {
995  hit1 = pair.getHits().at(2);
996  missPlane = 1;
997  }
998 
999  // Which plane is missing?
1000  geo::WireID wireID0 = hit0->getHit().WireID();
1001  geo::WireID wireID1 = hit1->getHit().WireID();
1002 
1003  // Ok, recover the wireID expected in the third plane...
1004  geo::WireID wireIn(wireID0.Cryostat,wireID0.TPC,missPlane,0);
1005  geo::WireID wireID = NearestWireID(pair.getPosition(), wireIn);
1006 
1007  // There can be a round off issue so check the next wire as well
1008  bool wireStatus = m_channelStatus[wireID.Plane][wireID.Wire] < maxChanStatus && m_channelStatus[wireID.Plane][wireID.Wire] >= minChanStatus;
1009  bool wireOneStatus = m_channelStatus[wireID.Plane][wireID.Wire+1] < maxChanStatus && m_channelStatus[wireID.Plane][wireID.Wire+1] >= minChanStatus;
1010 
1011  // Make sure they are of at least the minimum status
1012  if(wireStatus || wireOneStatus)
1013  {
1014  // Sort out which is the wire we're dealing with
1015  if (!wireStatus) wireID.Wire += 1;
1016 
1017  // Want to refine position since we "know" the missing wire
1018  geo::WireIDIntersection widIntersect0;
1019 
1020  if (m_geometry->WireIDsIntersect(wireID0, wireID, widIntersect0))
1021  {
1022  geo::WireIDIntersection widIntersect1;
1023 
1024  if (m_geometry->WireIDsIntersect(wireID1, wireID, widIntersect1))
1025  {
1026  float newPosition[] = {pair.getPosition()[0],pair.getPosition()[1],pair.getPosition()[2]};
1027 
1028  newPosition[1] = (newPosition[1] + widIntersect0.y + widIntersect1.y) / 3.;
1029  newPosition[2] = (newPosition[2] + widIntersect0.z + widIntersect1.z - 2. * m_zPosOffset) / 3.;
1030 
1031  pairOut = pair;
1032  pairOut.setWireID(wireID);
1033  pairOut.setPosition(newPosition);
1034 
1037 
1040 
1041  result = true;
1042  }
1043  }
1044  }
1045 
1046  return result;
1047 }
1048 
1049 const reco::ClusterHit2D* StandardHit3DBuilder::FindBestMatchingHit(const Hit2DSet& hit2DSet, const reco::ClusterHit3D& pair, float pairDeltaTimeLimits) const
1050 {
1051  static const float minCharge(0.);
1052 
1053  const reco::ClusterHit2D* bestVHit(0);
1054 
1055  float pairAvePeakTime(pair.getAvePeakTime());
1056 
1057  // Idea is to loop through the input set of hits and look for the best combination
1058  for (const auto& hit2D : hit2DSet)
1059  {
1060  if (hit2D->getHit().Integral() < minCharge) continue;
1061 
1062  float hitVPeakTime(hit2D->getTimeTicks());
1063  float deltaPeakTime(pairAvePeakTime-hitVPeakTime);
1064 
1065  if (deltaPeakTime > pairDeltaTimeLimits) continue;
1066 
1067  if (deltaPeakTime < -pairDeltaTimeLimits) break;
1068 
1069  pairDeltaTimeLimits = fabs(deltaPeakTime);
1070  bestVHit = hit2D;
1071  }
1072 
1073  return bestVHit;
1074 }
1075 
1076 int StandardHit3DBuilder::FindNumberInRange(const Hit2DSet& hit2DSet, const reco::ClusterHit3D& pair, float range) const
1077 {
1078  static const float minCharge(0.);
1079 
1080  int numberInRange(0);
1081  float pairAvePeakTime(pair.getAvePeakTime());
1082 
1083  // Idea is to loop through the input set of hits and look for the best combination
1084  for (const auto& hit2D : hit2DSet)
1085  {
1086  if (hit2D->getHit().Integral() < minCharge) continue;
1087 
1088  float hitVPeakTime(hit2D->getTimeTicks());
1089  float deltaPeakTime(pairAvePeakTime-hitVPeakTime);
1090 
1091  if (deltaPeakTime > range) continue;
1092 
1093  if (deltaPeakTime < -range) break;
1094 
1095  numberInRange++;
1096  }
1097 
1098  return numberInRange;
1099 }
1100 
1101 geo::WireID StandardHit3DBuilder::NearestWireID(const float* position, const geo::WireID& wireIDIn) const
1102 {
1103  geo::WireID wireID = wireIDIn;
1104 
1105  // Embed the call to the geometry's services nearest wire id method in a try-catch block
1106  try
1107  {
1108  // Switch from NearestWireID to this method to avoid the roundoff error issues...
1109  double distanceToWire = m_geometry->Plane(wireIDIn).WireCoordinate(position);
1110 
1111  wireID.Wire = int(distanceToWire);
1112  }
1113  catch(std::exception& exc)
1114  {
1115  // This can happen, almost always because the coordinates are **just** out of range
1116  mf::LogWarning("Cluster3D") << "Exception caught finding nearest wire, position - " << exc.what() << std::endl;
1117 
1118  // Assume extremum for wire number depending on z coordinate
1119  if (position[2] < 0.5 * m_geometry->DetLength()) wireID.Wire = 0;
1120  else wireID.Wire = m_geometry->Nwires(wireIDIn.Plane) - 1;
1121  }
1122 
1123  return wireID;
1124 }
1125 
1126 //------------------------------------------------------------------------------------------------------------------------------------------
1128 {
1129  // Sort by "modified start time" of pulse
1130  return left->getHit().PeakTime() < right->getHit().PeakTime();
1131 }
1132 
1134 {
1135  return left->getHit().PeakTime() < right->getHit().PeakTime();
1136 }
1137 
1138 //------------------------------------------------------------------------------------------------------------------------------------------
1140  RecobHitToPtrMap& hitToPtrMap) const
1141 {
1145  art::Handle< std::vector<recob::Hit> > recobHitHandle;
1146  evt.getByLabel(m_hitFinderTag, recobHitHandle);
1147 
1148  if (!recobHitHandle.isValid()) return;
1149 
1150  cet::cpu_timer theClockMakeHits;
1151 
1152  if (m_enableMonitoring) theClockMakeHits.start();
1153 
1154  // We'll want to correct the hit times for the plane offsets
1155  // (note this is already taken care of when converting to position)
1156  std::map<geo::PlaneID, double> planeOffsetMap;
1157 
1158  // Reserve memory for the hit vector
1159  m_clusterHit2DMasterVec.reserve(recobHitHandle->size());
1160 
1161  // Initialize the plane to hit vector map
1162  for(size_t cryoIdx = 0; cryoIdx < m_geometry->Ncryostats(); cryoIdx++)
1163  {
1164  for(size_t tpcIdx = 0; tpcIdx < m_geometry->NTPC(); tpcIdx++)
1165  {
1166  m_planeToHitVectorMap[geo::PlaneID(cryoIdx,tpcIdx,0)] = HitVector();
1167  m_planeToHitVectorMap[geo::PlaneID(cryoIdx,tpcIdx,1)] = HitVector();
1168  m_planeToHitVectorMap[geo::PlaneID(cryoIdx,tpcIdx,2)] = HitVector();
1169 
1170  // What we want here are the relative offsets between the planes
1171  // Note that plane 0 is assumed the "first" plane and is the reference
1172  planeOffsetMap[geo::PlaneID(cryoIdx,tpcIdx,0)] = 0.;
1173  planeOffsetMap[geo::PlaneID(cryoIdx,tpcIdx,1)] = m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,1))
1174  - m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,0));
1175  planeOffsetMap[geo::PlaneID(cryoIdx,tpcIdx,2)] = m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,2))
1176  - m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,0));
1177 
1178  std::cout << "***> plane 0 offset: " << planeOffsetMap[geo::PlaneID(cryoIdx,tpcIdx,0)] << ", plane 1: " << planeOffsetMap[geo::PlaneID(cryoIdx,tpcIdx,1)] << ", plane 2: " << planeOffsetMap[geo::PlaneID(cryoIdx,tpcIdx,2)] << std::endl;
1179  std::cout << " Det prop plane 0: " << m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,0)) << ", plane 1: " << m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,1)) << ", plane 2: " << m_detector->GetXTicksOffset(geo::PlaneID(cryoIdx,tpcIdx,2)) << ", Trig: " << m_detector->TriggerOffset() << std::endl;
1180  }
1181  }
1182 
1183  // Cycle through the recob hits to build ClusterHit2D objects and insert
1184  // them into the map
1185  for (size_t cIdx = 0; cIdx < recobHitHandle->size(); cIdx++)
1186  {
1187  art::Ptr<recob::Hit> recobHit(recobHitHandle, cIdx);
1188 
1189  const geo::WireID& hitWireID(recobHit->WireID());
1190 
1191  double hitPeakTime(recobHit->PeakTime() - planeOffsetMap[recobHit->WireID().planeID()]);
1192  double xPosition(m_detector->ConvertTicksToX(recobHit->PeakTime(), hitWireID.Plane, hitWireID.TPC, hitWireID.Cryostat));
1193 
1194  m_clusterHit2DMasterVec.emplace_back(0, 0., 0., xPosition, hitPeakTime, *recobHit);
1195 
1196  m_planeToHitVectorMap[recobHit->WireID().planeID()].push_back(&m_clusterHit2DMasterVec.back());
1197  m_planeToWireToHitSetMap[recobHit->WireID().planeID()][recobHit->WireID().Wire].insert(&m_clusterHit2DMasterVec.back());
1198 
1199  const recob::Hit* recobHitPtr = recobHit.get();
1200  hitToPtrMap[recobHitPtr] = recobHit;
1201  }
1202 
1203  // Make a loop through to sort the recover hits in time order
1204  for(auto& hitVectorMap : m_planeToHitVectorMap)
1205  std::sort(hitVectorMap.second.begin(), hitVectorMap.second.end(), SetHitTimeOrder);
1206 
1207  if (m_enableMonitoring)
1208  {
1209  theClockMakeHits.stop();
1210 
1211  m_timeVector[COLLECTARTHITS] = theClockMakeHits.accumulated_real_time();
1212  }
1213 
1214  mf::LogDebug("Cluster3D") << ">>>>> Number of ART hits: " << m_clusterHit2DMasterVec.size() << std::endl;
1215 }
1216 
1217 //------------------------------------------------------------------------------------------------------------------------------------------
1218 
1220 } // namespace lar_cluster3d
bool makeHitTriplet(reco::ClusterHit3D &pairOut, const reco::ClusterHit3D &pairIn, const reco::ClusterHit2D *hit2) const
Make a 3D HitPair object by checking two hits.
PlaneID const & planeID() const
Definition: geo_types.h:355
void setWireID(const geo::WireID &wid) const
Definition: Cluster3D.cxx:155
double z
z position of intersection
Definition: geo_types.h:494
std::map< size_t, HitVector > HitVectorMap
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:112
#define DEFINE_ART_CLASS_TOOL(tool)
Definition: ToolMacros.h:45
PlaneGeo const & Plane(unsigned int const p, unsigned int const tpc=0, unsigned int const cstat=0) const
Returns the specified wire.
float getTimeTicks() const
Definition: Cluster3D.h:72
virtual int TriggerOffset() const =0
const detinfo::DetectorProperties * m_detector
intermediate_table::iterator iterator
std::vector< size_t > ChannelStatusVec
define data structure for keeping track of channel status
geo::WireID NearestWireID(const float *position, const geo::WireID &wireID) const
Jacket the calls to finding the nearest wire in order to intercept the exceptions if out of range...
bool SetHitTimeOrder(const reco::ClusterHit2D *left, const reco::ClusterHit2D *right)
Planes which measure V.
Definition: geo_types.h:77
geo::WireID WireID() const
Initial tdc tick for hit.
Definition: Hit.h:234
float RMS() const
RMS of the hit shape, in tick units.
Definition: Hit.h:221
What follows are several highly useful typedefs which we want to expose to the outside world...
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:250
art::InputTag m_hitFinderTag
Data members to follow.
int findGoodHitPairs(const reco::ClusterHit2D *, HitVector::iterator &, HitVector::iterator &, HitMatchPairVecMap &) const
bool operator()(const reco::ClusterHit2D *, const reco::ClusterHit2D *) const
std::vector< geo::WireID > ChannelToWire(raw::ChannelID_t const channel) const
Returns a list of wires connected to the specified TPC channel.
CryostatID_t Cryostat
Index of cryostat.
Definition: geo_types.h:130
void BuildHit3D(reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
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.
float Integral() const
Integral under the calibrated signal waveform of the hit, in tick x ADC units.
Definition: Hit.h:225
WireID_t Wire
Index of the wire within its plane.
Definition: geo_types.h:313
float getSigmaPeakTime() const
Definition: Cluster3D.h:152
unsigned int Ncryostats() const
Returns the number of cryostats in the detector.
std::map< geo::PlaneID, WireToHitSetMap > PlaneToWireToHitSetMap
geo::Length_t WirePitch(geo::PlaneID const &planeid) const
Returns the distance between two consecutive wires.
unsigned int Nwires(unsigned int p, unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wires in the specified plane.
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 Hit3DBuilder(const art::Event &evt, reco::HitPairList &hitPairList, RecobHitToPtrMap &) const override
Given a set of recob hits, run DBscan to form 3D clusters.
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:222
std::list< std::unique_ptr< reco::ClusterHit3D >> HitPairList
Definition: Cluster3D.h:319
std::vector< reco::ClusterHit2D * > HitVector
bool isValid() const
Definition: Handle.h:190
float getAvePeakTime() const
Definition: Cluster3D.h:150
std::vector< std::unique_ptr< reco::ClusterHit3D >> HitPairVector
Planes which measure U.
Definition: geo_types.h:76
unsigned int Nplanes(unsigned int tpc=0, unsigned int cstat=0) const
Returns the total number of wire planes in the specified TPC.
Int_t max
Definition: plot.C:27
std::map< geo::TPCID, PlaneToWireToHitSetMap > TPCToPlaneToWireToHitSetMap
std::map< geo::PlaneID, HitVector > PlaneToHitVectorMap
T get(std::string const &key) const
Definition: ParameterSet.h:231
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.
const recob::Hit & getHit() const
Definition: Cluster3D.h:73
geo::Length_t DetLength(geo::TPCID const &tpcid) const
Returns the length of the active volume of the specified TPC.
float getXPosition() const
Definition: Cluster3D.h:71
TimeValues
enumerate the possible values for time checking if monitoring timing
Definition: IHit3DBuilder.h:57
The geometry of one entire detector, as served by art.
Definition: Geometry.h:110
std::map< const recob::Hit *, art::Ptr< recob::Hit >> RecobHitToPtrMap
Defines a structure mapping art representation to internal.
Definition: IHit3DBuilder.h:44
PlaneID_t Plane
Index of the plane within its TPC.
Definition: geo_types.h:258
IHit3DBuilder interface class definiton.
Definition: IHit3DBuilder.h:26
const float * getPosition() const
Definition: Cluster3D.h:145
std::pair< const reco::ClusterHit2D *, reco::ClusterHit3D > HitMatchPair
This builds a list of candidate hit pairs from lists of hits on two planes.
Detector simulation of raw signals on wires.
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:85
void findGoodTriplets(HitMatchPairVecMap &, HitMatchPairVecMap &, reco::HitPairList &, bool=false) const
This algorithm takes lists of hit pairs and finds good triplets.
unsigned int NTPC(unsigned int cstat=0) const
Returns the total number of TPCs in the specified cryostat.
std::map< unsigned int, Hit2DSet > WireToHitSetMap
const lariov::ChannelStatusProvider * m_channelFilter
void initialize(size_t id, unsigned int statusBits, const float *position, float totalCharge, float avePeakTime, float deltaPeakTime, float sigmaPeakTime, float hitChiSquare, float docaToAxis, float arclenToPoca, const ClusterHit2DVec &hitVec, const std::vector< float > &hitDelTSigVec, const std::vector< geo::WireID > &wireIDVec)
Definition: Cluster3D.cxx:122
void CollectArtHits(const art::Event &evt, RecobHitToPtrMap &hitToPtrMap) const
Extract the ART hits and the ART hit-particle relationships.
Filters for channels, events, etc.
bool WireIDsIntersect(WireID const &wid1, WireID const &wid2, geo::Point_t &intersection) const
Computes the intersection between two wires.
virtual double ConvertTicksToX(double ticks, int p, int t, int c) const =0
constexpr auto const & left(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
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:219
Utility object to perform functions of association.
void setID(const size_t &id) const
Definition: Cluster3D.h:162
T const * get() const
Definition: Ptr.h:321
void setPosition(const float *pos) const
Definition: Cluster3D.h:169
bool getByLabel(std::string const &label, std::string const &productInstanceName, Handle< PROD > &result) const
Definition: DataViewImpl.h:344
std::vector< reco::ClusterHit2D > Hit2DVector
std::vector< ChannelStatusVec > ChannelStatusByPlaneVec
double y
y position of intersection
Definition: geo_types.h:493
MaybeLogger_< ELseverityLevel::ELsev_success, false > LogDebug
size_t BuildHitPairMap(PlaneToHitVectorMap &planeToHitVectorMap, reco::HitPairList &hitPairList) const
Given the ClusterHit2D objects, build the HitPairMap.
MaybeLogger_< ELseverityLevel::ELsev_warning, false > LogWarning
std::map< geo::WireID, HitMatchPairVec > HitMatchPairVecMap
std::vector< std::pair< HitVector::iterator, HitVector::iterator >> PlaneHitVectorItrPairVec
Given the ClusterHit2D objects, build the HitPairMap.
HLT enums.
StandardHit3DBuilder(fhicl::ParameterSet const &pset)
Constructor.
Planes which measure W (third view for Bo, MicroBooNE, etc).
Definition: geo_types.h:78
2D representation of charge deposited in the TDC/wire plane
Definition: Hit.h:49
bool makeDeadChannelPair(reco::ClusterHit3D &pairOut, const reco::ClusterHit3D &pair, size_t maxStatus=4, size_t minStatus=0, float minOverlap=0.2) const
Make a 3D HitPair object from a valid pair and a dead channel in the missing plane.
double WireCoordinate(Point const &point) const
Returns the coordinate of the point on the plane, in wire units.
Definition: PlaneGeo.h:725
TPCID_t TPC
Index of the TPC within its cryostat.
Definition: geo_types.h:203
float getTimeToExecute(IHit3DBuilder::TimeValues index) const override
If monitoring, recover the time to execute a particular function.
const ClusterHit2DVec & getHits() const
Definition: Cluster3D.h:156
std::set< const reco::ClusterHit2D *, Hit2DSetCompare > Hit2DSet
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
void BuildChannelStatusVec(PlaneToWireToHitSetMap &planeToWiretoHitSetMap) const
Create the internal channel status vector (assume will eventually be event-by-event) ...
unsigned getStatusBits() const
Definition: Cluster3D.h:68
virtual double GetXTicksOffset(int p, int t, int c) const =0
void setStatusBit(unsigned bits) const
Definition: Cluster3D.h:75