CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
RPCPointVsRecHit Class Reference

#include <RPCPointVsRecHit.h>

Inheritance diagram for RPCPointVsRecHit:
edm::EDAnalyzer

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 
void beginJob ()
 
void endJob ()
 
 RPCPointVsRecHit (const edm::ParameterSet &pset)
 
 ~RPCPointVsRecHit ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Attributes

DQMStoredbe_
 
RPCValidHistograms h_
 
edm::InputTag recHitLabel_
 
edm::InputTag refHitLabel_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 18 of file RPCPointVsRecHit.h.

Constructor & Destructor Documentation

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

Definition at line 16 of file RPCPointVsRecHit.cc.

References dbe_, edm::ParameterSet::getParameter(), and cppFunctionSkipper::operator.

17 {
18  refHitLabel_ = pset.getParameter<edm::InputTag>("refHit");
19  recHitLabel_ = pset.getParameter<edm::InputTag>("recHit");
20 
22  if ( !dbe_ )
23  {
24  edm::LogError("RPCPointVsRecHit") << "No DQMStore instance\n";
25  return;
26  }
27 
28  // Book MonitorElements
29  const std::string subDir = pset.getParameter<std::string>("subDir");
30  h_.bookHistograms(dbe_, subDir);
31 }
T getParameter(std::string const &) const
void bookHistograms(DQMStore *dbe, const std::string subDir)
edm::InputTag refHitLabel_
RPCValidHistograms h_
edm::InputTag recHitLabel_
RPCPointVsRecHit::~RPCPointVsRecHit ( )

Definition at line 33 of file RPCPointVsRecHit.cc.

34 {
35 }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 45 of file RPCPointVsRecHit.cc.

References dbe_, edm::EventSetup::get(), edm::Event::getByLabel(), RPCRoll::id(), match(), RPCDetId::region(), relativeConstraints::ring, RPCDetId::ring(), mathSSE::sqrt(), relativeConstraints::station, and RPCDetId::station().

