Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:25:27

0001 #include "SimMuon/MCTruth/interface/RPCHitAssociator.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 
0004 using namespace std;
0005 
0006 // Constructor
0007 RPCHitAssociator::Config::Config(const edm::ParameterSet &conf, edm::ConsumesCollector iC)
0008     : RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag")),
0009       // CrossingFrame used or not ?
0010       RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsTag")),
0011       RPCsimhitsXFTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag")),
0012       crossingframe(conf.getParameter<bool>("crossingframe")) {
0013   if (crossingframe) {
0014     RPCsimhitsXFToken_ = iC.consumes<CrossingFrame<PSimHit>>(RPCsimhitsXFTag);
0015   } else if (!RPCsimhitsTag.label().empty()) {
0016     RPCsimhitsToken_ = iC.consumes<edm::PSimHitContainer>(RPCsimhitsTag);
0017   }
0018 
0019   RPCdigisimlinkToken_ = iC.consumes<edm::DetSetVector<RPCDigiSimLink>>(RPCdigisimlinkTag);
0020 }
0021 
0022 RPCHitAssociator::RPCHitAssociator(const edm::Event &e, const Config &conf) : theConfig(conf) { initEvent(e); }
0023 
0024 void RPCHitAssociator::initEvent(const edm::Event &e)
0025 
0026 {
0027   if (theConfig.crossingframe) {
0028     LogTrace("RPCHitAssociator") << "getting CrossingFrame<PSimHit> collection - " << theConfig.RPCsimhitsXFTag;
0029     CrossingFrame<PSimHit> const &cf = e.get(theConfig.RPCsimhitsXFToken_);
0030 
0031     std::unique_ptr<MixCollection<PSimHit>> RPCsimhits(new MixCollection<PSimHit>(&cf));
0032     LogTrace("RPCHitAssociator") << "... size = " << RPCsimhits->size();
0033 
0034     //   MixCollection<PSimHit> & simHits = *hits;
0035 
0036     for (MixCollection<PSimHit>::MixItr hitItr = RPCsimhits->begin(); hitItr != RPCsimhits->end(); ++hitItr) {
0037       _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
0038     }
0039 
0040   } else if (!theConfig.RPCsimhitsTag.label().empty()) {
0041     LogTrace("RPCHitAssociator") << "getting PSimHit collection - " << theConfig.RPCsimhitsTag;
0042     edm::PSimHitContainer const &RPCsimhits = e.get(theConfig.RPCsimhitsToken_);
0043     LogTrace("RPCHitAssociator") << "... size = " << RPCsimhits.size();
0044 
0045     // arrange the hits by detUnit
0046     for (edm::PSimHitContainer::const_iterator hitItr = RPCsimhits.begin(); hitItr != RPCsimhits.end(); ++hitItr) {
0047       _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
0048     }
0049   }
0050 
0051   LogTrace("RPCHitAssociator") << "getting RPCDigiSimLink collection - " << theConfig.RPCdigisimlinkTag;
0052   _thelinkDigis = e.getHandle(theConfig.RPCdigisimlinkToken_);
0053 }
0054 // end of constructor
0055 
0056 std::vector<RPCHitAssociator::SimHitIdpr> RPCHitAssociator::associateRecHit(const TrackingRecHit &hit) const {
0057   std::vector<SimHitIdpr> matched;
0058 
0059   const TrackingRecHit *hitp = &hit;
0060   const RPCRecHit *rpcrechit = dynamic_cast<const RPCRecHit *>(hitp);
0061 
0062   if (rpcrechit) {
0063     RPCDetId rpcDetId = rpcrechit->rpcId();
0064     int fstrip = rpcrechit->firstClusterStrip();
0065     int cls = rpcrechit->clusterSize();
0066     int bx = rpcrechit->BunchX();
0067 
0068     for (int i = fstrip; i < fstrip + cls; ++i) {
0069       std::set<RPCDigiSimLink> links = findRPCDigiSimLink(rpcDetId.rawId(), i, bx);
0070 
0071       if (links.empty())
0072         LogTrace("RPCHitAssociator") << "*** WARNING in RPCHitAssociator::associateRecHit, RPCRecHit " << *rpcrechit
0073                                      << ", strip " << i << " has no associated RPCDigiSimLink !" << endl;
0074 
0075       for (std::set<RPCDigiSimLink>::iterator itlink = links.begin(); itlink != links.end(); ++itlink) {
0076         SimHitIdpr currentId(itlink->getTrackId(), itlink->getEventId());
0077         if (find(matched.begin(), matched.end(), currentId) == matched.end())
0078           matched.push_back(currentId);
0079       }
0080     }
0081 
0082   } else
0083     LogTrace("RPCHitAssociator") << "*** WARNING in RPCHitAssociator::associateRecHit, null "
0084                                     "dynamic_cast !";
0085 
0086   return matched;
0087 }
0088 
0089 std::set<RPCDigiSimLink> RPCHitAssociator::findRPCDigiSimLink(uint32_t rpcDetId, int strip, int bx) const {
0090   std::set<RPCDigiSimLink> links;
0091 
0092   for (edm::DetSetVector<RPCDigiSimLink>::const_iterator itlink = _thelinkDigis->begin();
0093        itlink != _thelinkDigis->end();
0094        itlink++) {
0095     for (edm::DetSet<RPCDigiSimLink>::const_iterator digi_iter = itlink->data.begin(); digi_iter != itlink->data.end();
0096          ++digi_iter) {
0097       uint32_t detid = digi_iter->getDetUnitId();
0098       int str = digi_iter->getStrip();
0099       int bunchx = digi_iter->getBx();
0100 
0101       if (detid == rpcDetId && str == strip && bunchx == bx) {
0102         links.insert(*digi_iter);
0103       }
0104     }
0105   }
0106 
0107   return links;
0108 }