CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
VertexTableProducer Class Reference

#include <PhysicsTools/VertexTableProducer/plugins/VertexTableProducer.cc>

Inheritance diagram for VertexTableProducer:
edm::stream::EDProducer<>

Public Member Functions

 VertexTableProducer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

const double dlenMin_
 
const double dlenSigMin_
 
const StringCutObjectSelector< reco::VertexgoodPvCut_
 
const std::string goodPvCutString_
 
const edm::EDGetTokenT< pat::PackedCandidateCollectionpfc_
 
const std::string pvName_
 
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvs_
 
const edm::EDGetTokenT< edm::ValueMap< float > > pvsScore_
 
const StringCutObjectSelector< reco::CandidatesvCut_
 
const std::string svDoc_
 
const std::string svName_
 
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svs_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 51 of file VertexTableProducer.cc.

Constructor & Destructor Documentation

◆ VertexTableProducer()

VertexTableProducer::VertexTableProducer ( const edm::ParameterSet params)
explicit

Definition at line 78 of file VertexTableProducer.cc.

79  : pvs_(consumes<std::vector<reco::Vertex>>(params.getParameter<edm::InputTag>("pvSrc"))),
80  pfc_(consumes<pat::PackedCandidateCollection>(params.getParameter<edm::InputTag>("pfcSrc"))),
81  pvsScore_(consumes<edm::ValueMap<float>>(params.getParameter<edm::InputTag>("pvSrc"))),
82  svs_(consumes<edm::View<reco::VertexCompositePtrCandidate>>(params.getParameter<edm::InputTag>("svSrc"))),
83  svCut_(params.getParameter<std::string>("svCut"), true),
84  goodPvCut_(params.getParameter<std::string>("goodPvCut"), true),
85  goodPvCutString_(params.getParameter<std::string>("goodPvCut")),
86  pvName_(params.getParameter<std::string>("pvName")),
87  svName_(params.getParameter<std::string>("svName")),
88  svDoc_(params.getParameter<std::string>("svDoc")),
89  dlenMin_(params.getParameter<double>("dlenMin")),
90  dlenSigMin_(params.getParameter<double>("dlenSigMin"))
91 
92 {
93  produces<nanoaod::FlatTable>("pv");
94  produces<nanoaod::FlatTable>("otherPVs");
95  produces<nanoaod::FlatTable>("svs");
96  produces<edm::PtrVector<reco::Candidate>>();
97 }
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvs_
const std::string pvName_
const std::string svDoc_
const edm::EDGetTokenT< edm::ValueMap< float > > pvsScore_
const std::string goodPvCutString_
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svs_
const std::string svName_
const StringCutObjectSelector< reco::Vertex > goodPvCut_
const edm::EDGetTokenT< pat::PackedCandidateCollection > pfc_
const StringCutObjectSelector< reco::Candidate > svCut_

Member Function Documentation

◆ fillDescriptions()

void VertexTableProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 223 of file VertexTableProducer.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and AlCaHLTBitMon_QueryRunRegistry::string.

