File indexing completed on 2025-06-04 22:36:16
0001 #include <memory>
0002
0003
0004 #include "FWCore/Framework/interface/Frameworkfwd.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/StreamID.h"
0012
0013 #include "DataFormats/VertexReco/interface/Vertex.h"
0014 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0015
0016 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0017
0018 #include "DataFormats/NanoAOD/interface/FlatTable.h"
0019 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0020 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
0021 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
0022 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
0023 #include "DataFormats/Common/interface/ValueMap.h"
0024
0025 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0026
0027
0028
0029
0030
0031 class HLTVertexTableProducer : public edm::stream::EDProducer<> {
0032 public:
0033 explicit HLTVertexTableProducer(const edm::ParameterSet&);
0034 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035
0036 private:
0037 void produce(edm::Event&, const edm::EventSetup&) override;
0038
0039
0040 const edm::EDGetTokenT<std::vector<reco::Vertex>> pvs_;
0041 const edm::EDGetTokenT<reco::PFCandidateCollection> pfc_;
0042 const edm::EDGetTokenT<edm::ValueMap<float>> pvsScore_;
0043 const StringCutObjectSelector<reco::Vertex> goodPvCut_;
0044 const std::string goodPvCutString_;
0045 const std::string pvName_;
0046 const double dlenMin_, dlenSigMin_;
0047 };
0048
0049
0050
0051
0052
0053 HLTVertexTableProducer::HLTVertexTableProducer(const edm::ParameterSet& params)
0054 : pvs_(consumes<std::vector<reco::Vertex>>(params.getParameter<edm::InputTag>("pvSrc"))),
0055 pfc_(consumes<reco::PFCandidateCollection>(params.getParameter<edm::InputTag>("pfSrc"))),
0056 pvsScore_(consumes<edm::ValueMap<float>>(params.getParameter<edm::InputTag>("pvSrc"))),
0057 goodPvCut_(params.getParameter<std::string>("goodPvCut"), true),
0058 goodPvCutString_(params.getParameter<std::string>("goodPvCut")),
0059 pvName_(params.getParameter<std::string>("pvName")),
0060 dlenMin_(params.getParameter<double>("dlenMin")),
0061 dlenSigMin_(params.getParameter<double>("dlenSigMin"))
0062
0063 {
0064 produces<nanoaod::FlatTable>("PV");
0065 produces<edm::PtrVector<reco::VertexCompositePtrCandidate>>();
0066 }
0067
0068
0069
0070
0071
0072
0073 void HLTVertexTableProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0074 using namespace edm;
0075
0076
0077 auto pvsIn = iEvent.getHandle(pvs_);
0078 if (!pvsIn.isValid()) {
0079 edm::LogWarning("HLTVertexTableProducer")
0080 << "Invalid handle for " << pvName_ << " in primary vertex input collection";
0081 return;
0082 }
0083 const auto& pvsScoreProd = iEvent.get(pvsScore_);
0084
0085
0086 auto pfcIn = iEvent.getHandle(pfc_);
0087 if (!pfcIn.isValid()) {
0088 edm::LogWarning("HLTVertexTableProducer")
0089 << "Invalid handle for " << pvName_ << " in PF candidate input collection";
0090 return;
0091 }
0092
0093 std::vector<float> v_ndof;
0094 std::vector<float> v_chi2;
0095 std::vector<float> v_x;
0096 std::vector<float> v_y;
0097 std::vector<float> v_z;
0098 std::vector<float> v_xError;
0099 std::vector<float> v_yError;
0100 std::vector<float> v_zError;
0101 std::vector<uint8_t> v_is_good;
0102 std::vector<uint8_t> v_nTracks;
0103 std::vector<float> v_pv_score;
0104 std::vector<float> v_pv_sumpt2;
0105 std::vector<float> v_pv_sumpx;
0106 std::vector<float> v_pv_sumpy;
0107
0108 for (size_t i = 0; i < (*pvsIn).size(); i++) {
0109 v_ndof.push_back((*pvsIn)[i].ndof());
0110 v_chi2.push_back((*pvsIn)[i].normalizedChi2());
0111 v_x.push_back((*pvsIn)[i].x());
0112 v_y.push_back((*pvsIn)[i].y());
0113 v_z.push_back((*pvsIn)[i].z());
0114 v_xError.push_back((*pvsIn)[i].xError());
0115 v_yError.push_back((*pvsIn)[i].yError());
0116 v_zError.push_back((*pvsIn)[i].zError());
0117 v_nTracks.push_back((*pvsIn)[i].nTracks());
0118 v_is_good.push_back(goodPvCut_((*pvsIn)[i]));
0119 v_pv_score.push_back(pvsScoreProd.get(pvsIn.id(), i));
0120
0121 float pv_sumpt2 = 0;
0122 float pv_sumpx = 0;
0123 float pv_sumpy = 0;
0124 for (const auto& obj : *pfcIn) {
0125
0126 if (obj.charge() == 0)
0127 continue;
0128 double dz = fabs(obj.trackRef()->dz((*pvsIn)[i].position()));
0129 bool include_pfc = false;
0130 if (dz < 0.2) {
0131 include_pfc = true;
0132 for (size_t j = 0; j < (*pvsIn).size() && j != i; j++) {
0133 double newdz = fabs(obj.trackRef()->dz((*pvsIn)[j].position()));
0134 if (newdz < dz) {
0135 include_pfc = false;
0136 break;
0137 }
0138 }
0139 }
0140 if (include_pfc) {
0141 float pfc_pt = obj.pt();
0142 pv_sumpt2 += pfc_pt * pfc_pt;
0143 pv_sumpx += obj.px();
0144 pv_sumpy += obj.py();
0145 }
0146 }
0147
0148 v_pv_sumpt2.push_back(pv_sumpt2);
0149 v_pv_sumpx.push_back(pv_sumpx);
0150 v_pv_sumpy.push_back(pv_sumpy);
0151 }
0152
0153
0154 auto pvTable = std::make_unique<nanoaod::FlatTable>((*pvsIn).size(), pvName_, true);
0155 pvTable->addColumn<float>("ndof", v_ndof, "primary vertex number of degrees of freedom", 8);
0156 pvTable->addColumn<float>("chi2", v_chi2, "primary vertex reduced chi2", 8);
0157 pvTable->addColumn<float>("x", v_x, "primary vertex x coordinate", 10);
0158 pvTable->addColumn<float>("y", v_y, "primary vertex y coordinate", 10);
0159 pvTable->addColumn<float>("z", v_z, "primary vertex z coordinate", 16);
0160 pvTable->addColumn<float>("xError", v_xError, "primary vertex error in x coordinate", 10);
0161 pvTable->addColumn<float>("yError", v_yError, "primary vertex error in y coordinate", 10);
0162 pvTable->addColumn<float>("zError", v_zError, "primary vertex error in z coordinate", 16);
0163 pvTable->addColumn<uint8_t>(
0164 "isGood", v_is_good, "wheter the primary vertex passes selection: " + goodPvCutString_ + ")");
0165 pvTable->addColumn<uint8_t>("nTracks", v_nTracks, "primary vertex number of associated tracks");
0166 pvTable->addColumn<float>("score", v_pv_score, "primary vertex score, i.e. sum pt2 of clustered objects", 8);
0167 pvTable->addColumn<float>(
0168 "sumpt2", v_pv_sumpt2, "sum pt2 of pf charged candidates within dz=0.2 for the main primary vertex", 10);
0169 pvTable->addColumn<float>(
0170 "sumpx", v_pv_sumpx, "sum px of pf charged candidates within dz=0.2 for the main primary vertex", 10);
0171 pvTable->addColumn<float>(
0172 "sumpy", v_pv_sumpy, "sum py of pf charged candidates within dz=0.2 for the main primary vertex", 10);
0173
0174 iEvent.put(std::move(pvTable), "PV");
0175 }
0176
0177
0178 void HLTVertexTableProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0179 edm::ParameterSetDescription desc;
0180
0181 desc.add<std::string>("pvName")->setComment("name of the flat table ouput");
0182 desc.add<edm::InputTag>("pvSrc")->setComment(
0183 "std::vector<reco::Vertex> and ValueMap<float> primary vertex input collections");
0184 desc.add<edm::InputTag>("pfSrc")->setComment("reco::PFCandidateCollection PF candidates input collections");
0185 desc.add<std::string>("goodPvCut")->setComment("selection on the primary vertex");
0186
0187 desc.add<double>("dlenMin")->setComment("minimum value of dl to select secondary vertex");
0188 desc.add<double>("dlenSigMin")->setComment("minimum value of dl significance to select secondary vertex");
0189
0190 descriptions.addWithDefaultLabel(desc);
0191 }
0192
0193
0194 DEFINE_FWK_MODULE(HLTVertexTableProducer);