CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Types | Public Member Functions | Public Attributes | Private Types | Private Member Functions | Private Attributes
TrackerHitAssociator Class Reference

#include <TrackerHitAssociator.h>

Classes

struct  Config
 

Public Types

typedef std::map
< simHitCollectionID,
std::vector< PSimHit > > 
simhit_collectionMap
 
typedef std::map< unsigned int,
std::vector< PSimHit > > 
simhit_map
 
typedef std::pair
< simHitCollectionID, unsigned
int > 
simhitAddr
 
typedef std::pair< unsigned
int, unsigned int > 
simHitCollectionID
 

Public Member Functions

void associateCluster (const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< PSimHit > &simhit) const
 
std::vector< SimHitIdprassociateGSMatchedRecHit (const SiTrackerGSMatchedRecHit2D *gsmrechit) const
 
std::vector< SimHitIdprassociateGSRecHit (const SiTrackerGSRecHit2D *gsrechit) const
 
std::vector< PSimHitassociateHit (const TrackingRecHit &thit) const
 
std::vector< SimHitIdprassociateHitId (const TrackingRecHit &thit) const
 
void associateHitId (const TrackingRecHit &thit, std::vector< SimHitIdpr > &simhitid, std::vector< simhitAddr > *simhitCFPos=0) const
 
std::vector< SimHitIdprassociateMatchedRecHit (const SiStripMatchedRecHit2D *matchedrechit, std::vector< simhitAddr > *simhitCFPos=0) const
 
std::vector< PSimHitassociateMultiRecHit (const SiTrackerMultiRecHit *multirechit) const
 
std::vector< SimHitIdprassociateMultiRecHitId (const SiTrackerMultiRecHit *multirechit, std::vector< simhitAddr > *simhitCFPos=0) const
 
void associatePixelRecHit (const SiPixelRecHit *pixelrechit, std::vector< SimHitIdpr > &simhitid, std::vector< simhitAddr > *simhitCFPos=0) const
 
std::vector< SimHitIdprassociateProjectedRecHit (const ProjectedSiStripRecHit2D *projectedrechit, std::vector< simhitAddr > *simhitCFPos=0) const
 
