File indexing completed on 2024-09-07 04:35:27
0001 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0002 #include "CondCore/Utilities/interface/PayloadInspector.h"
0003 #include "CondCore/CondDB/interface/Time.h"
0004 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0005 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0006 #include "CondCore/EcalPlugins/plugins/EcalDrawUtils.h"
0007
0008
0009 #include "CondFormats/EcalObjects/interface/EcalLaserAPDPNRatios.h"
0010
0011 #include "TH2F.h"
0012 #include "TCanvas.h"
0013 #include "TStyle.h"
0014 #include "TLine.h"
0015 #include "TLatex.h"
0016
0017 #include <memory>
0018 #include <sstream>
0019
0020 namespace {
0021 enum { kEBChannels = 61200, kEEChannels = 14648 };
0022 enum {
0023 MIN_IETA = 1,
0024 MIN_IPHI = 1,
0025 MAX_IETA = 85,
0026 MAX_IPHI = 360,
0027 EBhistEtaMax = 171
0028 };
0029 enum {
0030 IX_MIN = 1,
0031 IY_MIN = 1,
0032 IX_MAX = 100,
0033 IY_MAX = 100,
0034 EEhistXMax = 220
0035 };
0036
0037
0038
0039
0040
0041
0042
0043
0044 class EcalLaserAPDPNRatiosEBMap : public cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios> {
0045 public:
0046 EcalLaserAPDPNRatiosEBMap()
0047 : cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios>("ECAL Barrel APDPNRatios - map ",
0048 "iphi",
0049 MAX_IPHI,
0050 MIN_IPHI,
0051 MAX_IPHI + 1,
0052 "ieta",
0053 EBhistEtaMax,
0054 -MAX_IETA,
0055 MAX_IETA + 1) {
0056 Base::setSingleIov(true);
0057 }
0058
0059
0060 bool fill() override {
0061 auto tag = PlotBase::getTag<0>();
0062 for (auto const& iov : tag.iovs) {
0063 std::shared_ptr<EcalLaserAPDPNRatios> payload = Base::fetchPayload(std::get<1>(iov));
0064 if (payload.get()) {
0065
0066 for (int iphi = 1; iphi < 361; iphi++)
0067 fillWithValue(iphi, 0, 1);
0068
0069 for (int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
0070 uint32_t rawid = EBDetId::unhashIndex(cellid);
0071 float p2 = (payload->getLaserMap())[rawid].p2;
0072
0073 fillWithValue((EBDetId(rawid)).iphi(), (EBDetId(rawid)).ieta(), p2);
0074 }
0075 }
0076 }
0077
0078 return true;
0079 }
0080 };
0081
0082 class EcalLaserAPDPNRatiosEEMap : public cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios> {
0083 private:
0084 int EEhistSplit = 20;
0085
0086 public:
0087 EcalLaserAPDPNRatiosEEMap()
0088 : cond::payloadInspector::Histogram2D<EcalLaserAPDPNRatios>("ECAL Endcap APDPNRatios - map ",
0089 "ix",
0090 EEhistXMax,
0091 IX_MIN,
0092 EEhistXMax + 1,
0093 "iy",
0094 IY_MAX,
0095 IY_MIN,
0096 IY_MAX + 1) {
0097 Base::setSingleIov(true);
0098 }
0099
0100 bool fill() override {
0101 auto tag = PlotBase::getTag<0>();
0102 for (auto const& iov : tag.iovs) {
0103 std::shared_ptr<EcalLaserAPDPNRatios> payload = Base::fetchPayload(std::get<1>(iov));
0104 if (payload.get()) {
0105
0106 for (int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
0107 for (int iy = IY_MIN; iy < IY_MAX + 1; iy++)
0108 fillWithValue(ix, iy, 0);
0109
0110 for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {
0111 if (!EEDetId::validHashIndex(cellid))
0112 continue;
0113 uint32_t rawid = EEDetId::unhashIndex(cellid);
0114 float p2 = (payload->getLaserMap())[rawid].p2;
0115 EEDetId myEEId(rawid);
0116 if (myEEId.zside() == -1)
0117 fillWithValue(myEEId.ix(), myEEId.iy(), p2);
0118 else
0119 fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), p2);
0120 }
0121 }
0122 }
0123 return true;
0124 }
0125 };
0126
0127
0128
0129
0130 class EcalLaserAPDPNRatiosPlot : public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios> {
0131 public:
0132 EcalLaserAPDPNRatiosPlot()
0133 : cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios>("ECAL Laser APDPNRatios - map ") {
0134 setSingleIov(true);
0135 }
0136
0137 bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0138 TH2F** barrel = new TH2F*[3];
0139 TH2F** endc_p = new TH2F*[3];
0140 TH2F** endc_m = new TH2F*[3];
0141 float pEBmin[3], pEEmin[3];
0142
0143 for (int i = 0; i < 3; i++) {
0144 barrel[i] =
0145 new TH2F(Form("EBp%i", i), Form("EB p%i", i + 1), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
0146 endc_p[i] =
0147 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);
0148 endc_m[i] =
0149 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);
0150 pEBmin[i] = 10.;
0151 pEEmin[i] = 10.;
0152 }
0153
0154 auto iov = iovs.front();
0155 std::shared_ptr<EcalLaserAPDPNRatios> payload = fetchPayload(std::get<1>(iov));
0156 unsigned long IOV = std::get<0>(iov);
0157 int run = 0;
0158 if (IOV < 4294967296)
0159 run = std::get<0>(iov);
0160 else
0161 run = IOV >> 32;
0162 if (payload.get()) {
0163
0164 for (int cellid = 0; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
0165 uint32_t rawid = EBDetId::unhashIndex(cellid);
0166 Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
0167 Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
0168 if (eta > 0.)
0169 eta = eta - 0.5;
0170 else
0171 eta = eta + 0.5;
0172 float p1 = (payload->getLaserMap())[rawid].p1;
0173 if (p1 < pEBmin[0])
0174 pEBmin[0] = p1;
0175 float p2 = (payload->getLaserMap())[rawid].p2;
0176 if (p2 < pEBmin[1])
0177 pEBmin[1] = p2;
0178 float p3 = (payload->getLaserMap())[rawid].p3;
0179 if (p3 < pEBmin[2])
0180 pEBmin[2] = p3;
0181 barrel[0]->Fill(phi, eta, p1);
0182 barrel[1]->Fill(phi, eta, p2);
0183 barrel[2]->Fill(phi, eta, p3);
0184 }
0185
0186
0187 for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {
0188 if (!EEDetId::validHashIndex(cellid))
0189 continue;
0190 uint32_t rawid = EEDetId::unhashIndex(cellid);
0191 EEDetId myEEId(rawid);
0192 float p1 = (payload->getLaserMap())[rawid].p1;
0193 if (p1 < pEEmin[0])
0194 pEEmin[0] = p1;
0195 float p2 = (payload->getLaserMap())[rawid].p2;
0196 if (p2 < pEEmin[1])
0197 pEEmin[1] = p2;
0198 float p3 = (payload->getLaserMap())[rawid].p3;
0199 if (p3 < pEEmin[2])
0200 pEEmin[2] = p3;
0201 if (myEEId.zside() == 1) {
0202 endc_p[0]->Fill(myEEId.ix(), myEEId.iy(), p1);
0203 endc_p[1]->Fill(myEEId.ix(), myEEId.iy(), p2);
0204 endc_p[2]->Fill(myEEId.ix(), myEEId.iy(), p3);
0205 } else {
0206 endc_m[0]->Fill(myEEId.ix(), myEEId.iy(), p1);
0207 endc_m[1]->Fill(myEEId.ix(), myEEId.iy(), p2);
0208 endc_m[2]->Fill(myEEId.ix(), myEEId.iy(), p3);
0209 }
0210 }
0211 }
0212 else
0213 return false;
0214
0215 gStyle->SetPalette(1);
0216 gStyle->SetOptStat(0);
0217 TCanvas canvas("CC map", "CC map", 2800, 2600);
0218 TLatex t1;
0219 t1.SetNDC();
0220 t1.SetTextAlign(26);
0221 t1.SetTextSize(0.05);
0222 if (IOV < 4294967296)
0223 t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APD/PN, IOV %i", run));
0224 else {
0225 time_t t = run;
0226 char buf[256];
0227 struct tm lt;
0228 localtime_r(&t, <);
0229 strftime(buf, sizeof(buf), "%F %R:%S", <);
0230 buf[sizeof(buf) - 1] = 0;
0231 t1.DrawLatex(0.5, 0.96, Form("Ecal Laser APD/PN, IOV %s", buf));
0232 }
0233
0234 float xmi[3] = {0.0, 0.26, 0.74};
0235 float xma[3] = {0.26, 0.74, 1.00};
0236 TPad*** pad = new TPad**[3];
0237 for (int i = 0; i < 3; i++) {
0238 pad[i] = new TPad*[3];
0239 for (int obj = 0; obj < 3; obj++) {
0240 float yma = 0.94 - (0.32 * i);
0241 float ymi = yma - 0.28;
0242 pad[i][obj] = new TPad(Form("p_%i_%i", obj, i), Form("p_%i_%i", obj, i), xmi[obj], ymi, xma[obj], yma);
0243 pad[i][obj]->Draw();
0244 }
0245 }
0246
0247 for (int i = 0; i < 3; i++) {
0248
0249
0250 float xmin = pEBmin[i] * 10.;
0251 int min = (int)xmin;
0252 pEBmin[i] = (float)min / 10.;
0253
0254 xmin = pEEmin[i] * 10.;
0255 min = (int)xmin;
0256 pEEmin[i] = (float)min / 10.;
0257
0258 pad[i][0]->cd();
0259 DrawEE(endc_m[i], pEEmin[i], 1.1);
0260 pad[i][1]->cd();
0261 DrawEB(barrel[i], pEBmin[i], 1.1);
0262 pad[i][2]->cd();
0263 DrawEE(endc_p[i], pEEmin[i], 1.1);
0264 }
0265
0266 std::string ImageName(m_imageFileName);
0267 canvas.SaveAs(ImageName.c_str());
0268 return true;
0269 }
0270 };
0271
0272
0273
0274
0275 template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags, int method>
0276 class EcalLaserAPDPNRatiosBase : public cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios, nIOVs, ntags> {
0277 public:
0278 EcalLaserAPDPNRatiosBase()
0279 : cond::payloadInspector::PlotImage<EcalLaserAPDPNRatios, nIOVs, ntags>("ECAL Laser APDPNRatios difference") {}
0280
0281 bool fill() override {
0282 TH2F** barrel = new TH2F*[3];
0283 TH2F** endc_p = new TH2F*[3];
0284 TH2F** endc_m = new TH2F*[3];
0285 float pEBmin[3], pEEmin[3], pEBmax[3], pEEmax[3];
0286 for (int i = 0; i < 3; i++) {
0287 barrel[i] =
0288 new TH2F(Form("EBp%i", i), Form("EB p%i", i + 1), MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
0289 endc_p[i] =
0290 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);
0291 endc_m[i] =
0292 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);
0293 pEBmin[i] = 10.;
0294 pEEmin[i] = 10.;
0295 pEBmax[i] = -10.;
0296 pEEmax[i] = -10.;
0297 }
0298 unsigned int run[2] = {0, 0};
0299 std::string l_tagname[2];
0300 unsigned long IOV = 0;
0301 float pEB[3][kEBChannels], pEE[3][kEEChannels];
0302 auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
0303 l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
0304 auto firstiov = iovs.front();
0305 IOV = std::get<0>(firstiov);
0306 if (IOV < 4294967296)
0307 run[0] = std::get<0>(firstiov);
0308 else
0309 run[0] = IOV >> 32;
0310 std::tuple<cond::Time_t, cond::Hash> lastiov;
0311 if (ntags == 2) {
0312 auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
0313 l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
0314 lastiov = tag2iovs.front();
0315 } else {
0316 lastiov = iovs.back();
0317 l_tagname[1] = l_tagname[0];
0318 }
0319 IOV = std::get<0>(lastiov);
0320 if (IOV < 4294967296)
0321 run[1] = std::get<0>(lastiov);
0322 else
0323 run[1] = IOV >> 32;
0324
0325 for (int irun = 0; irun < nIOVs; irun++) {
0326 std::shared_ptr<EcalLaserAPDPNRatios> payload;
0327 if (irun == 0) {
0328 payload = this->fetchPayload(std::get<1>(firstiov));
0329 } else {
0330 payload = this->fetchPayload(std::get<1>(lastiov));
0331 }
0332 if (payload.get()) {
0333
0334 for (int cellid = 0; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
0335 uint32_t rawid = EBDetId::unhashIndex(cellid);
0336 if (irun == 0) {
0337 pEB[0][cellid] = (payload->getLaserMap())[rawid].p1;
0338 pEB[1][cellid] = (payload->getLaserMap())[rawid].p2;
0339 pEB[2][cellid] = (payload->getLaserMap())[rawid].p3;
0340 } else {
0341 Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
0342 Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
0343 if (eta > 0.)
0344 eta = eta - 0.5;
0345 else
0346 eta = eta + 0.5;
0347 double dr;
0348 if (method == 0)
0349 dr = (payload->getLaserMap())[rawid].p1 - pEB[0][cellid];
0350 else {
0351 if (pEB[0][cellid] == 0.) {
0352 if ((payload->getLaserMap())[rawid].p1 == 0.)
0353 dr = 1.;
0354 else
0355 dr = 9999.;
0356 } else
0357 dr = (payload->getLaserMap())[rawid].p1 / pEB[0][cellid];
0358 }
0359 if (dr < pEBmin[0])
0360 pEBmin[0] = dr;
0361 if (dr > pEBmax[0])
0362 pEBmax[0] = dr;
0363 barrel[0]->Fill(phi, eta, dr);
0364 if (method == 0)
0365 dr = (payload->getLaserMap())[rawid].p2 - pEB[1][cellid];
0366 else {
0367 if (pEB[1][cellid] == 0.) {
0368 if ((payload->getLaserMap())[rawid].p2 == 0.)
0369 dr = 1.;
0370 else
0371 dr = 9999.;
0372 } else
0373 dr = (payload->getLaserMap())[rawid].p2 / pEB[1][cellid];
0374 }
0375 if (dr < pEBmin[1])
0376 pEBmin[1] = dr;
0377 if (dr > pEBmax[1])
0378 pEBmax[1] = dr;
0379 barrel[1]->Fill(phi, eta, dr);
0380 if (method == 0)
0381 dr = (payload->getLaserMap())[rawid].p3 - pEB[2][cellid];
0382 else {
0383 if (pEB[2][cellid] == 0.) {
0384 if ((payload->getLaserMap())[rawid].p3 == 0.)
0385 dr = 1.;
0386 else
0387 dr = 9999.;
0388 } else
0389 dr = (payload->getLaserMap())[rawid].p3 / pEB[2][cellid];
0390 }
0391 if (dr < pEBmin[2])
0392 pEBmin[2] = dr;
0393 if (dr > pEBmax[2])
0394 pEBmax[2] = dr;
0395 barrel[2]->Fill(phi, eta, dr);
0396 }
0397 }
0398
0399
0400 for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {
0401 if (!EEDetId::validHashIndex(cellid))
0402 continue;
0403 uint32_t rawid = EEDetId::unhashIndex(cellid);
0404 EEDetId myEEId(rawid);
0405 if (irun == 0) {
0406 pEE[0][cellid] = (payload->getLaserMap())[rawid].p1;
0407 pEE[1][cellid] = (payload->getLaserMap())[rawid].p2;
0408 pEE[2][cellid] = (payload->getLaserMap())[rawid].p3;
0409 } else {
0410 double dr1, dr2, dr3;
0411 if (method == 0)
0412 dr1 = (payload->getLaserMap())[rawid].p1 - pEE[0][cellid];
0413 else {
0414 if (pEE[0][cellid] == 0.) {
0415 if ((payload->getLaserMap())[rawid].p1 == 0.)
0416 dr1 = 1.;
0417 else
0418 dr1 = 9999.;
0419 } else
0420 dr1 = (payload->getLaserMap())[rawid].p1 / pEE[0][cellid];
0421 }
0422 if (dr1 < pEEmin[0])
0423 pEEmin[0] = dr1;
0424 if (dr1 > pEEmax[0])
0425 pEEmax[0] = dr1;
0426 if (method == 0)
0427 dr2 = (payload->getLaserMap())[rawid].p2 - pEE[1][cellid];
0428 else {
0429 if (pEE[1][cellid] == 0.) {
0430 if ((payload->getLaserMap())[rawid].p2 == 0.)
0431 dr2 = 1.;
0432 else
0433 dr2 = 9999.;
0434 } else
0435 dr2 = (payload->getLaserMap())[rawid].p2 / pEE[1][cellid];
0436 }
0437 if (dr2 < pEEmin[1])
0438 pEEmin[1] = dr2;
0439 if (dr2 > pEEmax[1])
0440 pEEmax[1] = dr2;
0441 if (method == 0)
0442 dr3 = (payload->getLaserMap())[rawid].p3 - pEE[2][cellid];
0443 else {
0444 if (pEE[0][cellid] == 0.) {
0445 if ((payload->getLaserMap())[rawid].p3 == 0.)
0446 dr3 = 1.;
0447 else
0448 dr3 = 9999.;
0449 } else
0450 dr3 = (payload->getLaserMap())[rawid].p3 / pEE[2][cellid];
0451 }
0452 if (dr3 < pEEmin[2])
0453 pEEmin[2] = dr3;
0454 if (dr3 > pEEmax[2])
0455 pEEmax[2] = dr3;
0456 if (myEEId.zside() == 1) {
0457 endc_p[0]->Fill(myEEId.ix(), myEEId.iy(), dr1);
0458 endc_p[1]->Fill(myEEId.ix(), myEEId.iy(), dr2);
0459 endc_p[2]->Fill(myEEId.ix(), myEEId.iy(), dr3);
0460 } else {
0461 endc_m[0]->Fill(myEEId.ix(), myEEId.iy(), dr1);
0462 endc_m[1]->Fill(myEEId.ix(), myEEId.iy(), dr2);
0463 endc_m[2]->Fill(myEEId.ix(), myEEId.iy(), dr3);
0464 }
0465 }
0466 }
0467 }
0468 else
0469 return false;
0470 }
0471
0472 gStyle->SetPalette(1);
0473 gStyle->SetOptStat(0);
0474 TCanvas canvas("CC map", "CC map", 2800, 2600);
0475 TLatex t1;
0476 t1.SetNDC();
0477 t1.SetTextAlign(26);
0478 int len = l_tagname[0].length() + l_tagname[1].length();
0479 std::string dr[2] = {"-", "/"};
0480 if (IOV < 4294967296) {
0481 if (ntags == 2) {
0482 if (len < 80) {
0483 t1.SetTextSize(0.02);
0484 t1.DrawLatex(0.5,
0485 0.96,
0486 Form("%s IOV %i %s %s IOV %i",
0487 l_tagname[1].c_str(),
0488 run[1],
0489 dr[method].c_str(),
0490 l_tagname[0].c_str(),
0491 run[0]));
0492 } else {
0493 t1.SetTextSize(0.03);
0494 t1.DrawLatex(0.5, 0.96, Form("Ecal LaserAPDPNRatios, IOV %i %s %i", run[1], dr[method].c_str(), run[0]));
0495 }
0496 } else {
0497 t1.SetTextSize(0.03);
0498 t1.DrawLatex(0.5, 0.96, Form("%s, IOV %i %s %i", l_tagname[0].c_str(), run[1], dr[method].c_str(), run[0]));
0499 }
0500 } else {
0501 time_t t = run[0];
0502 char buf0[256], buf1[256];
0503 struct tm lt;
0504 localtime_r(&t, <);
0505 strftime(buf0, sizeof(buf0), "%F %R:%S", <);
0506 buf0[sizeof(buf0) - 1] = 0;
0507 t = run[1];
0508 localtime_r(&t, <);
0509 strftime(buf1, sizeof(buf1), "%F %R:%S", <);
0510 buf1[sizeof(buf1) - 1] = 0;
0511 if (ntags == 2) {
0512 if (len < 80) {
0513 t1.SetTextSize(0.02);
0514 t1.DrawLatex(0.5,
0515 0.96,
0516 Form("%s IOV %i %s %s IOV %i",
0517 l_tagname[1].c_str(),
0518 run[1],
0519 dr[method].c_str(),
0520 l_tagname[0].c_str(),
0521 run[0]));
0522 } else {
0523 t1.SetTextSize(0.03);
0524 t1.DrawLatex(0.5, 0.96, Form("Ecal LaserAPDPNRatios, IOV %i %s %i", run[1], dr[method].c_str(), run[0]));
0525 }
0526 } else {
0527 t1.SetTextSize(0.03);
0528 t1.DrawLatex(0.5, 0.96, Form("%s, IOV %i %s %i", l_tagname[0].c_str(), run[1], dr[method].c_str(), run[0]));
0529 }
0530 }
0531 float xmi[3] = {0.0, 0.24, 0.76};
0532 float xma[3] = {0.24, 0.76, 1.00};
0533 TPad*** pad = new TPad**[3];
0534 for (int i = 0; i < 3; i++) {
0535 pad[i] = new TPad*[3];
0536 for (int obj = 0; obj < 3; obj++) {
0537 float yma = 0.94 - (0.32 * i);
0538 float ymi = yma - 0.28;
0539 pad[i][obj] = new TPad(Form("p_%i_%i", obj, i), Form("p_%i_%i", obj, i), xmi[obj], ymi, xma[obj], yma);
0540 pad[i][obj]->Draw();
0541 }
0542 }
0543
0544 for (int i = 0; i < 3; i++) {
0545
0546
0547 float xmin = (pEBmin[i] - 0.009) * 100.;
0548 int min = (int)xmin;
0549 pEBmin[i] = (float)min / 100.;
0550 float xmax = (pEBmax[i] + 0.009) * 100.;
0551 int max = (int)xmax;
0552 pEBmax[i] = (float)max / 100.;
0553
0554 xmin = (pEEmin[i] + 0.009) * 100.;
0555 min = (int)xmin;
0556 pEEmin[i] = (float)min / 100.;
0557 xmax = (pEEmax[i] + 0.009) * 100.;
0558 max = (int)xmax;
0559 pEEmax[i] = (float)max / 100.;
0560
0561 pad[i][0]->cd();
0562 DrawEE(endc_m[i], pEEmin[i], pEEmax[i]);
0563 endc_m[i]->GetZaxis()->SetLabelSize(0.02);
0564 pad[i][1]->cd();
0565 DrawEB(barrel[i], pEBmin[i], pEBmax[i]);
0566 pad[i][2]->cd();
0567 DrawEE(endc_p[i], pEEmin[i], pEEmax[i]);
0568 endc_p[i]->GetZaxis()->SetLabelSize(0.02);
0569 }
0570
0571 std::string ImageName(this->m_imageFileName);
0572 canvas.SaveAs(ImageName.c_str());
0573 return true;
0574 }
0575 };
0576 using EcalLaserAPDPNRatiosDiffOneTag = EcalLaserAPDPNRatiosBase<cond::payloadInspector::SINGLE_IOV, 1, 0>;
0577 using EcalLaserAPDPNRatiosDiffTwoTags = EcalLaserAPDPNRatiosBase<cond::payloadInspector::SINGLE_IOV, 2, 0>;
0578 using EcalLaserAPDPNRatiosRatioOneTag = EcalLaserAPDPNRatiosBase<cond::payloadInspector::SINGLE_IOV, 1, 1>;
0579 using EcalLaserAPDPNRatiosRatioTwoTags = EcalLaserAPDPNRatiosBase<cond::payloadInspector::SINGLE_IOV, 2, 1>;
0580
0581 }
0582
0583
0584 PAYLOAD_INSPECTOR_MODULE(EcalLaserAPDPNRatios) {
0585 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosEBMap);
0586 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosEEMap);
0587 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosPlot);
0588 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosDiffOneTag);
0589 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosDiffTwoTags);
0590 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosRatioOneTag);
0591 PAYLOAD_INSPECTOR_CLASS(EcalLaserAPDPNRatiosRatioTwoTags);
0592 }