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  myEvent_(e),
28  doPixel_( true ),
29  doStrip_( true ),
30  doTrackAssoc_( false ) {
31  trackerContainers.clear();
32  //
33  // Take by default all tracker SimHits
34  //
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  // Step A: Get Inputs
49  //MAKE THIS PRIVATE MEMBERS
50  // edm::Handle<CrossingFrame<PSimHit> > cf_simhit;
51  // std::vector<const CrossingFrame<PSimHit> *> cf_simhitvec;
52 
53  for(uint32_t i = 0; i< trackerContainers.size();i++){
55  cf_simhitvec.push_back(cf_simhit.product());
56  }
57 
58  std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(new MixCollection<PSimHit>(cf_simhitvec));
59  TrackerHits = (*allTrackerHits);
60 
61  //Loop on PSimHit
62  SimHitMap.clear();
63 
65  for (isim=allTrackerHits->begin(); isim!= allTrackerHits->end();isim++) {
66  SimHitMap[(*isim).detUnitId()].push_back((*isim));
67  }
68 
69  if(doStrip_) e.getByLabel("simSiStripDigis", stripdigisimlink);
70  if(doPixel_) e.getByLabel("simSiPixelDigis", pixeldigisimlink);
71 
72 }
73 
74 //
75 // Constructor with configurables
76 //
78  myEvent_(e),
79  doPixel_( conf.getParameter<bool>("associatePixel") ),
80  doStrip_( conf.getParameter<bool>("associateStrip") ),
81  doTrackAssoc_( conf.getParameter<bool>("associateRecoTracks") ){
82 
83  //if track association there is no need to acces the CrossingFrame
84  if(!doTrackAssoc_) {
85 
86  trackerContainers.clear();
87  trackerContainers = conf.getParameter<std::vector<std::string> >("ROUList");
88 
89  // Step A: Get Inputs
90  // edm::Handle<CrossingFrame<PSimHit> > cf_simhit;
91  // std::vector<const CrossingFrame<PSimHit> *> cf_simhitvec;
92  for(uint32_t i = 0; i< trackerContainers.size();i++){
94  cf_simhitvec.push_back(cf_simhit.product());
95  }
96 
97  // std::cout << "SIMHITVEC SIZE = " << cf_simhitvec.size() << std::endl;
98 
99  // TrackerHits = new MixCollection<PSimHit>(cf_simhitvec);
100  //std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(TrackerHits);
101  // std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(new MixCollection<PSimHit>(cf_simhitvec));
102  std::auto_ptr<MixCollection<PSimHit> > allTrackerHits(new MixCollection<PSimHit>(cf_simhitvec));
103  TrackerHits = (*allTrackerHits);
104 
105  //Loop on PSimHit
106  SimHitMap.clear();
107 
109  for (isim=allTrackerHits->begin(); isim!= allTrackerHits->end();isim++) {
110  SimHitMap[(*isim).detUnitId()].push_back((*isim));
111  }
112 
113  }
114 
115  if(doStrip_) e.getByLabel("simSiStripDigis", stripdigisimlink);
116  if(doPixel_) e.getByLabel("simSiPixelDigis", pixeldigisimlink);
117 
118 }
119 
120 std::vector<PSimHit> TrackerHitAssociator::associateHit(const TrackingRecHit & thit)
121 {
122 
123  //check in case of TTRH
124  if(const TransientTrackingRecHit * ttrh = dynamic_cast<const TransientTrackingRecHit *>(&thit)) {
125  std::cout << "calling associateHit for TransientTRH" << std::endl;
126  return associateHit(*ttrh->hit());
127  }
128 
129  //vector with the matched SimHit
130  std::vector<PSimHit> result;
131  // std::vector<PSimHit> result_old;
132  //initialize vectors!
133  simtrackid.clear();
134  simhitCFPos.clear();
135  StripHits = false;
136 
137  //get the Detector type of the rechit
138  DetId detid= thit.geographicalId();
139  uint32_t detID = detid.rawId();
140 
141  // cout << "Associator ---> get Detid " << detID << endl;
142  //check we are in the strip tracker
143  if(detid.subdetId() == StripSubdetector::TIB ||
144  detid.subdetId() == StripSubdetector::TOB ||
145  detid.subdetId() == StripSubdetector::TID ||
146  detid.subdetId() == StripSubdetector::TEC)
147  {
148  //check if it is a simple SiStripRecHit2D
149  if(const SiStripRecHit2D * rechit =
150  dynamic_cast<const SiStripRecHit2D *>(&thit))
151  {
152  //std::cout << "associate to hit2D" << std::endl;
154  }
155  else if(const SiStripRecHit1D * rechit =
156  dynamic_cast<const SiStripRecHit1D *>(&thit)) //check if it is a SiStripRecHit1D
157  {
158  //std::cout << "associate to hit1D" << std::endl;
160  }
161  else if(const SiStripMatchedRecHit2D * rechit =
162  dynamic_cast<const SiStripMatchedRecHit2D *>(&thit)) //check if it is a matched SiStripMatchedRecHit2D
163  {
164  //std::cout << "associate to matched" << std::endl;
166  }
167  else if(const ProjectedSiStripRecHit2D * rechit =
168  dynamic_cast<const ProjectedSiStripRecHit2D *>(&thit)) //check if it is a ProjectedSiStripRecHit2D
169  {
170  //std::cout << "associate to projectedHit" << std::endl;
172  detid = rechit->originalHit().geographicalId();
173  detID = detid.rawId();
174  }
175  else {
176  //std::cout << "associate to Invalid strip hit??" << std::endl;
177  //throw cms::Exception("Unknown RecHit Type") << "TrackerHitAssociator failed first casting of " << typeid(thit).name() << " type ";
178  }
179 
180  }
181  //check we are in the pixel tracker
182  if( (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelBarrel ||
183  (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelEndcap)
184  {
185  if(const SiPixelRecHit * rechit = dynamic_cast<const SiPixelRecHit *>(&thit))
186  {
187  //std::cout << "associate to pixelHit" << std::endl;
189  }
190  }
191  //check if these are GSRecHits (from FastSim)
192 
193  if(const SiTrackerGSRecHit2D * rechit = dynamic_cast<const SiTrackerGSRecHit2D *>(&thit))
194  {
195  simtrackid = associateGSRecHit(rechit);
196  }
197  if (const SiTrackerMultiRecHit * rechit = dynamic_cast<const SiTrackerMultiRecHit *>(&thit)){
198  return associateMultiRecHit(rechit);
199  }
200 
201  //check if these are GSMatchedRecHits (from FastSim)
202  if(const SiTrackerGSMatchedRecHit2D * rechit = dynamic_cast<const SiTrackerGSMatchedRecHit2D *>(&thit))
203  {
205  }
206 
207  //
208  //Save the SimHits in a vector. for the macthed hits both the rphi and stereo simhits are saved.
209  //
210 
211  if(StripHits){
212  //USE THIS FOR STRIPS
213  // std::cout << "NEW SIZE = " << simhitCFPos.size() << std::endl;
214 
215  for(size_t i=0; i<simhitCFPos.size(); i++){
216  //std::cout << "NEW CFPOS " << simhitCFPos[i] << endl;
217  //std::cout << "NEW LOCALPOS " << TrackerHits.getObject(simhitCFPos[i]).localPosition() << endl;
218  result.push_back( TrackerHits.getObject(simhitCFPos[i]));
219  }
220  }else {
221 
222  //now get the SimHit from the trackid
223  vector<PSimHit> simHit;
224  std::map<unsigned int, std::vector<PSimHit> >::const_iterator it = SimHitMap.find(detID);
225  simHit.clear();
226  if (it!= SimHitMap.end()){
227  simHit = it->second;
228  vector<PSimHit>::const_iterator simHitIter = simHit.begin();
229  vector<PSimHit>::const_iterator simHitIterEnd = simHit.end();
230  for (;simHitIter != simHitIterEnd; ++simHitIter) {
231  const PSimHit ihit = *simHitIter;
232  unsigned int simHitid = ihit.trackId();
233  EncodedEventId simHiteid = ihit.eventId();
234 
235  for(size_t i=0; i<simtrackid.size();i++){
236  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second){
237  // cout << "Associator ---> ID" << ihit.trackId() << " Simhit x= " << ihit.localPosition().x()
238  // << " y= " << ihit.localPosition().y() << " z= " << ihit.localPosition().x() << endl;
239  result.push_back(ihit);
240  }
241  }
242  }
243  }else{
245  std::map<unsigned int, std::vector<PSimHit> >::const_iterator itrphi =
246  SimHitMap.find(detID+2);//iterator to the simhit in the rphi module
247  std::map<unsigned int, std::vector<PSimHit> >::const_iterator itster =
248  SimHitMap.find(detID+1);//iterator to the simhit in the stereo module
249  if (itrphi!= SimHitMap.end()&&itster!=SimHitMap.end()){
250  simHit = itrphi->second;
251  simHit.insert(simHit.end(),(itster->second).begin(),(itster->second).end());
252  vector<PSimHit>::const_iterator simHitIter = simHit.begin();
253  vector<PSimHit>::const_iterator simHitIterEnd = simHit.end();
254  for (;simHitIter != simHitIterEnd; ++simHitIter) {
255  const PSimHit ihit = *simHitIter;
256  unsigned int simHitid = ihit.trackId();
257  EncodedEventId simHiteid = ihit.eventId();
258 
259  //===>>>>>change here!!!!
260  // for(size_t i=0; i<simtrackid.size();i++){
261  // cout << " GluedDet Associator --> check sihit id's = " << simHitid <<"; compared id's = "<< simtrackid[i] <<endl;
262  // if(simHitid == simtrackid[i]){ //exclude the geant particles. they all have the same id
263  for(size_t i=0; i<simtrackid.size();i++){
264  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second){
265  // cout << "GluedDet Associator ---> ID" << ihit.trackId() << " Simhit x= " << ihit.localPosition().x()
266  // << " y= " << ihit.localPosition().y() << " z= " << ihit.localPosition().x() << endl;
267  result.push_back(ihit);
268  }
269  }
270  }
271  }
272  }
273  }
274 
275 
276  return result;
277 }
278 
279 //std::vector<unsigned int> TrackerHitAssociator::associateHitId(const TrackingRecHit & thit)
280 std::vector< SimHitIdpr > TrackerHitAssociator::associateHitId(const TrackingRecHit & thit)
281 {
282  std::vector< SimHitIdpr > simhitid;
283  associateHitId(thit, simhitid);
284  return simhitid;
285 }
286 
287 void TrackerHitAssociator::associateHitId(const TrackingRecHit & thit, std::vector< SimHitIdpr > & simtkid)
288 {
289 
290  //check in case of TTRH
291  if(const TransientTrackingRecHit * ttrh = dynamic_cast<const TransientTrackingRecHit *>(&thit)) {
292  associateHitId(*ttrh->hit(), simtkid);
293  }
294  else{
295  simtkid.clear();
296  simhitCFPos.clear();
297  //vector with the associated Simtkid
298  //vector with the CF position of the associated simhits
299 
300  //get the Detector type of the rechit
301  DetId detid= thit.geographicalId();
302  //apparently not used
303  // uint32_t detID = detid.rawId();
304  if (const SiTrackerMultiRecHit * rechit = dynamic_cast<const SiTrackerMultiRecHit *>(&thit)){
305  simtkid=associateMultiRecHitId(rechit);
306  }
307 
308  //cout << "Associator ---> get Detid " << detID << endl;
309  //check we are in the strip tracker
310  if(detid.subdetId() == StripSubdetector::TIB ||
311  detid.subdetId() == StripSubdetector::TOB ||
312  detid.subdetId() == StripSubdetector::TID ||
313  detid.subdetId() == StripSubdetector::TEC)
314  {
315  //check if it is a simple SiStripRecHit2D
316  if(const SiStripRecHit2D * rechit =
317  dynamic_cast<const SiStripRecHit2D *>(&thit))
318  {
319  associateSimpleRecHit(rechit, simtkid);
320  }
321  //check if it is a matched SiStripMatchedRecHit2D
322  else if(const SiStripRecHit1D * rechit =
323  dynamic_cast<const SiStripRecHit1D *>(&thit))
324  {
325  associateSiStripRecHit1D(rechit,simtkid);
326  }
327  //check if it is a matched SiStripMatchedRecHit2D
328  else if(const SiStripMatchedRecHit2D * rechit =
329  dynamic_cast<const SiStripMatchedRecHit2D *>(&thit))
330  {
331  simtkid = associateMatchedRecHit(rechit);
332  }
333  //check if it is a ProjectedSiStripRecHit2D
334  else if(const ProjectedSiStripRecHit2D * rechit =
335  dynamic_cast<const ProjectedSiStripRecHit2D *>(&thit))
336  {
337  simtkid = associateProjectedRecHit(rechit);
338  }
339  else{
340  //std::cout << "associate to invalid" << std::endl;
341  //throw cms::Exception("Unknown RecHit Type") << "TrackerHitAssociator failed second casting of " << typeid(thit).name() << " type ";
342  }
343  }
344  //check we are in the pixel tracker
345  else if( (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelBarrel ||
346  (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelEndcap)
347  {
348  if(const SiPixelRecHit * rechit = dynamic_cast<const SiPixelRecHit *>(&thit))
349  {
350  associatePixelRecHit(rechit,simtkid );
351  }
352  }
353  //check if these are GSRecHits (from FastSim)
354  if(const SiTrackerGSRecHit2D * rechit = dynamic_cast<const SiTrackerGSRecHit2D *>(&thit))
355  {
356  simtkid = associateGSRecHit(rechit);
357  }
358  if(const SiTrackerGSMatchedRecHit2D * rechit = dynamic_cast<const SiTrackerGSMatchedRecHit2D *>(&thit))
359  {
360  simtkid = associateGSMatchedRecHit(rechit);
361  }
362 
363  }
364 }
365 
366 
367 void TrackerHitAssociator::associateSimpleRecHit(const SiStripRecHit2D * simplerechit, std::vector<SimHitIdpr>& simtrackid)
368 {
369  const SiStripCluster* clust = 0;
370  if(simplerechit->cluster().isNonnull())
371  {
372  clust=&(*simplerechit->cluster());
373  }
374  else if(simplerechit->cluster_regional().isNonnull())
375  {
376  clust=&(*simplerechit->cluster_regional());
377  }
378 
380 }
381 
382 
383 //This function could be templated to avoid to repeat the same code twice??
384 void TrackerHitAssociator::associateSiStripRecHit1D(const SiStripRecHit1D * simplerechit, std::vector<SimHitIdpr>& simtrackid)
385 {
386  const SiStripCluster* clust = 0;
387  if(simplerechit->cluster().isNonnull())
388  {
389  clust=&(*simplerechit->cluster());
390  }
391  else if(simplerechit->cluster_regional().isNonnull())
392  {
393  clust=&(*simplerechit->cluster_regional());
394  }
396 }
397 
399  const uint32_t& detID,
400  std::vector<SimHitIdpr>& theSimtrackid, std::vector<PSimHit>& simhit)
401 {
402 // Caller needs to clear simhit before calling this function
403 
404  //initialize vector
405  theSimtrackid.clear();
406  //initialize class vector
407  simhitCFPos.clear();
408 
409  associateSimpleRecHitCluster(clust, detID, theSimtrackid);
410 
411  for(size_t i=0; i<simhitCFPos.size(); i++){
412  simhit.push_back(TrackerHits.getObject(simhitCFPos[i]));
413  }
414 }
415 
417  const uint32_t& detID,
418  std::vector<PSimHit>& simhit)
419 {
420 // Caller needs to clear simhit before calling this function
421 
422  //initialize class vectors
423  simtrackid.clear();
424  simhitCFPos.clear();
425 
427 
428  for(size_t i=0; i<simhitCFPos.size(); i++){
429  simhit.push_back(TrackerHits.getObject(simhitCFPos[i]));
430  }
431 }
432 
434  const uint32_t& detID,
435  std::vector<SimHitIdpr>& simtrackid){
436  // std::cout <<"ASSOCIATE SIMPLE RECHIT" << std::endl;
437  StripHits =true;
438 
439  //DetId detid= clust->geographicalId();
440  // uint32_t detID = detid.rawId();
441 
442  //to store temporary charge information
443  std::vector<SimHitIdpr> cache_simtrackid;
444  cache_simtrackid.clear();
445 
446  std::map<SimHitIdpr, vector<float> > temp_simtrackid;
447  temp_simtrackid.clear();
448 
450  if(isearch != stripdigisimlink->end()) { //if it is not empty
451  //link_detset is a structure, link_detset.data is a std::vector<StripDigiSimLink>
452  //edm::DetSet<StripDigiSimLink> link_detset = (*stripdigisimlink)[detID];
453  edm::DetSet<StripDigiSimLink> link_detset = (*isearch);
454 
455  if(clust!=0){//the cluster is valid
456 
457  // const edm::Ref<edm::DetSetVector<SiStripCluster>, SiStripCluster, edm::refhelper::FindForDetSetVector<SiStripCluster> > clust=simplerechit->cluster();
458 
459  //float chg;
460  int clusiz = clust->amplitudes().size();
461  int first = clust->firstStrip();
462  int last = first + clusiz;
463 
464 // std::cout << "CLUSTERSIZE " << clusiz << " first strip = " << first << " last strip = " << last-1 << std::endl;
465 // std::cout << " detID = " << detID << " DETSET size = " << link_detset.data.size() << std::endl;
466  //use a vector
467  std::vector<SimHitIdpr> idcachev;
468  std::vector<int> CFposcachev;
469  for(edm::DetSet<StripDigiSimLink>::const_iterator linkiter = link_detset.data.begin(); linkiter != link_detset.data.end(); linkiter++){
470  //StripDigiSimLink link = *linkiter;
471 
472  if( (int)(linkiter->channel()) >= first && (int)(linkiter->channel()) < last ){
473 
474  //check this digisimlink
475 // printf("%s%4d%s%8d%s%3d%s%8.4f\n", "CHANNEL = ", linkiter->channel(), " TrackID = ", linkiter->SimTrackId(),
476 // " Process = ", TrackerHits.getObject(linkiter->CFposition()-1).processType(), " fraction = ", linkiter->fraction());
477  /*
478  std::cout << "CHECKING CHANNEL = " << linkiter->channel() << std::endl;
479  std::cout << "TrackID = " << linkiter->SimTrackId() << std::endl;
480  std::cout << "Position = " << linkiter->CFposition() << std::endl;
481  std::cout << " POS -1 = " << TrackerHits.getObject(linkiter->CFposition()-1).localPosition() << std::endl;
482  std::cout << " Process = " << TrackerHits.getObject(linkiter->CFposition()-1).processType() << std::endl;
483  std::cout << " fraction = " << linkiter->fraction() << std::endl;
484  */
485 
486  SimHitIdpr currentId(linkiter->SimTrackId(), linkiter->eventId());
487 
488  //create a vector with the list of SimTrack ID's of the tracks that contributed to the RecHit
489  //write the id only once in the vector
490 
491  if(find(idcachev.begin(),idcachev.end(),currentId ) == idcachev.end()){
492  /*
493  std::cout << " Adding track id = " << currentId.first
494  << " Event id = " << currentId.second.event()
495  << " Bunch Xing = " << currentId.second.bunchCrossing()
496  << std::endl;
497  */
498  idcachev.push_back(currentId);
499  simtrackid.push_back(currentId);
500  }
501 
502  //create a vector that contains all the position (in the MixCollection) of the SimHits that contributed to the RecHit
503  //write position only once
504  int currentCFPos = linkiter->CFposition()-1;
505  if(find(CFposcachev.begin(),CFposcachev.end(),currentCFPos ) == CFposcachev.end()){
506  /*
507  std::cout << "CHECKING CHANNEL = " << linkiter->channel() << std::endl;
508  std::cout << "\tTrackID = " << linkiter->SimTrackId() << "\tCFPos = " << currentCFPos << std::endl;
509  std::cout << "\tLocal Pos = " << TrackerHits.getObject(currentCFPos).localPosition()
510  << "\tProcess = " << TrackerHits.getObject(currentCFPos).processType() << std::endl;
511  */
512  CFposcachev.push_back(currentCFPos);
513  simhitCFPos.push_back(currentCFPos);
514  // simhitassoc.push_back( TrackerHits.getObject(currentCFPos));
515  }
516  }
517  }
518  }
519  else {
520  edm::LogError("TrackerHitAssociator")<<"no cluster reference attached";
521  }
522  }
523 }
524 
525 std::vector<SimHitIdpr> TrackerHitAssociator::associateMatchedRecHit(const SiStripMatchedRecHit2D * matchedrechit)
526 {
527 
528  StripHits = true;
529 
530 
531  vector<SimHitIdpr> matched_mono;
532  vector<SimHitIdpr> matched_st;
533  matched_mono.clear();
534  matched_st.clear();
535 
536  const SiStripRecHit2D mono = matchedrechit->monoHit();
537  const SiStripRecHit2D st = matchedrechit->stereoHit();
538  //associate the two simple hits separately
539  associateSimpleRecHit(&mono,matched_mono );
540  associateSimpleRecHit(&st, matched_st );
541 
542  //save in a vector all the simtrack-id's that are common to mono and stereo hits
543  if(!matched_mono.empty() && !matched_st.empty()){
544  simtrackid.clear(); //final result vector
545  // std::vector<unsigned int> idcachev;
546  std::vector<SimHitIdpr> idcachev;
547  //for(vector<unsigned int>::iterator mhit=matched_mono.begin(); mhit != matched_mono.end(); mhit++){
548  for(vector<SimHitIdpr>::iterator mhit=matched_mono.begin(); mhit != matched_mono.end(); mhit++){
549  //save only once the ID
550  if(find(idcachev.begin(), idcachev.end(),(*mhit)) == idcachev.end()) {
551  idcachev.push_back(*mhit);
552  //save if the stereoID matched the monoID
553  if(find(matched_st.begin(), matched_st.end(),(*mhit))!=matched_st.end()){
554  simtrackid.push_back(*mhit);
555  //std::cout << "matched case: saved ID " << (*mhit) << std::endl;
556  }
557  }
558  }
559  }
560  return simtrackid;
561 }
562 
563 
564 std::vector<SimHitIdpr> TrackerHitAssociator::associateProjectedRecHit(const ProjectedSiStripRecHit2D * projectedrechit)
565 {
566  StripHits = true;
567 
568 
569  //projectedRecHit is a "matched" rechit with only one component
570 
571  vector<SimHitIdpr> matched_mono;
572  matched_mono.clear();
573 
574  const SiStripRecHit2D mono = projectedrechit->originalHit();
575  associateSimpleRecHit(&mono, matched_mono);
576  return matched_mono;
577 }
578 
579 //std::vector<unsigned int> TrackerHitAssociator::associatePixelRecHit(const SiPixelRecHit * pixelrechit)
580 void TrackerHitAssociator::associatePixelRecHit(const SiPixelRecHit * pixelrechit, std::vector<SimHitIdpr> & simtrackid)
581 {
582  StripHits = false;
583 
584  //
585  // Pixel associator
586  //
587  DetId detid= pixelrechit->geographicalId();
588  uint32_t detID = detid.rawId();
589 
591  if(isearch != pixeldigisimlink->end()) { //if it is not empty
592  // edm::DetSet<PixelDigiSimLink> link_detset = (*pixeldigisimlink)[detID];
593  edm::DetSet<PixelDigiSimLink> link_detset = (*isearch);
594  // edm::Ref< edm::DetSetVector<SiPixelCluster>, SiPixelCluster> const& cluster = pixelrechit->cluster();
595  SiPixelRecHit::ClusterRef const& cluster = pixelrechit->cluster();
596 
597  //check the reference is valid
598 
599  if(!(cluster.isNull())){//if the cluster is valid
600 
601  int minPixelRow = (*cluster).minPixelRow();
602  int maxPixelRow = (*cluster).maxPixelRow();
603  int minPixelCol = (*cluster).minPixelCol();
604  int maxPixelCol = (*cluster).maxPixelCol();
605  //std::cout << " Cluster minRow " << minPixelRow << " maxRow " << maxPixelRow << std::endl;
606  //std::cout << " Cluster minCol " << minPixelCol << " maxCol " << maxPixelCol << std::endl;
607  edm::DetSet<PixelDigiSimLink>::const_iterator linkiter = link_detset.data.begin();
608  int dsl = 0;
609  // std::vector<unsigned int> idcachev;
610  std::vector<SimHitIdpr> idcachev;
611  for( ; linkiter != link_detset.data.end(); linkiter++) {
612  dsl++;
613  std::pair<int,int> pixel_coord = PixelDigi::channelToPixel(linkiter->channel());
614  //std::cout << " " << dsl << ") Digi link: row " << pixel_coord.first << " col " << pixel_coord.second << std::endl;
615  if( pixel_coord.first <= maxPixelRow &&
616  pixel_coord.first >= minPixelRow &&
617  pixel_coord.second <= maxPixelCol &&
618  pixel_coord.second >= minPixelCol ) {
619  //std::cout << " !-> trackid " << linkiter->SimTrackId() << endl;
620  //std::cout << " fraction " << linkiter->fraction() << endl;
621  SimHitIdpr currentId(linkiter->SimTrackId(), linkiter->eventId());
622  // if(find(idcachev.begin(),idcachev.end(),linkiter->SimTrackId()) == idcachev.end()){
623  if(find(idcachev.begin(),idcachev.end(),currentId) == idcachev.end()){
624  // simtrackid.push_back(linkiter->SimTrackId());
625  //idcachev.push_back(linkiter->SimTrackId());
626  simtrackid.push_back(currentId);
627  idcachev.push_back(currentId);
628  }
629  }
630  }
631  }
632  else{
633  edm::LogError("TrackerHitAssociator")<<"no Pixel cluster reference attached";
634 
635  }
636  }
637 
638 
639 }
640 
641 std::vector<SimHitIdpr> TrackerHitAssociator::associateGSRecHit(const SiTrackerGSRecHit2D * gsrechit)
642 {
643  StripHits = false;
644  //GSRecHit is the FastSimulation RecHit that contains the TrackId already
645 
646  vector<SimHitIdpr> simtrackid;
647  simtrackid.clear();
648  SimHitIdpr currentId(gsrechit->simtrackId(), EncodedEventId(gsrechit->eeId()));
649  simtrackid.push_back(currentId);
650  return simtrackid;
651 }
652 
653 std::vector<PSimHit> TrackerHitAssociator::associateMultiRecHit(const SiTrackerMultiRecHit * multirechit){
654  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
655  // std::vector<PSimHit> assimhits;
656  int size=multirechit->weights().size(), idmostprobable=0;
657 
658  for (int i=0; i<size; i++){
659  if(multirechit->weight(i)>multirechit->weight(idmostprobable)) idmostprobable=i;
660  }
661 
662  // std::vector<PSimHit> assimhits = associateHit(**mostpropable);
663  // assimhits.insert(assimhits.end(), asstocurrent.begin(), asstocurrent.end());
664  //std::cout << "Returning " << assimhits.size() << " simhits" << std::endl;
665  return associateHit(*componenthits[idmostprobable]);
666 }
667 
668 std::vector<SimHitIdpr> TrackerHitAssociator::associateMultiRecHitId(const SiTrackerMultiRecHit * multirechit){
669  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
670  int size=multirechit->weights().size(), idmostprobable=0;
671 
672  for (int i=0; i<size; i++){
673  if(multirechit->weight(i)>multirechit->weight(idmostprobable)) idmostprobable=i;
674  }
675 
676  return associateHitId(*componenthits[idmostprobable]);
677 }
678 
680 {
681  StripHits = false;
682  //GSRecHit is the FastSimulation RecHit that contains the TrackId already
683 
684  vector<SimHitIdpr> simtrackid;
685  simtrackid.clear();
686  SimHitIdpr currentId(gsmrechit->simtrackId(), EncodedEventId(gsmrechit->eeId()));
687  simtrackid.push_back(currentId);
688  return simtrackid;
689 }
690 
ClusterRef cluster() const
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
MixCollection< PSimHit > TrackerHits
edm::Handle< edm::DetSetVector< StripDigiSimLink > > stripdigisimlink
std::vector< SimHitIdpr > associateGSRecHit(const SiTrackerGSRecHit2D *gsrechit)
ClusterRegionalRef cluster_regional() const
TrackerHitAssociator(const edm::Event &e)
void associateSimpleRecHit(const SiStripRecHit2D *simplerechit, std::vector< SimHitIdpr > &simhitid)
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
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit)
void associateSiStripRecHit1D(const SiStripRecHit1D *simplerechit, std::vector< SimHitIdpr > &simhitid)
std::vector< SimHitIdpr > simtrackid
float weight(unsigned int i) const
const int & simtrackId() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
U second(std::pair< T, U > const &p)
iterator end()
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
std::vector< const CrossingFrame< PSimHit > * > cf_simhitvec
bool isNull() const
Checks for null.
Definition: Ref.h:247
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
void associateSimpleRecHitCluster(const SiStripCluster *clust, const uint32_t &detID, std::vector< SimHitIdpr > &theSimtrackid, std::vector< PSimHit > &simhit)
std::vector< float > const & weights() const
std::vector< int > simhitCFPos
tuple result
Definition: query.py:137
ClusterRef cluster() const
EncodedEventId eventId() const
Definition: PSimHit.h:105
#define end
Definition: vmac.h:38
std::vector< SimHitIdpr > associateMultiRecHitId(const SiTrackerMultiRecHit *multirechit)
bool first
Definition: L1TdeRCT.cc:94
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
std::vector< PSimHit > associateMultiRecHit(const SiTrackerMultiRecHit *multirechit)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
tuple conf
Definition: dbtoconf.py:185
edm::Handle< CrossingFrame< PSimHit > > cf_simhit
Definition: DetId.h:20
iterator begin()
const T & getObject(unsigned int ip) const
Definition: MixCollection.h:30
const uint32_t & eeId() const
std::pair< uint32_t, EncodedEventId > SimHitIdpr
SiStripRecHit2D stereoHit() const
const uint32_t & eeId() const
ClusterRef cluster() const
Definition: SiPixelRecHit.h:41
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:41
T const * product() const
Definition: Handle.h:74
SiStripRecHit2D monoHit() const
collection_type data
Definition: DetSet.h:79
#define begin
Definition: vmac.h:31
std::vector< PSimHit > associateHit(const TrackingRecHit &thit)
unsigned int trackId() const
Definition: PSimHit.h:102
tuple cout
Definition: gather_cfg.py:121
DetId geographicalId() const
collection_type::const_iterator const_iterator
Definition: DetSet.h:34
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:106
std::vector< SimHitIdpr > associateProjectedRecHit(const ProjectedSiStripRecHit2D *projectedrechit)
const SiStripRecHit2D & originalHit() const
tuple size
Write out results.
std::vector< SimHitIdpr > associateMatchedRecHit(const SiStripMatchedRecHit2D *matchedrechit)
const std::vector< uint8_t > & amplitudes() const
void associatePixelRecHit(const SiPixelRecHit *pixelrechit, std::vector< SimHitIdpr > &simhitid)
Our base class.
Definition: SiPixelRecHit.h:22
std::vector< SimHitIdpr > associateGSMatchedRecHit(const SiTrackerGSMatchedRecHit2D *gsmrechit)