CMS 3D CMS Logo

HGCalSimHitValidation.cc
Go to the documentation of this file.
1 // system include files
2 #include <cmath>
3 #include <iostream>
4 #include <fstream>
5 #include <vector>
6 #include <map>
7 #include <string>
8 
14 
21 
24 
35 
42 
47 
48 #include "CLHEP/Geometry/Point3D.h"
49 #include "CLHEP/Geometry/Transform3D.h"
50 #include "CLHEP/Geometry/Vector3D.h"
51 #include "CLHEP/Units/GlobalSystemOfUnits.h"
52 #include "CLHEP/Units/GlobalPhysicalConstants.h"
53 
55 public:
56  struct energysum {
58  etotal = 0;
59  for (int i = 0; i < 6; ++i)
60  eTime[i] = 0.;
61  }
62  double eTime[6], etotal;
63  };
64 
65  struct hitsinfo {
67  x = y = z = phi = eta = 0.0;
68  cell = cell2 = sector = sector2 = type = layer = 0;
69  }
70  double x, y, z, phi, eta;
71  int cell, cell2, sector, sector2, type, layer;
72  };
73 
75  ~HGCalSimHitValidation() override {}
76 
77  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
78 
79 protected:
80  void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
81  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
82  void analyze(const edm::Event&, const edm::EventSetup&) override;
83 
84 private:
85  void analyzeHits(std::vector<PCaloHit>& hits);
86  void fillOccupancyMap(std::map<int, int>& OccupancyMap, int layer);
87  void fillHitsInfo(std::pair<hitsinfo, energysum> hit_, unsigned int itimeslice, double esum);
89 
90  // ----------member data ---------------------------
94  std::vector<double> times_;
99  unsigned int layers_;
101  std::map<uint32_t, HepGeom::Transform3D> transMap_;
102 
103  std::vector<MonitorElement*> HitOccupancy_Plus_, HitOccupancy_Minus_;
104  std::vector<MonitorElement*> EtaPhi_Plus_, EtaPhi_Minus_;
106  static const unsigned int maxTime_ = 6;
107  std::vector<MonitorElement*> energy_[maxTime_];
108  unsigned int nTimes_;
109 };
110 
112  : nameDetector_(iConfig.getParameter<std::string>("DetectorName")),
113  caloHitSource_(iConfig.getParameter<std::string>("CaloHitSource")),
114  times_(iConfig.getParameter<std::vector<double> >("TimeSlices")),
115  verbosity_(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
116  testNumber_(iConfig.getUntrackedParameter<bool>("TestNumber", true)),
117  symmDet_(true),
118  firstLayer_(1) {
119  heRebuild_ = (nameDetector_ == "HCal") ? true : false;
120  tok_hepMC_ = consumes<edm::HepMCProduct>(edm::InputTag("generatorSmeared"));
121  tok_hits_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", caloHitSource_));
122  nTimes_ = (times_.size() > maxTime_) ? maxTime_ : times_.size();
123 }
124 
127  std::vector<double> times = {25.0, 1000.0};
128  desc.add<std::string>("DetectorName", "HGCalEESensitive");
129  desc.add<std::string>("CaloHitSource", "HGCHitsEE");
130  desc.add<std::vector<double> >("TimeSlices", times);
131  desc.addUntracked<int>("Verbosity", 0);
132  desc.addUntracked<bool>("TestNumber", true);
133  descriptions.add("hgcalSimHitValidationEE", desc);
134 }
135 
137  //Generator input
138  if (verbosity_ > 0) {
140  iEvent.getByToken(tok_hepMC_, evtMC);
141  if (!evtMC.isValid()) {
142  edm::LogVerbatim("HGCalValidation") << "no HepMCProduct found";
143  } else {
144  const HepMC::GenEvent* myGenEvent = evtMC->GetEvent();
145  unsigned int k(0);
146  for (HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end();
147  ++p, ++k) {
148  edm::LogVerbatim("HGCalValidation") << "Particle[" << k << "] with pt " << (*p)->momentum().perp() << " eta "
149  << (*p)->momentum().eta() << " phi " << (*p)->momentum().phi();
150  }
151  }
152  }
153 
154  //Now the hits
155  edm::Handle<edm::PCaloHitContainer> theCaloHitContainers;
156  iEvent.getByToken(tok_hits_, theCaloHitContainers);
157  if (theCaloHitContainers.isValid()) {
158  if (verbosity_ > 0)
159  edm::LogVerbatim("HGCalValidation") << " PcalohitItr = " << theCaloHitContainers->size();
160  std::vector<PCaloHit> caloHits;
161  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
162  if (heRebuild_ && testNumber_) {
163  for (unsigned int i = 0; i < caloHits.size(); ++i) {
164  unsigned int id_ = caloHits[i].id();
166  if (hid.subdet() != int(HcalEndcap))
167  hid = HcalDetId(HcalEmpty, hid.ieta(), hid.iphi(), hid.depth());
168  caloHits[i].setID(hid.rawId());
169  if (verbosity_ > 0)
170  edm::LogVerbatim("HGCalValidation") << "Hit[" << i << "] " << hid;
171  }
172  }
173  analyzeHits(caloHits);
174  } else if (verbosity_ > 0) {
175  edm::LogVerbatim("HGCalValidation") << "PCaloHitContainer does not exist!";
176  }
177 }
178 
179 void HGCalSimHitValidation::analyzeHits(std::vector<PCaloHit>& hits) {
180  std::map<int, int> OccupancyMap_plus, OccupancyMap_minus;
181  OccupancyMap_plus.clear();
182  OccupancyMap_minus.clear();
183 
184  std::map<uint32_t, std::pair<hitsinfo, energysum> > map_hits;
185  map_hits.clear();
186 
187  if (verbosity_ > 0)
188  edm::LogVerbatim("HGCalValidation") << nameDetector_ << " with " << hits.size() << " PcaloHit elements";
189  unsigned int nused(0);
190  for (unsigned int i = 0; i < hits.size(); i++) {
191  double energy = hits[i].energy();
192  double time = hits[i].time();
193  uint32_t id_ = hits[i].id();
194  int cell, sector, subsector(0), layer, zside;
195  int subdet(0), cell2(0), type(0);
196  if (heRebuild_) {
197  HcalDetId detId = HcalDetId(id_);
198  subdet = detId.subdet();
199  if (subdet != static_cast<int>(HcalEndcap))
200  continue;
201  cell = detId.ietaAbs();
202  sector = detId.iphi();
203  subsector = 1;
204  layer = detId.depth();
205  zside = detId.zside();
206  } else if ((hgcons_->geomMode() == HGCalGeometryMode::Hexagon8) ||
208  HGCSiliconDetId detId = HGCSiliconDetId(id_);
209  subdet = ForwardEmpty;
210  cell = detId.cellU();
211  cell2 = detId.cellV();
212  sector = detId.waferU();
213  subsector = detId.waferV();
214  type = detId.type();
215  layer = detId.layer();
216  zside = detId.zside();
217  } else if (hgcons_->geomMode() == HGCalGeometryMode::Square) {
218  HGCalTestNumbering::unpackSquareIndex(id_, zside, layer, sector, subsector, cell);
219  } else if (hgcons_->geomMode() == HGCalGeometryMode::Trapezoid) {
221  subdet = ForwardEmpty;
222  sector = detId.ietaAbs();
223  cell = detId.iphi();
224  subsector = 1;
225  type = detId.type();
226  layer = detId.layer();
227  zside = detId.zside();
228  } else {
229  HGCalTestNumbering::unpackHexagonIndex(id_, subdet, zside, layer, sector, type, cell);
230  }
231  nused++;
232  if (verbosity_ > 1)
233  edm::LogVerbatim("HGCalValidation")
234  << "Detector " << nameDetector_ << " zside = " << zside << " sector|wafer = " << sector << ":" << subsector
235  << " type = " << type << " layer = " << layer << " cell = " << cell << ":" << cell2 << " energy = " << energy
236  << " energyem = " << hits[i].energyEM() << " energyhad = " << hits[i].energyHad() << " time = " << time;
237 
238  HepGeom::Point3D<float> gcoord;
239  if (heRebuild_) {
240  std::pair<double, double> etaphi = hcons_->getEtaPhi(subdet, zside * cell, sector);
241  double rz = hcons_->getRZ(subdet, zside * cell, layer);
242  if (verbosity_ > 2)
243  edm::LogVerbatim("HGCalValidation") << "i/p " << subdet << ":" << zside << ":" << cell << ":" << sector << ":"
244  << layer << " o/p " << etaphi.first << ":" << etaphi.second << ":" << rz;
245  gcoord = HepGeom::Point3D<float>(rz * cos(etaphi.second) / cosh(etaphi.first),
246  rz * sin(etaphi.second) / cosh(etaphi.first),
247  rz * tanh(etaphi.first));
248  } else if (hgcons_->geomMode() == HGCalGeometryMode::Square) {
249  std::pair<float, float> xy = hgcons_->locateCell(cell, layer, subsector, false);
250  const HepGeom::Point3D<float> lcoord(xy.first, xy.second, 0);
251  int subs = (symmDet_ ? 0 : subsector);
252  id_ = HGCalTestNumbering::packSquareIndex(zside, layer, sector, subs, 0);
253  gcoord = (transMap_[id_] * lcoord);
254  } else {
255  std::pair<float, float> xy;
258  xy = hgcons_->locateCell(layer, sector, subsector, cell, cell2, false, true);
259  } else if (hgcons_->geomMode() == HGCalGeometryMode::Trapezoid) {
260  xy = hgcons_->locateCellTrap(layer, sector, cell, false);
261  } else {
262  xy = hgcons_->locateCell(cell, layer, sector, false);
263  }
264  double zp = hgcons_->waferZ(layer, false);
265  if (zside < 0)
266  zp = -zp;
267  float xp = (zp < 0) ? -xy.first : xy.first;
268  gcoord = HepGeom::Point3D<float>(xp, xy.second, zp);
269  }
270  double tof = (gcoord.mag() * CLHEP::mm) / CLHEP::c_light;
271  if (verbosity_ > 1)
272  edm::LogVerbatim("HGCalValidation")
273  << std::hex << id_ << std::dec << " global coordinate " << gcoord << " time " << time << ":" << tof;
274  time -= tof;
275 
276  energysum esum;
277  hitsinfo hinfo;
278  if (map_hits.count(id_) != 0) {
279  hinfo = map_hits[id_].first;
280  esum = map_hits[id_].second;
281  } else {
282  hinfo.x = gcoord.x();
283  hinfo.y = gcoord.y();
284  hinfo.z = gcoord.z();
285  hinfo.sector = sector;
286  hinfo.sector2 = subsector;
287  hinfo.cell = cell;
288  hinfo.cell2 = cell2;
289  hinfo.type = type;
290  hinfo.layer = layer - firstLayer_;
291  hinfo.phi = gcoord.getPhi();
292  hinfo.eta = gcoord.getEta();
293  }
294  esum.etotal += energy;
295  for (unsigned int k = 0; k < nTimes_; ++k) {
296  if (time > 0 && time < times_[k])
297  esum.eTime[k] += energy;
298  }
299 
300  if (verbosity_ > 1)
301  edm::LogVerbatim("HGCalValidation") << " ----------------------- gx = " << hinfo.x << " gy = " << hinfo.y
302  << " gz = " << hinfo.z << " phi = " << hinfo.phi << " eta = " << hinfo.eta;
303  map_hits[id_] = std::pair<hitsinfo, energysum>(hinfo, esum);
304  }
305  if (verbosity_ > 0)
306  edm::LogVerbatim("HGCalValidation") << nameDetector_ << " with " << map_hits.size()
307  << " detector elements being hit";
308 
309  std::map<uint32_t, std::pair<hitsinfo, energysum> >::iterator itr;
310  for (itr = map_hits.begin(); itr != map_hits.end(); ++itr) {
311  hitsinfo hinfo = (*itr).second.first;
312  energysum esum = (*itr).second.second;
313  int layer = hinfo.layer;
314  double eta = hinfo.eta;
315 
316  for (unsigned int itimeslice = 0; itimeslice < nTimes_; itimeslice++) {
317  fillHitsInfo((*itr).second, itimeslice, esum.eTime[itimeslice]);
318  }
319 
320  if (eta > 0.0)
321  fillOccupancyMap(OccupancyMap_plus, layer);
322  else
323  fillOccupancyMap(OccupancyMap_minus, layer);
324  }
325  if (verbosity_ > 0)
326  edm::LogVerbatim("HGCalValidation") << "With map:used:total " << hits.size() << "|" << nused << "|"
327  << map_hits.size() << " hits";
328 
329  for (auto const& itr : OccupancyMap_plus) {
330  int layer = itr.first;
331  int occupancy = itr.second;
332  HitOccupancy_Plus_.at(layer)->Fill(occupancy);
333  }
334  for (auto const& itr : OccupancyMap_minus) {
335  int layer = itr.first;
336  int occupancy = itr.second;
337  HitOccupancy_Minus_.at(layer)->Fill(occupancy);
338  }
339 }
340 
341 void HGCalSimHitValidation::fillOccupancyMap(std::map<int, int>& OccupancyMap, int layer) {
342  if (OccupancyMap.find(layer) != OccupancyMap.end()) {
343  ++OccupancyMap[layer];
344  } else {
345  OccupancyMap[layer] = 1;
346  }
347 }
348 
349 void HGCalSimHitValidation::fillHitsInfo(std::pair<hitsinfo, energysum> hits, unsigned int itimeslice, double esum) {
350  unsigned int ilayer = hits.first.layer;
351  if (ilayer < layers_) {
352  energy_[itimeslice].at(ilayer)->Fill(esum);
353  if (itimeslice == 0) {
354  EtaPhi_Plus_.at(ilayer)->Fill(hits.first.eta, hits.first.phi);
355  EtaPhi_Minus_.at(ilayer)->Fill(hits.first.eta, hits.first.phi);
356  }
357  } else {
358  if (verbosity_ > 0)
359  edm::LogVerbatim("HGCalValidation")
360  << "Problematic Hit for " << nameDetector_ << " at sector " << hits.first.sector << ":" << hits.first.sector2
361  << " layer " << hits.first.layer << " cell " << hits.first.cell << ":" << hits.first.cell2 << " energy "
362  << hits.second.etotal;
363  }
364 }
365 
367  if (verbosity_ > 0)
368  edm::LogVerbatim("HGCalValidation") << "Initialize HGCalDDDConstants for " << nameDetector_ << " : " << hgcons_;
369 
370  if (hgcons_->geomMode() == HGCalGeometryMode::Square) {
371  const DDCompactView& cview = *ddViewH;
372  std::string attribute = "Volume";
374 
375  DDSpecificsMatchesValueFilter filter{DDValue(attribute, value, 0)};
376  DDFilteredView fv(cview, filter);
377  bool dodet = fv.firstChild();
378 
379  while (dodet) {
380  const DDSolid& sol = fv.logicalPart().solid();
381  const std::string& name = sol.name().fullname();
382  int isd = (name.find(nameDetector_) == std::string::npos) ? -1 : 1;
383  if (isd > 0) {
384  std::vector<int> copy = fv.copyNumbers();
385  int nsiz = static_cast<int>(copy.size());
386  int lay = (nsiz > 0) ? copy[nsiz - 1] : -1;
387  int sec = (nsiz > 1) ? copy[nsiz - 2] : -1;
388  int zp = (nsiz > 3) ? copy[nsiz - 4] : -1;
389  if (zp != 1)
390  zp = -1;
391  const DDTrap& trp = static_cast<DDTrap>(sol);
392  int subs = (trp.alpha1() > 0 ? 1 : 0);
393  symmDet_ = (trp.alpha1() == 0 ? true : false);
394  uint32_t id = HGCalTestNumbering::packSquareIndex(zp, lay, sec, subs, 0);
395  DD3Vector x, y, z;
396  fv.rotation().GetComponents(x, y, z);
397  const CLHEP::HepRep3x3 rotation(x.X(), y.X(), z.X(), x.Y(), y.Y(), z.Y(), x.Z(), y.Z(), z.Z());
398  const CLHEP::HepRotation hr(rotation);
399  const CLHEP::Hep3Vector h3v(fv.translation().X(), fv.translation().Y(), fv.translation().Z());
400  const HepGeom::Transform3D ht3d(hr, h3v);
401  transMap_.insert(std::make_pair(id, ht3d));
402  if (verbosity_ > 2)
403  edm::LogVerbatim("HGCalValidation") << HGCalDetId(id) << " Transform using " << h3v << " and " << hr;
404  }
405  dodet = fv.next();
406  }
407  if (verbosity_ > 0)
408  edm::LogVerbatim("HGCalValidation") << "Finds " << transMap_.size() << " elements and SymmDet_ = " << symmDet_;
409  }
410  return true;
411 }
412 
413 // ------------ method called when starting to processes a run ------------
415  if (heRebuild_) {
417  iSetup.get<HcalRecNumberingRecord>().get(pHRNDC);
418  hcons_ = &(*pHRNDC);
419  layers_ = hcons_->getMaxDepth(1);
420  } else {
422  iSetup.get<IdealGeometryRecord>().get(nameDetector_, pHGDC);
423  hgcons_ = &(*pHGDC);
424  layers_ = hgcons_->layers(false);
427  iSetup.get<IdealGeometryRecord>().get(pDD);
428  defineGeometry(pDD);
429  }
430  if (verbosity_ > 0)
431  edm::LogVerbatim("HGCalValidation") << nameDetector_ << " defined with " << layers_ << " Layers with first at "
432  << firstLayer_;
433 }
434 
436  iB.setCurrentFolder("HGCAL/HGCalSimHitsV/" + nameDetector_);
437 
438  std::ostringstream histoname;
439  for (unsigned int il = 0; il < layers_; ++il) {
440  int ilayer = firstLayer_ + static_cast<int>(il);
441  auto istr1 = std::to_string(ilayer);
442  while (istr1.size() < 2) {
443  istr1.insert(0, "0");
444  }
445  histoname.str("");
446  histoname << "HitOccupancy_Plus_layer_" << istr1;
447  HitOccupancy_Plus_.push_back(iB.book1D(histoname.str().c_str(), "HitOccupancy_Plus", 501, -0.5, 500.5));
448  histoname.str("");
449  histoname << "HitOccupancy_Minus_layer_" << istr1;
450  HitOccupancy_Minus_.push_back(iB.book1D(histoname.str().c_str(), "HitOccupancy_Minus", 501, -0.5, 500.5));
451 
452  histoname.str("");
453  histoname << "EtaPhi_Plus_"
454  << "layer_" << istr1;
455  EtaPhi_Plus_.push_back(iB.book2D(histoname.str().c_str(), "Occupancy", 31, 1.45, 3.0, 72, -CLHEP::pi, CLHEP::pi));
456  histoname.str("");
457  histoname << "EtaPhi_Minus_"
458  << "layer_" << istr1;
459  EtaPhi_Minus_.push_back(
460  iB.book2D(histoname.str().c_str(), "Occupancy", 31, -3.0, -1.45, 72, -CLHEP::pi, CLHEP::pi));
461 
462  for (unsigned int itimeslice = 0; itimeslice < nTimes_; itimeslice++) {
463  histoname.str("");
464  histoname << "energy_time_" << itimeslice << "_layer_" << istr1;
465  energy_[itimeslice].push_back(iB.book1D(histoname.str().c_str(), "energy_", 100, 0, 0.1));
466  }
467  }
468 
469  MeanHitOccupancy_Plus_ = iB.book1D("MeanHitOccupancy_Plus", "MeanHitOccupancy_Plus", layers_, 0.5, layers_ + 0.5);
470  MeanHitOccupancy_Minus_ = iB.book1D("MeanHitOccupancy_Minus", "MeanHitOccupancy_Minus", layers_, 0.5, layers_ + 0.5);
471 }
472 
474 //define this as a plug-in
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::pair< T, T > etaphi(T x, T y, T z)
Definition: FastMath.h:162
MonitorElement * MeanHitOccupancy_Plus_
type
Definition: HCALResponse.h:21
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const N & name() const
Definition: DDBase.h:59
void fillOccupancyMap(std::map< int, int > &OccupancyMap, int layer)
std::pair< double, double > getEtaPhi(const int &subdet, const int &ieta, const int &iphi) const
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:138
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
std::vector< MonitorElement * > HitOccupancy_Plus_
void analyze(const edm::Event &, const edm::EventSetup &) override
ProductID id() const
Definition: HandleBase.cc:13
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
nav_type copyNumbers() const
return the stack of copy numbers
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalDetId.h:141
int waferU() const
int cellV() const
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
edm::EDGetTokenT< edm::HepMCProduct > tok_hepMC_
const DDRotationMatrix & rotation() const
The absolute rotation of the current node.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
int type() const
get the type
void analyzeHits(std::vector< PCaloHit > &hits)
int zside() const
get the z-side of the cell (1/-1)
int zside(DetId const &)
std::vector< MonitorElement * > HitOccupancy_Minus_
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
int cellU() const
get the cell #&#39;s in u,v or in x,y
A DDSolid represents the shape of a part.
Definition: DDSolid.h:39
MonitorElement * MeanHitOccupancy_Minus_
const Double_t pi
bool defineGeometry(edm::ESTransientHandle< DDCompactView > &ddViewH)
int depth() const
get the tower depth
Definition: HcalDetId.h:164
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int layers(bool reco) const
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DD3Vector
A DD Translation is currently implemented with Root Vector3D.
Definition: DDTranslation.h:6
bool next()
set current node to the next node in the filtered tree
const HcalDDDRecConstants * hcons_
int type() const
get the type
susybsm::HSCParticleRef hr
Definition: classes.h:26
int layer() const
get the layer #
std::pair< float, float > locateCellTrap(int lay, int ieta, int iphi, bool reco) const
double getRZ(const int &subdet, const int &ieta, const int &depth) const
int waferV() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static const unsigned int maxTime_
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
int ieta() const
get the cell ieta
Definition: HcalDetId.h:155
const std::string fullname() const
Definition: DDName.h:40
Interface to a Trapezoid.
Definition: DDSolid.h:77
int iphi() const
get the phi index
HGCalGeometryMode::GeometryMode geomMode() const
Definition: value.py:1
std::vector< MonitorElement * > EtaPhi_Plus_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
HGCalSimHitValidation(const edm::ParameterSet &)
double waferZ(int layer, bool reco) const
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:148
int iphi() const
get the cell iphi
Definition: HcalDetId.h:157
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
int getMaxDepth(const int &type) const
static void unpackSquareIndex(const uint32_t &idx, int &z, int &lay, int &sec, int &subsec, int &cell)
double alpha1(void) const
Angle with respect to the y axis from the centre of the side at y=-pDy1 to the centre at y=+pDy1 of t...
Definition: DDSolid.cc:143
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< MonitorElement * > energy_[maxTime_]
int layer() const
get the layer #
int zside() const
get the z-side of the cell (1/-1)
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Definition: DQMStore.cc:266
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
T get() const
Definition: EventSetup.h:73
std::map< uint32_t, HepGeom::Transform3D > transMap_
bool firstChild()
set the current node to the first child ...
const HGCalDDDConstants * hgcons_
int firstLayer() const
std::vector< MonitorElement * > EtaPhi_Minus_
const DDTranslation & translation() const
The absolute translation of the current node.
DetId relabel(const uint32_t testId) const
static uint32_t packSquareIndex(int z, int lay, int sec, int subsec, int cell)
void fillHitsInfo(std::pair< hitsinfo, energysum > hit_, unsigned int itimeslice, double esum)
std::vector< double > times_
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
Definition: Run.h:45