CMS 3D CMS Logo

EcalLaserAPDPNRatiosRef_PayloadInspector.cc
Go to the documentation of this file.
8 // the data format of the condition to be inspected
10 
11 #include "TH2F.h"
12 #include "TCanvas.h"
13 #include "TStyle.h"
14 #include "TLine.h"
15 #include "TLatex.h"
16 
17 #include <memory>
18 #include <sstream>
19 
20 namespace {
21  enum {kEBChannels = 61200, kEEChannels = 14648};
22  enum {MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360, EBhistEtaMax = 171}; // barrel lower and upper bounds on eta and phi
23  enum {IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100, EEhistXMax = 220}; // endcaps lower and upper bounds on x and y
24 
25  /*******************************************************
26 
27  2d histogram of ECAL barrel Laser APDPN Ratios Ref of 1 IOV
28 
29  *******************************************************/
30 
31  // inherit from one of the predefined plot class: Histogram2D
32  class EcalLaserAPDPNRatiosRefEBMap : public cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatiosRef> {
33 
34  public:
35  EcalLaserAPDPNRatiosRefEBMap() : cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatiosRef>("ECAL Barrel Laser APDPN Ratios Ref- map ",
36  "iphi", MAX_IPHI, MIN_IPHI, MAX_IPHI + 1,
37  "ieta", EBhistEtaMax, -MAX_IETA, MAX_IETA + 1) {
38  Base::setSingleIov( true );
39  }
40 
41  // Histogram2D::fill (virtual) needs be overridden - the implementation should use fillWithValue
42  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
43 
44  for (auto const & iov: iovs) {
45  std::shared_ptr<EcalLaserAPDPNRatiosRef> payload = Base::fetchPayload( std::get<1>(iov) );
46  if( payload.get() ){
47  // looping over the EB channels, via the dense-index, mapped into EBDetId's
48  if (payload->barrelItems().empty()) return false;
49  // set to -1 for ieta 0 (no crystal)
50  for(int iphi = MIN_IPHI; iphi < MAX_IPHI+1; iphi++) fillWithValue(iphi, 0, -1);
51 
52  for(int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
53  uint32_t rawid = EBDetId::unhashIndex(cellid);
54 
55  // check the existence of ECAL LaserAPDPNRatiosRef, for a given ECAL barrel channel
56  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
57  if (value_ptr == payload->end())
58  continue; // cell absent from payload
59 
60  float weight = (float)(*value_ptr);
61 
62  // fill the Histogram2D here
63  fillWithValue( (EBDetId(rawid)).iphi() , (EBDetId(rawid)).ieta(), weight);
64  }// loop over cellid
65  }// if payload.get()
66  }// loop over IOV's (1 in this case)
67 
68  return true;
69 
70  }//fill method
71  };
72 
73  /*******************************************************
74 
75  2d histogram of ECAL EndCaps Laser APDPN Ratios Ref of 1 IOV
76 
77  *******************************************************/
78 
79  class EcalLaserAPDPNRatiosRefEEMap : public cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatiosRef> {
80 
81  private:
82  int EEhistSplit = 20;
83 
84  public:
85  EcalLaserAPDPNRatiosRefEEMap() : cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatiosRef>( "ECAL Endcap Laser APDPN Ratios Ref- map ",
86  "ix", EEhistXMax, IX_MIN, EEhistXMax + 1,
87  "iy", IY_MAX, IY_MIN, IY_MAX + 1) {
88  Base::setSingleIov( true );
89  }
90 
91  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
92 
93  for (auto const & iov: iovs) {
94  std::shared_ptr<EcalLaserAPDPNRatiosRef> payload = Base::fetchPayload( std::get<1>(iov) );
95  if( payload.get() ){
96  if (payload->endcapItems().empty()) return false;
97 
98  // set to -1 everywhwere
99  for(int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
100  for(int iy = IY_MAX; iy < IY_MAX + 1; iy++)
101  fillWithValue(ix, iy, -1);
102 
103  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid){ // loop on EE cells
104  if (EEDetId::validHashIndex(cellid)){
105  uint32_t rawid = EEDetId::unhashIndex(cellid);
106  EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
107  if (value_ptr == payload->end())
108  continue; // cell absent from payload
109 
110  float weight = (float)(*value_ptr);
111  EEDetId myEEId(rawid);
112  if(myEEId.zside() == -1)
113  fillWithValue(myEEId.ix(), myEEId.iy(), weight);
114  else
115  fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), weight);
116  } // validDetId
117  } // loop over cellid
118 
119  } // payload
120  } // loop over IOV's (1 in this case)
121  return true;
122  }// fill method
123 
124  };
125 
126 
127  /*************************************************
128  2d plot of Ecal Laser APDPN Ratios Ref of 1 IOV
129  *************************************************/
130  class EcalLaserAPDPNRatiosRefPlot : public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatiosRef> {
131 
132  public:
133  EcalLaserAPDPNRatiosRefPlot() : cond::payloadInspector::PlotImage<EcalLaserAPDPNRatiosRef>("Ecal Laser APDPN Ratios Ref - map ") {
134  setSingleIov( true );
135  }
136 
137  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
138  TH2F* barrel = new TH2F("EB", "Laser APDPN Ratios Ref EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
139  TH2F* endc_p = new TH2F("EE+", "Laser APDPN Ratios Ref EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
140  TH2F* endc_m = new TH2F("EE-", "Laser APDPN Ratios Ref EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
141 
142  auto iov = iovs.front();
143  std::shared_ptr<EcalLaserAPDPNRatiosRef> payload = fetchPayload( std::get<1>(iov) );
144  unsigned int run = std::get<0>(iov);
145 
146  if( payload.get() ){
147 
148  if (payload->barrelItems().empty())
149  return false;
150 
151  fillEBMap_SingleIOV<EcalLaserAPDPNRatiosRef>(payload, barrel);
152 
153 
154  if (payload->endcapItems().empty())
155  return false;
156 
157  fillEEMap_SingleIOV<EcalLaserAPDPNRatiosRef>(payload, endc_m, endc_p);
158 
159  } // payload
160 
161  gStyle->SetPalette(1);
162  gStyle->SetOptStat(0);
163  TCanvas canvas("CC map","CC map", 1600, 450);
164  TLatex t1;
165  t1.SetNDC();
166  t1.SetTextAlign(26);
167  t1.SetTextSize(0.04);
168  t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APDPN Ratios Ref, IOV %i", run));
169 
170  float xmi[3] = {0.0 , 0.24, 0.76};
171  float xma[3] = {0.24, 0.76, 1.00};
172  TPad** pad = new TPad*;
173  for (int obj = 0; obj < 3; obj++) {
174  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
175  pad[obj]->Draw();
176  }
177  // EcalDrawMaps ICMap;
178  pad[0]->cd();
179  // ICMap.DrawEE(endc_m, 0., 2.);
180  DrawEE(endc_m, 0., 2.5);
181  pad[1]->cd();
182  // ICMap.DrawEB(barrel, 0., 2.);
183  DrawEB(barrel, 0., 2.5);
184  pad[2]->cd();
185  // ICMap.DrawEE(endc_p, 0., 2.);
186  DrawEE(endc_p, 0., 2.5);
187 
188  std::string ImageName(m_imageFileName);
189  canvas.SaveAs(ImageName.c_str());
190  return true;
191  }// fill method
192  };
193 
194  /*****************************************************************
195  2d plot of Ecal Laser APDPN Ratios Ref difference between 2 IOVs
196  *****************************************************************/
197  class EcalLaserAPDPNRatiosRefDiff : public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatiosRef> {
198 
199  public:
200  EcalLaserAPDPNRatiosRefDiff() : cond::payloadInspector::PlotImage<EcalLaserAPDPNRatiosRef>("Ecal Laser APDPN Ratios Ref difference ") {
201  setSingleIov(false);
202  }
203 
204  bool fill( const std::vector<std::tuple<cond::Time_t,cond::Hash> >& iovs ) override{
205  TH2F* barrel = new TH2F("EB", "Laser APDPN Ratios Ref EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
206  TH2F* endc_p = new TH2F("EE+", "Laser APDPN Ratios Ref EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
207  TH2F* endc_m = new TH2F("EE-", "Laser APDPN Ratios Ref EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
208  float pEBmin, pEEmin, pEBmax, pEEmax;
209  pEBmin = 10.;
210  pEEmin = 10.;
211  pEBmax = -10.;
212  pEEmax = -10.;
213 
214  unsigned int run[2], irun = 0;
215  float pEB[kEBChannels], pEE[kEEChannels];
216  for ( auto const & iov: iovs) {
217  std::shared_ptr<EcalLaserAPDPNRatiosRef> payload = fetchPayload( std::get<1>(iov) );
218  run[irun] = std::get<0>(iov);
219 
220  if( payload.get() ){
221 
222  if (payload->barrelItems().empty())
223  return false;
224 
225  fillEBMap_DiffIOV<EcalLaserAPDPNRatiosRef>(payload, barrel, irun, pEB, pEBmin, pEBmax);
226 
227  if (payload->endcapItems().empty())
228  return false;
229 
230  fillEEMap_DiffIOV<EcalLaserAPDPNRatiosRef>(payload, endc_m, endc_p, irun, pEE, pEEmin, pEEmax);
231 
232  }// payload
233  irun++;
234 
235  }// loop over IOVs
236 
237  gStyle->SetPalette(1);
238  gStyle->SetOptStat(0);
239  TCanvas canvas("CC map","CC map", 1600, 450);
240  TLatex t1;
241  t1.SetNDC();
242  t1.SetTextAlign(26);
243  t1.SetTextSize(0.04);
244  t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APDPN Ratios Ref Diff, IOV %i - %i", run[1], run[0]));
245 
246  float xmi[3] = {0.0 , 0.24, 0.76};
247  float xma[3] = {0.24, 0.76, 1.00};
248  TPad** pad = new TPad*;
249 
250  for (int obj = 0; obj < 3; obj++) {
251  pad[obj] = new TPad(Form("p_%i", obj),Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
252  pad[obj]->Draw();
253  }
254 
255  pad[0]->cd();
256  DrawEE(endc_m, pEEmin, pEEmax);
257  pad[1]->cd();
258  DrawEB(barrel, pEBmin, pEBmax);
259  pad[2]->cd();
260  DrawEE(endc_p, pEEmin, pEEmax);
261 
262  std::string ImageName(m_imageFileName);
263  canvas.SaveAs(ImageName.c_str());
264  return true;
265  }// fill method
266  };
267 
268 
269 /*******************************************************
270  2d plot of Ecal Laser APDPN Ratios Ref Summary of 1 IOV
271  *******************************************************/
272 class EcalLaserAPDPNRatiosRefSummaryPlot: public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatiosRef>{
273  public:
274  EcalLaserAPDPNRatiosRefSummaryPlot():
275  cond::payloadInspector::PlotImage<EcalLaserAPDPNRatiosRef>("Ecal Laser APDPN RatiosRef Summary- map "){
276  setSingleIov(true);
277  }
278 
279  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs)override {
280  auto iov=iovs.front();
281  std::shared_ptr <EcalLaserAPDPNRatiosRef> payload = fetchPayload(std::get<1> (iov));
282  unsigned int run=std::get<0> (iov);
283  TH2F* align;
284  int NbRows;
285 
286  if(payload.get()){
287  NbRows=2;
288  align=new TH2F("", "", 0, 0, 0, 0, 0, 0);
289 
290 
291  float mean_x_EB=0.0f;
292  float mean_x_EE=0.0f;
293 
294  float rms_EB=0.0f;
295  float rms_EE=0.0f;
296 
297  int num_x_EB=0;
298  int num_x_EE=0;
299 
300 
301  payload->summary(mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
302  fillTableWithSummary(align, "Ecal Laser APDPN Ratios Ref", mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
303 
304  }else
305  return false;
306 
307  gStyle->SetPalette(1);
308  gStyle->SetOptStat(0);
309  TCanvas canvas("CC map", "CC map", 1000, 1000);
310  TLatex t1;
311  t1.SetNDC();
312  t1.SetTextAlign(26);
313  t1.SetTextSize(0.04);
314  t1.SetTextColor(2);
315  t1.DrawLatex(0.5, 0.96,Form("Ecal Laser APDPN Ratios Ref Summary, IOV %i", run));
316 
317 
318  TPad* pad = new TPad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
319  pad->Draw();
320  pad->cd();
321  align->Draw("TEXT");
322 
323  drawTable(NbRows,4);
324 
325 
326  std::string ImageName(m_imageFileName);
327  canvas.SaveAs(ImageName.c_str());
328 
329  return true;
330  }
331 };
332 
333 
334 }
335 
336 // Register the classes as boost python plugin
338  PAYLOAD_INSPECTOR_CLASS( EcalLaserAPDPNRatiosRefEBMap );
339  PAYLOAD_INSPECTOR_CLASS( EcalLaserAPDPNRatiosRefEEMap );
340  PAYLOAD_INSPECTOR_CLASS( EcalLaserAPDPNRatiosRefPlot );
341  PAYLOAD_INSPECTOR_CLASS( EcalLaserAPDPNRatiosRefDiff );
342  PAYLOAD_INSPECTOR_CLASS( EcalLaserAPDPNRatiosRefSummaryPlot );
343 }
void fillTableWithSummary(TH2F *&align, std::string title, const float &mean_x_EB, const float &rms_EB, const int &num_x_EB, const float &mean_x_EE, const float &rms_EE, const int &num_x_EE)
Definition: weight.py:1
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:65
static bool validHashIndex(int i)
Definition: EEDetId.h:239
void DrawEE(TH2F *endc, float min, float max)
Definition: EcalDrawUtils.h:29
void DrawEB(TH2F *ebmap, float min, float max)
Definition: EcalDrawUtils.h:4
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
static const int MIN_HASH
Definition: EBDetId.h:149
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
std::vector< Item >::const_iterator const_iterator
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:110
Definition: plugin.cc:24
def canvas(sub, attr)
Definition: svgfig.py:482
void drawTable(int nbRows, int nbColumns)