Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-06-07 22:22:31

0001 // -*- C++ -*-
0002 // Package:    SiPixelPhase1CompareVertexSoA
0003 // Class:      SiPixelPhase1CompareVertexSoA
0004 //
0005 /**\class SiPixelPhase1CompareVertexSoA SiPixelPhase1CompareVertexSoA.cc 
0006 */
0007 //
0008 // Author: Suvankar Roy Chowdhury
0009 //
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 #include "DataFormats/Common/interface/Handle.h"
0017 // DQM Histograming
0018 #include "DQMServices/Core/interface/MonitorElement.h"
0019 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0020 #include "DQMServices/Core/interface/DQMStore.h"
0021 #include "CUDADataFormats/Vertex/interface/ZVertexHeterogeneous.h"
0022 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0023 
0024 class SiPixelPhase1CompareVertexSoA : public DQMEDAnalyzer {
0025 public:
0026   using IndToEdm = std::vector<uint16_t>;
0027   explicit SiPixelPhase1CompareVertexSoA(const edm::ParameterSet&);
0028   ~SiPixelPhase1CompareVertexSoA() override = default;
0029   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0030   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0031   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0032 
0033 private:
0034   const edm::EDGetTokenT<ZVertexHeterogeneous> tokenSoAVertexCPU_;
0035   const edm::EDGetTokenT<ZVertexHeterogeneous> tokenSoAVertexGPU_;
0036   const edm::EDGetTokenT<reco::BeamSpot> tokenBeamSpot_;
0037   const std::string topFolderName_;
0038   const float dzCut_;
0039   MonitorElement* hnVertex_;
0040   MonitorElement* hx_;
0041   MonitorElement* hy_;
0042   MonitorElement* hz_;
0043   MonitorElement* hchi2_;
0044   MonitorElement* hchi2oNdof_;
0045   MonitorElement* hptv2_;
0046   MonitorElement* hntrks_;
0047   MonitorElement* hxdiff_;
0048   MonitorElement* hydiff_;
0049   MonitorElement* hzdiff_;
0050 };
0051 
0052 //
0053 // constructors
0054 //
0055 
0056 SiPixelPhase1CompareVertexSoA::SiPixelPhase1CompareVertexSoA(const edm::ParameterSet& iConfig)
0057     : tokenSoAVertexCPU_(consumes<ZVertexHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelVertexSrcCPU"))),
0058       tokenSoAVertexGPU_(consumes<ZVertexHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelVertexSrcGPU"))),
0059       tokenBeamSpot_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
0060       topFolderName_(iConfig.getParameter<std::string>("topFolderName")),
0061       dzCut_(iConfig.getParameter<double>("dzCut")) {}
0062 
0063 //
0064 // -- Analyze
0065 //
0066 void SiPixelPhase1CompareVertexSoA::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0067   const auto& vsoaHandleCPU = iEvent.getHandle(tokenSoAVertexCPU_);
0068   const auto& vsoaHandleGPU = iEvent.getHandle(tokenSoAVertexGPU_);
0069   if (not vsoaHandleCPU or not vsoaHandleGPU) {
0070     edm::LogWarning out("SiPixelPhase1CompareTrackSoA");
0071     if (not vsoaHandleCPU) {
0072       out << "reference (cpu) tracks not found; ";
0073     }
0074     if (not vsoaHandleGPU) {
0075       out << "target (gpu) tracks not found; ";
0076     }
0077     out << "the comparison will not run.";
0078     return;
0079   }
0080 
0081   auto const& vsoaCPU = *vsoaHandleCPU->get();
0082   int nVerticesCPU = vsoaCPU.nvFinal;
0083   auto const& vsoaGPU = *vsoaHandleGPU->get();
0084   int nVerticesGPU = vsoaGPU.nvFinal;
0085 
0086   auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
0087   float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
0088   if (!bsHandle.isValid()) {
0089     edm::LogWarning("PixelVertexProducer") << "No beamspot found. returning vertexes with (0,0,Z) ";
0090   } else {
0091     const reco::BeamSpot& bs = *bsHandle;
0092     x0 = bs.x0();
0093     y0 = bs.y0();
0094     z0 = bs.z0();
0095     dxdz = bs.dxdz();
0096     dydz = bs.dydz();
0097   }
0098 
0099   for (int ivc = 0; ivc < nVerticesCPU; ivc++) {
0100     auto sic = vsoaCPU.sortInd[ivc];
0101     auto zc = vsoaCPU.zv[sic];
0102     auto xc = x0 + dxdz * zc;
0103     auto yc = y0 + dydz * zc;
0104     zc += z0;
0105 
0106     auto ndofCPU = vsoaCPU.ndof[sic];
0107     auto chi2CPU = vsoaCPU.chi2[sic];
0108 
0109     const int32_t notFound = -1;
0110     int32_t closestVtxidx = notFound;
0111     float mindz = dzCut_;
0112 
0113     for (int ivg = 0; ivg < nVerticesGPU; ivg++) {
0114       auto sig = vsoaGPU.sortInd[ivg];
0115       auto zgc = vsoaGPU.zv[sig] + z0;
0116       auto zDist = std::abs(zc - zgc);
0117       //insert some matching condition
0118       if (zDist > dzCut_)
0119         continue;
0120       if (mindz > zDist) {
0121         mindz = zDist;
0122         closestVtxidx = sig;
0123       }
0124     }
0125     if (closestVtxidx == notFound)
0126       continue;
0127 
0128     auto zg = vsoaGPU.zv[closestVtxidx];
0129     auto xg = x0 + dxdz * zg;
0130     auto yg = y0 + dydz * zg;
0131     zg += z0;
0132     auto ndofGPU = vsoaGPU.ndof[closestVtxidx];
0133     auto chi2GPU = vsoaGPU.chi2[closestVtxidx];
0134 
0135     hx_->Fill(xc - x0, xg - x0);
0136     hy_->Fill(yc - y0, yg - y0);
0137     hz_->Fill(zc, zg);
0138     hxdiff_->Fill(xc - xg);
0139     hydiff_->Fill(yc - yg);
0140     hzdiff_->Fill(zc - zg);
0141     hchi2_->Fill(chi2CPU, chi2GPU);
0142     hchi2oNdof_->Fill(chi2CPU / ndofCPU, chi2GPU / ndofGPU);
0143     hptv2_->Fill(vsoaCPU.ptv2[sic], vsoaGPU.ptv2[closestVtxidx]);
0144     hntrks_->Fill(ndofCPU + 1, ndofGPU + 1);
0145   }
0146   hnVertex_->Fill(nVerticesCPU, nVerticesGPU);
0147 }
0148 
0149 //
0150 // -- Book Histograms
0151 //
0152 void SiPixelPhase1CompareVertexSoA::bookHistograms(DQMStore::IBooker& ibooker,
0153                                                    edm::Run const& iRun,
0154                                                    edm::EventSetup const& iSetup) {
0155   ibooker.cd();
0156   ibooker.setCurrentFolder(topFolderName_);
0157 
0158   // FIXME: all the 2D correlation plots are quite heavy in terms of memory consumption, so a as soon as DQM supports either TH2I or THnSparse
0159   // these should be moved to a less resource consuming format
0160   hnVertex_ = ibooker.book2I("nVertex", "# of Vertex;CPU;GPU", 101, -0.5, 100.5, 101, -0.5, 100.5);
0161   hx_ = ibooker.book2I("vx", "Vertez x - Beamspot x;CPU;GPU", 50, -0.1, 0.1, 50, -0.1, 0.1);
0162   hy_ = ibooker.book2I("vy", "Vertez y - Beamspot y;CPU;GPU", 50, -0.1, 0.1, 50, -0.1, 0.1);
0163   hz_ = ibooker.book2I("vz", "Vertez z;CPU;GPU", 30, -30., 30., 30, -30., 30.);
0164   hchi2_ = ibooker.book2I("chi2", "Vertex chi-squared;CPU;GPU", 40, 0., 20., 40, 0., 20.);
0165   hchi2oNdof_ = ibooker.book2I("chi2oNdof", "Vertex chi-squared/Ndof;CPU;GPU", 40, 0., 20., 40, 0., 20.);
0166   hptv2_ = ibooker.book2I("ptsq", "Vertex p_T squared;CPU;GPU", 200, 0., 200., 200, 0., 200.);
0167   hntrks_ = ibooker.book2I("ntrk", "#tracks associated;CPU;GPU", 100, -0.5, 99.5, 100, -0.5, 99.5);
0168   hntrks_ = ibooker.book2I("ntrk", "#tracks associated;CPU;GPU", 100, -0.5, 99.5, 100, -0.5, 99.5);
0169   hxdiff_ = ibooker.book1D("vxdiff", ";Vertex x difference (CPU - GPU);#entries", 100, -0.001, 0.001);
0170   hydiff_ = ibooker.book1D("vydiff", ";Vertex y difference (CPU - GPU);#entries", 100, -0.001, 0.001);
0171   hzdiff_ = ibooker.book1D("vzdiff", ";Vertex z difference (CPU - GPU);#entries", 100, -2.5, 2.5);
0172 }
0173 
0174 void SiPixelPhase1CompareVertexSoA::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0175   // monitorpixelVertexSoA
0176   edm::ParameterSetDescription desc;
0177   desc.add<edm::InputTag>("pixelVertexSrcCPU", edm::InputTag("pixelVerticesSoA@cpu"));
0178   desc.add<edm::InputTag>("pixelVertexSrcGPU", edm::InputTag("pixelVerticesSoA@cuda"));
0179   desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
0180   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexCompareSoAGPUvsCPU");
0181   desc.add<double>("dzCut", 1.);
0182   descriptions.addWithDefaultLabel(desc);
0183 }
0184 DEFINE_FWK_MODULE(SiPixelPhase1CompareVertexSoA);