CMS 3D CMS Logo

HGCalHitValidation.cc
Go to the documentation of this file.
1 //
3 // Package: HGCalHitValidation
4 // Class: HGCalHitValidation
5 //
26 
34 
47 
53 
54 #include <cmath>
55 #include <memory>
56 #include <iostream>
57 #include <string>
58 #include <vector>
59 
60 //#define EDM_ML_DEBUG
61 
63 
64 public:
65 
66  explicit HGCalHitValidation( const edm::ParameterSet& );
67  ~HGCalHitValidation() override;
68  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
69 
70 protected:
71  typedef std::tuple<float,float,float,float> HGCHitTuple;
72 
73  void dqmBeginRun(edm::Run const&, edm::EventSetup const&) override;
74  void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
75  void analyze(const edm::Event&, const edm::EventSetup&) override;
76  void analyzeHGCalSimHit(edm::Handle<std::vector<PCaloHit>> const& simHits,
77  int idet, MonitorElement *hist,
78  std::map<unsigned int, HGCHitTuple>&);
79  template<class T1>
80  void analyzeHGCalRecHit(T1 const & theHits,
81  std::map<unsigned int, HGCHitTuple> const& hitRefs);
82 
83 private:
84  //HGC Geometry
85  std::vector<const HGCalDDDConstants*> hgcCons_;
86  std::vector<const HGCalGeometry*> hgcGeometry_;
90  std::vector<std::string> geometrySource_;
91  std::vector<int> ietaExcludeBH_;
92  bool ifHCAL_;
93 
102 
103  //histogram related stuff
107 
112 
116 
117 };
118 
119 
121 
122  geometrySource_ = cfg.getUntrackedParameter< std::vector<std::string> >("geometrySource");
123  eeSimHitToken_ = consumes<std::vector<PCaloHit>>(cfg.getParameter<edm::InputTag>("eeSimHitSource"));
124  fhSimHitToken_ = consumes<std::vector<PCaloHit>>(cfg.getParameter<edm::InputTag>("fhSimHitSource"));
125  bhSimHitToken_ = consumes<std::vector<PCaloHit>>(cfg.getParameter<edm::InputTag>("bhSimHitSource"));
126  eeRecHitToken_ = consumes<HGCeeRecHitCollection>(cfg.getParameter<edm::InputTag>("eeRecHitSource"));
127  fhRecHitToken_ = consumes<HGChefRecHitCollection>(cfg.getParameter<edm::InputTag>("fhRecHitSource"));
128  ietaExcludeBH_ = cfg.getParameter<std::vector<int> >("ietaExcludeBH");
129  ifHCAL_ = cfg.getParameter<bool>("ifHCAL");
130  if (ifHCAL_)
131  bhRecHitTokenh_ = consumes<HBHERecHitCollection>(cfg.getParameter<edm::InputTag>("bhRecHitSource"));
132  else
133  bhRecHitTokeng_ = consumes<HGChebRecHitCollection>(cfg.getParameter<edm::InputTag>("bhRecHitSource"));
134 
135 #ifdef EDM_ML_DEBUG
136  edm::LogInfo("HGCalValid") << "Exclude the following "
137  << ietaExcludeBH_.size()
138  << " ieta values from BH plots (BH "
139  << ifHCAL_ << ") ";
140  for (unsigned int k=0; k<ietaExcludeBH_.size(); ++k)
141  edm::LogInfo("HGCalValid") << " [" << k << "] " << ietaExcludeBH_[k];
142 #endif
143 }
144 
146 
148  //The following says we do not know what parameters are allowed so do no validation
149  // Please change this to state exactly what you do use, even if it is no parameters
151  desc.setUnknown();
152  descriptions.addDefault(desc);
153 }
154 
156  edm::Run const&,
157  edm::EventSetup const&) {
158 
159  iB.setCurrentFolder("HGCalSimHitsV/HitValidation");
160 
161  //initiating histograms
162  heedzVsZ = iB.book2D("heedzVsZ","",7200,-360,360,100,-0.1,0.1);
163  heedyVsY = iB.book2D("heedyVsY","",400,-200,200,100,-0.02,0.02);
164  heedxVsX = iB.book2D("heedxVsX","",400,-200,200,100,-0.02,0.02);
165  heeRecVsSimZ = iB.book2D("heeRecVsSimZ","",7200,-360,360,7200,-360,360);
166  heeRecVsSimY = iB.book2D("heeRecVsSimY","",400,-200,200,400,-200,200);
167  heeRecVsSimX = iB.book2D("heeRecVsSimX","",400,-200,200,400,-200,200);
168 
169  hefdzVsZ = iB.book2D("hefdzVsZ","",8200,-410,410,100,-0.1,0.1);
170  hefdyVsY = iB.book2D("hefdyVsY","",400,-200,200,100,-0.02,0.02);
171  hefdxVsX = iB.book2D("hefdxVsX","",400,-200,200,100,-0.02,0.02);
172  hefRecVsSimZ = iB.book2D("hefRecVsSimZ","",8200,-410,410,8200,-410,410);
173  hefRecVsSimY = iB.book2D("hefRecVsSimY","",400,-200,200,400,-200,200);
174  hefRecVsSimX = iB.book2D("hefRecVsSimX","",400,-200,200,400,-200,200);
175 
176  hebdzVsZ = iB.book2D("hebdzVsZ","",1080,-540,540,100,-1.0,1.0);
177  hebdPhiVsPhi = iB.book2D("hebdPhiVsPhi","",M_PI*100,-0.5,M_PI+0.5,200,-0.2,0.2);
178  hebdEtaVsEta = iB.book2D("hebdEtaVsEta","",1000,-5,5,200,-0.1,0.1);
179  hebRecVsSimZ = iB.book2D("hebRecVsSimZ","",1080,-540,540,1080,-540,540);
180  hebRecVsSimY = iB.book2D("hebRecVsSimY","",400,-200,200,400,-200,200);
181  hebRecVsSimX = iB.book2D("hebRecVsSimX","",400,-200,200,400,-200,200);
182 
183  heeEnRec = iB.book1D("heeEnRec","",1000,0,10);
184  heeEnSim = iB.book1D("heeEnSim","",1000,0,0.01);
185  heeEnSimRec = iB.book2D("heeEnSimRec","",200,0,0.002,200,0,0.2);
186 
187  hefEnRec = iB.book1D("hefEnRec","",1000,0,10);
188  hefEnSim = iB.book1D("hefEnSim","",1000,0,0.01);
189  hefEnSimRec = iB.book2D("hefEnSimRec","",200,0,0.001,200,0,0.5);
190 
191  hebEnRec = iB.book1D("hebEnRec","",1000,0,15);
192  hebEnSim = iB.book1D("hebEnSim","",1000,0,0.01);
193  hebEnSimRec = iB.book2D("hebEnSimRec","",200,0,0.02,200,0,4);
194 
195 }
196 
198  edm::EventSetup const& iSetup) {
199  //initiating hgc Geometry
200  for (size_t i=0; i<geometrySource_.size(); i++) {
201  if (geometrySource_[i].find("Hcal") != std::string::npos) {
203  iSetup.get<HcalSimNumberingRecord>().get(pHSNDC);
204  if (pHSNDC.isValid()) {
205  hcCons_ = pHSNDC.product();
206  hgcCons_.push_back(nullptr);
207  } else {
208  edm::LogWarning("HGCalValid") << "Cannot initiate HcalDDDSimConstants: "
209  << geometrySource_[i] << std::endl;
210  }
212  iSetup.get<HcalRecNumberingRecord>().get(pHRNDC);
213  if (pHRNDC.isValid()) {
214  hcConr_ = pHRNDC.product();
215  } else {
216  edm::LogWarning("HGCalValid") << "Cannot initiate HcalDDDRecConstants: "
217  << geometrySource_[i] << std::endl;
218  }
220  iSetup.get<CaloGeometryRecord>().get(caloG);
221  if (caloG.isValid()) {
222  const CaloGeometry* geo = caloG.product();
224  hgcGeometry_.push_back(nullptr);
225  } else {
226  edm::LogWarning("HGCalValid") << "Cannot initiate HcalGeometry for "
227  << geometrySource_[i] << std::endl;
228  }
229  } else {
231  iSetup.get<IdealGeometryRecord>().get(geometrySource_[i],hgcCons);
232  if (hgcCons.isValid()) {
233  hgcCons_.push_back(hgcCons.product());
234  } else {
235  edm::LogWarning("HGCalValid") << "Cannot initiate HGCalDDDConstants for "
236  << geometrySource_[i] << std::endl;
237  }
239  iSetup.get<IdealGeometryRecord>().get(geometrySource_[i],hgcGeom);
240  if(hgcGeom.isValid()) {
241  hgcGeometry_.push_back(hgcGeom.product());
242  } else {
243  edm::LogWarning("HGCalValid") << "Cannot initiate HGCalGeometry for "
244  << geometrySource_[i] << std::endl;
245  }
246  }
247  }
248 }
249 
251  std::map<unsigned int, HGCHitTuple> eeHitRefs, fhHitRefs, bhHitRefs;
252 
253  //Accesing ee simhits
255  iEvent.getByToken(eeSimHitToken_, eeSimHits);
256 
257  if (eeSimHits.isValid()) {
258  analyzeHGCalSimHit(eeSimHits, 0, heeEnSim, eeHitRefs);
259 #ifdef EDM_ML_DEBUG
260  for (std::map<unsigned int,HGCHitTuple>::iterator itr=eeHitRefs.begin();
261  itr != eeHitRefs.end(); ++itr) {
262  int idx = std::distance(eeHitRefs.begin(),itr);
263  edm::LogInfo("HGCalValid") << "EEHit[" << idx << "] " << std::hex
264  << itr->first << std::dec << "; Energy "
265  << std::get<0>(itr->second)
266  << "; Position (" << std::get<1>(itr->second)
267  << ", " << std::get<2>(itr->second) <<", "
268  << std::get<3>(itr->second) << ")" <<std::endl;
269  }
270 #endif
271  } else {
272  edm::LogWarning("HGCalValid") << "No EE SimHit Found " << std::endl;
273  }
274 
275  //Accesing fh simhits
277  iEvent.getByToken(fhSimHitToken_, fhSimHits);
278  if (fhSimHits.isValid()) {
279  analyzeHGCalSimHit(fhSimHits, 1, hefEnSim, fhHitRefs);
280 #ifdef EDM_ML_DEBUG
281  for (std::map<unsigned int,HGCHitTuple>::iterator itr=fhHitRefs.begin();
282  itr != fhHitRefs.end(); ++itr) {
283  int idx = std::distance(fhHitRefs.begin(),itr);
284  edm::LogInfo("HGCalValid") << "FHHit[" << idx << "] " << std::hex
285  << itr->first << std::dec << "; Energy "
286  << std::get<0>(itr->second) << "; Position ("
287  << std::get<1>(itr->second) << ", "
288  << std::get<2>(itr->second) <<", "
289  << std::get<3>(itr->second) << ")" <<std::endl;
290  }
291 #endif
292  } else {
293  edm::LogWarning("HGCalValid") << "No FH SimHit Found " << std::endl;
294  }
295 
296  //Accessing bh simhits
298  iEvent.getByToken(bhSimHitToken_, bhSimHits);
299  if (bhSimHits.isValid()) {
300  for (std::vector<PCaloHit>::const_iterator simHit = bhSimHits->begin();
301  simHit != bhSimHits->end(); ++simHit) {
302  int subdet, z, depth, eta, phi, lay;
303  HcalTestNumbering::unpackHcalIndex(simHit->id(), subdet, z, depth, eta, phi, lay);
304 
305  if (subdet == static_cast<int>(HcalEndcap)) {
306  HcalCellType::HcalCell cell = hcCons_->cell(subdet, z, lay, eta, phi);
307  double zp = cell.rz/10;
308 
309  HcalDDDRecConstants::HcalID idx = hcConr_->getHCID(subdet,eta,phi,lay,depth);
310  int sign = (z==0)?(-1):(1);
311  zp *= sign;
312  HcalDetId id = HcalDetId(HcalEndcap,sign*idx.eta,idx.phi,idx.depth);
313 
314  float energy = simHit->energy();
315  float energySum(energy);
316  if (bhHitRefs.count(id.rawId()) != 0) energySum += std::get<0>(bhHitRefs[id.rawId()]);
317  hebEnSim->Fill(energy);
318  if (std::find(ietaExcludeBH_.begin(),ietaExcludeBH_.end(),idx.eta) ==
319  ietaExcludeBH_.end()) {
320  bhHitRefs[id.rawId()] = std::make_tuple(energySum,cell.eta,cell.phi,zp);
321 #ifdef EDM_ML_DEBUG
322  edm::LogInfo("HGCalValid") << "Accept " << id << std::endl;
323  } else {
324  edm::LogInfo("HGCalValid") << "Reject " << id << std::endl;
325 #endif
326  }
327  }
328  }
329 #ifdef EDM_ML_DEBUG
330  for (std::map<unsigned int,HGCHitTuple>::iterator itr=bhHitRefs.begin();
331  itr != bhHitRefs.end(); ++itr) {
332  int idx = std::distance(bhHitRefs.begin(),itr);
333  edm::LogInfo("HGCalValid") << "BHHit[" << idx << "] " << std::hex
334  << itr->first << std::dec << "; Energy "
335  << std::get<0>(itr->second) << "; Position ("
336  << std::get<1>(itr->second) << ", "
337  << std::get<2>(itr->second) <<", "
338  << std::get<3>(itr->second) << ")" <<std::endl;
339  }
340 #endif
341  } else {
342  edm::LogWarning("HGCalValid") << "No BH SimHit Found " << std::endl;
343  }
344 
345  //accessing EE Rechit information
347  iEvent.getByToken(eeRecHitToken_, eeRecHit);
348  if (eeRecHit.isValid()) {
349  const HGCeeRecHitCollection* theHits = (eeRecHit.product());
350  for (auto it = theHits->begin(); it != theHits->end(); ++it) {
351  double energy = it->energy();
352  heeEnRec->Fill(energy);
353  std::map<unsigned int, HGCHitTuple>::const_iterator itr = eeHitRefs.find(it->id().rawId());
354  if (itr != eeHitRefs.end()) {
355  GlobalPoint xyz = hgcGeometry_[0]->getPosition(it->id());
356  heeRecVsSimX->Fill(std::get<1>(itr->second),xyz.x());
357  heeRecVsSimY->Fill(std::get<2>(itr->second),xyz.y());
358  heeRecVsSimZ->Fill(std::get<3>(itr->second),xyz.z());
359  heedxVsX->Fill(std::get<1>(itr->second),(xyz.x()-std::get<1>(itr->second)));
360  heedyVsY->Fill(std::get<2>(itr->second),(xyz.y()-std::get<2>(itr->second)));
361  heedzVsZ->Fill(std::get<3>(itr->second),(xyz.z()-std::get<3>(itr->second)));
362  heeEnSimRec->Fill(std::get<0>(itr->second),energy);
363 #ifdef EDM_ML_DEBUG
364  edm::LogInfo("HGCalValid") << "EEHit: " << std::hex << it->id().rawId()
365  << std::dec << " Sim ("
366  << std::get<0>(itr->second) << ", "
367  << std::get<1>(itr->second) << ", "
368  << std::get<2>(itr->second) << ", "
369  << std::get<3>(itr->second) << ") Rec ("
370  << energy << ", " << xyz.x() << ", "
371  << xyz.y() << ", " << xyz.z() << ")\n";
372 #endif
373  }
374  }
375  } else {
376  edm::LogWarning("HGCalValid") << "No EE RecHit Found " << std::endl;
377  }
378 
379  //accessing FH Rechit information
381  iEvent.getByToken(fhRecHitToken_, fhRecHit);
382  if (fhRecHit.isValid()) {
383  const HGChefRecHitCollection* theHits = (fhRecHit.product());
384  for (auto it = theHits->begin(); it!=theHits->end(); ++it) {
385  double energy = it->energy();
386  hefEnRec->Fill(energy);
387  std::map<unsigned int, HGCHitTuple>::const_iterator itr = fhHitRefs.find(it->id().rawId());
388  if (itr != fhHitRefs.end()) {
389  GlobalPoint xyz = hgcGeometry_[1]->getPosition(it->id());
390 
391  hefRecVsSimX->Fill(std::get<1>(itr->second),xyz.x());
392  hefRecVsSimY->Fill(std::get<2>(itr->second),xyz.y());
393  hefRecVsSimZ->Fill(std::get<3>(itr->second),xyz.z());
394  hefdxVsX->Fill(std::get<1>(itr->second),(xyz.x()-std::get<1>(itr->second)));
395  hefdyVsY->Fill(std::get<2>(itr->second),(xyz.y()-std::get<2>(itr->second)));
396  hefdzVsZ->Fill(std::get<3>(itr->second),(xyz.z()-std::get<3>(itr->second)));
397  hefEnSimRec->Fill(std::get<0>(itr->second),energy);
398 #ifdef EDM_ML_DEBUG
399  edm::LogInfo("HGCalValid") << "FHHit: " << std::hex << it->id().rawId()
400  << std::dec << " Sim ("
401  << std::get<0>(itr->second) << ", "
402  << std::get<1>(itr->second) << ", "
403  << std::get<2>(itr->second) << ", "
404  << std::get<3>(itr->second) << ") Rec ("
405  << energy << "," << xyz.x() << ", "
406  << xyz.y() << ", " << xyz.z() << ")\n";
407 #endif
408  }
409  }
410  } else {
411  edm::LogWarning("HGCalValid") << "No FH RecHit Found " << std::endl;
412  }
413 
414 
415  //accessing BH Rechit information
416  if (ifHCAL_) {
418  iEvent.getByToken(bhRecHitTokenh_, bhRecHit);
419  if (bhRecHit.isValid()) {
420  const HBHERecHitCollection* theHits = (bhRecHit.product());
421  analyzeHGCalRecHit(theHits, bhHitRefs);
422  } else {
423  edm::LogWarning("HGCalValid") << "No BH RecHit Found " << std::endl;
424  }
425  } else {
427  iEvent.getByToken(bhRecHitTokeng_, bhRecHit);
428  if (bhRecHit.isValid()) {
429  const HGChebRecHitCollection* theHits = (bhRecHit.product());
430  analyzeHGCalRecHit(theHits, bhHitRefs);
431  } else {
432  edm::LogWarning("HGCalValid") << "No BH RecHit Found " << std::endl;
433  }
434  }
435 }
436 
437 void HGCalHitValidation::analyzeHGCalSimHit(edm::Handle<std::vector<PCaloHit>> const& simHits,
438  int idet, MonitorElement *hist,
439  std::map<unsigned int, HGCHitTuple>& hitRefs) {
440 
441  const HGCalTopology &hTopo=hgcGeometry_[idet]->topology();
442  for (std::vector<PCaloHit>::const_iterator simHit = simHits->begin(); simHit != simHits->end(); ++simHit) {
443  int subdet, zside, layer, wafer, celltype, cell;
444  HGCalTestNumbering::unpackHexagonIndex(simHit->id(), subdet, zside, layer, wafer, celltype, cell);
445  std::pair<float, float> xy = hgcCons_[idet]->locateCell(cell,layer,wafer,false);
446  float zp = hgcCons_[idet]->waferZ(layer,false);
447  if (zside < 0) zp = -zp;
448  float xp = (zp<0) ? -xy.first/10 : xy.first/10;
449  float yp = xy.second/10.0;
450 
451  //skip this hit if after ganging it is not valid
452  std::pair<int,int> recoLayerCell=hgcCons_[idet]->simToReco(cell,layer,wafer,hTopo.detectorType());
453  cell = recoLayerCell.first;
454  layer = recoLayerCell.second;
455 
456  //skip this hit if after ganging it is not valid
457  if (layer<0 || cell<0) {
458  } else {
459  //assign the RECO DetId
460  HGCalDetId id = HGCalDetId((ForwardSubdetector)(subdet),zside,layer,celltype,wafer,cell);
461  float energy = simHit->energy();
462 
463  float energySum(energy);
464  if (hitRefs.count(id.rawId()) != 0) energySum += std::get<0>(hitRefs[id.rawId()]);
465  hitRefs[id.rawId()] = std::make_tuple(energySum,xp,yp,zp);
466  hist->Fill(energy);
467  }
468  }
469 }
470 
471 template<class T1>
472 void HGCalHitValidation::analyzeHGCalRecHit(T1 const & theHits,
473  std::map<unsigned int, HGCHitTuple> const& hitRefs) {
474  for (auto it = theHits->begin(); it!=theHits->end(); ++it) {
475  DetId id = it->id();
476  if (id.subdetId() == (int)(HcalEndcap)) {
477  double energy = it->energy();
478  hebEnRec->Fill(energy);
479  GlobalPoint xyz = hcGeometry_->getGeometry(id)->getPosition();
480 
481  std::map<unsigned int, HGCHitTuple>::const_iterator itr = hitRefs.find(id.rawId());
482  if (itr != hitRefs.end()) {
483  float ang3 = xyz.phi().value(); // returns the phi in radians
484  double fac = sinh(std::get<1>(itr->second));
485  double pT = std::get<3>(itr->second) / fac;
486  double xp = pT * cos(std::get<2>(itr->second));
487  double yp = pT * sin(std::get<2>(itr->second));
488  hebRecVsSimX->Fill(xp,xyz.x());
489  hebRecVsSimY->Fill(yp,xyz.y());
490  hebRecVsSimZ->Fill(std::get<3>(itr->second),xyz.z());
491  hebdEtaVsEta->Fill(std::get<1>(itr->second),(xyz.eta()-std::get<1>(itr->second)));
492  hebdPhiVsPhi->Fill(std::get<2>(itr->second),(ang3-std::get<2>(itr->second)));
493  hebdzVsZ->Fill(std::get<3>(itr->second),(xyz.z()-std::get<3>(itr->second)));
494  hebEnSimRec->Fill(std::get<0>(itr->second),energy);
495 
496 #ifdef EDM_ML_DEBUG
497  edm::LogInfo("HGCalValid") << "BHHit: " << std::hex << id.rawId()
498  << std::dec << " Sim ("
499  << std::get<0>(itr->second) << ", "
500  << std::get<1>(itr->second) << ", "
501  << std::get<2>(itr->second) << ", "
502  << std::get<3>(itr->second) << ") Rec ("
503  << energy << ", " << xyz.x() << ", "
504  << xyz.y() << ", " << xyz.z() << ")\n";
505 #endif
506  }
507  }
508  }
509 }
510 
511 //define this as a plug-in
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:44
edm::EDGetTokenT< std::vector< PCaloHit > > fhSimHitToken_
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * hebRecVsSimX
MonitorElement * hefRecVsSimY
MonitorElement * hefRecVsSimX
edm::InputTag bhSimHitSource
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * hebRecVsSimY
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
MonitorElement * hefEnSimRec
MonitorElement * heedzVsZ
bool detectorType() const
int zside(DetId const &)
MonitorElement * heeRecVsSimX
HcalCellType::HcalCell cell(const int &det, const int &zside, const int &depth, const int &etaR, const int &iphi) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
ForwardSubdetector
std::vector< const HGCalDDDConstants * > hgcCons_
MonitorElement * heeRecVsSimZ
const HcalDDDSimConstants * hcCons_
MonitorElement * heeRecVsSimY
MonitorElement * heeEnRec
void Fill(long long x)
std::vector< const HGCalGeometry * > hgcGeometry_
MonitorElement * hebEnRec
edm::EDGetTokenT< HGCeeRecHitCollection > eeRecHitToken_
std::tuple< float, float, float, float > HGCHitTuple
edm::EDGetTokenT< std::vector< PCaloHit > > eeSimHitToken_
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
MonitorElement * hebEnSim
MonitorElement * hefdyVsY
MonitorElement * heeEnSimRec
void analyzeHGCalRecHit(T1 const &theHits, std::map< unsigned int, HGCHitTuple > const &hitRefs)
std::vector< int > ietaExcludeBH_
MonitorElement * hefEnSim
HcalID getHCID(int subdet, int ieta, int iphi, int lay, int idepth) const
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
MonitorElement * hebRecVsSimZ
std::vector< std::string > geometrySource_
MonitorElement * hefEnRec
const CaloSubdetectorGeometry * hcGeometry_
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
MonitorElement * heedxVsX
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * heeEnSim
edm::EDGetTokenT< std::vector< PCaloHit > > bhSimHitToken_
bool isValid() const
Definition: HandleBase.h:74
edm::EDGetTokenT< HBHERecHitCollection > bhRecHitTokenh_
MonitorElement * hebdPhiVsPhi
edm::InputTag fhSimHitSource
#define M_PI
int k[5][pyjets_maxn]
const_iterator end() const
MonitorElement * hebEnSimRec
Definition: DetId.h:18
edm::InputTag eeSimHitSource
HGCalHitValidation(const edm::ParameterSet &)
MonitorElement * hefdxVsX
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
T const * product() const
Definition: Handle.h:81
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
void dqmBeginRun(edm::Run const &, edm::EventSetup const &) override
const T & get() const
Definition: EventSetup.h:59
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
T eta() const
Definition: PV3DBase.h:76
const HcalDDDRecConstants * hcConr_
edm::EDGetTokenT< HGChebRecHitCollection > bhRecHitTokeng_
void analyzeHGCalSimHit(edm::Handle< std::vector< PCaloHit >> const &simHits, int idet, MonitorElement *hist, std::map< unsigned int, HGCHitTuple > &)
MonitorElement * hefRecVsSimZ
MonitorElement * hefdzVsZ
double energySum(const DataFrame &df, int fs, int ls)
bool isValid() const
Definition: ESHandle.h:47
MonitorElement * hebdzVsZ
T x() const
Definition: PV3DBase.h:62
T const * product() const
Definition: ESHandle.h:86
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:43
MonitorElement * hebdEtaVsEta
edm::EDGetTokenT< HGChefRecHitCollection > fhRecHitToken_