CMS 3D CMS Logo

HGCalHitValidation.cc
Go to the documentation of this file.
1 //
3 // Package: HGCalHitValidation
4 // Class: HGCalHitValidation
5 //
20 
26 
38 
43 
44 #include <cmath>
45 #include <memory>
46 #include <iostream>
47 #include <string>
48 #include <vector>
49 
50 //#define EDM_ML_DEBUG
51 
53 public:
54  explicit HGCalHitValidation(const edm::ParameterSet&);
55  ~HGCalHitValidation() override;
56  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
57 
58 protected:
59  typedef std::tuple<float, float, float, float> HGCHitTuple;
60 
61  void dqmBeginRun(edm::Run const&, edm::EventSetup const&) override;
62  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
63  void analyze(const edm::Event&, const edm::EventSetup&) override;
64  void analyzeHGCalSimHit(edm::Handle<std::vector<PCaloHit>> const& simHits,
65  int idet,
67  std::map<unsigned int, HGCHitTuple>&);
68 
69 private:
70  //HGC Geometry
71  std::vector<const HGCalDDDConstants*> hgcCons_;
72  std::vector<const HGCalGeometry*> hgcGeometry_;
73  std::vector<std::string> geometrySource_;
74  std::vector<int> ietaExcludeBH_;
75 
83  std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord>> tok_ddd_;
84  std::vector<edm::ESGetToken<HGCalGeometry, IdealGeometryRecord>> tok_geom_;
85 
86  //histogram related stuff
90 
95 
99 };
100 
102  geometrySource_ = cfg.getParameter<std::vector<std::string>>("geometrySource");
103  eeSimHitToken_ = consumes<std::vector<PCaloHit>>(cfg.getParameter<edm::InputTag>("eeSimHitSource"));
104  fhSimHitToken_ = consumes<std::vector<PCaloHit>>(cfg.getParameter<edm::InputTag>("fhSimHitSource"));
105  bhSimHitToken_ = consumes<std::vector<PCaloHit>>(cfg.getParameter<edm::InputTag>("bhSimHitSource"));
106  eeRecHitToken_ = consumes<HGCeeRecHitCollection>(cfg.getParameter<edm::InputTag>("eeRecHitSource"));
107  fhRecHitToken_ = consumes<HGChefRecHitCollection>(cfg.getParameter<edm::InputTag>("fhRecHitSource"));
108  bhRecHitToken_ = consumes<HGChebRecHitCollection>(cfg.getParameter<edm::InputTag>("bhRecHitSource"));
109  ietaExcludeBH_ = cfg.getParameter<std::vector<int>>("ietaExcludeBH");
110 
111  for (const auto& name : geometrySource_) {
112  tok_ddd_.emplace_back(
113  esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag{"", name}));
114  tok_geom_.emplace_back(
115  esConsumes<HGCalGeometry, IdealGeometryRecord, edm::Transition::BeginRun>(edm::ESInputTag{"", name}));
116  }
117 
118 #ifdef EDM_ML_DEBUG
119  edm::LogInfo("HGCalValid") << "Exclude the following " << ietaExcludeBH_.size() << " ieta values from BH plots";
120  for (unsigned int k = 0; k < ietaExcludeBH_.size(); ++k)
121  edm::LogInfo("HGCalValid") << " [" << k << "] " << ietaExcludeBH_[k];
122 #endif
123 }
124 
126 
129  std::vector<std::string> source = {"HGCalEESensitive", "HGCalHESiliconSensitive", "HGCalHEScintillatorSensitive"};
130  desc.add<std::vector<std::string>>("geometrySource", source);
131  desc.add<edm::InputTag>("eeSimHitSource", edm::InputTag("g4SimHits", "HGCHitsEE"));
132  desc.add<edm::InputTag>("fhSimHitSource", edm::InputTag("g4SimHits", "HGCHitsHEfront"));
133  desc.add<edm::InputTag>("bhSimHitSource", edm::InputTag("g4SimHits", "HGCHitsHEback"));
134  desc.add<edm::InputTag>("eeRecHitSource", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
135  desc.add<edm::InputTag>("fhRecHitSource", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
136  desc.add<edm::InputTag>("bhRecHitSource", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
137  std::vector<int> dummy;
138  desc.add<std::vector<int>>("ietaExcludeBH", dummy);
139  descriptions.add("hgcalHitValidation", desc);
140 }
141 
143  iB.setCurrentFolder("HGCAL/HGCalSimHitsV/HitValidation");
144 
145  //initiating histograms
146  heedzVsZ = iB.book2D("heedzVsZ", "", 720, -360, 360, 100, -0.1, 0.1);
147  heedyVsY = iB.book2D("heedyVsY", "", 400, -200, 200, 100, -0.02, 0.02);
148  heedxVsX = iB.book2D("heedxVsX", "", 400, -200, 200, 100, -0.02, 0.02);
149  heeRecVsSimZ = iB.book2D("heeRecVsSimZ", "", 720, -360, 360, 720, -360, 360);
150  heeRecVsSimY = iB.book2D("heeRecVsSimY", "", 400, -200, 200, 400, -200, 200);
151  heeRecVsSimX = iB.book2D("heeRecVsSimX", "", 400, -200, 200, 400, -200, 200);
152 
153  hefdzVsZ = iB.book2D("hefdzVsZ", "", 820, -410, 410, 100, -0.1, 0.1);
154  hefdyVsY = iB.book2D("hefdyVsY", "", 400, -200, 200, 100, -0.02, 0.02);
155  hefdxVsX = iB.book2D("hefdxVsX", "", 400, -200, 200, 100, -0.02, 0.02);
156  hefRecVsSimZ = iB.book2D("hefRecVsSimZ", "", 820, -410, 410, 820, -410, 410);
157  hefRecVsSimY = iB.book2D("hefRecVsSimY", "", 400, -200, 200, 400, -200, 200);
158  hefRecVsSimX = iB.book2D("hefRecVsSimX", "", 400, -200, 200, 400, -200, 200);
159 
160  hebdzVsZ = iB.book2D("hebdzVsZ", "", 1080, -540, 540, 100, -1.0, 1.0);
161  hebdPhiVsPhi = iB.book2D("hebdPhiVsPhi", "", M_PI * 100, -0.5, M_PI + 0.5, 200, -0.2, 0.2);
162  hebdEtaVsEta = iB.book2D("hebdEtaVsEta", "", 1000, -5, 5, 200, -0.1, 0.1);
163  hebRecVsSimZ = iB.book2D("hebRecVsSimZ", "", 1080, -540, 540, 1080, -540, 540);
164  hebRecVsSimY = iB.book2D("hebRecVsSimY", "", 400, -200, 200, 400, -200, 200);
165  hebRecVsSimX = iB.book2D("hebRecVsSimX", "", 400, -200, 200, 400, -200, 200);
166 
167  heeEnRec = iB.book1D("heeEnRec", "", 1000, 0, 10);
168  heeEnSim = iB.book1D("heeEnSim", "", 1000, 0, 0.01);
169  heeEnSimRec = iB.book2D("heeEnSimRec", "", 200, 0, 0.002, 200, 0, 0.2);
170 
171  hefEnRec = iB.book1D("hefEnRec", "", 1000, 0, 10);
172  hefEnSim = iB.book1D("hefEnSim", "", 1000, 0, 0.01);
173  hefEnSimRec = iB.book2D("hefEnSimRec", "", 200, 0, 0.001, 200, 0, 0.5);
174 
175  hebEnRec = iB.book1D("hebEnRec", "", 1000, 0, 15);
176  hebEnSim = iB.book1D("hebEnSim", "", 1000, 0, 0.01);
177  hebEnSimRec = iB.book2D("hebEnSimRec", "", 200, 0, 0.02, 200, 0, 4);
178 }
179 
180 void HGCalHitValidation::dqmBeginRun(edm::Run const& iRun, edm::EventSetup const& iSetup) {
181  //initiating hgc Geometry
182  for (size_t i = 0; i < geometrySource_.size(); i++) {
184  if (hgcCons.isValid()) {
185  hgcCons_.push_back(hgcCons.product());
186  } else {
187  edm::LogWarning("HGCalValid") << "Cannot initiate HGCalDDDConstants for " << geometrySource_[i] << std::endl;
188  }
190  if (hgcGeom.isValid()) {
191  hgcGeometry_.push_back(hgcGeom.product());
192  } else {
193  edm::LogWarning("HGCalValid") << "Cannot initiate HGCalGeometry for " << geometrySource_[i] << std::endl;
194  }
195  }
196 }
197 
199  std::map<unsigned int, HGCHitTuple> eeHitRefs, fhHitRefs, bhHitRefs;
200 
201  //Accesing ee simhits
203  iEvent.getByToken(eeSimHitToken_, eeSimHits);
204 
205  if (eeSimHits.isValid()) {
206  analyzeHGCalSimHit(eeSimHits, 0, heeEnSim, eeHitRefs);
207 #ifdef EDM_ML_DEBUG
208  for (std::map<unsigned int, HGCHitTuple>::iterator itr = eeHitRefs.begin(); itr != eeHitRefs.end(); ++itr) {
209  int idx = std::distance(eeHitRefs.begin(), itr);
210  edm::LogInfo("HGCalValid") << "EEHit[" << idx << "] " << std::hex << itr->first << std::dec << "; Energy "
211  << std::get<0>(itr->second) << "; Position (" << std::get<1>(itr->second) << ", "
212  << std::get<2>(itr->second) << ", " << std::get<3>(itr->second) << ")";
213  }
214 #endif
215  } else {
216  edm::LogVerbatim("HGCalValid") << "No EE SimHit Found ";
217  }
218 
219  //Accesing fh simhits
221  iEvent.getByToken(fhSimHitToken_, fhSimHits);
222  if (fhSimHits.isValid()) {
223  analyzeHGCalSimHit(fhSimHits, 1, hefEnSim, fhHitRefs);
224 #ifdef EDM_ML_DEBUG
225  for (std::map<unsigned int, HGCHitTuple>::iterator itr = fhHitRefs.begin(); itr != fhHitRefs.end(); ++itr) {
226  int idx = std::distance(fhHitRefs.begin(), itr);
227  edm::LogInfo("HGCalValid") << "FHHit[" << idx << "] " << std::hex << itr->first << std::dec << "; Energy "
228  << std::get<0>(itr->second) << "; Position (" << std::get<1>(itr->second) << ", "
229  << std::get<2>(itr->second) << ", " << std::get<3>(itr->second) << ")";
230  }
231 #endif
232  } else {
233  edm::LogVerbatim("HGCalValid") << "No FH SimHit Found ";
234  }
235 
236  //Accessing bh simhits
238  iEvent.getByToken(bhSimHitToken_, bhSimHits);
239  if (bhSimHits.isValid()) {
240  analyzeHGCalSimHit(bhSimHits, 2, hebEnSim, bhHitRefs);
241 #ifdef EDM_ML_DEBUG
242  for (std::map<unsigned int, HGCHitTuple>::iterator itr = bhHitRefs.begin(); itr != bhHitRefs.end(); ++itr) {
243  int idx = std::distance(bhHitRefs.begin(), itr);
244  edm::LogInfo("HGCalValid") << "BHHit[" << idx << "] " << std::hex << itr->first << std::dec << "; Energy "
245  << std::get<0>(itr->second) << "; Position (" << std::get<1>(itr->second) << ", "
246  << std::get<2>(itr->second) << ", " << std::get<3>(itr->second) << ")";
247  }
248 #endif
249  } else {
250  edm::LogVerbatim("HGCalValid") << "No BH SimHit Found ";
251  }
252 
253  //accessing EE Rechit information
255  iEvent.getByToken(eeRecHitToken_, eeRecHit);
256  if (eeRecHit.isValid()) {
257  const HGCeeRecHitCollection* theHits = (eeRecHit.product());
258  for (auto it = theHits->begin(); it != theHits->end(); ++it) {
259  double energy = it->energy();
260  heeEnRec->Fill(energy);
261  std::map<unsigned int, HGCHitTuple>::const_iterator itr = eeHitRefs.find(it->id().rawId());
262  if (itr != eeHitRefs.end()) {
263  GlobalPoint xyz = hgcGeometry_[0]->getPosition(it->id());
264  heeRecVsSimX->Fill(std::get<1>(itr->second), xyz.x());
265  heeRecVsSimY->Fill(std::get<2>(itr->second), xyz.y());
266  heeRecVsSimZ->Fill(std::get<3>(itr->second), xyz.z());
267  heedxVsX->Fill(std::get<1>(itr->second), (xyz.x() - std::get<1>(itr->second)));
268  heedyVsY->Fill(std::get<2>(itr->second), (xyz.y() - std::get<2>(itr->second)));
269  heedzVsZ->Fill(std::get<3>(itr->second), (xyz.z() - std::get<3>(itr->second)));
270  heeEnSimRec->Fill(std::get<0>(itr->second), energy);
271 #ifdef EDM_ML_DEBUG
272  edm::LogInfo("HGCalValid") << "EEHit: " << std::hex << it->id().rawId() << std::dec << " Sim ("
273  << std::get<0>(itr->second) << ", " << std::get<1>(itr->second) << ", "
274  << std::get<2>(itr->second) << ", " << std::get<3>(itr->second) << ") Rec ("
275  << energy << ", " << xyz.x() << ", " << xyz.y() << ", " << xyz.z() << ")";
276 #endif
277  }
278  }
279  } else {
280  edm::LogVerbatim("HGCalValid") << "No EE RecHit Found ";
281  }
282 
283  //accessing FH Rechit information
285  iEvent.getByToken(fhRecHitToken_, fhRecHit);
286  if (fhRecHit.isValid()) {
287  const HGChefRecHitCollection* theHits = (fhRecHit.product());
288  for (auto it = theHits->begin(); it != theHits->end(); ++it) {
289  double energy = it->energy();
290  hefEnRec->Fill(energy);
291  std::map<unsigned int, HGCHitTuple>::const_iterator itr = fhHitRefs.find(it->id().rawId());
292  if (itr != fhHitRefs.end()) {
293  GlobalPoint xyz = hgcGeometry_[1]->getPosition(it->id());
294 
295  hefRecVsSimX->Fill(std::get<1>(itr->second), xyz.x());
296  hefRecVsSimY->Fill(std::get<2>(itr->second), xyz.y());
297  hefRecVsSimZ->Fill(std::get<3>(itr->second), xyz.z());
298  hefdxVsX->Fill(std::get<1>(itr->second), (xyz.x() - std::get<1>(itr->second)));
299  hefdyVsY->Fill(std::get<2>(itr->second), (xyz.y() - std::get<2>(itr->second)));
300  hefdzVsZ->Fill(std::get<3>(itr->second), (xyz.z() - std::get<3>(itr->second)));
301  hefEnSimRec->Fill(std::get<0>(itr->second), energy);
302 #ifdef EDM_ML_DEBUG
303  edm::LogInfo("HGCalValid") << "FHHit: " << std::hex << it->id().rawId() << std::dec << " Sim ("
304  << std::get<0>(itr->second) << ", " << std::get<1>(itr->second) << ", "
305  << std::get<2>(itr->second) << ", " << std::get<3>(itr->second) << ") Rec ("
306  << energy << "," << xyz.x() << ", " << xyz.y() << ", " << xyz.z() << ")";
307 #endif
308  }
309  }
310  } else {
311  edm::LogVerbatim("HGCalValid") << "No FH RecHit Found ";
312  }
313 
314  //accessing BH Rechit information
316  iEvent.getByToken(bhRecHitToken_, bhRecHit);
317  if (bhRecHit.isValid()) {
318  const HGChebRecHitCollection* theHits = (bhRecHit.product());
319  for (auto it = theHits->begin(); it != theHits->end(); ++it) {
320  double energy = it->energy();
321  hebEnRec->Fill(energy);
322  std::map<unsigned int, HGCHitTuple>::const_iterator itr = bhHitRefs.find(it->id().rawId());
323  GlobalPoint xyz = hgcGeometry_[2]->getPosition(it->id());
324  if (itr != bhHitRefs.end()) {
325  float ang3 = xyz.phi().value(); // returns the phi in radians
326  double fac = sinh(std::get<1>(itr->second));
327  double pT = std::get<3>(itr->second) / fac;
328  double xp = pT * cos(std::get<2>(itr->second));
329  double yp = pT * sin(std::get<2>(itr->second));
330  hebRecVsSimX->Fill(xp, xyz.x());
331  hebRecVsSimY->Fill(yp, xyz.y());
332  hebRecVsSimZ->Fill(std::get<3>(itr->second), xyz.z());
333  hebdEtaVsEta->Fill(std::get<1>(itr->second), (xyz.eta() - std::get<1>(itr->second)));
334  hebdPhiVsPhi->Fill(std::get<2>(itr->second), (ang3 - std::get<2>(itr->second)));
335  hebdzVsZ->Fill(std::get<3>(itr->second), (xyz.z() - std::get<3>(itr->second)));
336  hebEnSimRec->Fill(std::get<0>(itr->second), energy);
337 
338 #ifdef EDM_ML_DEBUG
339  edm::LogInfo("HGCalValid") << "BHHit: " << std::hex << it->id().rawId() << std::dec << " Sim ("
340  << std::get<0>(itr->second) << ", " << std::get<1>(itr->second) << ", "
341  << std::get<2>(itr->second) << ", " << std::get<3>(itr->second) << ") Rec ("
342  << energy << ", " << xyz.x() << ", " << xyz.y() << ", " << xyz.z() << ")\n";
343 #endif
344  }
345  }
346  } else {
347  edm::LogVerbatim("HGCalValid") << "No BH RecHit Found ";
348  }
349 }
350 
351 void HGCalHitValidation::analyzeHGCalSimHit(edm::Handle<std::vector<PCaloHit>> const& simHits,
352  int idet,
354  std::map<unsigned int, HGCHitTuple>& hitRefs) {
355  const HGCalTopology& hTopo = hgcGeometry_[idet]->topology();
356  for (std::vector<PCaloHit>::const_iterator simHit = simHits->begin(); simHit != simHits->end(); ++simHit) {
357  int subdet, zside, layer, wafer, celltype, cell;
358  HGCalTestNumbering::unpackHexagonIndex(simHit->id(), subdet, zside, layer, wafer, celltype, cell);
359  std::pair<float, float> xy = hgcCons_[idet]->locateCell(cell, layer, wafer, false);
360  float zp = hgcCons_[idet]->waferZ(layer, false);
361  if (zside < 0)
362  zp = -zp;
363  float xp = (zp < 0) ? -xy.first / 10 : xy.first / 10;
364  float yp = xy.second / 10.0;
365 
366  //skip this hit if after ganging it is not valid
367  std::pair<int, int> recoLayerCell = hgcCons_[idet]->simToReco(cell, layer, wafer, hTopo.detectorType());
368  cell = recoLayerCell.first;
369  layer = recoLayerCell.second;
370 
371  //skip this hit if after ganging it is not valid
372  if (layer < 0 || cell < 0) {
373  } else {
374  //assign the RECO DetId
375  HGCalDetId id = HGCalDetId((ForwardSubdetector)(subdet), zside, layer, celltype, wafer, cell);
376  float energy = simHit->energy();
377 
378  float energySum(energy);
379  if (hitRefs.count(id.rawId()) != 0)
380  energySum += std::get<0>(hitRefs[id.rawId()]);
381  hitRefs[id.rawId()] = std::make_tuple(energySum, xp, yp, zp);
382  hist->Fill(energy);
383  }
384  }
385 }
386 
387 //define this as a plug-in
ConfigurationDescriptions.h
edm::ESHandle::product
T const * product() const
Definition: ESHandle.h:86
dqm::impl::MonitorElement
Definition: MonitorElement.h:98
Handle.h
Geom::Phi::value
T1 value() const
Explicit access to value in case implicit conversion not OK.
Definition: Phi.h:75
mps_fire.i
i
Definition: mps_fire.py:428
edm::ESInputTag
Definition: ESInputTag.h:87
ESTransientHandle.h
edm::Handle::product
T const * product() const
Definition: Handle.h:70
ESHandle.h
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
HGCNumberingScheme.h
HGCalHitValidation::hebEnRec
MonitorElement * hebEnRec
Definition: HGCalHitValidation.cc:96
edm::Run
Definition: Run.h:45
ecaldqm::zside
int zside(DetId const &)
Definition: EcalDQMCommonUtils.cc:189
ForwardSubdetector
ForwardSubdetector
Definition: ForwardSubdetector.h:4
edm::EDGetTokenT
Definition: EDGetToken.h:33
HGCalGeometryMode.h
HGCalHitValidation::tok_ddd_
std::vector< edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > > tok_ddd_
Definition: HGCalHitValidation.cc:83
HGCalTestNumbering::unpackHexagonIndex
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
Definition: HGCalTestNumbering.cc:46
HGCalHitValidation::bhSimHitSource
edm::InputTag bhSimHitSource
Definition: HGCalHitValidation.cc:76
HGCalHitValidation::heedxVsX
MonitorElement * heedxVsX
Definition: HGCalHitValidation.cc:87
HGCalHitValidation::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: HGCalHitValidation.cc:127
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
HGCalHitValidation::hebRecVsSimX
MonitorElement * hebRecVsSimX
Definition: HGCalHitValidation.cc:93
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
dqm::implementation::NavigatorBase::setCurrentFolder
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
HGCalHitValidation::bookHistograms
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: HGCalHitValidation.cc:142
FastTrackerRecHitCombiner_cfi.simHits
simHits
Definition: FastTrackerRecHitCombiner_cfi.py:5
edm::SortedCollection
Definition: SortedCollection.h:49
HGCalHitValidation::heedyVsY
MonitorElement * heedyVsY
Definition: HGCalHitValidation.cc:87
DQMStore.h
HGCalHitValidation::eeRecHitToken_
edm::EDGetTokenT< HGCeeRecHitCollection > eeRecHitToken_
Definition: HGCalHitValidation.cc:80
HGCalHitValidation::heeEnSim
MonitorElement * heeEnSim
Definition: HGCalHitValidation.cc:98
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
HGCalHitValidation::HGCalHitValidation
HGCalHitValidation(const edm::ParameterSet &)
Definition: HGCalHitValidation.cc:101
HGCalHitValidation::hefEnSim
MonitorElement * hefEnSim
Definition: HGCalHitValidation.cc:97
HGCalHitValidation::heedzVsZ
MonitorElement * heedzVsZ
Definition: HGCalHitValidation.cc:87
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
ForwardSubdetector.h
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
funct::sin
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
CastorDataFrameFilter_impl::energySum
double energySum(const DataFrame &df, int fs, int ls)
Definition: CastorDataFrameFilter.cc:22
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
HGCalHitValidation::hefEnRec
MonitorElement * hefEnRec
Definition: HGCalHitValidation.cc:97
MakerMacros.h
HGCalHitValidation::fhRecHitToken_
edm::EDGetTokenT< HGChefRecHitCollection > fhRecHitToken_
Definition: HGCalHitValidation.cc:81
funct::cos
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
HGCalHitValidation::bhSimHitToken_
edm::EDGetTokenT< std::vector< PCaloHit > > bhSimHitToken_
Definition: HGCalHitValidation.cc:79
HGCRecHit.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
Service.h
rpcPointValidation_cfi.simHit
simHit
Definition: rpcPointValidation_cfi.py:24
PVValHelper::pT
Definition: PVValidationHelpers.h:71
HGCalTopology::detectorType
bool detectorType() const
Definition: HGCalTopology.h:119
dqm::impl::MonitorElement::Fill
void Fill(long long x)
Definition: MonitorElement.h:290
HGCalHitValidation::hefRecVsSimY
MonitorElement * hefRecVsSimY
Definition: HGCalHitValidation.cc:92
source
static const std::string source
Definition: EdmProvDump.cc:47
edm::ESHandle< HGCalDDDConstants >
HGCalHitValidation::~HGCalHitValidation
~HGCalHitValidation() override
Definition: HGCalHitValidation.cc:125
HGCalHitValidation::heeRecVsSimX
MonitorElement * heeRecVsSimX
Definition: HGCalHitValidation.cc:91
HCALHighEnergyHPDFilter_cfi.energy
energy
Definition: HCALHighEnergyHPDFilter_cfi.py:5
HGCalHitValidation::hebRecVsSimZ
MonitorElement * hebRecVsSimZ
Definition: HGCalHitValidation.cc:93
HGCalHitValidation::dqmBeginRun
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
Definition: HGCalHitValidation.cc:180
dqmdumpme.k
k
Definition: dqmdumpme.py:60
HGCalHitValidation::hebEnSimRec
MonitorElement * hebEnSimRec
Definition: HGCalHitValidation.cc:94
HGCalHitValidation::heeRecVsSimY
MonitorElement * heeRecVsSimY
Definition: HGCalHitValidation.cc:91
Point3DBase< float, GlobalTag >
EDGetToken.h
HGCalHitValidation::hefRecVsSimZ
MonitorElement * hefRecVsSimZ
Definition: HGCalHitValidation.cc:92
geometryCSVtoXML.xy
xy
Definition: geometryCSVtoXML.py:19
HGCalHitValidation::hebRecVsSimY
MonitorElement * hebRecVsSimY
Definition: HGCalHitValidation.cc:93
DQMEDAnalyzer.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
HGCalHitValidation::eeSimHitSource
edm::InputTag eeSimHitSource
Definition: HGCalHitValidation.cc:76
DQMEDAnalyzer
Definition: DQMEDAnalyzer.py:1
HGCalHitValidation::hebdEtaVsEta
MonitorElement * hebdEtaVsEta
Definition: HGCalHitValidation.cc:89
HGCalHitValidation::hefRecVsSimX
MonitorElement * hefRecVsSimX
Definition: HGCalHitValidation.cc:92
CaloSubdetectorGeometry.h
HGCRecHitCollections.h
HGCalGeometry.h
edm::ParameterSet
Definition: ParameterSet.h:47
HGCalHitValidation::fhSimHitSource
edm::InputTag fhSimHitSource
Definition: HGCalHitValidation.cc:76
Event.h
HGCalHitValidation::fhSimHitToken_
edm::EDGetTokenT< std::vector< PCaloHit > > fhSimHitToken_
Definition: HGCalHitValidation.cc:78
gpuVertexFinder::hist
__shared__ Hist hist
Definition: gpuClusterTracksDBSCAN.h:48
HGCalHitValidation::bhRecHitToken_
edm::EDGetTokenT< HGChebRecHitCollection > bhRecHitToken_
Definition: HGCalHitValidation.cc:82
HGCalHitValidation::hgcGeometry_
std::vector< const HGCalGeometry * > hgcGeometry_
Definition: HGCalHitValidation.cc:72
HGCalHitValidation::hefdxVsX
MonitorElement * hefdxVsX
Definition: HGCalHitValidation.cc:88
PCaloHit.h
PV3DBase::eta
T eta() const
Definition: PV3DBase.h:73
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
iEvent
int iEvent
Definition: GenABIO.cc:224
HGCalHitValidation::heeRecVsSimZ
MonitorElement * heeRecVsSimZ
Definition: HGCalHitValidation.cc:91
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
edm::EventSetup::getHandle
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:155
hgcalTriggerNtuples_cfi.fhSimHits
fhSimHits
Definition: hgcalTriggerNtuples_cfi.py:51
IdealGeometryRecord.h
edm::EventSetup
Definition: EventSetup.h:58
HGCalHitValidation::hefdyVsY
MonitorElement * hefdyVsY
Definition: HGCalHitValidation.cc:88
hgcalTriggerNtuples_cfi.eeSimHits
eeSimHits
Definition: hgcalTriggerNtuples_cfi.py:50
HGCalHitValidation::hgcCons_
std::vector< const HGCalDDDConstants * > hgcCons_
Definition: HGCalHitValidation.cc:71
edm::ESHandleBase::isValid
bool isValid() const
Definition: ESHandle.h:44
HGCalTopology
Definition: HGCalTopology.h:12
HGCalHitValidation
Definition: HGCalHitValidation.cc:52
looper.cfg
cfg
Definition: looper.py:296
HGCalHitValidation::heeEnRec
MonitorElement * heeEnRec
Definition: HGCalHitValidation.cc:98
HGCalHitValidation::tok_geom_
std::vector< edm::ESGetToken< HGCalGeometry, IdealGeometryRecord > > tok_geom_
Definition: HGCalHitValidation.cc:84
HGCalDetId
Definition: HGCalDetId.h:8
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
HGCalDetId.h
HGCalHitValidation::hebdzVsZ
MonitorElement * hebdzVsZ
Definition: HGCalHitValidation.cc:89
HGCalHitValidation::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: HGCalHitValidation.cc:198
bscTrigger_cfi.theHits
theHits
Definition: bscTrigger_cfi.py:18
HGCalHitValidation::heeEnSimRec
MonitorElement * heeEnSimRec
Definition: HGCalHitValidation.cc:94
Frameworkfwd.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
HGCalHitValidation::hefdzVsZ
MonitorElement * hefdzVsZ
Definition: HGCalHitValidation.cc:88
dqm::implementation::IBooker::book2D
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:177
HGCalHitValidation::hefEnSimRec
MonitorElement * hefEnSimRec
Definition: HGCalHitValidation.cc:94
CaloGeometry.h
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
HGCalHitValidation::geometrySource_
std::vector< std::string > geometrySource_
Definition: HGCalHitValidation.cc:73
Exception.h
PCaloHitContainer.h
dqm::implementation::IBooker
Definition: DQMStore.h:43
HGCalHitValidation::eeSimHitToken_
edm::EDGetTokenT< std::vector< PCaloHit > > eeSimHitToken_
Definition: HGCalHitValidation.cc:77
HGCalHitValidation::analyzeHGCalSimHit
void analyzeHGCalSimHit(edm::Handle< std::vector< PCaloHit >> const &simHits, int idet, MonitorElement *hist, std::map< unsigned int, HGCHitTuple > &)
Definition: HGCalHitValidation.cc:351
HGCalDDDConstants.h
ParameterSet.h
dummy
Definition: DummySelector.h:38
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
hgcalTriggerNtuples_cfi.bhSimHits
bhSimHits
Definition: hgcalTriggerNtuples_cfi.py:52
edm::Event
Definition: Event.h:73
HLT_FULL_cff.distance
distance
Definition: HLT_FULL_cff.py:7746
HGCalHitValidation::hebEnSim
MonitorElement * hebEnSim
Definition: HGCalHitValidation.cc:96
edm::Log
Definition: MessageLogger.h:70
TauDecayModes.dec
dec
Definition: TauDecayModes.py:142
edm::InputTag
Definition: InputTag.h:15
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
HGCalHitValidation::hebdPhiVsPhi
MonitorElement * hebdPhiVsPhi
Definition: HGCalHitValidation.cc:89
HGCalHitValidation::ietaExcludeBH_
std::vector< int > ietaExcludeBH_
Definition: HGCalHitValidation.cc:74
dqm::implementation::IBooker::book1D
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
HGCalTestNumbering.h
HGCalHitValidation::HGCHitTuple
std::tuple< float, float, float, float > HGCHitTuple
Definition: HGCalHitValidation.cc:59