CMS 3D CMS Logo

HGCalSimHitStudy.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 
10 
14 
25 
28 
31 
32 #include "CLHEP/Units/GlobalSystemOfUnits.h"
33 #include "CLHEP/Units/GlobalPhysicalConstants.h"
34 
35 #include "TH1D.h"
36 #include "TH2D.h"
37 
38 class HGCalSimHitStudy : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
39 public:
40  struct hitsinfo {
42  phi = eta = energy = time = 0.0;
43  layer = 0;
44  }
45  double phi, eta, energy, time;
46  int layer;
47  };
48 
49  explicit HGCalSimHitStudy(const edm::ParameterSet&);
50  ~HGCalSimHitStudy() override = default;
51  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
52 
53 protected:
54  void beginJob() override;
55  void beginRun(edm::Run const&, edm::EventSetup const&) override;
56  void analyze(edm::Event const&, edm::EventSetup const&) override;
57  void endRun(edm::Run const&, edm::EventSetup const&) override {}
58 
59 private:
60  void analyzeHits(int, const std::string&, const std::vector<PCaloHit>&);
61 
62  // ----------member data ---------------------------
63  const std::vector<std::string> nameDetectors_, caloHitSources_;
64  const double rmin_, rmax_, zmin_, zmax_;
65  const double etamin_, etamax_;
67  const bool ifNose_, ifLayer_;
68  const std::vector<edm::ESGetToken<HGCalDDDConstants, IdealGeometryRecord> > tok_hgcGeom_;
69  const std::vector<edm::EDGetTokenT<edm::PCaloHitContainer> > tok_hits_;
70  std::vector<const HGCalDDDConstants*> hgcons_;
71  std::vector<int> layers_, layerFront_;
72 
73  //histogram related stuff
74  std::vector<TH2D*> h_RZ_, h_EtaPhi_, h_EtFiZp_, h_EtFiZm_, h_XY_;
75  std::vector<TH1D*> h_E_, h_T_, h_LayerZp_, h_LayerZm_;
76  std::vector<TH1D*> h_W1_, h_W2_, h_C1_, h_C2_, h_Ly_;
77 };
78 
80  : nameDetectors_(iConfig.getParameter<std::vector<std::string> >("detectorNames")),
81  caloHitSources_(iConfig.getParameter<std::vector<std::string> >("caloHitSources")),
82  rmin_(iConfig.getUntrackedParameter<double>("rMin", 0.0)),
83  rmax_(iConfig.getUntrackedParameter<double>("rMax", 3000.0)),
84  zmin_(iConfig.getUntrackedParameter<double>("zMin", 3000.0)),
85  zmax_(iConfig.getUntrackedParameter<double>("zMax", 6000.0)),
86  etamin_(iConfig.getUntrackedParameter<double>("etaMin", 1.0)),
87  etamax_(iConfig.getUntrackedParameter<double>("etaMax", 3.0)),
88  nbinR_(iConfig.getUntrackedParameter<int>("nBinR", 300)),
89  nbinZ_(iConfig.getUntrackedParameter<int>("nBinZ", 300)),
90  nbinEta_(iConfig.getUntrackedParameter<int>("nBinEta", 200)),
91  nLayers_(iConfig.getUntrackedParameter<int>("layers", 50)),
92  verbosity_(iConfig.getUntrackedParameter<int>("verbosity", 0)),
93  ifNose_(iConfig.getUntrackedParameter<bool>("ifNose", false)),
94  ifLayer_(iConfig.getUntrackedParameter<bool>("ifLayer", false)),
95  tok_hgcGeom_{
97  [this](const std::string& name) {
98  return esConsumes<HGCalDDDConstants, IdealGeometryRecord, edm::Transition::BeginRun>(
99  edm::ESInputTag{"", name});
100  })},
101  tok_hits_{edm::vector_transform(caloHitSources_, [this](const std::string& source) {
102  return consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits", source));
103  })} {
104  usesResource(TFileService::kSharedResource);
105 }
106 
109  std::vector<std::string> names = {"HGCalEESensitive", "HGCalHESiliconSensitive", "HGCalHEScintillatorSensitive"};
110  std::vector<std::string> sources = {"HGCHitsEE", "HGCHitsHEfront", "HGCHitsHEback"};
111  desc.add<std::vector<std::string> >("detectorNames", names);
112  desc.add<std::vector<std::string> >("caloHitSources", sources);
113  desc.addUntracked<double>("rMin", 0.0);
114  desc.addUntracked<double>("rMax", 3000.0);
115  desc.addUntracked<double>("zMin", 3000.0);
116  desc.addUntracked<double>("zMax", 6000.0);
117  desc.addUntracked<double>("etaMin", 1.0);
118  desc.addUntracked<double>("etaMax", 3.0);
119  desc.addUntracked<int>("nBinR", 300);
120  desc.addUntracked<int>("nBinZ", 300);
121  desc.addUntracked<int>("nBinEta", 200);
122  desc.addUntracked<int>("layers", 50);
123  desc.addUntracked<int>("verbosity", 0);
124  desc.addUntracked<bool>("ifNose", false);
125  desc.addUntracked<bool>("ifLayer", false);
126  descriptions.add("hgcalSimHitStudy", desc);
127 }
128 
130  //Now the hits
131  for (unsigned int k = 0; k < tok_hits_.size(); ++k) {
132  const edm::Handle<edm::PCaloHitContainer>& theCaloHitContainers = iEvent.getHandle(tok_hits_[k]);
133  if (theCaloHitContainers.isValid()) {
134  if (verbosity_ > 0)
135  edm::LogVerbatim("HGCalValidation") << " PcalohitItr = " << theCaloHitContainers->size();
136  std::vector<PCaloHit> caloHits;
137  caloHits.insert(caloHits.end(), theCaloHitContainers->begin(), theCaloHitContainers->end());
138  analyzeHits(k, nameDetectors_[k], caloHits);
139  } else if (verbosity_ > 0) {
140  edm::LogVerbatim("HGCalValidation") << "PCaloHitContainer does not "
141  << "exist for " << nameDetectors_[k];
142  }
143  }
144 }
145 
146 void HGCalSimHitStudy::analyzeHits(int ih, std::string const& name, std::vector<PCaloHit> const& hits) {
147  if (verbosity_ > 0)
148  edm::LogVerbatim("HGCalValidation") << name << " with " << hits.size() << " PcaloHit elements";
149 
150  std::map<uint32_t, hitsinfo> map_hits;
151  map_hits.clear();
152 
153  unsigned int nused(0);
154  for (auto const& hit : hits) {
155  double energy = hit.energy();
156  double time = hit.time();
157  uint32_t id = hit.id();
158  int cell, sector, sector2(0), layer, zside;
159  int subdet(0), cell2(0), type(0);
160  HepGeom::Point3D<float> gcoord;
161  std::pair<float, float> xy;
162  if (ifNose_) {
163  HFNoseDetId detId = HFNoseDetId(id);
164  subdet = detId.subdetId();
165  cell = detId.cellU();
166  cell2 = detId.cellV();
167  sector = detId.waferU();
168  sector2 = detId.waferV();
169  type = detId.type();
170  layer = detId.layer();
171  zside = detId.zside();
172  xy = hgcons_[ih]->locateCell(zside, layer, sector, sector2, cell, cell2, false, true, false, false);
173  h_W2_[ih]->Fill(sector2);
174  h_C2_[ih]->Fill(cell2);
175  } else if (hgcons_[ih]->waferHexagon8()) {
176  HGCSiliconDetId detId = HGCSiliconDetId(id);
177  subdet = static_cast<int>(detId.det());
178  cell = detId.cellU();
179  cell2 = detId.cellV();
180  sector = detId.waferU();
181  sector2 = detId.waferV();
182  type = detId.type();
183  layer = detId.layer();
184  zside = detId.zside();
185  xy = hgcons_[ih]->locateCell(zside, layer, sector, sector2, cell, cell2, false, true, false, false);
186  h_W2_[ih]->Fill(sector2);
187  h_C2_[ih]->Fill(cell2);
188  } else if (hgcons_[ih]->tileTrapezoid()) {
190  subdet = static_cast<int>(detId.det());
191  sector = detId.ieta();
192  cell = detId.iphi();
193  type = detId.type();
194  layer = detId.layer();
195  zside = detId.zside();
196  xy = hgcons_[ih]->locateCellTrap(zside, layer, sector, cell, false, false);
197  } else {
198  edm::LogError("HGCalValidation") << "HGCalSimHitStudy: Wrong geometry mode " << hgcons_[ih]->geomMode();
199  continue;
200  }
201  double zp = hgcons_[ih]->waferZ(layer, false);
202  if (zside < 0)
203  zp = -zp;
204  double xp = (zp < 0) ? -xy.first : xy.first;
205  gcoord = HepGeom::Point3D<float>(xp, xy.second, zp);
206  if (verbosity_ > 2)
207  edm::LogVerbatim("HGCalValidation")
208  << "i/p " << subdet << ":" << zside << ":" << layer << ":" << sector << ":" << sector2 << ":" << cell << ":"
209  << cell2 << " o/p " << xy.first << ":" << xy.second << ":" << zp;
210  nused++;
211  double tof = (gcoord.mag() * CLHEP::mm) / CLHEP::c_light;
212  if (verbosity_ > 1)
213  edm::LogVerbatim("HGCalValidation")
214  << "Detector " << name << " zside = " << zside << " layer = " << layer << " type = " << type
215  << " wafer = " << sector << ":" << sector2 << " cell = " << cell << ":" << cell2 << " positon = " << gcoord
216  << " energy = " << energy << " time = " << time << ":" << tof;
217  time -= tof;
218  if (time < 0)
219  time = 0;
220  hitsinfo hinfo;
221  if (map_hits.count(id) != 0) {
222  hinfo = map_hits[id];
223  } else {
224  hinfo.layer = layer + layerFront_[ih];
225  hinfo.phi = gcoord.getPhi();
226  hinfo.eta = gcoord.getEta();
227  hinfo.time = time;
228  }
229  hinfo.energy += energy;
230  map_hits[id] = hinfo;
231 
232  //Fill in histograms
233  h_RZ_[0]->Fill(std::abs(gcoord.z()), gcoord.rho());
234  h_RZ_[ih + 1]->Fill(std::abs(gcoord.z()), gcoord.rho());
235  if (ifLayer_) {
236  if (hinfo.layer <= static_cast<int>(h_XY_.size()))
237  h_XY_[hinfo.layer - 1]->Fill(gcoord.x(), gcoord.y());
238  } else {
239  h_EtaPhi_[0]->Fill(std::abs(hinfo.eta), hinfo.phi);
240  h_EtaPhi_[ih + 1]->Fill(std::abs(hinfo.eta), hinfo.phi);
241  }
242  h_Ly_[ih]->Fill(layer);
243  h_W1_[ih]->Fill(sector);
244  h_C1_[ih]->Fill(cell);
245  }
246  if (verbosity_ > 0)
247  edm::LogVerbatim("HGCalValidation") << name << " with " << map_hits.size() << ":" << nused << " detector elements"
248  << " being hit";
249 
250  for (auto const& hit : map_hits) {
251  hitsinfo hinfo = hit.second;
252  if (verbosity_ > 1)
253  edm::LogVerbatim("HGCalValidation")
254  << " ---------------------- eta = " << hinfo.eta << " phi = " << hinfo.phi << " layer = " << hinfo.layer
255  << " E = " << hinfo.energy << " T = " << hinfo.time;
256  h_E_[0]->Fill(hinfo.energy);
257  h_E_[ih + 1]->Fill(hinfo.energy);
258  h_T_[0]->Fill(hinfo.time);
259  h_T_[ih + 1]->Fill(hinfo.time);
260  if (hinfo.eta > 0) {
261  if (!ifLayer_) {
262  h_EtFiZp_[0]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
263  h_EtFiZp_[ih + 1]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
264  }
265  h_LayerZp_[0]->Fill(hinfo.layer, hinfo.energy);
266  h_LayerZp_[ih + 1]->Fill(hinfo.layer, hinfo.energy);
267  } else {
268  if (!ifLayer_) {
269  h_EtFiZm_[0]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
270  h_EtFiZm_[ih + 1]->Fill(std::abs(hinfo.eta), hinfo.phi, hinfo.energy);
271  }
272  h_LayerZm_[0]->Fill(hinfo.layer, hinfo.energy);
273  h_LayerZm_[ih + 1]->Fill(hinfo.layer, hinfo.energy);
274  }
275  }
276 }
277 
278 // ------------ method called when starting to processes a run ------------
280  for (unsigned int k = 0; k < nameDetectors_.size(); ++k) {
281  hgcons_.emplace_back(&iSetup.getData(tok_hgcGeom_[k]));
282  layers_.emplace_back(hgcons_.back()->layers(false));
283  layerFront_.emplace_back(hgcons_.back()->firstLayer());
284  if (verbosity_ > 0)
285  edm::LogVerbatim("HGCalValidation") << nameDetectors_[k] << " defined with " << layers_.back() << " Layers with "
286  << layerFront_.back() << " in front";
287  }
288 }
289 
292 
293  std::ostringstream name, title;
294  for (unsigned int ih = 0; ih <= nameDetectors_.size(); ++ih) {
295  name.str("");
296  title.str("");
297  if (ih == 0) {
298  name << "RZ_AllDetectors";
299  title << "R vs Z for All Detectors";
300  } else {
301  name << "RZ_" << nameDetectors_[ih - 1];
302  title << "R vs Z for " << nameDetectors_[ih - 1];
303  }
304  h_RZ_.emplace_back(
305  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinZ_, zmin_, zmax_, nbinR_, rmin_, rmax_));
306  if (ifLayer_) {
307  if (ih == 0) {
308  for (int ly = 0; ly < nLayers_; ++ly) {
309  name.str("");
310  title.str("");
311  name << "XY_L" << (ly + 1);
312  title << "Y vs X at Layer " << (ly + 1);
313  h_XY_.emplace_back(
314  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinR_, -rmax_, rmax_, nbinR_, -rmax_, rmax_));
315  }
316  }
317  } else {
318  name.str("");
319  title.str("");
320  if (ih == 0) {
321  name << "EtaPhi_AllDetectors";
322  title << "#phi vs #eta for All Detectors";
323  } else {
324  name << "EtaPhi_" << nameDetectors_[ih - 1];
325  title << "#phi vs #eta for " << nameDetectors_[ih - 1];
326  }
327  h_EtaPhi_.emplace_back(
328  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinEta_, etamin_, etamax_, 200, -M_PI, M_PI));
329  name.str("");
330  title.str("");
331  if (ih == 0) {
332  name << "EtFiZp_AllDetectors";
333  title << "#phi vs #eta (+z) for All Detectors";
334  } else {
335  name << "EtFiZp_" << nameDetectors_[ih - 1];
336  title << "#phi vs #eta (+z) for " << nameDetectors_[ih - 1];
337  }
338  h_EtFiZp_.emplace_back(
339  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinEta_, etamin_, etamax_, 200, -M_PI, M_PI));
340  name.str("");
341  title.str("");
342  if (ih == 0) {
343  name << "EtFiZm_AllDetectors";
344  title << "#phi vs #eta (-z) for All Detectors";
345  } else {
346  name << "EtFiZm_" << nameDetectors_[ih - 1];
347  title << "#phi vs #eta (-z) for " << nameDetectors_[ih - 1];
348  }
349  h_EtFiZm_.emplace_back(
350  fs->make<TH2D>(name.str().c_str(), title.str().c_str(), nbinEta_, etamin_, etamax_, 200, -M_PI, M_PI));
351  }
352  name.str("");
353  title.str("");
354  if (ih == 0) {
355  name << "LayerZp_AllDetectors";
356  title << "Energy vs Layer (+z) for All Detectors";
357  } else {
358  name << "LayerZp_" << nameDetectors_[ih - 1];
359  title << "Energy vs Layer (+z) for " << nameDetectors_[ih - 1];
360  }
361  h_LayerZp_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 60, 0.0, 60.0));
362  name.str("");
363  title.str("");
364  if (ih == 0) {
365  name << "LayerZm_AllDetectors";
366  title << "Energy vs Layer (-z) for All Detectors";
367  } else {
368  name << "LayerZm_" << nameDetectors_[ih - 1];
369  title << "Energy vs Layer (-z) for " << nameDetectors_[ih - 1];
370  }
371  h_LayerZm_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 60, 0.0, 60.0));
372 
373  name.str("");
374  title.str("");
375  if (ih == 0) {
376  name << "E_AllDetectors";
377  title << "Energy Deposit for All Detectors";
378  } else {
379  name << "E_" << nameDetectors_[ih - 1];
380  title << "Energy Deposit for " << nameDetectors_[ih - 1];
381  }
382  h_E_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 1000, 0.0, 1.0));
383 
384  name.str("");
385  title.str("");
386  if (ih == 0) {
387  name << "T_AllDetectors";
388  title << "Time for All Detectors";
389  } else {
390  name << "T_" << nameDetectors_[ih - 1];
391  title << "Time for " << nameDetectors_[ih - 1];
392  }
393  h_T_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 1000, 0.0, 200.0));
394  }
395 
396  for (unsigned int ih = 0; ih < nameDetectors_.size(); ++ih) {
397  name.str("");
398  title.str("");
399  name << "LY_" << nameDetectors_[ih];
400  title << "Layer number for " << nameDetectors_[ih];
401  h_Ly_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 200, 0, 100));
402  if (nameDetectors_[ih] == "HGCalHEScintillatorSensitive") {
403  name.str("");
404  title.str("");
405  name << "IR_" << nameDetectors_[ih];
406  title << "Radius index for " << nameDetectors_[ih];
407  h_W1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 200, -50, 50));
408  name.str("");
409  title.str("");
410  name << "FI_" << nameDetectors_[ih];
411  title << "#phi index for " << nameDetectors_[ih];
412  h_C1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 720, 0, 360));
413  } else {
414  name.str("");
415  title.str("");
416  name << "WU_" << nameDetectors_[ih];
417  title << "u index of wafers for " << nameDetectors_[ih];
418  h_W1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 200, -50, 50));
419  name.str("");
420  title.str("");
421  name << "WV_" << nameDetectors_[ih];
422  title << "v index of wafers for " << nameDetectors_[ih];
423  h_W2_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 100, -50, 50));
424  name.str("");
425  title.str("");
426  name << "CU_" << nameDetectors_[ih];
427  title << "u index of cells for " << nameDetectors_[ih];
428  h_C1_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 100, 0, 50));
429  name.str("");
430  title.str("");
431  name << "CV_" << nameDetectors_[ih];
432  title << "v index of cells for " << nameDetectors_[ih];
433  h_C2_.emplace_back(fs->make<TH1D>(name.str().c_str(), title.str().c_str(), 100, 0, 50));
434  }
435  }
436 }
437 
439 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
void beginJob() override
const std::vector< std::string > nameDetectors_
Log< level::Info, true > LogVerbatim
std::vector< TH2D * > h_XY_
void analyze(edm::Event const &, edm::EventSetup const &) override
std::vector< TH1D * > h_LayerZm_
std::vector< TH2D * > h_RZ_
std::vector< TH1D * > h_T_
std::vector< TH1D * > h_W1_
void beginRun(edm::Run const &, edm::EventSetup const &) override
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
int layer() const
get the layer #
Definition: HFNoseDetId.h:57
std::vector< TH2D * > h_EtFiZm_
std::vector< int > layerFront_
std::vector< TH1D * > h_C1_
int waferU() const
Definition: HFNoseDetId.h:76
int cellU() const
get the cell #&#39;s in u,v or in x,y
const std::vector< edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > > tok_hgcGeom_
int zside(DetId const &)
int type() const
get/set the type
Log< level::Error, false > LogError
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
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
const std::string names[nVars_]
HGCalSimHitStudy(const edm::ParameterSet &)
const double etamin_
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
const double zmin_
int layer() const
get the layer #
int iEvent
Definition: GenABIO.cc:224
int iphi() const
get the phi index
~HGCalSimHitStudy() override=default
std::vector< TH1D * > h_LayerZp_
void analyzeHits(int, const std::string &, const std::vector< PCaloHit > &)
int layer() const
get the layer #
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< TH2D * > h_EtFiZp_
int waferU() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int cellU() const
get the cell #&#39;s in u,v or in x,y
Definition: HFNoseDetId.h:60
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< TH1D * > h_C2_
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
int zside() const
get the z-side of the cell (1/-1)
const double rmin_
int type() const
get the type
Definition: HFNoseDetId.h:51
#define M_PI
std::vector< int > layers_
unsigned int id
const double etamax_
int waferV() const
const double zmax_
const double rmax_
int cellV() const
Definition: HFNoseDetId.h:61
int zside() const
get the z-side of the cell (1/-1)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
int waferV() const
Definition: HFNoseDetId.h:79
int cellV() const
bool isValid() const
Definition: HandleBase.h:70
const std::vector< edm::EDGetTokenT< edm::PCaloHitContainer > > tok_hits_
std::vector< TH1D * > h_Ly_
std::vector< TH1D * > h_E_
int type() const
get the type
std::vector< const HGCalDDDConstants * > hgcons_
std::vector< TH1D * > h_W2_
static std::string const source
Definition: EdmProvDump.cc:49
const std::vector< std::string > caloHitSources_
std::vector< TH2D * > h_EtaPhi_
void endRun(edm::Run const &, edm::EventSetup const &) override
Definition: Run.h:45
int zside() const
get the z-side of the cell (1/-1)
Definition: HFNoseDetId.h:54