Back to home page

Project CMSSW displayed by LXR

 
 

    


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

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