CMS 3D CMS Logo

TemplatedVertexArbitrator.h
Go to the documentation of this file.
1 #ifndef TemplatedVertexArbitrator_h
2 #define TemplatedVertexArbitrator_h
3 #include <memory>
4 #include <set>
5 
6 
11 
12 
18 
25 
26 
32 
38 
39 
43 
45 
46 //#define VTXDEBUG
47 
48 const unsigned int nTracks(const reco::Vertex & sv) {return sv.nTracks();}
49 const unsigned int nTracks(const reco::VertexCompositePtrCandidate & sv) {return sv.numberOfSourceCandidatePtrs();}
50 
51 template <class InputContainer, class VTX>
53  public:
54  typedef std::vector<VTX> Product;
56 
57 
58  virtual void produce(edm::Event &event, const edm::EventSetup &es) override ;
59 
60  private:
61  bool trackFilter(const reco::TrackRef &track) const;
62 
67  std::unique_ptr<TrackVertexArbitration<VTX> > theArbitrator;
68 };
69 
70 
71 template <class InputContainer, class VTX>
73 {
74  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
75  token_secondaryVertex = consumes<Product>(params.getParameter<edm::InputTag>("secondaryVertices"));
76  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
77  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
78  produces<Product>();
79  theArbitrator.reset( new TrackVertexArbitration<VTX>(params) );
80 }
81 
82 template <class InputContainer, class VTX>
84 {
85  using namespace reco;
86 
88  event.getByToken(token_secondaryVertex, secondaryVertices);
89  Product theSecVertexColl = *(secondaryVertices.product());
90 
92  event.getByToken(token_primaryVertex, primaryVertices);
93 
94  auto recoVertices = std::make_unique<Product>();
95  if(primaryVertices->size()!=0){
96  const reco::Vertex &pv = (*primaryVertices)[0];
97 
99  event.getByToken(token_tracks, tracks);
100 
102  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
103  trackBuilder);
104 
106  event.getByToken(token_beamSpot,beamSpot);
107 
108  std::vector<TransientTrack> selectedTracks;
109  for(typename InputContainer::const_iterator track = tracks->begin();
110  track != tracks->end(); ++track) {
111  selectedTracks.push_back(tthelpers::buildTT(tracks,trackBuilder,track - tracks->begin()));
112  }
113 
114 
115  // const edm::RefVector< TrackCollection > tracksForArbitration= selectedTracks;
116  Product theRecoVertices = theArbitrator->trackVertexArbitrator(beamSpot, pv, selectedTracks,
117  theSecVertexColl);
118 
119  for(unsigned int ivtx=0; ivtx < theRecoVertices.size(); ivtx++){
120  if ( !(nTracks(theRecoVertices[ivtx]) > 1) ) continue;
121  recoVertices->push_back(theRecoVertices[ivtx]);
122  }
123 
124  }
125  event.put(std::move(recoVertices));
126 
127 
128 
129 }
130 
131 
132 
133 #endif
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
T getParameter(std::string const &) const
const unsigned int nTracks(const reco::Vertex &sv)
edm::EDGetTokenT< Product > token_secondaryVertex
virtual void produce(edm::Event &event, const edm::EventSetup &es) override
std::unique_ptr< TrackVertexArbitration< VTX > > theArbitrator
edm::EDGetTokenT< InputContainer > token_tracks
TemplatedVertexArbitrator(const edm::ParameterSet &params)
def pv(vc)
Definition: MetAnalyzer.py:6
virtual size_type numberOfSourceCandidatePtrs() const
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
T const * product() const
Definition: Handle.h:81
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:169
const T & get() const
Definition: EventSetup.h:55
bool trackFilter(const reco::TrackRef &track) const
fixed size matrix
reco::TransientTrack buildTT(edm::Handle< reco::TrackCollection > &tracks, edm::ESHandle< TransientTrackBuilder > &trackbuilder, unsigned int k)
Definition: TTHelpers.h:8
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1