CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Static Public Member Functions | Private Attributes
SiPixelCompareTrackSoAAlpaka< T > Class Template Reference
Inheritance diagram for SiPixelCompareTrackSoAAlpaka< T >:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Types

using PixelTrackSoA = TracksHost< T >
 
- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

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
 
 SiPixelCompareTrackSoAAlpaka (const edm::ParameterSet &)
 
 ~SiPixelCompareTrackSoAAlpaka () override=default
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 

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_tkAllHost_
 
MonitorElementhphi_z_tkAllHostMatched_
 
MonitorElementhphidiffMatched_
 
MonitorElementhpt_
 
MonitorElementhpt_eta_tkAllHost_
 
MonitorElementhpt_eta_tkAllHostMatched_
 
MonitorElementhptdiffMatched_
 
MonitorElementhptLogLog_
 
MonitorElementhquality_
 
MonitorElementhtip_
 
MonitorElementhtipdiffMatched_
 
MonitorElementhz_
 
MonitorElementhzdiffMatched_
 
const pixelTrack::Quality minQuality_
 
const edm::EDGetTokenT< PixelTrackSoAtokenSoATrackDevice_
 
const edm::EDGetTokenT< PixelTrackSoAtokenSoATrackHost_
 
const std::string topFolderName_
 
const bool useQualityCut_
 

Additional Inherited Members

- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

template<typename T>
class SiPixelCompareTrackSoAAlpaka< T >

Definition at line 59 of file SiPixelCompareTrackSoAAlpaka.cc.

Member Typedef Documentation

◆ PixelTrackSoA

template<typename T >
using SiPixelCompareTrackSoAAlpaka< T >::PixelTrackSoA = TracksHost<T>

Definition at line 61 of file SiPixelCompareTrackSoAAlpaka.cc.

Constructor & Destructor Documentation

◆ SiPixelCompareTrackSoAAlpaka()

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

Definition at line 116 of file SiPixelCompareTrackSoAAlpaka.cc.

117  : tokenSoATrackHost_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcHost"))),
118  tokenSoATrackDevice_(consumes<PixelTrackSoA>(iConfig.getParameter<edm::InputTag>("pixelTrackSrcDevice"))),
119  topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
120  useQualityCut_(iConfig.getParameter<bool>("useQualityCut")),
122  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 > tokenSoATrackHost_
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackDevice_

◆ ~SiPixelCompareTrackSoAAlpaka()

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

Member Function Documentation

◆ analyze()

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

Reimplemented from DQMEDAnalyzer.

Definition at line 128 of file SiPixelCompareTrackSoAAlpaka.cc.

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

