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
RPCRecHitValid Class Reference

#include <RPCRecHitValid.h>

Inheritance diagram for RPCRecHitValid:
edm::EDAnalyzer

Public Member Functions

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

Private Attributes

DQMStoredbe_
 
RPCValidHistograms h_
 
bool isStandAloneMode_
 
edm::InputTag recHitLabel_
 
std::string rootFileName_
 
edm::InputTag simHitLabel_
 

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)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 20 of file RPCRecHitValid.h.

Constructor & Destructor Documentation

RPCRecHitValid::RPCRecHitValid ( const edm::ParameterSet pset)

Definition at line 17 of file RPCRecHitValid.cc.

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

18 {
19  rootFileName_ = pset.getUntrackedParameter<string>("rootFileName", "");
20  simHitLabel_ = pset.getParameter<edm::InputTag>("simHit");
21  recHitLabel_ = pset.getParameter<edm::InputTag>("recHit");
22 
23  isStandAloneMode_ = pset.getUntrackedParameter<bool>("standAloneMode", false);
24 
26  if ( !dbe_ )
27  {
28  edm::LogError("RPCRecHitValid") << "No DQMStore instance\n";
29  return;
30  }
31 
32  // Book MonitorElements
33  const std::string subDir = pset.getParameter<std::string>("subDir");
34  h_.bookHistograms(dbe_, subDir);
35 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DQMStore * dbe_
void bookHistograms(DQMStore *dbe, const std::string subDir)
edm::InputTag recHitLabel_
RPCValidHistograms h_
edm::InputTag simHitLabel_
std::string rootFileName_
RPCRecHitValid::~RPCRecHitValid ( )

Definition at line 37 of file RPCRecHitValid.cc.

References dbe_, and DQMStore::save().

38 {
39  if ( dbe_ )
40  {
41  if ( !rootFileName_.empty() ) dbe_->save(rootFileName_);
42  }
43 }
DQMStore * dbe_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:1898
std::string rootFileName_

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 53 of file RPCRecHitValid.cc.

References abs, dbe_, edm::EventSetup::get(), edm::Event::getByLabel(), RPCRoll::id(), match(), pos, RPCDetId::region(), relativeConstraints::ring, RPCDetId::ring(), relativeConstraints::station, RPCDetId::station(), and GeomDet::toGlobal().

54 {
55  if ( !dbe_ )
56  {
57  edm::LogError("RPCRecHitValid") << "No DQMStore instance\n";
58  return;
59  }
60 
61  // Get the RPC Geometry
63  eventSetup.get<MuonGeometryRecord>().get(rpcGeom);
64 
65  // Retrieve SimHits from the event
67  if ( !event.getByLabel(simHitLabel_, simHitHandle) )
68  {
69  edm::LogInfo("RPCRecHitValid") << "Cannot find simHit collection\n";
70  return;
71  }
72 
73  // Retrieve RecHits from the event
75  if ( !event.getByLabel(recHitLabel_, recHitHandle) )
76  {
77  edm::LogInfo("RPCRecHitValid") << "Cannot find recHit collection\n";
78  return;
79  }
80 
81  typedef edm::PSimHitContainer::const_iterator SimHitIter;
82  typedef RPCRecHitCollection::const_iterator RecHitIter;
83 
84  // Loop over simHits, fill histograms which does not need associations
85  for ( SimHitIter simHitIter = simHitHandle->begin();
86  simHitIter != simHitHandle->end(); ++simHitIter )
87  {
88  if ( abs(simHitIter->particleType()) != 13 ) continue;
89 
90  const RPCDetId detId = static_cast<const RPCDetId>(simHitIter->detUnitId());
91  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
92  if ( !roll ) continue;
93 
94  const int region = roll->id().region();
95  const int ring = roll->id().ring();
96  //const int sector = roll->id().sector();
97  const int station = roll->id().station();
98  //const int layer = roll->id().layer();
99  //const int subSector = roll->id().subsector();
100 
101  if ( region == 0 )
102  {
103  h_.nRefHit_W->Fill(ring);
104  h_.nRefHit_WvsR->Fill(ring, station);
105  }
106  else
107  {
108  h_.nRefHit_D->Fill(region*station);
109  h_.nRefHit_DvsR->Fill(region*station, ring);
110  }
111 
112  const GlobalPoint pos = roll->toGlobal(simHitIter->localPosition());
113  }
114 
115  // Loop over recHits, fill histograms which does not need associations
116  for ( RecHitIter recHitIter = recHitHandle->begin();
117  recHitIter != recHitHandle->end(); ++recHitIter )
118  {
119  const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
120  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId()));
121  if ( !roll ) continue;
122 
123  const int region = roll->id().region();
124  const int ring = roll->id().ring();
125  //const int sector = roll->id().sector();
126  const int station = roll->id().station();
127  //const int layer = roll->id().layer();
128  //const int subSector = roll->id().subsector();
129 
130  h_.clusterSize->Fill(recHitIter->clusterSize());
131 
132  if ( region == 0 )
133  {
134  h_.nRecHit_W->Fill(ring);
135  h_.nRecHit_WvsR->Fill(ring, station);
136  }
137  else
138  {
139  h_.nRecHit_D->Fill(region*station);
140  h_.nRecHit_DvsR->Fill(ring, region*station);
141  }
142 
143  const GlobalPoint pos = roll->toGlobal(recHitIter->localPosition());
144  }
145 
146  // Start matching SimHits to RecHits
147  typedef std::map<SimHitIter, RecHitIter> SimToRecHitMap;
148  SimToRecHitMap simToRecHitMap;
149 
150  for ( SimHitIter simHitIter = simHitHandle->begin();
151  simHitIter != simHitHandle->end(); ++simHitIter )
152  {
153  if ( abs(simHitIter->particleType()) != 13 ) continue;
154 
155  const RPCDetId simDetId = static_cast<const RPCDetId>(simHitIter->detUnitId());
156  const RPCRoll* simRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(simDetId));
157  if ( !simRoll ) continue;
158 
159  const double simX = simHitIter->localPosition().x();
160 
161  for ( RecHitIter recHitIter = recHitHandle->begin();
162  recHitIter != recHitHandle->end(); ++recHitIter )
163  {
164  const RPCDetId recDetId = static_cast<const RPCDetId>(recHitIter->rpcId());
165  const RPCRoll* recRoll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(recDetId));
166  if ( !recRoll ) continue;
167 
168  if ( simDetId != recDetId ) continue;
169 
170  const double recX = recHitIter->localPosition().x();
171  const double newDx = fabs(recX - simX);
172 
173  // Associate SimHit to RecHit
174  SimToRecHitMap::const_iterator prevSimToReco = simToRecHitMap.find(simHitIter);
175  if ( prevSimToReco == simToRecHitMap.end() )
176  {
177  simToRecHitMap.insert(std::make_pair(simHitIter, recHitIter));
178  }
179  else
180  {
181  const double oldDx = fabs(prevSimToReco->second->localPosition().x() - simX);
182 
183  if ( newDx < oldDx )
184  {
185  simToRecHitMap[simHitIter] = recHitIter;
186  }
187  }
188  }
189  }
190 
191  // Now we have simHit-recHit mapping
192  // So we can fill up relavant histograms
193  for ( SimToRecHitMap::const_iterator match = simToRecHitMap.begin();
194  match != simToRecHitMap.end(); ++match )
195  {
196  SimHitIter simHitIter = match->first;
197  RecHitIter recHitIter = match->second;
198 
199  const RPCDetId detId = static_cast<const RPCDetId>(simHitIter->detUnitId());
200  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
201 
202  const int region = roll->id().region();
203  const int ring = roll->id().ring();
204  //const int sector = roll->id().sector();
205  const int station = roll->id().station();
206  //const int layer = roll->id().layer();
207  //const int subsector = roll->id().subsector();
208 
209  const double simX = simHitIter->localPosition().x();
210  const double recX = recHitIter->localPosition().x();
211  const double errX = recHitIter->localPositionError().xx();
212  const double dX = recX - simX;
213  const double pull = errX == 0 ? -999 : dX/errX;
214 
215  const GlobalPoint simPos = roll->toGlobal(simHitIter->localPosition());
216  const GlobalPoint recPos = roll->toGlobal(recHitIter->localPosition());
217 
218  if ( region == 0 )
219  {
220  h_.res_W->Fill(dX);
221  h_.pull_W->Fill(pull);
222  h_.nMatchedRefHit_W->Fill(ring);
223  h_.nMatchedRefHit_WvsR->Fill(ring, station);
224 
225  h_.res2_W->Fill(ring, dX);
226  h_.pull2_W->Fill(ring, pull);
227  }
228  else
229  {
230  h_.res_D->Fill(dX);
231  h_.pull_D->Fill(pull);
232  h_.nMatchedRefHit_D->Fill(region*station);
233  h_.nMatchedRefHit_DvsR->Fill(region*station, ring);
234 
235  h_.res2_D->Fill(region*station, dX);
236  h_.pull2_D->Fill(region*station, pull);
237  }
238 
239  }
240 
241  // Find Lost hits
242  for ( SimHitIter simHitIter = simHitHandle->begin();
243  simHitIter != simHitHandle->end(); ++simHitIter )
244  {
245  const RPCDetId detId = static_cast<const RPCDetId>(simHitIter->detUnitId());
246  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
247 
248  const int region = roll->id().region();
249  const int ring = roll->id().ring();
250  //const int sector = roll->id().sector();
251  const int station = roll->id().station();
252  //const int layer = roll->id().layer();
253  //const int subsector = roll->id().subsector();
254 
255  bool matched = false;
256  for ( SimToRecHitMap::const_iterator match = simToRecHitMap.begin();
257  match != simToRecHitMap.end(); ++match )
258  {
259  if ( simHitIter == match->first )
260  {
261  matched = true;
262  break;
263  }
264  }
265 
266  if ( !matched )
267  {
268  if ( region == 0 )
269  {
270  h_.nUnMatchedRefHit_W->Fill(ring);
271  h_.nUnMatchedRefHit_WvsR->Fill(ring, station);
272  }
273  else
274  {
275  h_.nUnMatchedRefHit_D->Fill(region*station);
276  h_.nUnMatchedRefHit_DvsR->Fill(region*station, ring);
277  }
278  }
279  }
280 
281  // Find Noisy hits
282  for ( RecHitIter recHitIter = recHitHandle->begin();
283  recHitIter != recHitHandle->end(); ++recHitIter )
284  {
285  const RPCDetId detId = static_cast<const RPCDetId>(recHitIter->rpcId());
286  const RPCRoll* roll = dynamic_cast<const RPCRoll*>(rpcGeom->roll(detId));
287 
288  const int region = roll->id().region();
289  const int ring = roll->id().ring();
290  //const int sector = roll->id().sector();
291  const int station = roll->id().station();
292  //const int layer = roll->id().layer();
293  //const int subsector = roll->id().subsector();
294 
295  bool matched = false;
296  for ( SimToRecHitMap::const_iterator match = simToRecHitMap.begin();
297  match != simToRecHitMap.end(); ++match )
298  {
299  if ( recHitIter == match->second )
300  {
301  matched = true;
302  break;
303  }
304  }
305 
306  if ( !matched )
307  {
308  if ( region == 0 )
309  {
310  h_.nUnMatchedRecHit_W->Fill(ring);
311  h_.nUnMatchedRecHit_WvsR->Fill(ring, station);
312  }
313  else
314  {
315  h_.nUnMatchedRecHit_D->Fill(region*station);
316  h_.nUnMatchedRecHit_DvsR->Fill(region*station, ring);
317  }
318 
319  const GlobalPoint pos = roll->toGlobal(recHitIter->localPosition());
320 // h_[HName::NoisyHitEta]->Fill(pos.eta());
321  }
322  }
323 }
DQMStore * dbe_
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
#define abs(x)
Definition: mlp_lapack.h:159
void Fill(long long x)
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:46
RPCDetId id() const
Definition: RPCRoll.cc:24
edm::InputTag recHitLabel_
int ring() const
Definition: RPCDetId.h:74
RPCValidHistograms h_
edm::InputTag simHitLabel_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
const T & get() const
Definition: EventSetup.h:55
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:65
int station() const
Definition: RPCDetId.h:98
void RPCRecHitValid::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 45 of file RPCRecHitValid.cc.

46 {
47 }
void RPCRecHitValid::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 49 of file RPCRecHitValid.cc.

50 {
51 }

Member Data Documentation

DQMStore* RPCRecHitValid::dbe_
private

Definition at line 33 of file RPCRecHitValid.h.

RPCValidHistograms RPCRecHitValid::h_
private

Definition at line 37 of file RPCRecHitValid.h.

bool RPCRecHitValid::isStandAloneMode_
private

Definition at line 35 of file RPCRecHitValid.h.

edm::InputTag RPCRecHitValid::recHitLabel_
private

Definition at line 31 of file RPCRecHitValid.h.

std::string RPCRecHitValid::rootFileName_
private

Definition at line 34 of file RPCRecHitValid.h.

edm::InputTag RPCRecHitValid::simHitLabel_
private

Definition at line 31 of file RPCRecHitValid.h.