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 
10 
16 
23 
29 
35 
39 
41 
42 //#define VTXDEBUG
43 
44 const unsigned int nTracks(const reco::Vertex &sv) { return sv.nTracks(); }
45 const unsigned int nTracks(const reco::VertexCompositePtrCandidate &sv) { return sv.numberOfSourceCandidatePtrs(); }
46 
47 template <class InputContainer, class VTX>
49 public:
50  typedef std::vector<VTX> Product;
52 
55  pdesc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
56  pdesc.add<edm::InputTag>("primaryVertices", edm::InputTag("offlinePrimaryVertices"));
58  pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
59  pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("vertexMerger"));
61  pdesc.add<edm::InputTag>("tracks", edm::InputTag("particleFlow"));
62  pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("candidateVertexMerger"));
63  } else {
64  pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
65  pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("vertexMerger"));
66  }
67  pdesc.add<double>("dLenFraction", 0.333);
68  pdesc.add<double>("dRCut", 0.4);
69  pdesc.add<double>("distCut", 0.04);
70  pdesc.add<double>("sigCut", 5.0);
71  pdesc.add<double>("fitterSigmacut", 3.0);
72  pdesc.add<double>("fitterTini", 256);
73  pdesc.add<double>("fitterRatio", 0.25);
74  pdesc.add<int>("trackMinLayers", 4);
75  pdesc.add<double>("trackMinPt", 0.4);
76  pdesc.add<int>("trackMinPixels", 1);
77  pdesc.add<double>("maxTimeSignificance", 3.5);
79  cdesc.add("trackVertexArbitratorDefault", pdesc);
81  cdesc.add("candidateVertexArbitratorDefault", pdesc);
82  } else {
83  cdesc.addDefault(pdesc);
84  }
85  }
86 
87  void produce(edm::Event &event, const edm::EventSetup &es) override;
88 
89 private:
90  bool trackFilter(const reco::TrackRef &track) const;
91 
96  std::unique_ptr<TrackVertexArbitration<VTX> > theArbitrator;
97 };
98 
99 template <class InputContainer, class VTX>
101  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
102  token_secondaryVertex = consumes<Product>(params.getParameter<edm::InputTag>("secondaryVertices"));
103  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
104  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
105  produces<Product>();
106  theArbitrator.reset(new TrackVertexArbitration<VTX>(params));
107 }
108 
109 template <class InputContainer, class VTX>
111  using namespace reco;
112 
114  event.getByToken(token_secondaryVertex, secondaryVertices);
115  Product theSecVertexColl = *(secondaryVertices.product());
116 
118  event.getByToken(token_primaryVertex, primaryVertices);
119 
120  auto recoVertices = std::make_unique<Product>();
121  if (!primaryVertices->empty()) {
122  const reco::Vertex &pv = (*primaryVertices)[0];
123 
125  event.getByToken(token_tracks, tracks);
126 
128  es.get<TransientTrackRecord>().get("TransientTrackBuilder", trackBuilder);
129 
131  event.getByToken(token_beamSpot, beamSpot);
132 
133  std::vector<TransientTrack> selectedTracks;
134  for (typename InputContainer::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
135  selectedTracks.push_back(tthelpers::buildTT(tracks, trackBuilder, track - tracks->begin()));
136  }
137 
138  // const edm::RefVector< TrackCollection > tracksForArbitration= selectedTracks;
139  Product theRecoVertices = theArbitrator->trackVertexArbitrator(beamSpot, pv, selectedTracks, theSecVertexColl);
140 
141  for (unsigned int ivtx = 0; ivtx < theRecoVertices.size(); ivtx++) {
142  if (!(nTracks(theRecoVertices[ivtx]) > 1))
143  continue;
144  recoVertices->push_back(theRecoVertices[ivtx]);
145  }
146  }
147  event.put(std::move(recoVertices));
148 }
149 
150 #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
void produce(edm::Event &event, const edm::EventSetup &es) override
size_type numberOfSourceCandidatePtrs() const override
std::unique_ptr< TrackVertexArbitration< VTX > > theArbitrator
static void fillDescriptions(edm::ConfigurationDescriptions &cdesc)
edm::EDGetTokenT< InputContainer > token_tracks
void addDefault(ParameterSetDescription const &psetDescription)
TemplatedVertexArbitrator(const edm::ParameterSet &params)
def pv(vc)
Definition: MetAnalyzer.py:7
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
T const * product() const
Definition: Handle.h:69
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:146
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool trackFilter(const reco::TrackRef &track) const
fixed size matrix
T get() const
Definition: EventSetup.h:73
reco::TransientTrack buildTT(edm::Handle< reco::TrackCollection > &tracks, edm::ESHandle< TransientTrackBuilder > &trackbuilder, unsigned int k)
Definition: TTHelpers.h:10
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1