Skip to content
Merged
Show file tree
Hide file tree
Changes from 4 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
31 changes: 31 additions & 0 deletions PhysicsTools/NanoAOD/plugins/VertexTableProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -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
//
Expand All @@ -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_;
Expand All @@ -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),
Expand Down Expand Up @@ -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);
Expand All @@ -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++) {
Copy link
Contributor

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

if ((*pfcIn)[i].charge() == 0) {
continue;
} // skip neutrals
double dz = fabs((*pfcIn)[i].dz((*pvsIn)[0].position()));
int include_pfc = 0;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

use boolean type

if (dz < 0.2)
include_pfc = 1;
for (size_t j = 1; j < (*pvsIn).size(); j++) {
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

you do not want to loop over the over pv if include_pfc is already false

double newdz = fabs((*pfcIn)[i].dz((*pvsIn)[j].position()));
if (newdz < dz)
include_pfc = 0;
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

you probably need to add a "break" statement here, to avoid running over all vertices

} // 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
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

should you be more precise and say "pf candidate" instead of tracks ?


auto otherPVsTable =
std::make_unique<nanoaod::FlatTable>((*pvsIn).size() > 4 ? 3 : (*pvsIn).size() - 1, "Other" + pvName_, false);
std::vector<float> pvsz;
Expand Down Expand Up @@ -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");
Expand Down
1 change: 1 addition & 0 deletions PhysicsTools/NanoAOD/python/nanoDQM_cfi.py
Original file line number Diff line number Diff line change
Expand Up @@ -622,6 +622,7 @@
Plot1D('npvs', 'npvs', 20, 0, 60, 'total number of reconstructed primary vertices'),
Plot1D('npvsGood', 'npvsGood', 20, 0, 60, 'total number of Good primary vertices'),
Plot1D('score', 'score', 20, 0, 300000, 'main primary vertex score, i.e. sum pt2 of clustered objects'),
Plot1D('sumpt2', 'sumpt2', 100, 0, 300000, 'main primary vertex sum pt2 of the charged pf candidates'),
Plot1D('x', 'x', 20, -0.3, 0.3, 'main primary vertex position x coordinate'),
Plot1D('y', 'y', 20, -0.3, 0.3, 'main primary vertex position y coordinate'),
Plot1D('z', 'z', 20, -20, 20, 'main primary vertex position z coordinate'),
Expand Down
1 change: 1 addition & 0 deletions PhysicsTools/NanoAOD/python/vertices_cff.py
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,7 @@
vertexTable = cms.EDProducer("VertexTableProducer",
pvSrc = cms.InputTag("offlineSlimmedPrimaryVertices"),
goodPvCut = cms.string("!isFake && ndof > 4 && abs(z) <= 24 && position.Rho <= 2"),
pfcSrc = cms.InputTag("packedPFCandidates"),
svSrc = cms.InputTag("linkedObjects", "vertices"),
svCut = cms.string(""), # careful: adding a cut here would make the collection matching inconsistent with the SV table
dlenMin = cms.double(0),
Expand Down