CMS 3D CMS Logo

EcalTPGPedestals_PayloadInspector.cc
Go to the documentation of this file.
7 
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.h"
12 #include "TCanvas.h"
13 #include "TLine.h"
14 #include "TStyle.h"
15 #include "TLatex.h"
16 #include "TPave.h"
17 #include "TPaveStats.h"
18 #include <string>
19 #include <fstream>
20 
21 namespace {
22  enum {kEBChannels = 61200, kEEChannels = 14648, kGains = 3, kSides = 2};
23  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360}; // barrel lower and upper bounds on eta and phi
24  enum {IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100}; // endcaps lower and upper bounds on x and y
25  int gainValues[kGains] = {12, 6, 1};
26 
27  /**************************************************
28  2d plot of ECAL TPGPedestals of 1 IOV
29  **************************************************/
30  class EcalTPGPedestalsPlot : public cond::payloadInspector::PlotImage<EcalTPGPedestals> {
31  public:
32  EcalTPGPedestalsPlot() : cond::payloadInspector::PlotImage<EcalTPGPedestals>("ECAL Gain Ratios - map ") {
33  setSingleIov( true );
34  }
35 
36  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
37  TH2F** barrel_m = new TH2F*[kGains];
38  TH2F** endc_p_m = new TH2F*[kGains];
39  TH2F** endc_m_m = new TH2F*[kGains];
40  float mEBmin[kGains], mEEmin[kGains], mEBmax[kGains], mEEmax[kGains];
41  for (int gainId = 0; gainId < kGains; gainId++) {
42  barrel_m[gainId] = new TH2F(Form("EBm%i", gainId), Form("EB mean_x%i ", gainValues[gainId]), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
43  endc_p_m[gainId] = new TH2F(Form("EE+m%i", gainId), Form("EE+ mean_x%i", gainValues[gainId]), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
44  endc_m_m[gainId] = new TH2F(Form("EE-m%i", gainId), Form("EE- mean_x%i", gainValues[gainId]), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
45  mEBmin[gainId] = 10.;
46  mEEmin[gainId] = 10.;
47  mEBmax[gainId] = -10.;
48  mEEmax[gainId] = -10.;
49  }
50 
51  // std::ofstream fout;
52  // fout.open("./bid.txt");
53  auto iov = iovs.front();
54  std::shared_ptr<EcalTPGPedestals> payload = fetchPayload( std::get<1>(iov) );
55  unsigned int run = std::get<0>(iov);
56  if( payload.get() ){
57  for (int sign=0; sign < kSides; sign++) {
58  int thesign = sign==1 ? 1:-1;
59 
60  for (int ieta = 0; ieta < MAX_IETA; ieta++) {
61  for (int iphi = 0; iphi < MAX_IPHI; iphi++) {
62  EBDetId id((ieta+1)*thesign, iphi+1);
63  float y = -1 - ieta;
64  if(sign == 1) y = ieta;
65  float val = (*payload)[id.rawId()].mean_x12;
66  barrel_m[0]->Fill(iphi, y, val);
67  if(val < mEBmin[0]) mEBmin[0] = val;
68  if(val > mEBmax[0]) mEBmax[0] = val;
69  val = (*payload)[id.rawId()].mean_x6;
70  barrel_m[1]->Fill(iphi, y, val);
71  if(val < mEBmin[1]) mEBmin[1] = val;
72  if(val > mEBmax[1]) mEBmax[1] = val;
73  val = (*payload)[id.rawId()].mean_x1;
74  barrel_m[2]->Fill(iphi, y, val);
75  if(val < mEBmin[2]) mEBmin[2] = val;
76  if(val > mEBmax[2]) mEBmax[2] = val;
77  } // iphi
78  } // ieta
79 
80  for (int ix = 0; ix < IX_MAX; ix++) {
81  for (int iy = 0; iy < IY_MAX; iy++) {
82  if (! EEDetId::validDetId(ix+1,iy+1,thesign)) continue;
83  EEDetId id(ix+1,iy+1,thesign);
84  float val = (*payload)[id.rawId()].mean_x12;
85  if (thesign==1) endc_p_m[0]->Fill(ix + 1, iy + 1, val);
86  else endc_m_m[0]->Fill(ix + 1, iy + 1, val);
87  if(val < mEEmin[0]) mEEmin[0] = val;
88  if(val > mEEmax[0]) mEEmax[0] = val;
89  val = (*payload)[id.rawId()].mean_x6;
90  if (thesign==1) endc_p_m[1]->Fill(ix + 1, iy + 1, val);
91  else endc_m_m[1]->Fill(ix + 1, iy + 1, val);
92  if(val < mEEmin[1]) mEEmin[1] = val;
93  if(val > mEEmax[1]) mEEmax[1] = val;
94  val = (*payload)[id.rawId()].mean_x1;
95  if (thesign==1) endc_p_m[2]->Fill(ix + 1, iy + 1, val);
96  else endc_m_m[2]->Fill(ix + 1, iy + 1, val);
97  if(val < mEEmin[2]) mEEmin[2] = val;
98  if(val > mEEmax[2]) mEEmax[2] = val;
99  // fout << " x " << ix << " y " << " val " << val << std::endl;
100  } // iy
101  } // ix
102  } // side
103  } // if payload.get()
104  else return false;
105  // std::cout << " min " << rEEmin[2] << " max " << rEEmax[2] << std::endl;
106  // fout.close();
107 
108  gStyle->SetPalette(1);
109  gStyle->SetOptStat(0);
110  TCanvas canvas("CC map","CC map",1200,900);
111  TLatex t1;
112  t1.SetNDC();
113  t1.SetTextAlign(26);
114  t1.SetTextSize(0.05);
115  t1.DrawLatex(0.5, 0.96, Form("Ecal Gain TPGPedestals, IOV %i", run));
116 
117  float xmi[3] = {0.0 , 0.22, 0.78};
118  float xma[3] = {0.22, 0.78, 1.00};
119  TPad*** pad = new TPad**[kGains];
120  for (int gId = 0; gId < kGains; gId++) {
121  pad[gId] = new TPad*[3];
122  for (int obj = 0; obj < 3; obj++) {
123  float yma = 0.94 - (0.32 * gId);
124  float ymi = yma - 0.30;
125  pad[gId][obj] = new TPad(Form("p_%i_%i", obj, gId),Form("p_%i_%i", obj, gId),
126  xmi[obj], ymi, xma[obj], yma);
127  pad[gId][obj]->Draw();
128  }
129  }
130 
131  for (int gId = 0; gId < kGains; gId++) {
132  pad[gId][0]->cd();
133  DrawEE(endc_m_m[gId], mEEmin[gId], mEEmax[gId]);
134  pad[gId][1]->cd();
135  DrawEB(barrel_m[gId], mEBmin[gId], mEBmax[gId]);
136  pad[gId][2]->cd();
137  DrawEE(endc_p_m[gId], mEEmin[gId], mEEmax[gId]);
138  }
139 
140  std::string ImageName(m_imageFileName);
141  canvas.SaveAs(ImageName.c_str());
142  return true;
143  }// fill method
144  };
145 
146  /******************************************************************
147  2d plot of ECAL TPGPedestals difference between 2 IOVs
148  ******************************************************************/
149  class EcalTPGPedestalsDiff : public cond::payloadInspector::PlotImage<EcalTPGPedestals> {
150 
151  public:
152  EcalTPGPedestalsDiff() : cond::payloadInspector::PlotImage<EcalTPGPedestals>("ECAL Gain Ratios difference") {
153  setSingleIov(false);
154  }
155 
156  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
157  TH2F** barrel_m = new TH2F*[kGains];
158  TH2F** endc_p_m = new TH2F*[kGains];
159  TH2F** endc_m_m = new TH2F*[kGains];
160  float mEBmin[kGains], mEEmin[kGains], mEBmax[kGains], mEEmax[kGains];
161  float mEB[kGains][kEBChannels], mEE[kGains][kEEChannels];
162  for (int gainId = 0; gainId < kGains; gainId++) {
163  barrel_m[gainId] = new TH2F(Form("EBm%i", gainId), Form("EB mean_x%i ", gainValues[gainId]), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
164  endc_p_m[gainId] = new TH2F(Form("EE+m%i", gainId), Form("EE+ mean_x%i", gainValues[gainId]), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
165  endc_m_m[gainId] = new TH2F(Form("EE-m%i", gainId), Form("EE- mean_x%i", gainValues[gainId]), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
166  mEBmin[gainId] = 10.;
167  mEEmin[gainId] = 10.;
168  mEBmax[gainId] = -10.;
169  mEEmax[gainId] = -10.;
170  }
171 
172  unsigned int run[2], irun = 0;
173  //float gEB[3][kEBChannels], gEE[3][kEEChannels];
174  for ( auto const & iov: iovs) {
175  std::shared_ptr<EcalTPGPedestals> payload = fetchPayload( std::get<1>(iov) );
176  run[irun] = std::get<0>(iov);
177  if( payload.get() ){
178  for (int sign=0; sign < kSides; sign++) {
179  int thesign = sign==1 ? 1:-1;
180 
181  for (int ieta = 0; ieta < MAX_IETA; ieta++) {
182  for (int iphi = 0; iphi < MAX_IPHI; iphi++) {
183  EBDetId id((ieta+1)*thesign, iphi+1);
184  int hashindex = id.hashedIndex();
185  float y = -1 - ieta;
186  if(sign == 1) y = ieta;
187  float val = (*payload)[id.rawId()].mean_x12;
188  if(irun == 0) {
189  mEB[0][hashindex] = val;
190  }
191  else {
192  float diff = val - mEB[0][hashindex];
193  barrel_m[0]->Fill(iphi, y, diff);
194  if(diff < mEBmin[0]) mEBmin[0] = diff;
195  if(diff > mEBmax[0]) mEBmax[0] = diff;
196  }
197  val = (*payload)[id.rawId()].mean_x6;
198  if(irun == 0) {
199  mEB[1][hashindex] = val;
200  }
201  else {
202  float diff = val - mEB[1][hashindex];
203  barrel_m[1]->Fill(iphi, y, diff);
204  if(diff < mEBmin[1]) mEBmin[1] = diff;
205  if(diff > mEBmax[1]) mEBmax[1] = diff;
206  }
207  val = (*payload)[id.rawId()].mean_x1;
208  if(irun == 0) {
209  mEB[2][hashindex] = val;
210  }
211  else {
212  float diff = val - mEB[2][hashindex];
213  barrel_m[2]->Fill(iphi, y, diff);
214  if(diff < mEBmin[2]) mEBmin[2] = diff;
215  if(diff > mEBmax[2]) mEBmax[2] = diff;
216  }
217  } // iphi
218  } // ieta
219 
220  for (int ix = 0; ix < IX_MAX; ix++) {
221  for (int iy = 0; iy < IY_MAX; iy++) {
222  if (! EEDetId::validDetId(ix+1,iy+1,thesign)) continue;
223  EEDetId id(ix+1,iy+1,thesign);
224  int hashindex = id.hashedIndex();
225  float val = (*payload)[id.rawId()].mean_x12;
226  if(irun == 0) {
227  mEE[0][hashindex] = val;
228  }
229  else {
230  float diff = val - mEE[0][hashindex];
231  if (thesign==1) endc_p_m[0]->Fill(ix + 1, iy + 1, diff);
232  else endc_m_m[0]->Fill(ix + 1, iy + 1, diff);
233  if(diff < mEEmin[0]) mEEmin[0] = diff;
234  if(diff > mEEmax[0]) mEEmax[0] = diff;
235  }
236  val = (*payload)[id.rawId()].mean_x6;
237  if(irun == 0) {
238  mEE[1][hashindex] = val;
239  }
240  else {
241  float diff = val - mEE[1][hashindex];
242  if (thesign==1) endc_p_m[1]->Fill(ix + 1, iy + 1, diff);
243  else endc_m_m[1]->Fill(ix + 1, iy + 1, diff);
244  if(diff < mEEmin[1]) mEEmin[1] = diff;
245  if(diff > mEEmax[1]) mEEmax[1] = diff;
246  }
247  val = (*payload)[id.rawId()].mean_x1;
248  if(irun == 0) {
249  mEE[2][hashindex] = val;
250  }
251  else {
252  float diff = val - mEE[2][hashindex];
253  if (thesign==1) endc_p_m[2]->Fill(ix + 1, iy + 1, diff);
254  else endc_m_m[2]->Fill(ix + 1, iy + 1, diff);
255  if(diff < mEEmin[2]) mEEmin[2] = diff;
256  if(diff > mEEmax[2]) mEEmax[2] = diff;
257  }
258  // fout << " x " << ix << " y " << " diff " << diff << std::endl;
259  } // iy
260  } // ix
261  } // side
262  } // if payload.get()
263  else return false;
264  irun++;
265  } // loop over IOVs
266 
267  gStyle->SetPalette(1);
268  gStyle->SetOptStat(0);
269  TCanvas canvas("CC map","CC map",1200, 900);
270  TLatex t1;
271  t1.SetNDC();
272  t1.SetTextAlign(26);
273  t1.SetTextSize(0.05);
274  t1.DrawLatex(0.5, 0.96, Form("Ecal TPGPedestals, IOV %i - %i", run[1], run[0]));
275 
276  float xmi[3] = {0.0 , 0.22, 0.78};
277  float xma[3] = {0.22, 0.78, 1.00};
278  TPad*** pad = new TPad**[kGains];
279  for (int gId = 0; gId < kGains; gId++) {
280  pad[gId] = new TPad*[3];
281  for (int obj = 0; obj < 3; obj++) {
282  float yma = 0.94 - (0.32 * gId);
283  float ymi = yma - 0.30;
284  pad[gId][obj] = new TPad(Form("p_%i_%i", obj, gId),Form("p_%i_%i", obj, gId),
285  xmi[obj], ymi, xma[obj], yma);
286  pad[gId][obj]->Draw();
287  }
288  }
289 
290  for (int gId = 0; gId < kGains; gId++) {
291  pad[gId][0]->cd();
292  DrawEE(endc_m_m[gId], mEEmin[gId], mEEmax[gId]);
293  pad[gId][1]->cd();
294  DrawEB(barrel_m[gId], mEBmin[gId], mEBmax[gId]);
295  pad[gId][2]->cd();
296  DrawEE(endc_p_m[gId], mEEmin[gId], mEEmax[gId]);
297  }
298 
299  std::string ImageName(m_imageFileName);
300  canvas.SaveAs(ImageName.c_str());
301  return true;
302  }// fill method
303  };
304 
305 } // close namespace
306 
307  // Register the classes as boost python plugin
309  PAYLOAD_INSPECTOR_CLASS(EcalTPGPedestalsPlot);
310  PAYLOAD_INSPECTOR_CLASS(EcalTPGPedestalsDiff);
311 }
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
static const int kSides
const Int_t gainValues[kGains]
void DrawEE(TH2F *endc, float min, float max)
Definition: EcalDrawUtils.h:29
void DrawEB(TH2F *ebmap, float min, float max)
Definition: EcalDrawUtils.h:4
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
virtual bool fill(const std::vector< std::tuple< cond::Time_t, cond::Hash > > &iovs)=0
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
static bool validDetId(int crystal_ix, int crystal_iy, int iz)
Definition: EEDetId.h:248
constexpr int gainId(sample_type sample)
get the gainId (2 bits)
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:482
const Int_t kGains