CMS 3D CMS Logo

EcalFloatCondObjectContainer_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, kRMS = 5};
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 FloatCondObjectContainer of 1 IOV
25  *****************************************************/
26  class EcalFloatCondObjectContainerPlot : public cond::payloadInspector::PlotImage<EcalFloatCondObjectContainer> {
27 
28  public:
29  EcalFloatCondObjectContainerPlot() : cond::payloadInspector::PlotImage<EcalFloatCondObjectContainer>("ECAL FloatCondObjectContainer - 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", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
35  TH2F* endc_p = new TH2F("EE+","EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
36  TH2F* endc_m = new TH2F("EE-","EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
37  double EBmean = 0., EBrms = 0., EEmean = 0., EErms = 0.;
38  int EBtot = 0, EEtot = 0;
39 
40  auto iov = iovs.front();
41  std::shared_ptr<EcalFloatCondObjectContainer> payload = fetchPayload( std::get<1>(iov) );
42  unsigned int run = std::get<0>(iov);
43  if( payload.get() ){
44  for (int ieta = -MAX_IETA; ieta <= MAX_IETA; ieta++) {
45  Double_t eta = (Double_t)ieta;
46  if(ieta == 0) continue;
47  else if(ieta > 0.) eta = eta - 0.5; // 0.5 to 84.5
48  else eta = eta + 0.5; // -84.5 to -0.5
49  for (int iphi = 1; iphi <= MAX_IPHI; iphi++) {
50  Double_t phi = (Double_t)iphi - 0.5;
51  EBDetId id(ieta, iphi);
52  double val = (*payload)[id.rawId()];
53  barrel->Fill(phi, eta, val);
54  EBmean = EBmean + val;
55  EBrms = EBrms + val * val;
56  EBtot++;
57  }
58  }
59 
60  for (int sign = 0; sign < kSides; sign++) {
61  int thesign = sign==1 ? 1:-1;
62  for (int ix = 1; ix <= IX_MAX; ix++) {
63  for (int iy = 1; iy <= IY_MAX; iy++) {
64  if (! EEDetId::validDetId(ix, iy, thesign)) continue;
65  EEDetId id(ix, iy, thesign);
66  double val = (*payload)[id.rawId()];
67  EEmean = EEmean + val;
68  EErms = EErms + val * val;
69  EEtot++;
70  if (thesign==1)
71  endc_p->Fill(ix, iy, val);
72  else
73  endc_m->Fill(ix, iy, val);
74  }// iy
75  } // ix
76  } // side
77  } // payload
78  double vt =(double)EBtot;
79  EBmean = EBmean / vt;
80  EBrms = EBrms / vt - (EBmean * EBmean);
81  EBrms = sqrt(EBrms);
82  if(EBrms == 0.) EBrms = 0.001;
83  double pEBmin = EBmean - kRMS * EBrms;
84  double pEBmax = EBmean + kRMS * EBrms;
85  // std::cout << " mean " << EBmean << " rms " << EBrms << " entries " << EBtot << " min " << pEBmin << " max " << pEBmax << std::endl;
86  vt =(double)EEtot;
87  EEmean = EEmean / vt;
88  EErms = EErms / vt - (EEmean * EEmean);
89  EErms = sqrt(EErms);
90  if(EErms == 0.) EErms = 0.001;
91  double pEEmin = EEmean - kRMS * EErms;
92  double pEEmax = EEmean + kRMS * EErms;
93  // std::cout << " mean " << EEmean << " rms " << EErms << " entries " << EEtot << " min " << pEEmin << " max " << pEEmax << std::endl;
94 
95  gStyle->SetPalette(1);
96  gStyle->SetOptStat(0);
97  TCanvas canvas("CC map","CC map", 1600, 450);
98  TLatex t1;
99  t1.SetNDC();
100  t1.SetTextAlign(26);
101  t1.SetTextSize(0.05);
102  t1.DrawLatex(0.5, 0.96, Form("Ecal FloatCondObjectContainer, IOV %i", run));
103 
104  float xmi[3] = {0.0 , 0.24, 0.76};
105  float xma[3] = {0.24, 0.76, 1.00};
106  TPad** pad = new TPad*;
107  for (int obj = 0; obj < 3; obj++) {
108  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
109  pad[obj]->Draw();
110  }
111 
112  pad[0]->cd();
113  DrawEE(endc_m, pEEmin, pEEmax);
114  pad[1]->cd();
115  DrawEB(barrel, pEBmin, pEBmax);
116  pad[2]->cd();
117  DrawEE(endc_p, pEEmin, pEEmax);
118 
119  std::string ImageName(m_imageFileName);
120  canvas.SaveAs(ImageName.c_str());
121  return true;
122  }// fill method
123  };
124 
125  /************************************************************************
126  2d plot of ECAL FloatCondObjectContainer difference between 2 IOVs
127  ************************************************************************/
128  class EcalFloatCondObjectContainerDiff : public cond::payloadInspector::PlotImage<EcalFloatCondObjectContainer> {
129 
130  public:
131  EcalFloatCondObjectContainerDiff() : cond::payloadInspector::PlotImage<EcalFloatCondObjectContainer>("ECAL FloatCondObjectContainer difference") {
132  setSingleIov(false);
133  }
134 
135  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
136  TH2F* barrel = new TH2F("EB","EB difference", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
137  TH2F* endc_p = new TH2F("EE+","EE+ difference", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
138  TH2F* endc_m = new TH2F("EE-","EE- difference", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
139  double EBmean = 0., EBrms = 0., EEmean = 0., EErms = 0.;
140  int EBtot = 0, EEtot = 0;
141 
142  // unsigned int run[2] = {0, 0}, irun = 0;
143  unsigned int run[2], irun = 0;
144  float vEB[kEBChannels], vEE[kEEChannels];
145  for ( auto const & iov: iovs) {
146  std::shared_ptr<EcalFloatCondObjectContainer> payload = fetchPayload( std::get<1>(iov) );
147  run[irun] = std::get<0>(iov);
148  if( payload.get() ){
149  for (int ieta = -MAX_IETA; ieta <= MAX_IETA; ieta++) {
150  Double_t eta = (Double_t)ieta;
151  if(ieta == 0) continue;
152  else if(ieta > 0.) eta = eta - 0.5; // 0.5 to 84.5
153  else eta = eta + 0.5; // -84.5 to -0.5
154  for (int iphi = 1; iphi <= MAX_IPHI; iphi++) {
155  Double_t phi = (Double_t)iphi - 0.5;
156  EBDetId id(ieta, iphi);
157  int channel = id.hashedIndex();
158  double val = (*payload)[id.rawId()];
159  if(irun == 0) vEB[channel] = val;
160  else {
161  double diff = val - vEB[channel];
162  barrel->Fill(phi, eta, diff);
163  EBmean = EBmean + diff;
164  EBrms = EBrms + diff * diff;
165  EBtot++;
166  // std::cout << " entry " << EBtot << " mean " << EBmean << " rms " << EBrms << std::endl;
167  }
168  }
169  }
170 
171  for (int sign = 0; sign < kSides; sign++) {
172  int thesign = sign==1 ? 1:-1;
173  for (int ix = 1; ix <= IX_MAX; ix++) {
174  for (int iy = 1; iy <= IY_MAX; iy++) {
175  if (! EEDetId::validDetId(ix, iy, thesign)) continue;
176  EEDetId id(ix, iy, thesign);
177  int channel = id.hashedIndex();
178  double val = (*payload)[id.rawId()];
179  if(irun == 0) vEE[channel] = val;
180  else {
181  double diff = val - vEE[channel];
182  EEmean = EEmean + diff;
183  EErms = EErms + diff * diff;
184  EEtot++;
185  if (thesign==1)
186  endc_p->Fill(ix, iy, diff);
187  else
188  endc_m->Fill(ix, iy, diff);
189  }
190  }// iy
191  } // ix
192  } // side
193  } // payload
194  else return false;
195  irun++;
196  } // loop over IOVs
197 
198  double vt =(double)EBtot;
199  EBmean = EBmean / vt;
200  EBrms = EBrms / vt - (EBmean * EBmean);
201  EBrms = sqrt(EBrms);
202  if(EBrms == 0.) EBrms = 0.001;
203  double pEBmin = EBmean - kRMS * EBrms;
204  double pEBmax = EBmean + kRMS * EBrms;
205  // std::cout << " mean " << EBmean << " rms " << EBrms << " entries " << EBtot << " min " << pEBmin << " max " << pEBmax << std::endl;
206  vt =(double)EEtot;
207  EEmean = EEmean / vt;
208  EErms = EErms / vt - (EEmean * EEmean);
209  EErms = sqrt(EErms);
210  if(EErms == 0.) EErms = 0.001;
211  double pEEmin = EEmean - kRMS * EErms;
212  double pEEmax = EEmean + kRMS * EErms;
213  // std::cout << " mean " << EEmean << " rms " << EErms << " entries " << EEtot << " min " << pEEmin << " max " << pEEmax << std::endl;
214 
215  gStyle->SetPalette(1);
216  gStyle->SetOptStat(0);
217  TCanvas canvas("CC map","CC map", 1600, 450);
218  TLatex t1;
219  t1.SetNDC();
220  t1.SetTextAlign(26);
221  t1.SetTextSize(0.05);
222  t1.DrawLatex(0.5, 0.96, Form("Ecal FloatCondObjectContainer, IOV %i - %i", run[1], run[0]));
223 
224  float xmi[3] = {0.0 , 0.24, 0.76};
225  float xma[3] = {0.24, 0.76, 1.00};
226  TPad** pad = new TPad*;
227  for (int obj = 0; obj < 3; obj++) {
228  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
229  pad[obj]->Draw();
230  }
231 
232  pad[0]->cd();
233  DrawEE(endc_m, pEEmin, pEEmax);
234  pad[1]->cd();
235  DrawEB(barrel, pEBmin, pEBmax);
236  pad[2]->cd();
237  DrawEE(endc_p, pEEmin, pEEmax);
238 
239  std::string ImageName(m_imageFileName);
240  canvas.SaveAs(ImageName.c_str());
241  return true;
242  }// fill method
243  };
244 
245 } // close namespace
246 
247 // Register the classes as boost python plugin
249  PAYLOAD_INSPECTOR_CLASS(EcalFloatCondObjectContainerPlot);
250  PAYLOAD_INSPECTOR_CLASS(EcalFloatCondObjectContainerDiff);
251 }
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
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:482