Back to home page

Project CMSSW displayed by LXR

 
 

    


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/Phase1PixelMaps.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("Phase1PixelMaps testing", "[Phase1PixelMaps]") {
0011   //_____________________________________________________________
0012   SECTION("Check barrel plotting") {
0013     gStyle->SetOptStat(0);
0014     Phase1PixelMaps theMap("");
0015     TCanvas c = TCanvas("c", "c", 1200, 1200);
0016     theMap.bookBarrelHistograms("mytest", "test", "z-axis");
0017     theMap.drawBarrelMaps("mytest", c, "colz0");
0018     theMap.beautifyAllHistograms();
0019     c.SaveAs("Phase1PixelMaps_Barrel.png");
0020     REQUIRE(true);
0021   }
0022 
0023   //_____________________________________________________________
0024   SECTION("Check endcap plotting") {
0025     gStyle->SetOptStat(0);
0026     Phase1PixelMaps theMap("");
0027     TCanvas c = TCanvas("c", "c", 1200, 800);
0028     theMap.bookForwardHistograms("mytest", "test", "z-axis");
0029     theMap.drawForwardMaps("mytest", c, "colz0");
0030     theMap.beautifyAllHistograms();
0031     c.SaveAs("Phase1PixelMaps_Endcap.png");
0032     REQUIRE(true);
0033   }
0034 
0035   //_____________________________________________________________
0036   SECTION("Check summary plotting") {
0037     gStyle->SetOptStat(0);
0038     Phase1PixelMaps theMap("textAL");  // needed to not show the axis
0039     TCanvas c = TCanvas("c", "c", 1200, 800);
0040     theMap.bookBarrelHistograms("mytest", "test", "z-axis");
0041     theMap.bookForwardHistograms("mytest", "test", "z-axis");
0042     theMap.beautifyAllHistograms();
0043     theMap.drawSummaryMaps("mytest", c);
0044     c.SaveAs("Phase1PixelMaps_Summary.png");
0045     REQUIRE(true);
0046   }
0047 
0048   //_____________________________________________________________
0049   SECTION("Check summary filling") {
0050     gStyle->SetOptStat(0);
0051     Phase1PixelMaps theMap("COLZA L");
0052     TCanvas c = TCanvas("c", "c", 1200, 800);
0053     theMap.bookBarrelHistograms("mytest", "test", "z-axis");
0054     theMap.bookForwardHistograms("mytest", "test", "z-axis");
0055 
0056     SiPixelDetInfoFileReader reader_ = SiPixelDetInfoFileReader(edm::FileInPath(k_geo).fullPath());
0057     const auto& detIds = reader_.getAllDetIds();
0058     int count = 0;
0059     for (const auto& it : detIds) {
0060       count++;
0061       int subid = DetId(it).subdetId();
0062       if (subid == PixelSubdetector::PixelBarrel) {
0063         theMap.fillBarrelBin("mytest", it, count);
0064       } else if (subid == PixelSubdetector::PixelEndcap) {
0065         theMap.fillForwardBin("mytest", it, count);
0066       }
0067     }
0068     theMap.beautifyAllHistograms();
0069     theMap.drawSummaryMaps("mytest", c);
0070     c.SaveAs("Phase1PixelMaps_Summary_Filled.png");
0071     REQUIRE(true);
0072   }
0073 
0074   //_____________________________________________________________
0075   SECTION("Check summary filling V2") {
0076     gStyle->SetOptStat(0);
0077     gStyle->SetPalette(kRainBow);
0078     Phase1PixelMaps theMap("COLZA L");
0079 
0080     TCanvas c = TCanvas("c", "c", 1200, 800);
0081     theMap.book("mytest", "module counts", "module counts");
0082 
0083     SiPixelDetInfoFileReader reader_ = SiPixelDetInfoFileReader(edm::FileInPath(k_geo).fullPath());
0084     const auto& detIds = reader_.getAllDetIds();
0085     int count = 0;
0086     for (const auto& it : detIds) {
0087       count++;
0088       theMap.fill("mytest", it, count);
0089     }
0090 
0091     theMap.beautifyAllHistograms();
0092     theMap.drawSummaryMaps("mytest", c);
0093     c.SaveAs("Phase1PixelMaps_Summary_Filled_V2.png");
0094     REQUIRE(true);
0095   }
0096 
0097   //_____________________________________________________________
0098   SECTION("Check summary filling V3") {
0099     gStyle->SetOptStat(0);
0100     gStyle->SetPalette(kBlackBody);
0101     Phase1PixelMaps theMap("COLZA L");
0102 
0103     TCanvas c = TCanvas("c", "c", 1200, 800);
0104     theMap.book("mytest", "module counts", "module counts");
0105 
0106     SiPixelDetInfoFileReader reader_ = SiPixelDetInfoFileReader(edm::FileInPath(k_geo).fullPath());
0107     const auto& detIds = reader_.getAllDetIds();
0108     int count = 0;
0109     for (const auto& it : detIds) {
0110       count++;
0111       theMap.fill("mytest", it, count);
0112     }
0113 
0114     theMap.setNoRescale();
0115     theMap.beautifyAllHistograms();
0116     theMap.drawSummaryMaps("mytest", c);
0117     theMap.setBarrelScale("mytest", std::make_pair(0., 100.));
0118     theMap.setForwardScale("mytest", std::make_pair(0., 20.));
0119     c.SaveAs("Phase1PixelMaps_Summary_Filled_V3.png");
0120     REQUIRE(true);
0121   }
0122 }