Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-11-18 03:07:17

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