Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:17

0001 #include "DQMServices/Core/interface/DQMEDHarvester.h"
0002 #include "DQMServices/Core/interface/DQMStore.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/ServiceRegistry/interface/Service.h"
0009 #include "FWCore/Framework/interface/MakerMacros.h"
0010 
0011 // for string manipulations
0012 #include <fmt/printf.h>
0013 
0014 class SiPixelTrackComparisonHarvester : public DQMEDHarvester {
0015 public:
0016   explicit SiPixelTrackComparisonHarvester(const edm::ParameterSet&);
0017   ~SiPixelTrackComparisonHarvester() override = default;
0018   void dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) override;
0019   void project2DalongDiagonal(MonitorElement* input2D, DQMStore::IBooker& ibooker);
0020   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0021 
0022 private:
0023   // ----------member data ---------------------------
0024   const std::string topFolder_;
0025 };
0026 
0027 SiPixelTrackComparisonHarvester::SiPixelTrackComparisonHarvester(const edm::ParameterSet& iConfig)
0028     : topFolder_(iConfig.getParameter<std::string>("topFolderName")) {}
0029 
0030 void SiPixelTrackComparisonHarvester::dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) {
0031   MonitorElement* hpt_eta_tkAllCPU = igetter.get(topFolder_ + "/ptetatrkAllCPU");
0032   MonitorElement* hpt_eta_tkAllCPUmatched = igetter.get(topFolder_ + "/ptetatrkAllCPUmatched");
0033   MonitorElement* hphi_z_tkAllCPU = igetter.get(topFolder_ + "/phiztrkAllCPU");
0034   MonitorElement* hphi_z_tkAllCPUmatched = igetter.get(topFolder_ + "/phiztrkAllCPUmatched");
0035 
0036   if (hpt_eta_tkAllCPU == nullptr or hpt_eta_tkAllCPUmatched == nullptr or hphi_z_tkAllCPU == nullptr or
0037       hphi_z_tkAllCPUmatched == nullptr) {
0038     edm::LogError("SiPixelTrackComparisonHarvester")
0039         << "MEs needed for this module are not found in the input file. Skipping.";
0040     return;
0041   }
0042 
0043   ibooker.cd();
0044   ibooker.setCurrentFolder(topFolder_);
0045   MonitorElement* hpt_eta_matchRatio = ibooker.book2D(
0046       "matchingeff_pt_eta", "Efficiency of track matching; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
0047   MonitorElement* hphi_z_matchRatio = ibooker.book2D(
0048       "matchingeff_phi_z", "Efficiency of track matching; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
0049 
0050   hpt_eta_matchRatio->divide(hpt_eta_tkAllCPUmatched, hpt_eta_tkAllCPU, 1., 1., "B");
0051   hphi_z_matchRatio->divide(hphi_z_tkAllCPUmatched, hphi_z_tkAllCPU, 1., 1., "B");
0052 
0053   // now create the 1D projection from the 2D histograms
0054   std::vector<std::string> listOfMEsToProject = {"nTracks",
0055                                                  "nLooseAndAboveTracks",
0056                                                  "nLooseAndAboveTracks_matched",
0057                                                  "nRecHits",
0058                                                  "nLayers",
0059                                                  "nChi2ndof",
0060                                                  "charge",
0061                                                  "pt",
0062                                                  "eta",
0063                                                  "phi",
0064                                                  "z",
0065                                                  "tip"};
0066   for (const auto& me : listOfMEsToProject) {
0067     MonitorElement* input2D = igetter.get(topFolder_ + "/" + me);
0068     this->project2DalongDiagonal(input2D, ibooker);
0069   }
0070 }
0071 
0072 void SiPixelTrackComparisonHarvester::project2DalongDiagonal(MonitorElement* input2D, DQMStore::IBooker& ibooker) {
0073   if (input2D == nullptr) {
0074     edm::LogError("SiPixelTrackComparisonHarvester")
0075         << "MEs needed for diagonal projection are not found in the input file. Skipping.";
0076     return;
0077   }
0078 
0079   ibooker.cd();
0080   ibooker.setCurrentFolder(topFolder_ + "/projectedDifferences");
0081   const auto& h_name = fmt::sprintf("%s_proj", input2D->getName());
0082   const auto& h_title = fmt::sprintf(";%s CPU -GPU difference", input2D->getTitle());
0083   const auto& span = (input2D->getAxisMax() - input2D->getAxisMin());
0084   const auto& b_w = span / input2D->getNbinsX();
0085   const auto& nbins = ((input2D->getNbinsX() % 2) == 0) ? input2D->getNbinsX() + 1 : input2D->getNbinsX();
0086 
0087   MonitorElement* diagonalized = ibooker.book1D(h_name, h_title, nbins, -span / 2., span / 2.);
0088 
0089   // collect all the entry on each diagonal of the 2D histogram
0090   for (int i = 1; i <= input2D->getNbinsX(); i++) {
0091     for (int j = 1; j <= input2D->getNbinsY(); j++) {
0092       diagonalized->Fill((i - j) * b_w, input2D->getBinContent(i, j));
0093     }
0094   }
0095 
0096   // zero the error on the bin as it's sort of meaningless for the way we fill it
0097   // by collecting the entry on each diagonal
0098   for (int bin = 1; bin <= diagonalized->getNbinsX(); bin++) {
0099     diagonalized->setBinError(bin, 0.f);
0100   }
0101 }
0102 
0103 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0104 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0105 void SiPixelTrackComparisonHarvester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0106   edm::ParameterSetDescription desc;
0107   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareGPUvsCPU/");
0108   descriptions.addWithDefaultLabel(desc);
0109 }
0110 
0111 //define this as a plug-in
0112 DEFINE_FWK_MODULE(SiPixelTrackComparisonHarvester);