CMS 3D CMS Logo

EcalTPGWeightGroup_PayloadInspector.cc
Go to the documentation of this file.
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.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 
21  enum {kEBTotalTowers = 2448, kEETotalTowers = 1584};
22  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 17, MAX_IPHI = 72}; // barrel lower and upper bounds on eta and phi
23  enum {IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100, EEhistXMax = 220}; // endcaps lower and upper bounds on x and y
24 /***********************************************
25  2d plot of EcalTPGWeightGroup of 1 IOV
26 ************************************************/
27  class EcalTPGWeightGroupPlot : public cond::payloadInspector::PlotImage<EcalTPGWeightGroup> {
28 
29  public:
30  EcalTPGWeightGroupPlot() : cond::payloadInspector::PlotImage<EcalTPGWeightGroup>("EcalTPGWeightGroup - map ") {
31  setSingleIov(true);
32  }
33 
34  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
35 
36  uint32_t minEB=0;
37  uint32_t maxEB=2;
38  uint32_t minEE=0;
39  uint32_t maxEE=2;
40 
41  TH2F* barrel = new TH2F("EB","EB Tower Status", 72, 0, 72, 34, -17, 17);
42  TH2F* endc_p = new TH2F("EE+","EE+ Tower Status",22, 0, 22, 22, 0, 22);
43  TH2F* endc_m = new TH2F("EE-","EE- Tower Status",22, 0, 22, 22, 0, 22);
44 
45  auto iov = iovs.front();
46  std::shared_ptr<EcalTPGWeightGroup> payload = fetchPayload( std::get<1>(iov) );
47  unsigned int run = std::get<0>(iov);
48 
49  if( payload.get() ){
50 
51 
52  const EcalTPGWeightGroup::EcalTPGGroupsMap & map=(*payload).getMap();
54 
55  for(it = map.begin() ; it != map.end() ; it++){
57  EcalTrigTowerDetId ttId_eb(rawid);
58 
59  //std::cout<<(*it).second << std::endl;
60 
61  if(ttId_eb.subDet()==1){
62  //barrel
63  int ieta = ttId_eb.ieta();
64  if(ieta < 0) ieta--; // -1 to -17
65  int iphi = ttId_eb.iphi() - 1; // 0 to 71
66 
67  if(minEB > (*it).second)
68  minEB = (*it).second;
69 
70  if(maxEB < (*it).second)
71  maxEB = (*it).second;
72 
73  barrel->Fill(iphi, ieta,(*it).second);
74  }
75 
76 
77  if(EcalScDetId::validHashIndex((*it).first)){
78  EcalScDetId rawid_ee = EcalScDetId::unhashIndex((*it).first);
79  EcalScDetId ttId_ee(rawid_ee);
80 
81  //endcaps
82 
83  int ix = ttId_ee.ix();
84  int iy = ttId_ee.iy();
85  int zside = ttId_ee.zside();
86 
87 
88  if(minEE > (*it).second)
89  minEE = (*it).second;
90 
91  if(maxEE < (*it).second)
92  maxEE = (*it).second;
93 
94  if(zside == 1){
95  endc_p->Fill(ix, iy, (*it).second);
96  }else{
97  endc_m->Fill(ix, iy, (*it).second);
98  }
99 
100  }
101 
102  }
103 
104 
105  } // payload
106 
107 
108  TCanvas canvas("CC map","CC map",800,800);
109  TLatex t1;
110  t1.SetNDC();
111  t1.SetTextAlign(26);
112  t1.SetTextSize(0.05);
113  t1.DrawLatex(0.5, 0.96, Form("Ecal TPG WeightGroup, IOV %i", run));
114 
115  //TPad* padb = new TPad("padb","padb", 0., 0.55, 1., 1.);
116  TPad* padb = new TPad("padb","padb", 0., 0.45, 1., 0.9);
117  padb->Draw();
118  TPad* padem = new TPad("padem","padem", 0., 0., 0.45, 0.45);
119  padem->Draw();
120  TPad* padep = new TPad("padep","padep", 0.55, 0., 1., 0.45);
121  padep->Draw();
122 
123  TLine* l = new TLine(0., 0., 0., 0.);
124  l->SetLineWidth(1);
125  padb->cd();
126  barrel->SetStats(false);
127  barrel->SetMaximum(maxEB);
128  barrel->SetMinimum(minEB);
129  barrel->Draw("colz");
130 
131  //barrel->Draw("col");
132  for(int i = 0; i <17; i++) {
133  Double_t x = 4.+ (i * 4);
134  l = new TLine(x, -17., x, 17.);
135  l->Draw();
136  }
137 
138  l = new TLine(0., 0., 72., 0.);
139  l->Draw();
140 
141 
142  padem->cd();
143  DrawEE_Tower(endc_m,l,minEE,maxEE);
144 
145 
146  padep->cd();
147  DrawEE_Tower(endc_p,l,minEE,maxEE);
148 
149 
150  std::string ImageName(m_imageFileName);
151  canvas.SaveAs(ImageName.c_str());
152  return true;
153  }// fill method
154 
155 };
156 
157 
158 }
159 
160 
161 // Register the classes as boost python plugin
163  PAYLOAD_INSPECTOR_CLASS(EcalTPGWeightGroupPlot);
164 }
static EcalTrigTowerDetId detIdFromDenseIndex(uint32_t di)
std::map< uint32_t, uint32_t >::const_iterator EcalTPGGroupsMapItr
Definition: EcalTPGGroups.h:23
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
static EcalScDetId unhashIndex(int hi)
Definition: EcalScDetId.h:118
int zside(DetId const &)
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
virtual bool fill(const std::vector< std::tuple< cond::Time_t, cond::Hash > > &iovs)=0
void DrawEE_Tower(TH2F *endc, TLine *l, double minScale, double maxScale)
Definition: EcalDrawUtils.h:71
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
static bool validHashIndex(int hi)
Definition: EcalScDetId.h:140
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:481
std::map< uint32_t, uint32_t > EcalTPGGroupsMap
Definition: EcalTPGGroups.h:22