CMS 3D CMS Logo

RPCHitAssociator.cc
Go to the documentation of this file.
2 
3 using namespace std;
4 
5 // Constructor
7  : RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag")),
8  // CrossingFrame used or not ?
9  crossingframe(conf.getParameter<bool>("crossingframe")),
10  RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsTag")),
11  RPCsimhitsXFTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag")) {
12  if (crossingframe) {
14  } else if (!RPCsimhitsTag.label().empty()) {
16  }
17 
19 }
20 
22  const edm::EventSetup &eventSetup,
23  const edm::ParameterSet &conf)
24  : RPCdigisimlinkTag(conf.getParameter<edm::InputTag>("RPCdigisimlinkTag")),
25  // CrossingFrame used or not ?
26  crossingframe(conf.getParameter<bool>("crossingframe")),
27  RPCsimhitsTag(conf.getParameter<edm::InputTag>("RPCsimhitsTag")),
28  RPCsimhitsXFTag(conf.getParameter<edm::InputTag>("RPCsimhitsXFTag")) {
29  initEvent(e, eventSetup);
30 }
31 
33 
34 {
35  if (crossingframe) {
37  LogTrace("RPCHitAssociator") << "getting CrossingFrame<PSimHit> collection - " << RPCsimhitsXFTag;
39 
40  std::unique_ptr<MixCollection<PSimHit>> RPCsimhits(new MixCollection<PSimHit>(cf.product()));
41  LogTrace("RPCHitAssociator") << "... size = " << RPCsimhits->size();
42 
43  // MixCollection<PSimHit> & simHits = *hits;
44 
45  for (MixCollection<PSimHit>::MixItr hitItr = RPCsimhits->begin(); hitItr != RPCsimhits->end(); ++hitItr) {
46  _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
47  }
48 
49  } else if (!RPCsimhitsTag.label().empty()) {
51  LogTrace("RPCHitAssociator") << "getting PSimHit collection - " << RPCsimhitsTag;
52  e.getByLabel(RPCsimhitsTag, RPCsimhits);
53  LogTrace("RPCHitAssociator") << "... size = " << RPCsimhits->size();
54 
55  // arrange the hits by detUnit
56  for (edm::PSimHitContainer::const_iterator hitItr = RPCsimhits->begin(); hitItr != RPCsimhits->end(); ++hitItr) {
57  _SimHitMap[hitItr->detUnitId()].push_back(*hitItr);
58  }
59  }
60 
62  LogTrace("RPCHitAssociator") << "getting RPCDigiSimLink collection - " << RPCdigisimlinkTag;
63  e.getByLabel(RPCdigisimlinkTag, thelinkDigis);
64  _thelinkDigis = thelinkDigis;
65 }
66 // end of constructor
67 
68 std::vector<RPCHitAssociator::SimHitIdpr> RPCHitAssociator::associateRecHit(const TrackingRecHit &hit) const {
69  std::vector<SimHitIdpr> matched;
70 
71  const TrackingRecHit *hitp = &hit;
72  const RPCRecHit *rpcrechit = dynamic_cast<const RPCRecHit *>(hitp);
73 
74  if (rpcrechit) {
75  RPCDetId rpcDetId = rpcrechit->rpcId();
76  int fstrip = rpcrechit->firstClusterStrip();
77  int cls = rpcrechit->clusterSize();
78  int bx = rpcrechit->BunchX();
79 
80  for (int i = fstrip; i < fstrip + cls; ++i) {
81  std::set<RPCDigiSimLink> links = findRPCDigiSimLink(rpcDetId.rawId(), i, bx);
82 
83  if (links.empty())
84  LogTrace("RPCHitAssociator") << "*** WARNING in RPCHitAssociator::associateRecHit, RPCRecHit " << *rpcrechit
85  << ", strip " << i << " has no associated RPCDigiSimLink !" << endl;
86 
87  for (std::set<RPCDigiSimLink>::iterator itlink = links.begin(); itlink != links.end(); ++itlink) {
88  SimHitIdpr currentId(itlink->getTrackId(), itlink->getEventId());
89  if (find(matched.begin(), matched.end(), currentId) == matched.end())
90  matched.push_back(currentId);
91  }
92  }
93 
94  } else
95  LogTrace("RPCHitAssociator") << "*** WARNING in RPCHitAssociator::associateRecHit, null "
96  "dynamic_cast !";
97 
98  return matched;
99 }
100 
101 std::set<RPCDigiSimLink> RPCHitAssociator::findRPCDigiSimLink(uint32_t rpcDetId, int strip, int bx) const {
102  std::set<RPCDigiSimLink> links;
103 
105  itlink != _thelinkDigis->end();
106  itlink++) {
107  for (edm::DetSet<RPCDigiSimLink>::const_iterator digi_iter = itlink->data.begin(); digi_iter != itlink->data.end();
108  ++digi_iter) {
109  uint32_t detid = digi_iter->getDetUnitId();
110  int str = digi_iter->getStrip();
111  int bunchx = digi_iter->getBx();
112 
113  if (detid == rpcDetId && str == strip && bunchx == bx) {
114  links.insert(*digi_iter);
115  }
116  }
117  }
118 
119  return links;
120 }
int firstClusterStrip() const
Definition: RPCRecHit.h:111
RPCHitAssociator(const edm::ParameterSet &, edm::ConsumesCollector &&ic)
edm::InputTag RPCdigisimlinkTag
edm::EDGetTokenT< CrossingFrame< PSimHit > > RPCsimhitsXFToken_
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
void initEvent(const edm::Event &, const edm::EventSetup &)
edm::Handle< edm::DetSetVector< RPCDigiSimLink > > _thelinkDigis
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< SimHitIdpr > associateRecHit(const TrackingRecHit &hit) const
int clusterSize() const
Definition: RPCRecHit.h:115
RPCDetId rpcId() const
Return the rpcId.
Definition: RPCRecHit.h:103
std::map< unsigned int, edm::PSimHitContainer > _SimHitMap
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:480
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:361
#define LogTrace(id)
edm::InputTag RPCsimhitsTag
T const * product() const
Definition: Handle.h:74
edm::EDGetTokenT< edm::PSimHitContainer > RPCsimhitsToken_
edm::InputTag RPCsimhitsXFTag
std::pair< uint32_t, EncodedEventId > SimHitIdpr
std::string const & label() const
Definition: InputTag.h:36
std::set< RPCDigiSimLink > findRPCDigiSimLink(uint32_t rpcDetId, int strip, int bx) const
HLT enums.
int BunchX() const
Definition: RPCRecHit.h:107
std::vector< PSimHit > PSimHitContainer
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:346
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:104
#define str(s)