CMS 3D CMS Logo

EcalLaserAPDPNRatios_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"
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 {
23  MIN_IETA = 1,
24  MIN_IPHI = 1,
25  MAX_IETA = 85,
26  MAX_IPHI = 360,
27  EBhistEtaMax = 171
28  }; // barrel lower and upper bounds on eta and phi
29  enum {
30  IX_MIN = 1,
31  IY_MIN = 1,
32  IX_MAX = 100,
33  IY_MAX = 100,
34  EEhistXMax = 220
35  }; // endcaps lower and upper bounds on x and y
36 
37  /*******************************************************
38 
39  2d histogram of ECAL barrel APDPNRatios of 1 IOV
40 
41  *******************************************************/
42 
43  // inherit from one of the predefined plot class: Histogram2D
44  class EcalLaserAPDPNRatiosEBMap : public cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios> {
45  public:
46  EcalLaserAPDPNRatiosEBMap()
47  : cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios>("ECAL Barrel APDPNRatios - map ",
48  "iphi",
49  MAX_IPHI,
50  MIN_IPHI,
51  MAX_IPHI + 1,
52  "ieta",
54  -MAX_IETA,
55  MAX_IETA + 1) {
56  Base::setSingleIov(true);
57  }
58 
59  // Histogram2D::fill (virtual) needs be overridden - the implementation should use fillWithValue
60  bool fill() override {
61  auto tag = PlotBase::getTag<0>();
62  for (auto const& iov : tag.iovs) {
63  std::shared_ptr<EcalLaserAPDPNRatios> payload = Base::fetchPayload(std::get<1>(iov));
64  if (payload.get()) {
65  // set to -1 for ieta 0 (no crystal)
66  for (int iphi = 1; iphi < 361; iphi++)
67  fillWithValue(iphi, 0, -1);
68 
69  for (int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
70  uint32_t rawid = EBDetId::unhashIndex(cellid);
71  float p2 = (payload->getLaserMap())[rawid].p2;
72  // fill the Histogram2D here
73  fillWithValue((EBDetId(rawid)).iphi(), (EBDetId(rawid)).ieta(), p2);
74  } // loop over cellid
75  } // if payload.get()
76  } // loop over IOV's (1 in this case)
77 
78  return true;
79  } // fill method
80  };
81 
82  class EcalLaserAPDPNRatiosEEMap : public cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios> {
83  private:
84  int EEhistSplit = 20;
85 
86  public:
87  EcalLaserAPDPNRatiosEEMap()
88  : cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios>("ECAL Endcap APDPNRatios - map ",
89  "ix",
90  EEhistXMax,
91  IX_MIN,
92  EEhistXMax + 1,
93  "iy",
94  IY_MAX,
95  IY_MIN,
96  IY_MAX + 1) {
97  Base::setSingleIov(true);
98  }
99 
100  bool fill() override {
101  auto tag = PlotBase::getTag<0>();
102  for (auto const& iov : tag.iovs) {
103  std::shared_ptr<EcalLaserAPDPNRatios> payload = Base::fetchPayload(std::get<1>(iov));
104  if (payload.get()) {
105  // set to -1 everywhwere
106  for (int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
107  for (int iy = IY_MAX; iy < IY_MAX + 1; iy++)
108  fillWithValue(ix, iy, -1);
109 
110  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {
111  if (!EEDetId::validHashIndex(cellid))
112  continue;
113  uint32_t rawid = EEDetId::unhashIndex(cellid);
114  float p2 = (payload->getLaserMap())[rawid].p2;
115  EEDetId myEEId(rawid);
116  if (myEEId.zside() == -1)
117  fillWithValue(myEEId.ix(), myEEId.iy(), p2);
118  else
119  fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), p2);
120  } // loop over cellid
121  } // payload
122  } // loop over IOV's (1 in this case)
123  return true;
124  } // fill method
125  };
126 
127  /*************************************************
128  2d plot of ECAL IntercalibConstants of 1 IOV
129  *************************************************/
130  class EcalLaserAPDPNRatiosPlot : public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios> {
131  public:
132  EcalLaserAPDPNRatiosPlot()
133  : cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios>("ECAL Laser APDPNRatios - 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*[3];
139  TH2F** endc_p = new TH2F*[3];
140  TH2F** endc_m = new TH2F*[3];
141  float pEBmin[3], pEEmin[3];
142 
143  for (int i = 0; i < 3; i++) {
144  barrel[i] =
145  new TH2F(Form("EBp%i", i), Form("EB p%i", i + 1), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
146  endc_p[i] =
147  new TH2F(Form("EE+p%i", i), Form("EE+ p%i", i + 1), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
148  endc_m[i] =
149  new TH2F(Form("EE-p%i", i), Form("EE- p%i", i + 1), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
150  pEBmin[i] = 10.;
151  pEEmin[i] = 10.;
152  }
153 
154  auto iov = iovs.front();
155  std::shared_ptr<EcalLaserAPDPNRatios> payload = fetchPayload(std::get<1>(iov));
156  unsigned long IOV = std::get<0>(iov);
157  int run = 0;
158  if (IOV < 4294967296)
159  run = std::get<0>(iov);
160  else // time type IOV
161  run = IOV >> 32;
162  if (payload.get()) {
163  // looping over the EB channels, via the dense-index, mapped into EBDetId's
164  for (int cellid = 0; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) { // loop on EB cells
165  uint32_t rawid = EBDetId::unhashIndex(cellid);
166  Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
167  Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
168  if (eta > 0.)
169  eta = eta - 0.5; // 0.5 to 84.5
170  else
171  eta = eta + 0.5; // -84.5 to -0.5
172  float p1 = (payload->getLaserMap())[rawid].p1;
173  if (p1 < pEBmin[0])
174  pEBmin[0] = p1;
175  float p2 = (payload->getLaserMap())[rawid].p2;
176  if (p2 < pEBmin[1])
177  pEBmin[1] = p2;
178  float p3 = (payload->getLaserMap())[rawid].p3;
179  if (p3 < pEBmin[2])
180  pEBmin[2] = p3;
181  barrel[0]->Fill(phi, eta, p1);
182  barrel[1]->Fill(phi, eta, p2);
183  barrel[2]->Fill(phi, eta, p3);
184  } // loop over cellid
185 
186  // looping over the EE channels
187  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {
188  if (!EEDetId::validHashIndex(cellid))
189  continue;
190  uint32_t rawid = EEDetId::unhashIndex(cellid);
191  EEDetId myEEId(rawid);
192  float p1 = (payload->getLaserMap())[rawid].p1;
193  if (p1 < pEEmin[0])
194  pEEmin[0] = p1;
195  float p2 = (payload->getLaserMap())[rawid].p2;
196  if (p2 < pEEmin[1])
197  pEEmin[1] = p2;
198  float p3 = (payload->getLaserMap())[rawid].p3;
199  if (p3 < pEEmin[2])
200  pEEmin[2] = p3;
201  if (myEEId.zside() == 1) {
202  endc_p[0]->Fill(myEEId.ix(), myEEId.iy(), p1);
203  endc_p[1]->Fill(myEEId.ix(), myEEId.iy(), p2);
204  endc_p[2]->Fill(myEEId.ix(), myEEId.iy(), p3);
205  } else {
206  endc_m[0]->Fill(myEEId.ix(), myEEId.iy(), p1);
207  endc_m[1]->Fill(myEEId.ix(), myEEId.iy(), p2);
208  endc_m[2]->Fill(myEEId.ix(), myEEId.iy(), p3);
209  }
210  } // validDetId
211  } // if payload.get()
212  else
213  return false;
214 
215  gStyle->SetPalette(1);
216  gStyle->SetOptStat(0);
217  TCanvas canvas("CC map", "CC map", 2800, 2600);
218  TLatex t1;
219  t1.SetNDC();
220  t1.SetTextAlign(26);
221  t1.SetTextSize(0.05);
222  if (IOV < 4294967296)
223  t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APD/PN, IOV %i", run));
224  else { // time type IOV
225  time_t t = run;
226  char buf[256];
227  struct tm lt;
228  localtime_r(&t, &lt);
229  strftime(buf, sizeof(buf), "%F %R:%S", &lt);
230  buf[sizeof(buf) - 1] = 0;
231  t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APD/PN, IOV %s", buf));
232  }
233 
234  float xmi[3] = {0.0, 0.26, 0.74};
235  float xma[3] = {0.26, 0.74, 1.00};
236  TPad*** pad = new TPad**[3];
237  for (int i = 0; i < 3; i++) {
238  pad[i] = new TPad*[3];
239  for (int obj = 0; obj < 3; obj++) {
240  float yma = 0.94 - (0.32 * i);
241  float ymi = yma - 0.28;
242  pad[i][obj] = new TPad(Form("p_%i_%i", obj, i), Form("p_%i_%i", obj, i), xmi[obj], ymi, xma[obj], yma);
243  pad[i][obj]->Draw();
244  }
245  }
246 
247  for (int i = 0; i < 3; i++) {
248  // compute histo limits with some rounding
249  // std::cout << " before " << pEBmin[i];
250  float xmin = pEBmin[i] * 10.;
251  int min = (int)xmin;
252  pEBmin[i] = (float)min / 10.;
253  // std::cout << " after " << pEBmin[i] << std::endl << " before " << pEEmin[i];
254  xmin = pEEmin[i] * 10.;
255  min = (int)xmin;
256  pEEmin[i] = (float)min / 10.;
257  // std::cout << " after " << pEEmin[i] << std::endl;
258  pad[i][0]->cd();
259  DrawEE(endc_m[i], pEEmin[i], 1.1);
260  pad[i][1]->cd();
261  DrawEB(barrel[i], pEBmin[i], 1.1);
262  pad[i][2]->cd();
263  DrawEE(endc_p[i], pEEmin[i], 1.1);
264  }
265 
266  std::string ImageName(m_imageFileName);
267  canvas.SaveAs(ImageName.c_str());
268  return true;
269  } // fill method
270  };
271 
272  /*****************************************************************
273  2d plot of ECAL IntercalibConstants difference between 2 IOVs
274  ******************************************************************/
275  class EcalLaserAPDPNRatiosDiff : public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios> {
276  public:
277  EcalLaserAPDPNRatiosDiff()
278  : cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios>("ECAL Laser APDPNRatios difference") {
279  setSingleIov(false);
280  }
281 
282  bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
283  TH2F** barrel = new TH2F*[3];
284  TH2F** endc_p = new TH2F*[3];
285  TH2F** endc_m = new TH2F*[3];
286  float pEBmin[3], pEEmin[3], pEBmax[3], pEEmax[3];
287 
288  for (int i = 0; i < 3; i++) {
289  barrel[i] =
290  new TH2F(Form("EBp%i", i), Form("EB p%i", i + 1), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
291  endc_p[i] =
292  new TH2F(Form("EE+p%i", i), Form("EE+ p%i", i + 1), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
293  endc_m[i] =
294  new TH2F(Form("EE-p%i", i), Form("EE- p%i", i + 1), IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
295  pEBmin[i] = 10.;
296  pEEmin[i] = 10.;
297  pEBmax[i] = -10.;
298  pEEmax[i] = -10.;
299  }
300  unsigned int run[2] = {0, 0}, irun = 0;
301  unsigned long IOV = 0;
302  float pEB[3][kEBChannels], pEE[3][kEEChannels];
303  for (auto const& iov : iovs) {
304  std::shared_ptr<EcalLaserAPDPNRatios> payload = fetchPayload(std::get<1>(iov));
305  IOV = std::get<0>(iov);
306  if (IOV < 4294967296)
307  run[irun] = std::get<0>(iov);
308  else // time type IOV
309  run[irun] = IOV >> 32;
310  if (payload.get()) {
311  // looping over the EB channels, via the dense-index, mapped into EBDetId's
312  for (int cellid = 0; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) { // loop on EB cells
313  uint32_t rawid = EBDetId::unhashIndex(cellid);
314  if (irun == 0) {
315  pEB[0][cellid] = (payload->getLaserMap())[rawid].p1;
316  pEB[1][cellid] = (payload->getLaserMap())[rawid].p2;
317  pEB[2][cellid] = (payload->getLaserMap())[rawid].p3;
318  } else {
319  Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
320  Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
321  if (eta > 0.)
322  eta = eta - 0.5; // 0.5 to 84.5
323  else
324  eta = eta + 0.5; // -84.5 to -0.5
325  double diff = (payload->getLaserMap())[rawid].p1 - pEB[0][cellid];
326  if (diff < pEBmin[0])
327  pEBmin[0] = diff;
328  if (diff > pEBmax[0])
329  pEBmax[0] = diff;
330  barrel[0]->Fill(phi, eta, diff);
331  diff = (payload->getLaserMap())[rawid].p2 - pEB[1][cellid];
332  if (diff < pEBmin[1])
333  pEBmin[1] = diff;
334  if (diff > pEBmax[1])
335  pEBmax[1] = diff;
336  barrel[1]->Fill(phi, eta, diff);
337  diff = (payload->getLaserMap())[rawid].p3 - pEB[2][cellid];
338  if (diff < pEBmin[2])
339  pEBmin[2] = diff;
340  if (diff > pEBmax[2])
341  pEBmax[2] = diff;
342  barrel[2]->Fill(phi, eta, diff);
343  }
344  } // loop over cellid
345 
346  // looping over the EE channels
347  for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {
348  if (!EEDetId::validHashIndex(cellid))
349  continue;
350  uint32_t rawid = EEDetId::unhashIndex(cellid);
351  EEDetId myEEId(rawid);
352  if (irun == 0) {
353  pEE[0][cellid] = (payload->getLaserMap())[rawid].p1;
354  pEE[1][cellid] = (payload->getLaserMap())[rawid].p2;
355  pEE[2][cellid] = (payload->getLaserMap())[rawid].p3;
356  } else {
357  double diff1 = (payload->getLaserMap())[rawid].p1 - pEE[0][cellid];
358  if (diff1 < pEEmin[0])
359  pEEmin[0] = diff1;
360  if (diff1 > pEEmax[0])
361  pEEmax[0] = diff1;
362  double diff2 = (payload->getLaserMap())[rawid].p2 - pEE[1][cellid];
363  if (diff2 < pEEmin[1])
364  pEEmin[1] = diff2;
365  if (diff2 > pEEmax[1])
366  pEEmax[1] = diff2;
367  double diff3 = (payload->getLaserMap())[rawid].p3 - pEE[2][cellid];
368  if (diff3 < pEEmin[2])
369  pEEmin[2] = diff3;
370  if (diff3 > pEEmax[2])
371  pEEmax[2] = diff3;
372  if (myEEId.zside() == 1) {
373  endc_p[0]->Fill(myEEId.ix(), myEEId.iy(), diff1);
374  endc_p[1]->Fill(myEEId.ix(), myEEId.iy(), diff2);
375  endc_p[2]->Fill(myEEId.ix(), myEEId.iy(), diff3);
376  } else {
377  endc_m[0]->Fill(myEEId.ix(), myEEId.iy(), diff1);
378  endc_m[1]->Fill(myEEId.ix(), myEEId.iy(), diff2);
379  endc_m[2]->Fill(myEEId.ix(), myEEId.iy(), diff3);
380  }
381  }
382  } // loop over cellid
383  } // if payload.get()
384  else
385  return false;
386  irun++;
387  } // loop over IOVs
388 
389  gStyle->SetPalette(1);
390  gStyle->SetOptStat(0);
391  TCanvas canvas("CC map", "CC map", 2800, 2600);
392  TLatex t1;
393  t1.SetNDC();
394  t1.SetTextAlign(26);
395  t1.SetTextSize(0.05);
396  if (IOV < 4294967296) {
397  t1.SetTextSize(0.05);
398  t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APD/PN, IOV %i - %i", run[1], run[0]));
399  } else { // time type IOV
400  time_t t = run[0];
401  char buf0[256], buf1[256];
402  struct tm lt;
403  localtime_r(&t, &lt);
404  strftime(buf0, sizeof(buf0), "%F %R:%S", &lt);
405  buf0[sizeof(buf0) - 1] = 0;
406  t = run[1];
407  localtime_r(&t, &lt);
408  strftime(buf1, sizeof(buf1), "%F %R:%S", &lt);
409  buf1[sizeof(buf1) - 1] = 0;
410  t1.SetTextSize(0.015);
411  t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APD/PN, IOV %s - %s", buf1, buf0));
412  }
413 
414  float xmi[3] = {0.0, 0.24, 0.76};
415  float xma[3] = {0.24, 0.76, 1.00};
416  TPad*** pad = new TPad**[3];
417  for (int i = 0; i < 3; i++) {
418  pad[i] = new TPad*[3];
419  for (int obj = 0; obj < 3; obj++) {
420  float yma = 0.94 - (0.32 * i);
421  float ymi = yma - 0.28;
422  pad[i][obj] = new TPad(Form("p_%i_%i", obj, i), Form("p_%i_%i", obj, i), xmi[obj], ymi, xma[obj], yma);
423  pad[i][obj]->Draw();
424  }
425  }
426 
427  for (int i = 0; i < 3; i++) {
428  // compute histo limits with some rounding
429  // std::cout << " before min " << pEBmin[i] << " max " << pEBmax[i];
430  float xmin = (pEBmin[i] - 0.009) * 100.;
431  int min = (int)xmin;
432  pEBmin[i] = (float)min / 100.;
433  float xmax = (pEBmax[i] + 0.009) * 100.;
434  int max = (int)xmax;
435  pEBmax[i] = (float)max / 100.;
436  // std::cout << " after min " << pEBmin[i] << " max " << pEBmax[i] << std::endl << " before min " << pEEmin[i] << " max " << pEEmax[i];
437  xmin = (pEEmin[i] + 0.009) * 100.;
438  min = (int)xmin;
439  pEEmin[i] = (float)min / 100.;
440  xmax = (pEEmax[i] + 0.009) * 100.;
441  max = (int)xmax;
442  pEEmax[i] = (float)max / 100.;
443  // std::cout << " after min " << pEEmin[i] << " max " << pEEmax[i]<< std::endl;
444  pad[i][0]->cd();
445  DrawEE(endc_m[i], pEEmin[i], pEEmax[i]);
446  endc_m[i]->GetZaxis()->SetLabelSize(0.02);
447  pad[i][1]->cd();
448  DrawEB(barrel[i], pEBmin[i], pEBmax[i]);
449  pad[i][2]->cd();
450  DrawEE(endc_p[i], pEEmin[i], pEEmax[i]);
451  endc_p[i]->GetZaxis()->SetLabelSize(0.02);
452  }
453 
454  std::string ImageName(m_imageFileName);
455  canvas.SaveAs(ImageName.c_str());
456  return true;
457  } // fill method
458  };
459 
460 } // namespace
461 
462 // Register the classes as boost python plugin
464  PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosEBMap);
465  PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosEEMap);
466  PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosPlot);
467  PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosDiff);
468 }
change_name.diff
diff
Definition: change_name.py:13
svgfig.canvas
def canvas(*sub, **attr)
Definition: svgfig.py:482
mps_fire.i
i
Definition: mps_fire.py:428
Reference_intrackfit_cff.barrel
list barrel
Definition: Reference_intrackfit_cff.py:37
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
kEEChannels
Definition: EcalFloatCondObjectContainerUtils.h:9
MAX_IPHI
Definition: EcalFloatCondObjectContainerUtils.h:14
IX_MAX
Definition: EcalFloatCondObjectContainerUtils.h:20
IY_MAX
Definition: EcalFloatCondObjectContainerUtils.h:21
min
T min(T a, T b)
Definition: MathUtil.h:58
EBDetId
Definition: EBDetId.h:17
EEDetId::unhashIndex
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:65
PayloadInspector.h
EBDetId.h
EEDetId.h
EBDetId::unhashIndex
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:110
cond::payloadInspector::Plot3D< PayloadType, float, float, float, IOV_M, 1 >::fetchPayload
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
Definition: PayloadInspector.h:515
PAYLOAD_INSPECTOR_CLASS
#define PAYLOAD_INSPECTOR_CLASS(CLASS_NAME)
Definition: PayloadInspectorModule.h:10
PayloadInspectorModule.h
DrawEE
void DrawEE(TH2F *endc, float min, float max)
Definition: EcalDrawUtils.h:29
LEDCalibrationChannels.iphi
iphi
Definition: LEDCalibrationChannels.py:64
RandomServiceHelper.t1
t1
Definition: RandomServiceHelper.py:256
EEDetId::kSizeForDenseIndexing
Definition: EEDetId.h:329
GlobalPosition_Frontier_DevDB_cff.tag
tag
Definition: GlobalPosition_Frontier_DevDB_cff.py:11
DrawEB
void DrawEB(TH2F *ebmap, float min, float max)
Definition: EcalDrawUtils.h:4
PVValHelper::eta
Definition: PVValidationHelpers.h:69
MIN_IETA
Definition: EcalFloatCondObjectContainerUtils.h:11
p2
double p2[4]
Definition: TauolaWrapper.h:90
jets_cff.payload
payload
Definition: jets_cff.py:32
PAYLOAD_INSPECTOR_MODULE
#define PAYLOAD_INSPECTOR_MODULE(PAYLOAD_TYPENAME)
Definition: PayloadInspectorModule.h:8
getGTfromDQMFile.obj
obj
Definition: getGTfromDQMFile.py:32
EEDetId
Definition: EEDetId.h:14
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
cond
Definition: plugin.cc:23
cond::payloadInspector::Histogram2D
Definition: PayloadInspector.h:771
LEDCalibrationChannels.ieta
ieta
Definition: LEDCalibrationChannels.py:63
Time.h
IY_MIN
Definition: EcalFloatCondObjectContainerUtils.h:19
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
cond::payloadInspector::Histogram2D::fillWithValue
void fillWithValue(float xvalue, float yvalue, float weight=1)
Definition: PayloadInspector.h:810
cond::payloadInspector::PlotImage::fetchPayload
std::shared_ptr< PayloadType > fetchPayload(const cond::Hash &payloadHash)
Definition: PayloadInspector.h:863
createfilelist.int
int
Definition: createfilelist.py:10
MIN_IPHI
Definition: EcalFloatCondObjectContainerUtils.h:12
EEhistXMax
Definition: EcalFloatCondObjectContainerUtils.h:22
EEDetId::validHashIndex
static bool validHashIndex(int i)
Definition: EEDetId.h:239
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
p1
double p1[4]
Definition: TauolaWrapper.h:89
EcalDrawUtils.h
EcalLaserAPDPNRatios.h
visDQMUpload.buf
buf
Definition: visDQMUpload.py:154
EBhistEtaMax
Definition: EcalFloatCondObjectContainerUtils.h:15
DDAxes::phi
writedatasetfile.run
run
Definition: writedatasetfile.py:27
cond::payloadInspector::PlotImpl::fill
virtual bool fill()=0
TrackerOfflineValidation_Dqm_cff.xmax
xmax
Definition: TrackerOfflineValidation_Dqm_cff.py:11
cond::payloadInspector::PlotImage
Definition: PayloadInspector.h:852
EBDetId::kSizeForDenseIndexing
Definition: EBDetId.h:155
p3
double p3[4]
Definition: TauolaWrapper.h:91
EBDetId::MIN_HASH
static const int MIN_HASH
Definition: EBDetId.h:149
EcalLaserAPDPNRatios
Definition: EcalLaserAPDPNRatios.h:14
kEBChannels
Definition: EcalFloatCondObjectContainerUtils.h:9
TrackerOfflineValidation_Dqm_cff.xmin
xmin
Definition: TrackerOfflineValidation_Dqm_cff.py:10
IX_MIN
Definition: EcalFloatCondObjectContainerUtils.h:18
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
cond::payloadInspector::Histogram2D::fill
bool fill() override
Definition: PayloadInspector.h:820
MAX_IETA
Definition: EcalFloatCondObjectContainerUtils.h:13