CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
SiPixelPhase1CompareTrackSoA Class Reference

#include <SiPixelPhase1CompareTrackSoA.cc>

Inheritance diagram for SiPixelPhase1CompareTrackSoA:

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
 
 SiPixelPhase1CompareTrackSoA (const edm::ParameterSet &)
 
 ~SiPixelPhase1CompareTrackSoA () override=default
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

const float dr2cut_
 
MonitorElementhchi2_
 
MonitorElementhChi2VsEta_
 
MonitorElementhChi2VsPhi_
 
MonitorElementheta_
 
MonitorElementhetadiffMatched_
 
MonitorElementhnHits_
 
MonitorElementhnHitsVsEta_
 
MonitorElementhnHitsVsPhi_
 
MonitorElementhnLayers_
 
MonitorElementhnLayersVsEta_
 
MonitorElementhnLayersVsPhi_
 
MonitorElementhnLooseAndAboveTracks_
 
MonitorElementhnLooseAndAboveTracks_matched_
 
MonitorElementhnTracks_
 
MonitorElementhphi_
 
MonitorElementhphi_z_tkAllCPU_
 
MonitorElementhphi_z_tkAllCPUMatched_
 
MonitorElementhphidiffMatched_
 
MonitorElementhpt_
 
MonitorElementhpt_eta_tkAllCPU_
 
MonitorElementhpt_eta_tkAllCPUMatched_
 
MonitorElementhptdiffMatched_
 
MonitorElementhptLogLog_
 
MonitorElementhquality_
 
MonitorElementhtip_
 
MonitorElementhz_
 
MonitorElementhzdiffMatched_
 
const pixelTrack::Quality minQuality_
 
const edm::EDGetTokenT< PixelTrackHeterogeneoustokenSoATrackCPU_
 
const edm::EDGetTokenT< PixelTrackHeterogeneoustokenSoATrackGPU_
 
const std::string topFolderName_
 
const bool useQualityCut_
 

Detailed Description

Definition at line 65 of file SiPixelPhase1CompareTrackSoA.cc.

Constructor & Destructor Documentation

◆ SiPixelPhase1CompareTrackSoA()

SiPixelPhase1CompareTrackSoA::SiPixelPhase1CompareTrackSoA ( const edm::ParameterSet iConfig)
explicit

Definition at line 115 of file SiPixelPhase1CompareTrackSoA.cc.

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")),
121  dr2cut_(iConfig.getParameter<double>("deltaR2cut")) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Quality qualityByName(std::string const &name)
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrackGPU_
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrackCPU_

◆ ~SiPixelPhase1CompareTrackSoA()

SiPixelPhase1CompareTrackSoA::~SiPixelPhase1CompareTrackSoA ( )
overridedefault

Member Function Documentation

◆ analyze()

void SiPixelPhase1CompareTrackSoA::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 126 of file SiPixelPhase1CompareTrackSoA.cc.

References reco::deltaPhi(), reco::deltaR2(), dr2cut_, dqm::impl::MonitorElement::Fill(), hchi2_, heta_, hetadiffMatched_, hnHits_, hnLayers_, hnLooseAndAboveTracks_, hnLooseAndAboveTracks_matched_, hnTracks_, hphi_, hphi_z_tkAllCPU_, hphi_z_tkAllCPUMatched_, hphidiffMatched_, hpt_, hpt_eta_tkAllCPU_, hpt_eta_tkAllCPUMatched_, hptdiffMatched_, hptLogLog_, htip_, hz_, hzdiffMatched_, iEvent, minQuality_, notFound, or, MillePedeFileConverter_cfg::out, tokenSoATrackCPU_, tokenSoATrackGPU_, and useQualityCut_.

