CMS 3D CMS Logo

SiPixelMonitorVertexSoAAlpaka.cc
Go to the documentation of this file.
1 // -*- C++ -*-
3 // Package: SiPixelMonitorVertexSoAAlpaka
4 // Class: SiPixelMonitorVertexSoAAlpaka
5 //
8 //
9 // Author: Suvankar Roy Chowdhury
10 //
20 // DQM Histograming
26 
28 public:
29  using IndToEdm = std::vector<uint16_t>;
31  ~SiPixelMonitorVertexSoAAlpaka() override = default;
32  void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
33  void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
34  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
35 
36 private:
48 };
49 
50 //
51 // constructors
52 //
53 
55  : tokenSoAVertex_(consumes<ZVertexHost>(iConfig.getParameter<edm::InputTag>("pixelVertexSrc"))),
56  tokenBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
57  topFolderName_(iConfig.getParameter<std::string>("topFolderName")) {}
58 
59 //
60 // -- Analyze
61 //
63  const auto& vsoaHandle = iEvent.getHandle(tokenSoAVertex_);
64  if (!vsoaHandle.isValid()) {
65  edm::LogWarning("SiPixelMonitorVertexSoAAlpaka") << "No Vertex SoA found \n returning!" << std::endl;
66  return;
67  }
68 
69  auto const& vsoa = *vsoaHandle;
70  int nVertices = vsoa.view().nvFinal();
71  auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
72  float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
73  if (!bsHandle.isValid()) {
74  edm::LogWarning("SiPixelMonitorVertexSoAAlpaka") << "No beamspot found. returning vertexes with (0,0,Z) ";
75  } else {
76  const reco::BeamSpot& bs = *bsHandle;
77  x0 = bs.x0();
78  y0 = bs.y0();
79  z0 = bs.z0();
80  dxdz = bs.dxdz();
81  dydz = bs.dydz();
82  }
83 
84  for (int iv = 0; iv < nVertices; iv++) {
85  auto si = vsoa.view()[iv].sortInd();
86  auto z = vsoa.view()[si].zv();
87  auto x = x0 + dxdz * z;
88  auto y = y0 + dydz * z;
89 
90  z += z0;
91  hx->Fill(x);
92  hy->Fill(y);
93  hz->Fill(z);
94  auto ndof = vsoa.view()[si].ndof();
95  hchi2->Fill(vsoa.view()[si].chi2());
96  hchi2oNdof->Fill(vsoa.view()[si].chi2() / ndof);
97  hptv2->Fill(vsoa.view()[si].ptv2());
98  hntrks->Fill(ndof + 1);
99  }
101 }
102 
103 //
104 // -- Book Histograms
105 //
107  edm::Run const& iRun,
108  edm::EventSetup const& iSetup) {
109  //std::string top_folder = ""//
110  ibooker.cd();
112  hnVertex = ibooker.book1D("nVertex", ";# of Vertices;#entries", 101, -0.5, 100.5);
113  hx = ibooker.book1D("vx", ";Vertex x;#entries", 10, -5., 5.);
114  hy = ibooker.book1D("vy", ";Vertex y;#entries", 10, -5., 5.);
115  hz = ibooker.book1D("vz", ";Vertex z;#entries", 30, -30., 30);
116  hchi2 = ibooker.book1D("chi2", ";Vertex chi-squared;#entries", 40, 0., 20.);
117  hchi2oNdof = ibooker.book1D("chi2oNdof", ";Vertex chi-squared/Ndof;#entries", 40, 0., 20.);
118  hptv2 = ibooker.book1D("ptsq", ";Vertex #sum (p_{T})^{2};#entries", 200, 0., 200.);
119  hntrks = ibooker.book1D("ntrk", ";#tracks associated;#entries", 100, -0.5, 99.5);
120 }
121 
123  // monitorpixelVertexSoA
125  desc.add<edm::InputTag>("pixelVertexSrc", edm::InputTag("pixelVerticesAlpaka"));
126  desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
127  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexAlpaka");
128  descriptions.addWithDefaultLabel(desc);
129 }
130 
float dydz
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< reco::BeamSpot > tokenBeamSpot_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
float dxdz
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
SiPixelMonitorVertexSoAAlpaka(const edm::ParameterSet &)
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
fixed size matrix
HLT enums.
const edm::EDGetTokenT< ZVertexHost > tokenSoAVertex_
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
~SiPixelMonitorVertexSoAAlpaka() override=default
Definition: Run.h:45