46 {
47  if ( !dbe_ )
48  {
49  edm::LogError("RPCPointVsRecHit") << "No DQMStore instance\n";
50  return;
51  }
52 
53  // Get the RPC Geometry
55  eventSetup.get<MuonGeometryRecord>().get(rpcGeom);
56 
57  // Retrieve RefHits from the event
59  if ( !event.getByLabel(refHitLabel_, refHitHandle) )
60  {
61  edm::LogInfo("RPCPointVsRecHit") << "Cannot find reference hit collection\n";
62  return;
63  }
64 
65  // Retrieve RecHits from the event
67  if ( !event.getByLabel(recHitLabel_, recHitHandle) )
68  {
69  edm::LogInfo("RPCPointVsRecHit") << "Cannot find recHit collection\n";
70  return;
71  }
72 
73  typedef RPCRecHitCollection::const_iterator RecHitIter;
74 
75  // Loop over refHits, fill histograms which does not need associations
76  int nRefHitBarrel = 0, nRefHitEndcap = 0;
77  for ( RecHitIter refHitIter = refHitHandle->begin();
78  refHitIter != refHitHandle->end(); ++refHitIter )
79  {
80  const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
81  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
82  if ( !roll ) continue;
83 
84  const int region = roll->id().region();
85  const int ring = roll->id().ring();
86  //const int sector = roll->id().sector();
87  const int station = roll->id().station();
88  //const int layer = roll->id().layer();
89  //const int subSector = roll->id().subsector();
90 
91  if ( region == 0 )
92  {
96  }
97  else
98  {
99  h_.refHitOccupancyEndcap_disk->Fill(region*station);
100  h_.refHitOccupancyEndcap_disk_ring->Fill(region*station, ring);
101  }
102 
103  }
104  h_.nRefHitBarrel->Fill(nRefHitBarrel);
105  h_.nRefHitEndcap->Fill(nRefHitEndcap);
106 
107  // Loop over recHits, fill histograms which does not need associations
108  int sumClusterSizeBarrel = 0, sumClusterSizeEndcap = 0;
109  int nRecHitBarrel = 0, nRecHitEndcap = 0;
110  for ( RecHitIter recHitIter = recHitHandle->begin();
111  recHitIter != recHitHandle->end(); ++recHitIter )
112  {
113  const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
114  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
115  if ( !roll ) continue;
116 
117  const int region = roll->id().region();
118  const int ring = roll->id().ring();
119  //const int sector = roll->id().sector();
120  const int station = (roll->id().station());
121  //const int layer = roll->id().layer();
122  //const int subSector = roll->id().subsector();
123 
124  h_.clusterSize->Fill(recHitIter->clusterSize());
125 
126  if ( region == 0 )
127  {
128  ++nRecHitBarrel;
129  sumClusterSizeBarrel += recHitIter->clusterSize();
130  h_.clusterSizeBarrel->Fill(recHitIter->clusterSize());
134  }
135  else
136  {
137  ++nRecHitEndcap;
138  sumClusterSizeEndcap += recHitIter->clusterSize();
139  h_.clusterSizeEndcap->Fill(recHitIter->clusterSize());
140  h_.recHitOccupancyEndcap_disk->Fill(region*station);
141  h_.recHitOccupancyEndcap_disk_ring->Fill(ring, region*station);
142  }
143 
144  }
145  const double nRecHit = nRecHitBarrel+nRecHitEndcap;
146  h_.nRecHitBarrel->Fill(nRecHitBarrel);
147  h_.nRecHitEndcap->Fill(nRecHitEndcap);
148  if ( nRecHit > 0 )
149  {
150  const int sumClusterSize = sumClusterSizeBarrel+sumClusterSizeEndcap;
151  h_.avgClusterSize->Fill(double(sumClusterSize)/nRecHit);
152 
153  if ( nRecHitBarrel > 0 )
154  {
155  h_.avgClusterSizeBarrel->Fill(double(sumClusterSizeBarrel)/nRecHitBarrel);
156  }
157  if ( nRecHitEndcap > 0 )
158  {
159  h_.avgClusterSizeEndcap->Fill(double(sumClusterSizeEndcap)/nRecHitEndcap);
160  }
161  }
162 
163  // Start matching RefHits to RecHits
164  typedef std::map<RecHitIter, RecHitIter> RecToRecHitMap;
165  RecToRecHitMap refToRecHitMap;
166 
167  for ( RecHitIter refHitIter = refHitHandle->begin();
168  refHitIter != refHitHandle->end(); ++refHitIter )
169  {
170  const RPCDetId refDetId = static_cast<const RPCDetId>(refHitIter->rpcId());
171  const RPCRoll* refRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(refDetId));
172  if ( !refRoll ) continue;
173 
174  const double refX = refHitIter->localPosition().x();
175 
176  for ( RecHitIter recHitIter = recHitHandle->begin();
177  recHitIter != recHitHandle->end(); ++recHitIter )
178  {
179  const RPCDetId recDetId = static_cast<const RPCDetId>(recHitIter->rpcId());
180  const RPCRoll* recRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(recDetId));
181  if ( !recRoll ) continue;
182 
183  if ( refDetId != recDetId ) continue;
184 
185  const double recX = recHitIter->localPosition().x();
186  const double newDx = fabs(recX - refX);
187 
188  // Associate RefHit to RecHit
189  RecToRecHitMap::const_iterator prevRefToReco = refToRecHitMap.find(refHitIter);
190  if ( prevRefToReco == refToRecHitMap.end() )
191  {
192  refToRecHitMap.insert(std::make_pair(refHitIter, recHitIter));
193  }
194  else
195  {
196  const double oldDx = fabs(prevRefToReco->second->localPosition().x() - refX);
197 
198  if ( newDx < oldDx )
199  {
200  refToRecHitMap[refHitIter] = recHitIter;
201  }
202  }
203  }
204  }
205 
206  // Now we have refHit-recHit mapping
207  // So we can fill up relavant histograms
208  for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
209  match != refToRecHitMap.end(); ++match )
210  {
211  RecHitIter refHitIter = match->first;
212  RecHitIter recHitIter = match->second;
213 
214  const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
215  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
216 
217  const int region = roll->id().region();
218  const int ring = roll->id().ring();
219  //const int sector = roll->id().sector();
220  const int station = roll->id().station();
221  //const int layer = roll->id().layer();
222  //const int subsector = roll->id().subsector();
223 
224  const double refX = refHitIter->localPosition().x();
225  const double recX = recHitIter->localPosition().x();
226  const double errX = sqrt(recHitIter->localPositionError().xx());
227  const double dX = recX - refX;
228  const double pull = errX == 0 ? -999 : dX/errX;
229 
230  //const GlobalPoint refPos = roll->toGlobal(refHitIter->localPosition());
231  //const GlobalPoint recPos = roll->toGlobal(recHitIter->localPosition());
232 
233  if ( region == 0 )
234  {
235  h_.resBarrel->Fill(dX);
236  h_.pullBarrel->Fill(pull);
240 
241  h_.res_wheel_res->Fill(ring, dX);
242  h_.res_station_res->Fill(station, dX);
243  h_.pull_wheel_pull->Fill(ring, pull);
244  h_.pull_station_pull->Fill(station, pull);
245  }
246  else
247  {
248  h_.resEndcap->Fill(dX);
249  h_.pullEndcap->Fill(pull);
250  h_.matchOccupancyEndcap_disk->Fill(region*station);
251  h_.matchOccupancyEndcap_disk_ring->Fill(region*station, ring);
252 
253  h_.res_disk_res->Fill(region*station, dX);
254  h_.res_ring_res->Fill(ring, dX);
255  h_.pull_disk_pull->Fill(region*station, pull);
256  h_.pull_ring_pull->Fill(ring, pull);
257  }
258  }
259 
260 /*
261  // Find Lost hits
262  for ( RecHitIter refHitIter = refHitHandle->begin();
263  refHitIter != refHitHandle->end(); ++refHitIter )
264  {
265  const RPCDetId detId = static_cast<const RPCDetId>(refHitIter->rpcId());
266  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
267 
268  const int region = roll->id().region();
269  const int ring = roll->id().ring();
270  //const int sector = roll->id().sector();
271  const int station = roll->id().station();
272  //const int layer = roll->id().layer();
273  //const int subsector = roll->id().subsector();
274 
275  bool matched = false;
276  for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
277  match != refToRecHitMap.end(); ++match )
278  {
279  if ( refHitIter == match->first )
280  {
281  matched = true;
282  break;
283  }
284  }
285 
286  if ( !matched )
287  {
288  if ( region == 0 )
289  {
290  h_.nUrefHitOccupancyBarrel_wheel->Fill(ring);
291  h_.nUrefHitOccupancyBarrel_wheel_ring->Fill(ring, station);
292  }
293  else
294  {
295  h_.nUnMatchedRefHit_disk->Fill(region*station);
296  h_.nUnMatchedRefHit_disk_ring->Fill(region*station, ring);
297  }
298  }
299  }
300 */
301 
302  // Find Noisy hits
303  for ( RecHitIter recHitIter = recHitHandle->begin();
304  recHitIter != recHitHandle->end(); ++recHitIter )
305  {
306  const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
307  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
308 
309  const int region = roll->id().region();
310  const int ring = roll->id().ring();
311  //const int sector = roll->id().sector();
312  const int station = roll->id().station();
313  //const int layer = roll->id().layer();
314  //const int subsector = roll->id().subsector();
315 
316  bool matched = false;
317  for ( RecToRecHitMap::const_iterator match = refToRecHitMap.begin();
318  match != refToRecHitMap.end(); ++match )
319  {
320  if ( recHitIter == match->second )
321  {
322  matched = true;
323  break;
324  }
325  }
326 
327  if ( !matched )
328  {
329  if ( region == 0 )
330  {
333  h_.umOccupancyBarrel_wheel_station->Fill(ring, station);
334  }
335  else
336  {
337  h_.umOccupancyEndcap_disk->Fill(region*station);
338  h_.umOccupancyEndcap_disk_ring->Fill(region*station, ring);
339  }
340 
341  //const GlobalPoint pos = roll->toGlobal(recHitIter->localPosition());
342  //h_[HName::NoisyHitEta]->Fill(pos.eta());
343  }
344  }
345 }
void Fill(long long x)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
RPCDetId id() const
Definition: RPCRoll.cc:24
int ring() const
Definition: RPCDetId.h:76
T sqrt(T t)
Definition: SSEVec.h:46
edm::InputTag refHitLabel_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
const T & get() const
Definition: EventSetup.h:55
RPCValidHistograms h_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
edm::InputTag recHitLabel_
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:67
int station() const
Definition: RPCDetId.h:100
void RPCPointVsRecHit::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 37 of file RPCPointVsRecHit.cc.

38 {
39 }
void RPCPointVsRecHit::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 41 of file RPCPointVsRecHit.cc.

42 {
43 }

Member Data Documentation

DQMStore* RPCPointVsRecHit::dbe_
private

Definition at line 31 of file RPCPointVsRecHit.h.

RPCValidHistograms RPCPointVsRecHit::h_
private

Definition at line 33 of file RPCPointVsRecHit.h.

edm::InputTag RPCPointVsRecHit::recHitLabel_
private

Definition at line 29 of file RPCPointVsRecHit.h.

edm::InputTag RPCPointVsRecHit::refHitLabel_
private

Definition at line 29 of file RPCPointVsRecHit.h.