Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-08 02:18:42

0001 // -*- C++ -*-
0002 //
0003 // Package:    PhysicsTools/NanoAOD
0004 // Class:      PVertexBPHTable
0005 //
0006 /*
0007  Simple primary vertex table
0008  Description: [one line class summary]
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  gkaratha
0015 //         Created:  Mon, 28 Aug 2024 09:26:39 GMT
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0024 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0025 #include "DataFormats/Common/interface/ValueMap.h"
0026 #include "DataFormats/NanoAOD/interface/FlatTable.h"
0027 #include "DataFormats/VertexReco/interface/Vertex.h"
0028 #include "FWCore/Framework/interface/Event.h"
0029 #include "FWCore/Framework/interface/Frameworkfwd.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 #include "FWCore/Framework/interface/stream/EDProducer.h"
0032 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0033 #include "FWCore/Utilities/interface/StreamID.h"
0034 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
0035 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
0036 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0037 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0038 
0039 //
0040 // class declaration
0041 //
0042 
0043 class PVertexBPHTable : public edm::stream::EDProducer<> {
0044 public:
0045   explicit PVertexBPHTable(const edm::ParameterSet&);
0046   ~PVertexBPHTable() override;
0047 
0048 private:
0049   void beginStream(edm::StreamID) override;
0050   void produce(edm::Event&, const edm::EventSetup&) override;
0051   void endStream() override;
0052 
0053   // ----------member data ---------------------------
0054 
0055   const edm::EDGetTokenT<std::vector<reco::Vertex>> pvs_;
0056   const edm::EDGetTokenT<edm::ValueMap<float>> pvsScore_;
0057   const StringCutObjectSelector<reco::Vertex> goodPvCut_;
0058   const std::string pvName_;
0059 };
0060 
0061 //
0062 // constructors and destructor
0063 //
0064 PVertexBPHTable::PVertexBPHTable(const edm::ParameterSet& params)
0065     : pvs_(consumes<std::vector<reco::Vertex>>(params.getParameter<edm::InputTag>("pvSrc"))),
0066       pvsScore_(consumes<edm::ValueMap<float>>(params.getParameter<edm::InputTag>("pvSrc"))),
0067       goodPvCut_(params.getParameter<std::string>("goodPvCut"), true),
0068       pvName_(params.getParameter<std::string>("pvName"))
0069 
0070 {
0071   produces<nanoaod::FlatTable>("pv");
0072   produces<edm::PtrVector<reco::Candidate>>();
0073 }
0074 
0075 PVertexBPHTable::~PVertexBPHTable() {
0076   // do anything here that needs to be done at destruction time
0077   // (e.g. close files, deallocate resources etc.)
0078 }
0079 
0080 //
0081 // member functions
0082 //
0083 
0084 // ------------ method called to produce the data  ------------
0085 
0086 void PVertexBPHTable::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0087   using namespace edm;
0088   const auto& pvsScoreProd = iEvent.get(pvsScore_);
0089   auto pvsCol = iEvent.getHandle(pvs_);
0090 
0091   auto selCandPv = std::make_unique<PtrVector<reco::Candidate>>();
0092   std::vector<float> pvscore, chi2, covXX, covYY, covZZ, covXY, covXZ, covYZ, vx, vy, vz, pt, eta, phi, mass, ndof;
0093   std::vector<int> charge, ntracks;
0094 
0095   size_t i = 0;
0096   for (const auto& pv : *pvsCol) {
0097     if (!goodPvCut_(pv)) {
0098       i++;
0099       continue;
0100     }
0101     pvscore.push_back(pvsScoreProd.get(pvsCol.id(), i));
0102     ntracks.push_back(pv.tracksSize());
0103     chi2.push_back(pv.chi2());
0104     covXX.push_back(pv.covariance(0, 0));
0105     covYY.push_back(pv.covariance(1, 1));
0106     covZZ.push_back(pv.covariance(2, 2));
0107     covXY.push_back(pv.covariance(0, 1));
0108     covXZ.push_back(pv.covariance(0, 2));
0109     covYZ.push_back(pv.covariance(1, 2));
0110     vx.push_back(pv.x());
0111     vy.push_back(pv.y());
0112     vz.push_back(pv.z());
0113     pt.push_back(pv.p4().pt());
0114     eta.push_back(pv.p4().eta());
0115     phi.push_back(pv.p4().phi());
0116     mass.push_back(pv.p4().M());
0117     ndof.push_back(pv.ndof());
0118     i++;
0119   }
0120   auto table = std::make_unique<nanoaod::FlatTable>(pvscore.size(), pvName_, false, false);
0121   table->addColumn<float>("score", pvscore, "", 10);
0122   table->addColumn<float>("vx", vx, "", 10);
0123   table->addColumn<float>("vy", vy, "", 10);
0124   table->addColumn<float>("vz", vz, "", 10);
0125   table->addColumn<float>("pt", pt, "", 10);
0126   table->addColumn<float>("eta", eta, "", 10);
0127   table->addColumn<float>("phi", phi, "", 10);
0128   table->addColumn<float>("mass", mass, "", 10);
0129   table->addColumn<float>("chi2", chi2, "", 10);
0130   table->addColumn<float>("ndof", ndof, "", 10);
0131   table->addColumn<float>("covXX", covXX, "", 10);
0132   table->addColumn<float>("covYY", covYY, "", 10);
0133   table->addColumn<float>("covZZ", covZZ, "", 10);
0134   table->addColumn<float>("covXY", covXY, "", 10);
0135   table->addColumn<float>("covXZ", covXZ, "", 10);
0136   table->addColumn<float>("covYZ", covYZ, "", 10);
0137   table->addColumn<uint8_t>("ntracks", ntracks, "");
0138 
0139   iEvent.put(std::move(table), "pv");
0140 }
0141 
0142 // ------------ method called once each stream before processing any runs, lumis
0143 // or events  ------------
0144 void PVertexBPHTable::beginStream(edm::StreamID) {}
0145 
0146 // ------------ method called once each stream after processing all runs, lumis
0147 // and events  ------------
0148 void PVertexBPHTable::endStream() {}
0149 
0150 // define this as a plug-in
0151 DEFINE_FWK_MODULE(PVertexBPHTable);