CMS 3D CMS Logo

RPCHitAssociator.cc

Go to the documentation of this file.
00001 #include "SimMuon/MCTruth/interface/RPCHitAssociator.h"
00002 
00003 using namespace std;
00004 
00005 typedef std::pair<uint32_t, EncodedEventId> SimHitIdpr;
00006 
00007 // Constructor
00008 RPCHitAssociator::RPCHitAssociator(const edm::Event& e, const edm::EventSetup& eventSetup, const edm::ParameterSet& conf):
00009 
00010   RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag")),
00011   RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag"))
00012 
00013 {
00014   
00015   edm::Handle<CrossingFrame<PSimHit> > cf;
00016   LogTrace("RPCHitAssociator") <<"getting CrossingFrame<PSimHit> collection - "<<RPCsimhitsTag;
00017   e.getByLabel(RPCsimhitsTag, cf);
00018   
00019   std::auto_ptr<MixCollection<PSimHit> > 
00020     hits( new MixCollection<PSimHit>(cf.product()) );
00021   MixCollection<PSimHit> & simHits = *hits;
00022 
00023   for(MixCollection<PSimHit>::MixItr hitItr = simHits.begin();
00024       hitItr != simHits.end(); ++hitItr) 
00025   {
00026     _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
00027   }
00028 
00029   edm::Handle< edm::DetSetVector<RPCDigiSimLink> > thelinkDigis;
00030   LogTrace("RPCHitAssociator") <<"getting RPCDigiSimLink collection - "<<RPCdigisimlinkTag;
00031   e.getByLabel(RPCdigisimlinkTag, thelinkDigis);
00032   _thelinkDigis = thelinkDigis;
00033 }
00034 // end of constructor
00035 
00036 std::vector<SimHitIdpr> RPCHitAssociator::associateRecHit(const TrackingRecHit & hit) {
00037  
00038   std::vector<SimHitIdpr> matched;
00039   matched.clear();
00040 
00041   const TrackingRecHit * hitp = &hit;
00042   const RPCRecHit * rpcrechit = dynamic_cast<const RPCRecHit *>(hitp);
00043 
00044   RPCDetId rpcDetId = rpcrechit->rpcId();
00045   int fstrip = rpcrechit->firstClusterStrip();
00046   int cls = rpcrechit->clusterSize();
00047   int bx = rpcrechit->BunchX();
00048 
00049   for(int i = fstrip; i < fstrip+cls; ++i) {
00050     std::set<RPCDigiSimLink> links = findRPCDigiSimLink(rpcDetId.rawId(),i,bx);
00051     for(std::set<RPCDigiSimLink>::iterator itlink = links.begin(); itlink != links.end(); ++itlink) {
00052       SimHitIdpr currentId(itlink->getTrackId(),itlink->getEventId());
00053       if(find(matched.begin(),matched.end(),currentId ) == matched.end())
00054         matched.push_back(currentId);
00055     }
00056   }
00057 
00058   return  matched;
00059 }
00060   
00061 std::set<RPCDigiSimLink>  RPCHitAssociator::findRPCDigiSimLink(uint32_t rpcDetId, int strip, int bx){
00062 
00063   std::set<RPCDigiSimLink> links;
00064 
00065   for (edm::DetSetVector<RPCDigiSimLink>::const_iterator itlink = _thelinkDigis->begin(); itlink != _thelinkDigis->end(); itlink++){
00066     for(edm::DetSet<RPCDigiSimLink>::const_iterator digi_iter=itlink->data.begin();digi_iter != itlink->data.end();++digi_iter){
00067 
00068       uint32_t detid = digi_iter->getDetUnitId();
00069       int str = digi_iter->getStrip();
00070       int bunchx = digi_iter->getBx();
00071 
00072       if(detid == rpcDetId && str == strip && bunchx == bx){
00073         links.insert(*digi_iter);
00074       }
00075 
00076     }
00077   }
00078 
00079   return links;
00080 }
00081 
00082 

Generated on Tue Jun 9 17:47:41 2009 for CMSSW by  doxygen 1.5.4