File indexing completed on 2024-04-06 12:01:44
0001 #ifndef CONDCORE_SIPIXELPLUGINS_PIXELREGIONCONTAINERS_H
0002 #define CONDCORE_SIPIXELPLUGINS_PIXELREGIONCONTAINERS_H
0003
0004 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0005 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0006 #include "DataFormats/TrackerCommon/interface/PixelBarrelName.h"
0007 #include "DataFormats/TrackerCommon/interface/PixelEndcapName.h"
0008 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0009 #include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
0010 #include "FWCore/ParameterSet/interface/FileInPath.h"
0011 #include <boost/range/adaptor/indexed.hpp>
0012 #include <cstdlib>
0013 #include "TH1.h"
0014 #include "TLatex.h"
0015
0016 namespace PixelRegions {
0017
0018 inline std::string itoa(int i) {
0019 char temp[20];
0020 sprintf(temp, "%d", i);
0021 return ((std::string)temp);
0022 }
0023
0024
0025
0026
0027
0028
0029
0030 enum PixelId {
0031
0032 L1 = 1100, L2 = 1200, L3 = 1300, L4 = 1400,
0033
0034 Rm1l = 2111, Rm1u = 2112, Rm2l = 2121, Rm2u = 2122, Rm3l = 2131, Rm3u = 2132,
0035
0036 Rp1l = 2211, Rp1u = 2212, Rp2l = 2221, Rp2u = 2222, Rp3l = 2231, Rp3u = 2232,
0037
0038 Ph2EmR1 = 3101, Ph2EmR2 = 3102, Ph2EmR3 = 3103, Ph2EmR4 = 3104,
0039 Ph2EpR1 = 3201, Ph2EpR2 = 3202, Ph2EpR3 = 3203, Ph2EpR4 = 3204,
0040 Ph2FmR1 = 4101, Ph2FmR2 = 4102, Ph2FmR3 = 4103, Ph2FmR4 = 4104, Ph2FmR5 = 4105,
0041 Ph2FpR1 = 4201, Ph2FpR2 = 4202, Ph2FpR3 = 4203, Ph2FpR4 = 4204, Ph2FpR5 = 4205,
0042 End = 99999
0043 };
0044
0045 const std::vector<PixelId> PixelIDs = {
0046
0047 PixelId::L1, PixelId::L2, PixelId::L3, PixelId::L4,
0048
0049 PixelId::Rm1l, PixelId::Rm1u, PixelId::Rm2l, PixelId::Rm2u, PixelId::Rm3l, PixelId::Rm3u,
0050
0051 PixelId::Rp1l, PixelId::Rp1u, PixelId::Rp2l, PixelId::Rp2u, PixelId::Rp3l, PixelId::Rp3u,
0052
0053 PixelId::Ph2EmR1, PixelId::Ph2EmR2, PixelId::Ph2EmR3, PixelId::Ph2EmR4,
0054 PixelId::Ph2EpR1, PixelId::Ph2EpR2, PixelId::Ph2EpR3, PixelId::Ph2EpR4,
0055 PixelId::Ph2FmR1, PixelId::Ph2FmR2, PixelId::Ph2FmR3, PixelId::Ph2FmR4, PixelId::Ph2FmR5,
0056 PixelId::Ph2FpR1, PixelId::Ph2FpR2, PixelId::Ph2FpR3, PixelId::Ph2FpR4, PixelId::Ph2FpR5,
0057 PixelId::End
0058 };
0059
0060
0061
0062 const std::vector<std::string> IDlabels = {"Barrel Pixel L1",
0063 "Barrel Pixel L2",
0064 "Barrel Pixel L3",
0065 "Barrel Pixel L4",
0066 "FPIX(-) Disk 1 inner ring",
0067 "FPIX(-) Disk 1 outer ring",
0068 "FPIX(-) Disk 2 inner ring",
0069 "FPIX(-) Disk 2 outer ring",
0070 "FPIX(-) Disk 3 inner ring",
0071 "FPIX(-) Disk 3 outer ring",
0072 "FPIX(+) Disk 1 inner ring",
0073 "FPIX(+) Disk 1 outer ring",
0074 "FPIX(+) Disk 2 inner ring",
0075 "FPIX(+) Disk 2 outer ring",
0076 "FPIX(+) Disk 3 inner ring",
0077 "FPIX(+) Disk 3 outer ring",
0078 "EPix(-) Ring 1",
0079 "EPix(-) Ring 2",
0080 "EPix(-) Ring 3",
0081 "EPix(-) Ring 4",
0082 "EPix(+) Ring 1",
0083 "EPix(+) Ring 2",
0084 "EPix(+) Ring 3",
0085 "EPix(+) Ring 4",
0086 "FPix(-) Ring 1",
0087 "FPix(-) Ring 2",
0088 "FPix(-) Ring 3",
0089 "FPix(-) Ring 4",
0090 "FPix(-) Ring 5",
0091 "FPix(+) Ring 1",
0092 "FPix(+) Ring 2",
0093 "FPix(+) Ring 3",
0094 "FPix(+) Ring 4",
0095 "FPix(+) Ring 5",
0096 "END"};
0097
0098
0099 [[maybe_unused]] static const std::vector<std::string> getIDLabels(const SiPixelPI::phase& ph, bool isBarrel) {
0100 std::vector<std::string> out;
0101 for (const auto& label : IDlabels) {
0102 if (isBarrel) {
0103 if (label.find("Barrel") != std::string::npos) {
0104 out.push_back(label);
0105 }
0106 } else {
0107 if (ph == SiPixelPI::phase::two) {
0108 if (label.find("Ring") != std::string::npos) {
0109 out.push_back(label);
0110 }
0111 } else {
0112 if (label.find("ring") != std::string::npos) {
0113 out.push_back(label);
0114 }
0115 }
0116 }
0117 }
0118 return out;
0119 }
0120
0121
0122 static const PixelId calculateBPixID(const unsigned int layer) {
0123
0124 PixelId bpixLayer = static_cast<PixelId>(1000 + 100 * layer);
0125 return bpixLayer;
0126 }
0127
0128
0129 static const PixelId calculateFPixID(const SiPixelPI::phase& ph,
0130 const unsigned int side,
0131 const unsigned int disk,
0132 const unsigned int ring) {
0133
0134 using namespace SiPixelPI;
0135 unsigned int prefix(2000);
0136 unsigned int disk_(0);
0137 if (ph > phase::one) {
0138 if (disk < 9) {
0139 prefix += 1000;
0140 } else {
0141 prefix += 2000;
0142 }
0143 } else {
0144 disk_ = disk;
0145 }
0146 PixelId fpixRing = static_cast<PixelId>(prefix + 100 * side + 10 * disk_ + ring);
0147 return fpixRing;
0148 }
0149
0150
0151 static const PixelId detIdToPixelId(const unsigned int detid,
0152 const TrackerTopology* trackTopo,
0153 const SiPixelPI::phase& ph) {
0154 using namespace SiPixelPI;
0155 DetId detId = DetId(detid);
0156 unsigned int subid = detId.subdetId();
0157 unsigned int pixid = 0;
0158 if (subid == PixelSubdetector::PixelBarrel) {
0159 int layer = trackTopo->pxbLayer(detId);
0160 pixid = calculateBPixID(layer);
0161 } else if (subid == PixelSubdetector::PixelEndcap) {
0162 int side = trackTopo->pxfSide(detId);
0163 int disk = trackTopo->pxfDisk(detId);
0164 int ring(0);
0165 switch (ph) {
0166 case phase::zero:
0167 ring = SiPixelPI::ring(detid, *trackTopo, false);
0168 break;
0169 case phase::one:
0170 ring = SiPixelPI::ring(detid, *trackTopo, true);
0171 break;
0172 case phase::two:
0173
0174
0175 ring = trackTopo->pxfBlade(detId);
0176 break;
0177 default:
0178 throw cms::Exception("LogicalError") << " there is not such phase as " << ph;
0179 }
0180 pixid = calculateFPixID(ph, side, disk, ring);
0181 }
0182 PixelId pixID = static_cast<PixelId>(pixid);
0183 return pixID;
0184 }
0185
0186
0187 [[maybe_unused]] static const std::vector<uint32_t> attachedDets(const PixelRegions::PixelId theId,
0188 const TrackerTopology* trackTopo,
0189 const SiPixelPI::phase& ph) {
0190 using namespace SiPixelPI;
0191 std::vector<uint32_t> out = {};
0192 edm::FileInPath m_fp;
0193
0194 switch (ph) {
0195 case phase::zero:
0196
0197 m_fp = edm::FileInPath("CalibTracker/SiPixelESProducers/data/PixelSkimmedGeometry.txt");
0198 break;
0199 case phase::one:
0200
0201 m_fp = edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt");
0202 break;
0203 case phase::two:
0204 m_fp = edm::FileInPath("SLHCUpgradeSimulations/Geometry/data/PhaseII/Tilted/PixelSkimmedGeometryT14.txt");
0205 break;
0206 default:
0207 throw cms::Exception("LogicalError") << " there is not such phase as " << ph;
0208 }
0209
0210 SiPixelDetInfoFileReader pxlreader(m_fp.fullPath());
0211 const std::vector<uint32_t>& pxldetids = pxlreader.getAllDetIds();
0212 for (const auto& d : pxldetids) {
0213 auto ID = detIdToPixelId(d, trackTopo, ph);
0214 if (ID == theId) {
0215 out.push_back(d);
0216 }
0217 }
0218
0219
0220
0221
0222 if (out.empty()) {
0223 out.push_back(0xFFFFFFFF);
0224 }
0225
0226 COUT << "ID:" << theId << " ";
0227 for (const auto& entry : out) {
0228 COUT << entry << ",";
0229 }
0230 COUT << std::endl;
0231
0232 return out;
0233 }
0234
0235
0236
0237
0238 class PixelRegionContainers {
0239 public:
0240 PixelRegionContainers(const TrackerTopology* t_topo, const SiPixelPI::phase& ph)
0241 : m_trackerTopo(t_topo), m_Phase(ph) {
0242
0243 m_isLog = false;
0244 if (m_trackerTopo) {
0245 m_isTrackerTopologySet = true;
0246 } else {
0247 m_isTrackerTopologySet = false;
0248 }
0249 }
0250
0251 ~PixelRegionContainers() {}
0252
0253
0254 void setTheTopo(const TrackerTopology* t_topo) {
0255 m_trackerTopo = t_topo;
0256 m_isTrackerTopologySet = true;
0257 }
0258
0259
0260 const TrackerTopology* getTheTopo() { return m_trackerTopo; }
0261
0262
0263 void bookAll(std::string title_label,
0264 std::string x_label,
0265 std::string y_label,
0266 const int nbins,
0267 const float xmin,
0268 const float xmax) {
0269 using namespace SiPixelPI;
0270 for (const auto& pixelId : PixelIDs | boost::adaptors::indexed(0)) {
0271 if (m_Phase == phase::two &&
0272 pixelId.value() > PixelId::L4 &&
0273 pixelId.value() < PixelId::Ph2EmR1
0274 ) {
0275 continue;
0276 }
0277
0278 m_theMap[pixelId.value()] = std::make_shared<TH1F>((title_label + itoa(pixelId.value())).c_str(),
0279 Form("%s %s;%s;%s",
0280 (IDlabels.at(pixelId.index())).c_str(),
0281 title_label.c_str(),
0282 x_label.c_str(),
0283 y_label.c_str()),
0284 nbins,
0285 xmin,
0286 xmax);
0287 }
0288 }
0289
0290
0291 void fill(const unsigned int detid, const float value) {
0292
0293 assert(m_trackerTopo);
0294
0295
0296 PixelRegions::PixelId myId = PixelRegions::detIdToPixelId(detid, m_trackerTopo, m_Phase);
0297 if (m_theMap.find(myId) != m_theMap.end()) {
0298 m_theMap[myId]->Fill(value);
0299 } else {
0300 throw cms::Exception("LogicError") << detid << " :=> " << myId << " is not a recongnized PixelId enumerator! \n"
0301 << m_trackerTopo->print(detid);
0302 }
0303 }
0304
0305
0306 void draw(TCanvas& canv, bool isBarrel, const char* option = "bar2", bool isPhase1Comparison = false) {
0307 using namespace SiPixelPI;
0308 if (isBarrel) {
0309 if (canv.GetListOfPrimitives()->GetSize() == 0) {
0310
0311 canv.Divide(2, 2);
0312 }
0313 for (int j = 1; j <= 4; j++) {
0314 if (!m_isLog) {
0315 canv.cd(j);
0316 } else {
0317 canv.cd(j)->SetLogy();
0318 }
0319 if ((j == 4) && (m_Phase < 1) && !isPhase1Comparison) {
0320 m_theMap.at(PixelIDs[j - 1])->Draw("AXIS");
0321 TLatex t2;
0322 t2.SetTextAlign(22);
0323 t2.SetTextSize(0.1);
0324 t2.SetTextAngle(45);
0325 t2.SetTextFont(61);
0326 t2.SetTextColor(kBlack);
0327 t2.DrawLatexNDC(0.5, 0.5, "Not in Phase-0!");
0328 } else {
0329 m_theMap.at(PixelIDs[j - 1])->Draw(option);
0330 }
0331 }
0332 } else {
0333
0334
0335 const std::array<int, 18> phase2Pattern = {{1, 2, 3, 4,
0336 6, 7, 8, 9,
0337 11, 12, 13, 14, 15,
0338 16, 17, 18, 19, 20}};
0339
0340
0341 if (canv.GetListOfPrimitives()->GetSize() == 0) {
0342 canv.Divide(m_Phase == phase::two ? 5 : 4, m_Phase == phase::two ? 4 : 3);
0343 }
0344 const int maxPads = (m_Phase == phase::two) ? 18 : 12;
0345 for (int j = 1; j <= maxPads; j++) {
0346 unsigned int mapIndex = m_Phase == 2 ? j + 15 : j + 3;
0347 if (!m_isLog) {
0348 canv.cd((m_Phase == phase::two) ? phase2Pattern[j - 1] : j);
0349
0350 } else {
0351 canv.cd(j)->SetLogy();
0352 }
0353 if ((j % 6 == 5 || j % 6 == 0) && (m_Phase < 1) && !isPhase1Comparison) {
0354 m_theMap.at(PixelIDs[mapIndex])->Draw("AXIS");
0355 TLatex t2;
0356 t2.SetTextAlign(22);
0357 t2.SetTextSize(0.1);
0358 t2.SetTextAngle(45);
0359 t2.SetTextFont(61);
0360 t2.SetTextColor(kBlack);
0361 t2.DrawLatexNDC(0.5, 0.5, "Not in Phase-0!");
0362 } else {
0363 m_theMap.at(PixelIDs[mapIndex])->Draw(option);
0364 }
0365 }
0366 }
0367 }
0368
0369
0370 void beautify(const int linecolor = kBlack, const int fillcolor = kRed) {
0371 for (const auto& plot : m_theMap) {
0372 plot.second->SetTitle("");
0373 if (!m_isLog) {
0374 plot.second->GetYaxis()->SetRangeUser(0., plot.second->GetMaximum() * 1.30);
0375 } else {
0376 plot.second->GetYaxis()->SetRangeUser(0.1, plot.second->GetMaximum() * 100.);
0377 }
0378 plot.second->SetLineColor(linecolor);
0379 if (fillcolor > 0) {
0380 plot.second->SetFillColor(fillcolor);
0381 }
0382 plot.second->SetMarkerStyle(20);
0383 plot.second->SetMarkerSize(1);
0384 SiPixelPI::makeNicePlotStyle(plot.second.get());
0385 plot.second->SetStats(true);
0386 }
0387 }
0388
0389
0390 void setLogScale() { m_isLog = true; }
0391
0392
0393 void stats(int index = 0) {
0394 for (const auto& plot : m_theMap) {
0395 TPaveStats* st = (TPaveStats*)plot.second->FindObject("stats");
0396 if (st) {
0397 st->SetTextSize(0.03);
0398 st->SetLineColor(10);
0399 if (plot.second->GetFillColor() != 0) {
0400 st->SetTextColor(plot.second->GetFillColor());
0401 } else {
0402 st->SetTextColor(plot.second->GetLineColor());
0403 }
0404 SiPixelPI::adjustStats(st, 0.13, 0.85 - index * 0.08, 0.36, 0.93 - index * 0.08);
0405 }
0406 }
0407 }
0408
0409
0410 std::shared_ptr<TH1F>& getHistoFromMap(const PixelRegions::PixelId& theId) {
0411 auto it = m_theMap.find(theId);
0412 if (it != m_theMap.end()) {
0413 return it->second;
0414 } else {
0415 throw cms::Exception("LogicError")
0416 << "PixelRegionContainer::getHistoFromMap(): No histogram is available for PixelId: " << theId << "\n";
0417 }
0418 }
0419
0420
0421 void rescaleMax(PixelRegionContainers& the2ndContainer) {
0422 for (const auto& plot : m_theMap) {
0423 auto thePixId = plot.first;
0424 auto extrema = SiPixelPI::getExtrema((plot.second).get(), the2ndContainer.getHistoFromMap(thePixId).get());
0425 plot.second->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
0426 the2ndContainer.getHistoFromMap(thePixId)->GetYaxis()->SetRangeUser(extrema.first, extrema.second);
0427 }
0428 }
0429
0430 private:
0431 const TrackerTopology* m_trackerTopo;
0432 bool m_isTrackerTopologySet;
0433 SiPixelPI::phase m_Phase;
0434 std::map<PixelId, std::shared_ptr<TH1F>> m_theMap;
0435 int m_nbins;
0436 float m_xim, m_xmax;
0437 bool m_isLog;
0438 };
0439 }
0440 #endif