CMS 3D CMS Logo

HGCalSimHitsClient.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <fstream>
3 #include <unistd.h>
4 #include <vector>
5 
8 
11 
19 
22 
24 private:
25  //member data
27  const int nTimes_, verbosity_;
29  unsigned int layers_;
30 
31 public:
32  explicit HGCalSimHitsClient(const edm::ParameterSet &);
33  ~HGCalSimHitsClient() override = default;
34 
35  void beginRun(const edm::Run &run, const edm::EventSetup &c) override;
36  void dqmEndJob(DQMStore::IBooker &ib, DQMStore::IGetter &ig) override;
38  int simHitsEndjob(const std::vector<MonitorElement *> &hgcalMEs);
39 };
40 
42  : nameDetector_(iConfig.getParameter<std::string>("DetectorName")),
43  nTimes_(iConfig.getParameter<int>("TimeSlices")),
44  verbosity_(iConfig.getUntrackedParameter<int>("Verbosity", 0)),
46  edm::ESInputTag{"", nameDetector_})) {}
47 
49  const HGCalDDDConstants *hgcons = &iSetup.getData(tok_hgcal_);
50  layers_ = hgcons->layers(true);
51  if (verbosity_ > 0)
52  edm::LogVerbatim("HGCalValidation") << "Initialize HGCalSimHitsClient for " << nameDetector_ << " : " << layers_;
53 }
54 
56 
58  ig.setCurrentFolder("/");
59  if (verbosity_ > 0)
60  edm::LogVerbatim("HGCalValidation") << " runClient";
61  std::vector<MonitorElement *> hgcalMEs;
62  std::vector<std::string> fullDirPath = ig.getSubdirs();
63 
64  for (unsigned int i = 0; i < fullDirPath.size(); i++) {
65  if (verbosity_ > 0)
66  edm::LogVerbatim("HGCalValidation") << "fullPath: " << fullDirPath.at(i);
67  ig.setCurrentFolder(fullDirPath.at(i));
68  std::vector<std::string> fullSubDirPath = ig.getSubdirs();
69 
70  for (unsigned int j = 0; j < fullSubDirPath.size(); j++) {
71  if (verbosity_ > 0)
72  edm::LogVerbatim("HGCalValidation") << "fullSubPath: " << fullSubDirPath.at(j);
73  std::string nameDirectory = "HGCAL/HGCalSimHitsV/" + nameDetector_;
74 
75  if (strcmp(fullSubDirPath.at(j).c_str(), nameDirectory.c_str()) == 0) {
76  hgcalMEs = ig.getContents(fullSubDirPath.at(j));
77  if (verbosity_ > 0)
78  edm::LogVerbatim("HGCalValidation") << "hgcalMES size : " << hgcalMEs.size();
79  if (!simHitsEndjob(hgcalMEs))
80  edm::LogWarning("HGCalValidation") << "\nError in SimhitsEndjob!";
81  }
82  }
83  }
84 }
85 
86 int HGCalSimHitsClient::simHitsEndjob(const std::vector<MonitorElement *> &hgcalMEs) {
87  std::vector<MonitorElement *> energy_[6];
88  std::vector<MonitorElement *> EtaPhi_Plus_, EtaPhi_Minus_;
89  std::vector<MonitorElement *> HitOccupancy_Plus_, HitOccupancy_Minus_;
90  MonitorElement *MeanHitOccupancy_Plus_, *MeanHitOccupancy_Minus_;
91 
92  std::ostringstream name;
93  double nevent;
94  int nbinsx, nbinsy;
95  for (unsigned int ilayer = 0; ilayer < layers_; ilayer++) {
96  for (int itimeslice = 0; itimeslice < nTimes_; itimeslice++) {
97  //Energy
98  name.str("");
99  name << "energy_time_" << itimeslice << "_layer_" << ilayer;
100  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
101  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
102  energy_[itimeslice].push_back(hgcalMEs[ih]);
103  }
104  }
105  //normalization
106  nevent = energy_[itimeslice].at(ilayer)->getEntries();
107  nbinsx = energy_[itimeslice].at(ilayer)->getNbinsX();
108  for (int i = 1; i <= nbinsx; i++) {
109  double binValue = energy_[itimeslice].at(ilayer)->getBinContent(i) / nevent;
110  energy_[itimeslice].at(ilayer)->setBinContent(i, binValue);
111  }
112  }
113 
114  //EtaPhi 2d plots
115  name.str("");
116  name << "EtaPhi_Plus_"
117  << "layer_" << ilayer;
118  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
119  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
120  EtaPhi_Plus_.push_back(hgcalMEs[ih]);
121  }
122  }
123 
124  name.str("");
125  name << "EtaPhi_Minus_"
126  << "layer_" << ilayer;
127  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
128  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
129  EtaPhi_Minus_.push_back(hgcalMEs[ih]);
130  }
131  }
132  //normalization EtaPhi
133  nevent = EtaPhi_Plus_.at(ilayer)->getEntries();
134  nbinsx = EtaPhi_Plus_.at(ilayer)->getNbinsX();
135  nbinsy = EtaPhi_Plus_.at(ilayer)->getNbinsY();
136  for (int i = 1; i <= nbinsx; ++i) {
137  for (int j = 1; j <= nbinsy; ++j) {
138  double binValue = EtaPhi_Plus_.at(ilayer)->getBinContent(i, j) / nevent;
139  EtaPhi_Plus_.at(ilayer)->setBinContent(i, j, binValue);
140  }
141  }
142 
143  nevent = EtaPhi_Minus_.at(ilayer)->getEntries();
144  nbinsx = EtaPhi_Minus_.at(ilayer)->getNbinsX();
145  nbinsy = EtaPhi_Plus_.at(ilayer)->getNbinsY();
146  for (int i = 1; i <= nbinsx; ++i) {
147  for (int j = 1; j <= nbinsy; ++j) {
148  double binValue = EtaPhi_Minus_.at(ilayer)->getBinContent(i, j) / nevent;
149  EtaPhi_Minus_.at(ilayer)->setBinContent(i, j, binValue);
150  }
151  }
152 
153  //HitOccupancy
154  name.str("");
155  name << "HitOccupancy_Plus_layer_" << ilayer;
156  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
157  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
158  HitOccupancy_Plus_.push_back(hgcalMEs[ih]);
159  }
160  }
161  name.str("");
162  name << "HitOccupancy_Minus_layer_" << ilayer;
163  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
164  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
165  HitOccupancy_Minus_.push_back(hgcalMEs[ih]);
166  }
167  }
168 
169  nevent = HitOccupancy_Plus_.at(ilayer)->getEntries();
170  nbinsx = HitOccupancy_Plus_.at(ilayer)->getNbinsX();
171  for (int i = 1; i <= nbinsx; ++i) {
172  double binValue = HitOccupancy_Plus_.at(ilayer)->getBinContent(i) / nevent;
173  HitOccupancy_Plus_.at(ilayer)->setBinContent(i, binValue);
174  }
175 
176  nevent = HitOccupancy_Minus_.at(ilayer)->getEntries();
177  nbinsx = HitOccupancy_Minus_.at(ilayer)->getNbinsX();
178  for (int i = 1; i <= nbinsx; ++i) {
179  double binValue = HitOccupancy_Minus_.at(ilayer)->getBinContent(i) / nevent;
180  HitOccupancy_Minus_.at(ilayer)->setBinContent(i, binValue);
181  }
182 
183  } //loop over layers
184 
185  name.str("");
186  name << "MeanHitOccupancy_Plus";
187  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
188  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
189  MeanHitOccupancy_Plus_ = hgcalMEs[ih];
190  for (int ilayer = 0; ilayer < static_cast<int>(layers_); ++ilayer) {
191  double meanVal = HitOccupancy_Plus_.at(ilayer)->getMean();
192  MeanHitOccupancy_Plus_->setBinContent(ilayer + 1, meanVal);
193  }
194  break;
195  }
196  }
197 
198  name.str("");
199  name << "MeanHitOccupancy_Minus";
200  for (unsigned int ih = 0; ih < hgcalMEs.size(); ih++) {
201  if (strcmp(hgcalMEs[ih]->getName().c_str(), name.str().c_str()) == 0) {
202  MeanHitOccupancy_Minus_ = hgcalMEs[ih];
203  for (int ilayer = 0; ilayer < static_cast<int>(layers_); ++ilayer) {
204  double meanVal = HitOccupancy_Minus_.at(ilayer)->getMean();
205  MeanHitOccupancy_Minus_->setBinContent(ilayer + 1, meanVal);
206  }
207  break;
208  }
209  }
210 
211  return 1;
212 }
213 
215 
Log< level::Info, true > LogVerbatim
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
virtual void runClient_(DQMStore::IBooker &ib, DQMStore::IGetter &ig)
int simHitsEndjob(const std::vector< MonitorElement *> &hgcalMEs)
void dqmEndJob(DQMStore::IBooker &ib, DQMStore::IGetter &ig) override
void beginRun(const edm::Run &run, const edm::EventSetup &c) override
int nevent
Definition: AMPTWrapper.h:84
Transition
Definition: Transition.h:12
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int layers(bool reco) const
bool getData(T &iHolder) const
Definition: EventSetup.h:122
HGCalSimHitsClient(const edm::ParameterSet &)
const edm::ESGetToken< HGCalDDDConstants, IdealGeometryRecord > tok_hgcal_
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
virtual double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
HLT enums.
std::string getName(const G4String &)
Definition: ForwardName.cc:3
~HGCalSimHitsClient() override=default
const std::string nameDetector_
Log< level::Warning, false > LogWarning
Definition: Run.h:45
ib
Definition: cuy.py:661
virtual std::vector< dqm::harvesting::MonitorElement * > getContents(std::string const &path) const
Definition: DQMStore.cc:610
virtual DQM_DEPRECATED std::vector< std::string > getSubdirs() const
Definition: DQMStore.cc:717