CMS 3D CMS Logo

PATSecondaryVertexSlimmer.cc
Go to the documentation of this file.
1 #include <string>
2 #include <memory>
3 
4 // user include files
8 
12 
19 
20 
21 namespace pat {
23  public:
25  ~PATSecondaryVertexSlimmer() override;
26 
27  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
28  private:
34  };
35 }
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> >(iConfig.getParameter<edm::InputTag>("packedPFCandidates"))),
42  map2_(mayConsume<edm::Association<pat::PackedCandidateCollection> >(iConfig.existsAs<edm::InputTag>("lostTracksCandidates") ? iConfig.getParameter<edm::InputTag>("lostTracksCandidates") : edm::InputTag("lostTracks") ))
43 {
44  produces< reco::VertexCompositePtrCandidateCollection >();
45 }
46 
48 
50 
51  auto outPtr = std::make_unique<reco::VertexCompositePtrCandidateCollection>();
52 
54  iEvent.getByToken(src_, candVertices);
55 
57  iEvent.getByToken(map_,pf2pc);
58 
59  // if reco::VertexCompositePtrCandidate secondary vertices are present
60  if( candVertices.isValid() )
61  {
62  outPtr->reserve(candVertices->size());
63  for (unsigned int i = 0, n = candVertices->size(); i < n; ++i) {
64 
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 
72  if((*pf2pc)[*it].isNonnull() && (*pf2pc)[*it]->numberOfHits() > 0)
73  v.addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc)[*it]) ));
74  }
75 
76  outPtr->push_back( v );
77  }
78  }
79  // otherwise fallback to reco::Vertex secondary vertices
80  else
81  {
83  iEvent.getByToken(srcLegacy_, vertices);
84 
85  if( vertices.isValid() ){
86 
88  iEvent.getByToken(map2_,pf2pc2);
89 
90  outPtr->reserve(vertices->size());
91  for (unsigned int i = 0, n = vertices->size(); i < n; ++i) {
92  const reco::Vertex &v = (*vertices)[i];
93  outPtr->push_back(reco::VertexCompositePtrCandidate(0,v.p4(),v.position(), v.error(), v.chi2(), v.ndof()));
94 
95  for(reco::Vertex::trackRef_iterator it=v.tracks_begin(); it != v.tracks_end(); it++) {
96  if(v.trackWeight(*it)>0.5) {
97  if((*pf2pc)[*it].isNonnull() && (*pf2pc)[*it]->numberOfHits() > 0) {
98  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc)[*it]) ));
99  }
100  else {
101  if((*pf2pc2)[*it].isNonnull()) {
102  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc2)[*it]) ));
103  }
104  else { edm::LogError("PATSecondaryVertexSlimmer") << "HELPME" << std::endl;}
105  }
106  }
107  }
108  }
109  }else {
112  iEvent.getByToken( srcV0s_, srcV0s );
113 
115  iEvent.getByToken(map2_,pf2pc2);
116 
117  outPtr->reserve(srcV0s->size());
118  for (unsigned int i = 0, n = srcV0s->size(); i < n; ++i) {
119  const reco::VertexCompositeCandidate &v = (*srcV0s)[i];
120  outPtr->push_back(reco::VertexCompositePtrCandidate(0,v.p4(),v.vertex(), v.vertexCovariance(), v.vertexChi2(), v.vertexNdof()));
121 
122  for(size_t dIdx=0; dIdx< v.numberOfDaughters() ; dIdx++){
123  reco::TrackRef trackRef = (dynamic_cast<const reco::RecoChargedCandidate*>(v.daughter(dIdx)))->track();
124  if((*pf2pc)[trackRef].isNonnull() && (*pf2pc)[trackRef]->numberOfHits() > 0) {
125  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc)[trackRef]) ));
126  }
127  else {
128  if((*pf2pc2)[trackRef].isNonnull()) {
129  outPtr->back().addDaughter(reco::CandidatePtr(edm::refToPtr((*pf2pc2)[trackRef]) ));
130  }
131  else { edm::LogError("PATSecondaryVertexSlimmer") << "HELPME" << std::endl;}
132  }
133  }
134  }
135 
136  } // if reco::Vertex
137  } // if Candidate
138 
139  iEvent.put(std::move(outPtr));
140 }
141 
const edm::EDGetTokenT< std::vector< reco::Vertex > > srcLegacy_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
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
double vertexChi2() const override
chi-squares
double vertexCovariance(int i, int j) const override
(i, j)-th element of error matrix, i, j = 0, ... 2
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > map2_
std::vector< pat::PackedCandidate > PackedCandidateCollection
const edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > map_
PATSecondaryVertexSlimmer(const edm::ParameterSet &)
const Point & position() const
position
Definition: Vertex.h:109
std::vector< VertexCompositePtrCandidate > VertexCompositePtrCandidateCollection
collection of Candidate objects
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
virtual const daughters & daughterPtrVector() const
references to daughtes
const edm::EDGetTokenT< reco::VertexCompositeCandidateCollection > srcV0s_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:230
const Point & vertex() const override
vertex position (overwritten by PF...)
double chi2() const
chi-squares
Definition: Vertex.h:98
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
float trackWeight(const TREF &r) const
returns the weight with which a Track has contributed to the vertex-fit.
Definition: Vertex.h:81
bool isValid() const
Definition: HandleBase.h:74
virtual void clearDaughters()
clear daughter references
double ndof() const
Definition: Vertex.h:105
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:135
Error error() const
return SMatrix
Definition: Vertex.h:139
fixed size matrix
HLT enums.
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
void addDaughter(const CandidatePtr &)
add a daughter via a reference
size_type numberOfDaughters() const override
number of daughters
const edm::EDGetTokenT< reco::VertexCompositePtrCandidateCollection > src_
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
def move(src, dest)
Definition: eostools.py:511