128  {
129  using helper = TracksUtilities<T>;
130  const auto& tsoaHandleHost = iEvent.getHandle(tokenSoATrackHost_);
131  const auto& tsoaHandleDevice = iEvent.getHandle(tokenSoATrackDevice_);
132  if (not tsoaHandleHost or not tsoaHandleDevice) {
133  edm::LogWarning out("SiPixelCompareTrackSoAAlpaka");
134  if (not tsoaHandleHost) {
135  out << "reference (cpu) tracks not found; ";
136  }
137  if (not tsoaHandleDevice) {
138  out << "target (gpu) tracks not found; ";
139  }
140  out << "the comparison will not run.";
141  return;
142  }
143 
144  auto const& tsoaHost = *tsoaHandleHost;
145  auto const& tsoaDevice = *tsoaHandleDevice;
146  auto maxTracksHost = tsoaHost.view().metadata().size(); //this should be same for both?
147  auto maxTracksDevice = tsoaDevice.view().metadata().size(); //this should be same for both?
148  auto const* qualityHost = tsoaHost.view().quality();
149  auto const* qualityDevice = tsoaDevice.view().quality();
150  int32_t nTracksHost = 0;
151  int32_t nTracksDevice = 0;
152  int32_t nLooseAndAboveTracksHost = 0;
153  int32_t nLooseAndAboveTracksHost_matchedDevice = 0;
154  int32_t nLooseAndAboveTracksDevice = 0;
155 
156  //Loop over Device tracks and store the indices of the loose tracks. Whats happens if useQualityCut_ is false?
157  std::vector<int32_t> looseTrkidxDevice;
158  for (int32_t jt = 0; jt < maxTracksDevice; ++jt) {
159  if (helper::nHits(tsoaDevice.view(), jt) == 0)
160  break; // this is a guard
161  if (!(tsoaDevice.view()[jt].pt() > 0.))
162  continue;
163  nTracksDevice++;
164  if (useQualityCut_ && qualityDevice[jt] < minQuality_)
165  continue;
166  nLooseAndAboveTracksDevice++;
167  looseTrkidxDevice.emplace_back(jt);
168  }
169 
170  //Now loop over Host tracks//nested loop for loose gPU tracks
171  for (int32_t it = 0; it < maxTracksHost; ++it) {
172  int nHitsHost = helper::nHits(tsoaHost.view(), it);
173 
174  if (nHitsHost == 0)
175  break; // this is a guard
176 
177  float ptHost = tsoaHost.view()[it].pt();
178  float etaHost = tsoaHost.view()[it].eta();
179  float phiHost = reco::phi(tsoaHost.view(), it);
180  float zipHost = reco::zip(tsoaHost.view(), it);
181  float tipHost = reco::tip(tsoaHost.view(), it);
182 
183  if (!(ptHost > 0.))
184  continue;
185  nTracksHost++;
186  if (useQualityCut_ && qualityHost[it] < minQuality_)
187  continue;
188  nLooseAndAboveTracksHost++;
189  //Now loop over loose Device trk and find the closest in DeltaR//do we need pt cut?
190  const int32_t notFound = -1;
191  int32_t closestTkidx = notFound;
192  float mindr2 = dr2cut_;
193 
194  for (auto gid : looseTrkidxDevice) {
195  float etaDevice = tsoaDevice.view()[gid].eta();
196  float phiDevice = reco::phi(tsoaDevice.view(), gid);
197  float dr2 = reco::deltaR2(etaHost, phiHost, etaDevice, phiDevice);
198  if (dr2 > dr2cut_)
199  continue; // this is arbitrary
200  if (mindr2 > dr2) {
201  mindr2 = dr2;
202  closestTkidx = gid;
203  }
204  }
205 
206  hpt_eta_tkAllHost_->Fill(etaHost, ptHost); //all Host tk
207  hphi_z_tkAllHost_->Fill(phiHost, zipHost);
208  if (closestTkidx == notFound)
209  continue;
210  nLooseAndAboveTracksHost_matchedDevice++;
211 
212  hchi2_->Fill(tsoaHost.view()[it].chi2(), tsoaDevice.view()[closestTkidx].chi2());
213  hCharge_->Fill(reco::charge(tsoaHost.view(), it), reco::charge(tsoaDevice.view(), closestTkidx));
214  hnHits_->Fill(helper::nHits(tsoaHost.view(), it), helper::nHits(tsoaDevice.view(), closestTkidx));
215  hnLayers_->Fill(tsoaHost.view()[it].nLayers(), tsoaDevice.view()[closestTkidx].nLayers());
216  hpt_->Fill(tsoaHost.view()[it].pt(), tsoaDevice.view()[closestTkidx].pt());
217  hptLogLog_->Fill(tsoaHost.view()[it].pt(), tsoaDevice.view()[closestTkidx].pt());
218  heta_->Fill(etaHost, tsoaDevice.view()[closestTkidx].eta());
219  hphi_->Fill(phiHost, reco::phi(tsoaDevice.view(), closestTkidx));
220  hz_->Fill(zipHost, reco::zip(tsoaDevice.view(), closestTkidx));
221  htip_->Fill(tipHost, reco::tip(tsoaDevice.view(), closestTkidx));
222  hptdiffMatched_->Fill(ptHost - tsoaDevice.view()[closestTkidx].pt());
223  hCurvdiffMatched_->Fill((reco::charge(tsoaHost.view(), it) / tsoaHost.view()[it].pt()) -
224  (reco::charge(tsoaDevice.view(), closestTkidx) / tsoaDevice.view()[closestTkidx].pt()));
225  hetadiffMatched_->Fill(etaHost - tsoaDevice.view()[closestTkidx].eta());
226  hphidiffMatched_->Fill(reco::deltaPhi(phiHost, reco::phi(tsoaDevice.view(), closestTkidx)));
227  hzdiffMatched_->Fill(zipHost - reco::zip(tsoaDevice.view(), closestTkidx));
228  htipdiffMatched_->Fill(tipHost - reco::tip(tsoaDevice.view(), closestTkidx));
229  hpt_eta_tkAllHostMatched_->Fill(etaHost, tsoaHost.view()[it].pt()); //matched to gpu
230  hphi_z_tkAllHostMatched_->Fill(etaHost, zipHost);
231  }
232  hnTracks_->Fill(nTracksHost, nTracksDevice);
233  hnLooseAndAboveTracks_->Fill(nLooseAndAboveTracksHost, nLooseAndAboveTracksDevice);
234  hnLooseAndAboveTracks_matched_->Fill(nLooseAndAboveTracksHost, nLooseAndAboveTracksHost_matchedDevice);
235 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Definition: helper.py:1
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float charge(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:73
void Fill(long long x)
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
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
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float phi(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:80
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackHost_
static const GlobalPoint notFound(0, 0, 0)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackDevice_

◆ bookHistograms()

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

Implements DQMEDAnalyzer.

Definition at line 241 of file SiPixelCompareTrackSoAAlpaka.cc.

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

243  {
244  iBook.cd();
245  iBook.setCurrentFolder(topFolderName_);
246 
247  // clang-format off
248  std::string toRep = "Number of tracks";
249  // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports THnSparse
250  // these should be moved to a less resource consuming format
251  hnTracks_ = iBook.book2I("nTracks", fmt::format("{} per event; Host; Device",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
252  hnLooseAndAboveTracks_ = iBook.book2I("nLooseAndAboveTracks", fmt::format("{} (quality #geq loose) per event; Host; Device",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
253  hnLooseAndAboveTracks_matched_ = iBook.book2I("nLooseAndAboveTracks_matched", fmt::format("{} (quality #geq loose) per event; Host; Device",toRep), 501, -0.5, 500.5, 501, -0.5, 500.5);
254 
255  toRep = "Number of all RecHits per track (quality #geq loose)";
256  hnHits_ = iBook.book2I("nRecHits", fmt::format("{};Host;Device",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
257 
258  toRep = "Number of all layers per track (quality #geq loose)";
259  hnLayers_ = iBook.book2I("nLayers", fmt::format("{};Host;Device",toRep), 15, -0.5, 14.5, 15, -0.5, 14.5);
260 
261  toRep = "Track (quality #geq loose) #chi^{2}/ndof";
262  hchi2_ = iBook.book2I("nChi2ndof", fmt::format("{};Host;Device",toRep), 40, 0., 20., 40, 0., 20.);
263 
264  toRep = "Track (quality #geq loose) charge";
265  hCharge_ = iBook.book2I("charge",fmt::format("{};Host;Device",toRep),3, -1.5, 1.5, 3, -1.5, 1.5);
266 
267  hpt_ = iBook.book2I("pt", "Track (quality #geq loose) p_{T} [GeV];Host;Device", 200, 0., 200., 200, 0., 200.);
268  hptLogLog_ = make2DIfLog(iBook, true, true, "ptLogLog", "Track (quality #geq loose) p_{T} [GeV];Host;Device", 200, log10(0.5), log10(200.), 200, log10(0.5), log10(200.));
269  heta_ = iBook.book2I("eta", "Track (quality #geq loose) #eta;Host;Device", 30, -3., 3., 30, -3., 3.);
270  hphi_ = iBook.book2I("phi", "Track (quality #geq loose) #phi;Host;Device", 30, -M_PI, M_PI, 30, -M_PI, M_PI);
271  hz_ = iBook.book2I("z", "Track (quality #geq loose) z [cm];Host;Device", 30, -30., 30., 30, -30., 30.);
272  htip_ = iBook.book2I("tip", "Track (quality #geq loose) TIP [cm];Host;Device", 100, -0.5, 0.5, 100, -0.5, 0.5);
273  //1D difference plots
274  hptdiffMatched_ = iBook.book1D("ptdiffmatched", " p_{T} diff [GeV] between matched tracks; #Delta p_{T} [GeV]", 60, -30., 30.);
275  hCurvdiffMatched_ = iBook.book1D("curvdiffmatched", "q/p_{T} diff [GeV] between matched tracks; #Delta q/p_{T} [GeV]", 60, -30., 30.);
276  hetadiffMatched_ = iBook.book1D("etadiffmatched", " #eta diff between matched tracks; #Delta #eta", 160, -0.04 ,0.04);
277  hphidiffMatched_ = iBook.book1D("phidiffmatched", " #phi diff between matched tracks; #Delta #phi", 160, -0.04 ,0.04);
278  hzdiffMatched_ = iBook.book1D("zdiffmatched", " z diff between matched tracks; #Delta z [cm]", 300, -1.5, 1.5);
279  htipdiffMatched_ = iBook.book1D("tipdiffmatched", " TIP diff between matched tracks; #Delta TIP [cm]", 300, -1.5, 1.5);
280  //2D plots for eff
281  hpt_eta_tkAllHost_ = iBook.book2I("ptetatrkAllHost", "Track (quality #geq loose) on Host; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
282  hpt_eta_tkAllHostMatched_ = iBook.book2I("ptetatrkAllHostmatched", "Track (quality #geq loose) on Host matched to Device track; #eta; p_{T} [GeV];", 30, -M_PI, M_PI, 200, 0., 200.);
283 
284  hphi_z_tkAllHost_ = iBook.book2I("phiztrkAllHost", "Track (quality #geq loose) on Host; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
285  hphi_z_tkAllHostMatched_ = iBook.book2I("phiztrkAllHostmatched", "Track (quality #geq loose) on Host; #phi; z [cm];", 30, -M_PI, M_PI, 30, -30., 30.);
286 
287 }
#define M_PI

◆ fillDescriptions()

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

Definition at line 290 of file SiPixelCompareTrackSoAAlpaka.cc.

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

290  {
291  // monitorpixelTrackSoA
293  desc.add<edm::InputTag>("pixelTrackSrcHost", edm::InputTag("pixelTracksAlpakaSerial"));
294  desc.add<edm::InputTag>("pixelTrackSrcDevice", edm::InputTag("pixelTracksAlpaka"));
295  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackCompareDeviceVSHost");
296  desc.add<bool>("useQualityCut", true);
297  desc.add<std::string>("minQuality", "loose");
298  desc.add<double>("deltaR2cut", 0.04);
299  descriptions.addWithDefaultLabel(desc);
300 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

Member Data Documentation

◆ dr2cut_

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

Definition at line 75 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hCharge_

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

Definition at line 85 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hchi2_

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

Definition at line 86 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hChi2VsEta_

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

Definition at line 88 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hChi2VsPhi_

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

Definition at line 87 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hCurvdiffMatched_

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

Definition at line 98 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ heta_

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

Definition at line 91 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hetadiffMatched_

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

Definition at line 99 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnHits_

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

Definition at line 79 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnHitsVsEta_

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

Definition at line 81 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnHitsVsPhi_

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

Definition at line 80 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnLayers_

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

Definition at line 82 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnLayersVsEta_

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

Definition at line 84 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnLayersVsPhi_

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

Definition at line 83 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnLooseAndAboveTracks_

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

Definition at line 77 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnLooseAndAboveTracks_matched_

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

Definition at line 78 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hnTracks_

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

Definition at line 76 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hphi_

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

Definition at line 92 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hphi_z_tkAllHost_

template<typename T >
MonitorElement* SiPixelCompareTrackSoAAlpaka< T >::hphi_z_tkAllHost_
private

Definition at line 107 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hphi_z_tkAllHostMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoAAlpaka< T >::hphi_z_tkAllHostMatched_
private

Definition at line 108 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hphidiffMatched_

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

Definition at line 100 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hpt_

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

Definition at line 89 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hpt_eta_tkAllHost_

template<typename T >
MonitorElement* SiPixelCompareTrackSoAAlpaka< T >::hpt_eta_tkAllHost_
private

Definition at line 105 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hpt_eta_tkAllHostMatched_

template<typename T >
MonitorElement* SiPixelCompareTrackSoAAlpaka< T >::hpt_eta_tkAllHostMatched_
private

Definition at line 106 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hptdiffMatched_

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

Definition at line 97 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hptLogLog_

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

Definition at line 90 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hquality_

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

Definition at line 95 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ htip_

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

Definition at line 94 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ htipdiffMatched_

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

Definition at line 102 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hz_

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

Definition at line 93 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ hzdiffMatched_

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

Definition at line 101 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ minQuality_

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

Definition at line 74 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ tokenSoATrackDevice_

template<typename T >
const edm::EDGetTokenT<PixelTrackSoA> SiPixelCompareTrackSoAAlpaka< T >::tokenSoATrackDevice_
private

Definition at line 71 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ tokenSoATrackHost_

template<typename T >
const edm::EDGetTokenT<PixelTrackSoA> SiPixelCompareTrackSoAAlpaka< T >::tokenSoATrackHost_
private

Definition at line 70 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ topFolderName_

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

Definition at line 72 of file SiPixelCompareTrackSoAAlpaka.cc.

◆ useQualityCut_

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

Definition at line 73 of file SiPixelCompareTrackSoAAlpaka.cc.