CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCHitAssociator.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 
5 
6 
7 // Constructor
10  RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag")),
11  // CrossingFrame used or not ?
12  crossingframe(conf.getParameter<bool>("crossingframe")),
13  RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsTag")),
14  RPCsimhitsXFTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag"))
15 {
16  if (crossingframe){
18  } else if (!RPCsimhitsTag.label().empty()) {
20  }
21 
23 }
24 
26  RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag")),
27  // CrossingFrame used or not ?
28  crossingframe(conf.getParameter<bool>("crossingframe")),
29  RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsTag")),
30  RPCsimhitsXFTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag"))
31 {
32  initEvent(e,eventSetup);
33 }
34 
35 
37 
38 
39 {
40  if (crossingframe) {
41 
43  LogTrace("RPCHitAssociator") <<"getting CrossingFrame<PSimHit> collection - "<<RPCsimhitsXFTag;
45 
46  std::auto_ptr<MixCollection<PSimHit> >
47  RPCsimhits( new MixCollection<PSimHit>(cf.product()) );
48  LogTrace("RPCHitAssociator") <<"... size = "<<RPCsimhits->size();
49 
50  // MixCollection<PSimHit> & simHits = *hits;
51 
52  for(MixCollection<PSimHit>::MixItr hitItr = RPCsimhits->begin();
53  hitItr != RPCsimhits->end(); ++hitItr)
54  {
55  _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
56  }
57 
58  } else if (!RPCsimhitsTag.label().empty()) {
60  LogTrace("RPCHitAssociator") <<"getting PSimHit collection - "<<RPCsimhitsTag;
61  e.getByLabel(RPCsimhitsTag, RPCsimhits);
62  LogTrace("RPCHitAssociator") <<"... size = "<<RPCsimhits->size();
63 
64  // arrange the hits by detUnit
65  for(edm::PSimHitContainer::const_iterator hitItr = RPCsimhits->begin();
66  hitItr != RPCsimhits->end(); ++hitItr)
67  {
68  _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
69  }
70  }
71 
73  LogTrace("RPCHitAssociator") <<"getting RPCDigiSimLink collection - "<<RPCdigisimlinkTag;
74  e.getByLabel(RPCdigisimlinkTag, thelinkDigis);
75  _thelinkDigis = thelinkDigis;
76 }
77 // end of constructor
78 
79 std::vector<RPCHitAssociator::SimHitIdpr> RPCHitAssociator::associateRecHit(const TrackingRecHit & hit) const {
80 
81  std::vector<SimHitIdpr> matched;
82 
83  const TrackingRecHit * hitp = &hit;
84  const RPCRecHit * rpcrechit = dynamic_cast<const RPCRecHit *>(hitp);
85 
86  if (rpcrechit) {
87 
88  RPCDetId rpcDetId = rpcrechit->rpcId();
89  int fstrip = rpcrechit->firstClusterStrip();
90  int cls = rpcrechit->clusterSize();
91  int bx = rpcrechit->BunchX();
92 
93  for(int i = fstrip; i < fstrip+cls; ++i) {
94  std::set<RPCDigiSimLink> links = findRPCDigiSimLink(rpcDetId.rawId(),i,bx);
95 
96  if (links.empty()) edm::LogWarning("RPCHitAssociator")
97  <<"*** WARNING in RPCHitAssociator::associateRecHit, RPCRecHit "<<*rpcrechit<<", strip "<<i<<" has no associated RPCDigiSimLink !"<<endl;
98 
99  for(std::set<RPCDigiSimLink>::iterator itlink = links.begin(); itlink != links.end(); ++itlink) {
100  SimHitIdpr currentId(itlink->getTrackId(),itlink->getEventId());
101  if(find(matched.begin(),matched.end(),currentId ) == matched.end())
102  matched.push_back(currentId);
103  }
104  }
105 
106  } else edm::LogWarning("RPCHitAssociator")<<"*** WARNING in RPCHitAssociator::associateRecHit, null dynamic_cast !";
107 
108  return matched;
109 }
110 
111 std::set<RPCDigiSimLink> RPCHitAssociator::findRPCDigiSimLink(uint32_t rpcDetId, int strip, int bx) const {
112 
113  std::set<RPCDigiSimLink> links;
114 
116  for(edm::DetSet<RPCDigiSimLink>::const_iterator digi_iter=itlink->data.begin();digi_iter != itlink->data.end();++digi_iter){
117 
118  uint32_t detid = digi_iter->getDetUnitId();
119  int str = digi_iter->getStrip();
120  int bunchx = digi_iter->getBx();
121 
122  if(detid == rpcDetId && str == strip && bunchx == bx){
123  links.insert(*digi_iter);
124  }
125 
126  }
127  }
128 
129  return links;
130 }
131 
132 
int firstClusterStrip() const
Definition: RPCRecHit.h:105
int i
Definition: DBlmapReader.cc:9
RPCHitAssociator(const edm::ParameterSet &, edm::ConsumesCollector &&ic)
edm::InputTag RPCdigisimlinkTag
edm::EDGetTokenT< CrossingFrame< PSimHit > > RPCsimhitsXFToken_
void initEvent(const edm::Event &, const edm::EventSetup &)
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
std::vector< SimHitIdpr > associateRecHit(const TrackingRecHit &hit) const
int clusterSize() const
Definition: RPCRecHit.h:109
RPCDetId rpcId() const
Return the rpcId.
Definition: RPCRecHit.h:97
std::map< unsigned int, edm::PSimHitContainer > _SimHitMap
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:365
#define LogTrace(id)
tuple conf
Definition: dbtoconf.py:185
edm::InputTag RPCsimhitsTag
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< edm::PSimHitContainer > RPCsimhitsToken_
edm::InputTag RPCsimhitsXFTag
std::pair< uint32_t, EncodedEventId > SimHitIdpr
std::string const & label() const
Definition: InputTag.h:43
std::set< RPCDigiSimLink > findRPCDigiSimLink(uint32_t rpcDetId, int strip, int bx) const
int BunchX() const
Definition: RPCRecHit.h:101
std::vector< PSimHit > PSimHitContainer
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:350
edm::EDGetTokenT< edm::DetSetVector< RPCDigiSimLink > > RPCdigisimlinkToken_
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:108
edm::Handle< edm::DetSetVector< RPCDigiSimLink > > _thelinkDigis