#include <RPCHitAssociator.h>
Public Types | |
typedef edm::DetSetVector < RPCDigiSimLink > | RPCDigiSimLinks |
typedef std::pair< uint32_t, EncodedEventId > | SimHitIdpr |
Public Member Functions | |
std::vector< SimHitIdpr > | associateRecHit (const TrackingRecHit &hit) |
std::set< RPCDigiSimLink > | findRPCDigiSimLink (uint32_t rpcDetId, int strip, int bx) |
RPCHitAssociator (const edm::Event &, const edm::EventSetup &, const edm::ParameterSet &) | |
~RPCHitAssociator () | |
Private Attributes | |
std::map< unsigned int, edm::PSimHitContainer > | _SimHitMap |
edm::Handle< edm::DetSetVector < RPCDigiSimLink > > | _thelinkDigis |
bool | crossingframe |
edm::InputTag | RPCdigisimlinkTag |
edm::InputTag | RPCsimhitsTag |
edm::InputTag | RPCsimhitsXFTag |
Definition at line 31 of file RPCHitAssociator.h.
Definition at line 34 of file RPCHitAssociator.h.
typedef std::pair<uint32_t, EncodedEventId> RPCHitAssociator::SimHitIdpr |
Definition at line 35 of file RPCHitAssociator.h.
RPCHitAssociator::RPCHitAssociator | ( | const edm::Event & | e, |
const edm::EventSetup & | eventSetup, | ||
const edm::ParameterSet & | conf | ||
) |
Definition at line 6 of file RPCHitAssociator.cc.
References _SimHitMap, _thelinkDigis, MixCollection< T >::begin(), crossingframe, edm::Event::getByLabel(), edm::InputTag::label(), LogTrace, edm::Handle< T >::product(), RPCdigisimlinkTag, RPCsimhitsTag, and RPCsimhitsXFTag.
: RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag")), // CrossingFrame used or not ? crossingframe(conf.getParameter<bool>("crossingframe")), RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsTag")), RPCsimhitsXFTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag")) { if (crossingframe) { edm::Handle<CrossingFrame<PSimHit> > cf; LogTrace("RPCHitAssociator") <<"getting CrossingFrame<PSimHit> collection - "<<RPCsimhitsXFTag; e.getByLabel(RPCsimhitsXFTag, cf); std::auto_ptr<MixCollection<PSimHit> > RPCsimhits( new MixCollection<PSimHit>(cf.product()) ); LogTrace("RPCHitAssociator") <<"... size = "<<RPCsimhits->size(); // MixCollection<PSimHit> & simHits = *hits; for(MixCollection<PSimHit>::MixItr hitItr = RPCsimhits->begin(); hitItr != RPCsimhits->end(); ++hitItr) { _SimHitMap[hitItr->detUnitId()].push_back(*hitItr); } } else if (!RPCsimhitsTag.label().empty()) { edm::Handle<edm::PSimHitContainer> RPCsimhits; LogTrace("RPCHitAssociator") <<"getting PSimHit collection - "<<RPCsimhitsTag; e.getByLabel(RPCsimhitsTag, RPCsimhits); LogTrace("RPCHitAssociator") <<"... size = "<<RPCsimhits->size(); // arrange the hits by detUnit for(edm::PSimHitContainer::const_iterator hitItr = RPCsimhits->begin(); hitItr != RPCsimhits->end(); ++hitItr) { _SimHitMap[hitItr->detUnitId()].push_back(*hitItr); } } edm::Handle< edm::DetSetVector<RPCDigiSimLink> > thelinkDigis; LogTrace("RPCHitAssociator") <<"getting RPCDigiSimLink collection - "<<RPCdigisimlinkTag; e.getByLabel(RPCdigisimlinkTag, thelinkDigis); _thelinkDigis = thelinkDigis; }
RPCHitAssociator::~RPCHitAssociator | ( | ) | [inline] |
Definition at line 41 of file RPCHitAssociator.h.
{}
std::vector< RPCHitAssociator::SimHitIdpr > RPCHitAssociator::associateRecHit | ( | const TrackingRecHit & | hit | ) |
Definition at line 54 of file RPCHitAssociator.cc.
References RPCRecHit::BunchX(), RPCRecHit::clusterSize(), spr::find(), findRPCDigiSimLink(), RPCRecHit::firstClusterStrip(), i, DetId::rawId(), and RPCRecHit::rpcId().
Referenced by MuonAssociatorByHits::getMatchedIds().
{ std::vector<SimHitIdpr> matched; const TrackingRecHit * hitp = &hit; const RPCRecHit * rpcrechit = dynamic_cast<const RPCRecHit *>(hitp); if (rpcrechit) { RPCDetId rpcDetId = rpcrechit->rpcId(); int fstrip = rpcrechit->firstClusterStrip(); int cls = rpcrechit->clusterSize(); int bx = rpcrechit->BunchX(); for(int i = fstrip; i < fstrip+cls; ++i) { std::set<RPCDigiSimLink> links = findRPCDigiSimLink(rpcDetId.rawId(),i,bx); if (links.empty()) edm::LogWarning("RPCHitAssociator") <<"*** WARNING in RPCHitAssociator::associateRecHit, RPCRecHit "<<*rpcrechit<<", strip "<<i<<" has no associated RPCDigiSimLink !"<<endl; for(std::set<RPCDigiSimLink>::iterator itlink = links.begin(); itlink != links.end(); ++itlink) { SimHitIdpr currentId(itlink->getTrackId(),itlink->getEventId()); if(find(matched.begin(),matched.end(),currentId ) == matched.end()) matched.push_back(currentId); } } } else edm::LogWarning("RPCHitAssociator")<<"*** WARNING in RPCHitAssociator::associateRecHit, null dynamic_cast !"; return matched; }
std::set< RPCDigiSimLink > RPCHitAssociator::findRPCDigiSimLink | ( | uint32_t | rpcDetId, |
int | strip, | ||
int | bx | ||
) |
Definition at line 86 of file RPCHitAssociator.cc.
References _thelinkDigis, and cond::rpcobgas::detid.
Referenced by associateRecHit().
{ std::set<RPCDigiSimLink> links; for (edm::DetSetVector<RPCDigiSimLink>::const_iterator itlink = _thelinkDigis->begin(); itlink != _thelinkDigis->end(); itlink++){ for(edm::DetSet<RPCDigiSimLink>::const_iterator digi_iter=itlink->data.begin();digi_iter != itlink->data.end();++digi_iter){ uint32_t detid = digi_iter->getDetUnitId(); int str = digi_iter->getStrip(); int bunchx = digi_iter->getBx(); if(detid == rpcDetId && str == strip && bunchx == bx){ links.insert(*digi_iter); } } } return links; }
std::map<unsigned int, edm::PSimHitContainer> RPCHitAssociator::_SimHitMap [private] |
Definition at line 56 of file RPCHitAssociator.h.
Referenced by RPCHitAssociator().
Definition at line 49 of file RPCHitAssociator.h.
Referenced by findRPCDigiSimLink(), and RPCHitAssociator().
bool RPCHitAssociator::crossingframe [private] |
Definition at line 52 of file RPCHitAssociator.h.
Referenced by RPCHitAssociator().
Definition at line 50 of file RPCHitAssociator.h.
Referenced by RPCHitAssociator().
edm::InputTag RPCHitAssociator::RPCsimhitsTag [private] |
Definition at line 53 of file RPCHitAssociator.h.
Referenced by RPCHitAssociator().
Definition at line 54 of file RPCHitAssociator.h.
Referenced by RPCHitAssociator().