126  {
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 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrackGPU_
const edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrackCPU_
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
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
static const GlobalPoint notFound(0, 0, 0)

◆ bookHistograms()

void SiPixelPhase1CompareTrackSoA::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
override

Definition at line 225 of file SiPixelPhase1CompareTrackSoA.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2I(), dqm::implementation::NavigatorBase::cd(), hchi2_, heta_, hetadiffMatched_, hnHits_, hnLayers_, hnLooseAndAboveTracks_, hnLooseAndAboveTracks_matched_, hnTracks_, hphi_, hphi_z_tkAllCPU_, hphi_z_tkAllCPUMatched_, hphidiffMatched_, hpt_, hpt_eta_tkAllCPU_, hpt_eta_tkAllCPUMatched_, hptdiffMatched_, hptLogLog_, htip_, hz_, hzdiffMatched_, M_PI, dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, and topFolderName_.

227  {
228  iBook.cd();
229  iBook.setCurrentFolder(topFolderName_);
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 }
#define M_PI

◆ fillDescriptions()

void SiPixelPhase1CompareTrackSoA::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 268 of file SiPixelPhase1CompareTrackSoA.cc.

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

268  {
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 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ dr2cut_

const float SiPixelPhase1CompareTrackSoA::dr2cut_
private

Definition at line 79 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze().

◆ hchi2_

MonitorElement* SiPixelPhase1CompareTrackSoA::hchi2_
private

Definition at line 89 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hChi2VsEta_

MonitorElement* SiPixelPhase1CompareTrackSoA::hChi2VsEta_
private

Definition at line 91 of file SiPixelPhase1CompareTrackSoA.cc.

◆ hChi2VsPhi_

MonitorElement* SiPixelPhase1CompareTrackSoA::hChi2VsPhi_
private

Definition at line 90 of file SiPixelPhase1CompareTrackSoA.cc.

◆ heta_

MonitorElement* SiPixelPhase1CompareTrackSoA::heta_
private

Definition at line 94 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hetadiffMatched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hetadiffMatched_
private

Definition at line 101 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hnHits_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnHits_
private

Definition at line 83 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hnHitsVsEta_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnHitsVsEta_
private

Definition at line 85 of file SiPixelPhase1CompareTrackSoA.cc.

◆ hnHitsVsPhi_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnHitsVsPhi_
private

Definition at line 84 of file SiPixelPhase1CompareTrackSoA.cc.

◆ hnLayers_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnLayers_
private

Definition at line 86 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hnLayersVsEta_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnLayersVsEta_
private

Definition at line 88 of file SiPixelPhase1CompareTrackSoA.cc.

◆ hnLayersVsPhi_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnLayersVsPhi_
private

Definition at line 87 of file SiPixelPhase1CompareTrackSoA.cc.

◆ hnLooseAndAboveTracks_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnLooseAndAboveTracks_
private

Definition at line 81 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hnLooseAndAboveTracks_matched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnLooseAndAboveTracks_matched_
private

Definition at line 82 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hnTracks_

MonitorElement* SiPixelPhase1CompareTrackSoA::hnTracks_
private

Definition at line 80 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hphi_

MonitorElement* SiPixelPhase1CompareTrackSoA::hphi_
private

Definition at line 95 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hphi_z_tkAllCPU_

MonitorElement* SiPixelPhase1CompareTrackSoA::hphi_z_tkAllCPU_
private

Definition at line 107 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hphi_z_tkAllCPUMatched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hphi_z_tkAllCPUMatched_
private

Definition at line 108 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hphidiffMatched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hphidiffMatched_
private

Definition at line 102 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hpt_

MonitorElement* SiPixelPhase1CompareTrackSoA::hpt_
private

Definition at line 92 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hpt_eta_tkAllCPU_

MonitorElement* SiPixelPhase1CompareTrackSoA::hpt_eta_tkAllCPU_
private

Definition at line 105 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hpt_eta_tkAllCPUMatched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hpt_eta_tkAllCPUMatched_
private

Definition at line 106 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hptdiffMatched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hptdiffMatched_
private

Definition at line 100 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hptLogLog_

MonitorElement* SiPixelPhase1CompareTrackSoA::hptLogLog_
private

Definition at line 93 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hquality_

MonitorElement* SiPixelPhase1CompareTrackSoA::hquality_
private

Definition at line 98 of file SiPixelPhase1CompareTrackSoA.cc.

◆ htip_

MonitorElement* SiPixelPhase1CompareTrackSoA::htip_
private

Definition at line 97 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hz_

MonitorElement* SiPixelPhase1CompareTrackSoA::hz_
private

Definition at line 96 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ hzdiffMatched_

MonitorElement* SiPixelPhase1CompareTrackSoA::hzdiffMatched_
private

Definition at line 103 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

◆ minQuality_

const pixelTrack::Quality SiPixelPhase1CompareTrackSoA::minQuality_
private

Definition at line 78 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze().

◆ tokenSoATrackCPU_

const edm::EDGetTokenT<PixelTrackHeterogeneous> SiPixelPhase1CompareTrackSoA::tokenSoATrackCPU_
private

Definition at line 74 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze().

◆ tokenSoATrackGPU_

const edm::EDGetTokenT<PixelTrackHeterogeneous> SiPixelPhase1CompareTrackSoA::tokenSoATrackGPU_
private

Definition at line 75 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze().

◆ topFolderName_

const std::string SiPixelPhase1CompareTrackSoA::topFolderName_
private

Definition at line 76 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by bookHistograms().

◆ useQualityCut_

const bool SiPixelPhase1CompareTrackSoA::useQualityCut_
private

Definition at line 77 of file SiPixelPhase1CompareTrackSoA.cc.

Referenced by analyze().