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