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