CMS 3D CMS Logo

EcalTPGPhysicsConst_PayloadInspector.cc
Go to the documentation of this file.
7 
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.h" // a 2-D histogram with four bytes per cell (float)
12 #include "TCanvas.h"
13 #include "TLine.h"
14 #include "TStyle.h"
15 #include "TLatex.h"//write mathematical equations.
16 #include "TPave.h"
17 #include "TPaveStats.h"
18 #include <string>
19 #include <fstream>
20 
21 namespace {
22 /*****************************************
23  2d plot of Ecal TPG Physics Const of 1 IOV
24  ******************************************/
25 class EcalTPGPhysicsConstPlot: public cond::payloadInspector::PlotImage<EcalTPGPhysicsConst>{
26  public:
27  EcalTPGPhysicsConstPlot():
28  cond::payloadInspector::PlotImage<EcalTPGPhysicsConst>("ECAL TPG Physics Constant - map "){
29  setSingleIov(true);
30  }
31 
32  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
33  auto iov=iovs.front();
34  std::shared_ptr <EcalTPGPhysicsConst> payload = fetchPayload(std::get<1> (iov));
35  unsigned int run=std::get<0> (iov);
36  TH2F* align;
37  int NbRows;
38 
39  if(payload.get()){
40  EcalTPGPhysicsConstMap map=(*payload).getMap();
41  NbRows=map.size();
42 
43 
44  align =new TH2F("TPGPhysicsConstant","mapKey EtSat ttf_threshold_Low ttf_threshold_High FG_lowThreshold FG_highThreshold FG_lowRatio FG_highRatio",
45  8, 0, 8, NbRows, 0, NbRows);
46 
47  double row = NbRows - 0.5;
48  for (std::map<uint32_t, EcalTPGPhysicsConst::Item>::const_iterator it = map.begin();
49  it != map.end();it++) {
50 
51  uint32_t mapKey=it->first;
52  EcalTPGPhysicsConst::Item item=it->second;
53 
54  align->Fill(0.5, row, mapKey);
55  align->Fill(1.5, row, item.EtSat);
56  align->Fill(2.5, row, item.ttf_threshold_Low);
57  align->Fill(3.5, row, item.ttf_threshold_High);
58  align->Fill(4.5, row, item.FG_lowThreshold);
59  align->Fill(5.5, row, item.FG_highThreshold);
60  align->Fill(6.5, row, item.FG_lowRatio);
61  align->Fill(7.5, row, item.FG_highRatio);
62 
63  row = row - 1.;
64  }
65 
66 
67  }else
68  return false;
69 
70 
71  gStyle->SetPalette(1);
72  gStyle->SetOptStat(0);
73  TCanvas canvas("CC map", "CC map", 1000, 1000);
74  TLatex t1;
75  t1.SetNDC();
76  t1.SetTextAlign(26);
77  t1.SetTextSize(0.05);
78  t1.SetTextColor(2);
79  t1.DrawLatex(0.5, 0.96,Form("ECAL TPG Physics Constant, IOV %i", run));
80 
81  TPad* pad = new TPad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
82  pad->Draw();
83  pad->cd();
84  align->Draw("TEXT");
85 
86  drawTable(NbRows,8);
87 
88  align->GetXaxis()->SetTickLength(0.);
89  align->GetXaxis()->SetLabelSize(0.);
90  align->GetYaxis()->SetTickLength(0.);
91  align->GetYaxis()->SetLabelSize(0.);
92 
93  std::string ImageName(m_imageFileName);
94  canvas.SaveAs(ImageName.c_str());
95 
96  return true;
97  }
98 
99 };
100 
101 }
102 // Register the classes as boost python plugin
104  PAYLOAD_INSPECTOR_CLASS(EcalTPGPhysicsConstPlot);
105 }
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
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:481
std::map< uint32_t, EcalTPGPhysicsConst::Item > EcalTPGPhysicsConstMap
void drawTable(int nbRows, int nbColumns)