CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalDigisClient.cc
Go to the documentation of this file.
2 // -*- C++ -*-
3 //
4 // Package: HcalDigisClient
5 // Class: HcalDigisClient
6 //
14 //
15 // Original Author: Ali Fahim,22 R-013,+41227672649,
16 // Created: Wed Mar 23 22:54:28 CET 2011
17 //
18 //
19 
20 
21 // system include files
22 
24  outputFile_ = iConfig.getUntrackedParameter<std::string > ("outputFile", "HcalDigisClient.root");
25  dirName_ = iConfig.getParameter<std::string > ("DQMDirName");
26  msm_ = new std::map<std::string, MonitorElement*>();
27 
28 
29 }
30 
31 
33  delete msm_;
34 }
35 
36 
38 
39  std::string strtmp;
40  HistLim ietaLim(85, -42.5, 42.5);
41 
42  for (int depth = 1; depth <= 4; depth++) {
43  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth" + str(depth) + "_" + subdetopt;
44  book1D(ib,strtmp, ietaLim);
45  }
46 
47 }
48 
51  std::vector<MonitorElement*> hcalMEs;
52  // Since out folders are fixed to three, we can just go over these three folders
53  // i.e., CaloTowersV/CaloTowersTask, HcalRecHitsV/HcalRecHitTask, NoiseRatesV/NoiseRatesTask.
54  std::vector<std::string> fullPathHLTFolders = ig.getSubdirs();
55  for (unsigned int i = 0; i < fullPathHLTFolders.size(); i++) {
56  ig.setCurrentFolder(fullPathHLTFolders[i]);
57  std::vector<std::string> fullSubPathHLTFolders = ig.getSubdirs();
58  for (unsigned int j = 0; j < fullSubPathHLTFolders.size(); j++) {
59  if (strcmp(fullSubPathHLTFolders[j].c_str(), "HcalDigisV/HcalDigiTask") == 0) {
60  hcalMEs = ig.getContents(fullSubPathHLTFolders[j]);
61  if (!HcalDigisEndjob(hcalMEs, "HB"))
62  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HB";
63  if (!HcalDigisEndjob(hcalMEs, "HE"))
64  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HE";
65  if (!HcalDigisEndjob(hcalMEs, "HO"))
66  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HO";
67  if (!HcalDigisEndjob(hcalMEs, "HF"))
68  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HF"; }
69  }
70  }
71 }
72 
73 int HcalDigisClient::HcalDigisEndjob(const std::vector<MonitorElement*> &hcalMEs, std::string subdet_) {
74 
75  using namespace std;
76  string strtmp;
77 
78 
79  MonitorElement * nevtot(0);
80  MonitorElement * ieta_iphi_occupancy_map1(0);
81  MonitorElement * ieta_iphi_occupancy_map2(0);
82  MonitorElement * ieta_iphi_occupancy_map3(0);
83  MonitorElement * ieta_iphi_occupancy_map4(0);
84 
85 
86  // std::cout << " Number of histos " << hcalMEs.size() << std::endl;
87 
88  for (unsigned int ih = 0; ih < hcalMEs.size(); ih++) {
89 
90  if (hcalMEs[ih]->getName() == "nevtot") nevtot = hcalMEs[ih];
91 
92  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth1_" + subdet_;
93  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map1 = hcalMEs[ih];
94  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth2_" + subdet_;
95  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map2 = hcalMEs[ih];
96  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth3_" + subdet_;
97  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map3 = hcalMEs[ih];
98  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth4_" + subdet_;
99  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map4 = hcalMEs[ih];
100 
101  }
102 
103  if (nevtot == 0 ||
104  ieta_iphi_occupancy_map1 == 0 ||
105  ieta_iphi_occupancy_map2 == 0 ||
106  ieta_iphi_occupancy_map3 == 0 ||
107  ieta_iphi_occupancy_map4 == 0
108  ) {
109  edm::LogError("HcalDigisClient") << "No nevtot or maps histo found...";
110  return 0;
111  }
112 
113  int ev = nevtot->getEntries();
114  if(ev <= 0) {
115  edm::LogError("HcalDigisClient") << "normalization factor <= 0!";
116  return 0;
117  }
118 
119  float fev = (float) nevtot->getEntries();
120 
121  float sumphi_1, sumphi_2, sumphi_3, sumphi_4;
122  float phi_factor;
123  float cnorm;
124 
125  int nx = ieta_iphi_occupancy_map1->getNbinsX();
126  int ny = ieta_iphi_occupancy_map1->getNbinsY();
127 
128  for (int i = 1; i <= nx; i++) {
129  for (int j = 1; j <= ny; j++) {
130 
131  // occupancies
132 
133  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth1_" + subdet_;
134  cnorm = ieta_iphi_occupancy_map1->getBinContent(i, j) / fev;
135  ieta_iphi_occupancy_map1->setBinContent(i, j, cnorm);
136 
137  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth2_" + subdet_;
138  cnorm = ieta_iphi_occupancy_map2->getBinContent(i, j) / fev;
139  ieta_iphi_occupancy_map2->setBinContent(i, j, cnorm);
140 
141  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth3_" + subdet_;
142  cnorm = ieta_iphi_occupancy_map3->getBinContent(i, j) / fev;
143  ieta_iphi_occupancy_map3->setBinContent(i, j, cnorm);
144 
145  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth4_" + subdet_;
146  cnorm = ieta_iphi_occupancy_map4->getBinContent(i, j) / fev;
147  ieta_iphi_occupancy_map4->setBinContent(i, j, cnorm);
148 
149  }
150  }
151 
152  for (int i = 1; i <= 82; i++) {
153 
154  int ieta = i - 42; // -41 -1, 0 40
155  if (ieta >= 0) ieta += 1; // -41 -1, 1 41 - to make it detector-like
156 
157  if (ieta >= -20 && ieta <= 20) {
158  phi_factor = 72.;
159  } else {
160  if (ieta >= 40 || ieta <= -40)
161  phi_factor = 18.;
162  else
163  phi_factor = 36.;
164  }
165 
166  sumphi_1 = 0.;
167  sumphi_2 = 0.;
168  sumphi_3 = 0.;
169  sumphi_4 = 0.;
170 
171  for (int iphi = 1; iphi <= 72; iphi++) {
172  sumphi_1 += ieta_iphi_occupancy_map1->getTH1()->GetBinContent(ieta_iphi_occupancy_map1->getTH1()->FindFixBin(double(ieta),double(iphi)));
173  sumphi_2 += ieta_iphi_occupancy_map2->getTH1()->GetBinContent(ieta_iphi_occupancy_map2->getTH1()->FindFixBin(double(ieta),double(iphi)));
174  sumphi_3 += ieta_iphi_occupancy_map3->getTH1()->GetBinContent(ieta_iphi_occupancy_map3->getTH1()->FindFixBin(double(ieta),double(iphi)));
175  sumphi_4 += ieta_iphi_occupancy_map4->getTH1()->GetBinContent(ieta_iphi_occupancy_map4->getTH1()->FindFixBin(double(ieta),double(iphi)));
176  }
177 
178  //REMOVED (JRD) if (ieta >= 0) ieta -= 1; // -41 -1, 0 40 - to bring back to strtmp num !!!
179  double deta = double(ieta);
180 
181  // occupancies vs ieta
182  cnorm = sumphi_1 / phi_factor;
183  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth1_" + subdet_;
184  fill1D(strtmp, deta, cnorm);
185 
186  cnorm = sumphi_2 / phi_factor;
187  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth2_" + subdet_;
188  fill1D(strtmp, deta, cnorm);
189 
190  cnorm = sumphi_3 / phi_factor;
191  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth3_" + subdet_;
192  fill1D(strtmp, deta, cnorm);
193 
194  cnorm = sumphi_4 / phi_factor;
195  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth4_" + subdet_;
196  fill1D(strtmp, deta, cnorm);
197 
198  } // end of i-loop
199 
200  return 1;
201 }
202 
204  if (!msm_->count(name)) return NULL;
205  else return msm_->find(name)->second;
206 }
207 
209  std::stringstream out;
210  out << x;
211  return out.str();
212 }
213 
214 double HcalDigisClient::integralMETH2D(MonitorElement* ME, int i0, int i1, int j0, int j1) {
215  double sum(0);
216  for (int i = i0; i <= i1; i++) {
217  for (int j = j0; j <= j1; j++) {
218  sum += ME->getBinContent(i, j);
219  }
220  }
221 
222  return sum;
223 }
224 
226  int nx = ME->getNbinsX();
227  int ny = ME->getNbinsY();
228 
229  double content(0);
230  for (int i = 1; i <= nx; i++) {
231  for (int j = 1; j <= ny; j++) {
232  content = ME->getBinContent(i, j);
233  content *= s;
234  ME->setBinContent(i, j, content);
235  }
236  }
237 }
238 
239 //define this as a plug-in
241 
std::vector< MonitorElement * > getContents(Args &&...args)
Definition: DQMStore.h:197
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void setBinContent(int binx, double content)
set content of bin (1-D)
int ib
Definition: cuy.py:660
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::string dirName_
void scaleMETH2D(MonitorElement *ME, double s)
virtual void runClient(DQMStore::IBooker &ib, DQMStore::IGetter &ig)
#define NULL
Definition: scimark2.h:8
bool ev
int HcalDigisEndjob(const std::vector< MonitorElement * > &hcalMEs, std::string subdet_)
double getEntries(void) const
get # of entries
HcalDigisClient(const edm::ParameterSet &)
int getNbinsY(void) const
get # of bins in Y-axis
std::string outputFile_
Definition: ME.h:11
void fill1D(std::string name, double X, double weight=1)
MonitorElement * monitor(std::string name)
int j
Definition: DBlmapReader.cc:9
TH1 * getTH1(void) const
std::string str(int x)
void booking(DQMStore::IBooker &ib, std::string subdetopt)
void book1D(DQMStore::IBooker &ib, std::string name, int n, double min, double max)
std::map< std::string, MonitorElement * > * msm_
double integralMETH2D(MonitorElement *ME, int i0, int i1, int j0, int j1)
std::vector< std::string > getSubdirs(void)
Definition: DQMStore.cc:305
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:329
double getBinContent(int binx) const
get content of bin (1-D)
int getNbinsX(void) const
get # of bins in X-axis