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_
 
MonitorElementhCurvdiffMatched_
 
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_
 
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 124 of file SiPixelCompareTrackSoA.cc.

125  : tokenSoATrackCPU_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcCPU"))),
126  tokenSoATrackGPU_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcGPU"))),
127  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
128  useQualityCut_(iConfig.getParameter<bool>("useQualityCut")),
130  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 136 of file SiPixelCompareTrackSoA.cc.

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

136  {
137  using helper = TracksUtilities<T>;
138  const auto& tsoaHandleCPU = iEvent.getHandle(tokenSoATrackCPU_);
139  const auto& tsoaHandleGPU = iEvent.getHandle(tokenSoATrackGPU_);
140  if (not tsoaHandleCPU or not tsoaHandleGPU) {
141  edm::LogWarning out("SiPixelCompareTrackSoA");
142  if (not tsoaHandleCPU) {
143  out << "reference (cpu) tracks not found; ";
144  }
145  if (not tsoaHandleGPU) {
146  out << "target (gpu) tracks not found; ";
147  }
148  out << "the comparison will not run.";
149  return;
150  }
151 
152  auto const& tsoaCPU = *tsoaHandleCPU;
153  auto const& tsoaGPU = *tsoaHandleGPU;
154  auto maxTracksCPU = tsoaCPU.view().metadata().size(); //this should be same for both?
155  auto maxTracksGPU = tsoaGPU.view().metadata().size(); //this should be same for both?
156  auto const* qualityCPU = tsoaCPU.view().quality();
157  auto const* qualityGPU = tsoaGPU.view().quality();
158  int32_t nTracksCPU = 0;
159  int32_t nTracksGPU = 0;
160  int32_t nLooseAndAboveTracksCPU = 0;
161  int32_t nLooseAndAboveTracksCPU_matchedGPU = 0;
162  int32_t nLooseAndAboveTracksGPU = 0;
163 
164  //Loop over GPU tracks and store the indices of the loose tracks. Whats happens if useQualityCut_ is false?
165  std::vector<int32_t> looseTrkidxGPU;
166  for (int32_t jt = 0; jt < maxTracksGPU; ++jt) {
167  if (helper::nHits(tsoaGPU.view(), jt) == 0)
168  break; // this is a guard
169  if (!(tsoaGPU.view()[jt].pt() > 0.))
170  continue;
171  nTracksGPU++;
172  if (useQualityCut_ && qualityGPU[jt] < minQuality_)
173  continue;
174  nLooseAndAboveTracksGPU++;
175  looseTrkidxGPU.emplace_back(jt);
176  }
177 
178  //Now loop over CPU tracks//nested loop for loose gPU tracks
179  for (int32_t it = 0; it < maxTracksCPU; ++it) {
180  int nHitsCPU = helper::nHits(tsoaCPU.view(), it);
181 
182  if (nHitsCPU == 0)
183  break; // this is a guard
184 
185  float ptCPU = tsoaCPU.view()[it].pt();
186  float etaCPU = tsoaCPU.view()[it].eta();
187  float phiCPU = helper::phi(tsoaCPU.view(), it);
188  float zipCPU = helper::zip(tsoaCPU.view(), it);
189  float tipCPU = helper::tip(tsoaCPU.view(), it);
190 
191  if (!(ptCPU > 0.))
192  continue;
193  nTracksCPU++;
194  if (useQualityCut_ && qualityCPU[it] < minQuality_)
195  continue;
196  nLooseAndAboveTracksCPU++;
197  //Now loop over loose GPU trk and find the closest in DeltaR//do we need pt cut?
198  const int32_t notFound = -1;
199  int32_t closestTkidx = notFound;
200  float mindr2 = dr2cut_;
201 
202  for (auto gid : looseTrkidxGPU) {
203  float etaGPU = tsoaGPU.view()[gid].eta();
204  float phiGPU = helper::phi(tsoaGPU.view(), gid);
205  float dr2 = reco::deltaR2(etaCPU, phiCPU, etaGPU, phiGPU);
206  if (dr2 > dr2cut_)
207  continue; // this is arbitrary
208  if (mindr2 > dr2) {
209  mindr2 = dr2;
210  closestTkidx = gid;
211  }
212  }
213 
214  hpt_eta_tkAllCPU_->Fill(etaCPU, ptCPU); //all CPU tk
215  hphi_z_tkAllCPU_->Fill(phiCPU, zipCPU);
216  if (closestTkidx == notFound)
217  continue;
218  nLooseAndAboveTracksCPU_matchedGPU++;
219 
220  hchi2_->Fill(tsoaCPU.view()[it].chi2(), tsoaGPU.view()[closestTkidx].chi2());
221  hCharge_->Fill(helper::charge(tsoaCPU.view(), it), helper::charge(tsoaGPU.view(), closestTkidx));
222  hnHits_->Fill(helper::nHits(tsoaCPU.view(), it), helper::nHits(tsoaGPU.view(), closestTkidx));
223  hnLayers_->Fill(tsoaCPU.view()[it].nLayers(), tsoaGPU.view()[closestTkidx].nLayers());
224  hpt_->Fill(ptCPU, tsoaGPU.view()[closestTkidx].pt());
225  hptLogLog_->Fill(ptCPU, tsoaGPU.view()[closestTkidx].pt());
226  heta_->Fill(etaCPU, tsoaGPU.view()[closestTkidx].eta());
227  hphi_->Fill(phiCPU, helper::phi(tsoaGPU.view(), closestTkidx));
228  hz_->Fill(zipCPU, helper::zip(tsoaGPU.view(), closestTkidx));
229  htip_->Fill(tipCPU, helper::tip(tsoaGPU.view(), closestTkidx));
230  hptdiffMatched_->Fill(ptCPU - tsoaGPU.view()[closestTkidx].pt());
231  hCurvdiffMatched_->Fill((helper::charge(tsoaCPU.view(), it) / tsoaCPU.view()[it].pt()) -
232  (helper::charge(tsoaGPU.view(), closestTkidx) / tsoaGPU.view()[closestTkidx].pt()));
233  hetadiffMatched_->Fill(etaCPU - tsoaGPU.view()[closestTkidx].eta());
234  hphidiffMatched_->Fill(reco::deltaPhi(phiCPU, helper::phi(tsoaGPU.view(), closestTkidx)));
235  hzdiffMatched_->Fill(zipCPU - helper::zip(tsoaGPU.view(), closestTkidx));
236  htipdiffMatched_->Fill(tipCPU - helper::tip(tsoaGPU.view(), closestTkidx));
237  hpt_eta_tkAllCPUMatched_->Fill(etaCPU, tsoaCPU.view()[it].pt()); //matched to gpu
238  hphi_z_tkAllCPUMatched_->Fill(etaCPU, zipCPU);
239  }
240  hnTracks_->Fill(nTracksCPU, nTracksGPU);
241  hnLooseAndAboveTracks_->Fill(nLooseAndAboveTracksCPU, nLooseAndAboveTracksGPU);
242  hnLooseAndAboveTracks_matched_->Fill(nLooseAndAboveTracksCPU, nLooseAndAboveTracksCPU_matchedGPU);
243 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Definition: helper.py:1
MonitorElement * hphi_z_tkAllCPUMatched_
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackGPU_
MonitorElement * hnLooseAndAboveTracks_matched_
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
OutputIterator zip(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, InputIterator2 last2, OutputIterator result, Compare comp)
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 * hpt_eta_tkAllCPUMatched_
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 249 of file SiPixelCompareTrackSoA.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::book2I(), dqm::implementation::NavigatorBase::cd(), M_PI, dqm::implementation::NavigatorBase::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

251  {
252  iBook.cd();
253  iBook.setCurrentFolder(topFolderName_);
254 
255  // clang-format off
256  std::string toRep = "Number of tracks";
257  // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports THnSparse
258  // these should be moved to a less resource consuming format
259  hnTracks_ = iBook.book2I("nTracks", fmt::sprintf("%s per event; CPU; GPU",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
260  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);
261  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);
262 
263  toRep = "Number of all RecHits per track (quality #geq loose)";
264  hnHits_ = iBook.book2I("nRecHits", fmt::sprintf("%s;CPU;GPU",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
265 
266  toRep = "Number of all layers per track (quality #geq loose)";
267  hnLayers_ = iBook.book2I("nLayers", fmt::sprintf("%s;CPU;GPU",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
268 
269  toRep = "Track (quality #geq loose) #chi^{2}/ndof";
270  hchi2_ = iBook.book2I("nChi2ndof", fmt::sprintf("%s;CPU;GPU",toRep), 40, 0., 20., 40, 0., 20.);
271 
272  toRep = "Track (quality #geq loose) charge";
273  hCharge_ = iBook.book2I("charge",fmt::sprintf("%s;CPU;GPU",toRep),3, -1.5, 1.5, 3, -1.5, 1.5);
274 
275  hpt_ = iBook.book2I("pt", "Track (quality #geq loose) p_{T} [GeV];CPU;GPU", 200, 0., 200., 200, 0., 200.);
276  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.));
277  heta_ = iBook.book2I("eta", "Track (quality #geq loose) #eta;CPU;GPU", 30, -3., 3., 30, -3., 3.);
278  hphi_ = iBook.book2I("phi", "Track (quality #geq loose) #phi;CPU;GPU", 30, -M_PI, M_PI, 30, -M_PI, M_PI);
279  hz_ = iBook.book2I("z", "Track (quality #geq loose) z [cm];CPU;GPU", 30, -30., 30., 30, -30., 30.);
280  htip_ = iBook.book2I("tip", "Track (quality #geq loose) TIP [cm];CPU;GPU", 100, -0.5, 0.5, 100, -0.5, 0.5);
281  //1D difference plots
282  hptdiffMatched_ = iBook.book1D("ptdiffmatched", " p_{T} diff [GeV] between matched tracks; #Delta p_{T} [GeV]", 60, -30., 30.);
283  hCurvdiffMatched_ = iBook.book1D("curvdiffmatched", "q/p_{T} diff [GeV] between matched tracks; #Delta q/p_{T} [GeV]", 60, -30., 30.);
284  hetadiffMatched_ = iBook.book1D("etadiffmatched", " #eta diff between matched tracks; #Delta #eta", 160, -0.04 ,0.04);
285  hphidiffMatched_ = iBook.book1D("phidiffmatched", " #phi diff between matched tracks; #Delta #phi", 160, -0.04 ,0.04);
286  hzdiffMatched_ = iBook.book1D("zdiffmatched", " z diff between matched tracks; #Delta z [cm]", 300, -1.5, 1.5);
287  htipdiffMatched_ = iBook.book1D("tipdiffmatched", " TIP diff between matched tracks; #Delta TIP [cm]", 300, -1.5, 1.5);
288  //2D plots for eff
289  hpt_eta_tkAllCPU_ = iBook.book2I("ptetatrkAllCPU", "Track (quality #geq loose) on CPU; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
290  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.);
291 
292  hphi_z_tkAllCPU_ = iBook.book2I("phiztrkAllCPU", "Track (quality #geq loose) on CPU; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
293  hphi_z_tkAllCPUMatched_ = iBook.book2I("phiztrkAllCPUmatched", "Track (quality #geq loose) on CPU; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
294 
295 }
MonitorElement * hphi_z_tkAllCPUMatched_
const std::string topFolderName_
MonitorElement * hnLooseAndAboveTracks_matched_
#define M_PI
MonitorElement * hpt_eta_tkAllCPUMatched_
MonitorElement * hnLooseAndAboveTracks_

◆ fillDescriptions()

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

Definition at line 298 of file SiPixelCompareTrackSoA.cc.

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

298  {
299  // monitorpixelTrackSoA
301  desc.add<edm::InputTag>("pixelTrackSrcCPU", edm::InputTag("pixelTracksSoA@cpu"));
302  desc.add<edm::InputTag>("pixelTrackSrcGPU", edm::InputTag("pixelTracksSoA@cuda"));
303  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareGPUvsCPU");
304  desc.add<bool>("useQualityCut", true);
305  desc.add<std::string>("minQuality", "loose");
306  desc.add<double>("deltaR2cut", 0.04);
307  descriptions.addWithDefaultLabel(desc);
308 }
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 93 of file SiPixelCompareTrackSoA.cc.

◆ hchi2_

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

Definition at line 94 of file SiPixelCompareTrackSoA.cc.

◆ hChi2VsEta_

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

Definition at line 96 of file SiPixelCompareTrackSoA.cc.

◆ hChi2VsPhi_

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

Definition at line 95 of file SiPixelCompareTrackSoA.cc.

◆ hCurvdiffMatched_

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

Definition at line 106 of file SiPixelCompareTrackSoA.cc.

◆ heta_

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

Definition at line 99 of file SiPixelCompareTrackSoA.cc.

◆ hetadiffMatched_

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

Definition at line 107 of file SiPixelCompareTrackSoA.cc.

◆ hnHits_

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

Definition at line 87 of file SiPixelCompareTrackSoA.cc.

◆ hnHitsVsEta_

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

Definition at line 89 of file SiPixelCompareTrackSoA.cc.

◆ hnHitsVsPhi_

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

Definition at line 88 of file SiPixelCompareTrackSoA.cc.

◆ hnLayers_

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

Definition at line 90 of file SiPixelCompareTrackSoA.cc.

◆ hnLayersVsEta_

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

Definition at line 92 of file SiPixelCompareTrackSoA.cc.

◆ hnLayersVsPhi_

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

Definition at line 91 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 100 of file SiPixelCompareTrackSoA.cc.

◆ hphi_z_tkAllCPU_

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

Definition at line 115 of file SiPixelCompareTrackSoA.cc.

◆ hphi_z_tkAllCPUMatched_

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

Definition at line 116 of file SiPixelCompareTrackSoA.cc.

◆ hphidiffMatched_

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

Definition at line 108 of file SiPixelCompareTrackSoA.cc.

◆ hpt_

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

Definition at line 97 of file SiPixelCompareTrackSoA.cc.

◆ hpt_eta_tkAllCPU_

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

Definition at line 113 of file SiPixelCompareTrackSoA.cc.

◆ hpt_eta_tkAllCPUMatched_

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

Definition at line 114 of file SiPixelCompareTrackSoA.cc.

◆ hptdiffMatched_

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

Definition at line 105 of file SiPixelCompareTrackSoA.cc.

◆ hptLogLog_

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

Definition at line 98 of file SiPixelCompareTrackSoA.cc.

◆ hquality_

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

Definition at line 103 of file SiPixelCompareTrackSoA.cc.

◆ htip_

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

Definition at line 102 of file SiPixelCompareTrackSoA.cc.

◆ htipdiffMatched_

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

Definition at line 110 of file SiPixelCompareTrackSoA.cc.

◆ hz_

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

Definition at line 101 of file SiPixelCompareTrackSoA.cc.

◆ hzdiffMatched_

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

Definition at line 109 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.