CMS 3D CMS Logo

ESPedestalClient.cc
Go to the documentation of this file.
2 
5 
6 #include <memory>
7 #include <iostream>
8 #include <fstream>
9 #include <string>
10 
11 using namespace std;
12 
14  ESClient(ps),
15  fitPedestal_(ps.getUntrackedParameter<bool>("fitPedestal")),
16  fg_(new TF1("fg", "gaus")),
17  senZ_(),
18  senP_(),
19  senX_(),
20  senY_()
21 {
22  for (int i=0; i<2; i++)
23  for (int j=0; j<2; j++)
24  for (int k=0; k<40; k++)
25  for (int m=0; m<40; m++) {
26  hPed_[i][j][k][m] = nullptr;
27  hTotN_[i][j][k][m] = nullptr;
28  }
29 
30  std::string lutPath(ps.getUntrackedParameter<edm::FileInPath>("LookupTable").fullPath());
31 
32  // read in look-up table
33  int nLines, iz, ip, ix, iy, fed, kchip, pace, bundle, fiber, optorx;
34  ifstream file(lutPath);
35 
36  if( file.is_open() ) {
37 
38  file >> nLines;
39 
40  for (int i=0; i<nLines; ++i) {
41  file>> iz >> ip >> ix >> iy >> fed >> kchip >> pace >> bundle >> fiber >> optorx;
42 
43  senZ_.push_back(iz);
44  senP_.push_back(ip);
45  senX_.push_back(ix);
46  senY_.push_back(iy);
47  }
48 
49  file.close();
50  } else {
51  cout<<"ESPedestalClient : Look up table file can not be found in "<<lutPath<<endl;
52  }
53 }
54 
56  delete fg_;
57 }
58 
60 
61  if ( debug_ ) cout << "ESPedestalClient: endJob" << endl;
62 
63  // Preform pedestal fit
64  char hname[300];
65  int iz = 0;
66  if (fitPedestal_) {
67 
68  if ( verbose_ ) cout<<"ESPedestalClient: Fit Pedestal"<<endl;
69 
70  for (unsigned i=0; i<senZ_.size(); ++i) {
71 
72  iz = (senZ_[i]==1) ? 0:1;
73 
74  for (int is=0; is<32; ++is) {
75 
76  string dirname = prefixME_ + "/ESPedestalTask/";
77  sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is+1);
78  MonitorElement *meFit = _igetter.get(dirname+hname);
79 
80  if (meFit==nullptr) continue;
81  TH1F *rootHisto = meFit->getTH1F();
82 
83  rootHisto->Fit(fg_, "Q", "", 500, 1800);
84  rootHisto->Fit(fg_, "RQ", "", fg_->GetParameter(1)-2.*fg_->GetParameter(2),fg_->GetParameter(1)+2.*fg_->GetParameter(2));
85  hPed_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, (int)(fg_->GetParameter(1)+0.5));
86  hTotN_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, fg_->GetParameter(2));
87  }
88  }
89 
90  } else {
91 
92  if ( verbose_ ) cout<<"ESPedestalClient: Use Histogram Mean"<<endl;
93 
94  for (unsigned i=0; i<senZ_.size(); ++i) {
95 
96  iz = (senZ_[i]==1) ? 0:1;
97 
98  for (int is=0; is<32; ++is) {
99 
100  string dirname = prefixME_ + "/ESPedestalTask/";
101  sprintf(hname, "ADC Z %d P %d X %d Y %d Str %d", senZ_[i], senP_[i], senX_[i], senY_[i], is+1);
102  MonitorElement *meMean = _igetter.get(dirname+hname);
103 
104  if (meMean==nullptr) continue;
105  TH1F *rootHisto = meMean->getTH1F();
106 
107  hPed_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, (int)(rootHisto->GetMean()+0.5));
108  hTotN_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1]->setBinContent(is+1, rootHisto->GetRMS());
109 
110  }
111  }
112  }
113 }
114 
116 
117  // define histograms
118  _ibooker.setCurrentFolder(prefixME_+"/ESPedestalClient");
119 
120  char hname[300];
121  for (unsigned i=0; i<senZ_.size(); ++i) {
122 
123  int iz = (senZ_[i]==1) ? 0:1;
124 
125  sprintf(hname, "Ped Z %d P %d X %d Y %d", senZ_[i], senP_[i], senX_[i], senY_[i]);
126  hPed_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1] = _ibooker.book1D(hname, hname, 32, 0, 32);
127 
128  sprintf(hname, "Total Noise Z %d P %d X %d Y %d", senZ_[i], senP_[i], senX_[i], senY_[i]);
129  hTotN_[iz][senP_[i]-1][senX_[i]-1][senY_[i]-1] = _ibooker.book1D(hname, hname, 32, 0, 32);
130  }
131 
132 }
std::vector< int > senP_
T getUntrackedParameter(std::string const &, T const &) const
void setBinContent(int binx, double content)
set content of bin (1-D)
TH1F * getTH1F() const
MonitorElement * hTotN_[2][2][40][40]
~ESPedestalClient() override
std::vector< int > senX_
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::string prefixME_
Definition: ESClient.h:31
bool debug_
Definition: ESClient.h:34
void endJobAnalyze(DQMStore::IGetter &) override
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
std::vector< int > senZ_
MonitorElement * get(std::string const &path)
Definition: DQMStore.cc:303
int k[5][pyjets_maxn]
std::vector< int > senY_
bool verbose_
Definition: ESClient.h:33
std::string fullPath() const
Definition: FileInPath.cc:197
void book(DQMStore::IBooker &) override
ESPedestalClient(const edm::ParameterSet &)
MonitorElement * hPed_[2][2][40][40]