CMS 3D CMS Logo

SiPixelTrackComparisonHarvester.cc
Go to the documentation of this file.
10 
11 // for string manipulations
12 #include <fmt/printf.h>
13 
15 public:
17  ~SiPixelTrackComparisonHarvester() override = default;
18  void dqmEndJob(DQMStore::IBooker& ibooker, DQMStore::IGetter& igetter) override;
20  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
21 
22 private:
23  // ----------member data ---------------------------
25 };
26 
28  : topFolder_(iConfig.getParameter<std::string>("topFolderName")) {}
29 
31  MonitorElement* hpt_eta_tkAllReference = igetter.get(topFolder_ + "/ptetatrkAllReference");
32  if (hpt_eta_tkAllReference == nullptr) {
33  edm::LogError("SiPixelTrackComparisonHarvester")
34  << "MonitorElement not found: " << topFolder_ << "/ptetatrkAllReference. Skipping.";
35  return;
36  }
37 
38  MonitorElement* hpt_eta_tkAllReferencematched = igetter.get(topFolder_ + "/ptetatrkAllReferencematched");
39  if (hpt_eta_tkAllReferencematched == nullptr) {
40  edm::LogError("SiPixelTrackComparisonHarvester")
41  << "MonitorElement not found: " << topFolder_ << "/ptetatrkAllReferencematched. Skipping.";
42  return;
43  }
44 
45  MonitorElement* hphi_z_tkAllReference = igetter.get(topFolder_ + "/phiztrkAllReference");
46  if (hphi_z_tkAllReference == nullptr) {
47  edm::LogError("SiPixelTrackComparisonHarvester")
48  << "MonitorElement not found: " << topFolder_ << "/phiztrkAllReference. Skipping.";
49  return;
50  }
51 
52  MonitorElement* hphi_z_tkAllReferencematched = igetter.get(topFolder_ + "/phiztrkAllReferencematched");
53  if (hphi_z_tkAllReferencematched == nullptr) {
54  edm::LogError("SiPixelTrackComparisonHarvester")
55  << "MonitorElement not found: " << topFolder_ << "/phiztrkAllReferencematched. Skipping.";
56  return;
57  }
58 
59  ibooker.cd();
61  MonitorElement* hpt_eta_matchRatio = ibooker.book2D(
62  "matchingeff_pt_eta", "Efficiency of track matching; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
63  MonitorElement* hphi_z_matchRatio = ibooker.book2D(
64  "matchingeff_phi_z", "Efficiency of track matching; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
65 
66  hpt_eta_matchRatio->divide(hpt_eta_tkAllReferencematched, hpt_eta_tkAllReference, 1., 1., "B");
67  hphi_z_matchRatio->divide(hphi_z_tkAllReferencematched, hphi_z_tkAllReference, 1., 1., "B");
68 
69  // now create the 1D projection from the 2D histograms
70  std::vector<std::string> listOfMEsToProject = {"nTracks",
71  "nLooseAndAboveTracks",
72  "nLooseAndAboveTracks_matched",
73  "nRecHits",
74  "nLayers",
75  "nChi2ndof",
76  "charge",
77  "pt",
78  "curvature",
79  "eta",
80  "phi",
81  "z",
82  "tip"};
83  for (const auto& me : listOfMEsToProject) {
84  MonitorElement* input2D = igetter.get(topFolder_ + "/" + me);
85  edm::LogPrint("SiPixelTrackComparisonHarvester") << "processing " << topFolder_ + "/" + me;
86  this->project2DalongDiagonal(input2D, ibooker);
87  }
88 }
89 
91  if (input2D == nullptr) {
92  edm::LogError("SiPixelTrackComparisonHarvester")
93  << "ME needed for diagonal projection is not found in the input file at" << topFolder_ << ". Skipping.";
94  return;
95  }
96 
97  ibooker.cd();
98  ibooker.setCurrentFolder(topFolder_ + "/projectedDifferences");
99  const auto& h_name = fmt::sprintf("%s_proj", input2D->getName());
100  const auto& h_title = fmt::sprintf(";%s CPU -GPU difference", input2D->getTitle());
101  const auto& span = (input2D->getAxisMax() - input2D->getAxisMin());
102  const auto& b_w = span / input2D->getNbinsX();
103  const auto& nbins = ((input2D->getNbinsX() % 2) == 0) ? input2D->getNbinsX() + 1 : input2D->getNbinsX();
104 
105  MonitorElement* diagonalized = ibooker.book1D(h_name, h_title, nbins, -span / 2., span / 2.);
106 
107  // collect all the entry on each diagonal of the 2D histogram
108  for (int i = 1; i <= input2D->getNbinsX(); i++) {
109  for (int j = 1; j <= input2D->getNbinsY(); j++) {
110  diagonalized->Fill((i - j) * b_w, input2D->getBinContent(i, j));
111  }
112  }
113 
114  // zero the error on the bin as it's sort of meaningless for the way we fill it
115  // by collecting the entry on each diagonal
116  for (int bin = 1; bin <= diagonalized->getNbinsX(); bin++) {
117  diagonalized->setBinError(bin, 0.f);
118  }
119 }
120 
125  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareGPUvsCPU/");
126  descriptions.addWithDefaultLabel(desc);
127 }
128 
129 //define this as a plug-in
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
SiPixelTrackComparisonHarvester(const edm::ParameterSet &)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
Log< level::Error, false > LogError
void Fill(long long x)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void dqmEndJob(DQMStore::IBooker &ibooker, DQMStore::IGetter &igetter) override
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void project2DalongDiagonal(MonitorElement *input2D, DQMStore::IBooker &ibooker)
Log< level::Warning, true > LogPrint
virtual int getNbinsY() const
get # of bins in Y-axis
#define M_PI
virtual std::string getTitle() const
get MonitorElement title
~SiPixelTrackComparisonHarvester() override=default
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:712
virtual double getAxisMin(int axis=1) const
const std::string & getName() const
get name of ME
virtual double getAxisMax(int axis=1) const
virtual int getNbinsX() const
get # of bins in X-axis
virtual void setBinError(int binx, double error)
set uncertainty on content of bin (1-D)
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
virtual void divide(const MonitorElement *, const MonitorElement *, double, double, const char *)
Replace entries with results of dividing num by denom.
virtual double getBinContent(int binx) const
get content of bin (1-D)