CMS 3D CMS Logo

EcalTPGFineGrainEBGroup_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 
20 namespace {
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 
24  /***********************************************
25  2d plot of EcalTPGFineGrainEBGroup of 1 IOV
26  ************************************************/
27  class EcalTPGFineGrainEBGroupPlot : public cond::payloadInspector::PlotImage<EcalTPGFineGrainEBGroup> {
28 
29  public:
30  EcalTPGFineGrainEBGroupPlot() : cond::payloadInspector::PlotImage<EcalTPGFineGrainEBGroup>("EcalTPGFineGrainEBGroup - 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("EB","Ecal TPGFineGrain EB Group", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
36  int EBcount = 0;
37  double minEB=0, maxEB=1;
38 
39 
40  auto iov = iovs.front();
41  std::shared_ptr<EcalTPGFineGrainEBGroup> payload = fetchPayload( std::get<1>(iov) );
42  unsigned int run = std::get<0>(iov);
43  if( payload.get() ){
44 
45  const EcalTPGFineGrainEBGroup::EcalTPGGroupsMap &towerMap = (*payload).getMap();
46 
48  for(it = towerMap.begin(); it != towerMap.end(); ++it) {
49  //EcalTrigTowerDetId ttId = EcalTrigTowerDetId::detIdFromDenseIndex((*it).first);
50  EcalTrigTowerDetId ttId((*it).first);
51  int ieta = ttId.ieta();
52  //ieta--;
53  if(ieta > 0) ieta--; // -1 to -17
54  int iphi = ttId.iphi() - 1; // 0 to 71
55  // std::cout << " sub det " << ttId.subDet() << " phi " << iphi << " eta " << ieta << std::endl;
56  // ieta goes from -18 to -2 and 1 to 17. Change it to -17/-1 and 0/16
57 
58  //std::cout <<(*it).first<<std::endl;
59  //std::cout << " ieta " << ieta << " phi " << iphi << " value " << (*it).second << std::endl;
60 
61 
62  if(ttId.subDet() == 1) { // barrel
63 
64  barrel->Fill(iphi, ieta,(*it).second);
65 
66  if(maxEB<(*it).second)
67  maxEB=(*it).second;
68  if(minEB>(*it).second)
69  minEB=(*it).second;
70 
71  EBcount++;
72  }
73 
74  }
75  } // payload
76 
77  gStyle->SetPalette(1);
78  gStyle->SetOptStat(0);
79  // TCanvas canvas("CC map","CC map", 1600, 450);
80  Double_t w = 1400;
81  Double_t h = 1200;
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 TPGFine GrainEBGroup, IOV %i", run));
90 
91  TPad** pad = new TPad*;
92  for (int obj = 0; obj < 1; obj++) {
93  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), 0.0, 0.04, 1.0, 0.94);
94  pad[obj]->Draw();
95  }
96  t1.SetTextSize(0.03);
97  t1.DrawLatex(0.2, 0.88, Form("%i towers", EBcount));
98  // canvas.cd();
99  pad[0]->cd();
100  //barrel->SetStats(0);
101  barrel->SetMinimum(minEB);
102  barrel->SetMaximum(maxEB);
103  barrel->Draw("colz");
104  TLine* l = new TLine(0., 0., 0., 0.);
105  l->SetLineWidth(1);
106  for(int i = 0; i < MAX_IETA; i++) {
107  Double_t x = 4.+ (i * 4);
108  l = new TLine(x, -MAX_IETA, x, MAX_IETA);
109  l->Draw();
110  }
111  l = new TLine(0., 0., 72., 0.);
112  l->Draw();
113 
114  std::string ImageName(m_imageFileName);
115  canvas.SaveAs(ImageName.c_str());
116  return true;
117  }// fill method
118  };
119 
120 }
121 
122 
123 // Register the classes as boost python plugin
125  PAYLOAD_INSPECTOR_CLASS(EcalTPGFineGrainEBGroupPlot);
126 }
const double w
Definition: UKUtility.cc:23
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
std::map< uint32_t, uint32_t >::const_iterator EcalTPGGroupsMapItr
Definition: EcalTPGGroups.h:23
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
unsigned ttId(DetId const &)
#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)
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:481
std::map< uint32_t, uint32_t > EcalTPGGroupsMap
Definition: EcalTPGGroups.h:22