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