CMS 3D CMS Logo

EcalTPGFineGrainEBIdMap_PayloadInspector.cc
Go to the documentation of this file.
8 
9 // the data format of the condition to be inspected
11 
12 #include "TH2F.h"
13 #include "TCanvas.h"
14 #include "TStyle.h"
15 #include "TLine.h"
16 #include "TLatex.h"
17 
18 #include <string>
19 #include <cstring>
20 
21 namespace {
22  enum {TEMPLATESAMPLES=5};
23  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 17, MAX_IPHI = 72}; // barrel lower and upper bounds on eta and phi
24 
25  /***********************************************
26  2d plot of EcalTPGFineGrainEBIdMap of 1 IOV
27  ************************************************/
28  class EcalTPGFineGrainEBIdMapPlot : public cond::payloadInspector::PlotImage<EcalTPGFineGrainEBIdMap> {
29 
30  public:
31  EcalTPGFineGrainEBIdMapPlot() : cond::payloadInspector::PlotImage<EcalTPGFineGrainEBIdMap>("Ecal TPGFineGrainEBIdMap - map ") {
32  setSingleIov(true);
33  }
34 
35  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
36  TH2F** barrel = new TH2F*[TEMPLATESAMPLES];
37  double pEBmin[TEMPLATESAMPLES], pEBmax[TEMPLATESAMPLES];
38  std::string text[TEMPLATESAMPLES]= {"ThresholdETLow","ThresholdETHigh","RatioLow","RatioHigh","LUT"};
39  int EBcnt = 0;
40 
41  for(int s = 0; s < TEMPLATESAMPLES; ++s) {
42  char *y = new char[text[s].length() + 1];
43  std::strcpy(y, text[s].c_str());
44  barrel[s]=new TH2F("EB",y, MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
45 
46  }
47 
48  uint32_t ThresholdETLow =0;
49  uint32_t ThresholdETHigh =0;
50  uint32_t RatioLow =0;
51  uint32_t RatioHigh =0;
52  uint32_t LUT =0;
53 
54  auto iov = iovs.front();
55  std::shared_ptr<EcalTPGFineGrainEBIdMap> payload = fetchPayload( std::get<1>(iov) );
56  unsigned int run = std::get<0>(iov);
57  if( payload.get() ){
58 
59  const std::map<uint32_t, EcalTPGFineGrainConstEB> &towerMap = (*payload).getMap();
60  std::map<uint32_t, EcalTPGFineGrainConstEB>::const_iterator it = towerMap.begin();
61 
62 
63  for(int iphi = MIN_IPHI-1; iphi <= MAX_IPHI; iphi++){
64  for(int ieta= -1*MAX_IETA; ieta < MAX_IETA; ieta++){
65 
66 
67 
68  //if(ieta > 0) ieta--;
69 
70 
71  EcalTPGFineGrainConstEB fg=(*it).second;
72 
73  fg.getValues(ThresholdETLow,ThresholdETHigh,RatioLow,RatioHigh,LUT);
74 
75  barrel[0]->Fill(iphi, ieta,ThresholdETLow);
76  barrel[1]->Fill(iphi, ieta,ThresholdETHigh);
77  barrel[2]->Fill(iphi, ieta,RatioLow);
78  barrel[3]->Fill(iphi, ieta,RatioHigh);
79  barrel[4]->Fill(iphi, ieta,LUT);
80 
81 
82  if(ThresholdETLow<pEBmin[0])pEBmin[0]=ThresholdETLow;
83  if(ThresholdETHigh<pEBmin[1])pEBmin[1]=ThresholdETHigh;
84  if(RatioLow<pEBmin[2])pEBmin[2]=RatioLow;
85  if(RatioHigh<pEBmin[3])pEBmin[3]=RatioHigh;
86  if(LUT<pEBmin[4])pEBmin[4]=LUT;
87 
88  if(ThresholdETLow>pEBmax[0])pEBmax[0]=ThresholdETLow;
89  if(ThresholdETHigh>pEBmax[1])pEBmax[1]=ThresholdETHigh;
90  if(RatioLow>pEBmax[2])pEBmax[2]=RatioLow;
91  if(RatioHigh>pEBmax[3])pEBmax[3]=RatioHigh;
92  if(LUT>pEBmax[4])pEBmax[4]=LUT;
93 
94 
95 
96 
97  EBcnt++;
98  }
99  }
100 
101 
102 
103 
104  } // payload
105 
106  gStyle->SetPalette(1);
107  gStyle->SetOptStat(0);
108  TCanvas canvas("CC map","CC map", 1600, 2800);
109  TLatex t1;
110  t1.SetNDC();
111  t1.SetTextAlign(26);
112  t1.SetTextSize(0.04);
113  t1.DrawLatex(0.5, 0.96, Form("Ecal TPGFine Grain EBIdMap, IOV %i", run));
114 
115  float xmi = 0.24;
116  float xma = 0.76;
117 
118  TPad** pad = new TPad*[TEMPLATESAMPLES];
119  for (int s = 0; s < TEMPLATESAMPLES; s++) {
120  float yma = 0.94 - (0.16 * s);
121  float ymi = yma - 0.14;
122  //pad[s] = new TPad(text[s],Form("Towers %i", EBcnt),xmi, ymi, xma, yma);
123  char *y = new char[text[s].length() + 1];
124  std::strcpy(y, text[s].c_str());
125 
126  pad[s] = new TPad(Form("Towers %i", EBcnt),y,xmi, ymi, xma, yma);
127  pad[s]->Draw();
128  }
129 
130 
131  for(int s = 0; s < TEMPLATESAMPLES; s++) {
132  pad[s]->cd();
133 
134  if(pEBmin[s] == pEBmax[s]) { // same values everywhere!..
135  pEBmin[s] = pEBmin[s] - 1.e-06;
136  pEBmax[s] = pEBmax[s] + 1.e-06;
137  }
138 
139 
140  barrel[s]->SetMaximum(pEBmax[s]);
141  barrel[s]->SetMinimum(pEBmin[s]);
142  barrel[s]->Draw("colz");
143 
144  TLine* l = new TLine(0., 0., 0., 0.);
145  l->SetLineWidth(1);
146  for(int i = 0; i < MAX_IETA; i++) {
147  Double_t x = 4.+ (i * 4);
148  l = new TLine(x, -MAX_IETA, x, MAX_IETA);
149  l->Draw();
150  }
151 
152  }
153 
154 
155  std::string ImageName(m_imageFileName);
156  canvas.SaveAs(ImageName.c_str());
157 
158  return true;
159  }// fill method
160  };
161 
162 }
163 
164 
165 // Register the classes as boost python plugin
167  PAYLOAD_INSPECTOR_CLASS(EcalTPGFineGrainEBIdMapPlot);
168 }
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
virtual bool fill(const std::vector< std::tuple< cond::Time_t, cond::Hash > > &iovs)=0
void getValues(uint32_t &ThresholdETLow, uint32_t &ThresholdETHigh, uint32_t &RatioLow, uint32_t &RatioHigh, uint32_t &LUT) const
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:482
std::vector< unsigned short int > LUT
Definition: DTTracoLUTs.h:31