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