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;
115  const auto& pvsScoreProd = iEvent.get(pvsScore_);
116  auto pvsIn = iEvent.getHandle(pvs_);
117 
118  auto pvTable = std::make_unique<nanoaod::FlatTable>(1, pvName_, true);
119  pvTable->addColumnValue<float>("ndof", (*pvsIn)[0].ndof(), "main primary vertex number of degree of freedom", 8);
120  pvTable->addColumnValue<float>("x", (*pvsIn)[0].position().x(), "main primary vertex position x coordinate", 10);
121  pvTable->addColumnValue<float>("y", (*pvsIn)[0].position().y(), "main primary vertex position y coordinate", 10);
122  pvTable->addColumnValue<float>("z", (*pvsIn)[0].position().z(), "main primary vertex position z coordinate", 16);
123  pvTable->addColumnValue<float>("chi2", (*pvsIn)[0].normalizedChi2(), "main primary vertex reduced chi2", 8);
124  int goodPVs = 0;
125  for (const auto& pv : *pvsIn)
126  if (goodPvCut_(pv))
127  goodPVs++;
128  pvTable->addColumnValue<int>("npvs", pvsIn->size(), "total number of reconstructed primary vertices");
129  pvTable->addColumnValue<int>(
130  "npvsGood", goodPVs, "number of good reconstructed primary vertices. selection:" + goodPvCutString_);
131  pvTable->addColumnValue<float>(
132  "score", pvsScoreProd.get(pvsIn.id(), 0), "main primary vertex score, i.e. sum pt2 of clustered objects", 8);
133 
134  auto otherPVsTable =
135  std::make_unique<nanoaod::FlatTable>((*pvsIn).size() > 4 ? 3 : (*pvsIn).size() - 1, "Other" + pvName_, false);
136  std::vector<float> pvsz;
137  std::vector<float> pvscores;
138  for (size_t i = 1; i < (*pvsIn).size() && i < 4; i++) {
139  pvsz.push_back((*pvsIn)[i].position().z());
140  pvscores.push_back(pvsScoreProd.get(pvsIn.id(), i));
141  }
142  otherPVsTable->addColumn<float>("z", pvsz, "Z position of other primary vertices, excluding the main PV", 8);
143  otherPVsTable->addColumn<float>("score", pvscores, "scores of other primary vertices, excluding the main PV", 8);
144 
145  const auto& svsProd = iEvent.get(svs_);
146  auto selCandSv = std::make_unique<PtrVector<reco::Candidate>>();
147  std::vector<float> dlen, dlenSig, pAngle, dxy, dxySig;
148  std::vector<int> charge;
149  VertexDistance3D vdist;
150  VertexDistanceXY vdistXY;
151 
152  size_t i = 0;
153  const auto& PV0 = pvsIn->front();
154  for (const auto& sv : svsProd) {
155  if (svCut_(sv)) {
156  Measurement1D dl =
157  vdist.distance(PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
158  if (dl.value() > dlenMin_ and dl.significance() > dlenSigMin_) {
159  dlen.push_back(dl.value());
160  dlenSig.push_back(dl.significance());
161  edm::Ptr<reco::Candidate> c = svsProd.ptrAt(i);
162  selCandSv->push_back(c);
163  double dx = (PV0.x() - sv.vx()), dy = (PV0.y() - sv.vy()), dz = (PV0.z() - sv.vz());
164  double pdotv = (dx * sv.px() + dy * sv.py() + dz * sv.pz()) / sv.p() / sqrt(dx * dx + dy * dy + dz * dz);
165  pAngle.push_back(std::acos(pdotv));
166  Measurement1D d2d = vdistXY.distance(
167  PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
168  dxy.push_back(d2d.value());
169  dxySig.push_back(d2d.significance());
170 
171  int sum_charge = 0;
172  for (unsigned int id = 0; id < sv.numberOfDaughters(); ++id) {
173  const reco::Candidate* daughter = sv.daughter(id);
174  sum_charge += daughter->charge();
175  }
176  charge.push_back(sum_charge);
177  }
178  }
179  i++;
180  }
181 
182  auto svsTable = std::make_unique<nanoaod::FlatTable>(selCandSv->size(), svName_, false);
183  svsTable->setDoc(svDoc_);
184  // For SV we fill from here only stuff that cannot be created with the SimpleFlatTableProducer
185  svsTable->addColumn<float>("dlen", dlen, "decay length in cm", 10);
186  svsTable->addColumn<float>("dlenSig", dlenSig, "decay length significance", 10);
187  svsTable->addColumn<float>("dxy", dxy, "2D decay length in cm", 10);
188  svsTable->addColumn<float>("dxySig", dxySig, "2D decay length significance", 10);
189  svsTable->addColumn<float>("pAngle", pAngle, "pointing angle, i.e. acos(p_SV * (SV - PV)) ", 10);
190  svsTable->addColumn<int>("charge", charge, "sum of the charge of the SV tracks", 10);
191 
192  iEvent.put(std::move(pvTable), "pv");
193  iEvent.put(std::move(otherPVsTable), "otherPVs");
194  iEvent.put(std::move(svsTable), "svs");
195  iEvent.put(std::move(selCandSv));
196 }
197 
198 // ------------ method called once each stream before processing any runs, lumis or events ------------
200 
201 // ------------ method called once each stream after processing all runs, lumis and events ------------
203 
204 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
207 
208  desc.add<edm::InputTag>("pvSrc")->setComment(
209  "std::vector<reco::Vertex> and ValueMap<float> primary vertex input collections");
210  desc.add<std::string>("goodPvCut")->setComment("selection on the primary vertex");
211  desc.add<edm::InputTag>("svSrc")->setComment(
212  "reco::VertexCompositePtrCandidate compatible secondary vertex input collection");
213  desc.add<std::string>("svCut")->setComment("selection on the secondary vertex");
214 
215  desc.add<double>("dlenMin")->setComment("minimum value of dl to select secondary vertex");
216  desc.add<double>("dlenSigMin")->setComment("minimum value of dl significance to select secondary vertex");
217 
218  desc.add<std::string>("pvName")->setComment("name of the flat table ouput");
219  desc.add<std::string>("svName")->setComment("name of the flat table ouput");
220  desc.add<std::string>("svDoc")->setComment("a few words of documentation");
221 
222  descriptions.addWithDefaultLabel(desc);
223 }
224 
225 //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
const edm::EDGetTokenT< edm::ValueMap< float > > pvsScore_
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
void beginStream(edm::StreamID) override
const StringCutObjectSelector< reco::Vertex > goodPvCut_
const StringCutObjectSelector< reco::Candidate > svCut_
def move(src, dest)
Definition: eostools.py:511