CMS 3D CMS Logo

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

#include <SiPixelCompareTracks.cc>

Inheritance diagram for SiPixelCompareTracks< 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
 
template<typename U , typename V >
void analyzeSeparate (U tokenRef, V tokenTar, const edm::Event &iEvent)
 
void bookHistograms (DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
 
 SiPixelCompareTracks (const edm::ParameterSet &)
 
 ~SiPixelCompareTracks () 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_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< PixelTrackSoAtokenSoATrackReference_
 
const edm::EDGetTokenT< PixelTrackSoAtokenSoATrackTarget_
 
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 SiPixelCompareTracks< T >

Definition at line 72 of file SiPixelCompareTracks.cc.

Member Typedef Documentation

◆ PixelTrackSoA

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

Definition at line 74 of file SiPixelCompareTracks.cc.

Constructor & Destructor Documentation

◆ SiPixelCompareTracks()

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

Definition at line 134 of file SiPixelCompareTracks.cc.

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")),
140  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 > tokenSoATrackReference_
const pixelTrack::Quality minQuality_
const std::string topFolderName_
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackTarget_

◆ ~SiPixelCompareTracks()

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

Member Function Documentation

◆ analyze()

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

Reimplemented from DQMEDAnalyzer.

Definition at line 262 of file SiPixelCompareTracks.cc.

References iEvent.

262  {
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
266 }
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackReference_
void analyzeSeparate(U tokenRef, V tokenTar, const edm::Event &iEvent)
int iEvent
Definition: GenABIO.cc:224
const edm::EDGetTokenT< PixelTrackSoA > tokenSoATrackTarget_

◆ analyzeSeparate()

template<typename T >
template<typename U , typename V >
void SiPixelCompareTracks< T >::analyzeSeparate ( tokenRef,
tokenTar,
const edm::Event iEvent 
)

Definition at line 144 of file SiPixelCompareTracks.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().

144  {
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 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
MonitorElement * hnLooseAndAboveTracks_
Definition: helper.py:1
MonitorElement * hptdiffMatched_
MonitorElement * hphi_z_tkAllRef_
ALPAKA_FN_HOST_ACC ALPAKA_FN_INLINE constexpr float zip(ConstView const &tracks, int32_t i)
Definition: TracksSoA.h:90
MonitorElement * hphi_z_tkAllRefMatched_
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_
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
MonitorElement * hpt_eta_tkAllRef_
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
MonitorElement * hzdiffMatched_
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_
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits

◆ bookHistograms()

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

Implements DQMEDAnalyzer.

Definition at line 272 of file SiPixelCompareTracks.cc.

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

274  {
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 }
MonitorElement * hnLooseAndAboveTracks_
MonitorElement * hptdiffMatched_
MonitorElement * hphi_z_tkAllRef_
MonitorElement * hphi_z_tkAllRefMatched_
MonitorElement * hCurvdiffMatched_
MonitorElement * hnTracks_
MonitorElement * hpt_eta_tkAllRefMatched_
MonitorElement * hnLooseAndAboveTracks_matched_
MonitorElement * hphidiffMatched_
MonitorElement * hpt_eta_tkAllRef_
MonitorElement * hzdiffMatched_
#define M_PI
MonitorElement * hetadiffMatched_
MonitorElement * htipdiffMatched_
const std::string topFolderName_

◆ fillDescriptions()

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

Definition at line 321 of file SiPixelCompareTracks.cc.

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

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

Member Data Documentation

◆ dr2cut_

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

Definition at line 93 of file SiPixelCompareTracks.cc.

◆ hCharge_

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

Definition at line 103 of file SiPixelCompareTracks.cc.

◆ hchi2_

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

Definition at line 104 of file SiPixelCompareTracks.cc.

◆ hChi2VsEta_

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

Definition at line 106 of file SiPixelCompareTracks.cc.

◆ hChi2VsPhi_

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

Definition at line 105 of file SiPixelCompareTracks.cc.

◆ hCurvdiffMatched_

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

Definition at line 116 of file SiPixelCompareTracks.cc.

◆ heta_

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

Definition at line 109 of file SiPixelCompareTracks.cc.

◆ hetadiffMatched_

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

Definition at line 117 of file SiPixelCompareTracks.cc.

◆ hnHits_

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

Definition at line 97 of file SiPixelCompareTracks.cc.

◆ hnHitsVsEta_

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

Definition at line 99 of file SiPixelCompareTracks.cc.

◆ hnHitsVsPhi_

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

Definition at line 98 of file SiPixelCompareTracks.cc.

◆ hnLayers_

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

Definition at line 100 of file SiPixelCompareTracks.cc.

◆ hnLayersVsEta_

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

Definition at line 102 of file SiPixelCompareTracks.cc.

◆ hnLayersVsPhi_

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

Definition at line 101 of file SiPixelCompareTracks.cc.

◆ hnLooseAndAboveTracks_

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

Definition at line 95 of file SiPixelCompareTracks.cc.

◆ hnLooseAndAboveTracks_matched_

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

Definition at line 96 of file SiPixelCompareTracks.cc.

◆ hnTracks_

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

Definition at line 94 of file SiPixelCompareTracks.cc.

◆ hphi_

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

Definition at line 110 of file SiPixelCompareTracks.cc.

◆ hphi_z_tkAllRef_

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

Definition at line 125 of file SiPixelCompareTracks.cc.

◆ hphi_z_tkAllRefMatched_

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

Definition at line 126 of file SiPixelCompareTracks.cc.

◆ hphidiffMatched_

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

Definition at line 118 of file SiPixelCompareTracks.cc.

◆ hpt_

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

Definition at line 107 of file SiPixelCompareTracks.cc.

◆ hpt_eta_tkAllRef_

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

Definition at line 123 of file SiPixelCompareTracks.cc.

◆ hpt_eta_tkAllRefMatched_

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

Definition at line 124 of file SiPixelCompareTracks.cc.

◆ hptdiffMatched_

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

Definition at line 115 of file SiPixelCompareTracks.cc.

◆ hptLogLog_

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

Definition at line 108 of file SiPixelCompareTracks.cc.

◆ hquality_

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

Definition at line 113 of file SiPixelCompareTracks.cc.

◆ htip_

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

Definition at line 112 of file SiPixelCompareTracks.cc.

◆ htipdiffMatched_

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

Definition at line 120 of file SiPixelCompareTracks.cc.

◆ hz_

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

Definition at line 111 of file SiPixelCompareTracks.cc.

◆ hzdiffMatched_

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

Definition at line 119 of file SiPixelCompareTracks.cc.

◆ minQuality_

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

Definition at line 92 of file SiPixelCompareTracks.cc.

◆ tokenSoATrackReference_

template<typename T >
const edm::EDGetTokenT<PixelTrackSoA> SiPixelCompareTracks< T >::tokenSoATrackReference_
private

Definition at line 88 of file SiPixelCompareTracks.cc.

◆ tokenSoATrackTarget_

template<typename T >
const edm::EDGetTokenT<PixelTrackSoA> SiPixelCompareTracks< T >::tokenSoATrackTarget_
private

Definition at line 89 of file SiPixelCompareTracks.cc.

◆ topFolderName_

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

Definition at line 90 of file SiPixelCompareTracks.cc.

◆ useQualityCut_

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

Definition at line 91 of file SiPixelCompareTracks.cc.