CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Attributes
SiPixelCompareTrackSoA< T > Class Template Reference

#include <SiPixelCompareTrackSoA.cc>

Inheritance diagram for SiPixelCompareTrackSoA< T >:

Public Types

using PixelTrackSoA = TrackSoAHeterogeneousHost< T >
 

Public Member Functions

void analyze (const edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
void bookHistograms (DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
 
 SiPixelCompareTrackSoA (const edm::ParameterSet &)
 
 ~SiPixelCompareTrackSoA () override=default
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

const float dr2cut_
 
MonitorElementhCharge_
 
MonitorElementhchi2_
 
MonitorElementhChi2VsEta_
 
MonitorElementhChi2VsPhi_
 
MonitorElementhCurvature_
 
MonitorElementhCurvdiffMatched_
 
MonitorElementhDeltaNLooseAndAboveTracks_
 
MonitorElementhDeltaNLooseAndAboveTracks_matched_
 
MonitorElementhDeltaNTracks_
 
MonitorElementheta_
 
MonitorElementhetadiffMatched_
 
MonitorElementhnHits_
 
MonitorElementhnHitsVsEta_
 
MonitorElementhnHitsVsPhi_
 
MonitorElementhnLayers_
 
MonitorElementhnLayersVsEta_
 
MonitorElementhnLayersVsPhi_
 
MonitorElementhnLooseAndAboveTracks_
 
MonitorElementhnLooseAndAboveTracks_matched_
 
MonitorElementhnTracks_
 
MonitorElementhphi_
 
MonitorElementhphi_z_tkAllRef_
 
MonitorElementhphi_z_tkAllRefMatched_
 
MonitorElementhphidiffMatched_
 
MonitorElementhpt_
 
MonitorElementhpt_eta_tkAllRef_
 
MonitorElementhpt_eta_tkAllRefMatched_
 
MonitorElementhptdiffMatched_
 
MonitorElementhptLogLog_
 
MonitorElementhquality_
 
MonitorElementhtip_
 
MonitorElementhtipdiffMatched_
 
MonitorElementhz_
 
MonitorElementhzdiffMatched_
 
const pixelTrack::Quality minQuality_
 
const edm::EDGetTokenT< PixelTrackSoAtokenSoATrackCPU_
 
const edm::EDGetTokenT< PixelTrackSoAtokenSoATrackGPU_
 
const std::string topFolderName_
 
const bool useQualityCut_
 

Detailed Description

template<typename T>
class SiPixelCompareTrackSoA< T >

Definition at line 67 of file SiPixelCompareTrackSoA.cc.

Member Typedef Documentation

◆ PixelTrackSoA

Definition at line 69 of file SiPixelCompareTrackSoA.cc.

Constructor & Destructor Documentation

◆ SiPixelCompareTrackSoA()

template<typename T >
SiPixelCompareTrackSoA< T >::SiPixelCompareTrackSoA ( const edm::ParameterSet iConfig)
explicit

Definition at line 128 of file SiPixelCompareTrackSoA.cc.

129  : tokenSoATrackCPU_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcCPU"))),
130  tokenSoATrackGPU_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcGPU"))),
131  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
132  useQualityCut_(iConfig.getParameter<bool>("useQualityCut")),
134  dr2cut_(iConfig.getParameter<double>("deltaR2cut")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Quality qualityByName(std::string const &name)
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackGPU_
const std::string topFolderName_
const pixelTrack::Quality minQuality_
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackCPU_

◆ ~SiPixelCompareTrackSoA()

template<typename T >
SiPixelCompareTrackSoA< T >::~SiPixelCompareTrackSoA ( )
overridedefault

Member Function Documentation

◆ analyze()

template<typename T >
void SiPixelCompareTrackSoA< T >::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 140 of file SiPixelCompareTrackSoA.cc.

References ALCARECOTkAlJpsiMuMu_cff::charge, reco::deltaPhi(), reco::deltaR2(), PV3DBase< T, PVType, FrameType >::eta(), iEvent, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, SiStripPI::min, nHits, notFound, or, MillePedeFileConverter_cfg::out, qcdUeDQM_cfi::tip, and reco::zip().

140  {
141  using helper = TracksUtilities<T>;
142  const auto& tsoaHandleCPU = iEvent.getHandle(tokenSoATrackCPU_);
143  const auto& tsoaHandleGPU = iEvent.getHandle(tokenSoATrackGPU_);
144  if (not tsoaHandleCPU or not tsoaHandleGPU) {
145  edm::LogWarning out("SiPixelCompareTrackSoA");
146  if (not tsoaHandleCPU) {
147  out << "reference (cpu) tracks not found; ";
148  }
149  if (not tsoaHandleGPU) {
150  out << "target (gpu) tracks not found; ";
151  }
152  out << "the comparison will not run.";
153  return;
154  }
155 
156  auto const& tsoaCPU = *tsoaHandleCPU;
157  auto const& tsoaGPU = *tsoaHandleGPU;
158  auto maxTracksCPU = tsoaCPU.view().metadata().size(); //this should be same for both?
159  auto maxTracksGPU = tsoaGPU.view().metadata().size(); //this should be same for both?
160  auto const* qualityCPU = tsoaCPU.view().quality();
161  auto const* qualityGPU = tsoaGPU.view().quality();
162  int32_t nTracksCPU = 0;
163  int32_t nTracksGPU = 0;
164  int32_t nLooseAndAboveTracksCPU = 0;
165  int32_t nLooseAndAboveTracksCPU_matchedGPU = 0;
166  int32_t nLooseAndAboveTracksGPU = 0;
167 
168  //Loop over GPU tracks and store the indices of the loose tracks. Whats happens if useQualityCut_ is false?
169  std::vector<int32_t> looseTrkidxGPU;
170  for (int32_t jt = 0; jt < maxTracksGPU; ++jt) {
171  if (helper::nHits(tsoaGPU.view(), jt) == 0)
172  break; // this is a guard
173  if (!(tsoaGPU.view()[jt].pt() > 0.))
174  continue;
175  nTracksGPU++;
176  if (useQualityCut_ && qualityGPU[jt] < minQuality_)
177  continue;
178  nLooseAndAboveTracksGPU++;
179  looseTrkidxGPU.emplace_back(jt);
180  }
181 
182  //Now loop over CPU tracks//nested loop for loose gPU tracks
183  for (int32_t it = 0; it < maxTracksCPU; ++it) {
184  int nHitsCPU = helper::nHits(tsoaCPU.view(), it);
185 
186  if (nHitsCPU == 0)
187  break; // this is a guard
188 
189  float ptCPU = tsoaCPU.view()[it].pt();
190  float etaCPU = tsoaCPU.view()[it].eta();
191  float phiCPU = helper::phi(tsoaCPU.view(), it);
192  float zipCPU = helper::zip(tsoaCPU.view(), it);
193  float tipCPU = helper::tip(tsoaCPU.view(), it);
194  auto qCPU = helper::charge(tsoaCPU.view(), it);
195 
196  if (!(ptCPU > 0.))
197  continue;
198  nTracksCPU++;
199  if (useQualityCut_ && qualityCPU[it] < minQuality_)
200  continue;
201  nLooseAndAboveTracksCPU++;
202  //Now loop over loose GPU trk and find the closest in DeltaR//do we need pt cut?
203  const int32_t notFound = -1;
204  int32_t closestTkidx = notFound;
205  float mindr2 = dr2cut_;
206 
207  for (auto gid : looseTrkidxGPU) {
208  float etaGPU = tsoaGPU.view()[gid].eta();
209  float phiGPU = helper::phi(tsoaGPU.view(), gid);
210  float dr2 = reco::deltaR2(etaCPU, phiCPU, etaGPU, phiGPU);
211  if (dr2 > dr2cut_)
212  continue; // this is arbitrary
213  if (mindr2 > dr2) {
214  mindr2 = dr2;
215  closestTkidx = gid;
216  }
217  }
218 
219  hpt_eta_tkAllRef_->Fill(etaCPU, ptCPU); //all CPU tk
220  hphi_z_tkAllRef_->Fill(phiCPU, zipCPU);
221  if (closestTkidx == notFound)
222  continue;
223  nLooseAndAboveTracksCPU_matchedGPU++;
224 
225  hchi2_->Fill(tsoaCPU.view()[it].chi2(), tsoaGPU.view()[closestTkidx].chi2());
226  hCharge_->Fill(qCPU, helper::charge(tsoaGPU.view(), closestTkidx));
227  hnHits_->Fill(helper::nHits(tsoaCPU.view(), it), helper::nHits(tsoaGPU.view(), closestTkidx));
228  hnLayers_->Fill(tsoaCPU.view()[it].nLayers(), tsoaGPU.view()[closestTkidx].nLayers());
229  hpt_->Fill(ptCPU, tsoaGPU.view()[closestTkidx].pt());
230  hCurvature_->Fill(qCPU / ptCPU, helper::charge(tsoaGPU.view(), closestTkidx) / tsoaGPU.view()[closestTkidx].pt());
231  hptLogLog_->Fill(ptCPU, tsoaGPU.view()[closestTkidx].pt());
232  heta_->Fill(etaCPU, tsoaGPU.view()[closestTkidx].eta());
233  hphi_->Fill(phiCPU, helper::phi(tsoaGPU.view(), closestTkidx));
234  hz_->Fill(zipCPU, helper::zip(tsoaGPU.view(), closestTkidx));
235  htip_->Fill(tipCPU, helper::tip(tsoaGPU.view(), closestTkidx));
236  hptdiffMatched_->Fill(ptCPU - tsoaGPU.view()[closestTkidx].pt());
237  hCurvdiffMatched_->Fill((helper::charge(tsoaCPU.view(), it) / tsoaCPU.view()[it].pt()) -
238  (helper::charge(tsoaGPU.view(), closestTkidx) / tsoaGPU.view()[closestTkidx].pt()));
239  hetadiffMatched_->Fill(etaCPU - tsoaGPU.view()[closestTkidx].eta());
240  hphidiffMatched_->Fill(reco::deltaPhi(phiCPU, helper::phi(tsoaGPU.view(), closestTkidx)));
241  hzdiffMatched_->Fill(zipCPU - helper::zip(tsoaGPU.view(), closestTkidx));
242  htipdiffMatched_->Fill(tipCPU - helper::tip(tsoaGPU.view(), closestTkidx));
243  hpt_eta_tkAllRefMatched_->Fill(etaCPU, tsoaCPU.view()[it].pt()); //matched to gpu
244  hphi_z_tkAllRefMatched_->Fill(etaCPU, zipCPU);
245  }
246 
247  // Define a lambda function for filling the histograms
248  auto fillHistogram = [](auto& histogram, auto xValue, auto yValue) { histogram->Fill(xValue, yValue); };
249 
250  // Define a lambda for filling delta histograms
251  auto fillDeltaHistogram = [](auto& histogram, int cpuValue, int gpuValue) {
252  histogram->Fill(std::min(cpuValue, 1000), std::clamp(gpuValue - cpuValue, -100, 100));
253  };
254 
255  // Fill the histograms
256  fillHistogram(hnTracks_, nTracksCPU, nTracksGPU);
257  fillHistogram(hnLooseAndAboveTracks_, nLooseAndAboveTracksCPU, nLooseAndAboveTracksGPU);
258  fillHistogram(hnLooseAndAboveTracks_matched_, nLooseAndAboveTracksCPU, nLooseAndAboveTracksCPU_matchedGPU);
259 
260  fillDeltaHistogram(hDeltaNTracks_, nTracksCPU, nTracksGPU);
261  fillDeltaHistogram(hDeltaNLooseAndAboveTracks_, nLooseAndAboveTracksCPU, nLooseAndAboveTracksGPU);
262  fillDeltaHistogram(hDeltaNLooseAndAboveTracks_matched_, nLooseAndAboveTracksCPU, nLooseAndAboveTracksCPU_matchedGPU);
263 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Definition: helper.py:1
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackGPU_
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
MonitorElement * hDeltaNLooseAndAboveTracks_matched_
MonitorElement * hnLooseAndAboveTracks_matched_
void Fill(long long x)
MonitorElement * hpt_eta_tkAllRefMatched_
MonitorElement * hphi_z_tkAllRefMatched_
int iEvent
Definition: GenABIO.cc:224
MonitorElement * hDeltaNLooseAndAboveTracks_
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
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 * hnLooseAndAboveTracks_
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackCPU_

◆ bookHistograms()

template<typename T >
void SiPixelCompareTrackSoA< T >::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
override

Definition at line 269 of file SiPixelCompareTrackSoA.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2I(), dqm::implementation::NavigatorBase::cd(), ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), M_PI, Skims_PA_cff::name, dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, runGCPTkAlMap::title, multiplicitycorr_cfi::xBins, multiplicitycorr_cfi::xMax, photonAnalyzer_cfi::xMin, multiplicitycorr_cfi::yBins, multiplicitycorr_cfi::yMax, and photonAnalyzer_cfi::yMin.

271  {
272  iBook.cd();
273  iBook.setCurrentFolder(topFolderName_);
274 
275  // Define a helper function for booking histograms
276  std::string toRep = "Number of tracks";
277  auto bookTracksTH2I = [&](const std::string& name,
278  const std::string& title,
279  int xBins,
280  double xMin,
281  double xMax,
282  int yBins,
283  double yMin,
284  double yMax) {
285  return iBook.book2I(name, fmt::sprintf(title, toRep), xBins, xMin, xMax, yBins, yMin, yMax);
286  };
287 
288  // Define common parameters for different histogram types
289  constexpr int xBins = 501;
290  constexpr double xMin = -0.5;
291  constexpr double xMax = 1001.5;
292 
293  constexpr int dXBins = 1001;
294  constexpr double dXMin = -0.5;
295  constexpr double dXMax = 1000.5;
296 
297  constexpr int dYBins = 201;
298  constexpr double dYMin = -100.5;
299  constexpr double dYMax = 100.5;
300 
301  // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports THnSparse
302  // these should be moved to a less resource consuming format
303 
304  // Book histograms using the helper function
305  // clang-format off
306  hnTracks_ = bookTracksTH2I("nTracks", "%s per event; Reference; Target", xBins, xMin, xMax, xBins, xMin, xMax);
307  hnLooseAndAboveTracks_ = bookTracksTH2I("nLooseAndAboveTracks", "%s (quality #geq loose) per event; Reference; Target", xBins, xMin, xMax, xBins, xMin, xMax);
308  hnLooseAndAboveTracks_matched_ = bookTracksTH2I("nLooseAndAboveTracks_matched", "%s (quality #geq loose) per event; Reference; Target", xBins, xMin, xMax, xBins, xMin, xMax);
309 
310  hDeltaNTracks_ = bookTracksTH2I("deltaNTracks", "%s per event; Reference; Target - Reference", dXBins, dXMin, dXMax, dYBins, dYMin, dYMax);
311  hDeltaNLooseAndAboveTracks_ = bookTracksTH2I("deltaNLooseAndAboveTracks", "%s (quality #geq loose) per event; Reference; Target - Reference", dXBins, dXMin, dXMax, dYBins, dYMin, dYMax);
312  hDeltaNLooseAndAboveTracks_matched_ = bookTracksTH2I("deltaNLooseAndAboveTracks_matched", "%s (quality #geq loose) per event; Reference; Target - Reference", dXBins, dXMin, dXMax, dYBins, dYMin, dYMax);
313 
314  toRep = "Number of all RecHits per track (quality #geq loose)";
315  hnHits_ = iBook.book2I("nRecHits", fmt::sprintf("%s;CPU;GPU",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
316 
317  toRep = "Number of all layers per track (quality #geq loose)";
318  hnLayers_ = iBook.book2I("nLayers", fmt::sprintf("%s;CPU;GPU",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
319 
320  toRep = "Track (quality #geq loose) #chi^{2}/ndof";
321  hchi2_ = iBook.book2I("nChi2ndof", fmt::sprintf("%s;CPU;GPU",toRep), 40, 0., 20., 40, 0., 20.);
322 
323  toRep = "Track (quality #geq loose) charge";
324  hCharge_ = iBook.book2I("charge",fmt::sprintf("%s;CPU;GPU",toRep),3, -1.5, 1.5, 3, -1.5, 1.5);
325 
326  hpt_ = iBook.book2I("pt", "Track (quality #geq loose) p_{T} [GeV];CPU;GPU", 200, 0., 200., 200, 0., 200.);
327  hCurvature_ = iBook.book2I("curvature", "Track (quality #geq loose) q/p_{T} [GeV^{-1}];CPU;GPU", 60,- 3., 3., 60, -3., 3. );
328  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.));
329  heta_ = iBook.book2I("eta", "Track (quality #geq loose) #eta;CPU;GPU", 30, -3., 3., 30, -3., 3.);
330  hphi_ = iBook.book2I("phi", "Track (quality #geq loose) #phi;CPU;GPU", 30, -M_PI, M_PI, 30, -M_PI, M_PI);
331  hz_ = iBook.book2I("z", "Track (quality #geq loose) z [cm];CPU;GPU", 30, -30., 30., 30, -30., 30.);
332  htip_ = iBook.book2I("tip", "Track (quality #geq loose) TIP [cm];CPU;GPU", 100, -0.5, 0.5, 100, -0.5, 0.5);
333  //1D difference plots
334  hptdiffMatched_ = iBook.book1D("ptdiffmatched", " p_{T} diff [GeV] between matched tracks; #Delta p_{T} [GeV]", 61, -30.5, 30.5);
335  hCurvdiffMatched_ = iBook.book1D("curvdiffmatched", "q/p_{T} diff [GeV^{-1}] between matched tracks; #Delta q/p_{T} [GeV^{-1}]", 61, -3.05, 3.05);
336  hetadiffMatched_ = iBook.book1D("etadiffmatched", " #eta diff between matched tracks; #Delta #eta", 161, -0.045 ,0.045);
337  hphidiffMatched_ = iBook.book1D("phidiffmatched", " #phi diff between matched tracks; #Delta #phi", 161, -0.045 ,0.045);
338  hzdiffMatched_ = iBook.book1D("zdiffmatched", " z diff between matched tracks; #Delta z [cm]", 301, -1.55, 1.55);
339  htipdiffMatched_ = iBook.book1D("tipdiffmatched", " TIP diff between matched tracks; #Delta TIP [cm]", 301, -1.55, 1.55);
340  //2D plots for eff
341  hpt_eta_tkAllRef_ = iBook.book2I("ptetatrkAllReference", "Track (quality #geq loose) on CPU; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
342  hpt_eta_tkAllRefMatched_ = iBook.book2I("ptetatrkAllReferencematched", "Track (quality #geq loose) on CPU matched to GPU track; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
343 
344  hphi_z_tkAllRef_ = iBook.book2I("phiztrkAllReference", "Track (quality #geq loose) on CPU; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
345  hphi_z_tkAllRefMatched_ = iBook.book2I("phiztrkAllReferencematched", "Track (quality #geq loose) on CPU; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
346 
347 }
const std::string topFolderName_
MonitorElement * hDeltaNLooseAndAboveTracks_matched_
MonitorElement * hnLooseAndAboveTracks_matched_
MonitorElement * hpt_eta_tkAllRefMatched_
MonitorElement * hphi_z_tkAllRefMatched_
MonitorElement * hDeltaNLooseAndAboveTracks_
#define M_PI
MonitorElement * hnLooseAndAboveTracks_

◆ fillDescriptions()

template<typename T >
void SiPixelCompareTrackSoA< T >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 350 of file SiPixelCompareTrackSoA.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

350  {
351  // monitorpixelTrackSoA
353  desc.add<edm::InputTag>("pixelTrackSrcCPU", edm::InputTag("pixelTracksSoA@cpu"));
354  desc.add<edm::InputTag>("pixelTrackSrcGPU", edm::InputTag("pixelTracksSoA@cuda"));
355  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareGPUvsCPU");
356  desc.add<bool>("useQualityCut", true);
357  desc.add<std::string>("minQuality", "loose");
358  desc.add<double>("deltaR2cut", 0.02 * 0.02)->setComment("deltaR2 cut between track on CPU and GPU");
359  descriptions.addWithDefaultLabel(desc);
360 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ dr2cut_

template<typename T >
const float SiPixelCompareTrackSoA< T >::dr2cut_
private

Definition at line 83 of file SiPixelCompareTrackSoA.cc.

◆ hCharge_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hCharge_
private

Definition at line 96 of file SiPixelCompareTrackSoA.cc.

◆ hchi2_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hchi2_
private

Definition at line 97 of file SiPixelCompareTrackSoA.cc.

◆ hChi2VsEta_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hChi2VsEta_
private

Definition at line 99 of file SiPixelCompareTrackSoA.cc.

◆ hChi2VsPhi_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hChi2VsPhi_
private

Definition at line 98 of file SiPixelCompareTrackSoA.cc.

◆ hCurvature_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hCurvature_
private

Definition at line 101 of file SiPixelCompareTrackSoA.cc.

◆ hCurvdiffMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hCurvdiffMatched_
private

Definition at line 110 of file SiPixelCompareTrackSoA.cc.

◆ hDeltaNLooseAndAboveTracks_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hDeltaNLooseAndAboveTracks_
private

Definition at line 88 of file SiPixelCompareTrackSoA.cc.

◆ hDeltaNLooseAndAboveTracks_matched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hDeltaNLooseAndAboveTracks_matched_
private

Definition at line 89 of file SiPixelCompareTrackSoA.cc.

◆ hDeltaNTracks_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hDeltaNTracks_
private

Definition at line 87 of file SiPixelCompareTrackSoA.cc.

◆ heta_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::heta_
private

Definition at line 103 of file SiPixelCompareTrackSoA.cc.

◆ hetadiffMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hetadiffMatched_
private

Definition at line 111 of file SiPixelCompareTrackSoA.cc.

◆ hnHits_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnHits_
private

Definition at line 90 of file SiPixelCompareTrackSoA.cc.

◆ hnHitsVsEta_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnHitsVsEta_
private

Definition at line 92 of file SiPixelCompareTrackSoA.cc.

◆ hnHitsVsPhi_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnHitsVsPhi_
private

Definition at line 91 of file SiPixelCompareTrackSoA.cc.

◆ hnLayers_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnLayers_
private

Definition at line 93 of file SiPixelCompareTrackSoA.cc.

◆ hnLayersVsEta_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnLayersVsEta_
private

Definition at line 95 of file SiPixelCompareTrackSoA.cc.

◆ hnLayersVsPhi_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnLayersVsPhi_
private

Definition at line 94 of file SiPixelCompareTrackSoA.cc.

◆ hnLooseAndAboveTracks_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnLooseAndAboveTracks_
private

Definition at line 85 of file SiPixelCompareTrackSoA.cc.

◆ hnLooseAndAboveTracks_matched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnLooseAndAboveTracks_matched_
private

Definition at line 86 of file SiPixelCompareTrackSoA.cc.

◆ hnTracks_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hnTracks_
private

Definition at line 84 of file SiPixelCompareTrackSoA.cc.

◆ hphi_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hphi_
private

Definition at line 104 of file SiPixelCompareTrackSoA.cc.

◆ hphi_z_tkAllRef_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hphi_z_tkAllRef_
private

Definition at line 119 of file SiPixelCompareTrackSoA.cc.

◆ hphi_z_tkAllRefMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hphi_z_tkAllRefMatched_
private

Definition at line 120 of file SiPixelCompareTrackSoA.cc.

◆ hphidiffMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hphidiffMatched_
private

Definition at line 112 of file SiPixelCompareTrackSoA.cc.

◆ hpt_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hpt_
private

Definition at line 100 of file SiPixelCompareTrackSoA.cc.

◆ hpt_eta_tkAllRef_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hpt_eta_tkAllRef_
private

Definition at line 117 of file SiPixelCompareTrackSoA.cc.

◆ hpt_eta_tkAllRefMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hpt_eta_tkAllRefMatched_
private

Definition at line 118 of file SiPixelCompareTrackSoA.cc.

◆ hptdiffMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hptdiffMatched_
private

Definition at line 109 of file SiPixelCompareTrackSoA.cc.

◆ hptLogLog_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hptLogLog_
private

Definition at line 102 of file SiPixelCompareTrackSoA.cc.

◆ hquality_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hquality_
private

Definition at line 107 of file SiPixelCompareTrackSoA.cc.

◆ htip_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::htip_
private

Definition at line 106 of file SiPixelCompareTrackSoA.cc.

◆ htipdiffMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::htipdiffMatched_
private

Definition at line 114 of file SiPixelCompareTrackSoA.cc.

◆ hz_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hz_
private

Definition at line 105 of file SiPixelCompareTrackSoA.cc.

◆ hzdiffMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoA< T >::hzdiffMatched_
private

Definition at line 113 of file SiPixelCompareTrackSoA.cc.

◆ minQuality_

template<typename T >
const pixelTrack::Quality SiPixelCompareTrackSoA< T >::minQuality_
private

Definition at line 82 of file SiPixelCompareTrackSoA.cc.

◆ tokenSoATrackCPU_

template<typename T >
const edm::EDGetTokenT<PixelTrackSoA> SiPixelCompareTrackSoA< T >::tokenSoATrackCPU_
private

Definition at line 78 of file SiPixelCompareTrackSoA.cc.

◆ tokenSoATrackGPU_

template<typename T >
const edm::EDGetTokenT<PixelTrackSoA> SiPixelCompareTrackSoA< T >::tokenSoATrackGPU_
private

Definition at line 79 of file SiPixelCompareTrackSoA.cc.

◆ topFolderName_

template<typename T >
const std::string SiPixelCompareTrackSoA< T >::topFolderName_
private

Definition at line 80 of file SiPixelCompareTrackSoA.cc.

◆ useQualityCut_

template<typename T >
const bool SiPixelCompareTrackSoA< T >::useQualityCut_
private

Definition at line 81 of file SiPixelCompareTrackSoA.cc.