CMS 3D CMS Logo

EcalPulseShapes_PayloadInspector.cc
Go to the documentation of this file.
6 
7 // the data format of the condition to be inspected
9 
10 #include "TH2F.h"
11 #include "TProfile.h"
12 #include "TCanvas.h"
13 #include "TStyle.h"
14 #include "TLine.h"
15 #include "TLatex.h"
16 
17 #include <string>
18 
19 namespace {
20  enum {kEBChannels = 61200, kEEChannels = 14648, kSides = 2, kRMS = 3, TEMPLATESAMPLES = 12};
21  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360}; // barrel lower and upper bounds on eta and phi
22  enum {IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100}; // endcaps lower and upper bounds on x and y
23 
24  /*****************************************************
25  2d plot of ECAL PulseShapes of 1 IOV
26  *****************************************************/
27  class EcalPulseShapesPlot : public cond::payloadInspector::PlotImage<EcalPulseShapes> {
28 
29  public:
30  EcalPulseShapesPlot() : cond::payloadInspector::PlotImage<EcalPulseShapes>("ECAL PulseShapes - map ") {
31  setSingleIov(true);
32  }
33 
34  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
35  TH2F** barrel = new TH2F*[TEMPLATESAMPLES];
36  TH2F** endc_p = new TH2F*[TEMPLATESAMPLES];
37  TH2F** endc_m = new TH2F*[TEMPLATESAMPLES];
38  double EBmean[TEMPLATESAMPLES], EBrms[TEMPLATESAMPLES], EEmean[TEMPLATESAMPLES], EErms[TEMPLATESAMPLES];
39  int EBtot[TEMPLATESAMPLES], EEtot[TEMPLATESAMPLES];
40  double pEBmin[TEMPLATESAMPLES], pEBmax[TEMPLATESAMPLES], pEEmin[TEMPLATESAMPLES], pEEmax[TEMPLATESAMPLES];
41  for(int s = 0; s < TEMPLATESAMPLES; ++s) {
42  barrel[s] = new TH2F(Form("EBs%i",s),Form("sample %i EB",s), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
43  endc_p[s] = new TH2F(Form("EE+s%i",s),Form("sample %i EE+",s), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
44  endc_m[s] = new TH2F(Form("EE-s%i",s),Form("sample %i EE-",s), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
45  EBmean[s] = 0.;
46  EBrms[s] = 0.;
47  EEmean[s] = 0.;
48  EErms[s] = 0.;
49  EBtot[s] = 0;
50  EEtot[s] = 0;
51  /* pEBmin[s] = 10.;
52  pEBmax[s] = -10.;
53  pEEmin[s] = 10.;
54  pEEmax[s] = -10.;*/
55  }
56 
57  auto iov = iovs.front();
58  std::shared_ptr<EcalPulseShapes> payload = fetchPayload( std::get<1>(iov) );
59  unsigned int run = std::get<0>(iov);
60  if( payload.get() ){
61  // int chan = 0;
62  for (int ieta = -MAX_IETA; ieta <= MAX_IETA; ieta++) {
63  Double_t eta = (Double_t)ieta;
64  if(ieta == 0) continue;
65  else if(ieta > 0.) eta = eta - 0.5; // 0.5 to 84.5
66  else eta = eta + 0.5; // -84.5 to -0.5
67  for (int iphi = 1; iphi <= MAX_IPHI; iphi++) {
68  Double_t phi = (Double_t)iphi - 0.5;
69  EBDetId id(ieta, iphi);
70  for(int s = 0; s < TEMPLATESAMPLES; s++) {
71  double val = (*payload)[id.rawId()].pdfval[s];
72  barrel[s]->Fill(phi, eta, val);
73  EBmean[s] = EBmean[s] + val;
74  EBrms[s] = EBrms[s] + val * val;
75  EBtot[s]++;
76  // if(val < pEBmin[s]) pEBmin[s] = val;
77  // if(val > pEBmax[s]) pEBmax[s] = val;
78  }
79  }
80  }
81 
82  for (int sign = 0; sign < kSides; sign++) {
83  int thesign = sign==1 ? 1:-1;
84  for (int ix = 1; ix <= IX_MAX; ix++) {
85  for (int iy = 1; iy <= IY_MAX; iy++) {
86  if (! EEDetId::validDetId(ix, iy, thesign)) continue;
87  EEDetId id(ix, iy, thesign);
88  for(int s = 0; s < TEMPLATESAMPLES; s++) {
89  double val = (*payload)[id.rawId()].pdfval[s];
90  EEmean[s] = EEmean[s] + val;
91  EErms[s] = EErms[s] + val * val;
92  EEtot[s]++;
93  // if(val < pEEmin[s]) pEEmin[s] = val;
94  // if(val > pEEmax[s]) pEEmax[s] = val;
95  if (thesign==1)
96  endc_p[s]->Fill(ix, iy, val);
97  else
98  endc_m[s]->Fill(ix, iy, val);
99  }
100  }// iy
101  } // ix
102  } // side
103  } // payload
104 
105  for(int s = 0; s < TEMPLATESAMPLES; s++) {
106  std::cout << "EB sample " << s << " mean " << EBmean[s] << " rms " << EBrms[s] << " entries " << EBtot[s]
107  << " EE mean " << EEmean[s] << " rms " << EErms[s] << " entries " << EEtot[s] << std::endl;
108  double vt =(double)EBtot[s];
109  EBmean[s] = EBmean[s] / vt;
110  EBrms[s] = EBrms[s] / vt - (EBmean[s] * EBmean[s]);
111  if(EBrms[s] > 0) EBrms[s] = sqrt(EBrms[s]);
112  else EBrms[s] = 1.e-06;
113  pEBmin[s] = EBmean[s] - kRMS * EBrms[s];
114  pEBmax[s] = EBmean[s] + kRMS * EBrms[s];
115  std::cout << "EB sample " << s << " mean " << EBmean[s] << " rms " << EBrms[s] << " entries " << EBtot[s] << " min " << pEBmin[s] << " max " << pEBmax[s] << std::endl;
116  // if(pEBmin[s] <= 0.) pEBmin[s] = 1.e-06;
117  vt =(double)EEtot[s];
118  EEmean[s] = EEmean[s] / vt;
119  EErms[s] = EErms[s] / vt - (EEmean[s] * EEmean[s]);
120  if(EErms[s] > 0) EErms[s] = sqrt(EErms[s]);
121  else EErms[s] = 1.e-06;
122  pEEmin[s] = EEmean[s] - kRMS * EErms[s];
123  pEEmax[s] = EEmean[s] + kRMS * EErms[s];
124  std::cout << "EE sample " << s << " mean " << EEmean[s] << " rms " << EErms[s] << " entries " << EEtot[s] << " min " << pEEmin[s] << " max " << pEEmax[s] << std::endl;
125  // if(pEEmin[s] <= 0.) pEEmin[s] = 1.e-06;
126  }
127 
128  // for(int s = 0; s < TEMPLATESAMPLES; s++)
129  // std::cout << " sample " << s << " EB min " << pEBmin[s] << " max " << pEBmax[s] << " EE min " << pEEmin[s] << " max " << pEEmax[s] << std::endl;
130  gStyle->SetPalette(1);
131  gStyle->SetOptStat(0);
132  TCanvas canvas("CC map","CC map", 1600, 5600);
133  TLatex t1;
134  t1.SetNDC();
135  t1.SetTextAlign(26);
136  t1.SetTextSize(0.05);
137  t1.DrawLatex(0.5, 0.98, Form("Ecal PulseShapes, IOV %i", run));
138 
139  float xmi[3] = {0.0 , 0.24, 0.76};
140  float xma[3] = {0.24, 0.76, 1.00};
141  TPad*** pad = new TPad**[6];
142  for (int s = 0; s < TEMPLATESAMPLES; s++) {
143  pad[s] = new TPad*[3];
144  for (int obj = 0; obj < 3; obj++) {
145  float yma = 0.96 - (0.08 * s);
146  float ymi = yma - 0.07;
147  pad[s][obj] = new TPad(Form("p_%i_%i", obj, s),Form("p_%i_%i", obj, s),
148  xmi[obj], ymi, xma[obj], yma);
149  pad[s][obj]->Draw();
150  }
151  }
152 
153  int ipad=0;
154  for(int s = 0; s < TEMPLATESAMPLES; s++) {
155  pad[ipad][0]->cd();
156  if(pEBmin[s] == pEBmax[s]) { // same values everywhere!..
157  pEBmin[s] = pEBmin[s] - 1.e-06;
158  pEBmax[s] = pEBmax[s] + 1.e-06;
159  }
160  if(pEEmin[s] == pEEmax[s]) {
161  pEEmin[s] = pEEmin[s] - 1.e-06;
162  pEEmax[s] = pEEmax[s] + 1.e-06;
163  }
164  DrawEE(endc_m[s], pEEmin[s], pEEmax[s]);
165  pad[ipad][1]->cd();
166  DrawEB(barrel[s], pEBmin[s], pEBmax[s]);
167  pad[ipad][2]->cd();
168  DrawEE(endc_p[s], pEEmin[s], pEEmax[s]);
169  ipad++;
170  }
171 
172  std::string ImageName(m_imageFileName);
173  canvas.SaveAs(ImageName.c_str());
174 
175  return true;
176  }// fill method
177  };
178 
179  /********************************************
180  profile of ECAL PulseShapes for 1 IOV
181  ********************************************/
182  class EcalPulseShapesProfile : public cond::payloadInspector::PlotImage<EcalPulseShapes> {
183 
184  public:
185  EcalPulseShapesProfile() : cond::payloadInspector::PlotImage<EcalPulseShapes>("ECAL PulseShapes - Profile ") {
186  setSingleIov(true);
187  }
188 
189  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
190  TProfile* barrel = new TProfile("EB", "EB profile", TEMPLATESAMPLES, 0, TEMPLATESAMPLES);
191  TProfile* endcap = new TProfile("EE", "EE profile", TEMPLATESAMPLES, 0, TEMPLATESAMPLES);
192 
193  auto iov = iovs.front();
194  std::shared_ptr<EcalPulseShapes> payload = fetchPayload( std::get<1>(iov) );
195  unsigned int run = std::get<0>(iov);
196  if( payload.get() ){
197  // int chan = 0;
198  for (int ieta = -MAX_IETA; ieta <= MAX_IETA; ieta++) {
199  Double_t eta = (Double_t)ieta;
200  if(ieta == 0) continue;
201  else if(ieta > 0.) eta = eta - 0.5; // 0.5 to 84.5
202  else eta = eta + 0.5; // -84.5 to -0.5
203  for (int iphi = 1; iphi <= MAX_IPHI; iphi++) {
204  //Double_t phi = (Double_t)iphi - 0.5;
205  EBDetId id(ieta, iphi);
206  EcalPulseShape pulse = (*payload)[id.rawId()];
207  for(int s = 0; s < TEMPLATESAMPLES; s++) {
208  double val = pulse.pdfval[s];
209  barrel->Fill(s, val);
210  }
211  }
212  }
213 
214  for (int sign = 0; sign < kSides; sign++) {
215  int thesign = sign==1 ? 1:-1;
216  for (int ix = 1; ix <= IX_MAX; ix++) {
217  for (int iy = 1; iy <= IY_MAX; iy++) {
218  if (! EEDetId::validDetId(ix, iy, thesign)) continue;
219  EEDetId id(ix, iy, thesign);
220  EcalPulseShape pulse = (*payload)[id.rawId()];
221  for(int s = 0; s < TEMPLATESAMPLES; s++) {
222  double val = pulse.pdfval[s];
223  endcap->Fill(s, val);
224  }
225  }// iy
226  } // ix
227  } // side
228  } // payload
229 
230 
231  gStyle->SetPalette(1);
232  gStyle->SetOptStat(0);
233  TCanvas canvas("CC map","CC map", 500, 1000);
234  TLatex t1;
235  t1.SetNDC();
236  t1.SetTextAlign(26);
237  t1.SetTextSize(0.05);
238  t1.DrawLatex(0.5, 0.98, Form("Ecal PulseShapes, IOV %i", run));
239 
240  TPad** pad = new TPad*[2];
241  for (int s = 0; s < 2; s++) {
242  float yma = 0.96 - (0.47 * s);
243  float ymi = yma - 0.45;
244  pad[s] = new TPad(Form("p_%i", s),Form("p_%i", s), 0.0, ymi, 1.0, yma);
245  pad[s]->Draw();
246  }
247 
248  pad[0]->cd();
249  barrel->Draw("");
250  pad[1]->cd();
251  endcap->Draw("");
252 
253  std::string ImageName(m_imageFileName);
254  canvas.SaveAs(ImageName.c_str());
255  return true;
256  }// fill method
257 
258  };
259 
260 } // close namespace
261 
262 // Register the classes as boost python plugin
264  PAYLOAD_INSPECTOR_CLASS(EcalPulseShapesPlot);
265  PAYLOAD_INSPECTOR_CLASS(EcalPulseShapesProfile);
266 }
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
static const int kSides
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)
T sqrt(T t)
Definition: SSEVec.h:18
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
float pdfval[TEMPLATESAMPLES]
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:481