File indexing completed on 2024-09-07 04:35:29
0001 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0002 #include "CondCore/Utilities/interface/PayloadInspector.h"
0003 #include "DataFormats/EcalDetId/interface/EBDetId.h"
0004 #include "DataFormats/EcalDetId/interface/EEDetId.h"
0005 #include "CondCore/EcalPlugins/plugins/EcalDrawUtils.h"
0006
0007
0008 #include "CondFormats/EcalObjects/interface/EcalWeightXtalGroups.h"
0009
0010 #include "TH2F.h"
0011 #include "TCanvas.h"
0012 #include "TStyle.h"
0013 #include "TLine.h"
0014 #include "TLatex.h"
0015
0016 #include <string>
0017 #include <memory>
0018 #include <array>
0019
0020 namespace {
0021
0022 enum { kEBChannels = 61200, kEEChannels = 14648, kSides = 2, kRMS = 5, TEMPLATESAMPLES = 12 };
0023 enum { MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 85, MAX_IPHI = 360 };
0024 enum { IX_MIN = 1, IY_MIN = 1, IX_MAX = 100, IY_MAX = 100 };
0025
0026
0027
0028
0029 class EcalWeightXtalGroupsPlot : public cond::payloadInspector::PlotImage<EcalWeightXtalGroups> {
0030 public:
0031 EcalWeightXtalGroupsPlot()
0032 : cond::payloadInspector::PlotImage<EcalWeightXtalGroups>("Ecal Weight Xtal Groups - map ") {
0033 setSingleIov(true);
0034 }
0035
0036 bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0037 TH2F* barrel = new TH2F("EB", "mean EB", MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
0038 TH2F* endc_p = new TH2F("EE+", "mean EE+", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
0039 TH2F* endc_m = new TH2F("EE-", "mean EE-", IX_MAX, IX_MIN, IX_MAX + 1, IY_MAX, IY_MIN, IY_MAX + 1);
0040
0041 auto iov = iovs.front();
0042 std::shared_ptr<EcalWeightXtalGroups> payload = fetchPayload(std::get<1>(iov));
0043 unsigned int run = std::get<0>(iov);
0044
0045 if (payload.get()) {
0046 if (payload->barrelItems().empty())
0047 return false;
0048
0049 for (int cellid = EBDetId::MIN_HASH; cellid < EBDetId::kSizeForDenseIndexing; ++cellid) {
0050 uint32_t rawid = EBDetId::unhashIndex(cellid);
0051
0052 std::vector<EcalXtalGroupId>::const_iterator value_ptr = payload->find(rawid);
0053
0054 if (value_ptr == payload->end())
0055 continue;
0056
0057 unsigned int weight = (unsigned int)((*value_ptr).id());
0058 Double_t phi = (Double_t)(EBDetId(rawid)).iphi() - 0.5;
0059 Double_t eta = (Double_t)(EBDetId(rawid)).ieta();
0060
0061 if (eta > 0.)
0062 eta = eta - 0.5;
0063 else
0064 eta = eta + 0.5;
0065
0066 barrel->Fill(phi, eta, weight);
0067 }
0068
0069 if (payload->endcapItems().empty())
0070 return false;
0071
0072
0073 for (int iz = -1; iz < 2; iz = iz + 2)
0074 for (int iy = IY_MIN; iy < IY_MAX + IY_MIN; iy++)
0075 for (int ix = IX_MIN; ix < IX_MAX + IX_MIN; ix++)
0076 if (EEDetId::validDetId(ix, iy, iz)) {
0077 EEDetId myEEId = EEDetId(ix, iy, iz, EEDetId::XYMODE);
0078 uint32_t rawid = myEEId.rawId();
0079
0080 std::vector<EcalXtalGroupId>::const_iterator value_ptr = payload->find(rawid);
0081
0082 if (value_ptr == payload->end())
0083 continue;
0084
0085 unsigned int weight = (unsigned int)((*value_ptr).id());
0086
0087 if (iz == 1)
0088 endc_p->Fill(ix, iy, weight);
0089 else
0090 endc_m->Fill(ix, iy, weight);
0091 }
0092 }
0093
0094 gStyle->SetPalette(1);
0095 gStyle->SetOptStat(0);
0096 TCanvas canvas("CC map", "CC map", 1600, 450);
0097 TLatex t1;
0098 t1.SetNDC();
0099 t1.SetTextAlign(26);
0100 t1.SetTextSize(0.05);
0101 t1.DrawLatex(0.5, 0.96, Form("Ecal Weight Xtal Groups, IOV %i", run));
0102
0103 float xmi[3] = {0.0, 0.24, 0.76};
0104 float xma[3] = {0.24, 0.76, 1.00};
0105 std::array<std::unique_ptr<TPad>, 3> pad;
0106 for (int obj = 0; obj < 3; obj++) {
0107 pad[obj] = std::make_unique<TPad>(Form("p_%i", obj), Form("p_%i", obj), xmi[obj], 0.0, xma[obj], 0.94);
0108 pad[obj]->Draw();
0109 }
0110
0111 pad[0]->cd();
0112
0113 DrawEE(endc_m, 0., 2.5);
0114 pad[1]->cd();
0115
0116 DrawEB(barrel, 0., 2.5);
0117 pad[2]->cd();
0118
0119 DrawEE(endc_p, 0., 2.5);
0120
0121 std::string ImageName(m_imageFileName);
0122 canvas.SaveAs(ImageName.c_str());
0123 return true;
0124 }
0125 };
0126
0127 }
0128
0129
0130 PAYLOAD_INSPECTOR_MODULE(EcalWeightXtalGroups) { PAYLOAD_INSPECTOR_CLASS(EcalWeightXtalGroupsPlot); }