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_tkAllCPU = igetter.get(topFolder_ + "/ptetatrkAllCPU");
32  MonitorElement* hpt_eta_tkAllCPUmatched = igetter.get(topFolder_ + "/ptetatrkAllCPUmatched");
33  MonitorElement* hphi_z_tkAllCPU = igetter.get(topFolder_ + "/phiztrkAllCPU");
34  MonitorElement* hphi_z_tkAllCPUmatched = igetter.get(topFolder_ + "/phiztrkAllCPUmatched");
35 
36  if (hpt_eta_tkAllCPU == nullptr or hpt_eta_tkAllCPUmatched == nullptr or hphi_z_tkAllCPU == nullptr or
37  hphi_z_tkAllCPUmatched == nullptr) {
38  edm::LogError("SiPixelTrackComparisonHarvester")
39  << "MEs needed for this module are not found in the input file. Skipping.";
40  return;
41  }
42 
43  ibooker.cd();
45  MonitorElement* hpt_eta_matchRatio = ibooker.book2D(
46  "matchingeff_pt_eta", "Efficiency of track matching; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
47  MonitorElement* hphi_z_matchRatio = ibooker.book2D(
48  "matchingeff_phi_z", "Efficiency of track matching; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
49 
50  hpt_eta_matchRatio->divide(hpt_eta_tkAllCPUmatched, hpt_eta_tkAllCPU, 1., 1., "B");
51  hphi_z_matchRatio->divide(hphi_z_tkAllCPUmatched, hphi_z_tkAllCPU, 1., 1., "B");
52 
53  // now create the 1D projection from the 2D histograms
54  std::vector<std::string> listOfMEsToProject = {"nTracks",
55  "nLooseAndAboveTracks",
56  "nLooseAndAboveTracks_matched",
57  "nRecHits",
58  "nLayers",
59  "nChi2ndof",
60  "charge",
61  "pt",
62  "eta",
63  "phi",
64  "z",
65  "tip"};
66  for (const auto& me : listOfMEsToProject) {
67  MonitorElement* input2D = igetter.get(topFolder_ + "/" + me);
68  this->project2DalongDiagonal(input2D, ibooker);
69  }
70 }
71 
73  if (input2D == nullptr) {
74  edm::LogError("SiPixelTrackComparisonHarvester")
75  << "MEs needed for diagonal projection are not found in the input file. Skipping.";
76  return;
77  }
78 
79  ibooker.cd();
80  ibooker.setCurrentFolder(topFolder_ + "/projectedDifferences");
81  const auto& h_name = fmt::sprintf("%s_proj", input2D->getName());
82  const auto& h_title = fmt::sprintf(";%s CPU -GPU difference", input2D->getTitle());
83  const auto& span = (input2D->getAxisMax() - input2D->getAxisMin());
84  const auto& b_w = span / input2D->getNbinsX();
85  const auto& nbins = ((input2D->getNbinsX() % 2) == 0) ? input2D->getNbinsX() + 1 : input2D->getNbinsX();
86 
87  MonitorElement* diagonalized = ibooker.book1D(h_name, h_title, nbins, -span / 2., span / 2.);
88 
89  // collect all the entry on each diagonal of the 2D histogram
90  for (int i = 1; i <= input2D->getNbinsX(); i++) {
91  for (int j = 1; j <= input2D->getNbinsY(); j++) {
92  diagonalized->Fill((i - j) * b_w, input2D->getBinContent(i, j));
93  }
94  }
95 
96  // zero the error on the bin as it's sort of meaningless for the way we fill it
97  // by collecting the entry on each diagonal
98  for (int bin = 1; bin <= diagonalized->getNbinsX(); bin++) {
99  diagonalized->setBinError(bin, 0.f);
100  }
101 }
102 
107  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareGPUvsCPU/");
108  descriptions.addWithDefaultLabel(desc);
109 }
110 
111 //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)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
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)
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:212
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:697
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)