CMS 3D CMS Logo

EcalSimPulseShape_PayloadInspector.cc
Go to the documentation of this file.
4 
5 // the data format of the condition to be inspected
7 
8 #include "TProfile.h"
9 #include "TCanvas.h"
10 #include "TStyle.h"
11 #include "TLine.h"
12 #include "TLatex.h"
13 #include "TMarker.h"
14 
15 #include <string>
16 
17 namespace {
18  /********************************************
19  profile of ECAL SimPulseShape for 1 IOV
20  ********************************************/
21  class EcalSimPulseShapeProfile : public cond::payloadInspector::PlotImage<EcalSimPulseShape> {
22  public:
23  EcalSimPulseShapeProfile() : cond::payloadInspector::PlotImage<EcalSimPulseShape>("ECAL SimPulseShape - Profile ") {
24  setSingleIov(true);
25  }
26 
27  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> > &iovs) override {
28  auto iov = iovs.front();
29  std::shared_ptr<EcalSimPulseShape> payload = fetchPayload(std::get<1>(iov));
30  unsigned int run = std::get<0>(iov);
31  TProfile *barrel, *endcap, *apd;
32  int EBnbin = 0, EEnbin = 0, APDnbin = 0;
33  double EBxmax, EExmax, APDxmax, EBth, EEth, APDth;
34  if (payload.get()) {
35  EBth = (*payload).barrel_thresh;
36  EEth = (*payload).endcap_thresh;
37  APDth = (*payload).apd_thresh;
38  double time = (*payload).time_interval;
39  std::vector<double> EBshape = (*payload).barrel_shape;
40  std::vector<double> EEshape = (*payload).endcap_shape;
41  std::vector<double> APDshape = (*payload).apd_shape;
42  EBnbin = EBshape.size();
43  EBxmax = EBnbin * time;
44  EEnbin = EBshape.size();
45  EExmax = EEnbin * time;
46  APDnbin = APDshape.size();
47  APDxmax = APDnbin * time;
48 
49  /* std::cout << "_thresh barrel " << EBth << " endcap " << EEth << " apd " << APDth << std::endl
50  << " time interval " << time << std::endl
51  << " shape size barrel " << EBnbin << " endcap " << EEnbin << " apd " << APDnbin
52  << std::endl; */
53  barrel = new TProfile("EBshape", "", EBnbin, 0, EBxmax);
54  endcap = new TProfile("EEshape", "", EEnbin, 0, EExmax);
55  apd = new TProfile("APDshape", "", APDnbin, 0, APDxmax);
56  for (int s = 0; s < EBnbin; s++) {
57  double val = EBshape[s];
58  barrel->Fill(s, val);
59  }
60  for (int s = 0; s < EEnbin; s++) {
61  double val = EEshape[s];
62  endcap->Fill(s, val);
63  }
64  for (int s = 0; s < APDnbin; s++) {
65  double val = APDshape[s];
66  apd->Fill(s, val);
67  }
68  } // if payload.get()
69  else
70  return false;
71 
72  // gStyle->SetPalette(1);
73  gStyle->SetOptStat(0);
74  TCanvas canvas("ESPS", "ESPS", 1000, 500);
75  TLatex t1;
76  t1.SetNDC();
77  t1.SetTextAlign(26);
78  t1.SetTextSize(0.05);
79  t1.DrawLatex(0.5, 0.96, Form("Sim Pulse Shape, IOV %i", run));
80 
81  if (EBnbin == EEnbin && EBnbin == APDnbin) {
82  TPad *pad = new TPad("p_0", "p_0", 0.0, 0.0, 1.0, 0.95);
83  pad->Draw();
84  pad->cd();
85  barrel->SetXTitle("time (ns)");
86  barrel->SetYTitle("normalized amplitude (ADC#)");
87  barrel->SetMarkerColor(kBlack);
88  barrel->SetMarkerStyle(24);
89  // barrel->SetMarkerSize(0.5);
90  barrel->Draw("P");
91  TMarker *EBMarker = new TMarker(0.58, 0.85, 24);
92  EBMarker->SetNDC();
93  EBMarker->SetMarkerSize(1.0);
94  EBMarker->SetMarkerColor(kBlack);
95  EBMarker->Draw();
96  t1.SetTextAlign(12);
97  t1.DrawLatex(0.59, 0.85, Form("EB pulse, threshold %f", EBth));
98 
99  endcap->SetMarkerColor(kRed);
100  endcap->SetMarkerStyle(24);
101  // endcap->SetMarkerSize(0.5);
102  endcap->Draw("PSAME");
103  TMarker *EEMarker = new TMarker(0.58, 0.78, 24);
104  EEMarker->SetNDC();
105  EEMarker->SetMarkerSize(1.0);
106  EEMarker->SetMarkerColor(kRed);
107  EEMarker->Draw();
108  t1.SetTextColor(kRed);
109  t1.DrawLatex(0.59, 0.78, Form("EE pulse, threshold %f", EEth));
110 
111  apd->SetMarkerColor(kBlue);
112  apd->SetMarkerStyle(24);
113  // apd->SetMarkerSize(0.5);
114  apd->Draw("PSAME");
115  TMarker *APDMarker = new TMarker(0.58, 0.71, 24);
116  APDMarker->SetNDC();
117  APDMarker->SetMarkerSize(1.0);
118  APDMarker->SetMarkerColor(kBlue);
119  APDMarker->Draw();
120  t1.SetTextColor(kBlue);
121  t1.DrawLatex(0.59, 0.71, Form("APD pulse, threshold %f", APDth));
122  } else {
123  canvas.SetCanvasSize(1000, 1000);
124  TPad **pad = new TPad *[3];
125  for (int s = 0; s < 3; s++) {
126  float yma = 0.94 - (0.31 * s);
127  float ymi = yma - 0.29;
128  pad[s] = new TPad(Form("p_%i", s), Form("p_%i", s), 0.0, ymi, 1.0, yma);
129  pad[s]->Draw();
130  }
131  pad[0]->cd();
132  barrel->Draw("P");
133  barrel->SetXTitle("time (ns)");
134  barrel->SetYTitle("normalized amplitude (ADC#)");
135  barrel->SetMarkerColor(kBlack);
136  barrel->SetMarkerStyle(24);
137  TMarker *EBMarker = new TMarker(0.58, 0.80, 24);
138  EBMarker->SetNDC();
139  EBMarker->Draw();
140  EBMarker->SetMarkerSize(1.0);
141  EBMarker->SetMarkerColor(kBlack);
142  t1.SetTextAlign(12);
143  t1.DrawLatex(0.59, 0.80, Form("EB pulse, threshold %f", EBth));
144 
145  pad[1]->cd();
146  endcap->SetMarkerStyle(24);
147  endcap->SetMarkerColor(kRed);
148  endcap->Draw("P");
149  endcap->SetXTitle("time (ns)");
150  endcap->SetYTitle("normalized amplitude (ADC#)");
151  TMarker *EEMarker = new TMarker(0.58, 0.8, 24);
152  EEMarker->SetNDC();
153  EEMarker->Draw();
154  EEMarker->SetMarkerSize(1.0);
155  EEMarker->SetMarkerColor(kRed);
156  t1.SetTextColor(kRed);
157  t1.DrawLatex(0.59, 0.80, Form("EE pulse, threshold %f", EEth));
158 
159  pad[2]->cd();
160  apd->SetMarkerStyle(24);
161  apd->Draw("P");
162  apd->SetMarkerColor(kBlue);
163  apd->SetXTitle("time (ns)");
164  apd->SetYTitle("normalized amplitude (ADC#)");
165  TMarker *APDMarker = new TMarker(0.58, 0.8, 24);
166  APDMarker->SetNDC();
167  APDMarker->Draw();
168  APDMarker->SetMarkerSize(1.0);
169  APDMarker->SetMarkerColor(kBlue);
170  t1.SetTextColor(kBlue);
171  t1.DrawLatex(0.59, 0.80, Form("APD pulse, threshold %f", APDth));
172  }
173  std::string ImageName(m_imageFileName);
174  canvas.SaveAs(ImageName.c_str());
175  return true;
176  } // fill method
177  };
178 
179 } // namespace
180 
181 // Register the classes as boost python plugin
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
def canvas(sub, attr)
Definition: svgfig.py:482
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)