CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
SiPixelPhase1MonitorTrackSoA Class Reference

bookLayer More...

#include <SiPixelPhase1MonitorTrackSoA.cc>

Inheritance diagram for SiPixelPhase1MonitorTrackSoA:

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

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

MonitorElementhchi2
 
MonitorElementhChi2VsEta
 
MonitorElementhChi2VsPhi
 
MonitorElementheta
 
MonitorElementhnHits
 
MonitorElementhnHitsVsEta
 
MonitorElementhnHitsVsPhi
 
MonitorElementhnLayers
 
MonitorElementhnLayersVsEta
 
MonitorElementhnLayersVsPhi
 
MonitorElementhnLooseAndAboveTracks
 
MonitorElementhnTracks
 
MonitorElementhphi
 
MonitorElementhpt
 
MonitorElementhquality
 
MonitorElementhtip
 
MonitorElementhz
 
pixelTrack::Quality minQuality_
 
edm::EDGetTokenT
< PixelTrackHeterogeneous
tokenSoATrack_
 
std::string topFolderName_
 
bool useQualityCut_
 

Detailed Description

bookLayer

Definition at line 28 of file SiPixelPhase1MonitorTrackSoA.cc.

Constructor & Destructor Documentation

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

Definition at line 64 of file SiPixelPhase1MonitorTrackSoA.cc.

References edm::ParameterSet::getParameter(), minQuality_, pixelTrack::qualityByName(), AlCaHLTBitMon_QueryRunRegistry::string, tokenSoATrack_, topFolderName_, and useQualityCut_.

64  {
65  tokenSoATrack_ = consumes<PixelTrackHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelTrackSrc"));
66  topFolderName_ = iConfig.getParameter<std::string>("topFolderName"); //"SiPixelHeterogeneous/PixelTrackSoA";
67  useQualityCut_ = iConfig.getParameter<bool>("useQualityCut");
69 }
Quality qualityByName(std::string const &name)
edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrack_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
SiPixelPhase1MonitorTrackSoA::~SiPixelPhase1MonitorTrackSoA ( )
overridedefault

Member Function Documentation

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

Definition at line 74 of file SiPixelPhase1MonitorTrackSoA.cc.

References HLT_FULL_cff::chi2, PVValHelper::eta, dqm::impl::MonitorElement::Fill(), edm::Event::getHandle(), hchi2, hChi2VsEta, hChi2VsPhi, heta, hnHits, hnHitsVsEta, hnHitsVsPhi, hnLayers, hnLayersVsEta, hnLayersVsPhi, hnLooseAndAboveTracks, hnTracks, hphi, hpt, hquality, htip, hz, HLT_FULL_cff::maxTracks, minQuality_, nHits, mkfit::Config::nLayers, BeamSpotPI::nTracks, phi, DiDispStaMuonMonitor_cfi::pt, quality, HLT_FULL_cff::tip, tokenSoATrack_, useQualityCut_, and ComparisonHelper::zip().

