CMS 3D CMS Logo

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