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>(
121  "ndof", (*pvsIn)[0].ndof(), "main primary vertex number of degree of freedom", nanoaod::FlatTable::FloatColumn, 8);
122  pvTable->addColumnValue<float>(
123  "x", (*pvsIn)[0].position().x(), "main primary vertex position x coordinate", nanoaod::FlatTable::FloatColumn, 10);
124  pvTable->addColumnValue<float>(
125  "y", (*pvsIn)[0].position().y(), "main primary vertex position y coordinate", nanoaod::FlatTable::FloatColumn, 10);
126  pvTable->addColumnValue<float>(
127  "z", (*pvsIn)[0].position().z(), "main primary vertex position z coordinate", nanoaod::FlatTable::FloatColumn, 16);
128  pvTable->addColumnValue<float>(
129  "chi2", (*pvsIn)[0].normalizedChi2(), "main primary vertex reduced chi2", nanoaod::FlatTable::FloatColumn, 8);
130  int goodPVs = 0;
131  for (const auto& pv : *pvsIn)
132  if (goodPvCut_(pv))
133  goodPVs++;
134  pvTable->addColumnValue<int>(
135  "npvs", (*pvsIn).size(), "total number of reconstructed primary vertices", nanoaod::FlatTable::IntColumn);
136  pvTable->addColumnValue<int>("npvsGood",
137  goodPVs,
138  "number of good reconstructed primary vertices. selection:" + goodPvCutString_,
140  pvTable->addColumnValue<float>("score",
141  (*pvsScoreIn).get(pvsIn.id(), 0),
142  "main primary vertex score, i.e. sum pt2 of clustered objects",
144  8);
145 
146  auto otherPVsTable =
147  std::make_unique<nanoaod::FlatTable>((*pvsIn).size() > 4 ? 3 : (*pvsIn).size() - 1, "Other" + pvName_, false);
148  std::vector<float> pvsz;
149  for (size_t i = 1; i < (*pvsIn).size() && i < 4; i++)
150  pvsz.push_back((*pvsIn)[i - 1].position().z());
151  otherPVsTable->addColumn<float>(
152  "z", pvsz, "Z position of other primary vertices, excluding the main PV", nanoaod::FlatTable::FloatColumn, 8);
153 
155  iEvent.getByToken(svs_, svsIn);
156  auto selCandSv = std::make_unique<PtrVector<reco::Candidate>>();
157  std::vector<float> dlen, dlenSig, pAngle, dxy, dxySig;
158  VertexDistance3D vdist;
159  VertexDistanceXY vdistXY;
160 
161  size_t i = 0;
162  const auto& PV0 = pvsIn->front();
163  for (const auto& sv : *svsIn) {
164  if (svCut_(sv)) {
165  Measurement1D dl =
166  vdist.distance(PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
167  if (dl.value() > dlenMin_ and dl.significance() > dlenSigMin_) {
168  dlen.push_back(dl.value());
169  dlenSig.push_back(dl.significance());
170  edm::Ptr<reco::Candidate> c = svsIn->ptrAt(i);
171  selCandSv->push_back(c);
172  double dx = (PV0.x() - sv.vx()), dy = (PV0.y() - sv.vy()), dz = (PV0.z() - sv.vz());
173  double pdotv = (dx * sv.px() + dy * sv.py() + dz * sv.pz()) / sv.p() / sqrt(dx * dx + dy * dy + dz * dz);
174  pAngle.push_back(std::acos(pdotv));
175  Measurement1D d2d = vdistXY.distance(
176  PV0, VertexState(RecoVertex::convertPos(sv.position()), RecoVertex::convertError(sv.error())));
177  dxy.push_back(d2d.value());
178  dxySig.push_back(d2d.significance());
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", nanoaod::FlatTable::FloatColumn, 10);
187  svsTable->addColumn<float>("dlenSig", dlenSig, "decay length significance", nanoaod::FlatTable::FloatColumn, 10);
188  svsTable->addColumn<float>("dxy", dxy, "2D decay length in cm", nanoaod::FlatTable::FloatColumn, 10);
189  svsTable->addColumn<float>("dxySig", dxySig, "2D decay length significance", nanoaod::FlatTable::FloatColumn, 10);
190  svsTable->addColumn<float>(
191  "pAngle", pAngle, "pointing angle, i.e. acos(p_SV * (SV - PV)) ", nanoaod::FlatTable::FloatColumn, 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
nanoaod::FlatTable::FloatColumn
Definition: FlatTable.h:39
edm::StreamID
Definition: StreamID.h:30
VertexTableProducer::goodPvCut_
const StringCutObjectSelector< reco::Vertex > goodPvCut_
Definition: VertexTableProducer.cc:71
Measurement1D
Definition: Measurement1D.h:11
mps_fire.i
i
Definition: mps_fire.py:355
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
sistrip::View
View
Definition: ConstantsForView.h:26
CalibrationSummaryClient_cfi.params
params
Definition: CalibrationSummaryClient_cfi.py:14
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
Measurement1D::value
double value() const
Definition: Measurement1D.h:25
VertexTableProducer::VertexTableProducer
VertexTableProducer(const edm::ParameterSet &)
Definition: VertexTableProducer.cc:82
VertexDistance3D::distance
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
Definition: VertexDistance3D.cc:17
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
VertexTableProducer::pvs_
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvs_
Definition: VertexTableProducer.cc:67
EDProducer.h
VertexTableProducer::svCut_
const StringCutObjectSelector< reco::Candidate > svCut_
Definition: VertexTableProducer.cc:70
ConvertToFromReco.h
VertexDistance3D.h
VertexTableProducer::svs_
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svs_
Definition: VertexTableProducer.cc:69
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle
Definition: AssociativeIterator.h:50
VertexTableProducer::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: VertexTableProducer.cc:206
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
VertexDistanceXY.h
VertexTableProducer::beginStream
void beginStream(edm::StreamID) override
Definition: VertexTableProducer.cc:200
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
pfDeepBoostedJetPreprocessParams_cfi.sv
sv
Definition: pfDeepBoostedJetPreprocessParams_cfi.py:226
VertexCompositePtrCandidate.h
VertexState.h
VertexDistance3D
Definition: VertexDistance3D.h:13
VertexTableProducer::dlenSigMin_
const double dlenSigMin_
Definition: VertexTableProducer.cc:76
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
Vertex.h
VertexTableProducer::~VertexTableProducer
~VertexTableProducer() override
Definition: VertexTableProducer.cc:102
VertexTableProducer::svName_
const std::string svName_
Definition: VertexTableProducer.cc:74
funct::true
true
Definition: Factorize.h:173
Measurement1D::significance
double significance() const
Definition: Measurement1D.h:29
HLT_2018_cff.InputTag
InputTag
Definition: HLT_2018_cff.py:79016
edm::ParameterSet
Definition: ParameterSet.h:36
RecoVertex::convertError
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
Event.h
VertexTableProducer::pvName_
const std::string pvName_
Definition: VertexTableProducer.cc:73
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::ParameterSetDescription::setUnknown
void setUnknown()
Definition: ParameterSetDescription.cc:39
PVValHelper::dy
Definition: PVValidationHelpers.h:49
edm::EventSetup
Definition: EventSetup.h:57
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
VertexDistanceXY
Definition: VertexDistanceXY.h:11
VertexTableProducer::dlenMin_
const double dlenMin_
Definition: VertexTableProducer.cc:76
VertexTableProducer::endStream
void endStream() override
Definition: VertexTableProducer.cc:203
FlatTable.h
edm::Ptr< reco::Candidate >
ValueMap.h
VertexTableProducer::goodPvCutString_
const std::string goodPvCutString_
Definition: VertexTableProducer.cc:72
nanoaod::FlatTable::IntColumn
Definition: FlatTable.h:40
AlcaRecoSelection_cff.goodPVs
goodPVs
Good Primary Vertex Selection.
Definition: AlcaRecoSelection_cff.py:54
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
StringCutObjectSelector.h
HltBtagValidation_cff.Vertex
Vertex
Definition: HltBtagValidation_cff.py:32
VertexDistanceXY::distance
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
Definition: VertexDistanceXY.cc:19
PVValHelper::dxy
Definition: PVValidationHelpers.h:47
PVValHelper::dz
Definition: PVValidationHelpers.h:50
Frameworkfwd.h
VertexState
Definition: VertexState.h:13
StringCutObjectSelector< reco::Candidate >
RecoVertex::convertPos
reco::Vertex::Point convertPos(const GlobalPoint &p)
Definition: ConvertToFromReco.h:7
ParameterSet.h
VertexTableProducer::pvsScore_
const edm::EDGetTokenT< edm::ValueMap< float > > pvsScore_
Definition: VertexTableProducer.cc:68
edm::Event
Definition: Event.h:73
edm::ConfigurationDescriptions::addDefault
void addDefault(ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:99
StreamID.h
VertexTableProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: VertexTableProducer.cc:113
PVValHelper::dx
Definition: PVValidationHelpers.h:48
VertexTableProducer::svDoc_
const std::string svDoc_
Definition: VertexTableProducer.cc:75
VertexTableProducer
Definition: VertexTableProducer.cc:48