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 
44 //
45 // class declaration
46 //
47 
49 public:
50  explicit VertexTableProducer(const edm::ParameterSet&);
51  ~VertexTableProducer() override;
52 
53  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
54 
55 private:
56  void beginStream(edm::StreamID) override;
57  void produce(edm::Event&, const edm::EventSetup&) override;
58  void endStream() override;
59 
60  //virtual void beginRun(edm::Run const&, edm::EventSetup const&) override;
61  //virtual void endRun(edm::Run const&, edm::EventSetup const&) override;
62  //virtual void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
63  //virtual void endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&) override;
64 
65  // ----------member data ---------------------------
66 
76  const double dlenMin_, dlenSigMin_;
77 };
78 
79 //
80 // constructors and destructor
81 //
83  : pvs_(consumes<std::vector<reco::Vertex>>(params.getParameter<edm::InputTag>("pvSrc"))),
84  pvsScore_(consumes<edm::ValueMap<float>>(params.getParameter<edm::InputTag>("pvSrc"))),
85  svs_(consumes<edm::View<reco::VertexCompositePtrCandidate>>(params.getParameter<edm::InputTag>("svSrc"))),
86  svCut_(params.getParameter<std::string>("svCut"), true),
87  goodPvCut_(params.getParameter<std::string>("goodPvCut"), true),
88  goodPvCutString_(params.getParameter<std::string>("goodPvCut")),
89  pvName_(params.getParameter<std::string>("pvName")),
90  svName_(params.getParameter<std::string>("svName")),
91  svDoc_(params.getParameter<std::string>("svDoc")),
92  dlenMin_(params.getParameter<double>("dlenMin")),
93  dlenSigMin_(params.getParameter<double>("dlenSigMin"))
94 
95 {
96  produces<nanoaod::FlatTable>("pv");
97  produces<nanoaod::FlatTable>("otherPVs");
98  produces<nanoaod::FlatTable>("svs");
99  produces<edm::PtrVector<reco::Candidate>>();
100 }
101 
103  // do anything here that needs to be done at destruction time
104  // (e.g. close files, deallocate resources etc.)
105 }
106 
107 //
108 // member functions
109 //
110 
111 // ------------ method called to produce the data ------------
112 
114  using namespace edm;
117  iEvent.getByToken(pvs_, pvsIn);
118  iEvent.getByToken(pvsScore_, pvsScoreIn);
119  auto pvTable = std::make_unique<nanoaod::FlatTable>(1, pvName_, true);
120  pvTable->addColumnValue<float>("ndof", (*pvsIn)[0].ndof(), "main primary vertex number of degree of freedom", 8);
121  pvTable->addColumnValue<float>("x", (*pvsIn)[0].position().x(), "main primary vertex position x coordinate", 10);
122  pvTable->addColumnValue<float>("y", (*pvsIn)[0].position().y(), "main primary vertex position y coordinate", 10);
123  pvTable->addColumnValue<float>("z", (*pvsIn)[0].position().z(), "main primary vertex position z coordinate", 16);
124  pvTable->addColumnValue<float>("chi2", (*pvsIn)[0].normalizedChi2(), "main primary vertex reduced chi2", 8);
125  int goodPVs = 0;
126  for (const auto& pv : *pvsIn)
127  if (goodPvCut_(pv))
128  goodPVs++;
129  pvTable->addColumnValue<int>("npvs", pvsIn->size(), "total number of reconstructed primary vertices");
130  pvTable->addColumnValue<int>(
131  "npvsGood", goodPVs, "number of good reconstructed primary vertices. selection:" + goodPvCutString_);
132  pvTable->addColumnValue<float>(
133  "score", pvsScoreIn->get(pvsIn.id(), 0), "main primary vertex score, i.e. sum pt2 of clustered objects", 8);
134 
135  auto otherPVsTable =
136  std::make_unique<nanoaod::FlatTable>((*pvsIn).size() > 4 ? 3 : (*pvsIn).size() - 1, "Other" + pvName_, false);
137  std::vector<float> pvsz;
138  std::vector<float> pvscores;
139  for (size_t i = 1; i < (*pvsIn).size() && i < 4; i++) {
140  pvsz.push_back((*pvsIn)[i].position().z());
141  pvscores.push_back((*pvsScoreIn).get(pvsIn.id(), i));
142  }
143  otherPVsTable->addColumn<float>("z", pvsz, "Z position of other primary vertices, excluding the main PV", 8);
144  otherPVsTable->addColumn<float>("score", pvscores, "scores of other primary vertices, excluding the main PV", 8);
145 
147  iEvent.getByToken(svs_, svsIn);
148  auto selCandSv = std::make_unique<PtrVector<reco::Candidate>>();
149  std::vector<float> dlen, dlenSig, pAngle, dxy, dxySig;
150  std::vector<int> charge;
151  VertexDistance3D vdist;
152  VertexDistanceXY vdistXY;
153 
154  size_t i = 0;
155  const auto& PV0 = pvsIn->front();
156  for (const auto& sv : *svsIn) {
157  if (svCut_(sv)) {
158  Measurement1D dl =
159  vdist.distance(PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
160  if (dl.value() > dlenMin_ and dl.significance() > dlenSigMin_) {
161  dlen.push_back(dl.value());
162  dlenSig.push_back(dl.significance());
163  edm::Ptr<reco::Candidate> c = svsIn->ptrAt(i);
164  selCandSv->push_back(c);
165  double dx = (PV0.x() - sv.vx()), dy = (PV0.y() - sv.vy()), dz = (PV0.z() - sv.vz());
166  double pdotv = (dx * sv.px() + dy * sv.py() + dz * sv.pz()) / sv.p() / sqrt(dx * dx + dy * dy + dz * dz);
167  pAngle.push_back(std::acos(pdotv));
168  Measurement1D d2d = vdistXY.distance(
169  PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
170  dxy.push_back(d2d.value());
171  dxySig.push_back(d2d.significance());
172 
173  int sum_charge = 0;
174  for (unsigned int id = 0; id < sv.numberOfDaughters(); ++id) {
175  const reco::Candidate* daughter = sv.daughter(id);
176  sum_charge += daughter->charge();
177  }
178  charge.push_back(sum_charge);
179  }
180  }
181  i++;
182  }
183 
184  auto svsTable = std::make_unique<nanoaod::FlatTable>(selCandSv->size(), svName_, false);
185  // For SV we fill from here only stuff that cannot be created with the SimpleFlatTableProducer
186  svsTable->addColumn<float>("dlen", dlen, "decay length in cm", 10);
187  svsTable->addColumn<float>("dlenSig", dlenSig, "decay length significance", 10);
188  svsTable->addColumn<float>("dxy", dxy, "2D decay length in cm", 10);
189  svsTable->addColumn<float>("dxySig", dxySig, "2D decay length significance", 10);
190  svsTable->addColumn<float>("pAngle", pAngle, "pointing angle, i.e. acos(p_SV * (SV - PV)) ", 10);
191  svsTable->addColumn<int>("charge", charge, "sum of the charge of the SV tracks", 10);
192 
193  iEvent.put(std::move(pvTable), "pv");
194  iEvent.put(std::move(otherPVsTable), "otherPVs");
195  iEvent.put(std::move(svsTable), "svs");
196  iEvent.put(std::move(selCandSv));
197 }
198 
199 // ------------ method called once each stream before processing any runs, lumis or events ------------
201 
202 // ------------ method called once each stream after processing all runs, lumis and events ------------
204 
205 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
207  //The following says we do not know what parameters are allowed so do no validation
208  // Please change this to state exactly what you do use, even if it is no parameters
210  desc.setUnknown();
211  descriptions.addDefault(desc);
212 }
213 
214 //define this as a plug-in
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_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
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.
void addDefault(ParameterSetDescription const &psetDescription)
VertexTableProducer(const edm::ParameterSet &)
T sqrt(T t)
Definition: SSEVec.h:19
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svs_
def pv(vc)
Definition: MetAnalyzer.py:7
#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
void beginStream(edm::StreamID) override
const StringCutObjectSelector< reco::Vertex > goodPvCut_
const StringCutObjectSelector< reco::Candidate > svCut_
def move(src, dest)
Definition: eostools.py:511