Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:37

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 // the data format of the condition to be inspected
0009 #include "CondFormats/EcalObjects/interface/EcalTPGLinearizationConst.h"
0010 
0011 #include "TH2F.h"
0012 #include "TCanvas.h"
0013 #include "TLine.h"
0014 #include "TStyle.h"
0015 #include "TLatex.h"
0016 #include "TPave.h"
0017 #include "TPaveStats.h"
0018 #include <string>
0019 #include <fstream>
0020 
0021 namespace {
0022   enum { kEBChannels = 61200, kEEChannels = 14648, kGains = 3, kSides = 2 };
0023   enum { MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360 };  // barrel lower and upper bounds on eta and phi
0024   enum { IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100 };         // endcaps lower and upper bounds on x and y
0025   int gainValues[kGains] = {12, 6, 1};
0026 
0027   /**************************************************
0028      2d plot of ECAL TPGLinearizationConst of 1 IOV
0029   **************************************************/
0030   class EcalTPGLinearizationConstPlot : public cond::payloadInspector::PlotImage<EcalTPGLinearizationConst> {
0031   public:
0032     EcalTPGLinearizationConstPlot()
0033         : cond::payloadInspector::PlotImage<EcalTPGLinearizationConst>("ECAL Gain Ratios - map ") {
0034       setSingleIov(true);
0035     }
0036 
0037     bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0038       TH2F** barrel_m = new TH2F*[kGains];
0039       TH2F** endc_p_m = new TH2F*[kGains];
0040       TH2F** endc_m_m = new TH2F*[kGains];
0041       TH2F** barrel_r = new TH2F*[kGains];
0042       TH2F** endc_p_r = new TH2F*[kGains];
0043       TH2F** endc_m_r = new TH2F*[kGains];
0044       float mEBmin[kGains], mEEmin[kGains], mEBmax[kGains], mEEmax[kGains], rEBmin[kGains], rEEmin[kGains],
0045           rEBmax[kGains], rEEmax[kGains];
0046       for (int gainId = 0; gainId < kGains; gainId++) {
0047         barrel_m[gainId] = new TH2F(Form("EBm%i", gainId),
0048                                     Form("EB mult_x%i ", gainValues[gainId]),
0049                                     MAX_IPHI,
0050                                     0,
0051                                     MAX_IPHI,
0052                                     2 * MAX_IETA,
0053                                     -MAX_IETA,
0054                                     MAX_IETA);
0055         endc_p_m[gainId] = new TH2F(Form("EE+m%i", gainId),
0056                                     Form("EE+ mult_x%i", gainValues[gainId]),
0057                                     IX_MAX,
0058                                     IX_MIN,
0059                                     IX_MAX + 1,
0060                                     IY_MAX,
0061                                     IY_MIN,
0062                                     IY_MAX + 1);
0063         endc_m_m[gainId] = new TH2F(Form("EE-m%i", gainId),
0064                                     Form("EE- mult_x%i", gainValues[gainId]),
0065                                     IX_MAX,
0066                                     IX_MIN,
0067                                     IX_MAX + 1,
0068                                     IY_MAX,
0069                                     IY_MIN,
0070                                     IY_MAX + 1);
0071         barrel_r[gainId] = new TH2F(Form("EBr%i", gainId),
0072                                     Form("EB shift_x%i", gainValues[gainId]),
0073                                     MAX_IPHI,
0074                                     0,
0075                                     MAX_IPHI,
0076                                     2 * MAX_IETA,
0077                                     -MAX_IETA,
0078                                     MAX_IETA);
0079         endc_p_r[gainId] = new TH2F(Form("EE+r%i", gainId),
0080                                     Form("EE+ shift_x%i", gainValues[gainId]),
0081                                     IX_MAX,
0082                                     IX_MIN,
0083                                     IX_MAX + 1,
0084                                     IY_MAX,
0085                                     IY_MIN,
0086                                     IY_MAX + 1);
0087         endc_m_r[gainId] = new TH2F(Form("EE-r%i", gainId),
0088                                     Form("EE- shift_x%i", gainValues[gainId]),
0089                                     IX_MAX,
0090                                     IX_MIN,
0091                                     IX_MAX + 1,
0092                                     IY_MAX,
0093                                     IY_MIN,
0094                                     IY_MAX + 1);
0095         mEBmin[gainId] = 10.;
0096         mEEmin[gainId] = 10.;
0097         mEBmax[gainId] = -10.;
0098         mEEmax[gainId] = -10.;
0099         rEBmin[gainId] = 10.;
0100         rEEmin[gainId] = 10.;
0101         rEBmax[gainId] = -10.;
0102         rEEmax[gainId] = -10.;
0103       }
0104 
0105       //      std::ofstream fout;
0106       //      fout.open("./bid.txt");
0107       auto iov = iovs.front();
0108       std::shared_ptr<EcalTPGLinearizationConst> payload = fetchPayload(std::get<1>(iov));
0109       unsigned int run = std::get<0>(iov);
0110       if (payload.get()) {
0111         for (int sign = 0; sign < kSides; sign++) {
0112           int thesign = sign == 1 ? 1 : -1;
0113 
0114           for (int ieta = 0; ieta < MAX_IETA; ieta++) {
0115             for (int iphi = 0; iphi < MAX_IPHI; iphi++) {
0116               EBDetId id((ieta + 1) * thesign, iphi + 1);
0117               float y = -1 - ieta;
0118               if (sign == 1)
0119                 y = ieta;
0120               float val = (*payload)[id.rawId()].mult_x12;
0121               barrel_m[0]->Fill(iphi, y, val);
0122               if (val < mEBmin[0])
0123                 mEBmin[0] = val;
0124               if (val > mEBmax[0])
0125                 mEBmax[0] = val;
0126               val = (*payload)[id.rawId()].shift_x12;
0127               barrel_r[0]->Fill(iphi, y, val);
0128               if (val < rEBmin[0])
0129                 rEBmin[0] = val;
0130               if (val > rEBmax[0])
0131                 rEBmax[0] = val;
0132               val = (*payload)[id.rawId()].mult_x6;
0133               barrel_m[1]->Fill(iphi, y, val);
0134               if (val < mEBmin[1])
0135                 mEBmin[1] = val;
0136               if (val > mEBmax[1])
0137                 mEBmax[1] = val;
0138               val = (*payload)[id.rawId()].shift_x6;
0139               barrel_r[1]->Fill(iphi, y, val);
0140               if (val < rEBmin[1])
0141                 rEBmin[1] = val;
0142               if (val > rEBmax[1])
0143                 rEBmax[1] = val;
0144               val = (*payload)[id.rawId()].mult_x1;
0145               barrel_m[2]->Fill(iphi, y, val);
0146               if (val < mEBmin[2])
0147                 mEBmin[2] = val;
0148               if (val > mEBmax[2])
0149                 mEBmax[2] = val;
0150               val = (*payload)[id.rawId()].shift_x1;
0151               barrel_r[2]->Fill(iphi, y, val);
0152               if (val < rEBmin[2])
0153                 rEBmin[2] = val;
0154               if (val > rEBmax[2])
0155                 rEBmax[2] = val;
0156             }  // iphi
0157           }    // ieta
0158 
0159           for (int ix = 0; ix < IX_MAX; ix++) {
0160             for (int iy = 0; iy < IY_MAX; iy++) {
0161               if (!EEDetId::validDetId(ix + 1, iy + 1, thesign))
0162                 continue;
0163               EEDetId id(ix + 1, iy + 1, thesign);
0164               float val = (*payload)[id.rawId()].mult_x12;
0165               if (thesign == 1)
0166                 endc_p_m[0]->Fill(ix + 1, iy + 1, val);
0167               else
0168                 endc_m_m[0]->Fill(ix + 1, iy + 1, val);
0169               if (val < mEEmin[0])
0170                 mEEmin[0] = val;
0171               if (val > mEEmax[0])
0172                 mEEmax[0] = val;
0173               val = (*payload)[id.rawId()].shift_x12;
0174               if (thesign == 1)
0175                 endc_p_r[0]->Fill(ix + 1, iy + 1, val);
0176               else
0177                 endc_m_r[0]->Fill(ix + 1, iy + 1, val);
0178               if (val < rEEmin[0])
0179                 rEEmin[0] = val;
0180               if (val > rEEmax[0])
0181                 rEEmax[0] = val;
0182               val = (*payload)[id.rawId()].mult_x6;
0183               if (thesign == 1)
0184                 endc_p_m[1]->Fill(ix + 1, iy + 1, val);
0185               else
0186                 endc_m_m[1]->Fill(ix + 1, iy + 1, val);
0187               if (val < mEEmin[1])
0188                 mEEmin[1] = val;
0189               if (val > mEEmax[1])
0190                 mEEmax[1] = val;
0191               val = (*payload)[id.rawId()].shift_x6;
0192               if (thesign == 1)
0193                 endc_p_r[1]->Fill(ix + 1, iy + 1, val);
0194               else
0195                 endc_m_r[1]->Fill(ix + 1, iy + 1, val);
0196               if (val < rEEmin[1])
0197                 rEEmin[1] = val;
0198               if (val > rEEmax[1])
0199                 rEEmax[1] = val;
0200               val = (*payload)[id.rawId()].mult_x1;
0201               if (thesign == 1)
0202                 endc_p_m[2]->Fill(ix + 1, iy + 1, val);
0203               else
0204                 endc_m_m[2]->Fill(ix + 1, iy + 1, val);
0205               if (val < mEEmin[2])
0206                 mEEmin[2] = val;
0207               if (val > mEEmax[2])
0208                 mEEmax[2] = val;
0209               val = (*payload)[id.rawId()].shift_x1;
0210               if (thesign == 1)
0211                 endc_p_r[2]->Fill(ix + 1, iy + 1, val);
0212               else
0213                 endc_m_r[2]->Fill(ix + 1, iy + 1, val);
0214               if (val < rEEmin[2])
0215                 rEEmin[2] = val;
0216               if (val > rEEmax[2])
0217                 rEEmax[2] = val;
0218               //          fout << " x " << ix << " y " << " val " << val << std::endl;
0219             }  // iy
0220           }    // ix
0221         }      // side
0222       }        // if payload.get()
0223       else
0224         return false;
0225       //      std::cout << " min " << rEEmin[2] << " max " << rEEmax[2] << std::endl;
0226       //      fout.close();
0227 
0228       gStyle->SetPalette(1);
0229       gStyle->SetOptStat(0);
0230       TCanvas canvas("CC map", "CC map", 1200, 1800);
0231       TLatex t1;
0232       t1.SetNDC();
0233       t1.SetTextAlign(26);
0234       t1.SetTextSize(0.05);
0235       t1.DrawLatex(0.5, 0.96, Form("Ecal Gain TPGLinearizationConst, IOV %i", run));
0236 
0237       float xmi[3] = {0.0, 0.22, 0.78};
0238       float xma[3] = {0.22, 0.78, 1.00};
0239       TPad*** pad = new TPad**[6];
0240       for (int gId = 0; gId < 6; gId++) {
0241         pad[gId] = new TPad*[3];
0242         for (int obj = 0; obj < 3; obj++) {
0243           float yma = 0.94 - (0.16 * gId);
0244           float ymi = yma - 0.14;
0245           pad[gId][obj] = new TPad(Form("p_%i_%i", obj, gId), Form("p_%i_%i", obj, gId), xmi[obj], ymi, xma[obj], yma);
0246           pad[gId][obj]->Draw();
0247         }
0248       }
0249 
0250       for (int gId = 0; gId < kGains; gId++) {
0251         pad[gId][0]->cd();
0252         DrawEE(endc_m_m[gId], mEEmin[gId], mEEmax[gId]);
0253         pad[gId + 3][0]->cd();
0254         DrawEE(endc_m_r[gId], rEEmin[gId], rEEmax[gId]);
0255         pad[gId][1]->cd();
0256         DrawEB(barrel_m[gId], mEBmin[gId], mEBmax[gId]);
0257         pad[gId + 3][1]->cd();
0258         DrawEB(barrel_r[gId], rEBmin[gId], rEBmax[gId]);
0259         pad[gId][2]->cd();
0260         DrawEE(endc_p_m[gId], mEEmin[gId], mEEmax[gId]);
0261         pad[gId + 3][2]->cd();
0262         DrawEE(endc_p_r[gId], rEEmin[gId], rEEmax[gId]);
0263       }
0264 
0265       std::string ImageName(m_imageFileName);
0266       canvas.SaveAs(ImageName.c_str());
0267       return true;
0268     }  // fill method
0269   };
0270 
0271   /******************************************************************
0272      2d plot of ECAL TPGLinearizationConst difference between 2 IOVs
0273   ******************************************************************/
0274   template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags, int method>
0275   class EcalTPGLinearizationConstBase
0276       : public cond::payloadInspector::PlotImage<EcalTPGLinearizationConst, nIOVs, ntags> {
0277   public:
0278     EcalTPGLinearizationConstBase()
0279         : cond::payloadInspector::PlotImage<EcalTPGLinearizationConst, nIOVs, ntags>("ECAL Gain Ratios comparison") {}
0280 
0281     bool fill() override {
0282       TH2F** barrel_m = new TH2F*[kGains];
0283       TH2F** endc_p_m = new TH2F*[kGains];
0284       TH2F** endc_m_m = new TH2F*[kGains];
0285       TH2F** barrel_r = new TH2F*[kGains];
0286       TH2F** endc_p_r = new TH2F*[kGains];
0287       TH2F** endc_m_r = new TH2F*[kGains];
0288       float mEBmin[kGains], mEEmin[kGains], mEBmax[kGains], mEEmax[kGains], rEBmin[kGains], rEEmin[kGains],
0289           rEBmax[kGains], rEEmax[kGains];
0290       float mEB[kGains][kEBChannels], mEE[kGains][kEEChannels], rEB[kGains][kEBChannels], rEE[kGains][kEEChannels];
0291       for (int gainId = 0; gainId < kGains; gainId++) {
0292         barrel_m[gainId] = new TH2F(Form("EBm%i", gainId),
0293                                     Form("EB mult_x%i ", gainValues[gainId]),
0294                                     MAX_IPHI,
0295                                     0,
0296                                     MAX_IPHI,
0297                                     2 * MAX_IETA,
0298                                     -MAX_IETA,
0299                                     MAX_IETA);
0300         endc_p_m[gainId] = new TH2F(Form("EE+m%i", gainId),
0301                                     Form("EE+ mult_x%i", gainValues[gainId]),
0302                                     IX_MAX,
0303                                     IX_MIN,
0304                                     IX_MAX + 1,
0305                                     IY_MAX,
0306                                     IY_MIN,
0307                                     IY_MAX + 1);
0308         endc_m_m[gainId] = new TH2F(Form("EE-m%i", gainId),
0309                                     Form("EE- mult_x%i", gainValues[gainId]),
0310                                     IX_MAX,
0311                                     IX_MIN,
0312                                     IX_MAX + 1,
0313                                     IY_MAX,
0314                                     IY_MIN,
0315                                     IY_MAX + 1);
0316         barrel_r[gainId] = new TH2F(Form("EBr%i", gainId),
0317                                     Form("EB shift_x%i", gainValues[gainId]),
0318                                     MAX_IPHI,
0319                                     0,
0320                                     MAX_IPHI,
0321                                     2 * MAX_IETA,
0322                                     -MAX_IETA,
0323                                     MAX_IETA);
0324         endc_p_r[gainId] = new TH2F(Form("EE+r%i", gainId),
0325                                     Form("EE+ shift_x%i", gainValues[gainId]),
0326                                     IX_MAX,
0327                                     IX_MIN,
0328                                     IX_MAX + 1,
0329                                     IY_MAX,
0330                                     IY_MIN,
0331                                     IY_MAX + 1);
0332         endc_m_r[gainId] = new TH2F(Form("EE-r%i", gainId),
0333                                     Form("EE- shift_x%i", gainValues[gainId]),
0334                                     IX_MAX,
0335                                     IX_MIN,
0336                                     IX_MAX + 1,
0337                                     IY_MAX,
0338                                     IY_MIN,
0339                                     IY_MAX + 1);
0340         mEBmin[gainId] = 10.;
0341         mEEmin[gainId] = 10.;
0342         mEBmax[gainId] = -10.;
0343         mEEmax[gainId] = -10.;
0344         rEBmin[gainId] = 10.;
0345         rEEmin[gainId] = 10.;
0346         rEBmax[gainId] = -10.;
0347         rEEmax[gainId] = -10.;
0348       }
0349 
0350       unsigned int run[2];
0351       //float gEB[3][kEBChannels], gEE[3][kEEChannels];
0352       std::string l_tagname[2];
0353       auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
0354       l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
0355       auto firstiov = iovs.front();
0356       run[0] = std::get<0>(firstiov);
0357       std::tuple<cond::Time_t, cond::Hash> lastiov;
0358       if (ntags == 2) {
0359         auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
0360         l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
0361         lastiov = tag2iovs.front();
0362       } else {
0363         lastiov = iovs.back();
0364         l_tagname[1] = l_tagname[0];
0365       }
0366       run[1] = std::get<0>(lastiov);
0367       for (int irun = 0; irun < nIOVs; irun++) {
0368         std::shared_ptr<EcalTPGLinearizationConst> payload;
0369         if (irun == 0) {
0370           payload = this->fetchPayload(std::get<1>(firstiov));
0371         } else {
0372           payload = this->fetchPayload(std::get<1>(lastiov));
0373         }
0374         if (payload.get()) {
0375           float dr;
0376           for (int sign = 0; sign < kSides; sign++) {
0377             int thesign = sign == 1 ? 1 : -1;
0378 
0379             for (int ieta = 0; ieta < MAX_IETA; ieta++) {
0380               for (int iphi = 0; iphi < MAX_IPHI; iphi++) {
0381                 EBDetId id((ieta + 1) * thesign, iphi + 1);
0382                 int hashindex = id.hashedIndex();
0383                 float y = -1 - ieta;
0384                 if (sign == 1)
0385                   y = ieta;
0386                 float val = (*payload)[id.rawId()].mult_x12;
0387                 if (irun == 0) {
0388                   mEB[0][hashindex] = val;
0389                 } else {
0390                   if (method == 0)
0391                     dr = val - mEB[0][hashindex];
0392                   else {  // ratio
0393                     if (mEB[0][hashindex] == 0.) {
0394                       if (val == 0.)
0395                         dr = 1.;
0396                       else
0397                         dr = 9999.;
0398                     } else
0399                       dr = val / mEB[0][hashindex];
0400                   }
0401                   barrel_m[0]->Fill(iphi, y, dr);
0402                   if (dr < mEBmin[0])
0403                     mEBmin[0] = dr;
0404                   if (dr > mEBmax[0])
0405                     mEBmax[0] = dr;
0406                 }
0407                 val = (*payload)[id.rawId()].shift_x12;
0408                 if (irun == 0) {
0409                   rEB[0][hashindex] = val;
0410                 } else {
0411                   if (method == 0)
0412                     dr = val - rEB[0][hashindex];
0413                   else {  // ratio
0414                     if (rEB[0][hashindex] == 0.) {
0415                       if (val == 0.)
0416                         dr = 1.;
0417                       else
0418                         dr = 9999.;
0419                     } else
0420                       dr = val / rEB[0][hashindex];
0421                   }
0422                   barrel_r[0]->Fill(iphi, y, dr);
0423                   if (dr < rEBmin[0])
0424                     rEBmin[0] = dr;
0425                   if (dr > rEBmax[0])
0426                     rEBmax[0] = dr;
0427                 }
0428                 val = (*payload)[id.rawId()].mult_x6;
0429                 if (irun == 0) {
0430                   mEB[1][hashindex] = val;
0431                 } else {
0432                   if (method == 0)
0433                     dr = val - mEB[1][hashindex];
0434                   else {  // ratio
0435                     if (mEB[1][hashindex] == 0.) {
0436                       if (val == 0.)
0437                         dr = 1.;
0438                       else
0439                         dr = 9999.;
0440                     } else
0441                       dr = val / mEB[1][hashindex];
0442                   }
0443                   barrel_m[1]->Fill(iphi, y, dr);
0444                   if (dr < mEBmin[1])
0445                     mEBmin[1] = dr;
0446                   if (dr > mEBmax[1])
0447                     mEBmax[1] = dr;
0448                 }
0449                 val = (*payload)[id.rawId()].shift_x6;
0450                 if (irun == 0) {
0451                   rEB[1][hashindex] = val;
0452                 } else {
0453                   if (method == 0)
0454                     dr = val - rEB[1][hashindex];
0455                   else {  // ratio
0456                     if (rEB[1][hashindex] == 0.) {
0457                       if (val == 0.)
0458                         dr = 1.;
0459                       else
0460                         dr = 9999.;
0461                     } else
0462                       dr = val / rEB[1][hashindex];
0463                   }
0464                   barrel_r[1]->Fill(iphi, y, dr);
0465                   if (dr < rEBmin[1])
0466                     rEBmin[1] = dr;
0467                   if (dr > rEBmax[1])
0468                     rEBmax[1] = dr;
0469                 }
0470                 val = (*payload)[id.rawId()].mult_x1;
0471                 if (irun == 0) {
0472                   mEB[2][hashindex] = val;
0473                 } else {
0474                   if (method == 0)
0475                     dr = val - mEB[2][hashindex];
0476                   else {  // ratio
0477                     if (mEB[2][hashindex] == 0.) {
0478                       if (val == 0.)
0479                         dr = 1.;
0480                       else
0481                         dr = 9999.;
0482                     } else
0483                       dr = val / mEB[2][hashindex];
0484                   }
0485                   barrel_m[2]->Fill(iphi, y, dr);
0486                   if (dr < mEBmin[2])
0487                     mEBmin[2] = dr;
0488                   if (dr > mEBmax[2])
0489                     mEBmax[2] = dr;
0490                 }
0491                 val = (*payload)[id.rawId()].shift_x1;
0492                 if (irun == 0) {
0493                   rEB[2][hashindex] = val;
0494                 } else {
0495                   if (method == 0)
0496                     dr = val - rEB[2][hashindex];
0497                   else {  // ratio
0498                     if (rEB[2][hashindex] == 0.) {
0499                       if (val == 0.)
0500                         dr = 1.;
0501                       else
0502                         dr = 9999.;
0503                     } else
0504                       dr = val / rEB[2][hashindex];
0505                   }
0506                   barrel_r[2]->Fill(iphi, y, dr);
0507                   if (dr < rEBmin[2])
0508                     rEBmin[2] = dr;
0509                   if (dr > rEBmax[2])
0510                     rEBmax[2] = dr;
0511                 }
0512               }  // iphi
0513             }    // ieta
0514 
0515             for (int ix = 0; ix < IX_MAX; ix++) {
0516               for (int iy = 0; iy < IY_MAX; iy++) {
0517                 if (!EEDetId::validDetId(ix + 1, iy + 1, thesign))
0518                   continue;
0519                 EEDetId id(ix + 1, iy + 1, thesign);
0520                 int hashindex = id.hashedIndex();
0521                 float val = (*payload)[id.rawId()].mult_x12;
0522                 if (irun == 0) {
0523                   mEE[0][hashindex] = val;
0524                 } else {
0525                   if (method == 0)
0526                     dr = val - mEE[0][hashindex];
0527                   else {  // ratio
0528                     if (mEE[0][hashindex] == 0.) {
0529                       if (val == 0.)
0530                         dr = 1.;
0531                       else
0532                         dr = 9999.;
0533                     } else
0534                       dr = val / mEE[0][hashindex];
0535                   }
0536                   if (thesign == 1)
0537                     endc_p_m[0]->Fill(ix + 1, iy + 1, dr);
0538                   else
0539                     endc_m_m[0]->Fill(ix + 1, iy + 1, dr);
0540                   if (dr < mEEmin[0])
0541                     mEEmin[0] = dr;
0542                   if (dr > mEEmax[0])
0543                     mEEmax[0] = dr;
0544                 }
0545                 val = (*payload)[id.rawId()].shift_x12;
0546                 if (irun == 0) {
0547                   rEE[0][hashindex] = val;
0548                 } else {
0549                   if (method == 0)
0550                     dr = val - rEE[0][hashindex];
0551                   else {  // ratio
0552                     if (rEE[0][hashindex] == 0.) {
0553                       if (val == 0.)
0554                         dr = 1.;
0555                       else
0556                         dr = 9999.;
0557                     } else
0558                       dr = val / rEE[0][hashindex];
0559                   }
0560                   if (thesign == 1)
0561                     endc_p_r[0]->Fill(ix + 1, iy + 1, dr);
0562                   else
0563                     endc_m_r[0]->Fill(ix + 1, iy + 1, dr);
0564                   if (dr < rEEmin[0])
0565                     rEEmin[0] = dr;
0566                   if (dr > rEEmax[0])
0567                     rEEmax[0] = dr;
0568                 }
0569                 val = (*payload)[id.rawId()].mult_x6;
0570                 if (irun == 0) {
0571                   mEE[1][hashindex] = val;
0572                 } else {
0573                   if (method == 0)
0574                     dr = val - mEE[1][hashindex];
0575                   else {  // ratio
0576                     if (mEE[1][hashindex] == 0.) {
0577                       if (val == 0.)
0578                         dr = 1.;
0579                       else
0580                         dr = 9999.;
0581                     } else
0582                       dr = val / mEE[1][hashindex];
0583                   }
0584                   if (thesign == 1)
0585                     endc_p_m[1]->Fill(ix + 1, iy + 1, dr);
0586                   else
0587                     endc_m_m[1]->Fill(ix + 1, iy + 1, dr);
0588                   if (dr < mEEmin[1])
0589                     mEEmin[1] = dr;
0590                   if (dr > mEEmax[1])
0591                     mEEmax[1] = dr;
0592                 }
0593                 val = (*payload)[id.rawId()].shift_x6;
0594                 if (irun == 0) {
0595                   rEE[1][hashindex] = val;
0596                 } else {
0597                   if (method == 0)
0598                     dr = val - rEE[1][hashindex];
0599                   else {  // ratio
0600                     if (rEE[1][hashindex] == 0.) {
0601                       if (val == 0.)
0602                         dr = 1.;
0603                       else
0604                         dr = 9999.;
0605                     } else
0606                       dr = val / rEE[1][hashindex];
0607                   }
0608                   if (thesign == 1)
0609                     endc_p_r[1]->Fill(ix + 1, iy + 1, dr);
0610                   else
0611                     endc_m_r[1]->Fill(ix + 1, iy + 1, dr);
0612                   if (dr < rEEmin[1])
0613                     rEEmin[1] = dr;
0614                   if (dr > rEEmax[1])
0615                     rEEmax[1] = dr;
0616                 }
0617                 val = (*payload)[id.rawId()].mult_x1;
0618                 if (irun == 0) {
0619                   mEE[2][hashindex] = val;
0620                 } else {
0621                   if (method == 0)
0622                     dr = val - mEE[2][hashindex];
0623                   else {  // ratio
0624                     if (mEE[2][hashindex] == 0.) {
0625                       if (val == 0.)
0626                         dr = 1.;
0627                       else
0628                         dr = 9999.;
0629                     } else
0630                       dr = val / mEE[2][hashindex];
0631                   }
0632                   if (thesign == 1)
0633                     endc_p_m[2]->Fill(ix + 1, iy + 1, dr);
0634                   else
0635                     endc_m_m[2]->Fill(ix + 1, iy + 1, dr);
0636                   if (dr < mEEmin[2])
0637                     mEEmin[2] = dr;
0638                   if (dr > mEEmax[2])
0639                     mEEmax[2] = dr;
0640                 }
0641                 val = (*payload)[id.rawId()].shift_x1;
0642                 if (irun == 0) {
0643                   rEE[2][hashindex] = val;
0644                 } else {
0645                   if (method == 0)
0646                     dr = val - rEE[2][hashindex];
0647                   else {  // ratio
0648                     if (rEE[2][hashindex] == 0.) {
0649                       if (val == 0.)
0650                         dr = 1.;
0651                       else
0652                         dr = 9999.;
0653                     } else
0654                       dr = val / rEE[2][hashindex];
0655                   }
0656                   if (thesign == 1)
0657                     endc_p_r[2]->Fill(ix + 1, iy + 1, dr);
0658                   else
0659                     endc_m_r[2]->Fill(ix + 1, iy + 1, dr);
0660                   if (dr < rEEmin[2])
0661                     rEEmin[2] = dr;
0662                   if (dr > rEEmax[2])
0663                     rEEmax[2] = dr;
0664                 }
0665                 //        fout << " x " << ix << " y " << " dr " << dr << std::endl;
0666               }  // iy
0667             }    // ix
0668           }      // side
0669         }        //  if payload.get()
0670         else
0671           return false;
0672       }  // loop over IOVs
0673 
0674       gStyle->SetPalette(1);
0675       gStyle->SetOptStat(0);
0676       TCanvas canvas("CC map", "CC map", 1200, 1800);
0677       TLatex t1;
0678       t1.SetNDC();
0679       t1.SetTextAlign(26);
0680       int len = l_tagname[0].length() + l_tagname[1].length();
0681       std::string dr[2] = {"-", "/"};
0682       if (ntags == 2) {
0683         if (len < 70) {
0684           t1.SetTextSize(0.03);
0685           t1.DrawLatex(
0686               0.5,
0687               0.96,
0688               Form("%s %i %s %s %i", l_tagname[1].c_str(), run[1], dr[method].c_str(), l_tagname[0].c_str(), run[0]));
0689         } else {
0690           t1.SetTextSize(0.03);
0691           t1.DrawLatex(0.5, 0.96, Form("Ecal TPGLinearizationConst, IOV %i %s %i", run[1], dr[method].c_str(), run[0]));
0692         }
0693       } else {
0694         t1.SetTextSize(0.03);
0695         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]));
0696       }
0697 
0698       float xmi[3] = {0.0, 0.22, 0.78};
0699       float xma[3] = {0.22, 0.78, 1.00};
0700       TPad*** pad = new TPad**[6];
0701       for (int gId = 0; gId < 6; gId++) {
0702         pad[gId] = new TPad*[3];
0703         for (int obj = 0; obj < 3; obj++) {
0704           float yma = 0.94 - (0.16 * gId);
0705           float ymi = yma - 0.14;
0706           pad[gId][obj] = new TPad(Form("p_%i_%i", obj, gId), Form("p_%i_%i", obj, gId), xmi[obj], ymi, xma[obj], yma);
0707           pad[gId][obj]->Draw();
0708         }
0709       }
0710 
0711       for (int gId = 0; gId < kGains; gId++) {
0712         pad[gId][0]->cd();
0713         DrawEE(endc_m_m[gId], mEEmin[gId], mEEmax[gId]);
0714         pad[gId + 3][0]->cd();
0715         DrawEE(endc_m_r[gId], rEEmin[gId], rEEmax[gId]);
0716         pad[gId][1]->cd();
0717         DrawEB(barrel_m[gId], mEBmin[gId], mEBmax[gId]);
0718         pad[gId + 3][1]->cd();
0719         DrawEB(barrel_r[gId], rEBmin[gId], rEBmax[gId]);
0720         pad[gId][2]->cd();
0721         DrawEE(endc_p_m[gId], mEEmin[gId], mEEmax[gId]);
0722         pad[gId + 3][2]->cd();
0723         DrawEE(endc_p_r[gId], rEEmin[gId], rEEmax[gId]);
0724       }
0725 
0726       std::string ImageName(this->m_imageFileName);
0727       canvas.SaveAs(ImageName.c_str());
0728       return true;
0729     }  // fill method
0730   };   // class EcalPedestalsDiffBase
0731   using EcalTPGLinearizationConstDiffOneTag = EcalTPGLinearizationConstBase<cond::payloadInspector::SINGLE_IOV, 1, 0>;
0732   using EcalTPGLinearizationConstDiffTwoTags = EcalTPGLinearizationConstBase<cond::payloadInspector::SINGLE_IOV, 2, 0>;
0733   using EcalTPGLinearizationConstRatioOneTag = EcalTPGLinearizationConstBase<cond::payloadInspector::SINGLE_IOV, 1, 1>;
0734   using EcalTPGLinearizationConstRatioTwoTags = EcalTPGLinearizationConstBase<cond::payloadInspector::SINGLE_IOV, 2, 1>;
0735 
0736 }  // namespace
0737 
0738 // Register the classes as boost python plugin
0739 PAYLOAD_INSPECTOR_MODULE(EcalTPGLinearizationConst) {
0740   PAYLOAD_INSPECTOR_CLASS(EcalTPGLinearizationConstPlot);
0741   PAYLOAD_INSPECTOR_CLASS(EcalTPGLinearizationConstDiffOneTag);
0742   PAYLOAD_INSPECTOR_CLASS(EcalTPGLinearizationConstDiffTwoTags);
0743   PAYLOAD_INSPECTOR_CLASS(EcalTPGLinearizationConstRatioOneTag);
0744   PAYLOAD_INSPECTOR_CLASS(EcalTPGLinearizationConstRatioTwoTags);
0745 }