-
Notifications
You must be signed in to change notification settings - Fork 4.6k
Add PV sumPT2 of pf charged particles associated to PV in NanoAOD #43487
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Changes from 4 commits
eee28de
0a01e66
e7601d0
ba318e9
1c4c038
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
| Original file line number | Diff line number | Diff line change |
|---|---|---|
|
|
@@ -41,6 +41,9 @@ | |
| #include "RecoVertex/VertexPrimitives/interface/VertexState.h" | ||
| #include "DataFormats/Common/interface/ValueMap.h" | ||
|
|
||
| #include "DataFormats/PatCandidates/interface/PackedCandidate.h" | ||
| #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h" | ||
|
|
||
| // | ||
| // class declaration | ||
| // | ||
|
|
@@ -65,6 +68,7 @@ class VertexTableProducer : public edm::stream::EDProducer<> { | |
| // ----------member data --------------------------- | ||
|
|
||
| const edm::EDGetTokenT<std::vector<reco::Vertex>> pvs_; | ||
| const edm::EDGetTokenT<pat::PackedCandidateCollection> pfc_; | ||
| const edm::EDGetTokenT<edm::ValueMap<float>> pvsScore_; | ||
| const edm::EDGetTokenT<edm::View<reco::VertexCompositePtrCandidate>> svs_; | ||
| const StringCutObjectSelector<reco::Candidate> svCut_; | ||
|
|
@@ -81,6 +85,7 @@ class VertexTableProducer : public edm::stream::EDProducer<> { | |
| // | ||
| VertexTableProducer::VertexTableProducer(const edm::ParameterSet& params) | ||
| : pvs_(consumes<std::vector<reco::Vertex>>(params.getParameter<edm::InputTag>("pvSrc"))), | ||
| pfc_(consumes<pat::PackedCandidateCollection>(params.getParameter<edm::InputTag>("pfcSrc"))), | ||
| pvsScore_(consumes<edm::ValueMap<float>>(params.getParameter<edm::InputTag>("pvSrc"))), | ||
| svs_(consumes<edm::View<reco::VertexCompositePtrCandidate>>(params.getParameter<edm::InputTag>("svSrc"))), | ||
| svCut_(params.getParameter<std::string>("svCut"), true), | ||
|
|
@@ -115,6 +120,9 @@ void VertexTableProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSe | |
| const auto& pvsScoreProd = iEvent.get(pvsScore_); | ||
| auto pvsIn = iEvent.getHandle(pvs_); | ||
|
|
||
| //pf candidates collection | ||
| auto pfcIn = iEvent.getHandle(pfc_); | ||
|
|
||
| auto pvTable = std::make_unique<nanoaod::FlatTable>(1, pvName_, true); | ||
| pvTable->addColumnValue<float>("ndof", (*pvsIn)[0].ndof(), "main primary vertex number of degree of freedom", 8); | ||
| pvTable->addColumnValue<float>("x", (*pvsIn)[0].position().x(), "main primary vertex position x coordinate", 10); | ||
|
|
@@ -131,6 +139,28 @@ void VertexTableProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSe | |
| pvTable->addColumnValue<float>( | ||
| "score", pvsScoreProd.get(pvsIn.id(), 0), "main primary vertex score, i.e. sum pt2 of clustered objects", 8); | ||
|
|
||
| float pv_sumpt2 = 0.0; | ||
| for (size_t i = 0; i < (*pfcIn).size(); i++) { | ||
| if ((*pfcIn)[i].charge() == 0) { | ||
| continue; | ||
| } // skip neutrals | ||
| double dz = fabs((*pfcIn)[i].dz((*pvsIn)[0].position())); | ||
| int include_pfc = 0; | ||
|
||
| if (dz < 0.2) | ||
| include_pfc = 1; | ||
| for (size_t j = 1; j < (*pvsIn).size(); j++) { | ||
|
||
| double newdz = fabs((*pfcIn)[i].dz((*pvsIn)[j].position())); | ||
| if (newdz < dz) | ||
| include_pfc = 0; | ||
|
||
| } // this pf candidate belongs to other PV | ||
| if (include_pfc == 1) { | ||
| double pfc_pt = (*pfcIn)[i].pt(); | ||
| pv_sumpt2 += pfc_pt * pfc_pt; | ||
| } | ||
| } | ||
| pvTable->addColumnValue<float>( | ||
| "sumpt2", pv_sumpt2, "sum pt2 of tracks for the main primary vertex", 10); //precision from 8 to 10 | ||
|
||
|
|
||
| auto otherPVsTable = | ||
| std::make_unique<nanoaod::FlatTable>((*pvsIn).size() > 4 ? 3 : (*pvsIn).size() - 1, "Other" + pvName_, false); | ||
| std::vector<float> pvsz; | ||
|
|
@@ -207,6 +237,7 @@ void VertexTableProducer::fillDescriptions(edm::ConfigurationDescriptions& descr | |
|
|
||
| desc.add<edm::InputTag>("pvSrc")->setComment( | ||
| "std::vector<reco::Vertex> and ValueMap<float> primary vertex input collections"); | ||
| desc.add<edm::InputTag>("pfcSrc")->setComment("packedPFCandidates input collections"); | ||
| desc.add<std::string>("goodPvCut")->setComment("selection on the primary vertex"); | ||
| desc.add<edm::InputTag>("svSrc")->setComment( | ||
| "reco::VertexCompositePtrCandidate compatible secondary vertex input collection"); | ||
|
|
||
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
can you please use pfcand pointers (using auto for example) instead of index loops