CMS 3D CMS Logo

SiPixelCompareTracks.cc
Go to the documentation of this file.
1 // TODO: change file name to SiPixelCompareTracksSoA.cc when CUDA code is removed
2 
3 // -*- C++ -*-
4 // Package: SiPixelCompareTracks
5 // Class: SiPixelCompareTracks
6 //
9 //
10 // Author: Suvankar Roy Chowdhury
11 //
12 
13 // for string manipulations
14 #include <fmt/printf.h>
24 // DQM Histograming
28 // DataFormats
31 
32 namespace {
33  // same logic used for the MTV:
34  // cf https://github.com/cms-sw/cmssw/blob/master/Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc
36 
37  void setBinLog(TAxis* axis) {
38  int bins = axis->GetNbins();
39  float from = axis->GetXmin();
40  float to = axis->GetXmax();
41  float width = (to - from) / bins;
42  std::vector<float> new_bins(bins + 1, 0);
43  for (int i = 0; i <= bins; i++) {
44  new_bins[i] = TMath::Power(10, from + i * width);
45  }
46  axis->Set(bins, new_bins.data());
47  }
48 
49  void setBinLogX(TH1* h) {
50  TAxis* axis = h->GetXaxis();
51  setBinLog(axis);
52  }
53  void setBinLogY(TH1* h) {
54  TAxis* axis = h->GetYaxis();
55  setBinLog(axis);
56  }
57 
58  template <typename... Args>
59  dqm::reco::MonitorElement* make2DIfLog(DQMStore::IBooker& ibook, bool logx, bool logy, Args&&... args) {
60  auto h = std::make_unique<TH2I>(std::forward<Args>(args)...);
61  if (logx)
62  setBinLogX(h.get());
63  if (logy)
64  setBinLogY(h.get());
65  const auto& name = h->GetName();
66  return ibook.book2I(name, h.release());
67  }
68 } // namespace
69 
70 // TODO: change class name to SiPixelCompareTracksSoA when CUDA code is removed
71 template <typename T>
73 public:
75 
76  explicit SiPixelCompareTracks(const edm::ParameterSet&);
77  ~SiPixelCompareTracks() override = default;
78  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
79  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
80  // analyzeSeparate is templated to accept distinct types of SoAs
81  // The default use case is to use tracks from Alpaka reconstructed on CPU and GPU;
82  template <typename U, typename V>
83  void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event& iEvent);
84  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
85 
86 private:
87  // these two are both on Host but originally they have been produced on Host or on Device
91  const bool useQualityCut_;
93  const float dr2cut_;
114  //1D differences
121 
122  //for matching eff vs region: derive the ratio at harvesting
127 };
128 
129 //
130 // constructors
131 //
132 
133 template <typename T>
135  : tokenSoATrackReference_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackReferenceSoA"))),
136  tokenSoATrackTarget_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackTargetSoA"))),
137  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
138  useQualityCut_(iConfig.getParameter<bool>("useQualityCut")),
139  minQuality_(pixelTrack::qualityByName(iConfig.getParameter<std::string>("minQuality"))),
140  dr2cut_(iConfig.getParameter<double>("deltaR2cut")) {}
141 
142 template <typename T>
143 template <typename U, typename V>
145  using helper = TracksUtilities<T>;
146 
147  const auto& tsoaHandleRef = iEvent.getHandle(tokenRef);
148  const auto& tsoaHandleTar = iEvent.getHandle(tokenTar);
149 
150  if (not tsoaHandleRef or not tsoaHandleTar) {
151  edm::LogWarning out("SiPixelCompareTracks");
152  if (not tsoaHandleRef) {
153  out << "reference tracks not found; ";
154  }
155  if (not tsoaHandleTar) {
156  out << "target tracks not found; ";
157  }
158  out << "the comparison will not run.";
159  return;
160  }
161 
162  auto const& tsoaRef = *tsoaHandleRef;
163  auto const& tsoaTar = *tsoaHandleTar;
164 
165  auto maxTracksRef = tsoaRef.view().metadata().size(); //this should be same for both?
166  auto maxTracksTar = tsoaTar.view().metadata().size(); //this should be same for both?
167 
168  auto const* qualityRef = tsoaRef.view().quality();
169  auto const* qualityTar = tsoaTar.view().quality();
170 
171  int32_t nTracksRef = 0;
172  int32_t nTracksTar = 0;
173  int32_t nLooseAndAboveTracksRef = 0;
174  int32_t nLooseAndAboveTracksRef_matchedTar = 0;
175  int32_t nLooseAndAboveTracksTar = 0;
176 
177  //Loop over Tar tracks and store the indices of the loose tracks. Whats happens if useQualityCut_ is false?
178  std::vector<int32_t> looseTrkidxTar;
179  for (int32_t jt = 0; jt < maxTracksTar; ++jt) {
180  if (helper::nHits(tsoaTar.view(), jt) == 0)
181  break; // this is a guard
182  if (!(tsoaTar.view()[jt].pt() > 0.))
183  continue;
184  nTracksTar++;
185  if (useQualityCut_ && qualityTar[jt] < minQuality_)
186  continue;
187  nLooseAndAboveTracksTar++;
188  looseTrkidxTar.emplace_back(jt);
189  }
190 
191  //Now loop over Ref tracks//nested loop for loose gPU tracks
192  for (int32_t it = 0; it < maxTracksRef; ++it) {
193  int nHitsRef = helper::nHits(tsoaRef.view(), it);
194 
195  if (nHitsRef == 0)
196  break; // this is a guard
197 
198  float ptRef = tsoaRef.view()[it].pt();
199  float etaRef = tsoaRef.view()[it].eta();
200  float phiRef = reco::phi(tsoaRef.view(), it);
201  float zipRef = reco::zip(tsoaRef.view(), it);
202  float tipRef = reco::tip(tsoaRef.view(), it);
203 
204  if (!(ptRef > 0.))
205  continue;
206  nTracksRef++;
207  if (useQualityCut_ && qualityRef[it] < minQuality_)
208  continue;
209  nLooseAndAboveTracksRef++;
210  //Now loop over loose Tar trk and find the closest in DeltaR//do we need pt cut?
211  const int32_t notFound = -1;
212  int32_t closestTkidx = notFound;
213  float mindr2 = dr2cut_;
214 
215  for (auto gid : looseTrkidxTar) {
216  float etaTar = tsoaTar.view()[gid].eta();
217  float phiTar = reco::phi(tsoaTar.view(), gid);
218  float dr2 = reco::deltaR2(etaRef, phiRef, etaTar, phiTar);
219  if (dr2 > dr2cut_)
220  continue; // this is arbitrary
221  if (mindr2 > dr2) {
222  mindr2 = dr2;
223  closestTkidx = gid;
224  }
225  }
226 
227  hpt_eta_tkAllRef_->Fill(etaRef, ptRef); //all Ref tk
228  hphi_z_tkAllRef_->Fill(phiRef, zipRef);
229  if (closestTkidx == notFound)
230  continue;
231  nLooseAndAboveTracksRef_matchedTar++;
232 
233  hchi2_->Fill(tsoaRef.view()[it].chi2(), tsoaTar.view()[closestTkidx].chi2());
234  hCharge_->Fill(reco::charge(tsoaRef.view(), it), reco::charge(tsoaTar.view(), closestTkidx));
235  hnHits_->Fill(helper::nHits(tsoaRef.view(), it), helper::nHits(tsoaTar.view(), closestTkidx));
236  hnLayers_->Fill(tsoaRef.view()[it].nLayers(), tsoaTar.view()[closestTkidx].nLayers());
237  hpt_->Fill(ptRef, tsoaTar.view()[closestTkidx].pt());
238  hptLogLog_->Fill(ptRef, tsoaTar.view()[closestTkidx].pt());
239  heta_->Fill(etaRef, tsoaTar.view()[closestTkidx].eta());
240  hphi_->Fill(phiRef, reco::phi(tsoaTar.view(), closestTkidx));
241  hz_->Fill(zipRef, reco::zip(tsoaTar.view(), closestTkidx));
242  htip_->Fill(tipRef, reco::tip(tsoaTar.view(), closestTkidx));
243  hptdiffMatched_->Fill(ptRef - tsoaTar.view()[closestTkidx].pt());
244  hCurvdiffMatched_->Fill((reco::charge(tsoaRef.view(), it) / tsoaRef.view()[it].pt()) -
245  (reco::charge(tsoaTar.view(), closestTkidx) / tsoaTar.view()[closestTkidx].pt()));
246  hetadiffMatched_->Fill(etaRef - tsoaTar.view()[closestTkidx].eta());
247  hphidiffMatched_->Fill(reco::deltaPhi(phiRef, reco::phi(tsoaTar.view(), closestTkidx)));
248  hzdiffMatched_->Fill(zipRef - reco::zip(tsoaTar.view(), closestTkidx));
249  htipdiffMatched_->Fill(tipRef - reco::tip(tsoaTar.view(), closestTkidx));
250  hpt_eta_tkAllRefMatched_->Fill(etaRef, tsoaRef.view()[it].pt()); //matched to gpu
251  hphi_z_tkAllRefMatched_->Fill(etaRef, zipRef);
252  }
253  hnTracks_->Fill(nTracksRef, nTracksTar);
254  hnLooseAndAboveTracks_->Fill(nLooseAndAboveTracksRef, nLooseAndAboveTracksTar);
255  hnLooseAndAboveTracks_matched_->Fill(nLooseAndAboveTracksRef, nLooseAndAboveTracksRef_matchedTar);
256 }
257 
258 //
259 // -- Analyze
260 //
261 template <typename T>
263  // The default use case is to use vertices from Alpaka reconstructed on CPU and GPU;
264  // The function is left templated if any other cases need to be added
265  analyzeSeparate(tokenSoATrackReference_, tokenSoATrackTarget_, iEvent);
266 }
267 
268 //
269 // -- Book Histograms
270 //
271 template <typename T>
273  edm::Run const& iRun,
274  edm::EventSetup const& iSetup) {
275  iBook.cd();
276  iBook.setCurrentFolder(topFolderName_);
277 
278  // clang-format off
279  std::string toRep = "Number of tracks";
280  // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports THnSparse
281  // these should be moved to a less resource consuming format
282  hnTracks_ = iBook.book2I("nTracks", fmt::sprintf("%s per event; Reference; Target",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
283  hnLooseAndAboveTracks_ = iBook.book2I("nLooseAndAboveTracks", fmt::sprintf("%s (quality #geq loose) per event; Reference; Target",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
284  hnLooseAndAboveTracks_matched_ = iBook.book2I("nLooseAndAboveTracks_matched", fmt::sprintf("%s (quality #geq loose) per event; Reference; Target",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
285 
286  toRep = "Number of all RecHits per track (quality #geq loose)";
287  hnHits_ = iBook.book2I("nRecHits", fmt::sprintf("%s;Reference;Target",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
288 
289  toRep = "Number of all layers per track (quality #geq loose)";
290  hnLayers_ = iBook.book2I("nLayers", fmt::sprintf("%s;Reference;Target",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
291 
292  toRep = "Track (quality #geq loose) #chi^{2}/ndof";
293  hchi2_ = iBook.book2I("nChi2ndof", fmt::sprintf("%s;Reference;Target",toRep), 40, 0., 20., 40, 0., 20.);
294 
295  toRep = "Track (quality #geq loose) charge";
296  hCharge_ = iBook.book2I("charge",fmt::sprintf("%s;Reference;Target",toRep),3, -1.5, 1.5, 3, -1.5, 1.5);
297 
298  hpt_ = iBook.book2I("pt", "Track (quality #geq loose) p_{T} [GeV];Reference;Target", 200, 0., 200., 200, 0., 200.);
299  hptLogLog_ = make2DIfLog(iBook, true, true, "ptLogLog", "Track (quality #geq loose) p_{T} [GeV];Reference;Target", 200, log10(0.5), log10(200.), 200, log10(0.5), log10(200.));
300  heta_ = iBook.book2I("eta", "Track (quality #geq loose) #eta;Reference;Target", 30, -3., 3., 30, -3., 3.);
301  hphi_ = iBook.book2I("phi", "Track (quality #geq loose) #phi;Reference;Target", 30, -M_PI, M_PI, 30, -M_PI, M_PI);
302  hz_ = iBook.book2I("z", "Track (quality #geq loose) z [cm];Reference;Target", 30, -30., 30., 30, -30., 30.);
303  htip_ = iBook.book2I("tip", "Track (quality #geq loose) TIP [cm];Reference;Target", 100, -0.5, 0.5, 100, -0.5, 0.5);
304  //1D difference plots
305  hptdiffMatched_ = iBook.book1D("ptdiffmatched", " p_{T} diff [GeV] between matched tracks; #Delta p_{T} [GeV]", 60, -30., 30.);
306  hCurvdiffMatched_ = iBook.book1D("curvdiffmatched", "q/p_{T} diff [GeV] between matched tracks; #Delta q/p_{T} [GeV]", 60, -30., 30.);
307  hetadiffMatched_ = iBook.book1D("etadiffmatched", " #eta diff between matched tracks; #Delta #eta", 160, -0.04 ,0.04);
308  hphidiffMatched_ = iBook.book1D("phidiffmatched", " #phi diff between matched tracks; #Delta #phi", 160, -0.04 ,0.04);
309  hzdiffMatched_ = iBook.book1D("zdiffmatched", " z diff between matched tracks; #Delta z [cm]", 300, -1.5, 1.5);
310  htipdiffMatched_ = iBook.book1D("tipdiffmatched", " TIP diff between matched tracks; #Delta TIP [cm]", 300, -1.5, 1.5);
311  //2D plots for eff
312  hpt_eta_tkAllRef_ = iBook.book2I("ptetatrkAllReference", "Track (quality #geq loose) on Reference; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
313  hpt_eta_tkAllRefMatched_ = iBook.book2I("ptetatrkAllReferencematched", "Track (quality #geq loose) on Reference matched to Target track; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
314 
315  hphi_z_tkAllRef_ = iBook.book2I("phiztrkAllReference", "Track (quality #geq loose) on Reference; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
316  hphi_z_tkAllRefMatched_ = iBook.book2I("phiztrkAllReferencematched", "Track (quality #geq loose) on Reference; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
317 
318 }
319 
320 template<typename T>
322  // monitorpixelTrackSoA
324  desc.add<edm::InputTag>("pixelTrackReferenceSoA", edm::InputTag("pixelTracksAlpakaSerial"));
325  desc.add<edm::InputTag>("pixelTrackTargetSoA", edm::InputTag("pixelTracksAlpaka"));
326  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareDeviceVSHost");
327  desc.add<bool>("useQualityCut", true);
328  desc.add<std::string>("minQuality", "loose");
329  desc.add<double>("deltaR2cut", 0.04);
330  descriptions.addWithDefaultLabel(desc);
331 }
332 
333 // TODO: change module names to SiPixel*CompareTracksSoA when CUDA code is removed
334 
338 
342 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
MonitorElement * hnLooseAndAboveTracks_
Definition: helper.py:1
MonitorElement * hnLayersVsPhi_
MonitorElement * hptdiffMatched_
MonitorElement * hphi_z_tkAllRef_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
Quality qualityByName(std::string const &name)
SiPixelCompareTracks(const edm::ParameterSet &)
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackReference_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
MonitorElement * hphi_z_tkAllRefMatched_
void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event &iEvent)
T eta() const
Definition: PV3DBase.h:73
MonitorElement * hnLayersVsEta_
MonitorElement * hCurvdiffMatched_
MonitorElement * hnTracks_
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float charge(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:73
MonitorElement * hpt_eta_tkAllRefMatched_
MonitorElement * hnLooseAndAboveTracks_matched_
MonitorElement * hphidiffMatched_
dqm::reco::DQMStore DQMStore
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float tip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:85
int iEvent
Definition: GenABIO.cc:224
MonitorElement * hpt_eta_tkAllRef_
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t V
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
MonitorElement * hChi2VsEta_
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float phi(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:80
MonitorElement * hnHitsVsEta_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
MonitorElement * hzdiffMatched_
#define M_PI
~SiPixelCompareTracks() override=default
MonitorElement * hChi2VsPhi_
MonitorElement * hetadiffMatched_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static const GlobalPoint notFound(0, 0, 0)
const pixelTrack::Quality minQuality_
MonitorElement * htipdiffMatched_
HLT enums.
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
const std::string topFolderName_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
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:305
MonitorElement * hnHitsVsPhi_
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackTarget_
Definition: Run.h:45