CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ESPedestalClient.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <iostream>
3 #include <fstream>
4 #include <string>
5 #include <vector>
6 
9 
11 
12 #include <TH1F.h>
13 
14 using namespace edm;
15 using namespace std;
16 
18  fg(nullptr) {
19 
20  verbose_ = ps.getUntrackedParameter<bool>("verbose", true);
21  debug_ = ps.getUntrackedParameter<bool>("debug", true);
22  prefixME_ = ps.getUntrackedParameter<string>("prefixME", "EcalPreshower");
23  lookup_ = ps.getUntrackedParameter<FileInPath>("LookupTable");
24  enableCleanup_ = ps.getUntrackedParameter<bool>("enableCleanup", false);
25  fitPedestal_ = ps.getUntrackedParameter<bool>("fitPedestal", false);
26 
27  for (int i=0; i<2; i++)
28  for (int j=0; j<2; j++)
29  for (int k=0; k<40; k++)
30  for (int m=0; m<40; m++) {
31  hPed_[i][j][k][m] = 0;
32  hTotN_[i][j][k][m] = 0;
33  }
34 
35 }
36 
38  delete fg;
39 }
40 
42 
43  dqmStore_ = dqmStore;
44 
45  if ( debug_ ) cout << "ESPedestalClient: beginJob" << endl;
46 
47  ievt_ = 0;
48  jevt_ = 0;
49 
50 }
51 
53 
54  if ( debug_ ) cout << "ESPedestalClient: beginRun" << endl;
55 
56  jevt_ = 0;
57 
58  this->setup();
59 }
60 
62 
63  if ( debug_ ) cout << "ESPedestalClient: endJob, ievt = " << ievt_ << endl;
64 
65  // Preform pedestal fit
66  char hname[300];
67  int iz = 0;
68  if (fitPedestal_) {
69 
70  if ( verbose_ ) cout<<"ESPedestalClient: Fit Pedestal"<<endl;
71 
72  for (int i=0; i<nLines_; ++i) {
73 
74  iz = (senZ_[i]==1) ? 0:1;
75 
76  for (int is=0; is<32; ++is) {
77 
78  string dirname = prefixME_ + "/ESPedestalTask/";
79  sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is+1);
80  MonitorElement *meFit = dqmStore_->get(dirname+hname);
81 
82  if (meFit==0) continue;
83  TH1F *rootHisto = meFit->getTH1F();
84  rootHisto->Fit(fg, "Q", "", 500, 1800);
85  rootHisto->Fit(fg, "RQ", "", fg->GetParameter(1)-2.*fg->GetParameter(2),fg->GetParameter(1)+2.*fg->GetParameter(2));
86  hPed_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, (int)(fg->GetParameter(1)+0.5));
87  hTotN_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, fg->GetParameter(2));
88 
89  }
90  }
91 
92  } else {
93 
94  if ( verbose_ ) cout<<"ESPedestalClient: Use Histogram Mean"<<endl;
95 
96  for (int i=0; i<nLines_; ++i) {
97 
98  iz = (senZ_[i]==1) ? 0:1;
99 
100  for (int is=0; is<32; ++is) {
101 
102  string dirname = prefixME_ + "/ESPedestalTask/";
103  sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is+1);
104  MonitorElement *meMean = dqmStore_->get(dirname+hname);
105 
106  if (meMean==0) continue;
107  TH1F *rootHisto = meMean->getTH1F();
108 
109  hPed_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, (int)(rootHisto->GetMean()+0.5));
110  hTotN_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, rootHisto->GetRMS());
111 
112  }
113  }
114  }
115 
116  this->cleanup();
117 }
118 
120 
121  if ( debug_ ) cout << "ESPedestalClient: endRun, jevt = " << jevt_ << endl;
122 
123  this->cleanup();
124 }
125 
127 
128  // read in look-up table
129  int iz, ip, ix, iy, fed, kchip, pace, bundle, fiber, optorx;
130  ifstream file;
131 
132  file.open(lookup_.fullPath().c_str());
133  if( file.is_open() ) {
134 
135  file >> nLines_;
136 
137  for (int i=0; i<nLines_; ++i) {
138  file>> iz >> ip >> ix >> iy >> fed >> kchip >> pace >> bundle >> fiber >> optorx;
139 
140  senZ_[i] = iz;
141  senP_[i] = ip;
142  senX_[i] = ix;
143  senY_[i] = iy;
144  }
145 
146  } else {
147  cout<<"ESPedestalClient : Look up table file can not be found in "<<lookup_.fullPath().c_str()<<endl;
148  }
149 
150  // define histograms
151  dqmStore_->setCurrentFolder(prefixME_+"/ESPedestalClient");
152 
153  char hname[300];
154  for (int i=0; i<nLines_; ++i) {
155 
156  iz = (senZ_[i]==1) ? 0:1;
157 
158  sprintf(hname, "Ped Z %d P %d X %d Y %d", senZ_[i], senP_[i], senX_[i], senY_[i]);
159  hPed_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1] = dqmStore_->book1D(hname, hname, 32, 0, 32);
160 
161  sprintf(hname, "Total Noise Z %d P %d X %d Y %d", senZ_[i], senP_[i], senX_[i], senY_[i]);
162  hTotN_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1] = dqmStore_->book1D(hname, hname, 32, 0, 32);
163  }
164 
165  fg = new TF1("fg", "gaus");
166 
167 }
168 
170 
171  if( ! enableCleanup_ ) return;
172 
173  if ( debug_ ) cout << "ESPedestalClient: cleanup" << endl;
174 
175  for (int i=0; i<2; i++)
176  for (int j=0; j<2; j++)
177  for (int k=0; k<40; k++)
178  for (int m=0; m<40; m++) {
179  hPed_[i][j][k][m] = 0;
180  hTotN_[i][j][k][m] = 0;
181  }
182 
183 }
184 
186 
187  ievt_++;
188  jevt_++;
189 
190 }
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)
edm::FileInPath lookup_
void beginJob(DQMStore *dqmStore)
MonitorElement * hTotN_[2][2][40][40]
ESPedestalClient(const edm::ParameterSet &ps)
#define nullptr
virtual ~ESPedestalClient()
int j
Definition: DBlmapReader.cc:9
DQMStore * dqmStore_
int k[5][pyjets_maxn]
TH1F * getTH1F(void) const
std::string prefixME_
tuple cout
Definition: gather_cfg.py:121
std::string fullPath() const
Definition: FileInPath.cc:165
MonitorElement * hPed_[2][2][40][40]