74  {
75  const auto& tsoaHandle = iEvent.getHandle(tokenSoATrack_);
76  if (!tsoaHandle.isValid()) {
77  edm::LogWarning("SiPixelPhase1MonitorTrackSoA") << "No Track SoA found \n returning!" << std::endl;
78  return;
79  }
80 
81  auto const& tsoa = *((tsoaHandle.product())->get());
82  auto maxTracks = tsoa.stride();
83  auto const* quality = tsoa.qualityData();
84  int32_t nTracks = 0;
85  int32_t nLooseAndAboveTracks = 0;
86 
87  for (int32_t it = 0; it < maxTracks; ++it) {
88  auto nHits = tsoa.nHits(it);
89  auto nLayers = tsoa.nLayers(it);
90  if (nHits == 0)
91  break; // this is a guard
92  float pt = tsoa.pt(it);
93  if (!(pt > 0.))
94  continue;
95 
96  // fill the quality for all tracks
97  pixelTrack::Quality qual = tsoa.quality(it);
98  hquality->Fill(int(qual));
99  nTracks++;
100 
101  if (useQualityCut_ && quality[it] < minQuality_)
102  continue;
103 
104  // fill parameters only for quality >= loose
105  float chi2 = tsoa.chi2(it);
106  float phi = tsoa.phi(it);
107  float zip = tsoa.zip(it);
108  float eta = tsoa.eta(it);
109  float tip = tsoa.tip(it);
110 
111  hchi2->Fill(chi2);
112  hChi2VsPhi->Fill(phi, chi2);
113  hChi2VsEta->Fill(eta, chi2);
114  hnHits->Fill(nHits);
116  hnHitsVsPhi->Fill(phi, nHits);
117  hnHitsVsEta->Fill(eta, nHits);
118  hnLayersVsPhi->Fill(phi, nLayers);
119  hnLayersVsEta->Fill(eta, nLayers);
120  hpt->Fill(pt);
121  heta->Fill(eta);
122  hphi->Fill(phi);
123  hz->Fill(zip);
124  htip->Fill(tip);
125  nLooseAndAboveTracks++;
126  }
127  hnTracks->Fill(nTracks);
128  hnLooseAndAboveTracks->Fill(nLooseAndAboveTracks);
129 }
string quality
void Fill(long long x)
constexpr int nLayers
Definition: Config.h:73
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
OutputIterator zip(InputIterator1 first1, InputIterator1 last1, InputIterator2 first2, InputIterator2 last2, OutputIterator result, Compare comp)
edm::EDGetTokenT< PixelTrackHeterogeneous > tokenSoATrack_
caConstants::TupleMultiplicity const CAHitNtupletGeneratorKernelsGPU::HitToTuple const cms::cuda::AtomicPairCounter GPUCACell const *__restrict__ uint32_t const *__restrict__ gpuPixelDoublets::CellNeighborsVector const gpuPixelDoublets::CellTracksVector const GPUCACell::OuterHitOfCell const int32_t nHits
Log< level::Warning, false > LogWarning
void SiPixelPhase1MonitorTrackSoA::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  iSetup 
)
override

Definition at line 134 of file SiPixelPhase1MonitorTrackSoA.cc.

References dqm::implementation::IBooker::book1D(), dqm::implementation::IBooker::bookProfile(), dqm::implementation::NavigatorBase::cd(), hchi2, hChi2VsEta, hChi2VsPhi, heta, hnHits, hnHitsVsEta, hnHitsVsPhi, hnLayers, hnLayersVsEta, hnLayersVsPhi, hnLooseAndAboveTracks, hnTracks, hphi, hpt, hquality, htip, hz, mps_fire::i, M_PI, submitPVResolutionJobs::q, pixelTrack::qualityName, dqm::impl::MonitorElement::setBinLabel(), dqm::implementation::NavigatorBase::setCurrentFolder(), AlCaHLTBitMon_QueryRunRegistry::string, topFolderName_, and parallelization::uint().

