CMS 3D CMS Logo

SiPixelPhase1CompareTrackSoA.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 // Package: SiPixelPhase1CompareTrackSoA
3 // Class: SiPixelPhase1CompareTrackSoA
4 //
7 //
8 // Author: Suvankar Roy Chowdhury
9 //
19 // DQM Histograming
24 // for string manipulations
25 #include <fmt/printf.h>
26 
27 namespace {
28  // same logic used for the MTV:
29  // cf https://github.com/cms-sw/cmssw/blob/master/Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc
31 
32  void setBinLog(TAxis* axis) {
33  int bins = axis->GetNbins();
34  float from = axis->GetXmin();
35  float to = axis->GetXmax();
36  float width = (to - from) / bins;
37  std::vector<float> new_bins(bins + 1, 0);
38  for (int i = 0; i <= bins; i++) {
39  new_bins[i] = TMath::Power(10, from + i * width);
40  }
41  axis->Set(bins, new_bins.data());
42  }
43 
44  void setBinLogX(TH1* h) {
45  TAxis* axis = h->GetXaxis();
46  setBinLog(axis);
47  }
48  void setBinLogY(TH1* h) {
49  TAxis* axis = h->GetYaxis();
50  setBinLog(axis);
51  }
52 
53  template <typename... Args>
54  dqm::reco::MonitorElement* make2DIfLog(DQMStore::IBooker& ibook, bool logx, bool logy, Args&&... args) {
55  auto h = std::make_unique<TH2I>(std::forward<Args>(args)...);
56  if (logx)
57  setBinLogX(h.get());
58  if (logy)
59  setBinLogY(h.get());
60  const auto& name = h->GetName();
61  return ibook.book2I(name, h.release());
62  }
63 } // namespace
64 
66 public:
68  ~SiPixelPhase1CompareTrackSoA() override = default;
69  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
70  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
71  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
72 
73 private:
77  const bool useQualityCut_;
79  const float dr2cut_;
99  //1D differences
104  //for matching eff vs region:dervie ration at harvesting
109 };
110 
111 //
112 // constructors
113 //
114 
116  : tokenSoATrackCPU_(consumes<PixelTrackHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcCPU"))),
117  tokenSoATrackGPU_(consumes<PixelTrackHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcGPU"))),
118  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
119  useQualityCut_(iConfig.getParameter<bool>("useQualityCut")),
120  minQuality_(pixelTrack::qualityByName(iConfig.getParameter<std::string>("minQuality"))),
121  dr2cut_(iConfig.getParameter<double>("deltaR2cut")) {}
122 
123 //
124 // -- Analyze
125 //
127  const auto& tsoaHandleCPU = iEvent.getHandle(tokenSoATrackCPU_);
128  const auto& tsoaHandleGPU = iEvent.getHandle(tokenSoATrackGPU_);
129  if (not tsoaHandleCPU or not tsoaHandleGPU) {
130  edm::LogWarning out("SiPixelPhase1CompareTrackSoA");
131  if (not tsoaHandleCPU) {
132  out << "reference (cpu) tracks not found; ";
133  }
134  if (not tsoaHandleGPU) {
135  out << "target (gpu) tracks not found; ";
136  }
137  out << "the comparison will not run.";
138  return;
139  }
140 
141  auto const& tsoaCPU = *tsoaHandleCPU->get();
142  auto const& tsoaGPU = *tsoaHandleGPU->get();
143  auto maxTracksCPU = tsoaCPU.stride(); //this should be same for both?
144  auto maxTracksGPU = tsoaGPU.stride(); //this should be same for both?
145  auto const* qualityCPU = tsoaCPU.qualityData();
146  auto const* qualityGPU = tsoaGPU.qualityData();
147  int32_t nTracksCPU = 0;
148  int32_t nTracksGPU = 0;
149  int32_t nLooseAndAboveTracksCPU = 0;
150  int32_t nLooseAndAboveTracksCPU_matchedGPU = 0;
151  int32_t nLooseAndAboveTracksGPU = 0;
152 
153  //Loop over GPU tracks and store the indices of the loose tracks. Whats happens if useQualityCut_ is false?
154  std::vector<int32_t> looseTrkidxGPU;
155  for (int32_t jt = 0; jt < maxTracksGPU; ++jt) {
156  if (tsoaGPU.nHits(jt) == 0)
157  break; // this is a guard
158  if (!(tsoaGPU.pt(jt) > 0.))
159  continue;
160  nTracksGPU++;
161  if (useQualityCut_ && qualityGPU[jt] < minQuality_)
162  continue;
163  nLooseAndAboveTracksGPU++;
164  looseTrkidxGPU.emplace_back(jt);
165  }
166 
167  //Now loop over CPU tracks//nested loop for loose gPU tracks
168  for (int32_t it = 0; it < maxTracksCPU; ++it) {
169  if (tsoaCPU.nHits(it) == 0)
170  break; // this is a guard
171  if (!(tsoaCPU.pt(it) > 0.))
172  continue;
173  nTracksCPU++;
174  if (useQualityCut_ && qualityCPU[it] < minQuality_)
175  continue;
176  nLooseAndAboveTracksCPU++;
177  //Now loop over loose GPU trk and find the closest in DeltaR//do we need pt cut?
178  const int32_t notFound = -1;
179  int32_t closestTkidx = notFound;
180  float mindr2 = dr2cut_;
181  float etacpu = tsoaCPU.eta(it);
182  float phicpu = tsoaCPU.phi(it);
183  for (auto gid : looseTrkidxGPU) {
184  float etagpu = tsoaGPU.eta(gid);
185  float phigpu = tsoaGPU.phi(gid);
186  float dr2 = reco::deltaR2(etacpu, phicpu, etagpu, phigpu);
187  if (dr2 > dr2cut_)
188  continue; // this is arbitrary
189  if (mindr2 > dr2) {
190  mindr2 = dr2;
191  closestTkidx = gid;
192  }
193  }
194 
195  hpt_eta_tkAllCPU_->Fill(etacpu, tsoaCPU.pt(it)); //all CPU tk
196  hphi_z_tkAllCPU_->Fill(phicpu, tsoaCPU.zip(it));
197  if (closestTkidx == notFound)
198  continue;
199  nLooseAndAboveTracksCPU_matchedGPU++;
200 
201  hchi2_->Fill(tsoaCPU.chi2(it), tsoaGPU.chi2(closestTkidx));
202  hnHits_->Fill(tsoaCPU.nHits(it), tsoaGPU.nHits(closestTkidx));
203  hnLayers_->Fill(tsoaCPU.nLayers(it), tsoaGPU.nLayers(closestTkidx));
204  hpt_->Fill(tsoaCPU.pt(it), tsoaGPU.pt(closestTkidx));
205  hptLogLog_->Fill(tsoaCPU.pt(it), tsoaGPU.pt(closestTkidx));
206  heta_->Fill(etacpu, tsoaGPU.eta(closestTkidx));
207  hphi_->Fill(phicpu, tsoaGPU.phi(closestTkidx));
208  hz_->Fill(tsoaCPU.zip(it), tsoaGPU.zip(closestTkidx));
209  htip_->Fill(tsoaCPU.tip(it), tsoaGPU.tip(closestTkidx));
210  hptdiffMatched_->Fill(tsoaCPU.pt(it) - tsoaGPU.pt(closestTkidx));
211  hetadiffMatched_->Fill(etacpu - tsoaGPU.eta(closestTkidx));
212  hphidiffMatched_->Fill(reco::deltaPhi(phicpu, tsoaGPU.phi(closestTkidx)));
213  hzdiffMatched_->Fill(tsoaCPU.zip(it) - tsoaGPU.zip(closestTkidx));
214  hpt_eta_tkAllCPUMatched_->Fill(etacpu, tsoaCPU.pt(it)); //matched to gpu
215  hphi_z_tkAllCPUMatched_->Fill(phicpu, tsoaCPU.zip(it));
216  }
217  hnTracks_->Fill(nTracksCPU, nTracksGPU);
218  hnLooseAndAboveTracks_->Fill(nLooseAndAboveTracksCPU, nLooseAndAboveTracksGPU);
219  hnLooseAndAboveTracks_matched_->Fill(nLooseAndAboveTracksCPU, nLooseAndAboveTracksCPU_matchedGPU);
220 }
221 
222 //
223 // -- Book Histograms
224 //
226  edm::Run const& iRun,
227  edm::EventSetup const& iSetup) {
228  iBook.cd();
230 
231  // clang-format off
232  std::string toRep = "Number of tracks";
233  // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports either TH2I or THnSparse
234  // these should be moved to a less resource consuming format
235  hnTracks_ = iBook.book2I("nTracks", fmt::sprintf("%s per event; CPU; GPU",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
236  hnLooseAndAboveTracks_ = iBook.book2I("nLooseAndAboveTracks", fmt::sprintf("%s (quality #geq loose) per event; CPU; GPU",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
237  hnLooseAndAboveTracks_matched_ = iBook.book2I("nLooseAndAboveTracks_matched", fmt::sprintf("%s (quality #geq loose) per event; CPU; GPU",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
238 
239  toRep = "Number of all RecHits per track (quality #geq loose)";
240  hnHits_ = iBook.book2I("nRecHits", fmt::sprintf("%s;CPU;GPU",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
241 
242  toRep = "Number of all layers per track (quality #geq loose)";
243  hnLayers_ = iBook.book2I("nLayers", fmt::sprintf("%s;CPU;GPU",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
244 
245  toRep = "Track (quality #geq loose) #chi^{2}/ndof";
246  hchi2_ = iBook.book2I("nChi2ndof", fmt::sprintf("%s;CPU;GPU",toRep), 40, 0., 20., 40, 0., 20.);
247 
248  hpt_ = iBook.book2I("pt", "Track (quality #geq loose) p_{T} [GeV];CPU;GPU", 200, 0., 200., 200, 0., 200.);
249  hptLogLog_ = make2DIfLog(iBook, true, true, "ptLogLog", "Track (quality #geq loose) p_{T} [GeV];CPU;GPU", 200, log10(0.5), log10(200.), 200, log10(0.5), log10(200.));
250  heta_ = iBook.book2I("eta", "Track (quality #geq loose) #eta;CPU;GPU", 30, -3., 3., 30, -3., 3.);
251  hphi_ = iBook.book2I("phi", "Track (quality #geq loose) #phi;CPU;GPU", 30, -M_PI, M_PI, 30, -M_PI, M_PI);
252  hz_ = iBook.book2I("z", "Track (quality #geq loose) z [cm];CPU;GPU", 30, -30., 30., 30, -30., 30.);
253  htip_ = iBook.book2I("tip", "Track (quality #geq loose) TIP [cm];CPU;GPU", 100, -0.5, 0.5, 100, -0.5, 0.5);
254  //1D difference plots
255  hptdiffMatched_ = iBook.book1D("ptdiffmatched", " p_{T} diff [GeV] between matched tracks; #Delta p_{T} [GeV]", 60, -30., 30.);
256  hetadiffMatched_ = iBook.book1D("etadiffmatched", " #eta diff between matched tracks; #Delta #eta", 160, -0.04 ,0.04);
257  hphidiffMatched_ = iBook.book1D("phidiffmatched", " #phi diff between matched tracks; #Delta #phi", 160, -0.04 ,0.04);
258  hzdiffMatched_ = iBook.book1D("zdiffmatched", " z diff between matched tracks; #Delta z [cm]", 300, -1.5, 1.5);
259  //2D plots for eff
260  hpt_eta_tkAllCPU_ = iBook.book2I("ptetatrkAllCPU", "Track (quality #geq loose) on CPU; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
261  hpt_eta_tkAllCPUMatched_ = iBook.book2I("ptetatrkAllCPUmatched", "Track (quality #geq loose) on CPU matched to GPU track; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
262 
263  hphi_z_tkAllCPU_ = iBook.book2I("phiztrkAllCPU", "Track (quality #geq loose) on CPU; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
264  hphi_z_tkAllCPUMatched_ = iBook.book2I("phiztrkAllCPUmatched", "Track (quality #geq loose) on CPU; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
265 
266 }
267 
269  // monitorpixelTrackSoA
271  desc.add<edm::InputTag>("pixelTrackSrcCPU", edm::InputTag("pixelTracksSoA@cpu"));
272  desc.add<edm::InputTag>("pixelTrackSrcGPU", edm::InputTag("pixelTracksSoA@cuda"));
273  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareGPUvsCPU");
274  desc.add<bool>("useQualityCut", true);
275  desc.add<std::string>("minQuality", "loose");
276  desc.add<double>("deltaR2cut", 0.04);
277  descriptions.addWithDefaultLabel(desc);
278 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
Quality qualityByName(std::string const &name)
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrackGPU_
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrackCPU_
dqm::reco::DQMStore DQMStore
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
#define M_PI
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static const GlobalPoint notFound(0, 0, 0)
HLT enums.
~SiPixelPhase1CompareTrackSoA() override=default
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
MonitorElement * book2I(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:296
SiPixelPhase1CompareTrackSoA(const edm::ParameterSet &)
Definition: Run.h:45