CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
RPCPointVsRecHit Class Reference

#include <RPCPointVsRecHit.h>

Inheritance diagram for RPCPointVsRecHit:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 RPCPointVsRecHit (const edm::ParameterSet &pset)
 
 ~RPCPointVsRecHit () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Private Attributes

RPCValidHistograms h_
 
edm::EDGetTokenT< RPCRecHitCollectionrecHitToken_
 
edm::EDGetTokenT< RPCRecHitCollectionrefHitToken_
 
std::string subDir_
 

Detailed Description

Definition at line 19 of file RPCPointVsRecHit.h.

Constructor & Destructor Documentation

RPCPointVsRecHit::RPCPointVsRecHit ( const edm::ParameterSet pset)

Definition at line 15 of file RPCPointVsRecHit.cc.

References edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

16 {
17  refHitToken_ = consumes<RPCRecHitCollection>(pset.getParameter<edm::InputTag>("refHit"));
18  recHitToken_ = consumes<RPCRecHitCollection>(pset.getParameter<edm::InputTag>("recHit"));
19 
20  subDir_ = pset.getParameter<std::string>("subDir");
21 }
T getParameter(std::string const &) const
edm::EDGetTokenT< RPCRecHitCollection > refHitToken_
std::string subDir_
edm::EDGetTokenT< RPCRecHitCollection > recHitToken_
RPCPointVsRecHit::~RPCPointVsRecHit ( )
inlineoverride

Definition at line 23 of file RPCPointVsRecHit.h.

References analyze(), and bookHistograms().

23 {};

Member Function Documentation

void RPCPointVsRecHit::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
override

Definition at line 23 of file RPCPointVsRecHit.cc.

References fftjetvertexadder_cfi::errX, edm::EventSetup::get(), edm::Event::getByToken(), RPCRoll::id(), match(), RPCDetId::region(), relativeConstraints::ring, RPCDetId::ring(), RPCGeometry::roll(), mathSSE::sqrt(), relativeConstraints::station, and RPCDetId::station().

Referenced by ~RPCPointVsRecHit().

