Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:16

0001 // -*- C++ -*-
0002 ///bookLayer
0003 // Package:    SiPixelMonitorVertexSoAAlpaka
0004 // Class:      SiPixelMonitorVertexSoAAlpaka
0005 //
0006 /**\class SiPixelMonitorVertexSoAAlpaka SiPixelMonitorVertexSoAAlpaka.cc
0007 */
0008 //
0009 // Author: Suvankar Roy Chowdhury
0010 //
0011 #include "FWCore/Framework/interface/Frameworkfwd.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018 #include "FWCore/Utilities/interface/InputTag.h"
0019 #include "DataFormats/Common/interface/Handle.h"
0020 // DQM Histograming
0021 #include "DQMServices/Core/interface/MonitorElement.h"
0022 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0023 #include "DQMServices/Core/interface/DQMStore.h"
0024 #include "DataFormats/VertexSoA/interface/ZVertexHost.h"
0025 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0026 
0027 class SiPixelMonitorVertexSoAAlpaka : public DQMEDAnalyzer {
0028 public:
0029   using IndToEdm = std::vector<uint16_t>;
0030   explicit SiPixelMonitorVertexSoAAlpaka(const edm::ParameterSet&);
0031   ~SiPixelMonitorVertexSoAAlpaka() override = default;
0032   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0033   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0034   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035 
0036 private:
0037   const edm::EDGetTokenT<ZVertexHost> tokenSoAVertex_;
0038   const edm::EDGetTokenT<reco::BeamSpot> tokenBeamSpot_;
0039   std::string topFolderName_;
0040   MonitorElement* hnVertex;
0041   MonitorElement* hx;
0042   MonitorElement* hy;
0043   MonitorElement* hz;
0044   MonitorElement* hchi2;
0045   MonitorElement* hchi2oNdof;
0046   MonitorElement* hptv2;
0047   MonitorElement* hntrks;
0048 };
0049 
0050 //
0051 // constructors
0052 //
0053 
0054 SiPixelMonitorVertexSoAAlpaka::SiPixelMonitorVertexSoAAlpaka(const edm::ParameterSet& iConfig)
0055     : tokenSoAVertex_(consumes<ZVertexHost>(iConfig.getParameter<edm::InputTag>("pixelVertexSrc"))),
0056       tokenBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
0057       topFolderName_(iConfig.getParameter<std::string>("topFolderName")) {}
0058 
0059 //
0060 // -- Analyze
0061 //
0062 void SiPixelMonitorVertexSoAAlpaka::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0063   const auto& vsoaHandle = iEvent.getHandle(tokenSoAVertex_);
0064   if (!vsoaHandle.isValid()) {
0065     edm::LogWarning("SiPixelMonitorVertexSoAAlpaka") << "No Vertex SoA found \n returning!" << std::endl;
0066     return;
0067   }
0068 
0069   auto const& vsoa = *vsoaHandle;
0070   int nVertices = vsoa.view().nvFinal();
0071   auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
0072   float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
0073   if (!bsHandle.isValid()) {
0074     edm::LogWarning("SiPixelMonitorVertexSoAAlpaka") << "No beamspot found. returning vertexes with (0,0,Z) ";
0075   } else {
0076     const reco::BeamSpot& bs = *bsHandle;
0077     x0 = bs.x0();
0078     y0 = bs.y0();
0079     z0 = bs.z0();
0080     dxdz = bs.dxdz();
0081     dydz = bs.dydz();
0082   }
0083 
0084   for (int iv = 0; iv < nVertices; iv++) {
0085     auto si = vsoa.view()[iv].sortInd();
0086     auto z = vsoa.view()[si].zv();
0087     auto x = x0 + dxdz * z;
0088     auto y = y0 + dydz * z;
0089 
0090     z += z0;
0091     hx->Fill(x);
0092     hy->Fill(y);
0093     hz->Fill(z);
0094     auto ndof = vsoa.view()[si].ndof();
0095     hchi2->Fill(vsoa.view()[si].chi2());
0096     hchi2oNdof->Fill(vsoa.view()[si].chi2() / ndof);
0097     hptv2->Fill(vsoa.view()[si].ptv2());
0098     hntrks->Fill(ndof + 1);
0099   }
0100   hnVertex->Fill(nVertices);
0101 }
0102 
0103 //
0104 // -- Book Histograms
0105 //
0106 void SiPixelMonitorVertexSoAAlpaka::bookHistograms(DQMStore::IBooker& ibooker,
0107                                                    edm::Run const& iRun,
0108                                                    edm::EventSetup const& iSetup) {
0109   //std::string top_folder = ""//
0110   ibooker.cd();
0111   ibooker.setCurrentFolder(topFolderName_);
0112   hnVertex = ibooker.book1D("nVertex", ";# of Vertices;#entries", 101, -0.5, 100.5);
0113   hx = ibooker.book1D("vx", ";Vertex x;#entries", 10, -5., 5.);
0114   hy = ibooker.book1D("vy", ";Vertex y;#entries", 10, -5., 5.);
0115   hz = ibooker.book1D("vz", ";Vertex z;#entries", 30, -30., 30);
0116   hchi2 = ibooker.book1D("chi2", ";Vertex chi-squared;#entries", 40, 0., 20.);
0117   hchi2oNdof = ibooker.book1D("chi2oNdof", ";Vertex chi-squared/Ndof;#entries", 40, 0., 20.);
0118   hptv2 = ibooker.book1D("ptsq", ";Vertex #sum (p_{T})^{2};#entries", 200, 0., 200.);
0119   hntrks = ibooker.book1D("ntrk", ";#tracks associated;#entries", 100, -0.5, 99.5);
0120 }
0121 
0122 void SiPixelMonitorVertexSoAAlpaka::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0123   // monitorpixelVertexSoA
0124   edm::ParameterSetDescription desc;
0125   desc.add<edm::InputTag>("pixelVertexSrc", edm::InputTag("pixelVerticesAlpaka"));
0126   desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
0127   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexAlpaka");
0128   descriptions.addWithDefaultLabel(desc);
0129 }
0130 
0131 DEFINE_FWK_MODULE(SiPixelMonitorVertexSoAAlpaka);