CMS 3D CMS Logo

SiPixelPhase1MonitorVertexSoA.cc
Go to the documentation of this file.
1 // -*- C++ -*-
3 // Package: SiPixelPhase1MonitorVertexSoA
4 // Class: SiPixelPhase1MonitorVertexSoA
5 //
8 //
9 // Author: Suvankar Roy Chowdhury
10 //
20 // DQM Histograming
26 
28 public:
29  using IndToEdm = std::vector<uint16_t>;
31  ~SiPixelPhase1MonitorVertexSoA() 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<ZVertexHeterogeneous>(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 //
61 // -- Analyze
62 //
64  const auto& vsoaHandle = iEvent.getHandle(tokenSoAVertex_);
65  if (!vsoaHandle.isValid()) {
66  edm::LogWarning("SiPixelPhase1MonitorTrackSoA") << "No Vertex SoA found \n returning!" << std::endl;
67  return;
68  }
69 
70  auto const& vsoa = *((vsoaHandle.product())->get());
71  int nVertices = vsoa.nvFinal;
72  auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
73  float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
74  if (!bsHandle.isValid()) {
75  edm::LogWarning("PixelVertexProducer") << "No beamspot found. returning vertexes with (0,0,Z) ";
76  } else {
77  const reco::BeamSpot& bs = *bsHandle;
78  x0 = bs.x0();
79  y0 = bs.y0();
80  z0 = bs.z0();
81  dxdz = bs.dxdz();
82  dydz = bs.dydz();
83  }
84  for (int iv = 0; iv < nVertices; iv++) {
85  auto si = vsoa.sortInd[iv];
86  auto z = vsoa.zv[si];
87  auto x = x0 + dxdz * z;
88  auto y = y0 + dydz * z;
89  z += z0;
90  hx->Fill(x);
91  hy->Fill(y);
92  hz->Fill(z);
93  auto ndof = vsoa.ndof[si];
94  hchi2->Fill(vsoa.chi2[si]);
95  hchi2oNdof->Fill(vsoa.chi2[si] / ndof);
96  hptv2->Fill(vsoa.ptv2[si]);
97  hntrks->Fill(ndof + 1);
98  }
100 }
101 
102 //
103 // -- Book Histograms
104 //
106  edm::Run const& iRun,
107  edm::EventSetup const& iSetup) {
108  //std::string top_folder = ""//
109  ibooker.cd();
111  hnVertex = ibooker.book1D("nVertex", ";# of Vertex;#entries", 101, -0.5, 100.5);
112  hx = ibooker.book1D("vx", ";Vertex x;#entries", 10, -5., 5.);
113  hy = ibooker.book1D("vy", ";Vertex y;#entries", 10, -5., 5.);
114  hz = ibooker.book1D("vz", ";Vertex z;#entries", 30, -30., 30);
115  hchi2 = ibooker.book1D("chi2", ";Vertex chi-squared;#entries", 40, 0., 20.);
116  hchi2oNdof = ibooker.book1D("chi2oNdof", ";Vertex chi-squared/Ndof;#entries", 40, 0., 20.);
117  hptv2 = ibooker.book1D("ptsq", ";Vertex p_T squared;#entries", 200, 0., 200.);
118  hntrks = ibooker.book1D("ntrk", ";#tracks associated;#entries", 100, -0.5, 99.5);
119 }
120 
122  // monitorpixelVertexSoA
124  desc.add<edm::InputTag>("pixelVertexSrc", edm::InputTag("pixelVerticesSoA"));
125  desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
126  desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexSoA");
127  descriptions.addWithDefaultLabel(desc);
128 }
float dydz
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int32_t *__restrict__ iv
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
float dxdz
edm::EDGetTokenT< reco::BeamSpot > tokenBeamSpot_
~SiPixelPhase1MonitorVertexSoA() override=default
void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:224
SiPixelPhase1MonitorVertexSoA(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &iRun, edm::EventSetup const &iSetup) override
edm::EDGetTokenT< ZVertexHeterogeneous > 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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: Run.h:45