CMS 3D CMS Logo

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