CMS 3D CMS Logo

EcalTPGCrystalStatus_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 "TCanvas.h"
12 #include "TStyle.h"
13 #include "TLine.h"
14 #include "TLatex.h"
15 
16 #include <string>
17 
18 namespace {
19  enum {kEBChannels = 61200, kEEChannels = 14648, kSides = 2};
20  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360}; // barrel lower and upper bounds on eta and phi
21  enum {IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100}; // endcaps lower and upper bounds on x and y
22 
23  /***********************************************
24  2d plot of ECAL TPGCrystalStatus of 1 IOV
25  ************************************************/
26  class EcalTPGCrystalStatusPlot : public cond::payloadInspector::PlotImage<EcalTPGCrystalStatus> {
27 
28  public:
29  EcalTPGCrystalStatusPlot() : cond::payloadInspector::PlotImage<EcalTPGCrystalStatus>("ECAL TPGCrystalStatus - map ") {
30  setSingleIov(true);
31  }
32 
33  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
34  TH2F* barrel = new TH2F("EB","EB TPG Crystal Status", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
35  TH2F* endc_p = new TH2F("EE+","EE+ TPG Crystal Status", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
36  TH2F* endc_m = new TH2F("EE-","EE- TPG Crystal Status", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
37  int EBstat = 0, EEstat[2] = {0, 0};
38 
39  auto iov = iovs.front();
40  std::shared_ptr<EcalTPGCrystalStatus> payload = fetchPayload( std::get<1>(iov) );
41  unsigned int run = std::get<0>(iov);
42  if( payload.get() ){
43  for (int ieta = -MAX_IETA; ieta <= MAX_IETA; ieta++) {
44  Double_t eta = (Double_t)ieta;
45  if(ieta == 0) continue;
46  else if(ieta > 0.) eta = eta - 0.5; // 0.5 to 84.5
47  else eta = eta + 0.5; // -84.5 to -0.5
48  for (int iphi = 1; iphi <= MAX_IPHI; iphi++) {
49  Double_t phi = (Double_t)iphi - 0.5;
50  EBDetId id(ieta, iphi);
51  double val = (*payload)[id.rawId()].getStatusCode();
52  barrel->Fill(phi, eta, val);
53  if(val > 0) EBstat++;
54  }
55  }
56 
57  for (int sign = 0; sign < kSides; sign++) {
58  int thesign = sign==1 ? 1:-1;
59  for (int ix = 1; ix <= IX_MAX; ix++) {
60  for (int iy = 1; iy <= IY_MAX; iy++) {
61  if (! EEDetId::validDetId(ix, iy, thesign)) continue;
62  EEDetId id(ix, iy, thesign);
63  double val = (*payload)[id.rawId()].getStatusCode();
64  if (thesign==1) {
65  endc_p->Fill(ix, iy, val);
66  if(val > 0) EEstat[1]++;
67  }
68  else {
69  endc_m->Fill(ix, iy, val);
70  if(val > 0) EEstat[0]++;
71  }
72  }// iy
73  } // ix
74  } // side
75  } // payload
76 
77  gStyle->SetPalette(1);
78  gStyle->SetOptStat(0);
79  // TCanvas canvas("CC map","CC map", 1600, 450);
80  Double_t w = 1200;
81  Double_t h = 1400;
82  TCanvas canvas("c", "c", w, h);
83  canvas.SetWindowSize(w + (w - canvas.GetWw()), h + (h - canvas.GetWh()));
84 
85  TLatex t1;
86  t1.SetNDC();
87  t1.SetTextAlign(26);
88  t1.SetTextSize(0.05);
89  t1.DrawLatex(0.5, 0.96, Form("Ecal TPGCrystalStatus, IOV %i", run));
90 
91  // float xmi[3] = {0.0 , 0.24, 0.76};
92  // float xma[3] = {0.24, 0.76, 1.00};
93  float xmi[3] = {0.0, 0.0, 0.5};
94  float xma[3] = {1.0, 0.5, 1.0};
95  float ymi[3] = {0.47, 0.0, 0.0};
96  float yma[3] = {0.94, 0.47, 0.47};
97  TPad** pad = new TPad*;
98  for (int obj = 0; obj < 3; obj++) {
99  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], ymi[obj], xma[obj], yma[obj]);
100  pad[obj]->Draw();
101  }
102 
103  pad[0]->cd();
104  DrawEB(barrel, 0., 1.);
105  t1.DrawLatex(0.2, 0.94, Form("%i crystals", EBstat));
106  pad[1]->cd();
107  DrawEE(endc_m, 0., 1.);
108  t1.DrawLatex(0.15, 0.92, Form("%i crystals", EEstat[0]));
109  pad[2]->cd();
110  DrawEE(endc_p, 0., 1.);
111  t1.DrawLatex(0.15, 0.92, Form("%i crystals", EEstat[1]));
112 
113  std::string ImageName(m_imageFileName);
114  canvas.SaveAs(ImageName.c_str());
115  return true;
116  }// fill method
117  };
118 
119  /************************************************************************
120  2d plot of ECAL TPGCrystalStatus difference between 2 IOVs
121  ************************************************************************/
122  class EcalTPGCrystalStatusDiff : public cond::payloadInspector::PlotImage<EcalTPGCrystalStatus> {
123 
124  public:
125  EcalTPGCrystalStatusDiff() : cond::payloadInspector::PlotImage<EcalTPGCrystalStatus>("ECAL TPGCrystalStatus difference") {
126  setSingleIov(false);
127  }
128 
129  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
130  TH2F* barrel = new TH2F("EB","EB difference", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
131  TH2F* endc_p = new TH2F("EE+","EE+ difference", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
132  TH2F* endc_m = new TH2F("EE-","EE- difference", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
133  int EBstat = 0, EEstat[2] = {0, 0};
134 
135  unsigned int run[2] = {0, 0}, irun = 0;
136  float vEB[kEBChannels], vEE[kEEChannels];
137  for ( auto const & iov: iovs) {
138  std::shared_ptr<EcalTPGCrystalStatus> payload = fetchPayload( std::get<1>(iov) );
139  run[irun] = std::get<0>(iov);
140  if( payload.get() ){
141  for (int ieta = -MAX_IETA; ieta <= MAX_IETA; ieta++) {
142  Double_t eta = (Double_t)ieta;
143  if(ieta == 0) continue;
144  else if(ieta > 0.) eta = eta - 0.5; // 0.5 to 84.5
145  else eta = eta + 0.5; // -84.5 to -0.5
146  for (int iphi = 1; iphi <= MAX_IPHI; iphi++) {
147  Double_t phi = (Double_t)iphi - 0.5;
148  EBDetId id(ieta, iphi);
149  int channel = id.hashedIndex();
150  double val = (*payload)[id.rawId()].getStatusCode();
151  if(irun == 0) vEB[channel] = val;
152  else {
153  double diff = val - vEB[channel];
154  barrel->Fill(phi, eta, diff);
155  if(diff != 0) EBstat++;
156  // std::cout << " entry " << EBtot << " mean " << EBmean << " rms " << EBrms << std::endl;
157  }
158  }
159  }
160 
161  for (int sign = 0; sign < kSides; sign++) {
162  int thesign = sign==1 ? 1:-1;
163  for (int ix = 1; ix <= IX_MAX; ix++) {
164  for (int iy = 1; iy <= IY_MAX; iy++) {
165  if (! EEDetId::validDetId(ix, iy, thesign)) continue;
166  EEDetId id(ix, iy, thesign);
167  int channel = id.hashedIndex();
168  double val = (*payload)[id.rawId()].getStatusCode();
169  if(irun == 0) vEE[channel] = val;
170  else {
171  double diff = val - vEE[channel];
172  if (thesign==1) {
173  endc_p->Fill(ix, iy, diff);
174  if(diff != 0) EEstat[1]++;
175  }
176  else {
177  endc_m->Fill(ix, iy, diff);
178  if(diff != 0) EEstat[0]++;
179  }
180  }
181  }// iy
182  } // ix
183  } // side
184  } // payload
185  else return false;
186  irun++;
187  } // loop over IOVs
188 
189  gStyle->SetPalette(1);
190  gStyle->SetOptStat(0);
191  Double_t w = 1200;
192  Double_t h = 1400;
193  TCanvas canvas("c", "c", w, h);
194  canvas.SetWindowSize(w + (w - canvas.GetWw()), h + (h - canvas.GetWh()));
195 
196  TLatex t1;
197  t1.SetNDC();
198  t1.SetTextAlign(26);
199  t1.SetTextSize(0.05);
200  t1.DrawLatex(0.5, 0.96, Form("Ecal TPGCrystalStatus, IOV %i - %i", run[1], run[0]));
201 
202  // float xmi[3] = {0.0 , 0.24, 0.76};
203  // float xma[3] = {0.24, 0.76, 1.00};
204  float xmi[3] = {0.0, 0.0, 0.5};
205  float xma[3] = {1.0, 0.5, 1.0};
206  float ymi[3] = {0.47, 0.0, 0.0};
207  float yma[3] = {0.94, 0.47, 0.47};
208  TPad** pad = new TPad*;
209  for (int obj = 0; obj < 3; obj++) {
210  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], ymi[obj], xma[obj], yma[obj]);
211  pad[obj]->Draw();
212  }
213 
214  pad[0]->cd();
215  DrawEB(barrel, -1., 1.);
216  t1.DrawLatex(0.2, 0.94, Form("%i differences", EBstat));
217  pad[1]->cd();
218  DrawEE(endc_m, -1., 1.);
219  t1.DrawLatex(0.15, 0.92, Form("%i differences", EEstat[0]));
220  pad[2]->cd();
221  DrawEE(endc_p, -1., 1.);
222  t1.DrawLatex(0.15, 0.92, Form("%i differences", EEstat[1]));
223 
224  std::string ImageName(m_imageFileName);
225  canvas.SaveAs(ImageName.c_str());
226  return true;
227  }// fill method
228  };
229 
230 } // close namespace
231 
232 // Register the classes as boost python plugin
234  PAYLOAD_INSPECTOR_CLASS(EcalTPGCrystalStatusPlot);
235  PAYLOAD_INSPECTOR_CLASS(EcalTPGCrystalStatusDiff);
236 }
const double w
Definition: UKUtility.cc:23
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)
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
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:481