File indexing completed on 2024-04-06 12:01:36
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 "DataFormats/EcalDetId/interface/EcalSubdetector.h"
0006 #include "DataFormats/EcalDetId/interface/EcalTrigTowerDetId.h"
0007 #include "CondCore/EcalPlugins/plugins/EcalDrawUtils.h"
0008
0009
0010 #include "CondFormats/EcalObjects/interface/EcalTPGFineGrainEBIdMap.h"
0011
0012 #include "TH2F.h"
0013 #include "TCanvas.h"
0014 #include "TStyle.h"
0015 #include "TLine.h"
0016 #include "TLatex.h"
0017
0018 #include <string>
0019 #include <cstring>
0020
0021 namespace {
0022 enum { TEMPLATESAMPLES = 5 };
0023 enum { MIN_IETA = 1, MIN_IPHI = 1, MAX_IETA = 17, MAX_IPHI = 72 };
0024
0025
0026
0027
0028 class EcalTPGFineGrainEBIdMapPlot : public cond::payloadInspector::PlotImage<EcalTPGFineGrainEBIdMap> {
0029 public:
0030 EcalTPGFineGrainEBIdMapPlot()
0031 : cond::payloadInspector::PlotImage<EcalTPGFineGrainEBIdMap>("Ecal TPGFineGrainEBIdMap - map ") {
0032 setSingleIov(true);
0033 }
0034
0035 bool fill(const std::vector<std::tuple<cond::Time_t, cond::Hash> >& iovs) override {
0036 TH2F** barrel = new TH2F*[TEMPLATESAMPLES];
0037 double pEBmin[TEMPLATESAMPLES], pEBmax[TEMPLATESAMPLES];
0038 std::string text[TEMPLATESAMPLES] = {"ThresholdETLow", "ThresholdETHigh", "RatioLow", "RatioHigh", "LUT"};
0039 int EBcnt = 0;
0040
0041 for (int s = 0; s < TEMPLATESAMPLES; ++s) {
0042 char* y = new char[text[s].length() + 1];
0043 std::strcpy(y, text[s].c_str());
0044 barrel[s] = new TH2F("EB", y, MAX_IPHI, 0, MAX_IPHI, 2 * MAX_IETA, -MAX_IETA, MAX_IETA);
0045 }
0046
0047 uint32_t ThresholdETLow = 0;
0048 uint32_t ThresholdETHigh = 0;
0049 uint32_t RatioLow = 0;
0050 uint32_t RatioHigh = 0;
0051 uint32_t LUT = 0;
0052
0053 auto iov = iovs.front();
0054 std::shared_ptr<EcalTPGFineGrainEBIdMap> payload = fetchPayload(std::get<1>(iov));
0055 unsigned int run = std::get<0>(iov);
0056 if (payload.get()) {
0057 const std::map<uint32_t, EcalTPGFineGrainConstEB>& towerMap = (*payload).getMap();
0058 std::map<uint32_t, EcalTPGFineGrainConstEB>::const_iterator it = towerMap.begin();
0059
0060 for (int iphi = MIN_IPHI - 1; iphi <= MAX_IPHI; iphi++) {
0061 for (int ieta = -1 * MAX_IETA; ieta < MAX_IETA; ieta++) {
0062
0063
0064 EcalTPGFineGrainConstEB fg = (*it).second;
0065
0066 fg.getValues(ThresholdETLow, ThresholdETHigh, RatioLow, RatioHigh, LUT);
0067
0068 barrel[0]->Fill(iphi, ieta, ThresholdETLow);
0069 barrel[1]->Fill(iphi, ieta, ThresholdETHigh);
0070 barrel[2]->Fill(iphi, ieta, RatioLow);
0071 barrel[3]->Fill(iphi, ieta, RatioHigh);
0072 barrel[4]->Fill(iphi, ieta, LUT);
0073
0074 if (ThresholdETLow < pEBmin[0])
0075 pEBmin[0] = ThresholdETLow;
0076 if (ThresholdETHigh < pEBmin[1])
0077 pEBmin[1] = ThresholdETHigh;
0078 if (RatioLow < pEBmin[2])
0079 pEBmin[2] = RatioLow;
0080 if (RatioHigh < pEBmin[3])
0081 pEBmin[3] = RatioHigh;
0082 if (LUT < pEBmin[4])
0083 pEBmin[4] = LUT;
0084
0085 if (ThresholdETLow > pEBmax[0])
0086 pEBmax[0] = ThresholdETLow;
0087 if (ThresholdETHigh > pEBmax[1])
0088 pEBmax[1] = ThresholdETHigh;
0089 if (RatioLow > pEBmax[2])
0090 pEBmax[2] = RatioLow;
0091 if (RatioHigh > pEBmax[3])
0092 pEBmax[3] = RatioHigh;
0093 if (LUT > pEBmax[4])
0094 pEBmax[4] = LUT;
0095
0096 EBcnt++;
0097 }
0098 }
0099
0100 }
0101
0102 gStyle->SetPalette(1);
0103 gStyle->SetOptStat(0);
0104 TCanvas canvas("CC map", "CC map", 1600, 2800);
0105 TLatex t1;
0106 t1.SetNDC();
0107 t1.SetTextAlign(26);
0108 t1.SetTextSize(0.04);
0109 t1.DrawLatex(0.5, 0.96, Form("Ecal TPGFine Grain EBIdMap, IOV %i", run));
0110
0111 float xmi = 0.24;
0112 float xma = 0.76;
0113
0114 TPad** pad = new TPad*[TEMPLATESAMPLES];
0115 for (int s = 0; s < TEMPLATESAMPLES; s++) {
0116 float yma = 0.94 - (0.16 * s);
0117 float ymi = yma - 0.14;
0118
0119 char* y = new char[text[s].length() + 1];
0120 std::strcpy(y, text[s].c_str());
0121
0122 pad[s] = new TPad(Form("Towers %i", EBcnt), y, xmi, ymi, xma, yma);
0123 pad[s]->Draw();
0124 }
0125
0126 for (int s = 0; s < TEMPLATESAMPLES; s++) {
0127 pad[s]->cd();
0128
0129 if (pEBmin[s] == pEBmax[s]) {
0130 pEBmin[s] = pEBmin[s] - 1.e-06;
0131 pEBmax[s] = pEBmax[s] + 1.e-06;
0132 }
0133
0134 barrel[s]->SetMaximum(pEBmax[s]);
0135 barrel[s]->SetMinimum(pEBmin[s]);
0136 barrel[s]->Draw("colz");
0137
0138 TLine* l = new TLine(0., 0., 0., 0.);
0139 l->SetLineWidth(1);
0140 for (int i = 0; i < MAX_IETA; i++) {
0141 Double_t x = 4. + (i * 4);
0142 l = new TLine(x, -MAX_IETA, x, MAX_IETA);
0143 l->Draw();
0144 }
0145 }
0146
0147 std::string ImageName(m_imageFileName);
0148 canvas.SaveAs(ImageName.c_str());
0149
0150 return true;
0151 }
0152 };
0153
0154 }
0155
0156
0157 PAYLOAD_INSPECTOR_MODULE(EcalTPGFineGrainEBIdMap) { PAYLOAD_INSPECTOR_CLASS(EcalTPGFineGrainEBIdMapPlot); }