CMS 3D CMS Logo

TrackerHitAssociator.cc
Go to the documentation of this file.
1 // File: TrackerHitAssociator.cc
2 
3 #include <memory>
4 #include <string>
5 #include <vector>
6 
8 
10 
11 //--- for Geometry:
16 
18 
19 //for accumulate
20 #include <numeric>
21 #include <iostream>
22 
23 using namespace std;
24 using namespace edm;
25 
26 //
27 // Constructor for Config helper class, using default parameters
28 //
30  : doPixel_(true), doStrip_(true), useOTph2_(false), doTrackAssoc_(false), assocHitbySimTrack_(false) {
31  if (doStrip_) {
32  if (useOTph2_)
33  ph2OTrToken_ = iC.consumes<edm::DetSetVector<PixelDigiSimLink>>(edm::InputTag("simSiPixelDigis", "Tracker"));
34  else
35  stripToken_ = iC.consumes<edm::DetSetVector<StripDigiSimLink>>(edm::InputTag("simSiStripDigis"));
36  }
37  if (doPixel_) {
38  if (useOTph2_)
39  pixelToken_ = iC.consumes<edm::DetSetVector<PixelDigiSimLink>>(edm::InputTag("simSiPixelDigis", "Pixel"));
40  else
41  pixelToken_ = iC.consumes<edm::DetSetVector<PixelDigiSimLink>>(edm::InputTag("simSiPixelDigis"));
42  }
43  if (!doTrackAssoc_) {
44  std::vector<std::string> trackerContainers;
45  trackerContainers.reserve(12);
46  trackerContainers.emplace_back("g4SimHitsTrackerHitsTIBLowTof");
47  trackerContainers.emplace_back("g4SimHitsTrackerHitsTIBHighTof");
48  trackerContainers.emplace_back("g4SimHitsTrackerHitsTIDLowTof");
49  trackerContainers.emplace_back("g4SimHitsTrackerHitsTIDHighTof");
50  trackerContainers.emplace_back("g4SimHitsTrackerHitsTOBLowTof");
51  trackerContainers.emplace_back("g4SimHitsTrackerHitsTOBHighTof");
52  trackerContainers.emplace_back("g4SimHitsTrackerHitsTECLowTof");
53  trackerContainers.emplace_back("g4SimHitsTrackerHitsTECHighTof");
54  trackerContainers.emplace_back("g4SimHitsTrackerHitsPixelBarrelLowTof");
55  trackerContainers.emplace_back("g4SimHitsTrackerHitsPixelBarrelHighTof");
56  trackerContainers.emplace_back("g4SimHitsTrackerHitsPixelEndcapLowTof");
57  trackerContainers.emplace_back("g4SimHitsTrackerHitsPixelEndcapHighTof");
58  cfTokens_.reserve(trackerContainers.size());
59  simHitTokens_.reserve(trackerContainers.size());
60  for (auto const& trackerContainer : trackerContainers) {
61  cfTokens_.push_back(iC.consumes<CrossingFrame<PSimHit>>(edm::InputTag("mix", trackerContainer)));
62  simHitTokens_.push_back(iC.consumes<std::vector<PSimHit>>(edm::InputTag("g4SimHits", trackerContainer)));
63  }
64  }
65 }
66 
67 //
68 // Constructor for Config helper class, using configured parameters
69 //
71  : doPixel_(conf.getParameter<bool>("associatePixel")),
72  doStrip_(conf.getParameter<bool>("associateStrip")),
73  useOTph2_(conf.existsAs<bool>("usePhase2Tracker") ? conf.getParameter<bool>("usePhase2Tracker") : false),
74  //
75  doTrackAssoc_(conf.getParameter<bool>("associateRecoTracks")),
77  conf.existsAs<bool>("associateHitbySimTrack") ? conf.getParameter<bool>("associateHitbySimTrack") : false) {
78  if (doStrip_) {
79  if (useOTph2_)
80  ph2OTrToken_ =
81  iC.consumes<edm::DetSetVector<PixelDigiSimLink>>(conf.getParameter<edm::InputTag>("phase2TrackerSimLinkSrc"));
82  else
83  stripToken_ =
84  iC.consumes<edm::DetSetVector<StripDigiSimLink>>(conf.getParameter<edm::InputTag>("stripSimLinkSrc"));
85  }
86  if (doPixel_)
87  pixelToken_ = iC.consumes<edm::DetSetVector<PixelDigiSimLink>>(conf.getParameter<edm::InputTag>("pixelSimLinkSrc"));
88  if (!doTrackAssoc_) {
89  std::vector<std::string> trackerContainers(conf.getParameter<std::vector<std::string>>("ROUList"));
90  cfTokens_.reserve(trackerContainers.size());
91  simHitTokens_.reserve(trackerContainers.size());
92  for (auto const& trackerContainer : trackerContainers) {
93  cfTokens_.push_back(iC.consumes<CrossingFrame<PSimHit>>(edm::InputTag("mix", trackerContainer)));
94  simHitTokens_.push_back(iC.consumes<std::vector<PSimHit>>(edm::InputTag("g4SimHits", trackerContainer)));
95  }
96  }
97 }
98 
100  desc.setComment("auxilliary class to store information about recHit/simHit association");
101  desc.add<bool>("associatePixel", false);
102  desc.add<bool>("associateStrip", false);
103  desc.add<bool>("usePhase2Tracker", false);
104  desc.add<bool>("associateRecoTracks", false);
105  desc.add<bool>("associateHitbySimTrack", false);
106  desc.add<edm::InputTag>("phase2TrackerSimLinkSrc", edm::InputTag("simSiPixelDigis", "Tracker"));
107  desc.add<edm::InputTag>("stripSimLinkSrc", edm::InputTag("simSiStripDigis"));
108  desc.add<edm::InputTag>("pixelSimLinkSrc", edm::InputTag("simSiPixelDigis"));
109  desc.add<std::vector<std::string>>(
110  "ROUList", {"TrackerHitsTIBLowTof", "TrackerHitsTIBHighTof", "TrackerHitsTOBLowTof", "TrackerHitsTOBHighTof"});
111 }
112 
113 //
114 // Constructor supporting consumes interface
115 //
122  //if track association there is no need to access the input collections
123  if (!doTrackAssoc_) {
124  makeMaps(e, config);
125  }
126 
127  if (doStrip_) {
128  if (useOTph2_)
129  e.getByToken(config.ph2OTrToken_, ph2trackerdigisimlink);
130  else
131  e.getByToken(config.stripToken_, stripdigisimlink);
132  }
133  if (doPixel_)
134  e.getByToken(config.pixelToken_, pixeldigisimlink);
135 }
136 
138  // Step A: Get Inputs
139  // The collections are specified via ROUList in the configuration, and can
140  // be either crossing frames (e.g., mix/g4SimHitsTrackerHitsTIBLowTof)
141  // or just PSimHits (e.g., g4SimHits/TrackerHitsTIBLowTof)
142  if (assocHitbySimTrack_) {
143  for (auto const& cfToken : config.cfTokens_) {
145  int Nhits = 0;
146  if (theEvent.getByToken(cfToken, cf_simhit)) {
147  std::unique_ptr<MixCollection<PSimHit>> thisContainerHits(new MixCollection<PSimHit>(cf_simhit.product()));
148  for (auto const& isim : *thisContainerHits) {
149  DetId theDet(isim.detUnitId());
150  SimHitMap[theDet].push_back(isim);
151  ++Nhits;
152  }
153  LogDebug("TrkHitAssocTrace") << "simHits from crossing frames; map size = " << SimHitMap.size()
154  << ", Hit count = " << Nhits << std::endl;
155  }
156  }
157  for (auto const& simHitToken : config.simHitTokens_) {
159  int Nhits = 0;
160  if (theEvent.getByToken(simHitToken, simHits)) {
161  for (auto const& isim : *simHits) {
162  DetId theDet(isim.detUnitId());
163  SimHitMap[theDet].push_back(isim);
164  ++Nhits;
165  }
166  LogDebug("TrkHitAssocTrace") << "simHits from prompt collections; map size = " << SimHitMap.size()
167  << ", Hit count = " << Nhits << std::endl;
168  }
169  }
170  } else { // !assocHitbySimTrack_
171  const char* const highTag = "HighTof";
172  unsigned int tofBin;
174  subDetTofBin theSubDetTofBin;
175  unsigned int collectionIndex = 0;
176  for (auto const& cfToken : config.cfTokens_) {
177  collectionIndex++;
179  int Nhits = 0;
180  if (theEvent.getByToken(cfToken, cf_simhit)) {
181  std::unique_ptr<MixCollection<PSimHit>> thisContainerHits(new MixCollection<PSimHit>(cf_simhit.product()));
182  theEvent.labelsForToken(cfToken, labels);
183  if (std::strstr(labels.productInstance, highTag) != nullptr) {
184  tofBin = StripDigiSimLink::HighTof;
185  } else {
186  tofBin = StripDigiSimLink::LowTof;
187  }
188  for (auto const& isim : *thisContainerHits) {
189  DetId theDet(isim.detUnitId());
190  theSubDetTofBin = std::make_pair(theDet.subdetId(), tofBin);
191  SimHitCollMap[theSubDetTofBin] = collectionIndex;
192  SimHitMap[SimHitCollMap[theSubDetTofBin]].push_back(isim);
193  ++Nhits;
194  }
195  LogDebug("TrkHitAssocTrace") << "simHits from crossing frames " << collectionIndex << ": " << Nhits
196  << std::endl;
197  }
198  }
199  collectionIndex = 0;
200  for (auto const& simHitToken : config.simHitTokens_) {
201  collectionIndex++;
203  int Nhits = 0;
204  if (theEvent.getByToken(simHitToken, simHits)) {
205  theEvent.labelsForToken(simHitToken, labels);
206  if (std::strstr(labels.productInstance, highTag) != nullptr) {
207  tofBin = StripDigiSimLink::HighTof;
208  } else {
209  tofBin = StripDigiSimLink::LowTof;
210  }
211  for (auto const& isim : *simHits) {
212  DetId theDet(isim.detUnitId());
213  theSubDetTofBin = std::make_pair(theDet.subdetId(), tofBin);
214  SimHitCollMap[theSubDetTofBin] = collectionIndex;
215  SimHitMap[SimHitCollMap[theSubDetTofBin]].push_back(isim);
216  ++Nhits;
217  }
218  LogDebug("TrkHitAssocTrace") << "simHits from prompt collection " << collectionIndex << ": " << Nhits
219  << std::endl;
220  }
221  }
222  }
223 }
224 
225 std::vector<PSimHit> TrackerHitAssociator::associateHit(const TrackingRecHit& thit) const {
226  if (const SiTrackerMultiRecHit* rechit = dynamic_cast<const SiTrackerMultiRecHit*>(&thit)) {
227  return associateMultiRecHit(rechit);
228  }
229 
230  //vector with the matched SimHit
231  std::vector<PSimHit> result;
232 
233  if (doTrackAssoc_)
234  return result; // We don't want the SimHits for this RecHit
235 
236  // Vectors to contain lists of matched simTracks, simHits
237  std::vector<SimHitIdpr> simtrackid;
238  std::vector<simhitAddr> simhitCFPos;
239 
240  //get the Detector type of the rechit
241  DetId detid = thit.geographicalId();
242  uint32_t detID = detid.rawId();
243 
244  // Get the vectors of simtrackIDs and simHit addresses associated with this rechit
245  associateHitId(thit, simtrackid, &simhitCFPos);
246  LogDebug("TrkHitAssocTrace") << printDetBnchEvtTrk(detid, detID, simtrackid);
247 
248  // Get the vector of simHits associated with this rechit
249  if (!assocHitbySimTrack_ && !simhitCFPos.empty()) {
250  // We use the indices to the simHit collections taken
251  // from the DigiSimLinks and returned in simhitCFPos.
252  // simhitCFPos[i] contains the full address of the ith simhit:
253  // <collection index, simhit index>
254 
255  //check if the recHit is a SiStripMatchedRecHit2D
256  if (dynamic_cast<const SiStripMatchedRecHit2D*>(&thit)) {
257  for (auto const& theSimHitAddr : simhitCFPos) {
258  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
259  auto it = SimHitMap.find(theSimHitCollID);
260  if (it != SimHitMap.end()) {
261  unsigned int theSimHitIndex = theSimHitAddr.second;
262  if (theSimHitIndex < (it->second).size()) {
263  const PSimHit& theSimHit = (it->second)[theSimHitIndex];
264  // Try to remove ghosts by requiring a match to the simTrack also
265  unsigned int simHitid = theSimHit.trackId();
266  EncodedEventId simHiteid = theSimHit.eventId();
267  for (auto const& id : simtrackid) {
268  if (simHitid == id.first && simHiteid == id.second) {
269  result.push_back(theSimHit);
270  }
271  }
272  LogDebug("TrkHitAssocTrace") << "by CFpos, simHit detId = " << theSimHit.detUnitId() << " address = ("
273  << theSimHitAddr.first << ", " << theSimHitIndex
274  << "), process = " << theSimHit.processType() << " ("
275  << theSimHit.eventId().bunchCrossing() << ", " << theSimHit.eventId().event()
276  << ", " << theSimHit.trackId() << ")" << std::endl;
277  }
278  }
279  }
280  } else { // Not a SiStripMatchedRecHit2D
281  for (auto const& theSimHitAddr : simhitCFPos) {
282  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
283  auto it = SimHitMap.find(theSimHitCollID);
284  if (it != SimHitMap.end()) {
285  unsigned int theSimHitIndex = theSimHitAddr.second;
286  if (theSimHitIndex < (it->second).size()) {
287  result.push_back((it->second)[theSimHitIndex]);
288  LogDebug("TrkHitAssocTrace") << "by CFpos, simHit detId = " << (it->second)[theSimHitIndex].detUnitId()
289  << " address = (" << theSimHitCollID << ", " << theSimHitIndex
290  << "), process = " << (it->second)[theSimHitIndex].processType() << " ("
291  << (it->second)[theSimHitIndex].eventId().bunchCrossing() << ", "
292  << (it->second)[theSimHitIndex].eventId().event() << ", "
293  << (it->second)[theSimHitIndex].trackId() << ")" << std::endl;
294  }
295  }
296  }
297  }
298  return result;
299  } // if !assocHitbySimTrack
300 
301  // Get the SimHit from the trackid instead
302  auto it = SimHitMap.find(detID);
303  if (it != SimHitMap.end()) {
304  for (auto const& ihit : it->second) {
305  unsigned int simHitid = ihit.trackId();
306  EncodedEventId simHiteid = ihit.eventId();
307  for (auto id : simtrackid) {
308  if (simHitid == id.first && simHiteid == id.second) {
309  result.push_back(ihit);
310  LogDebug("TrkHitAssocTrace") << "by TrackID, simHit detId = " << ihit.detUnitId()
311  << ", process = " << ihit.processType() << " (" << ihit.eventId().bunchCrossing()
312  << ", " << ihit.eventId().event() << ", " << ihit.trackId() << ")" << std::endl;
313  break;
314  }
315  }
316  }
317 
318  } else {
320  auto itrphi = SimHitMap.find(detID + 2); //iterator to the simhit in the rphi module
321  auto itster = SimHitMap.find(detID + 1); //iterator to the simhit in the stereo module
322  if (itrphi != SimHitMap.end() && itster != SimHitMap.end()) {
323  std::vector<PSimHit> simHitVector = itrphi->second;
324  simHitVector.insert(simHitVector.end(), (itster->second).begin(), (itster->second).end());
325  for (auto const& ihit : simHitVector) {
326  unsigned int simHitid = ihit.trackId();
327  EncodedEventId simHiteid = ihit.eventId();
328  for (auto const& id : simtrackid) {
329  if (simHitid == id.first && simHiteid == id.second) {
330  result.push_back(ihit);
331  LogDebug("TrkHitAssocTrace") << "by TrackID, simHit detId = " << ihit.detUnitId()
332  << ", process = " << ihit.processType() << " ("
333  << ihit.eventId().bunchCrossing() << ", " << ihit.eventId().event() << ", "
334  << ihit.trackId() << ")" << std::endl;
335  break;
336  }
337  }
338  }
339  }
340  }
341 
342  return result;
343 }
344 
345 std::vector<SimHitIdpr> TrackerHitAssociator::associateHitId(const TrackingRecHit& thit) const {
346  std::vector<SimHitIdpr> simhitid;
347  associateHitId(thit, simhitid);
348  return simhitid;
349 }
350 
352  std::vector<SimHitIdpr>& simtkid,
353  std::vector<simhitAddr>* simhitCFPos) const {
354  simtkid.clear();
355 
356  if (const SiTrackerMultiRecHit* rechit = dynamic_cast<const SiTrackerMultiRecHit*>(&thit))
357  simtkid = associateMultiRecHitId(rechit, simhitCFPos);
358 
359  //check if it is a simple SiStripRecHit2D
360  if (const SiStripRecHit2D* rechit = dynamic_cast<const SiStripRecHit2D*>(&thit))
361  associateSiStripRecHit(rechit, simtkid, simhitCFPos);
362 
363  //check if it is a SiStripRecHit1D
364  else if (const SiStripRecHit1D* rechit = dynamic_cast<const SiStripRecHit1D*>(&thit))
365  associateSiStripRecHit(rechit, simtkid, simhitCFPos);
366 
367  //check if it is a SiStripMatchedRecHit2D
368  else if (const SiStripMatchedRecHit2D* rechit = dynamic_cast<const SiStripMatchedRecHit2D*>(&thit))
369  simtkid = associateMatchedRecHit(rechit, simhitCFPos);
370 
371  //check if it is a ProjectedSiStripRecHit2D
372  else if (const ProjectedSiStripRecHit2D* rechit = dynamic_cast<const ProjectedSiStripRecHit2D*>(&thit))
373  simtkid = associateProjectedRecHit(rechit, simhitCFPos);
374 
375  //check if it is a Phase2TrackerRecHit1D
376  else if (const Phase2TrackerRecHit1D* rechit = dynamic_cast<const Phase2TrackerRecHit1D*>(&thit))
377  associatePhase2TrackerRecHit(rechit, simtkid, simhitCFPos);
378 
379  //check if it is a SiPixelRecHit
380  else if (const SiPixelRecHit* rechit = dynamic_cast<const SiPixelRecHit*>(&thit))
381  associatePixelRecHit(rechit, simtkid, simhitCFPos);
382 
383  //check if these are GSRecHits (from FastSim)
384  if (trackerHitRTTI::isFast(thit))
385  simtkid = associateFastRecHit(static_cast<const FastTrackerRecHit*>(&thit));
386 }
387 
388 template <typename T>
389 inline void TrackerHitAssociator::associateSiStripRecHit(const T* simplerechit,
390  std::vector<SimHitIdpr>& simtrackid,
391  std::vector<simhitAddr>* simhitCFPos) const {
392  const SiStripCluster* clust = &(*simplerechit->cluster());
393  associateSimpleRecHitCluster(clust, simplerechit->geographicalId(), simtrackid, simhitCFPos);
394 }
395 
396 //
397 // Method for obtaining simTracks and simHits from a cluster
398 //
400  const DetId& detid,
401  std::vector<SimHitIdpr>& simtrackid,
402  std::vector<PSimHit>& simhit) const {
403  std::vector<simhitAddr> simhitCFPos;
404  associateSimpleRecHitCluster(clust, detid, simtrackid, &simhitCFPos);
405 
406  for (auto const& theSimHitAddr : simhitCFPos) {
407  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
408  auto it = SimHitMap.find(theSimHitCollID);
409 
410  if (it != SimHitMap.end()) {
411  unsigned int theSimHitIndex = theSimHitAddr.second;
412  if (theSimHitIndex < (it->second).size())
413  simhit.push_back((it->second)[theSimHitIndex]);
414  LogDebug("TrkHitAssocTrace") << "For cluster, simHit detId = " << (it->second)[theSimHitIndex].detUnitId()
415  << " address = (" << theSimHitCollID << ", " << theSimHitIndex
416  << "), process = " << (it->second)[theSimHitIndex].processType()
417  << " (bnch, evt, trk) = (" << (it->second)[theSimHitIndex].eventId().bunchCrossing()
418  << ", " << (it->second)[theSimHitIndex].eventId().event() << ", "
419  << (it->second)[theSimHitIndex].trackId() << ")" << std::endl;
420  }
421  }
422 }
423 
425  const DetId& detid,
426  std::vector<SimHitIdpr>& simtrackid,
427  std::vector<simhitAddr>* simhitCFPos) const {
428  uint32_t detID = detid.rawId();
429  auto isearch = stripdigisimlink->find(detID);
430  if (isearch != stripdigisimlink->end()) { //if it is not empty
431  auto link_detset = (*isearch);
432 
433  if (clust != nullptr) { //the cluster is valid
434  int clusiz = clust->amplitudes().size();
435  int first = clust->firstStrip();
436  int last = first + clusiz;
437 
438  LogDebug("TrkHitAssocDbg") << "Cluster size " << clusiz << " first strip = " << first
439  << " last strip = " << last - 1 << std::endl
440  << " detID = " << detID << " DETSET size = " << link_detset.data.size() << std::endl;
441  int channel;
442  for (const auto& linkiter : link_detset.data) {
443  channel = (int)(linkiter.channel());
444  if (channel >= first && channel < last) {
445  LogDebug("TrkHitAssocDbg") << "Channel = " << std::setw(4) << linkiter.channel()
446  << ", TrackID = " << std::setw(8) << linkiter.SimTrackId()
447  << ", tofBin = " << std::setw(3) << linkiter.TofBin()
448  << ", fraction = " << std::setw(8) << linkiter.fraction()
449  << ", Position = " << linkiter.CFposition() << std::endl;
450  SimHitIdpr currentId(linkiter.SimTrackId(), linkiter.eventId());
451 
452  //create a vector with the list of SimTrack ID's of the tracks that contributed to the RecHit
453  //write the id only once in the vector
454 
455  if (find(simtrackid.begin(), simtrackid.end(), currentId) == simtrackid.end()) {
456  LogDebug("TrkHitAssocDbg") << " Adding track id = " << currentId.first
457  << " Event id = " << currentId.second.event()
458  << " Bunch Xing = " << currentId.second.bunchCrossing() << std::endl;
459  simtrackid.push_back(currentId);
460  }
461 
462  if (simhitCFPos != nullptr) {
463  //create a vector that contains all the positions (in the MixCollection) of the SimHits that contributed to the RecHit
464  //write position only once
465  unsigned int currentCFPos = linkiter.CFposition();
466  unsigned int tofBin = linkiter.TofBin();
467  subDetTofBin theSubDetTofBin = std::make_pair(detid.subdetId(), tofBin);
468  auto it = SimHitCollMap.find(theSubDetTofBin);
469  if (it != SimHitCollMap.end()) {
470  simhitAddr currentAddr = std::make_pair(it->second, currentCFPos);
471  if (find(simhitCFPos->begin(), simhitCFPos->end(), currentAddr) == simhitCFPos->end()) {
472  simhitCFPos->push_back(currentAddr);
473  }
474  }
475  }
476  }
477  }
478  } else {
479  edm::LogError("TrackerHitAssociator") << "no cluster reference attached";
480  }
481  }
482 }
483 
484 std::vector<SimHitIdpr> TrackerHitAssociator::associateMatchedRecHit(const SiStripMatchedRecHit2D* matchedrechit,
485  std::vector<simhitAddr>* simhitCFPos) const {
486  std::vector<SimHitIdpr> matched_mono;
487  std::vector<SimHitIdpr> matched_st;
488 
489  const SiStripRecHit2D mono = matchedrechit->monoHit();
490  const SiStripRecHit2D st = matchedrechit->stereoHit();
491  //associate the two simple hits separately
492  associateSiStripRecHit(&mono, matched_mono, simhitCFPos);
493  associateSiStripRecHit(&st, matched_st, simhitCFPos);
494 
495  //save in a vector all the simtrack-id's that are common to mono and stereo hits
496  std::vector<SimHitIdpr> simtrackid;
497  if (!(matched_mono.empty() || matched_st.empty())) {
498  for (auto const& mhit : matched_mono) {
499  //save only once the ID
500  if (find(simtrackid.begin(), simtrackid.end(), mhit) == simtrackid.end()) {
501  //save if the stereoID matched the monoID
502  if (find(matched_st.begin(), matched_st.end(), mhit) != matched_st.end()) {
503  simtrackid.push_back(mhit);
504  }
505  }
506  }
507  }
508  return simtrackid;
509 }
510 
511 std::vector<SimHitIdpr> TrackerHitAssociator::associateProjectedRecHit(const ProjectedSiStripRecHit2D* projectedrechit,
512  std::vector<simhitAddr>* simhitCFPos) const {
513  //projectedRecHit is a "matched" rechit with only one component
514 
515  std::vector<SimHitIdpr> matched_mono;
516 
517  const SiStripRecHit2D mono = projectedrechit->originalHit();
518  associateSiStripRecHit(&mono, matched_mono, simhitCFPos);
519  return matched_mono;
520 }
521 
523  std::vector<SimHitIdpr>& simtrackid,
524  std::vector<simhitAddr>* simhitCFPos) const {
525  //
526  // Phase 2 outer tracker associator
527  //
528  DetId detid = rechit->geographicalId();
529  uint32_t detID = detid.rawId();
530 
531  auto isearch = ph2trackerdigisimlink->find(detID);
532  if (isearch != ph2trackerdigisimlink->end()) { //if it is not empty
533  auto link_detset = (*isearch);
534  Phase2TrackerRecHit1D::ClusterRef const& cluster = rechit->cluster();
535 
536  //check the reference is valid
537  if (!(cluster.isNull())) { //if the cluster is valid
538  int minRow = (*cluster).firstStrip();
539  int maxRow = (*cluster).firstStrip() + (*cluster).size();
540  int Col = (*cluster).column();
541  LogDebug("TrkHitAssocDbg") << " Cluster minRow " << minRow << " maxRow " << maxRow << " column " << Col
542  << std::endl;
543  int dsl = 0;
544  for (auto const& linkiter : link_detset.data) {
545  ++dsl;
546  std::pair<int, int> coord = Phase2TrackerDigi::channelToPixel(linkiter.channel());
547  LogDebug("TrkHitAssocDbg") << " " << dsl << ") Digi link: row " << coord.first << " col " << coord.second
548  << std::endl;
549  if (coord.first <= maxRow && coord.first >= minRow && coord.second == Col) {
550  LogDebug("TrkHitAssocDbg") << " !-> trackid " << linkiter.SimTrackId() << endl
551  << " fraction " << linkiter.fraction() << endl;
552  SimHitIdpr currentId(linkiter.SimTrackId(), linkiter.eventId());
553  if (find(simtrackid.begin(), simtrackid.end(), currentId) == simtrackid.end()) {
554  simtrackid.push_back(currentId);
555  }
556 
557  if (simhitCFPos != nullptr) {
558  //create a vector that contains all the positions (in the MixCollection) of the SimHits that contributed to the RecHit
559  //write position only once
560  unsigned int currentCFPos = linkiter.CFposition();
561  unsigned int tofBin = linkiter.TofBin();
562  subDetTofBin theSubDetTofBin = std::make_pair(detid.subdetId(), tofBin);
563  auto it = SimHitCollMap.find(theSubDetTofBin);
564  if (it != SimHitCollMap.end()) {
565  simhitAddr currentAddr = std::make_pair(it->second, currentCFPos);
566  if (find(simhitCFPos->begin(), simhitCFPos->end(), currentAddr) == simhitCFPos->end()) {
567  simhitCFPos->push_back(currentAddr);
568  }
569  }
570  }
571  }
572  } // end of simlink loop
573  } else {
574  edm::LogError("TrackerHitAssociator") << "no Phase2 outer tracker cluster reference attached";
575  }
576  }
577 }
578 
580  std::vector<SimHitIdpr>& simtrackid,
581  std::vector<simhitAddr>* simhitCFPos) const {
582  //
583  // Pixel associator
584  //
585  DetId detid = pixelrechit->geographicalId();
586  uint32_t detID = detid.rawId();
587 
588  auto isearch = pixeldigisimlink->find(detID);
589  if (isearch != pixeldigisimlink->end()) { //if it is not empty
590  auto link_detset = (*isearch);
591  SiPixelRecHit::ClusterRef const& cluster = pixelrechit->cluster();
592 
593  //check the reference is valid
594 
595  if (!(cluster.isNull())) { //if the cluster is valid
596 
597  int minPixelRow = (*cluster).minPixelRow();
598  int maxPixelRow = (*cluster).maxPixelRow();
599  int minPixelCol = (*cluster).minPixelCol();
600  int maxPixelCol = (*cluster).maxPixelCol();
601  LogDebug("TrkHitAssocDbg") << " Cluster minRow " << minPixelRow << " maxRow " << maxPixelRow << std::endl
602  << " Cluster minCol " << minPixelCol << " maxCol " << maxPixelCol << std::endl;
603  int dsl = 0;
604  for (auto const& linkiter : link_detset.data) {
605  ++dsl;
606  std::pair<int, int> pixel_coord = PixelDigi::channelToPixel(linkiter.channel());
607  LogDebug("TrkHitAssocDbg") << " " << dsl << ") Digi link: row " << pixel_coord.first << " col "
608  << pixel_coord.second << std::endl;
609  if (pixel_coord.first <= maxPixelRow && pixel_coord.first >= minPixelRow && pixel_coord.second <= maxPixelCol &&
610  pixel_coord.second >= minPixelCol) {
611  LogDebug("TrkHitAssocDbg") << " !-> trackid " << linkiter.SimTrackId() << endl
612  << " fraction " << linkiter.fraction() << endl;
613  SimHitIdpr currentId(linkiter.SimTrackId(), linkiter.eventId());
614  if (find(simtrackid.begin(), simtrackid.end(), currentId) == simtrackid.end()) {
615  simtrackid.push_back(currentId);
616  }
617 
618  if (simhitCFPos != nullptr) {
619  //create a vector that contains all the positions (in the MixCollection) of the SimHits that contributed to the RecHit
620  //write position only once
621  unsigned int currentCFPos = linkiter.CFposition();
622  unsigned int tofBin = linkiter.TofBin();
623  subDetTofBin theSubDetTofBin = std::make_pair(detid.subdetId(), tofBin);
624  auto it = SimHitCollMap.find(theSubDetTofBin);
625  if (it != SimHitCollMap.end()) {
626  simhitAddr currentAddr = std::make_pair(it->second, currentCFPos);
627  if (find(simhitCFPos->begin(), simhitCFPos->end(), currentAddr) == simhitCFPos->end()) {
628  simhitCFPos->push_back(currentAddr);
629  }
630  }
631  }
632  }
633  }
634  } else {
635  edm::LogError("TrackerHitAssociator") << "no Pixel cluster reference attached";
636  }
637  }
638 }
639 
640 std::vector<PSimHit> TrackerHitAssociator::associateMultiRecHit(const SiTrackerMultiRecHit* multirechit) const {
641  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
642  // std::vector<PSimHit> assimhits;
643  int size = multirechit->weights().size(), idmostprobable = 0;
644 
645  for (int i = 0; i < size; ++i) {
646  if (multirechit->weight(i) > multirechit->weight(idmostprobable))
647  idmostprobable = i;
648  }
649 
650  return associateHit(*componenthits[idmostprobable]);
651 }
652 
653 std::vector<SimHitIdpr> TrackerHitAssociator::associateMultiRecHitId(const SiTrackerMultiRecHit* multirechit,
654  std::vector<simhitAddr>* simhitCFPos) const {
655  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
656  int size = multirechit->weights().size(), idmostprobable = 0;
657 
658  for (int i = 0; i < size; ++i) {
659  if (multirechit->weight(i) > multirechit->weight(idmostprobable))
660  idmostprobable = i;
661  }
662 
663  std::vector<SimHitIdpr> simhitid;
664  associateHitId(*componenthits[idmostprobable], simhitid, simhitCFPos);
665  return simhitid;
666 }
667 
668 // fastsim
669 std::vector<SimHitIdpr> TrackerHitAssociator::associateFastRecHit(const FastTrackerRecHit* rechit) const {
670  vector<SimHitIdpr> simtrackid;
671  simtrackid.clear();
672  for (size_t index = 0, indexEnd = rechit->nSimTrackIds(); index < indexEnd; ++index) {
673  SimHitIdpr currentId(rechit->simTrackId(index), EncodedEventId(rechit->simTrackEventId(index)));
674  simtrackid.push_back(currentId);
675  }
676  return simtrackid;
677 }
678 
680  const uint32_t& detID,
681  std::vector<SimHitIdpr>& simtrackid) const {
682  std::stringstream message;
683  message << "recHit subdet, detID = " << detid.subdetId() << ", " << detID << ", (bnch, evt, trk) = ";
684  for (size_t i = 0; i < simtrackid.size(); ++i)
685  message << ", (" << simtrackid[i].second.bunchCrossing() << ", " << simtrackid[i].second.event() << ", "
686  << simtrackid[i].first << ")";
687  // message << std::endl;
688  return message.str();
689 }
ClusterRef cluster() const
Definition: SiPixelRecHit.h:47
size
Write out results.
std::string printDetBnchEvtTrk(const DetId &detid, const uint32_t &detID, std::vector< SimHitIdpr > &simtrackid) const
edm::Handle< edm::DetSetVector< PixelDigiSimLink > > pixeldigisimlink
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
uint16_t firstStrip() const
SiStripRecHit2D stereoHit() const
iterator find(det_id_type id)
Definition: DetSetVector.h:255
int event() const
get the contents of the subdetector field (should be protected?)
std::vector< float > const & weights() const
static void fillPSetDescription(edm::ParameterSetDescription &descriptions)
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:69
unsigned int detUnitId() const
Definition: PSimHit.h:99
std::vector< SimHitIdpr > associateFastRecHit(const FastTrackerRecHit *rechit) const
std::vector< edm::EDGetTokenT< CrossingFrame< PSimHit > > > cfTokens_
edm::Handle< edm::DetSetVector< StripDigiSimLink > > stripdigisimlink
ClusterRef cluster() const
void associatePixelRecHit(const SiPixelRecHit *pixelrechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=nullptr) const
std::pair< unsigned int, unsigned int > simhitAddr
void reserve(size_t s)
Definition: DetSetVector.h:145
edm::EDGetTokenT< edm::DetSetVector< PixelDigiSimLink > > ph2OTrToken_
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
void associateSiStripRecHit(const T *simplerechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=nullptr) const
T const * product() const
Definition: Handle.h:70
std::vector< SimHitIdpr > associateMatchedRecHit(const SiStripMatchedRecHit2D *matchedrechit, std::vector< simhitAddr > *simhitCFPos=nullptr) const
void associateCluster(const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< PSimHit > &simhit) const
void labelsForToken(EDGetToken const &iToken, ProductLabels &oLabels) const
Definition: Event.h:252
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:526
Definition: config.py:1
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::pair< unsigned int, unsigned int > subDetTofBin
edm::EDGetTokenT< edm::DetSetVector< PixelDigiSimLink > > pixelToken_
U second(std::pair< T, U > const &p)
void makeMaps(const edm::Event &theEvent, const Config &config)
std::vector< PSimHit > associateMultiRecHit(const SiTrackerMultiRecHit *multirechit) const
std::vector< PSimHit > associateHit(const TrackingRecHit &thit) const
SiStripCluster const & amplitudes() const
auto size() const
int bunchCrossing() const
get the detector field from this detid
unsigned short processType() const
Definition: PSimHit.h:131
virtual int32_t simTrackId(size_t i) const
std::vector< SimHitIdpr > associateMultiRecHitId(const SiTrackerMultiRecHit *multirechit, std::vector< simhitAddr > *simhitCFPos=nullptr) const
void associateSimpleRecHitCluster(const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=nullptr) const
unsigned int trackId() const
Definition: PSimHit.h:108
virtual size_t nSimTrackIds() const
EncodedEventId eventId() const
Definition: PSimHit.h:119
bool isFast(TrackingRecHit const &hit)
edm::Handle< edm::DetSetVector< PixelDigiSimLink > > ph2trackerdigisimlink
bool isNull() const
Checks for null.
Definition: Ref.h:229
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:316
edm::EDGetTokenT< edm::DetSetVector< StripDigiSimLink > > stripToken_
virtual int32_t simTrackEventId(size_t i) const
std::vector< edm::EDGetTokenT< std::vector< PSimHit > > > simHitTokens_
simhit_collectionMap SimHitCollMap
Definition: DetId.h:17
std::pair< uint32_t, EncodedEventId > SimHitIdpr
DetId geographicalId() const
TrackerHitAssociator(const edm::Event &e, const Config &config)
float weight(unsigned int i) const
static std::pair< unsigned int, unsigned int > channelToPixel(unsigned int ch)
HLT enums.
SiStripRecHit2D originalHit() const
std::vector< SimHitIdpr > associateProjectedRecHit(const ProjectedSiStripRecHit2D *projectedrechit, std::vector< simhitAddr > *simhitCFPos=nullptr) const
SiStripRecHit2D monoHit() const
void associatePhase2TrackerRecHit(const Phase2TrackerRecHit1D *rechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=nullptr) const
long double T
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
Our base class.
Definition: SiPixelRecHit.h:23
#define LogDebug(id)