24 {
25  // Get the RPC Geometry
27  eventSetup.get<MuonGeometryRecord>().get(rpcGeom);
28 
29  // Retrieve RefHits from the event
31  if ( !event.getByToken(refHitToken_, refHitHandle) )
32  {
33  edm::LogInfo("RPCPointVsRecHit") << "Cannot find reference hit collection\n";
34  return;
35  }
36 
37  // Retrieve RecHits from the event
39  if ( !event.getByToken(recHitToken_, recHitHandle) )
40  {
41  edm::LogInfo("RPCPointVsRecHit") << "Cannot find recHit collection\n";
42  return;
43  }
44 
45  typedef RPCRecHitCollection::const_iterator RecHitIter;
46 
47  // Loop over refHits, fill histograms which does not need associations
48  int nRefHitBarrel = 0, nRefHitEndcap = 0;
49  for ( RecHitIter refHitIter = refHitHandle->begin();
50  refHitIter != refHitHandle->end(); ++refHitIter )
51  {
52  const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
53  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
54  if ( !roll ) continue;
55 
56  const int region = roll->id().region();
57  const int ring = roll->id().ring();
58  //const int sector = roll->id().sector();
59  const int station = roll->id().station();
60  //const int layer = roll->id().layer();
61  //const int subSector = roll->id().subsector();
62 
63  if ( region == 0 )
64  {
68  }
69  else
70  {
71  h_.refHitOccupancyEndcap_disk->Fill(region*station);
72  h_.refHitOccupancyEndcap_disk_ring->Fill(region*station, ring);
73  }
74 
75  }
76  h_.nRefHitBarrel->Fill(nRefHitBarrel);
77  h_.nRefHitEndcap->Fill(nRefHitEndcap);
78 
79  // Loop over recHits, fill histograms which does not need associations
80  int sumClusterSizeBarrel = 0, sumClusterSizeEndcap = 0;
81  int nRecHitBarrel = 0, nRecHitEndcap = 0;
82  for ( RecHitIter recHitIter = recHitHandle->begin();
83  recHitIter != recHitHandle->end(); ++recHitIter )
84  {
85  const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
86  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
87  if ( !roll ) continue;
88 
89  const int region = roll->id().region();
90  const int ring = roll->id().ring();
91  //const int sector = roll->id().sector();
92  const int station = (roll->id().station());
93  //const int layer = roll->id().layer();
94  //const int subSector = roll->id().subsector();
95 
96  h_.clusterSize->Fill(recHitIter->clusterSize());
97 
98  if ( region == 0 )
99  {
100  ++nRecHitBarrel;
101  sumClusterSizeBarrel += recHitIter->clusterSize();
102  h_.clusterSizeBarrel->Fill(recHitIter->clusterSize());
106  }
107  else
108  {
109  ++nRecHitEndcap;
110  sumClusterSizeEndcap += recHitIter->clusterSize();
111  h_.clusterSizeEndcap->Fill(recHitIter->clusterSize());
112  h_.recHitOccupancyEndcap_disk->Fill(region*station);
113  h_.recHitOccupancyEndcap_disk_ring->Fill(ring, region*station);
114  }
115 
116  }
117  const double nRecHit = nRecHitBarrel+nRecHitEndcap;
118  h_.nRecHitBarrel->Fill(nRecHitBarrel);
119  h_.nRecHitEndcap->Fill(nRecHitEndcap);
120  if ( nRecHit > 0 )
121  {
122  const int sumClusterSize = sumClusterSizeBarrel+sumClusterSizeEndcap;
123  h_.avgClusterSize->Fill(double(sumClusterSize)/nRecHit);
124 
125  if ( nRecHitBarrel > 0 )
126  {
127  h_.avgClusterSizeBarrel->Fill(double(sumClusterSizeBarrel)/nRecHitBarrel);
128  }
129  if ( nRecHitEndcap > 0 )
130  {
131  h_.avgClusterSizeEndcap->Fill(double(sumClusterSizeEndcap)/nRecHitEndcap);
132  }
133  }
134 
135  // Start matching RefHits to RecHits
136  typedef std::map<RecHitIter, RecHitIter> RecToRecHitMap;
137  RecToRecHitMap refToRecHitMap;
138 
139  for ( RecHitIter refHitIter = refHitHandle->begin();
140  refHitIter != refHitHandle->end(); ++refHitIter )
141  {
142  const RPCDetId refDetId = static_cast<const RPCDetId>(refHitIter->rpcId());
143  const RPCRoll* refRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(refDetId));
144  if ( !refRoll ) continue;
145 
146  const double refX = refHitIter->localPosition().x();
147 
148  for ( RecHitIter recHitIter = recHitHandle->begin();
149  recHitIter != recHitHandle->end(); ++recHitIter )
150  {
151  const RPCDetId recDetId = static_cast<const RPCDetId>(recHitIter->rpcId());
152  const RPCRoll* recRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(recDetId));
153  if ( !recRoll ) continue;
154 
155  if ( refDetId != recDetId ) continue;
156 
157  const double recX = recHitIter->localPosition().x();
158  const double newDx = fabs(recX - refX);
159 
160  // Associate RefHit to RecHit
161  RecToRecHitMap::const_iterator prevRefToReco = refToRecHitMap.find(refHitIter);
162  if ( prevRefToReco == refToRecHitMap.end() )
163  {
164  refToRecHitMap.insert(std::make_pair(refHitIter, recHitIter));
165  }
166  else
167  {
168  const double oldDx = fabs(prevRefToReco->second->localPosition().x() - refX);
169 
170  if ( newDx < oldDx )
171  {
172  refToRecHitMap[refHitIter] = recHitIter;
173  }
174  }
175  }
176  }
177 
178  // Now we have refHit-recHit mapping
179  // So we can fill up relavant histograms
180  for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
181  match != refToRecHitMap.end(); ++match )
182  {
183  RecHitIter refHitIter = match->first;
184  RecHitIter recHitIter = match->second;
185 
186  const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
187  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
188 
189  const int region = roll->id().region();
190  const int ring = roll->id().ring();
191  //const int sector = roll->id().sector();
192  const int station = roll->id().station();
193  //const int layer = roll->id().layer();
194  //const int subsector = roll->id().subsector();
195 
196  const double refX = refHitIter->localPosition().x();
197  const double recX = recHitIter->localPosition().x();
198  const double errX = sqrt(recHitIter->localPositionError().xx());
199  const double dX = recX - refX;
200  const double pull = errX == 0 ? -999 : dX/errX;
201 
202  //const GlobalPoint refPos = roll->toGlobal(refHitIter->localPosition());
203  //const GlobalPoint recPos = roll->toGlobal(recHitIter->localPosition());
204 
205  if ( region == 0 )
206  {
207  h_.resBarrel->Fill(dX);
208  h_.pullBarrel->Fill(pull);
212 
213  h_.res_wheel_res->Fill(ring, dX);
214  h_.res_station_res->Fill(station, dX);
215  h_.pull_wheel_pull->Fill(ring, pull);
216  h_.pull_station_pull->Fill(station, pull);
217  }
218  else
219  {
220  h_.resEndcap->Fill(dX);
221  h_.pullEndcap->Fill(pull);
222  h_.matchOccupancyEndcap_disk->Fill(region*station);
223  h_.matchOccupancyEndcap_disk_ring->Fill(region*station, ring);
224 
225  h_.res_disk_res->Fill(region*station, dX);
226  h_.res_ring_res->Fill(ring, dX);
227  h_.pull_disk_pull->Fill(region*station, pull);
228  h_.pull_ring_pull->Fill(ring, pull);
229  }
230  }
231 
232 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
edm::EDGetTokenT< RPCRecHitCollection > refHitToken_
void Fill(long long x)
RPCDetId id() const
Definition: RPCRoll.cc:24
int ring() const
Definition: RPCDetId.h:72
T sqrt(T t)
Definition: SSEVec.h:18
RPCValidHistograms h_
T get() const
Definition: EventSetup.h:68
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
const RPCRoll * roll(RPCDetId id) const
Return a roll given its id.
Definition: RPCGeometry.cc:75
edm::EDGetTokenT< RPCRecHitCollection > recHitToken_
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:63
int station() const
Definition: RPCDetId.h:96
void RPCPointVsRecHit::bookHistograms ( DQMStore::IBooker booker,
edm::Run const &  run,
edm::EventSetup const &  eventSetup 
)
override

Definition at line 234 of file RPCPointVsRecHit.cc.

References DEFINE_FWK_MODULE.

Referenced by ~RPCPointVsRecHit().

236 {
237  // Book MonitorElements
238  h_.bookHistograms(booker, subDir_);
239 }
void bookHistograms(DQMStore::IBooker &booker, const std::string &subDir)
std::string subDir_
RPCValidHistograms h_

Member Data Documentation

RPCValidHistograms RPCPointVsRecHit::h_
private

Definition at line 32 of file RPCPointVsRecHit.h.

edm::EDGetTokenT<RPCRecHitCollection> RPCPointVsRecHit::recHitToken_
private

Definition at line 29 of file RPCPointVsRecHit.h.

edm::EDGetTokenT<RPCRecHitCollection> RPCPointVsRecHit::refHitToken_
private

Definition at line 29 of file RPCPointVsRecHit.h.

std::string RPCPointVsRecHit::subDir_
private

Definition at line 31 of file RPCPointVsRecHit.h.