Back to home page

Project CMSSW displayed by LXR

 
 

    


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