CMS 3D CMS Logo

EcalSimHitStudy.cc
Go to the documentation of this file.
3 
7 
13 
16 
22 
29 
33 
34 #include <TH1F.h>
35 #include <TH2F.h>
36 #include <cmath>
37 #include <iostream>
38 #include <fstream>
39 #include <map>
40 #include <memory>
41 #include <string>
42 #include <vector>
43 
44 //#define EDM_ML_DEBUG
45 
46 class EcalSimHitStudy : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
47 public:
49  ~EcalSimHitStudy() override {}
50  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
51 
52 protected:
53  void beginJob() override;
54  void analyze(edm::Event const&, edm::EventSetup const&) override;
55  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
56  void endRun(edm::Run const&, edm::EventSetup const&) override {}
57  void analyzeHits(std::vector<PCaloHit>&, int);
58 
59 private:
60  struct EcalHit {
61  uint16_t id;
62  double time, energy;
63  EcalHit(uint16_t i = 0, double t = 0, double e = 0) : id(i), time(t), energy(e) {}
64  };
65  static const int ndets_ = 2;
68  const std::vector<std::string> hitLab_;
69  const double maxEnergy_, tmax_, w0_;
70  const int selX_;
72  const std::vector<edm::EDGetTokenT<edm::PCaloHitContainer>> toks_calo_;
81 };
82 
85  g4Label_(ps.getUntrackedParameter<std::string>("ModuleLabel", "g4SimHits")),
86  hitLab_(ps.getUntrackedParameter<std::vector<std::string>>("CaloCollection")),
87  maxEnergy_(ps.getUntrackedParameter<double>("MaxEnergy", 200.0)),
88  tmax_(ps.getUntrackedParameter<double>("TimeCut", 100.0)),
89  w0_(ps.getUntrackedParameter<double>("W0", 4.7)),
90  selX_(ps.getUntrackedParameter<int>("SelectX", -1)),
91  tok_evt_(consumes<edm::HepMCProduct>(
92  edm::InputTag(ps.getUntrackedParameter<std::string>("SourceLabel", "VtxSmeared")))),
93  toks_calo_{edm::vector_transform(hitLab_, [this](const std::string& name) {
94  return consumes<edm::PCaloHitContainer>(edm::InputTag{g4Label_, name});
95  })} {
96  usesResource(TFileService::kSharedResource);
97 
98  edm::LogVerbatim("HitStudy") << "Module Label: " << g4Label_ << " Hits: " << hitLab_[0] << ", " << hitLab_[1]
99  << " MaxEnergy: " << maxEnergy_ << " Tmax: " << tmax_ << " Select " << selX_;
100 }
101 
104  std::vector<std::string> calonames;
105  desc.addUntracked<std::string>("ModuleLabel", "g4SimHits");
106  desc.addUntracked<std::vector<std::string>>("CaloCollection", calonames);
107  desc.addUntracked<std::string>("SourceLabel", "VtxSmeared");
108  desc.addUntracked<double>("MaxEnergy", 200.0);
109  desc.addUntracked<double>("TimeCut", 100.0);
110  desc.addUntracked<int>("SelectX", -1);
111  descriptions.add("EcalSimHitStudy", desc);
112 }
113 
116  if (!tfile.isAvailable())
117  throw cms::Exception("BadConfig") << "TFileService unavailable: "
118  << "please add it to config file";
119  char name[20], title[120];
120  sprintf(title, "Incident PT (GeV)");
121  ptInc_ = tfile->make<TH1F>("PtInc", title, 1000, 0., maxEnergy_);
122  ptInc_->GetXaxis()->SetTitle(title);
123  ptInc_->GetYaxis()->SetTitle("Events");
124  sprintf(title, "Incident Energy (GeV)");
125  eneInc_ = tfile->make<TH1F>("EneInc", title, 1000, 0., maxEnergy_);
126  eneInc_->GetXaxis()->SetTitle(title);
127  eneInc_->GetYaxis()->SetTitle("Events");
128  sprintf(title, "Incident #eta");
129  etaInc_ = tfile->make<TH1F>("EtaInc", title, 200, -5., 5.);
130  etaInc_->GetXaxis()->SetTitle(title);
131  etaInc_->GetYaxis()->SetTitle("Events");
132  sprintf(title, "Incident #phi");
133  phiInc_ = tfile->make<TH1F>("PhiInc", title, 200, -3.1415926, 3.1415926);
134  phiInc_->GetXaxis()->SetTitle(title);
135  phiInc_->GetYaxis()->SetTitle("Events");
136  std::string dets[ndets_] = {"EB", "EE"};
137  for (int i = 0; i < ndets_; i++) {
138  sprintf(name, "Hit%d", i);
139  sprintf(title, "Number of hits in %s", dets[i].c_str());
140  hit_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 1000, 0., 20000.);
141  hit_[i]->GetXaxis()->SetTitle(title);
142  hit_[i]->GetYaxis()->SetTitle("Events");
143  hit_[i]->Sumw2();
144  sprintf(name, "Time%d", i);
145  sprintf(title, "Time of the hit (ns) in %s", dets[i].c_str());
146  time_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 1000, 0., 1000.);
147  time_[i]->GetXaxis()->SetTitle(title);
148  time_[i]->GetYaxis()->SetTitle("Hits");
149  time_[i]->Sumw2();
150  sprintf(name, "TimeAll%d", i);
151  sprintf(title, "Hit time (ns) in %s (for first hit in the cell)", dets[i].c_str());
152  timeAll_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 1000, 0., 1000.);
153  timeAll_[i]->GetXaxis()->SetTitle(title);
154  timeAll_[i]->GetYaxis()->SetTitle("Hits");
155  timeAll_[i]->Sumw2();
156  sprintf(name, "Edep%d", i);
157  sprintf(title, "Energy deposit (GeV) in %s", dets[i].c_str());
158  edep_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 5000, 0., maxEnergy_);
159  edep_[i]->GetXaxis()->SetTitle(title);
160  edep_[i]->GetYaxis()->SetTitle("Hits");
161  edep_[i]->Sumw2();
162  sprintf(name, "EdepAll%d", i);
163  sprintf(title, "Total Energy deposit in the cell (GeV) in %s", dets[i].c_str());
164  edepAll_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 5000, 0., maxEnergy_);
165  edepAll_[i]->GetXaxis()->SetTitle(title);
166  edepAll_[i]->GetYaxis()->SetTitle("Hits");
167  edepAll_[i]->Sumw2();
168  sprintf(name, "EdepEM%d", i);
169  sprintf(title, "Energy deposit (GeV) by EM particles in %s", dets[i].c_str());
170  edepEM_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 5000, 0., maxEnergy_);
171  edepEM_[i]->GetXaxis()->SetTitle(title);
172  edepEM_[i]->GetYaxis()->SetTitle("Hits");
173  edepEM_[i]->Sumw2();
174  sprintf(name, "EdepHad%d", i);
175  sprintf(title, "Energy deposit (GeV) by hadrons in %s", dets[i].c_str());
176  edepHad_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 5000, 0., maxEnergy_);
177  edepHad_[i]->GetXaxis()->SetTitle(title);
178  edepHad_[i]->GetYaxis()->SetTitle("Hits");
179  edepHad_[i]->Sumw2();
180  sprintf(name, "Etot%d", i);
181  sprintf(title, "Total energy deposit (GeV) in %s", dets[i].c_str());
182  etot_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 5000, 0., maxEnergy_);
183  etot_[i]->GetXaxis()->SetTitle(title);
184  etot_[i]->GetYaxis()->SetTitle("Events");
185  etot_[i]->Sumw2();
186  sprintf(name, "EtotG%d", i);
187  sprintf(title, "Total energy deposit (GeV) in %s (t < 100 ns)", dets[i].c_str());
188  etotg_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 5000, 0., maxEnergy_);
189  etotg_[i]->GetXaxis()->SetTitle(title);
190  etotg_[i]->GetYaxis()->SetTitle("Events");
191  etotg_[i]->Sumw2();
192  sprintf(name, "r1by9%d", i);
193  sprintf(title, "E1/E9 in %s", dets[i].c_str());
194  r1by9_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 100, 0.0, 1.0);
195  r1by9_[i]->GetXaxis()->SetTitle(title);
196  r1by9_[i]->GetYaxis()->SetTitle("Events");
197  r1by9_[i]->Sumw2();
198  sprintf(name, "r1by25%d", i);
199  sprintf(title, "E1/E25 in %s", dets[i].c_str());
200  r1by25_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 100, 0.0, 1.0);
201  r1by25_[i]->GetXaxis()->SetTitle(title);
202  r1by25_[i]->GetYaxis()->SetTitle("Events");
203  r1by25_[i]->Sumw2();
204  sprintf(name, "r9by25%d", i);
205  sprintf(title, "E9/E25 in %s", dets[i].c_str());
206  r9by25_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 100, 0.0, 1.0);
207  r9by25_[i]->GetXaxis()->SetTitle(title);
208  r9by25_[i]->GetYaxis()->SetTitle("Events");
209  r9by25_[i]->Sumw2();
210  double ymax = (i == 0) ? 0.0005 : 0.005;
211  sprintf(name, "sEtaEta%d", i);
212  sprintf(title, "Cov(#eta,#eta) in %s", dets[i].c_str());
213  sEtaEta_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 1000, 0.0, ymax);
214  sEtaEta_[i]->GetXaxis()->SetTitle(title);
215  sEtaEta_[i]->GetYaxis()->SetTitle("Events");
216  sEtaEta_[i]->Sumw2();
217  sprintf(name, "sEtaPhi%d", i);
218  sprintf(title, "Cov(#eta,#phi) in %s", dets[i].c_str());
219  sEtaPhi_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 1000, 0.0, ymax);
220  sEtaPhi_[i]->GetXaxis()->SetTitle(title);
221  sEtaPhi_[i]->GetYaxis()->SetTitle("Events");
222  sEtaPhi_[i]->Sumw2();
223  ymax = (i == 0) ? 0.001 : 0.01;
224  sprintf(name, "sPhiPhi%d", i);
225  sprintf(title, "Cov(#phi,#phi) in %s", dets[i].c_str());
226  sPhiPhi_[i] = tfile->make<TH1F>(name, dets[i].c_str(), 1000, 0.0, ymax);
227  sPhiPhi_[i]->GetXaxis()->SetTitle(title);
228  sPhiPhi_[i]->GetYaxis()->SetTitle("Events");
229  sPhiPhi_[i]->Sumw2();
230  if (i == 0) {
231  sprintf(title, "%s+", dets[i].c_str());
232  poszp_[i] = tfile->make<TH2F>("poszp0", title, 100, 0, 100, 360, 0, 360);
233  poszp_[i]->GetXaxis()->SetTitle("i#eta");
234  poszp_[i]->GetYaxis()->SetTitle("i#phi");
235  sprintf(title, "%s-", dets[i].c_str());
236  poszn_[i] = tfile->make<TH2F>("poszn0", title, 100, 0, 100, 360, 0, 360);
237  poszn_[i]->GetXaxis()->SetTitle("i#eta");
238  poszn_[i]->GetYaxis()->SetTitle("i#phi");
239  } else {
240  sprintf(title, "%s+", dets[i].c_str());
241  poszp_[i] = tfile->make<TH2F>("poszp1", title, 100, -200, 200, 100, -200, 200);
242  poszp_[i]->GetXaxis()->SetTitle("x (cm)");
243  poszp_[i]->GetYaxis()->SetTitle("y (cm)");
244  sprintf(title, "%s-", dets[i].c_str());
245  poszn_[i] = tfile->make<TH2F>("poszn1", title, 100, -200, 200, 100, -200, 200);
246  poszn_[i]->GetXaxis()->SetTitle("x (cm)");
247  poszn_[i]->GetYaxis()->SetTitle("y (cm)");
248  }
249  poszp_[i]->GetYaxis()->SetTitleOffset(1.2);
250  poszp_[i]->Sumw2();
251  poszn_[i]->GetYaxis()->SetTitleOffset(1.2);
252  poszn_[i]->Sumw2();
253  }
254 }
255 
257 #ifdef EDM_ML_DEBUG
258  edm::LogVerbatim("HitStudy") << "Run = " << e.id().run() << " Event = " << e.id().event();
259 #endif
260  // get handles to calogeometry
261  geometry_ = &iS.getData(tokGeom_);
262 
263  double eInc = 0, etaInc = 0, phiInc = 0;
264  int type(-1);
266  e.getByToken(tok_evt_, EvtHandle);
267  if (EvtHandle.isValid()) {
268  const HepMC::GenEvent* myGenEvent = EvtHandle->GetEvent();
269 
270  HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
271  if (p != myGenEvent->particles_end()) {
272  eInc = (*p)->momentum().e();
273  etaInc = (*p)->momentum().eta();
274  phiInc = (*p)->momentum().phi();
275  }
276  double ptInc = eInc / std::cosh(etaInc);
277  ptInc_->Fill(ptInc);
278  eneInc_->Fill(eInc);
279  etaInc_->Fill(etaInc);
280  phiInc_->Fill(phiInc);
281 
282  if (std::abs(etaInc) < 1.46)
283  type = 0;
284  else if (std::abs(etaInc) > 1.49 && std::abs(etaInc) < 3.0)
285  type = 1;
286  }
287 
288  int typeMin = (type < 0) ? 0 : type;
289  int typeMax = (type < 0) ? 1 : type;
290  for (int type = typeMin; type <= typeMax; ++type) {
292  e.getByToken(toks_calo_[type], hitsCalo);
293  bool getHits = (hitsCalo.isValid());
294 #ifdef EDM_ML_DEBUG
295  edm::LogVerbatim("HitStudy") << "EcalSimHitStudy: Input flags Hits " << getHits << " with " << hitsCalo->size()
296  << " hits";
297 #endif
298  if (getHits) {
299  std::vector<PCaloHit> caloHits;
300  caloHits.insert(caloHits.end(), hitsCalo->begin(), hitsCalo->end());
301  if (!caloHits.empty())
302  analyzeHits(caloHits, type);
303  }
304  }
305 }
306 
307 void EcalSimHitStudy::analyzeHits(std::vector<PCaloHit>& hits, int indx) {
308  unsigned int nEC(0);
309  std::map<unsigned int, EcalHit> hitMap;
310  double etot(0), etotG(0);
311  for (auto hit : hits) {
312  double edep = hit.energy();
313  double time = hit.time();
314  unsigned int id_ = hit.id();
315  double edepEM = hit.energyEM();
316  double edepHad = hit.energyHad();
317  if (indx == 0) {
318  if ((hit.depth() == 1) || (hit.depth() == 2))
319  continue;
320  }
321  if (time <= tmax_) {
322  auto it = hitMap.find(id_);
323  if (it == hitMap.end()) {
324  uint16_t dep = hit.depth();
325  hitMap[id_] = EcalHit(dep, time, edep);
326  } else {
327  (it->second).energy += edep;
328  }
329  etotG += edep;
330  ++nEC;
331  }
332  time_[indx]->Fill(time);
333  edep_[indx]->Fill(edep);
334  edepEM_[indx]->Fill(edepEM);
335  edepHad_[indx]->Fill(edepHad);
336  etot += edep;
337  }
338 
339  double edepM(0);
340  unsigned int idM(0);
341  uint16_t depM(0);
342  for (auto it : hitMap) {
343  if (it.second.energy > edepM) {
344  idM = it.first;
345  edepM = it.second.energy;
346  depM = it.second.id;
347  }
348  }
349 
350  bool select(true);
351  if (selX_ >= 0) {
352  if ((depM & 0X4) != 0)
353  select = (selX_ > 0);
354  else
355  select = (selX_ == 0);
356  }
357 #ifdef EDM_ML_DEBUG
358  edm::LogVerbatim("HitStudy") << "EcalSimHitStudy::analyzeHits: Index " << indx << " Emax " << edepM << " IDMax "
359  << std::hex << idM << ":" << depM << std::dec << " Select " << select << ":" << selX_
360  << " Hits " << hits.size() << ":" << nEC << ":" << hitMap.size() << " ETotal " << etot
361  << ":" << etotG;
362 #endif
363  if (select) {
364  etot_[indx]->Fill(etot);
365  etotg_[indx]->Fill(etotG);
366  hit_[indx]->Fill(double(nEC));
367  for (auto it : hitMap) {
368  timeAll_[indx]->Fill((it.second).time);
369  edepAll_[indx]->Fill((it.second).energy);
370  DetId id(it.first);
371  if (indx == 0) {
372  if (EBDetId(id).zside() >= 0)
373  poszp_[indx]->Fill(EBDetId(id).ietaAbs(), EBDetId(id).iphi());
374  else
375  poszn_[indx]->Fill(EBDetId(id).ietaAbs(), EBDetId(id).iphi());
376  } else {
377  GlobalPoint gpos = geometry_->getGeometry(id)->getPosition();
378  if (EEDetId(id).zside() >= 0)
379  poszp_[indx]->Fill(gpos.x(), gpos.y());
380  else
381  poszn_[indx]->Fill(gpos.x(), gpos.y());
382  }
383  }
384 
385  math::XYZVector meanPosition(0.0, 0.0, 0.0);
386  std::vector<math::XYZVector> position;
387  std::vector<double> energy;
388  double e9(0), e25(0);
389  for (auto it : hitMap) {
390  DetId id(it.first);
391  int deta(99), dphi(99), dz(0);
392  if (indx == 0) {
393  deta = std::abs(EBDetId(id).ietaAbs() - EBDetId(idM).ietaAbs());
394  dphi = std::abs(EBDetId(id).iphi() - EBDetId(idM).iphi());
395  if (dphi > 180)
396  dphi = std::abs(dphi - 360);
397  dz = std::abs(EBDetId(id).zside() - EBDetId(idM).zside());
398  } else {
399  deta = std::abs(EEDetId(id).ix() - EEDetId(idM).ix());
400  dphi = std::abs(EEDetId(id).iy() - EEDetId(idM).iy());
401  dz = std::abs(EEDetId(id).zside() - EEDetId(idM).zside());
402  }
403  if (deta <= 1 && dphi <= 1 && dz < 1)
404  e9 += (it.second).energy;
405  if (deta <= 2 && dphi <= 2 && dz < 1) {
406  e25 += (it.second).energy;
407  GlobalPoint gpos = geometry_->getGeometry(id)->getPosition();
408  math::XYZVector pos(gpos.x(), gpos.y(), gpos.z());
409  meanPosition += (it.second).energy * pos;
410  position.push_back(pos);
411  energy.push_back((it.second).energy);
412  }
413  }
414  double r1by9 = (e9 > 0) ? (edepM / e9) : -1;
415  double r1by25 = (e25 > 0) ? (edepM / e25) : -1;
416  double r9by25 = (e25 > 0) ? (e9 / e25) : -1;
417 
418  meanPosition /= e25;
419  double denom(0), numEtaEta(0), numEtaPhi(0), numPhiPhi(0);
420  for (unsigned int k = 0; k < position.size(); ++k) {
421  double dEta = position[k].eta() - meanPosition.eta();
422  double dPhi = position[k].phi() - meanPosition.phi();
423  if (dPhi > +M_PI) {
424  dPhi = 2 * M_PI - dPhi;
425  }
426  if (dPhi < -M_PI) {
427  dPhi = 2 * M_PI + dPhi;
428  }
429 
430  double w = std::max(0.0, (w0_ + std::log(energy[k] / e25)));
431  denom += w;
432  numEtaEta += std::abs(w * dEta * dEta);
433  numEtaPhi += std::abs(w * dEta * dPhi);
434  numPhiPhi += std::abs(w * dPhi * dPhi);
435 #ifdef EDM_ML_DEBUG
436  edm::LogVerbatim("HitStudy") << "[" << k << "] dEta " << dEta << " dPhi " << dPhi << " Wt " << energy[k] / e25
437  << ":" << std::log(energy[k] / e25) << ":" << w;
438 #endif
439  }
440  double sEtaEta = (denom > 0) ? (numEtaEta / denom) : -1.0;
441  double sEtaPhi = (denom > 0) ? (numEtaPhi / denom) : -1.0;
442  double sPhiPhi = (denom > 0) ? (numPhiPhi / denom) : -1.0;
443 
444 #ifdef EDM_ML_DEBUG
445  edm::LogVerbatim("HitStudy") << "EcalSimHitStudy::Ratios " << r1by9 << " : " << r1by25 << " : " << r9by25
446  << " Covariances " << sEtaEta << " : " << sEtaPhi << " : " << sPhiPhi;
447 #endif
448  r1by9_[indx]->Fill(r1by9);
449  r1by25_[indx]->Fill(r1by25);
450  r9by25_[indx]->Fill(r9by25);
451  sEtaEta_[indx]->Fill(sEtaEta);
452  sEtaPhi_[indx]->Fill(sEtaPhi);
453  sPhiPhi_[indx]->Fill(sPhiPhi);
454  }
455 }
456 
457 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
Log< level::Info, true > LogVerbatim
TH1F * etotg_[ndets_]
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
TH1F * hit_[ndets_]
TH1F * time_[ndets_]
T w() const
T z() const
Definition: PV3DBase.h:61
const CaloGeometry * geometry_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
TH1F * r1by25_[ndets_]
const std::string g4Label_
int zside(DetId const &)
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
const edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
TH1F * sPhiPhi_[ndets_]
EcalSimHitStudy(const edm::ParameterSet &ps)
TH1F * edepEM_[ndets_]
void beginRun(edm::Run const &, edm::EventSetup const &) override
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
void analyze(edm::Event const &, edm::EventSetup const &) override
TH1F * edepHad_[ndets_]
Definition: tfile.py:1
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id.
Definition: CaloGeometry.cc:60
TH2F * poszp_[ndets_]
void beginJob() override
const std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > toks_calo_
void endRun(edm::Run const &, edm::EventSetup const &) override
const edm::ESGetToken< CaloGeometry, CaloGeometryRecord > tokGeom_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
void analyzeHits(std::vector< PCaloHit > &, int)
TH1F * sEtaPhi_[ndets_]
static const int ndets_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
#define M_PI
unsigned int id
Definition: DetId.h:17
TH1F * edep_[ndets_]
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
TH1F * sEtaEta_[ndets_]
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TH1F * edepAll_[ndets_]
bool isValid() const
Definition: HandleBase.h:70
HLT enums.
TH1F * timeAll_[ndets_]
const std::vector< std::string > hitLab_
static int position[264][3]
Definition: ReadPGInfo.cc:289
EcalHit(uint16_t i=0, double t=0, double e=0)
TH2F * poszn_[ndets_]
const double maxEnergy_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~EcalSimHitStudy() override
TH1F * r9by25_[ndets_]
TH1F * r1by9_[ndets_]
TH1F * etot_[ndets_]
const double tmax_
Definition: Run.h:45
const double w0_