CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
25 //
27  doPixel_( true ),
28  doStrip_( true ),
29  doTrackAssoc_( false ),
30  assocHitbySimTrack_( false ) {
31  //
32  // Take by default all tracker SimHits
33  //
34  vstring trackerContainers;
35  trackerContainers.push_back("g4SimHitsTrackerHitsTIBLowTof");
36  trackerContainers.push_back("g4SimHitsTrackerHitsTIBHighTof");
37  trackerContainers.push_back("g4SimHitsTrackerHitsTIDLowTof");
38  trackerContainers.push_back("g4SimHitsTrackerHitsTIDHighTof");
39  trackerContainers.push_back("g4SimHitsTrackerHitsTOBLowTof");
40  trackerContainers.push_back("g4SimHitsTrackerHitsTOBHighTof");
41  trackerContainers.push_back("g4SimHitsTrackerHitsTECLowTof");
42  trackerContainers.push_back("g4SimHitsTrackerHitsTECHighTof");
43  trackerContainers.push_back("g4SimHitsTrackerHitsPixelBarrelLowTof");
44  trackerContainers.push_back("g4SimHitsTrackerHitsPixelBarrelHighTof");
45  trackerContainers.push_back("g4SimHitsTrackerHitsPixelEndcapLowTof");
46  trackerContainers.push_back("g4SimHitsTrackerHitsPixelEndcapHighTof");
47 
48  makeMaps(e, trackerContainers);
49 
50  if(doStrip_) e.getByLabel("simSiStripDigis", stripdigisimlink);
51  if(doPixel_) e.getByLabel("simSiPixelDigis", pixeldigisimlink);
52 }
53 
55  doPixel_( conf.getParameter<bool>("associatePixel") ),
56  doStrip_( conf.getParameter<bool>("associateStrip") ),
57  doTrackAssoc_( conf.getParameter<bool>("associateRecoTracks") ) {
58 
59  assocHitbySimTrack_ = conf.existsAs<bool>("associateHitbySimTrack") ? conf.getParameter<bool>("associateHitbySimTrack") : false;
60 
61  if(doStrip_) iC.consumes<edm::DetSetVector<StripDigiSimLink> >(edm::InputTag("simSiStripDigis"));
62  if(doPixel_) iC.consumes<edm::DetSetVector<PixelDigiSimLink> >(edm::InputTag("simSiPixelDigis"));
63  }
64 
65 //
66 // Constructor with configurables
67 //
69  doPixel_( conf.getParameter<bool>("associatePixel") ),
70  doStrip_( conf.getParameter<bool>("associateStrip") ),
71  doTrackAssoc_( conf.getParameter<bool>("associateRecoTracks") ) {
72  assocHitbySimTrack_ = conf.existsAs<bool>("associateHitbySimTrack") ? conf.getParameter<bool>("associateHitbySimTrack") : false;
73 
74  //if track association there is no need to access the input collections
75  if(!doTrackAssoc_) {
76  vstring trackerContainers = conf.getParameter<std::vector<std::string> >("ROUList");
77  makeMaps(e, trackerContainers);
78  }
79 
80  if(doStrip_) e.getByLabel("simSiStripDigis", stripdigisimlink);
81  if(doPixel_) e.getByLabel("simSiPixelDigis", pixeldigisimlink);
82 }
83 
84 void TrackerHitAssociator::makeMaps(const edm::Event& theEvent, const vstring trackerContainers) {
85  // Step A: Get Inputs
86  // The collections are specified via ROUList in the configuration, and can
87  // be either crossing frames (e.g., mix/g4SimHitsTrackerHitsTIBLowTof)
88  // or just PSimHits (e.g., g4SimHits/TrackerHitsTIBLowTof)
89 
90  for(auto const& trackerContainer : trackerContainers) {
92  edm::InputTag tag_cf("mix", trackerContainer);
94  edm::InputTag tag_hits("g4SimHits", trackerContainer);
95  int Nhits = 0;
96  if (theEvent.getByLabel(tag_cf, cf_simhit)) {
97  std::auto_ptr<MixCollection<PSimHit> > thisContainerHits(new MixCollection<PSimHit>(cf_simhit.product()));
98  for (MixCollection<PSimHit>::iterator isim = thisContainerHits->begin();
99  isim != thisContainerHits->end(); isim++) {
100  DetId theDet((*isim).detUnitId());
101  if (assocHitbySimTrack_) {
102  SimHitMap[theDet].push_back((*isim));
103  } else {
104  unsigned int tofBin = StripDigiSimLink::LowTof;
105  if (trackerContainer.find(std::string("HighTof")) != std::string::npos) tofBin = StripDigiSimLink::HighTof;
106  simHitCollectionID theSimHitCollID = std::make_pair(theDet.subdetId(), tofBin);
107  SimHitCollMap[theSimHitCollID].push_back((*isim));
108  }
109  Nhits++;
110  }
111  // std::cout << "simHits from crossing frames; map size = " << SimHitCollMap.size() << ", Hit count = " << Nhits << std::endl;
112  } else {
113  theEvent.getByLabel(tag_hits, simHits);
114  for (std::vector<PSimHit>::const_iterator isim = simHits->begin();
115  isim != simHits->end(); isim++) {
116  DetId theDet((*isim).detUnitId());
117  if (assocHitbySimTrack_) {
118  SimHitMap[theDet].push_back((*isim));
119  } else {
120  unsigned int tofBin = StripDigiSimLink::LowTof;
121  if (trackerContainer.find(std::string("HighTof")) != std::string::npos) tofBin = StripDigiSimLink::HighTof;
122  simHitCollectionID theSimHitCollID = std::make_pair(theDet.subdetId(), tofBin);
123  SimHitCollMap[theSimHitCollID].push_back((*isim));
124  }
125  Nhits++;
126  }
127  // std::cout << "simHits from prompt collections; map size = " << SimHitCollMap.size() << ", Hit count = " << Nhits << std::endl;
128  }
129  }
130 }
131 
132 std::vector<PSimHit> TrackerHitAssociator::associateHit(const TrackingRecHit & thit) const
133 {
134 
135  if (const SiTrackerMultiRecHit * rechit = dynamic_cast<const SiTrackerMultiRecHit *>(&thit)){
136  return associateMultiRecHit(rechit);
137  }
138 
139  //vector with the matched SimHit
140  std::vector<PSimHit> result;
141 
142  if(doTrackAssoc_) return result;
143 
144  // Vectors to contain lists of matched simTracks, simHits
145  std::vector<SimHitIdpr> simtrackid;
146  std::vector<simhitAddr> simhitCFPos;
147 
148  //get the Detector type of the rechit
149  DetId detid= thit.geographicalId();
150  uint32_t detID = detid.rawId();
151 
152  // Get the vectors of simtrackIDs and simHit indices associated with this rechit
153  associateHitId(thit, simtrackid, &simhitCFPos);
154  // std::cout << "recHit subdet, detID = " << detid.subdetId() << ", " << detID << ", (bnch, evt, trk) = ";
155  // for (size_t i=0; i<simtrackid.size(); ++i)
156  // std::cout << ", (" << simtrackid[i].second.bunchCrossing() << ", "
157  // << simtrackid[i].second.event() << ", " << simtrackid[i].first << ")";
158  // std::cout << std::endl;
159 
160  // Get the vector of simHits associated with this rechit
161 
162  if (!assocHitbySimTrack_ && simhitCFPos.size() > 0) {
163  // We use the indices to the simHit collections taken
164  // from the DigiSimLinks and returned in simhitCFPos.
165  // simhitCFPos[i] contains the full address of the ith simhit:
166  // <collection, index> = <<subdet, tofBin>, index>
167 
168  //check if the recHit is a SiStripMatchedRecHit2D
169  bool isMatchedHit = false;
170  if(dynamic_cast<const SiStripMatchedRecHit2D *>(&thit))
171  isMatchedHit = true;
172 
173  for(size_t i=0; i<simhitCFPos.size(); i++) {
174  simhitAddr theSimHitAddr = simhitCFPos[i];
175  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
176  simhit_collectionMap::const_iterator it = SimHitCollMap.find(theSimHitCollID);
177  if (it!= SimHitCollMap.end()) {
178  unsigned int theSimHitIndex = theSimHitAddr.second;
179  if (theSimHitIndex < (it->second).size()) {
180  const PSimHit& theSimHit = (it->second)[theSimHitIndex];
181  if (isMatchedHit) {
182  // Try to remove ghosts by requiring a match to the simTrack also
183  unsigned int simHitid = theSimHit.trackId();
184  EncodedEventId simHiteid = theSimHit.eventId();
185  for(size_t i=0; i<simtrackid.size();i++) {
186  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second) {
187  result.push_back(theSimHit);
188  }
189  }
190  } else {
191  result.push_back(theSimHit);
192  }
193  // std::cout << "by CFpos, simHit detId = " << theSimHit.detUnitId() << " address = (" << (theSimHitAddr.first).first
194  // << ", " << (theSimHitAddr.first).second << ", " << theSimHitIndex
195  // << "), process = " << theSimHit.processType() << " (" << theSimHit.eventId().bunchCrossing()
196  // << ", " << theSimHit.eventId().event() << ", " << theSimHit.trackId() << ")" << std::endl;
197  }
198  }
199  }
200  return result;
201  }
202 
203  // Get the SimHit from the trackid instead
204  std::map<unsigned int, std::vector<PSimHit> >::const_iterator it = SimHitMap.find(detID);
205  if (it!= SimHitMap.end()) {
206  vector<PSimHit>::const_iterator simHitIter = (it->second).begin();
207  vector<PSimHit>::const_iterator simHitIterEnd = (it->second).end();
208  for (;simHitIter != simHitIterEnd; ++simHitIter) {
209  const PSimHit& ihit = *simHitIter;
210  unsigned int simHitid = ihit.trackId();
211  EncodedEventId simHiteid = ihit.eventId();
212  // std::cout << "by simTk, simHit, process = " << ihit.processType() << " (" << ihit.eventId().bunchCrossing()
213  // << ", " << ihit.eventId().event() << ", " << ihit.trackId() << ")";
214 
215  for(size_t i=0; i<simtrackid.size();i++) {
216  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second) {
217 // cout << "Associator ---> ID" << ihit.trackId() << " Simhit x= " << ihit.localPosition().x()
218 // << " y= " << ihit.localPosition().y() << " z= " << ihit.localPosition().x() << endl;
219  // std::cout << " matches";
220  result.push_back(ihit);
221  break;
222  }
223  }
224  // std::cout << std::endl;
225  }
226 
227  }else{
228 
230  std::map<unsigned int, std::vector<PSimHit> >::const_iterator itrphi =
231  SimHitMap.find(detID+2); //iterator to the simhit in the rphi module
232  std::map<unsigned int, std::vector<PSimHit> >::const_iterator itster =
233  SimHitMap.find(detID+1);//iterator to the simhit in the stereo module
234  if (itrphi!= SimHitMap.end()&&itster!=SimHitMap.end()) {
235  std::vector<PSimHit> simHitVector = itrphi->second;
236  simHitVector.insert(simHitVector.end(),(itster->second).begin(),(itster->second).end());
237  vector<PSimHit>::const_iterator simHitIter = simHitVector.begin();
238  vector<PSimHit>::const_iterator simHitIterEnd = simHitVector.end();
239  for (;simHitIter != simHitIterEnd; ++simHitIter) {
240  const PSimHit& ihit = *simHitIter;
241  unsigned int simHitid = ihit.trackId();
242  EncodedEventId simHiteid = ihit.eventId();
243  for(size_t i=0; i<simtrackid.size();i++) {
244  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second) {
245 // if(simHitid == simtrackid[i].first && simHiteid.bunchCrossing() == simtrackid[i].second.bunchCrossing()) {
246  // cout << "GluedDet Associator ---> ID" << ihit.trackId() << " Simhit x= " << ihit.localPosition().x()
247  // << " y= " << ihit.localPosition().y() << " z= " << ihit.localPosition().x() << endl;
248  result.push_back(ihit);
249  break;
250  }
251  }
252  }
253  }
254  }
255 
256  return result;
257 }
258 
259 std::vector< SimHitIdpr > TrackerHitAssociator::associateHitId(const TrackingRecHit & thit) const
260 {
261  std::vector< SimHitIdpr > simhitid;
262  associateHitId(thit, simhitid);
263  return simhitid;
264 }
265 
266 void TrackerHitAssociator::associateHitId(const TrackingRecHit & thit, std::vector< SimHitIdpr > & simtkid,
267  std::vector<simhitAddr>* simhitCFPos) const
268 {
269 
270  simtkid.clear();
271 
272  //get the Detector type of the rechit
273  DetId detid= thit.geographicalId();
274  if (const SiTrackerMultiRecHit * rechit = dynamic_cast<const SiTrackerMultiRecHit *>(&thit)){
275  simtkid=associateMultiRecHitId(rechit, simhitCFPos);
276  }
277 
278  //cout << "Associator ---> get Detid " << detID << endl;
279  //check we are in the strip tracker
280  if(detid.subdetId() == StripSubdetector::TIB ||
281  detid.subdetId() == StripSubdetector::TOB ||
282  detid.subdetId() == StripSubdetector::TID ||
283  detid.subdetId() == StripSubdetector::TEC)
284  {
285  //check if it is a simple SiStripRecHit2D
286  if(const SiStripRecHit2D * rechit =
287  dynamic_cast<const SiStripRecHit2D *>(&thit))
288  {
289  associateSiStripRecHit(rechit, simtkid, simhitCFPos);
290  }
291  //check if it is a SiStripRecHit1D
292  else if(const SiStripRecHit1D * rechit =
293  dynamic_cast<const SiStripRecHit1D *>(&thit))
294  {
295  associateSiStripRecHit(rechit, simtkid, simhitCFPos);
296  }
297  //check if it is a SiStripMatchedRecHit2D
298  else if(const SiStripMatchedRecHit2D * rechit =
299  dynamic_cast<const SiStripMatchedRecHit2D *>(&thit))
300  {
301  simtkid = associateMatchedRecHit(rechit, simhitCFPos);
302  }
303  //check if it is a ProjectedSiStripRecHit2D
304  else if(const ProjectedSiStripRecHit2D * rechit =
305  dynamic_cast<const ProjectedSiStripRecHit2D *>(&thit))
306  {
307  simtkid = associateProjectedRecHit(rechit, simhitCFPos);
308  }
309  else{
310  //std::cout << "associate to invalid" << std::endl;
311  //throw cms::Exception("Unknown RecHit Type") << "TrackerHitAssociator failed second casting of " << typeid(thit).name() << " type ";
312  }
313  }
314  //check we are in the pixel tracker
315  else if( (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelBarrel ||
316  (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelEndcap)
317  {
318  if(const SiPixelRecHit * rechit = dynamic_cast<const SiPixelRecHit *>(&thit))
319  {
320  associatePixelRecHit(rechit, simtkid, simhitCFPos);
321  }
322  }
323  //check if these are GSRecHits (from FastSim)
324  if(const SiTrackerGSRecHit2D * rechit = dynamic_cast<const SiTrackerGSRecHit2D *>(&thit))
325  {
326  simtkid = associateGSRecHit(rechit);
327  }
328  if(const SiTrackerGSMatchedRecHit2D * rechit = dynamic_cast<const SiTrackerGSMatchedRecHit2D *>(&thit))
329  {
330  simtkid = associateGSMatchedRecHit(rechit);
331  }
332 
333 }
334 
335 template<typename T>
336 void TrackerHitAssociator::associateSiStripRecHit(const T *simplerechit, std::vector<SimHitIdpr>& simtrackid, std::vector<simhitAddr>* simhitCFPos) const
337 {
338  const SiStripCluster* clust = &(*simplerechit->cluster());
339  associateSimpleRecHitCluster(clust, simplerechit->geographicalId(), simtrackid, simhitCFPos);
340 }
341 
342 //
343 // Method for obtaining simTracks and simHits from a cluster
344 //
346  const DetId& detid,
347  std::vector<SimHitIdpr>& simtrackid,
348  std::vector<PSimHit>& simhit) const {
349  std::vector<simhitAddr> simhitCFPos;
350  associateSimpleRecHitCluster(clust, detid, simtrackid, &simhitCFPos);
351 
352  for(size_t i=0; i<simhitCFPos.size(); i++) {
353  simhitAddr theSimHitAddr = simhitCFPos[i];
354  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
355  simhit_collectionMap::const_iterator it = SimHitCollMap.find(theSimHitCollID);
356  if (it!= SimHitCollMap.end()) {
357  unsigned int theSimHitIndex = theSimHitAddr.second;
358  if (theSimHitIndex < (it->second).size()) simhit.push_back((it->second)[theSimHitIndex]);
359  // const PSimHit& theSimHit = (it->second)[theSimHitIndex];
360  // std::cout << "For cluster, simHit detId = " << theSimHit.detUnitId() << " address = (" << (theSimHitAddr.first).first
361  // << ", " << (theSimHitAddr.first).second << ", " << theSimHitIndex
362  // << "), process = " << theSimHit.processType() << " (" << theSimHit.eventId().bunchCrossing()
363  // << ", " << theSimHit.eventId().event() << ", " << theSimHit.trackId() << ")" << std::endl;
364  }
365  }
366 }
367 
369  const DetId& detid,
370  std::vector<SimHitIdpr>& simtrackid,
371  std::vector<simhitAddr>* simhitCFPos) const {
372 
373  uint32_t detID = detid.rawId();
375  if(isearch != stripdigisimlink->end()) { //if it is not empty
376  edm::DetSet<StripDigiSimLink> link_detset = (*isearch);
377 
378  if(clust!=0){//the cluster is valid
379  int clusiz = clust->amplitudes().size();
380  int first = clust->firstStrip();
381  int last = first + clusiz;
382 
383 // std::cout << "CLUSTERSIZE " << clusiz << " first strip = " << first << " last strip = " << last-1 << std::endl;
384 // std::cout << " detID = " << detID << " DETSET size = " << link_detset.data.size() << std::endl;
385  //use a vector
386  std::vector<SimHitIdpr> idcachev;
387  std::vector<simhitAddr> CFposcachev;
388  for(edm::DetSet<StripDigiSimLink>::const_iterator linkiter = link_detset.data.begin(); linkiter != link_detset.data.end(); linkiter++){
389  if( (int)(linkiter->channel()) >= first && (int)(linkiter->channel()) < last ){
390 
391  //check this digisimlink
392 // printf("%s%4d%s%8d%s%3d%s%8.4f\n", "CHANNEL = ", linkiter->channel(), " TrackID = ", linkiter->SimTrackId(),
393 // " tofBin = ", linkiter->TofBin(), " fraction = ", linkiter->fraction());
394  /*
395  std::cout << "CHECKING CHANNEL = " << linkiter->channel() << std::endl;
396  std::cout << "TrackID = " << linkiter->SimTrackId() << std::endl;
397  std::cout << "Position = " << linkiter->CFposition() << std::endl;
398  std::cout << " fraction = " << linkiter->fraction() << std::endl;
399  */
400 
401  SimHitIdpr currentId(linkiter->SimTrackId(), linkiter->eventId());
402 
403  //create a vector with the list of SimTrack ID's of the tracks that contributed to the RecHit
404  //write the id only once in the vector
405 
406  if(find(idcachev.begin(),idcachev.end(),currentId ) == idcachev.end()){
407  /*
408  std::cout << " Adding track id = " << currentId.first
409  << " Event id = " << currentId.second.event()
410  << " Bunch Xing = " << currentId.second.bunchCrossing()
411  << std::endl;
412  */
413  idcachev.push_back(currentId);
414  simtrackid.push_back(currentId);
415  }
416 
417  if (simhitCFPos != 0) {
418  //create a vector that contains all the position (in the MixCollection) of the SimHits that contributed to the RecHit
419  //write position only once
420  unsigned int currentCFPos = linkiter->CFposition();
421  unsigned int tofBin = linkiter->TofBin();
422  simHitCollectionID theSimHitCollID = std::make_pair(detid.subdetId(), tofBin);
423  simhitAddr currentAddr = std::make_pair(theSimHitCollID, currentCFPos);
424 
425  if(find(CFposcachev.begin(), CFposcachev.end(), currentAddr ) == CFposcachev.end()) {
426 // std::cout << "CHECKING CHANNEL = " << linkiter->channel() << std::endl;
427 // std::cout << "\tTrackID = " << linkiter->SimTrackId() << "\tCFPos = " << currentCFPos <<"\ttofBin = " << tofBin << std::endl;
428 // simhit_collectionMap::const_iterator it = SimHitCollMap.find(theSimHitCollID);
429 // if (it!= SimHitCollMap.end()) {
430 // PSimHit theSimHit = it->second[currentCFPos];
431 // std::cout << "\tLocal Pos = " << theSimHit.localPosition()
432 // << "\tProcess = " << theSimHit.processType() << std::endl;
433 // }
434  CFposcachev.push_back(currentAddr);
435  simhitCFPos->push_back(currentAddr);
436  }
437  }
438  }
439  }
440  }
441  else {
442  edm::LogError("TrackerHitAssociator")<<"no cluster reference attached";
443  }
444  }
445 }
446 
447 std::vector<SimHitIdpr> TrackerHitAssociator::associateMatchedRecHit(const SiStripMatchedRecHit2D* matchedrechit, std::vector<simhitAddr>* simhitCFPos) const
448 {
449  std::vector<SimHitIdpr> matched_mono;
450  std::vector<SimHitIdpr> matched_st;
451 
452  const SiStripRecHit2D mono = matchedrechit->monoHit();
453  const SiStripRecHit2D st = matchedrechit->stereoHit();
454  //associate the two simple hits separately
455  associateSiStripRecHit(&mono, matched_mono, simhitCFPos);
456  associateSiStripRecHit(&st, matched_st, simhitCFPos);
457 
458  //save in a vector all the simtrack-id's that are common to mono and stereo hits
459  std::vector<SimHitIdpr> simtrackid;
460  if(!matched_mono.empty() && !matched_st.empty()){
461  std::vector<SimHitIdpr> idcachev;
462  for(vector<SimHitIdpr>::iterator mhit=matched_mono.begin(); mhit != matched_mono.end(); mhit++){
463  //save only once the ID
464  if(find(idcachev.begin(), idcachev.end(), (*mhit)) == idcachev.end()) {
465  idcachev.push_back(*mhit);
466  //save if the stereoID matched the monoID
467  if(find(matched_st.begin(), matched_st.end(), (*mhit)) != matched_st.end()) {
468  simtrackid.push_back(*mhit);
469  }
470  }
471  }
472  }
473 
474  return simtrackid;
475 }
476 
477 
478 std::vector<SimHitIdpr> TrackerHitAssociator::associateProjectedRecHit(const ProjectedSiStripRecHit2D * projectedrechit,
479  std::vector<simhitAddr>* simhitCFPos) const
480 {
481  //projectedRecHit is a "matched" rechit with only one component
482 
483  std::vector<SimHitIdpr> matched_mono;
484 
485  const SiStripRecHit2D mono = projectedrechit->originalHit();
486  associateSiStripRecHit(&mono, matched_mono, simhitCFPos);
487  return matched_mono;
488 }
489 
491  std::vector<SimHitIdpr> & simtrackid,
492  std::vector<simhitAddr>* simhitCFPos) const
493 {
494  //
495  // Pixel associator
496  //
497  DetId detid= pixelrechit->geographicalId();
498  uint32_t detID = detid.rawId();
499 
501  if(isearch != pixeldigisimlink->end()) { //if it is not empty
502  edm::DetSet<PixelDigiSimLink> link_detset = (*isearch);
503  SiPixelRecHit::ClusterRef const& cluster = pixelrechit->cluster();
504 
505  //check the reference is valid
506 
507  if(!(cluster.isNull())){//if the cluster is valid
508 
509  int minPixelRow = (*cluster).minPixelRow();
510  int maxPixelRow = (*cluster).maxPixelRow();
511  int minPixelCol = (*cluster).minPixelCol();
512  int maxPixelCol = (*cluster).maxPixelCol();
513  //std::cout << " Cluster minRow " << minPixelRow << " maxRow " << maxPixelRow << std::endl;
514  //std::cout << " Cluster minCol " << minPixelCol << " maxCol " << maxPixelCol << std::endl;
515  edm::DetSet<PixelDigiSimLink>::const_iterator linkiter = link_detset.data.begin();
516  int dsl = 0;
517  std::vector<SimHitIdpr> idcachev;
518  std::vector<simhitAddr> CFposcachev;
519  for( ; linkiter != link_detset.data.end(); linkiter++) {
520  dsl++;
521  std::pair<int,int> pixel_coord = PixelDigi::channelToPixel(linkiter->channel());
522  //std::cout << " " << dsl << ") Digi link: row " << pixel_coord.first << " col " << pixel_coord.second << std::endl;
523  if( pixel_coord.first <= maxPixelRow &&
524  pixel_coord.first >= minPixelRow &&
525  pixel_coord.second <= maxPixelCol &&
526  pixel_coord.second >= minPixelCol ) {
527  //std::cout << " !-> trackid " << linkiter->SimTrackId() << endl;
528  //std::cout << " fraction " << linkiter->fraction() << endl;
529  SimHitIdpr currentId(linkiter->SimTrackId(), linkiter->eventId());
530  if(find(idcachev.begin(),idcachev.end(),currentId) == idcachev.end()){
531  simtrackid.push_back(currentId);
532  idcachev.push_back(currentId);
533  }
534 
535  if (simhitCFPos != 0) {
536  //create a vector that contains all the position (in the MixCollection) of the SimHits that contributed to the RecHit
537  //write position only once
538  unsigned int currentCFPos = linkiter->CFposition();
539  unsigned int tofBin = linkiter->TofBin();
540  simHitCollectionID theSimHitCollID = std::make_pair(detid.subdetId(), tofBin);
541  simhitAddr currentAddr = std::make_pair(theSimHitCollID, currentCFPos);
542 
543  if(find(CFposcachev.begin(), CFposcachev.end(), currentAddr) == CFposcachev.end()) {
544  CFposcachev.push_back(currentAddr);
545  simhitCFPos->push_back(currentAddr);
546  }
547  }
548 
549  }
550  }
551  }
552  else{
553  edm::LogError("TrackerHitAssociator")<<"no Pixel cluster reference attached";
554 
555  }
556  }
557 }
558 
559 std::vector<SimHitIdpr> TrackerHitAssociator::associateGSRecHit(const SiTrackerGSRecHit2D * gsrechit) const
560 {
561  //GSRecHit is the FastSimulation RecHit that contains the TrackId already
562 
563  vector<SimHitIdpr> simtrackid;
564  simtrackid.clear();
565  SimHitIdpr currentId(gsrechit->simtrackId(), EncodedEventId(gsrechit->eeId()));
566  simtrackid.push_back(currentId);
567  return simtrackid;
568 }
569 
570 std::vector<PSimHit> TrackerHitAssociator::associateMultiRecHit(const SiTrackerMultiRecHit * multirechit) const{
571  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
572  // std::vector<PSimHit> assimhits;
573  int size=multirechit->weights().size(), idmostprobable=0;
574 
575  for (int i=0; i<size; i++){
576  if(multirechit->weight(i)>multirechit->weight(idmostprobable)) idmostprobable=i;
577  }
578 
579  return associateHit(*componenthits[idmostprobable]);
580 }
581 
582 std::vector<SimHitIdpr> TrackerHitAssociator::associateMultiRecHitId(const SiTrackerMultiRecHit * multirechit, std::vector<simhitAddr>* simhitCFPos) const{
583  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
584  int size=multirechit->weights().size(), idmostprobable=0;
585 
586  for (int i=0; i<size; i++){
587  if(multirechit->weight(i)>multirechit->weight(idmostprobable)) idmostprobable=i;
588  }
589 
590  std::vector< SimHitIdpr > simhitid;
591  associateHitId(*componenthits[idmostprobable], simhitid, simhitCFPos);
592  return simhitid;
593 }
594 
595 std::vector<SimHitIdpr> TrackerHitAssociator::associateGSMatchedRecHit(const SiTrackerGSMatchedRecHit2D * gsmrechit) const
596 {
597  //GSRecHit is the FastSimulation RecHit that contains the TrackId already
598 
599  vector<SimHitIdpr> simtrackid;
600  simtrackid.clear();
601  SimHitIdpr currentId(gsmrechit->simtrackId(), EncodedEventId(gsmrechit->eeId()));
602  simtrackid.push_back(currentId);
603  return simtrackid;
604 }
605 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:185
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< std::string > vstring
std::pair< simHitCollectionID, unsigned int > simhitAddr
void associatePixelRecHit(const SiPixelRecHit *pixelrechit, std::vector< SimHitIdpr > &simhitid, std::vector< simhitAddr > *simhitCFPos=0) const
uint16_t firstStrip() 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
void associateSimpleRecHitCluster(const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
float weight(unsigned int i) const
const int & simtrackId() 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)
std::vector< float > const & weights() const
tuple result
Definition: query.py:137
std::vector< SimHitIdpr > associateGSMatchedRecHit(const SiTrackerGSMatchedRecHit2D *gsmrechit) const
EncodedEventId eventId() const
Definition: PSimHit.h:105
#define end
Definition: vmac.h:37
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
bool isNull() const
Checks for null.
Definition: Ref.h:247
tuple conf
Definition: dbtoconf.py:185
TrackerHitAssociator(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
simhit_collectionMap SimHitCollMap
SiStripRecHit2D originalHit() const
Definition: DetId.h:18
void makeMaps(const edm::Event &theEvent, const vstring trackerContainers)
void associateCluster(const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< PSimHit > &simhit) const
const uint32_t & eeId() const
std::pair< uint32_t, EncodedEventId > SimHitIdpr
SiStripRecHit2D stereoHit() const
T const * product() const
Definition: Handle.h:81
const uint32_t & eeId() const
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:62
tuple simHits
Definition: trackerHits.py:16
SiStripRecHit2D monoHit() const
collection_type data
Definition: DetSet.h:78
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
#define begin
Definition: vmac.h:30
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
std::vector< SimHitIdpr > associateGSRecHit(const SiTrackerGSRecHit2D *gsrechit) const
DetId geographicalId() const
volatile std::atomic< bool > shutdown_flag false
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
long double T
std::vector< SimHitIdpr > associateMatchedRecHit(const SiStripMatchedRecHit2D *matchedrechit, std::vector< simhitAddr > *simhitCFPos=0) const
tuple size
Write out results.
void associateSiStripRecHit(const T *simplerechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
const std::vector< uint8_t > & amplitudes() const
std::pair< unsigned int, unsigned int > simHitCollectionID
Our base class.
Definition: SiPixelRecHit.h:23