CMS 3D CMS Logo

VertexTableProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PhysicsTools/NanoAOD
4 // Class: VertexTableProducer
5 //
13 //
14 // Original Author: Andrea Rizzi
15 // Created: Mon, 28 Aug 2017 09:26:39 GMT
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
31 
34 
36 
43 
46 
47 //
48 // class declaration
49 //
50 
52 public:
53  explicit VertexTableProducer(const edm::ParameterSet&);
54 
55  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
56 
57 private:
58  void produce(edm::Event&, const edm::EventSetup&) override;
59 
60  // ----------member data ---------------------------
61 
72  const double dlenMin_, dlenSigMin_;
73 };
74 
75 //
76 // constructors and destructor
77 //
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::VertexCompositePtrCandidate>>();
97 }
98 
99 //
100 // member functions
101 //
102 
103 // ------------ method called to produce the data ------------
104 
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::VertexCompositePtrCandidate>>();
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  selCandSv->push_back(svsProd.ptrAt(i));
186  double dx = (PV0.x() - sv.vx()), dy = (PV0.y() - sv.vy()), dz = (PV0.z() - sv.vz());
187  double pdotv = (dx * sv.px() + dy * sv.py() + dz * sv.pz()) / sv.p() / sqrt(dx * dx + dy * dy + dz * dz);
188  pAngle.push_back(std::acos(pdotv));
189  Measurement1D d2d = vdistXY.distance(
190  PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
191  dxy.push_back(d2d.value());
192  dxySig.push_back(d2d.significance());
193 
194  int sum_charge = 0;
195  for (unsigned int id = 0; id < sv.numberOfDaughters(); ++id) {
196  const reco::Candidate* daughter = sv.daughter(id);
197  sum_charge += daughter->charge();
198  }
199  charge.push_back(sum_charge);
200  }
201  }
202  i++;
203  }
204 
205  auto svsTable = std::make_unique<nanoaod::FlatTable>(selCandSv->size(), svName_, false);
206  svsTable->setDoc(svDoc_);
207  // For SV we fill from here only stuff that cannot be created with the SimpleFlatTableProducer
208  svsTable->addColumn<float>("dlen", dlen, "decay length in cm", 10);
209  svsTable->addColumn<float>("dlenSig", dlenSig, "decay length significance", 10);
210  svsTable->addColumn<float>("dxy", dxy, "2D decay length in cm", 10);
211  svsTable->addColumn<float>("dxySig", dxySig, "2D decay length significance", 10);
212  svsTable->addColumn<float>("pAngle", pAngle, "pointing angle, i.e. acos(p_SV * (SV - PV)) ", 10);
213  svsTable->addColumn<int16_t>("charge", charge, "sum of the charge of the SV tracks", 10);
214 
215  iEvent.put(std::move(pvTable), "pv");
216  iEvent.put(std::move(otherPVsTable), "otherPVs");
217  iEvent.put(std::move(svsTable), "svs");
218  iEvent.put(std::move(selCandSv));
219 }
220 
221 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
224 
225  desc.add<edm::InputTag>("pvSrc")->setComment(
226  "std::vector<reco::Vertex> and ValueMap<float> primary vertex input collections");
227  desc.add<edm::InputTag>("pfcSrc")->setComment("packedPFCandidates input collections");
228  desc.add<std::string>("goodPvCut")->setComment("selection on the primary vertex");
229  desc.add<edm::InputTag>("svSrc")->setComment(
230  "reco::VertexCompositePtrCandidate compatible secondary vertex input collection");
231  desc.add<std::string>("svCut")->setComment("selection on the secondary vertex");
232 
233  desc.add<double>("dlenMin")->setComment("minimum value of dl to select secondary vertex");
234  desc.add<double>("dlenSigMin")->setComment("minimum value of dl significance to select secondary vertex");
235 
236  desc.add<std::string>("pvName")->setComment("name of the flat table ouput");
237  desc.add<std::string>("svName")->setComment("name of the flat table ouput");
238  desc.add<std::string>("svDoc")->setComment("a few words of documentation");
239 
240  descriptions.addWithDefaultLabel(desc);
241 }
242 
243 //define this as a plug-in
reco::Vertex::Point convertPos(const GlobalPoint &p)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
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_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
std::vector< pat::PackedCandidate > PackedCandidateCollection
const edm::EDGetTokenT< edm::ValueMap< float > > pvsScore_
Definition: HeavyIon.h:7
const std::string goodPvCutString_
int iEvent
Definition: GenABIO.cc:224
goodPVs
Good Primary Vertex Selection.
VertexTableProducer(const edm::ParameterSet &)
T sqrt(T t)
Definition: SSEVec.h:19
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svs_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
virtual int charge() const =0
electric charge
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
void produce(edm::Event &, const edm::EventSetup &) override
double value() const
Definition: Measurement1D.h:25
double significance() const
Definition: Measurement1D.h:29
fixed size matrix
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