Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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 #include "CondCore/EcalPlugins/plugins/EcalFloatCondObjectContainerUtils.h"
0008 // the data format of the condition to be inspected
0009 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.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 #include <array>
0020 
0021 namespace {
0022   enum { kEBChannels = 61200, kEEChannels = 14648 };
0023   enum {
0024     MIN_IETA = 1,
0025     MIN_IPHI = 1,
0026     MAX_IETA = 85,
0027     MAX_IPHI = 360,
0028     EBhistEtaMax = 171
0029   };  // barrel lower and upper bounds on eta and phi
0030   enum {
0031     IX_MIN = 1,
0032     IY_MIN = 1,
0033     IX_MAX = 100,
0034     IY_MAX = 100,
0035     EEhistXMax = 220
0036   };  // endcaps lower and upper bounds on x and y
0037 
0038   /**************************************************************
0039      2d histogram of ECAL barrel Intercalib Constants of 1 IOV 
0040   ***************************************************************/
0041 
0042   // inherit from one of the predefined plot class: Histogram2D
0043   class EcalIntercalibConstantsEBMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstants> {
0044   public:
0045     EcalIntercalibConstantsEBMap()
0046         : cond::payloadInspector::Histogram2D<EcalIntercalibConstants>("ECAL Barrel Intercalib Constants - map ",
0047                                                                        "iphi",
0048                                                                        MAX_IPHI,
0049                                                                        MIN_IPHI,
0050                                                                        MAX_IPHI + 1,
0051                                                                        "ieta",
0052                                                                        EBhistEtaMax,
0053                                                                        -MAX_IETA,
0054                                                                        MAX_IETA + 1) {
0055       Base::setSingleIov(true);
0056     }
0057 
0058     // Histogram2D::fill (virtual) needs be overridden - the implementation should use fillWithValue
0059     bool fill() override {
0060       auto tag = PlotBase::getTag<0>();
0061       for (auto const& iov : tag.iovs) {
0062         std::shared_ptr<EcalIntercalibConstants> payload = Base::fetchPayload(std::get<1>(iov));
0063         if (payload.get()) {
0064           // looping over the EB channels, via the dense-index, mapped into EBDetId's
0065           if (payload->barrelItems().empty())
0066             return false;
0067           // set to 1 for ieta 0 (no crystal)
0068           for (int iphi = MIN_IPHI; iphi < MAX_IPHI + 1; iphi++)
0069             fillWithValue(iphi, 0, 1);
0070 
0071           for (int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
0072             uint32_t rawid = EBDetId::unhashIndex(cellid);
0073 
0074             // check the existence of ECAL Intercalib Constants, for a given ECAL barrel channel
0075             EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
0076             if (value_ptr == payload->end())
0077               continue;  // cell absent from payload
0078             float weight = (float)(*value_ptr);
0079 
0080             // fill the Histogram2D here
0081             fillWithValue((EBDetId(rawid)).iphi(), (EBDetId(rawid)).ieta(), weight);
0082           }  // loop over cellid
0083         }    // if payload.get()
0084       }      // loop over IOV's (1 in this case)
0085 
0086       return true;
0087     }  // fill method
0088   };
0089 
0090   class EcalIntercalibConstantsEEMap : public cond::payloadInspector::Histogram2D<EcalIntercalibConstants> {
0091   private:
0092     int EEhistSplit = 20;
0093 
0094   public:
0095     EcalIntercalibConstantsEEMap()
0096         : cond::payloadInspector::Histogram2D<EcalIntercalibConstants>("ECAL Endcap Intercalib Constants - map ",
0097                                                                        "ix",
0098                                                                        EEhistXMax,
0099                                                                        IX_MIN,
0100                                                                        EEhistXMax + 1,
0101                                                                        "iy",
0102                                                                        IY_MAX,
0103                                                                        IY_MIN,
0104                                                                        IY_MAX + 1) {
0105       Base::setSingleIov(true);
0106     }
0107 
0108     bool fill() override {
0109       auto tag = PlotBase::getTag<0>();
0110       for (auto const& iov : tag.iovs) {
0111         std::shared_ptr<EcalIntercalibConstants> payload = Base::fetchPayload(std::get<1>(iov));
0112         if (payload.get()) {
0113           if (payload->endcapItems().empty())
0114             return false;
0115 
0116           // set to 0 everywhwere
0117           for (int ix = IX_MIN; ix < EEhistXMax + 1; ix++)
0118             for (int iy = IY_MIN; iy < IY_MAX + 1; iy++)
0119               fillWithValue(ix, iy, 0);
0120 
0121           for (int cellid = 0; cellid < EEDetId::kSizeForDenseIndexing; ++cellid) {  // loop on EE cells
0122             if (EEDetId::validHashIndex(cellid)) {
0123               uint32_t rawid = EEDetId::unhashIndex(cellid);
0124               EcalFloatCondObjectContainer::const_iterator value_ptr = payload->find(rawid);
0125               if (value_ptr == payload->end())
0126                 continue;  // cell absent from payload
0127               float weight = (float)(*value_ptr);
0128               EEDetId myEEId(rawid);
0129               if (myEEId.zside() == -1)
0130                 fillWithValue(myEEId.ix(), myEEId.iy(), weight);
0131               else
0132                 fillWithValue(myEEId.ix() + IX_MAX + EEhistSplit, myEEId.iy(), weight);
0133             }  // validDetId
0134           }    // loop over cellid
0135         }      // payload
0136       }        // loop over IOV's (1 in this case)
0137       return true;
0138     }  // fill method
0139   };
0140 
0141   /*************************************************
0142      2d plot of ECAL IntercalibConstants of 1 IOV
0143   *************************************************/
0144   class EcalIntercalibConstantsPlot : public cond::payloadInspector::PlotImage<EcalIntercalibConstants> {
0145   public:
0146     EcalIntercalibConstantsPlot()
0147         : cond::payloadInspector::PlotImage<EcalIntercalibConstants>("ECAL Intercalib Constants - map ") {
0148       setSingleIov(true);
0149     }
0150 
0151     bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0152       TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
0153       TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
0154       TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
0155 
0156       auto iov = iovs.front();
0157       std::shared_ptr<EcalIntercalibConstants> payload = fetchPayload(std::get<1>(iov));
0158       unsigned int run = std::get<0>(iov);
0159 
0160       if (payload.get()) {
0161         if (payload->barrelItems().empty())
0162           return false;
0163 
0164         fillEBMap_SingleIOV<EcalIntercalibConstants>(payload, barrel);
0165 
0166         if (payload->endcapItems().empty())
0167           return false;
0168 
0169         fillEEMap_SingleIOV<EcalIntercalibConstants>(payload, endc_m, endc_p);
0170 
0171       }  // payload
0172 
0173       gStyle->SetPalette(1);
0174       gStyle->SetOptStat(0);
0175       TCanvas canvas("CC map", "CC map", 1600, 450);
0176       TLatex t1;
0177       t1.SetNDC();
0178       t1.SetTextAlign(26);
0179       t1.SetTextSize(0.05);
0180       t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants, IOV %i", run));
0181 
0182       float xmi[3] = {0.0, 0.24, 0.76};
0183       float xma[3] = {0.24, 0.76, 1.00};
0184       std::array<std::unique_ptr<TPad>, 3> pad;
0185       for (int obj = 0; obj < 3; obj++) {
0186         pad[obj] = std::make_unique<TPad>(Form("p_%i", obj), Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
0187         pad[obj]->Draw();
0188       }
0189       //      EcalDrawMaps ICMap;
0190       pad[0]->cd();
0191       //      ICMap.DrawEE(endc_m, 0., 2.);
0192       DrawEE(endc_m, 0., 2.5);
0193       pad[1]->cd();
0194       //      ICMap.DrawEB(barrel, 0., 2.);
0195       DrawEB(barrel, 0., 2.5);
0196       pad[2]->cd();
0197       //      ICMap.DrawEE(endc_p, 0., 2.);
0198       DrawEE(endc_p, 0., 2.5);
0199 
0200       std::string ImageName(m_imageFileName);
0201       canvas.SaveAs(ImageName.c_str());
0202       return true;
0203     }  // fill method
0204   };
0205 
0206   /**************************************************************************
0207      2d plot of ECAL IntercalibConstants difference or ratio between 2 IOVs
0208   ***************************************************************************/
0209   template <cond::payloadInspector::IOVMultiplicity nIOVs, int ntags, int method>
0210   class EcalIntercalibConstantsBase : public cond::payloadInspector::PlotImage<EcalIntercalibConstants, nIOVs, ntags> {
0211   public:
0212     EcalIntercalibConstantsBase()
0213         : cond::payloadInspector::PlotImage<EcalIntercalibConstants, nIOVs, ntags>(
0214               "ECAL Intercalib Constants comparison") {}
0215 
0216     bool fill() override {
0217       TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
0218       TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
0219       TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
0220       float pEBmin, pEEmin, pEBmax, pEEmax;
0221       pEBmin = 10.;
0222       pEEmin = 10.;
0223       pEBmax = -10.;
0224       pEEmax = -10.;
0225 
0226       unsigned int run[2];
0227       float pEB[kEBChannels], pEE[kEEChannels];
0228       std::string l_tagname[2];
0229       auto iovs = cond::payloadInspector::PlotBase::getTag<0>().iovs;
0230       l_tagname[0] = cond::payloadInspector::PlotBase::getTag<0>().name;
0231       auto firstiov = iovs.front();
0232       run[0] = std::get<0>(firstiov);
0233       std::tuple<cond::Time_t, cond::Hash> lastiov;
0234       if (ntags == 2) {
0235         auto tag2iovs = cond::payloadInspector::PlotBase::getTag<1>().iovs;
0236         l_tagname[1] = cond::payloadInspector::PlotBase::getTag<1>().name;
0237         lastiov = tag2iovs.front();
0238       } else {
0239         lastiov = iovs.back();
0240         l_tagname[1] = l_tagname[0];
0241       }
0242       run[1] = std::get<0>(lastiov);
0243       for (int irun = 0; irun < nIOVs; irun++) {
0244         std::shared_ptr<EcalIntercalibConstants> payload;
0245         if (irun == 0) {
0246           payload = this->fetchPayload(std::get<1>(firstiov));
0247         } else {
0248           payload = this->fetchPayload(std::get<1>(lastiov));
0249         }
0250         if (payload.get()) {
0251           if (payload->barrelItems().empty())
0252             return false;
0253           fillEBMap_TwoIOVs<EcalIntercalibConstants>(payload, barrel, irun, pEB, pEBmin, pEBmax, method);
0254 
0255           if (payload->endcapItems().empty())
0256             return false;
0257           fillEEMap_TwoIOVs<EcalIntercalibConstants>(payload, endc_m, endc_p, irun, pEE, pEEmin, pEEmax, method);
0258         }  // payload
0259       }    // loop over IOVs
0260 
0261       gStyle->SetPalette(1);
0262       gStyle->SetOptStat(0);
0263       TCanvas canvas("CC map", "CC map", 1600, 450);
0264       TLatex t1;
0265       t1.SetNDC();
0266       t1.SetTextAlign(26);
0267       int len = l_tagname[0].length() + l_tagname[1].length();
0268       std::string dr[2] = {"-", "/"};
0269       if (ntags == 2) {
0270         if (len < 170) {
0271           t1.SetTextSize(0.05);
0272           t1.DrawLatex(0.5,
0273                        0.96,
0274                        Form("%s IOV %i %s %s  IOV %i",
0275                             l_tagname[1].c_str(),
0276                             run[1],
0277                             dr[method].c_str(),
0278                             l_tagname[0].c_str(),
0279                             run[0]));
0280         } else {
0281           t1.SetTextSize(0.05);
0282           t1.DrawLatex(0.5, 0.96, Form("Ecal IntercalibConstants, IOV %i %s %i", run[1], dr[method].c_str(), run[0]));
0283         }
0284       } else {
0285         t1.SetTextSize(0.05);
0286         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]));
0287       }
0288       float xmi[3] = {0.0, 0.24, 0.76};
0289       float xma[3] = {0.24, 0.76, 1.00};
0290       std::array<std::unique_ptr<TPad>, 3> pad;
0291       for (int obj = 0; obj < 3; obj++) {
0292         pad[obj] = std::make_unique<TPad>(Form("p_%i", obj), Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
0293         pad[obj]->Draw();
0294       }
0295       pad[0]->cd();
0296       DrawEE(endc_m, pEEmin, pEEmax);
0297       pad[1]->cd();
0298       DrawEB(barrel, pEBmin, pEBmax);
0299       pad[2]->cd();
0300       DrawEE(endc_p, pEEmin, pEEmax);
0301 
0302       std::string ImageName(this->m_imageFileName);
0303       canvas.SaveAs(ImageName.c_str());
0304       return true;
0305     }  // fill method
0306   };   // class EcalIntercalibConstantsDiffBase
0307   using EcalIntercalibConstantsDiffOneTag = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 1, 0>;
0308   using EcalIntercalibConstantsDiffTwoTags = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 2, 0>;
0309   using EcalIntercalibConstantsRatioOneTag = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 1, 1>;
0310   using EcalIntercalibConstantsRatioTwoTags = EcalIntercalibConstantsBase<cond::payloadInspector::SINGLE_IOV, 2, 1>;
0311 
0312   /*********************************************************
0313      2d plot of Ecal Intercalib Constants Summary of 1 IOV
0314    *********************************************************/
0315   class EcalIntercalibConstantsSummaryPlot : public cond::payloadInspector::PlotImage<EcalIntercalibConstants> {
0316   public:
0317     EcalIntercalibConstantsSummaryPlot()
0318         : cond::payloadInspector::PlotImage<EcalIntercalibConstants>("Ecal Intercalib Constants Summary - map ") {
0319       setSingleIov(true);
0320     }
0321 
0322     bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0323       auto iov = iovs.front();
0324       std::shared_ptr<EcalIntercalibConstants> payload = fetchPayload(std::get<1>(iov));
0325       unsigned int run = std::get<0>(iov);
0326       TH2F* align;
0327       int NbRows;
0328 
0329       if (payload.get()) {
0330         NbRows = 2;
0331         align = new TH2F("", "", 0, 0, 0, 0, 0, 0);
0332 
0333         float mean_x_EB = 0.0f;
0334         float mean_x_EE = 0.0f;
0335 
0336         float rms_EB = 0.0f;
0337         float rms_EE = 0.0f;
0338 
0339         int num_x_EB = 0;
0340         int num_x_EE = 0;
0341 
0342         payload->summary(mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
0343         fillTableWithSummary(
0344             align, "Ecal Intercalib Constants", mean_x_EB, rms_EB, num_x_EB, mean_x_EE, rms_EE, num_x_EE);
0345       } else
0346         return false;
0347 
0348       gStyle->SetPalette(1);
0349       gStyle->SetOptStat(0);
0350       TCanvas canvas("CC map", "CC map", 1000, 1000);
0351       TLatex t1;
0352       t1.SetNDC();
0353       t1.SetTextAlign(26);
0354       t1.SetTextSize(0.04);
0355       t1.SetTextColor(2);
0356       t1.DrawLatex(0.5, 0.96, Form("Ecal Intercalib Constants Summary, IOV %i", run));
0357 
0358       TPad pad("pad", "pad", 0.0, 0.0, 1.0, 0.94);
0359       pad.Draw();
0360       pad.cd();
0361       align->Draw("TEXT");
0362       TLine* l = new TLine;
0363       l->SetLineWidth(1);
0364 
0365       drawTable(NbRows, 4);
0366 
0367       std::string ImageName(m_imageFileName);
0368       canvas.SaveAs(ImageName.c_str());
0369 
0370       return true;
0371     }
0372   };
0373 
0374 }  // namespace
0375 
0376 // Register the classes as boost python plugin
0377 PAYLOAD_INSPECTOR_MODULE(EcalIntercalibConstants) {
0378   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsEBMap);
0379   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsEEMap);
0380   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsPlot);
0381   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsDiffOneTag);
0382   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsDiffTwoTags);
0383   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsRatioOneTag);
0384   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsRatioTwoTags);
0385   PAYLOAD_INSPECTOR_CLASS(EcalIntercalibConstantsSummaryPlot);
0386 }