void associateSimpleRecHitCluster (const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
 
template<typename T >
void associateSiStripRecHit (const T *simplerechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
 
 TrackerHitAssociator (const edm::Event &e, const Config &config)
 
virtual ~TrackerHitAssociator ()
 

Public Attributes

simhit_collectionMap SimHitCollMap
 
simhit_map SimHitMap
 

Private Types

typedef std::vector< std::string > vstring
 

Private Member Functions

void makeMaps (const edm::Event &theEvent, const Config &config)
 

Private Attributes

bool assocHitbySimTrack_
 
bool doPixel_
 
bool doStrip_
 
bool doTrackAssoc_
 
edm::Handle< edm::DetSetVector
< PixelDigiSimLink > > 
pixeldigisimlink
 
edm::Handle< edm::DetSetVector
< StripDigiSimLink > > 
stripdigisimlink
 

Detailed Description

Definition at line 55 of file TrackerHitAssociator.h.

Member Typedef Documentation

Definition at line 104 of file TrackerHitAssociator.h.

typedef std::map<unsigned int, std::vector<PSimHit> > TrackerHitAssociator::simhit_map

Definition at line 102 of file TrackerHitAssociator.h.

typedef std::pair<simHitCollectionID, unsigned int> TrackerHitAssociator::simhitAddr

Definition at line 75 of file TrackerHitAssociator.h.

typedef std::pair<unsigned int, unsigned int> TrackerHitAssociator::simHitCollectionID

Definition at line 74 of file TrackerHitAssociator.h.

typedef std::vector<std::string> TrackerHitAssociator::vstring
private

Definition at line 108 of file TrackerHitAssociator.h.

Constructor & Destructor Documentation

TrackerHitAssociator::TrackerHitAssociator ( const edm::Event e,
const Config config 
)

Definition at line 83 of file TrackerHitAssociator.cc.

References doPixel_, doStrip_, doTrackAssoc_, edm::Event::getByToken(), makeMaps(), pixeldigisimlink, TrackerHitAssociator::Config::pixelToken_, stripdigisimlink, and TrackerHitAssociator::Config::stripToken_.

83  :
84  doPixel_(config.doPixel_),
85  doStrip_(config.doStrip_),
86  doTrackAssoc_(config.doTrackAssoc_),
87  assocHitbySimTrack_(config.assocHitbySimTrack_) {
88  //if track association there is no need to access the input collections
89  if(!doTrackAssoc_) {
90  makeMaps(e, config);
91  }
92 
93  if(doStrip_) e.getByToken(config.stripToken_, stripdigisimlink);
94  if(doPixel_) e.getByToken(config.pixelToken_, pixeldigisimlink);
95 }
edm::Handle< edm::DetSetVector< StripDigiSimLink > > stripdigisimlink
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
edm::Handle< edm::DetSetVector< PixelDigiSimLink > > pixeldigisimlink
void makeMaps(const edm::Event &theEvent, const Config &config)
virtual TrackerHitAssociator::~TrackerHitAssociator ( )
inlinevirtual

Definition at line 72 of file TrackerHitAssociator.h.

72 {}

Member Function Documentation

void TrackerHitAssociator::associateCluster ( const SiStripCluster clust,
const DetId detid,
std::vector< SimHitIdpr > &  simtrackid,
std::vector< PSimHit > &  simhit 
) const

Definition at line 403 of file TrackerHitAssociator.cc.

References associateSimpleRecHitCluster(), i, and SimHitCollMap.

406  {
407  std::vector<simhitAddr> simhitCFPos;
408  associateSimpleRecHitCluster(clust, detid, simtrackid, &simhitCFPos);
409 
410  for(size_t i=0, simhitCFPosSize = simhitCFPos.size(); i<simhitCFPosSize; ++i) {
411  simhitAddr theSimHitAddr = simhitCFPos[i];
412  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
413  simhit_collectionMap::const_iterator it = SimHitCollMap.find(theSimHitCollID);
414  if (it!= SimHitCollMap.end()) {
415  unsigned int theSimHitIndex = theSimHitAddr.second;
416  if (theSimHitIndex < (it->second).size()) simhit.push_back((it->second)[theSimHitIndex]);
417  // const PSimHit& theSimHit = (it->second)[theSimHitIndex];
418  // std::cout << "For cluster, simHit detId = " << theSimHit.detUnitId() << " address = (" << (theSimHitAddr.first).first
419  // << ", " << (theSimHitAddr.first).second << ", " << theSimHitIndex
420  // << "), process = " << theSimHit.processType() << " (" << theSimHit.eventId().bunchCrossing()
421  // << ", " << theSimHit.eventId().event() << ", " << theSimHit.trackId() << ")" << std::endl;
422  }
423  }
424 }
int i
Definition: DBlmapReader.cc:9
std::pair< simHitCollectionID, unsigned int > simhitAddr
void associateSimpleRecHitCluster(const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
simhit_collectionMap SimHitCollMap
std::pair< unsigned int, unsigned int > simHitCollectionID
std::vector< SimHitIdpr > TrackerHitAssociator::associateGSMatchedRecHit ( const SiTrackerGSMatchedRecHit2D gsmrechit) const

Definition at line 653 of file TrackerHitAssociator.cc.

References SiTrackerGSMatchedRecHit2D::eeId(), and SiTrackerGSMatchedRecHit2D::simtrackId().

Referenced by associateHitId().

654 {
655  //GSRecHit is the FastSimulation RecHit that contains the TrackId already
656 
657  vector<SimHitIdpr> simtrackid;
658  simtrackid.clear();
659  SimHitIdpr currentId(gsmrechit->simtrackId(), EncodedEventId(gsmrechit->eeId()));
660  simtrackid.push_back(currentId);
661  return simtrackid;
662 }
std::pair< uint32_t, EncodedEventId > SimHitIdpr
const uint32_t & eeId() const
std::vector< SimHitIdpr > TrackerHitAssociator::associateGSRecHit ( const SiTrackerGSRecHit2D gsrechit) const

Definition at line 617 of file TrackerHitAssociator.cc.

References SiTrackerGSRecHit2D::eeId(), and SiTrackerGSRecHit2D::simtrackId().

Referenced by associateHitId().

618 {
619  //GSRecHit is the FastSimulation RecHit that contains the TrackId already
620 
621  vector<SimHitIdpr> simtrackid;
622  simtrackid.clear();
623  SimHitIdpr currentId(gsrechit->simtrackId(), EncodedEventId(gsrechit->eeId()));
624  simtrackid.push_back(currentId);
625  return simtrackid;
626 }
const int & simtrackId() const
const uint32_t & eeId() const
std::pair< uint32_t, EncodedEventId > SimHitIdpr
std::vector< PSimHit > TrackerHitAssociator::associateHit ( const TrackingRecHit thit) const

Check if it's the gluedDet.

Definition at line 187 of file TrackerHitAssociator.cc.

References assocHitbySimTrack_, associateHitId(), associateMultiRecHit(), begin, cond::rpcobgas::detid, doTrackAssoc_, end, PSimHit::eventId(), plotBeamSpotDB::first, TrackingRecHit::geographicalId(), i, DetId::rawId(), query::result, edm::second(), SimHitCollMap, SimHitMap, and PSimHit::trackId().

Referenced by CkfDebugger::analyseCompatibleMeasurements(), SiPixelRecHitsValid::analyze(), TestHits::analyze(), TestSmoothHits::analyze(), TestTrackHits::analyze(), TestOutliers::analyze(), SiPixelErrorEstimation::analyze(), SiPixelTrackingRecHitsValid::analyze(), CkfDebugger::associated(), associateMultiRecHit(), CkfDebugger::correctTrajectory(), GlobalRecHitsAnalyzer::fillTrk(), GlobalRecHitsProducer::fillTrk(), spr::matchedSimTrack(), CkfDebugger::nextCorrectHits(), SiStripRecHitsValid::rechitanalysis(), SiStripTrackingRecHitsValid::rechitanalysis(), SiStripRecHitsValid::rechitanalysis_matched(), SiStripTrackingRecHitsValid::rechitanalysis_matched(), and CkfDebugger::testSeed().

188 {
189 
190  if (const SiTrackerMultiRecHit * rechit = dynamic_cast<const SiTrackerMultiRecHit *>(&thit)){
191  return associateMultiRecHit(rechit);
192  }
193 
194  //vector with the matched SimHit
195  std::vector<PSimHit> result;
196 
197  if(doTrackAssoc_) return result;
198 
199  // Vectors to contain lists of matched simTracks, simHits
200  std::vector<SimHitIdpr> simtrackid;
201  std::vector<simhitAddr> simhitCFPos;
202 
203  //get the Detector type of the rechit
204  DetId detid= thit.geographicalId();
205  uint32_t detID = detid.rawId();
206 
207  // Get the vectors of simtrackIDs and simHit indices associated with this rechit
208  associateHitId(thit, simtrackid, &simhitCFPos);
209  // std::cout << "recHit subdet, detID = " << detid.subdetId() << ", " << detID << ", (bnch, evt, trk) = ";
210  // for (size_t i=0; i<simtrackid.size(); ++i)
211  // std::cout << ", (" << simtrackid[i].second.bunchCrossing() << ", "
212  // << simtrackid[i].second.event() << ", " << simtrackid[i].first << ")";
213  // std::cout << std::endl;
214 
215  // Get the vector of simHits associated with this rechit
216 
217  if (!assocHitbySimTrack_ && simhitCFPos.size() > 0) {
218  // We use the indices to the simHit collections taken
219  // from the DigiSimLinks and returned in simhitCFPos.
220  // simhitCFPos[i] contains the full address of the ith simhit:
221  // <collection, index> = <<subdet, tofBin>, index>
222 
223  //check if the recHit is a SiStripMatchedRecHit2D
224  bool isMatchedHit = false;
225  if(dynamic_cast<const SiStripMatchedRecHit2D *>(&thit))
226  isMatchedHit = true;
227 
228  size_t simtrackidSize = simtrackid.size();
229  for(size_t i=0, simhitCFPosSize = simhitCFPos.size(); i<simhitCFPosSize; ++i) {
230  simhitAddr theSimHitAddr = simhitCFPos[i];
231  simHitCollectionID theSimHitCollID = theSimHitAddr.first;
232  simhit_collectionMap::const_iterator it = SimHitCollMap.find(theSimHitCollID);
233  if (it!= SimHitCollMap.end()) {
234  unsigned int theSimHitIndex = theSimHitAddr.second;
235  if (theSimHitIndex < (it->second).size()) {
236  const PSimHit& theSimHit = (it->second)[theSimHitIndex];
237  if (isMatchedHit) {
238  // Try to remove ghosts by requiring a match to the simTrack also
239  unsigned int simHitid = theSimHit.trackId();
240  EncodedEventId simHiteid = theSimHit.eventId();
241  for(size_t i=0; i<simtrackidSize;++i) {
242  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second) {
243  result.push_back(theSimHit);
244  }
245  }
246  } else {
247  result.push_back(theSimHit);
248  }
249  // std::cout << "by CFpos, simHit detId = " << theSimHit.detUnitId() << " address = (" << (theSimHitAddr.first).first
250  // << ", " << (theSimHitAddr.first).second << ", " << theSimHitIndex
251  // << "), process = " << theSimHit.processType() << " (" << theSimHit.eventId().bunchCrossing()
252  // << ", " << theSimHit.eventId().event() << ", " << theSimHit.trackId() << ")" << std::endl;
253  }
254  }
255  }
256  return result;
257  }
258 
259  // Get the SimHit from the trackid instead
260  std::map<unsigned int, std::vector<PSimHit> >::const_iterator it = SimHitMap.find(detID);
261  if (it!= SimHitMap.end()) {
262  vector<PSimHit>::const_iterator simHitIter = (it->second).begin();
263  vector<PSimHit>::const_iterator simHitIterEnd = (it->second).end();
264  size_t simtrackidSize = simtrackid.size();
265  for (;simHitIter != simHitIterEnd; ++simHitIter) {
266  const PSimHit& ihit = *simHitIter;
267  unsigned int simHitid = ihit.trackId();
268  EncodedEventId simHiteid = ihit.eventId();
269  // std::cout << "by simTk, simHit, process = " << ihit.processType() << " (" << ihit.eventId().bunchCrossing()
270  // << ", " << ihit.eventId().event() << ", " << ihit.trackId() << ")";
271 
272  for(size_t i=0; i<simtrackidSize;++i) {
273  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second) {
274 // cout << "Associator ---> ID" << ihit.trackId() << " Simhit x= " << ihit.localPosition().x()
275 // << " y= " << ihit.localPosition().y() << " z= " << ihit.localPosition().x() << endl;
276  // std::cout << " matches";
277  result.push_back(ihit);
278  break;
279  }
280  }
281  // std::cout << std::endl;
282  }
283 
284  }else{
285 
287  std::map<unsigned int, std::vector<PSimHit> >::const_iterator itrphi =
288  SimHitMap.find(detID+2); //iterator to the simhit in the rphi module
289  std::map<unsigned int, std::vector<PSimHit> >::const_iterator itster =
290  SimHitMap.find(detID+1);//iterator to the simhit in the stereo module
291  if (itrphi!= SimHitMap.end()&&itster!=SimHitMap.end()) {
292  std::vector<PSimHit> simHitVector = itrphi->second;
293  simHitVector.insert(simHitVector.end(),(itster->second).begin(),(itster->second).end());
294  vector<PSimHit>::const_iterator simHitIter = simHitVector.begin();
295  vector<PSimHit>::const_iterator simHitIterEnd = simHitVector.end();
296  size_t simtrackidSize = simtrackid.size();
297  for (;simHitIter != simHitIterEnd; ++simHitIter) {
298  const PSimHit& ihit = *simHitIter;
299  unsigned int simHitid = ihit.trackId();
300  EncodedEventId simHiteid = ihit.eventId();
301  for(size_t i=0; i<simtrackidSize; ++i) {
302  if(simHitid == simtrackid[i].first && simHiteid == simtrackid[i].second) {
303 // if(simHitid == simtrackid[i].first && simHiteid.bunchCrossing() == simtrackid[i].second.bunchCrossing()) {
304  // cout << "GluedDet Associator ---> ID" << ihit.trackId() << " Simhit x= " << ihit.localPosition().x()
305  // << " y= " << ihit.localPosition().y() << " z= " << ihit.localPosition().x() << endl;
306  result.push_back(ihit);
307  break;
308  }
309  }
310  }
311  }
312  }
313 
314  return result;
315 }
int i
Definition: DBlmapReader.cc:9
std::vector< PSimHit > associateMultiRecHit(const SiTrackerMultiRecHit *multirechit) const
std::pair< simHitCollectionID, unsigned int > simhitAddr
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
U second(std::pair< T, U > const &p)
tuple result
Definition: query.py:137
EncodedEventId eventId() const
Definition: PSimHit.h:105
#define end
Definition: vmac.h:37
simhit_collectionMap SimHitCollMap
Definition: DetId.h:18
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
#define begin
Definition: vmac.h:30
unsigned int trackId() const
Definition: PSimHit.h:102
DetId geographicalId() const
std::pair< unsigned int, unsigned int > simHitCollectionID
std::vector< SimHitIdpr > TrackerHitAssociator::associateHitId ( const TrackingRecHit thit) const
void TrackerHitAssociator::associateHitId ( const TrackingRecHit thit,
std::vector< SimHitIdpr > &  simhitid,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 324 of file TrackerHitAssociator.cc.

References associateGSMatchedRecHit(), associateGSRecHit(), associateMatchedRecHit(), associateMultiRecHitId(), associatePixelRecHit(), associateProjectedRecHit(), associateSiStripRecHit(), cond::rpcobgas::detid, TrackingRecHit::geographicalId(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, DetId::subdetId(), StripSubdetector::TEC, StripSubdetector::TIB, StripSubdetector::TID, and StripSubdetector::TOB.

326 {
327 
328  simtkid.clear();
329 
330  //get the Detector type of the rechit
331  DetId detid= thit.geographicalId();
332  if (const SiTrackerMultiRecHit * rechit = dynamic_cast<const SiTrackerMultiRecHit *>(&thit)){
333  simtkid=associateMultiRecHitId(rechit, simhitCFPos);
334  }
335 
336  //cout << "Associator ---> get Detid " << detID << endl;
337  //check we are in the strip tracker
338  if(detid.subdetId() == StripSubdetector::TIB ||
339  detid.subdetId() == StripSubdetector::TOB ||
340  detid.subdetId() == StripSubdetector::TID ||
341  detid.subdetId() == StripSubdetector::TEC)
342  {
343  //check if it is a simple SiStripRecHit2D
344  if(const SiStripRecHit2D * rechit =
345  dynamic_cast<const SiStripRecHit2D *>(&thit))
346  {
347  associateSiStripRecHit(rechit, simtkid, simhitCFPos);
348  }
349  //check if it is a SiStripRecHit1D
350  else if(const SiStripRecHit1D * rechit =
351  dynamic_cast<const SiStripRecHit1D *>(&thit))
352  {
353  associateSiStripRecHit(rechit, simtkid, simhitCFPos);
354  }
355  //check if it is a SiStripMatchedRecHit2D
356  else if(const SiStripMatchedRecHit2D * rechit =
357  dynamic_cast<const SiStripMatchedRecHit2D *>(&thit))
358  {
359  simtkid = associateMatchedRecHit(rechit, simhitCFPos);
360  }
361  //check if it is a ProjectedSiStripRecHit2D
362  else if(const ProjectedSiStripRecHit2D * rechit =
363  dynamic_cast<const ProjectedSiStripRecHit2D *>(&thit))
364  {
365  simtkid = associateProjectedRecHit(rechit, simhitCFPos);
366  }
367  else{
368  //std::cout << "associate to invalid" << std::endl;
369  //throw cms::Exception("Unknown RecHit Type") << "TrackerHitAssociator failed second casting of " << typeid(thit).name() << " type ";
370  }
371  }
372  //check we are in the pixel tracker
373  else if( (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelBarrel ||
374  (unsigned int)(detid.subdetId()) == PixelSubdetector::PixelEndcap)
375  {
376  if(const SiPixelRecHit * rechit = dynamic_cast<const SiPixelRecHit *>(&thit))
377  {
378  associatePixelRecHit(rechit, simtkid, simhitCFPos);
379  }
380  }
381  //check if these are GSRecHits (from FastSim)
382  if(const SiTrackerGSRecHit2D * rechit = dynamic_cast<const SiTrackerGSRecHit2D *>(&thit))
383  {
384  simtkid = associateGSRecHit(rechit);
385  }
386  if(const SiTrackerGSMatchedRecHit2D * rechit = dynamic_cast<const SiTrackerGSMatchedRecHit2D *>(&thit))
387  {
388  simtkid = associateGSMatchedRecHit(rechit);
389  }
390 
391 }
std::vector< SimHitIdpr > associateMultiRecHitId(const SiTrackerMultiRecHit *multirechit, std::vector< simhitAddr > *simhitCFPos=0) const
void associatePixelRecHit(const SiPixelRecHit *pixelrechit, std::vector< SimHitIdpr > &simhitid, std::vector< simhitAddr > *simhitCFPos=0) const
std::vector< SimHitIdpr > associateGSMatchedRecHit(const SiTrackerGSMatchedRecHit2D *gsmrechit) const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
Definition: DetId.h:18
std::vector< SimHitIdpr > associateProjectedRecHit(const ProjectedSiStripRecHit2D *projectedrechit, std::vector< simhitAddr > *simhitCFPos=0) const
std::vector< SimHitIdpr > associateGSRecHit(const SiTrackerGSRecHit2D *gsrechit) const
DetId geographicalId() const
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
Our base class.
Definition: SiPixelRecHit.h:23
std::vector< SimHitIdpr > TrackerHitAssociator::associateMatchedRecHit ( const SiStripMatchedRecHit2D matchedrechit,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 505 of file TrackerHitAssociator.cc.

References associateSiStripRecHit(), spr::find(), SiStripMatchedRecHit2D::monoHit(), and SiStripMatchedRecHit2D::stereoHit().

Referenced by associateHitId().

506 {
507  std::vector<SimHitIdpr> matched_mono;
508  std::vector<SimHitIdpr> matched_st;
509 
510  const SiStripRecHit2D mono = matchedrechit->monoHit();
511  const SiStripRecHit2D st = matchedrechit->stereoHit();
512  //associate the two simple hits separately
513  associateSiStripRecHit(&mono, matched_mono, simhitCFPos);
514  associateSiStripRecHit(&st, matched_st, simhitCFPos);
515 
516  //save in a vector all the simtrack-id's that are common to mono and stereo hits
517  std::vector<SimHitIdpr> simtrackid;
518  if(!matched_mono.empty() && !matched_st.empty()){
519  std::vector<SimHitIdpr> idcachev;
520  for(vector<SimHitIdpr>::iterator mhit=matched_mono.begin(), mhitEnd = matched_mono.end(); mhit != mhitEnd; ++mhit){
521  //save only once the ID
522  if(find(idcachev.begin(), idcachev.end(), (*mhit)) == idcachev.end()) {
523  idcachev.push_back(*mhit);
524  //save if the stereoID matched the monoID
525  if(find(matched_st.begin(), matched_st.end(), (*mhit)) != matched_st.end()) {
526  simtrackid.push_back(*mhit);
527  }
528  }
529  }
530  }
531 
532  return simtrackid;
533 }
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
SiStripRecHit2D stereoHit() const
SiStripRecHit2D monoHit() const
void associateSiStripRecHit(const T *simplerechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
std::vector< PSimHit > TrackerHitAssociator::associateMultiRecHit ( const SiTrackerMultiRecHit multirechit) const

Definition at line 628 of file TrackerHitAssociator.cc.

References associateHit(), i, SiTrackerMultiRecHit::recHits(), findQualityFiles::size, SiTrackerMultiRecHit::weight(), and SiTrackerMultiRecHit::weights().

Referenced by associateHit().

628  {
629  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
630  // std::vector<PSimHit> assimhits;
631  int size=multirechit->weights().size(), idmostprobable=0;
632 
633  for (int i=0; i<size; ++i){
634  if(multirechit->weight(i)>multirechit->weight(idmostprobable)) idmostprobable=i;
635  }
636 
637  return associateHit(*componenthits[idmostprobable]);
638 }
int i
Definition: DBlmapReader.cc:9
float weight(unsigned int i) const
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
std::vector< float > const & weights() const
std::vector< PSimHit > associateHit(const TrackingRecHit &thit) const
tuple size
Write out results.
std::vector< SimHitIdpr > TrackerHitAssociator::associateMultiRecHitId ( const SiTrackerMultiRecHit multirechit,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 640 of file TrackerHitAssociator.cc.

References associateHitId(), i, SiTrackerMultiRecHit::recHits(), findQualityFiles::size, SiTrackerMultiRecHit::weight(), and SiTrackerMultiRecHit::weights().

Referenced by associateHitId().

640  {
641  std::vector<const TrackingRecHit*> componenthits = multirechit->recHits();
642  int size=multirechit->weights().size(), idmostprobable=0;
643 
644  for (int i=0; i<size; ++i){
645  if(multirechit->weight(i)>multirechit->weight(idmostprobable)) idmostprobable=i;
646  }
647 
648  std::vector< SimHitIdpr > simhitid;
649  associateHitId(*componenthits[idmostprobable], simhitid, simhitCFPos);
650  return simhitid;
651 }
int i
Definition: DBlmapReader.cc:9
float weight(unsigned int i) const
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
std::vector< float > const & weights() const
std::vector< SimHitIdpr > associateHitId(const TrackingRecHit &thit) const
tuple size
Write out results.
void TrackerHitAssociator::associatePixelRecHit ( const SiPixelRecHit pixelrechit,
std::vector< SimHitIdpr > &  simhitid,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 548 of file TrackerHitAssociator.cc.

References PixelDigi::channelToPixel(), SiPixelRecHit::cluster(), edm::DetSet< T >::data, cond::rpcobgas::detid, spr::find(), TrackingRecHit::geographicalId(), edm::Ref< C, T, F >::isNull(), pixeldigisimlink, DetId::rawId(), and DetId::subdetId().

Referenced by associateHitId().

551 {
552  //
553  // Pixel associator
554  //
555  DetId detid= pixelrechit->geographicalId();
556  uint32_t detID = detid.rawId();
557 
559  if(isearch != pixeldigisimlink->end()) { //if it is not empty
560  edm::DetSet<PixelDigiSimLink> link_detset = (*isearch);
561  SiPixelRecHit::ClusterRef const& cluster = pixelrechit->cluster();
562 
563  //check the reference is valid
564 
565  if(!(cluster.isNull())){//if the cluster is valid
566 
567  int minPixelRow = (*cluster).minPixelRow();
568  int maxPixelRow = (*cluster).maxPixelRow();
569  int minPixelCol = (*cluster).minPixelCol();
570  int maxPixelCol = (*cluster).maxPixelCol();
571  //std::cout << " Cluster minRow " << minPixelRow << " maxRow " << maxPixelRow << std::endl;
572  //std::cout << " Cluster minCol " << minPixelCol << " maxCol " << maxPixelCol << std::endl;
573  edm::DetSet<PixelDigiSimLink>::const_iterator linkiter = link_detset.data.begin(), linkEnd = link_detset.data.end();
574  int dsl = 0;
575  std::vector<SimHitIdpr> idcachev;
576  std::vector<simhitAddr> CFposcachev;
577  for( ; linkiter != linkEnd; ++linkiter) {
578  ++dsl;
579  std::pair<int,int> pixel_coord = PixelDigi::channelToPixel(linkiter->channel());
580  //std::cout << " " << dsl << ") Digi link: row " << pixel_coord.first << " col " << pixel_coord.second << std::endl;
581  if( pixel_coord.first <= maxPixelRow &&
582  pixel_coord.first >= minPixelRow &&
583  pixel_coord.second <= maxPixelCol &&
584  pixel_coord.second >= minPixelCol ) {
585  //std::cout << " !-> trackid " << linkiter->SimTrackId() << endl;
586  //std::cout << " fraction " << linkiter->fraction() << endl;
587  SimHitIdpr currentId(linkiter->SimTrackId(), linkiter->eventId());
588  if(find(idcachev.begin(),idcachev.end(),currentId) == idcachev.end()){
589  simtrackid.push_back(currentId);
590  idcachev.push_back(currentId);
591  }
592 
593  if (simhitCFPos != 0) {
594  //create a vector that contains all the position (in the MixCollection) of the SimHits that contributed to the RecHit
595  //write position only once
596  unsigned int currentCFPos = linkiter->CFposition();
597  unsigned int tofBin = linkiter->TofBin();
598  simHitCollectionID theSimHitCollID = std::make_pair(detid.subdetId(), tofBin);
599  simhitAddr currentAddr = std::make_pair(theSimHitCollID, currentCFPos);
600 
601  if(find(CFposcachev.begin(), CFposcachev.end(), currentAddr) == CFposcachev.end()) {
602  CFposcachev.push_back(currentAddr);
603  simhitCFPos->push_back(currentAddr);
604  }
605  }
606 
607  }
608  }
609  }
610  else{
611  edm::LogError("TrackerHitAssociator")<<"no Pixel cluster reference attached";
612 
613  }
614  }
615 }
std::pair< simHitCollectionID, unsigned int > simhitAddr
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
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
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
Definition: DetId.h:18
std::pair< uint32_t, EncodedEventId > SimHitIdpr
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
static std::pair< int, int > channelToPixel(int ch)
Definition: PixelDigi.h:62
collection_type data
Definition: DetSet.h:78
DetId geographicalId() const
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
std::pair< unsigned int, unsigned int > simHitCollectionID
std::vector< SimHitIdpr > TrackerHitAssociator::associateProjectedRecHit ( const ProjectedSiStripRecHit2D projectedrechit,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 536 of file TrackerHitAssociator.cc.

References associateSiStripRecHit(), and ProjectedSiStripRecHit2D::originalHit().

Referenced by associateHitId().

538 {
539  //projectedRecHit is a "matched" rechit with only one component
540 
541  std::vector<SimHitIdpr> matched_mono;
542 
543  const SiStripRecHit2D mono = projectedrechit->originalHit();
544  associateSiStripRecHit(&mono, matched_mono, simhitCFPos);
545  return matched_mono;
546 }
SiStripRecHit2D originalHit() const
void associateSiStripRecHit(const T *simplerechit, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
void TrackerHitAssociator::associateSimpleRecHitCluster ( const SiStripCluster clust,
const DetId detid,
std::vector< SimHitIdpr > &  simtrackid,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 426 of file TrackerHitAssociator.cc.

References SiStripCluster::amplitudes(), edm::DetSet< T >::data, spr::find(), plotBeamSpotDB::first, SiStripCluster::firstStrip(), prof2calltree::last, DetId::rawId(), stripdigisimlink, and DetId::subdetId().

Referenced by associateCluster(), and associateSiStripRecHit().

429  {
430 
431  uint32_t detID = detid.rawId();
433  if(isearch != stripdigisimlink->end()) { //if it is not empty
434  edm::DetSet<StripDigiSimLink> link_detset = (*isearch);
435 
436  if(clust!=0){//the cluster is valid
437  int clusiz = clust->amplitudes().size();
438  int first = clust->firstStrip();
439  int last = first + clusiz;
440 
441 // std::cout << "CLUSTERSIZE " << clusiz << " first strip = " << first << " last strip = " << last-1 << std::endl;
442 // std::cout << " detID = " << detID << " DETSET size = " << link_detset.data.size() << std::endl;
443  //use a vector
444  std::vector<SimHitIdpr> idcachev;
445  std::vector<simhitAddr> CFposcachev;
446  for(edm::DetSet<StripDigiSimLink>::const_iterator linkiter = link_detset.data.begin(), linkEnd = link_detset.data.end(); linkiter != linkEnd; ++linkiter){
447  if( (int)(linkiter->channel()) >= first && (int)(linkiter->channel()) < last ){
448 
449  //check this digisimlink
450 // printf("%s%4d%s%8d%s%3d%s%8.4f\n", "CHANNEL = ", linkiter->channel(), " TrackID = ", linkiter->SimTrackId(),
451 // " tofBin = ", linkiter->TofBin(), " fraction = ", linkiter->fraction());
452  /*
453  std::cout << "CHECKING CHANNEL = " << linkiter->channel() << std::endl;
454  std::cout << "TrackID = " << linkiter->SimTrackId() << std::endl;
455  std::cout << "Position = " << linkiter->CFposition() << std::endl;
456  std::cout << " fraction = " << linkiter->fraction() << std::endl;
457  */
458 
459  SimHitIdpr currentId(linkiter->SimTrackId(), linkiter->eventId());
460 
461  //create a vector with the list of SimTrack ID's of the tracks that contributed to the RecHit
462  //write the id only once in the vector
463 
464  if(find(idcachev.begin(),idcachev.end(),currentId ) == idcachev.end()){
465  /*
466  std::cout << " Adding track id = " << currentId.first
467  << " Event id = " << currentId.second.event()
468  << " Bunch Xing = " << currentId.second.bunchCrossing()
469  << std::endl;
470  */
471  idcachev.push_back(currentId);
472  simtrackid.push_back(currentId);
473  }
474 
475  if (simhitCFPos != 0) {
476  //create a vector that contains all the position (in the MixCollection) of the SimHits that contributed to the RecHit
477  //write position only once
478  unsigned int currentCFPos = linkiter->CFposition();
479  unsigned int tofBin = linkiter->TofBin();
480  simHitCollectionID theSimHitCollID = std::make_pair(detid.subdetId(), tofBin);
481  simhitAddr currentAddr = std::make_pair(theSimHitCollID, 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_collectionMap::const_iterator it = SimHitCollMap.find(theSimHitCollID);
487 // if (it!= SimHitCollMap.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  else {
500  edm::LogError("TrackerHitAssociator")<<"no cluster reference attached";
501  }
502  }
503 }
edm::Handle< edm::DetSetVector< StripDigiSimLink > > stripdigisimlink
std::pair< simHitCollectionID, unsigned int > simhitAddr
uint16_t firstStrip() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::pair< uint32_t, EncodedEventId > SimHitIdpr
collection_type data
Definition: DetSet.h:78
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
const std::vector< uint8_t > & amplitudes() const
std::pair< unsigned int, unsigned int > simHitCollectionID
template<typename T >
void TrackerHitAssociator::associateSiStripRecHit ( const T simplerechit,
std::vector< SimHitIdpr > &  simtrackid,
std::vector< simhitAddr > *  simhitCFPos = 0 
) const

Definition at line 394 of file TrackerHitAssociator.cc.

References associateSimpleRecHitCluster().

Referenced by associateHitId(), associateMatchedRecHit(), and associateProjectedRecHit().

395 {
396  const SiStripCluster* clust = &(*simplerechit->cluster());
397  associateSimpleRecHitCluster(clust, simplerechit->geographicalId(), simtrackid, simhitCFPos);
398 }
void associateSimpleRecHitCluster(const SiStripCluster *clust, const DetId &detid, std::vector< SimHitIdpr > &simtrackid, std::vector< simhitAddr > *simhitCFPos=0) const
void TrackerHitAssociator::makeMaps ( const edm::Event theEvent,
const Config config 
)
private

Definition at line 97 of file TrackerHitAssociator.cc.

References assocHitbySimTrack_, TrackerHitAssociator::Config::cfTokens_, edm::Event::getByToken(), StripDigiSimLink::HighTof, HLT_25ns14e33_v1_cff::labels, edm::Event::labelsForToken(), StripDigiSimLink::LowTof, NULL, edm::Handle< T >::product(), edm::ProductLabels::productInstance, SimHitCollMap, SimHitMap, trackerHits::simHits, and TrackerHitAssociator::Config::simHitTokens_.

Referenced by TrackerHitAssociator().

97  {
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  const char* const highTag = "HighTof";
103  unsigned int tofBin;
105  if (assocHitbySimTrack_) {
106  for(auto const& cfToken : config.cfTokens_) {
108  //int Nhits = 0;
109  if (theEvent.getByToken(cfToken, cf_simhit)) {
110  std::unique_ptr<MixCollection<PSimHit> > thisContainerHits(new MixCollection<PSimHit>(cf_simhit.product()));
111  theEvent.labelsForToken(cfToken, labels);
112  if(std::strstr(labels.productInstance, highTag) != NULL) {
113  tofBin = StripDigiSimLink::HighTof;
114  } else {
115  tofBin = StripDigiSimLink::LowTof;
116  }
117  for (auto const& isim : *thisContainerHits) {
118  DetId theDet(isim.detUnitId());
119  SimHitMap[theDet].push_back(isim);
120  // ++Nhits;
121  }
122  // std::cout << "simHits from crossing frames; map size = " << SimHitCollMap.size() << ", Hit count = " << Nhits << std::endl;
123  }
124  }
125  for(auto const& simHitToken : config.simHitTokens_) {
127  //int Nhits = 0;
128  if(theEvent.getByToken(simHitToken, simHits)) {
129  theEvent.labelsForToken(simHitToken, labels);
130  if(std::strstr(labels.productInstance, highTag) != NULL) {
131  tofBin = StripDigiSimLink::HighTof;
132  } else {
133  tofBin = StripDigiSimLink::LowTof;
134  }
135  for (auto const& isim : *simHits) {
136  DetId theDet(isim.detUnitId());
137  SimHitMap[theDet].push_back(isim);
138  //++Nhits;
139  }
140  // std::cout << "simHits from prompt collections; map size = " << SimHitCollMap.size() << ", Hit count = " << Nhits << std::endl;
141  }
142  }
143  } else {
144  simHitCollectionID theSimHitCollID;
145  for(auto const& cfToken : config.cfTokens_) {
147  //int Nhits = 0;
148  if (theEvent.getByToken(cfToken, cf_simhit)) {
149  std::unique_ptr<MixCollection<PSimHit> > thisContainerHits(new MixCollection<PSimHit>(cf_simhit.product()));
150  theEvent.labelsForToken(cfToken, labels);
151  if(std::strstr(labels.productInstance, highTag) != NULL) {
152  tofBin = StripDigiSimLink::HighTof;
153  } else {
154  tofBin = StripDigiSimLink::LowTof;
155  }
156  for (auto const& isim : *thisContainerHits) {
157  DetId theDet(isim.detUnitId());
158  theSimHitCollID = std::make_pair(theDet.subdetId(), tofBin);
159  SimHitCollMap[theSimHitCollID].push_back(isim);
160  //++Nhits;
161  }
162  // std::cout << "simHits from crossing frames; map size = " << SimHitCollMap.size() << ", Hit count = " << Nhits << std::endl;
163  }
164  }
165  for(auto const& simHitToken : config.simHitTokens_) {
167  //int Nhits = 0;
168  if(theEvent.getByToken(simHitToken, simHits)) {
169  theEvent.labelsForToken(simHitToken, labels);
170  if(std::strstr(labels.productInstance, highTag) != NULL) {
171  tofBin = StripDigiSimLink::HighTof;
172  } else {
173  tofBin = StripDigiSimLink::LowTof;
174  }
175  for (auto const& isim : *simHits) {
176  DetId theDet(isim.detUnitId());
177  theSimHitCollID = std::make_pair(theDet.subdetId(), tofBin);
178  SimHitCollMap[theSimHitCollID].push_back(isim);
179  //++Nhits;
180  }
181  // std::cout << "simHits from prompt collections; map size = " << SimHitCollMap.size() << ", Hit count = " << Nhits << std::endl;
182  }
183  }
184  }
185 }
void labelsForToken(EDGetToken const &iToken, ProductLabels &oLabels) const
Definition: Event.h:221
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define NULL
Definition: scimark2.h:8
simhit_collectionMap SimHitCollMap
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
tuple simHits
Definition: trackerHits.py:16
char const * productInstance
Definition: ProductLabels.h:6
std::pair< unsigned int, unsigned int > simHitCollectionID

Member Data Documentation

bool TrackerHitAssociator::assocHitbySimTrack_
private

Definition at line 113 of file TrackerHitAssociator.h.

Referenced by associateHit(), and makeMaps().

bool TrackerHitAssociator::doPixel_
private

Definition at line 113 of file TrackerHitAssociator.h.

Referenced by TrackerHitAssociator().

bool TrackerHitAssociator::doStrip_
private

Definition at line 113 of file TrackerHitAssociator.h.

Referenced by TrackerHitAssociator().

bool TrackerHitAssociator::doTrackAssoc_
private

Definition at line 113 of file TrackerHitAssociator.h.

Referenced by associateHit(), and TrackerHitAssociator().

edm::Handle< edm::DetSetVector<PixelDigiSimLink> > TrackerHitAssociator::pixeldigisimlink
private

Definition at line 112 of file TrackerHitAssociator.h.

Referenced by associatePixelRecHit(), and TrackerHitAssociator().

simhit_collectionMap TrackerHitAssociator::SimHitCollMap

Definition at line 105 of file TrackerHitAssociator.h.

Referenced by associateCluster(), associateHit(), and makeMaps().

simhit_map TrackerHitAssociator::SimHitMap
edm::Handle< edm::DetSetVector<StripDigiSimLink> > TrackerHitAssociator::stripdigisimlink
private

Definition at line 111 of file TrackerHitAssociator.h.

Referenced by associateSimpleRecHitCluster(), and TrackerHitAssociator().