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 // $Id: HcalDigisClient.cc,v 1.3 2012/08/15 12:47:37 abdullin Exp $
18 //
19 //
20 
21 
22 // system include files
23 
26  outputFile_ = iConfig.getUntrackedParameter<std::string > ("outputFile", "HcalDigisClient.root");
27  dirName_ = iConfig.getParameter<std::string > ("DQMDirName");
28  if (!dbe_) edm::LogError("HcalDigisClient") << "unable to get DQMStore service, upshot is no client histograms will be made";
29  msm_ = new std::map<std::string, MonitorElement*>();
30  //if (iConfig.getUntrackedParameter<bool>("DQMStore", false)) if (dbe_) dbe_->setVerbose(0);
31 
32  // std::cout << "dirName: " << dirName_ << std::endl;
33  //dbe_->setCurrentFolder(dirName_);
34  dbe_->setCurrentFolder("HcalDigisV/HcalDigiTask");
35 
36  booking("HB");
37  booking("HE");
38  booking("HO");
39  booking("HF");
40 }
41 
43  using namespace edm;
44 
45 
46 }
47 
49 
50  std::string strtmp;
51  HistLim ietaLim(82, -41., 41.);
52 
53  for (int depth = 1; depth <= 4; depth++) {
54  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth" + str(depth) + "_" + subdetopt;
55  book1D(strtmp, ietaLim);
56  }
57 
58 }
59 
61  if (!dbe_) return; //we dont have the DQMStore so we cant do anything
63  std::vector<MonitorElement*> hcalMEs;
64  // Since out folders are fixed to three, we can just go over these three folders
65  // i.e., CaloTowersV/CaloTowersTask, HcalRecHitsV/HcalRecHitTask, NoiseRatesV/NoiseRatesTask.
66  std::vector<std::string> fullPathHLTFolders = dbe_->getSubdirs();
67  for (unsigned int i = 0; i < fullPathHLTFolders.size(); i++) {
68  dbe_->setCurrentFolder(fullPathHLTFolders[i]);
69  std::vector<std::string> fullSubPathHLTFolders = dbe_->getSubdirs();
70  for (unsigned int j = 0; j < fullSubPathHLTFolders.size(); j++) {
71  if (strcmp(fullSubPathHLTFolders[j].c_str(), "HcalDigisV/HcalDigiTask") == 0) {
72  hcalMEs = dbe_->getContents(fullSubPathHLTFolders[j]);
73  if (!HcalDigisEndjob(hcalMEs, "HB"))
74  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HB";
75  if (!HcalDigisEndjob(hcalMEs, "HE"))
76  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HE";
77  if (!HcalDigisEndjob(hcalMEs, "HO"))
78  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HO";
79  if (!HcalDigisEndjob(hcalMEs, "HF"))
80  edm::LogError("HcalDigisClient") << "Error in HcalDigisEndjob! HF"; }
81  }
82  }
83 }
84 
85 int HcalDigisClient::HcalDigisEndjob(const std::vector<MonitorElement*> &hcalMEs, std::string subdet_) {
86 
87  using namespace std;
88  string strtmp;
89 
90 
91  MonitorElement * nevtot(0);
92  MonitorElement * ieta_iphi_occupancy_map1(0);
93  MonitorElement * ieta_iphi_occupancy_map2(0);
94  MonitorElement * ieta_iphi_occupancy_map3(0);
95  MonitorElement * ieta_iphi_occupancy_map4(0);
96 
97 
98  std::cout << " Number of histos " << hcalMEs.size() << std::endl;
99 
100  for (unsigned int ih = 0; ih < hcalMEs.size(); ih++) {
101  if (hcalMEs[ih]->getName() == "nevtot") nevtot = hcalMEs[ih];
102 
103  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth1_" + subdet_;
104  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map1 = hcalMEs[ih];
105  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth2_" + subdet_;
106  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map2 = hcalMEs[ih];
107  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth3_" + subdet_;
108  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map3 = hcalMEs[ih];
109  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth4_" + subdet_;
110  if (hcalMEs[ih]->getName() == strtmp) ieta_iphi_occupancy_map4 = hcalMEs[ih];
111 
112  }//
113 
114  if (nevtot == 0 ||
115  ieta_iphi_occupancy_map1 == 0 ||
116  ieta_iphi_occupancy_map2 == 0 ||
117  ieta_iphi_occupancy_map3 == 0 ||
118  ieta_iphi_occupancy_map4 == 0
119  ) {
120  edm::LogError("HcalDigisClient") << "No nevtot or maps histo found...";
121  return 0;
122  }
123 
124  int ev = nevtot->getEntries();
125  if(ev <= 0) {
126  edm::LogError("HcalDigisClient") << "normalization factor <= 0!";
127  return 0;
128  }
129 
130  float fev = (float) nevtot->getEntries();
131 
132  int nx = ieta_iphi_occupancy_map1->getNbinsX();
133  int ny = ieta_iphi_occupancy_map1->getNbinsY();
134  float sumphi_1, sumphi_2, sumphi_3, sumphi_4;
135  float phi_factor;
136  float cnorm;
137 
138  for (int i = 1; i <= nx; i++) {
139  sumphi_1 = 0.;
140  sumphi_2 = 0.;
141  sumphi_3 = 0.;
142  sumphi_4 = 0.;
143 
144  for (int j = 1; j <= ny; j++) {
145 
146  // occupancies
147 
148  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth1_" + subdet_;
149  cnorm = ieta_iphi_occupancy_map1->getBinContent(i, j) / fev;
150  ieta_iphi_occupancy_map1->setBinContent(i, j, cnorm);
151  sumphi_1 += ieta_iphi_occupancy_map1->getBinContent(i, j);
152 
153  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth2_" + subdet_;
154  cnorm = ieta_iphi_occupancy_map2->getBinContent(i, j) / fev;
155  ieta_iphi_occupancy_map2->setBinContent(i, j, cnorm);
156  sumphi_2 += ieta_iphi_occupancy_map2->getBinContent(i, j);
157 
158  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth3_" + subdet_;
159  cnorm = ieta_iphi_occupancy_map3->getBinContent(i, j) / fev;
160  ieta_iphi_occupancy_map3->setBinContent(i, j, cnorm);
161  sumphi_3 += ieta_iphi_occupancy_map3->getBinContent(i, j);
162 
163  strtmp = "HcalDigiTask_ieta_iphi_occupancy_map_depth4_" + subdet_;
164  cnorm = ieta_iphi_occupancy_map4->getBinContent(i, j) / fev;
165  ieta_iphi_occupancy_map4->setBinContent(i, j, cnorm);
166  sumphi_4 += ieta_iphi_occupancy_map4->getBinContent(i, j);
167 
168  }
169 
170  int ieta = i - 42; // -41 -1, 0 40
171  if (ieta >= 0) ieta += 1; // -41 -1, 1 41 - to make it detector-like
172 
173  if (ieta >= -20 && ieta <= 20) {
174  phi_factor = 72.;
175  } else {
176  if (ieta >= 40 || ieta <= -40)
177  phi_factor = 18.;
178  else
179  phi_factor = 36.;
180  }
181 
182  if (ieta >= 0) ieta -= 1; // -41 -1, 0 40 - to bring back to strtmp num !!!
183  double deta = double(ieta);
184 
185  // occupancies vs ieta
186  cnorm = sumphi_1 / phi_factor;
187  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth1_" + subdet_;
188  fill1D(strtmp, deta, cnorm);
189 
190  cnorm = sumphi_2 / phi_factor;
191  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth2_" + subdet_;
192  fill1D(strtmp, deta, cnorm);
193 
194  cnorm = sumphi_3 / phi_factor;
195  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth3_" + subdet_;
196  fill1D(strtmp, deta, cnorm);
197 
198  cnorm = sumphi_4 / phi_factor;
199  strtmp = "HcalDigiTask_occupancy_vs_ieta_depth4_" + subdet_;
200  fill1D(strtmp, deta, cnorm);
201 
202  } // end of i-loop
203 
204  return 1;
205 }
206 
208  if (!msm_->count(name)) return NULL;
209  else return msm_->find(name)->second;
210 }
211 
213  std::stringstream out;
214  out << x;
215  return out.str();
216 }
217 
218 double HcalDigisClient::integralMETH2D(MonitorElement* ME, int i0, int i1, int j0, int j1) {
219  double sum(0);
220  for (int i = i0; i <= i1; i++) {
221  for (int j = j0; j <= j1; j++) {
222  sum += ME->getBinContent(i, j);
223  }
224  }
225 
226  return sum;
227 }
228 
230  int nx = ME->getNbinsX();
231  int ny = ME->getNbinsY();
232 
233  double content(0);
234  for (int i = 1; i <= nx; i++) {
235  for (int j = 1; j <= ny; j++) {
236  content = ME->getBinContent(i, j);
237  content *= s;
238  ME->setBinContent(i, j, content);
239  }
240  }
241 }
242 
243 //define this as a plug-in
245 
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)
std::vector< std::string > getSubdirs(void) const
Definition: DQMStore.cc:1424
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void analyze(const edm::Event &, const edm::EventSetup &)
std::string dirName_
void scaleMETH2D(MonitorElement *ME, double s)
#define NULL
Definition: scimark2.h:8
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)
int iEvent
Definition: GenABIO.cc:243
void booking(std::string subdetopt)
void book1D(std::string name, int n, double min, double max)
MonitorElement * monitor(std::string name)
int j
Definition: DBlmapReader.cc:9
std::vector< MonitorElement * > getContents(const std::string &path) const
Definition: DQMStore.cc:1502
std::string str(int x)
tuple out
Definition: dbtoconf.py:99
std::map< std::string, MonitorElement * > * msm_
double integralMETH2D(MonitorElement *ME, int i0, int i1, int j0, int j1)
double getBinContent(int binx) const
get content of bin (1-D)
int getNbinsX(void) const
get # of bins in X-axis
virtual void runClient()
tuple cout
Definition: gather_cfg.py:121
Definition: DDAxes.h:10
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434