File indexing completed on 2024-04-06 12:09:05
0001 #include "catch.hpp"
0002 #include <iostream>
0003 #include <numeric> // std::accumulate
0004 #include "TCanvas.h"
0005 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0006 #include "CalibTracker/SiPixelESProducers/interface/SiPixelDetInfoFileReader.h"
0007
0008 static const std::string k_geo = "SLHCUpgradeSimulations/Geometry/data/PhaseI/PixelSkimmedGeometry_phase1.txt";
0009
0010 TEST_CASE("Phase1PixelROCMaps testing", "[Phase1PixelROCMaps]") {
0011
0012 SECTION("Check barrel plotting") {
0013 gStyle->SetOptStat(0);
0014 Phase1PixelROCMaps theMap("");
0015 TCanvas c = TCanvas("c", "c", 1200, 1200);
0016 theMap.drawBarrelMaps(c, "testing barrel");
0017 c.SaveAs("Phase1PixelROCMaps_barrel.png");
0018 REQUIRE(theMap.getLayerMaps().size() == 4);
0019 }
0020
0021
0022 SECTION("Check endcap plotting") {
0023 gStyle->SetOptStat(0);
0024 Phase1PixelROCMaps theMap("");
0025 TCanvas c = TCanvas("c", "c", 1200, 600);
0026 theMap.drawForwardMaps(c, "testing endcaps");
0027 c.SaveAs("Phase1PixelROCMaps_endcap.png");
0028 REQUIRE(theMap.getRingMaps().size() == 2);
0029 }
0030
0031
0032 SECTION("Check whole detector plotting") {
0033 gStyle->SetOptStat(0);
0034 Phase1PixelROCMaps theMap("");
0035 TCanvas c = TCanvas("c", "c", 1200, 1600);
0036 theMap.drawMaps(c, "testing everything");
0037 c.SaveAs("Phase1PixelROCMaps_whole.png");
0038 REQUIRE(theMap.getLayerMaps().size() == 4);
0039 REQUIRE(theMap.getRingMaps().size() == 2);
0040 }
0041
0042
0043 SECTION("Check empty delta plotting") {
0044 gStyle->SetOptStat(0);
0045 Phase1PixelROCMaps theMap("#Delta", "#Delta");
0046 TCanvas c = TCanvas("c", "c", 1200, 1600);
0047 theMap.drawMaps(c, "testing empty #Delta");
0048 theMap.fillWholeModule(303042564, 1.);
0049 theMap.fillWholeModule(344912900, -1.);
0050 c.SaveAs("Phase1PixelROCMaps_emptyDelta.png");
0051 REQUIRE(theMap.getLayerMaps().size() == 4);
0052 REQUIRE(theMap.getRingMaps().size() == 2);
0053 }
0054
0055
0056 SECTION("Check filling whole modules") {
0057 Phase1PixelROCMaps theMap("");
0058 gStyle->SetOptStat(0);
0059 TCanvas c = TCanvas("c", "c", 1200, 1600);
0060 SiPixelDetInfoFileReader reader_ = SiPixelDetInfoFileReader(edm::FileInPath(k_geo).fullPath());
0061 const auto& detIds = reader_.getAllDetIds();
0062 for (const auto& it : detIds) {
0063 int subid = DetId(it).subdetId();
0064 if (subid == PixelSubdetector::PixelBarrel) {
0065 int module = theMap.findDetCoordinates(it).m_s_module;
0066 if (module % 2 == 0) {
0067 theMap.fillWholeModule(it, 1.);
0068 } else {
0069 theMap.fillWholeModule(it, -1.);
0070 }
0071 } else if (subid == PixelSubdetector::PixelEndcap) {
0072 int panel = theMap.findDetCoordinates(it).m_panel;
0073 if (panel % 2 == 0) {
0074 theMap.fillWholeModule(it, 1.);
0075 } else {
0076 theMap.fillWholeModule(it, -1.);
0077 }
0078 }
0079 }
0080 theMap.drawMaps(c, "testing whole modules");
0081 c.SaveAs("Phase1PixelROCMaps_fullModules.png");
0082
0083 int totalEntries = 0;
0084 const auto layerMaps = theMap.getLayerMaps();
0085 const auto ringMaps = theMap.getRingMaps();
0086
0087 int nBarrel = std::accumulate(layerMaps.begin(), layerMaps.end(), 0, [](int a, const std::shared_ptr<TH2> h) {
0088 return a += h.get()->GetEntries();
0089 });
0090
0091 int nForward = std::accumulate(ringMaps.begin(), ringMaps.end(), 0, [](int a, const std::shared_ptr<TH2> h) {
0092 return a += h.get()->GetEntries();
0093 });
0094
0095 totalEntries = nBarrel + nForward;
0096
0097 REQUIRE(totalEntries == (1856 * 16));
0098 }
0099
0100
0101 SECTION("Check filling in delta mode") {
0102 Phase1PixelROCMaps theMap("", "#Delta: flipped vs unflipped");
0103 gStyle->SetOptStat(0);
0104 TCanvas c = TCanvas("c", "c", 1200, 1600);
0105 SiPixelDetInfoFileReader reader_ = SiPixelDetInfoFileReader(edm::FileInPath(k_geo).fullPath());
0106 const auto& detIds = reader_.getAllDetIds();
0107 for (const auto& it : detIds) {
0108 bool isFlipped = theMap.findDetCoordinates(it).isFlipped();
0109 theMap.fillWholeModule(it, isFlipped ? 1. : -1);
0110 }
0111 theMap.drawMaps(c, "testing #Delta mode");
0112 c.SaveAs("Phase1PixelROCMaps_deltaMode.png");
0113
0114 int totalEntries = 0;
0115 const auto layerMaps = theMap.getLayerMaps();
0116 const auto ringMaps = theMap.getRingMaps();
0117
0118 int nBarrel = std::accumulate(layerMaps.begin(), layerMaps.end(), 0, [](int a, const std::shared_ptr<TH2> h) {
0119 return a += h.get()->GetEntries();
0120 });
0121
0122 int nForward = std::accumulate(ringMaps.begin(), ringMaps.end(), 0, [](int a, const std::shared_ptr<TH2> h) {
0123 return a += h.get()->GetEntries();
0124 });
0125
0126 totalEntries = nBarrel + nForward;
0127
0128 REQUIRE(totalEntries == (1856 * 16));
0129 }
0130
0131
0132 SECTION("Check filling selected ROCs") {
0133 Phase1PixelROCMaps theMap("");
0134 gStyle->SetOptStat(0);
0135 TCanvas c = TCanvas("c", "c", 1200, 1600);
0136 SiPixelDetInfoFileReader reader_ = SiPixelDetInfoFileReader(edm::FileInPath(k_geo).fullPath());
0137 const auto& detIds = reader_.getAllDetIds();
0138 for (const auto& it : detIds) {
0139 for (unsigned i = 0; i < 16; i++) {
0140 std::bitset<16> bad_rocs;
0141 bad_rocs.set(i);
0142 theMap.fillSelectedRocs(it, bad_rocs, i);
0143 }
0144
0145
0146 }
0147
0148 int totalEntries = 0;
0149 const auto layerMaps = theMap.getLayerMaps();
0150 const auto ringMaps = theMap.getRingMaps();
0151
0152 int nBarrel = std::accumulate(layerMaps.begin(), layerMaps.end(), 0, [](int a, const std::shared_ptr<TH2> h) {
0153 return a += h.get()->GetEntries();
0154 });
0155
0156 int nForward = std::accumulate(ringMaps.begin(), ringMaps.end(), 0, [](int a, const std::shared_ptr<TH2> h) {
0157 return a += h.get()->GetEntries();
0158 });
0159
0160 totalEntries = nBarrel + nForward;
0161
0162 theMap.drawMaps(c, "testing selected ROCs");
0163 c.SaveAs("Phase1PixelROCMaps_fullROCs.png");
0164
0165 REQUIRE(totalEntries == (1856 * 16));
0166 }
0167 }