136  {
137  iBook.cd();
138  iBook.setCurrentFolder(topFolderName_);
139 
140  // clang-format off
141  std::string toRep = "Number of tracks";
142  hnTracks = iBook.book1D("nTracks", fmt::sprintf(";%s per event;#events",toRep), 1001, -0.5, 1000.5);
143  hnLooseAndAboveTracks = iBook.book1D("nLooseAndAboveTracks", fmt::sprintf(";%s (quality #geq loose) per event;#events",toRep), 1001, -0.5, 1000.5);
144 
145  toRep = "Number of all RecHits per track (quality #geq loose)";
146  hnHits = iBook.book1D("nRecHits", fmt::sprintf(";%s;#tracks",toRep), 15, -0.5, 14.5);
147  hnHitsVsPhi = iBook.bookProfile("nHitsPerTrackVsPhi", fmt::sprintf("%s vs track #eta;Track #eta;%s",toRep,toRep), 30, -3., 3., 0., 15.);
148  hnHitsVsEta = iBook.bookProfile("nHitsPerTrackVsEta", fmt::sprintf("%s vs track #phi;Track #phi;%s",toRep,toRep), 30,-M_PI,-M_PI,0.,15.);
149 
150  toRep = "Number of all layers per track (quality #geq loose)";
151  hnLayers = iBook.book1D("nLayers", fmt::sprintf(";%s;#tracks",toRep), 15, -0.5, 14.5);
152  hnLayersVsPhi = iBook.bookProfile("nLayersPerTrackVsPhi", fmt::sprintf("%s vs track #eta;Track #eta;%s",toRep,toRep), 30, -3., 3., 0., 15.);
153  hnLayersVsEta = iBook.bookProfile("nLayersPerTrackVsEta", fmt::sprintf("%s vs track #phi;Track #phi;%s",toRep,toRep), 30,-M_PI,-M_PI,0.,15.);
154 
155  toRep = "Track (quality #geq loose) #chi^{2}/ndof";
156  hchi2 = iBook.book1D("nChi2ndof", fmt::sprintf(";%s;#tracks",toRep), 40, 0., 20.);
157  hChi2VsPhi = iBook.bookProfile("nChi2ndofVsPhi", fmt::sprintf("%s vs track #eta;Track #eta;%s",toRep,toRep), 30, -3., 3., 0., 20.);
158  hChi2VsEta = iBook.bookProfile("nChi2ndofVsEta", fmt::sprintf("%s vs track #phi;Track #phi;%s",toRep,toRep), 30, -M_PI, -M_PI, 0., 20.);
159  // clang-format on
160 
161  hpt = iBook.book1D("pt", ";Track (quality #geq loose) p_{T} [GeV];#tracks", 200, 0., 200.);
162  heta = iBook.book1D("eta", ";Track (quality #geq loose) #eta;#tracks", 30, -3., 3.);
163  hphi = iBook.book1D("phi", ";Track (quality #geq loose) #phi;#tracks", 30, -M_PI, M_PI);
164  hz = iBook.book1D("z", ";Track (quality #geq loose) z [cm];#tracks", 30, -30., 30.);
165  htip = iBook.book1D("tip", ";Track (quality #geq loose) TIP [cm];#tracks", 100, -0.5, 0.5);
166  hquality = iBook.book1D("quality", ";Track Quality;#tracks", 7, -0.5, 6.5);
167  uint i = 1;
168  for (const auto& q : pixelTrack::qualityName) {
169  hquality->setBinLabel(i, q, 1);
170  i++;
171  }
172 }
const std::string qualityName[qualitySize]
virtual void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
#define M_PI
void SiPixelPhase1MonitorTrackSoA::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 174 of file SiPixelPhase1MonitorTrackSoA.cc.

References edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, HLT_FULL_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

174  {
175  // monitorpixelTrackSoA
177  desc.add<edm::InputTag>("pixelTrackSrc", edm::InputTag("pixelTracksSoA"));
178  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackSoA");
179  desc.add<bool>("useQualityCut", true);
180  desc.add<std::string>("minQuality", "loose");
181  descriptions.addWithDefaultLabel(desc);
182 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ParameterDescriptionBase * add(U const &iLabel, T const &value)

Member Data Documentation

MonitorElement* SiPixelPhase1MonitorTrackSoA::hchi2
private

Definition at line 49 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hChi2VsEta
private

Definition at line 51 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hChi2VsPhi
private

Definition at line 50 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::heta
private

Definition at line 53 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnHits
private

Definition at line 43 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnHitsVsEta
private

Definition at line 45 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnHitsVsPhi
private

Definition at line 44 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnLayers
private

Definition at line 46 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnLayersVsEta
private

Definition at line 48 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnLayersVsPhi
private

Definition at line 47 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnLooseAndAboveTracks
private

Definition at line 42 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hnTracks
private

Definition at line 41 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hphi
private

Definition at line 54 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hpt
private

Definition at line 52 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hquality
private

Definition at line 57 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::htip
private

Definition at line 56 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

MonitorElement* SiPixelPhase1MonitorTrackSoA::hz
private

Definition at line 55 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and bookHistograms().

pixelTrack::Quality SiPixelPhase1MonitorTrackSoA::minQuality_
private

Definition at line 40 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and SiPixelPhase1MonitorTrackSoA().

edm::EDGetTokenT<PixelTrackHeterogeneous> SiPixelPhase1MonitorTrackSoA::tokenSoATrack_
private

Definition at line 37 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and SiPixelPhase1MonitorTrackSoA().

std::string SiPixelPhase1MonitorTrackSoA::topFolderName_
private

Definition at line 38 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by bookHistograms(), and SiPixelPhase1MonitorTrackSoA().

bool SiPixelPhase1MonitorTrackSoA::useQualityCut_
private

Definition at line 39 of file SiPixelPhase1MonitorTrackSoA.cc.

Referenced by analyze(), and SiPixelPhase1MonitorTrackSoA().