223  {
225 
226  desc.add<edm::InputTag>("pvSrc")->setComment(
227  "std::vector<reco::Vertex> and ValueMap<float> primary vertex input collections");
228  desc.add<edm::InputTag>("pfcSrc")->setComment("packedPFCandidates input collections");
229  desc.add<std::string>("goodPvCut")->setComment("selection on the primary vertex");
230  desc.add<edm::InputTag>("svSrc")->setComment(
231  "reco::VertexCompositePtrCandidate compatible secondary vertex input collection");
232  desc.add<std::string>("svCut")->setComment("selection on the secondary vertex");
233 
234  desc.add<double>("dlenMin")->setComment("minimum value of dl to select secondary vertex");
235  desc.add<double>("dlenSigMin")->setComment("minimum value of dl significance to select secondary vertex");
236 
237  desc.add<std::string>("pvName")->setComment("name of the flat table ouput");
238  desc.add<std::string>("svName")->setComment("name of the flat table ouput");
239  desc.add<std::string>("svDoc")->setComment("a few words of documentation");
240 
241  descriptions.addWithDefaultLabel(desc);
242 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ produce()

void VertexTableProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 105 of file VertexTableProducer.cc.

References HltBtagPostValidation_cff::c, ALCARECOTkAlJpsiMuMu_cff::charge, reco::Candidate::charge(), RecoVertex::convertError(), RecoVertex::convertPos(), VertexDistanceXY::distance(), VertexDistance3D::distance(), dlenMin_, dlenSigMin_, PVValHelper::dx, PVValHelper::dxy, PVValHelper::dy, PVValHelper::dz, goodPvCut_, goodPvCutString_, AlcaRecoSelection_cff::goodPVs, mps_fire::i, l1ctLayer2EG_cff::id, iEvent, dqmiolumiharvest::j, eostools::move(), getGTfromDQMFile::obj, pfc_, position, pvName_, pvs_, pvsScore_, Measurement1D::significance(), mathSSE::sqrt(), pfDeepBoostedJetPreprocessParams_cfi::sv, svCut_, svDoc_, svName_, svs_, Measurement1D::value(), and z.

105  {
106  using namespace edm;
107  const auto& pvsScoreProd = iEvent.get(pvsScore_);
108  auto pvsIn = iEvent.getHandle(pvs_);
109 
110  //pf candidates collection
111  auto pfcIn = iEvent.getHandle(pfc_);
112 
113  auto pvTable = std::make_unique<nanoaod::FlatTable>(1, pvName_, true);
114  pvTable->addColumnValue<float>("ndof", (*pvsIn)[0].ndof(), "main primary vertex number of degree of freedom", 8);
115  pvTable->addColumnValue<float>("x", (*pvsIn)[0].position().x(), "main primary vertex position x coordinate", 10);
116  pvTable->addColumnValue<float>("y", (*pvsIn)[0].position().y(), "main primary vertex position y coordinate", 10);
117  pvTable->addColumnValue<float>("z", (*pvsIn)[0].position().z(), "main primary vertex position z coordinate", 16);
118  pvTable->addColumnValue<float>("chi2", (*pvsIn)[0].normalizedChi2(), "main primary vertex reduced chi2", 8);
119  int goodPVs = 0;
120  for (const auto& pv : *pvsIn)
121  if (goodPvCut_(pv))
122  goodPVs++;
123  pvTable->addColumnValue<uint8_t>("npvs", pvsIn->size(), "total number of reconstructed primary vertices");
124  pvTable->addColumnValue<uint8_t>(
125  "npvsGood", goodPVs, "number of good reconstructed primary vertices. selection:" + goodPvCutString_);
126  pvTable->addColumnValue<float>(
127  "score", pvsScoreProd.get(pvsIn.id(), 0), "main primary vertex score, i.e. sum pt2 of clustered objects", 8);
128 
129  float pv_sumpt2 = 0.0, pv_sumpx = 0.0, pv_sumpy = 0.0;
130  for (const auto& obj : *pfcIn) {
131  if (obj.charge() == 0) {
132  continue;
133  } // skip neutrals
134  double dz = fabs(obj.dz((*pvsIn)[0].position()));
135  bool include_pfc = false;
136  if (dz < 0.2) {
137  include_pfc = true;
138  for (size_t j = 1; j < (*pvsIn).size(); j++) {
139  double newdz = fabs(obj.dz((*pvsIn)[j].position()));
140  if (newdz < dz) {
141  include_pfc = false;
142  break;
143  }
144  } // this pf candidate belongs to other PV
145  }
146  if (include_pfc) {
147  float pfc_pt = obj.pt();
148  pv_sumpt2 += pfc_pt * pfc_pt;
149  pv_sumpx += obj.px();
150  pv_sumpy += obj.py();
151  }
152  }
153  pvTable->addColumnValue<float>(
154  "sumpt2", pv_sumpt2, "sum pt2 of pf charged candidates for the main primary vertex", 10);
155  pvTable->addColumnValue<float>("sumpx", pv_sumpx, "sum px of pf charged candidates for the main primary vertex", 10);
156  pvTable->addColumnValue<float>("sumpy", pv_sumpy, "sum py of pf charged candidates for the main primary vertex", 10);
157 
158  auto otherPVsTable =
159  std::make_unique<nanoaod::FlatTable>((*pvsIn).size() > 4 ? 3 : (*pvsIn).size() - 1, "Other" + pvName_, false);
160  std::vector<float> pvsz;
161  std::vector<float> pvscores;
162  for (size_t i = 1; i < (*pvsIn).size() && i < 4; i++) {
163  pvsz.push_back((*pvsIn)[i].position().z());
164  pvscores.push_back(pvsScoreProd.get(pvsIn.id(), i));
165  }
166  otherPVsTable->addColumn<float>("z", pvsz, "Z position of other primary vertices, excluding the main PV", 8);
167  otherPVsTable->addColumn<float>("score", pvscores, "scores of other primary vertices, excluding the main PV", 8);
168 
169  const auto& svsProd = iEvent.get(svs_);
170  auto selCandSv = std::make_unique<PtrVector<reco::Candidate>>();
171  std::vector<float> dlen, dlenSig, pAngle, dxy, dxySig;
172  std::vector<int16_t> charge;
173  VertexDistance3D vdist;
174  VertexDistanceXY vdistXY;
175 
176  size_t i = 0;
177  const auto& PV0 = pvsIn->front();
178  for (const auto& sv : svsProd) {
179  if (svCut_(sv)) {
180  Measurement1D dl =
181  vdist.distance(PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
182  if (dl.value() > dlenMin_ and dl.significance() > dlenSigMin_) {
183  dlen.push_back(dl.value());
184  dlenSig.push_back(dl.significance());
185  edm::Ptr<reco::Candidate> c = svsProd.ptrAt(i);
186  selCandSv->push_back(c);
187  double dx = (PV0.x() - sv.vx()), dy = (PV0.y() - sv.vy()), dz = (PV0.z() - sv.vz());
188  double pdotv = (dx * sv.px() + dy * sv.py() + dz * sv.pz()) / sv.p() / sqrt(dx * dx + dy * dy + dz * dz);
189  pAngle.push_back(std::acos(pdotv));
190  Measurement1D d2d = vdistXY.distance(
191  PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
192  dxy.push_back(d2d.value());
193  dxySig.push_back(d2d.significance());
194 
195  int sum_charge = 0;
196  for (unsigned int id = 0; id < sv.numberOfDaughters(); ++id) {
197  const reco::Candidate* daughter = sv.daughter(id);
198  sum_charge += daughter->charge();
199  }
200  charge.push_back(sum_charge);
201  }
202  }
203  i++;
204  }
205 
206  auto svsTable = std::make_unique<nanoaod::FlatTable>(selCandSv->size(), svName_, false);
207  svsTable->setDoc(svDoc_);
208  // For SV we fill from here only stuff that cannot be created with the SimpleFlatTableProducer
209  svsTable->addColumn<float>("dlen", dlen, "decay length in cm", 10);
210  svsTable->addColumn<float>("dlenSig", dlenSig, "decay length significance", 10);
211  svsTable->addColumn<float>("dxy", dxy, "2D decay length in cm", 10);
212  svsTable->addColumn<float>("dxySig", dxySig, "2D decay length significance", 10);
213  svsTable->addColumn<float>("pAngle", pAngle, "pointing angle, i.e. acos(p_SV * (SV - PV)) ", 10);
214  svsTable->addColumn<int16_t>("charge", charge, "sum of the charge of the SV tracks", 10);
215 
216  iEvent.put(std::move(pvTable), "pv");
217  iEvent.put(std::move(otherPVsTable), "otherPVs");
218  iEvent.put(std::move(svsTable), "svs");
219  iEvent.put(std::move(selCandSv));
220 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvs_
const std::string pvName_
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
const std::string svDoc_
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
const edm::EDGetTokenT< edm::ValueMap< float > > pvsScore_
const std::string goodPvCutString_
int iEvent
Definition: GenABIO.cc:224
goodPVs
Good Primary Vertex Selection.
T sqrt(T t)
Definition: SSEVec.h:19
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svs_
virtual int charge() const =0
electric charge
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
double value() const
Definition: Measurement1D.h:25
double significance() const
Definition: Measurement1D.h:29
const std::string svName_
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:289
const StringCutObjectSelector< reco::Vertex > goodPvCut_
const edm::EDGetTokenT< pat::PackedCandidateCollection > pfc_
const StringCutObjectSelector< reco::Candidate > svCut_
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ dlenMin_

const double VertexTableProducer::dlenMin_
private

Definition at line 72 of file VertexTableProducer.cc.

Referenced by produce().

◆ dlenSigMin_

const double VertexTableProducer::dlenSigMin_
private

Definition at line 72 of file VertexTableProducer.cc.

Referenced by produce().

◆ goodPvCut_

const StringCutObjectSelector<reco::Vertex> VertexTableProducer::goodPvCut_
private

Definition at line 67 of file VertexTableProducer.cc.

Referenced by produce().

◆ goodPvCutString_

const std::string VertexTableProducer::goodPvCutString_
private

Definition at line 68 of file VertexTableProducer.cc.

Referenced by produce().

◆ pfc_

const edm::EDGetTokenT<pat::PackedCandidateCollection> VertexTableProducer::pfc_
private

Definition at line 63 of file VertexTableProducer.cc.

Referenced by produce().

◆ pvName_

const std::string VertexTableProducer::pvName_
private

Definition at line 69 of file VertexTableProducer.cc.

Referenced by produce().

◆ pvs_

const edm::EDGetTokenT<std::vector<reco::Vertex> > VertexTableProducer::pvs_
private

Definition at line 62 of file VertexTableProducer.cc.

Referenced by produce().

◆ pvsScore_

const edm::EDGetTokenT<edm::ValueMap<float> > VertexTableProducer::pvsScore_
private

Definition at line 64 of file VertexTableProducer.cc.

Referenced by produce().

◆ svCut_

const StringCutObjectSelector<reco::Candidate> VertexTableProducer::svCut_
private

Definition at line 66 of file VertexTableProducer.cc.

Referenced by produce().

◆ svDoc_

const std::string VertexTableProducer::svDoc_
private

Definition at line 71 of file VertexTableProducer.cc.

Referenced by produce().

◆ svName_

const std::string VertexTableProducer::svName_
private

Definition at line 70 of file VertexTableProducer.cc.

Referenced by produce().

◆ svs_

const edm::EDGetTokenT<edm::View<reco::VertexCompositePtrCandidate> > VertexTableProducer::svs_
private

Definition at line 65 of file VertexTableProducer.cc.

Referenced by produce().