CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PATSecondaryVertexSlimmer.cc
Go to the documentation of this file.
1 #include <string>
2 #include <memory>
3 
4 // user include files
8 
12 
19 
20 namespace pat {
22  public:
24  ~PATSecondaryVertexSlimmer() override;
25 
26  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
27 
28  private:
34  };
35 } // namespace pat
36 
38  : src_(consumes<reco::VertexCompositePtrCandidateCollection>(iConfig.getParameter<edm::InputTag>("src"))),
39  srcLegacy_(mayConsume<std::vector<reco::Vertex> >(iConfig.getParameter<edm::InputTag>("src"))),
40  srcV0s_(mayConsume<reco::VertexCompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("src"))),
41  map_(consumes<edm::Association<pat::PackedCandidateCollection> >(
42  iConfig.getParameter<edm::InputTag>("packedPFCandidates"))),
43  map2_(mayConsume<edm::Association<pat::PackedCandidateCollection> >(
44  iConfig.existsAs<edm::InputTag>("lostTracksCandidates")
45  ? iConfig.getParameter<edm::InputTag>("lostTracksCandidates")
46  : edm::InputTag("lostTracks"))) {
47  produces<reco::VertexCompositePtrCandidateCollection>();
48 }
49 
51 
53  auto outPtr = std::make_unique<reco::VertexCompositePtrCandidateCollection>();
54 
56  iEvent.getByToken(src_, candVertices);
57 
59  iEvent.getByToken(map_, pf2pc);
60 
61  // if reco::VertexCompositePtrCandidate secondary vertices are present
62  if (candVertices.isValid()) {
63  outPtr->reserve(candVertices->size());
64  for (unsigned int i = 0, n = candVertices->size(); i < n; ++i) {
65  reco::VertexCompositePtrCandidate v = (*candVertices)[i];
66 
67  std::vector<reco::CandidatePtr> daughters = v.daughterPtrVector();
68  v.clearDaughters();
69 
70  for (std::vector<reco::CandidatePtr>::const_iterator it = daughters.begin(); it != daughters.end(); ++it) {
71  if ((*pf2pc)[*it].isNonnull() && (*pf2pc)[*it]->numberOfHits() > 0)
72  v.addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc)[*it])));
73  }
74 
75  outPtr->push_back(v);
76  }
77  }
78  // otherwise fallback to reco::Vertex secondary vertices
79  else {
81  iEvent.getByToken(srcLegacy_, vertices);
82 
83  if (vertices.isValid()) {
85  iEvent.getByToken(map2_, pf2pc2);
86 
87  outPtr->reserve(vertices->size());
88  for (unsigned int i = 0, n = vertices->size(); i < n; ++i) {
89  const reco::Vertex& v = (*vertices)[i];
90  outPtr->push_back(reco::VertexCompositePtrCandidate(0, v.p4(), v.position(), v.error(), v.chi2(), v.ndof()));
91 
92  for (reco::Vertex::trackRef_iterator it = v.tracks_begin(); it != v.tracks_end(); it++) {
93  if (v.trackWeight(*it) > 0.5) {
94  if ((*pf2pc)[*it].isNonnull() && (*pf2pc)[*it]->numberOfHits() > 0) {
95  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc)[*it])));
96  } else {
97  if ((*pf2pc2)[*it].isNonnull()) {
98  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc2)[*it])));
99  } else {
100  edm::LogError("PATSecondaryVertexSlimmer") << "HELPME" << std::endl;
101  }
102  }
103  }
104  }
105  }
106  } else {
109  iEvent.getByToken(srcV0s_, srcV0s);
110 
112  iEvent.getByToken(map2_, pf2pc2);
113 
114  outPtr->reserve(srcV0s->size());
115  for (unsigned int i = 0, n = srcV0s->size(); i < n; ++i) {
116  const reco::VertexCompositeCandidate& v = (*srcV0s)[i];
117  outPtr->push_back(reco::VertexCompositePtrCandidate(
118  0, v.p4(), v.vertex(), v.vertexCovariance(), v.vertexChi2(), v.vertexNdof()));
119 
120  for (size_t dIdx = 0; dIdx < v.numberOfDaughters(); dIdx++) {
121  reco::TrackRef trackRef = (dynamic_cast<const reco::RecoChargedCandidate*>(v.daughter(dIdx)))->track();
122  if ((*pf2pc)[trackRef].isNonnull() && (*pf2pc)[trackRef]->numberOfHits() > 0) {
123  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc)[trackRef])));
124  } else {
125  if ((*pf2pc2)[trackRef].isNonnull()) {
126  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc2)[trackRef])));
127  } else {
128  edm::LogError("PATSecondaryVertexSlimmer") << "HELPME" << std::endl;
129  }
130  }
131  }
132  }
133 
134  } // if reco::Vertex
135  } // if Candidate
136 
137  iEvent.put(std::move(outPtr));
138 }
139 
const edm::EDGetTokenT< std::vector< reco::Vertex > > srcLegacy_
double vertexChi2() const override
chi-squares
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< VertexCompositeCandidate > VertexCompositeCandidateCollection
collection of Candidate objects
Ptr< typename C::value_type > refToPtr(Ref< C, typename C::value_type, refhelper::FindUsingAdvance< C, typename C::value_type > > const &ref)
Definition: RefToPtr.h:18
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > map2_
std::vector< pat::PackedCandidate > PackedCandidateCollection
const edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > map_
const Point & vertex() const override
vertex position (overwritten by PF...)
Log< level::Error, false > LogError
PATSecondaryVertexSlimmer(const edm::ParameterSet &)
const Point & position() const
position
Definition: Vertex.h:127
const LorentzVector & p4() const final
four-momentum Lorentz vector
std::vector< VertexCompositePtrCandidate > VertexCompositePtrCandidateCollection
collection of Candidate objects
virtual const daughters & daughterPtrVector() const
references to daughtes
const edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > srcV0s_
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
double chi2() const
chi-squares
Definition: Vertex.h:116
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.h:110
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:96
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.h:108
bool isValid() const
Definition: HandleBase.h:70
virtual void clearDaughters()
clear daughter references
double ndof() const
Definition: Vertex.h:123
math::XYZTLorentzVectorD p4(float mass=0.13957018, float minWeight=0.5) const
Returns the four momentum of the sum of the tracks, assuming the given mass for the decay products...
Definition: Vertex.cc:103
size_type numberOfDaughters() const override
number of daughters
Error error() const
return SMatrix
Definition: Vertex.h:163
void addDaughter(const CandidatePtr &)
add a daughter via a reference
const edm::EDGetTokenT< reco::VertexCompositePtrCandidateCollection > src_
double vertexCovariance(int i, int j) const override
(i, j)-th element of error matrix, i, j = 0, ... 2
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector&lt;TrackRef&gt;
Definition: Vertex.h:38