CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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  align = new TH2F("TPGPhysicsConstant",
44  "mapKey EtSat ttf_threshold_Low ttf_threshold_High FG_lowThreshold "
45  " FG_highThreshold FG_lowRatio FG_highRatio",
46  8,
47  0,
48  8,
49  NbRows,
50  0,
51  NbRows);
52 
53  double row = NbRows - 0.5;
54  for (std::map<uint32_t, EcalTPGPhysicsConst::Item>::const_iterator it = map.begin(); it != map.end(); it++) {
55  uint32_t mapKey = it->first;
56  EcalTPGPhysicsConst::Item item = it->second;
57 
58  align->Fill(0.5, row, mapKey);
59  align->Fill(1.5, row, item.EtSat);
60  align->Fill(2.5, row, item.ttf_threshold_Low);
61  align->Fill(3.5, row, item.ttf_threshold_High);
62  align->Fill(4.5, row, item.FG_lowThreshold);
63  align->Fill(5.5, row, item.FG_highThreshold);
64  align->Fill(6.5, row, item.FG_lowRatio);
65  align->Fill(7.5, row, item.FG_highRatio);
66 
67  row = row - 1.;
68  }
69 
70  } else
71  return false;
72 
73  gStyle->SetPalette(1);
74  gStyle->SetOptStat(0);
75  TCanvas canvas("CC map", "CC map", 1000, 1000);
76  TLatex t1;
77  t1.SetNDC();
78  t1.SetTextAlign(26);
79  t1.SetTextSize(0.05);
80  t1.SetTextColor(2);
81  t1.DrawLatex(0.5, 0.96, Form("ECAL TPG Physics Constant, IOV %i", run));
82 
83  TPad* pad = new TPad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
84  pad->Draw();
85  pad->cd();
86  align->Draw("TEXT");
87 
88  drawTable(NbRows, 8);
89 
90  align->GetXaxis()->SetTickLength(0.);
91  align->GetXaxis()->SetLabelSize(0.);
92  align->GetYaxis()->SetTickLength(0.);
93  align->GetYaxis()->SetLabelSize(0.);
94 
95  std::string ImageName(m_imageFileName);
96  canvas.SaveAs(ImageName.c_str());
97 
98  return true;
99  }
100  };
101 
102 } // namespace
103 // Register the classes as boost python plugin
def canvas
Definition: svgfig.py:482
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
std::map< uint32_t, EcalTPGPhysicsConst::Item > EcalTPGPhysicsConstMap
void drawTable(int nbRows, int nbColumns)
Definition: EcalDrawUtils.h:91
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)