CMS 3D CMS Logo

VertexArbitrators.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <set>
3 
8 
15 
22 
28 
34 
38 
40 
41 //#define VTXDEBUG
42 
43 inline const unsigned int nTracks(const reco::Vertex &sv) { return sv.nTracks(); }
44 inline const unsigned int nTracks(const reco::VertexCompositePtrCandidate &sv) {
45  return sv.numberOfSourceCandidatePtrs();
46 }
47 
48 template <class InputContainer, class VTX>
50 public:
51  typedef std::vector<VTX> Product;
53 
56  pdesc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
57  pdesc.add<edm::InputTag>("primaryVertices", edm::InputTag("offlinePrimaryVertices"));
59  pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
60  pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("vertexMerger"));
62  pdesc.add<edm::InputTag>("tracks", edm::InputTag("particleFlow"));
63  pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("candidateVertexMerger"));
64  } else {
65  pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
66  pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("vertexMerger"));
67  }
68  pdesc.add<double>("dLenFraction", 0.333);
69  pdesc.add<double>("dRCut", 0.4);
70  pdesc.add<double>("distCut", 0.04);
71  pdesc.add<double>("sigCut", 5.0);
72  pdesc.add<double>("fitterSigmacut", 3.0);
73  pdesc.add<double>("fitterTini", 256);
74  pdesc.add<double>("fitterRatio", 0.25);
75  pdesc.add<int>("trackMinLayers", 4);
76  pdesc.add<double>("trackMinPt", 0.4);
77  pdesc.add<int>("trackMinPixels", 1);
78  pdesc.add<double>("maxTimeSignificance", 3.5);
80  cdesc.add("trackVertexArbitratorDefault", pdesc);
82  cdesc.add("candidateVertexArbitratorDefault", pdesc);
83  } else {
84  cdesc.addDefault(pdesc);
85  }
86  }
87 
88  void produce(edm::Event &event, const edm::EventSetup &es) override;
89 
90 private:
91  bool trackFilter(const reco::TrackRef &track) const;
92 
98  std::unique_ptr<TrackVertexArbitration<VTX> > theArbitrator;
99 };
100 
101 template <class InputContainer, class VTX>
103  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
104  token_secondaryVertex = consumes<Product>(params.getParameter<edm::InputTag>("secondaryVertices"));
105  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
106  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
107  token_trackBuilder =
108  esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"));
109  produces<Product>();
110  theArbitrator.reset(new TrackVertexArbitration<VTX>(params));
111 }
112 
113 template <class InputContainer, class VTX>
115  using namespace reco;
116 
118  event.getByToken(token_secondaryVertex, secondaryVertices);
119  Product theSecVertexColl = *(secondaryVertices.product());
120 
122  event.getByToken(token_primaryVertex, primaryVertices);
123 
124  auto recoVertices = std::make_unique<Product>();
125  if (!primaryVertices->empty()) {
126  const reco::Vertex &pv = (*primaryVertices)[0];
127 
129  event.getByToken(token_tracks, tracks);
130 
131  edm::ESHandle<TransientTrackBuilder> trackBuilder = es.getHandle(token_trackBuilder);
132 
134  event.getByToken(token_beamSpot, beamSpot);
135 
136  std::vector<TransientTrack> selectedTracks;
137  for (typename InputContainer::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
138  selectedTracks.push_back(tthelpers::buildTT(tracks, trackBuilder, track - tracks->begin()));
139  }
140 
141  // const edm::RefVector< TrackCollection > tracksForArbitration= selectedTracks;
142  Product theRecoVertices = theArbitrator->trackVertexArbitrator(beamSpot, pv, selectedTracks, theSecVertexColl);
143 
144  for (unsigned int ivtx = 0; ivtx < theRecoVertices.size(); ivtx++) {
145  if (!(nTracks(theRecoVertices[ivtx]) > 1))
146  continue;
147  recoVertices->push_back(theRecoVertices[ivtx]);
148  }
149  }
150  event.put(std::move(recoVertices));
151 }
152 
156 
edm::EDGetTokenT< reco::BeamSpot > token_beamSpot
TemplatedVertexArbitrator< edm::View< reco::Candidate >, reco::VertexCompositePtrCandidate > CandidateVertexArbitrator
const unsigned int nTracks(const reco::Vertex &sv)
std::vector< VTX > Product
bool trackFilter(const reco::TrackRef &track) const
edm::EDGetTokenT< Product > token_secondaryVertex
void produce(edm::Event &event, const edm::EventSetup &es) 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
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
edm::EDGetTokenT< reco::VertexCollection > token_primaryVertex
auto const & tracks
cannot be loose
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > token_trackBuilder
reco::TransientTrack buildTT(edm::Handle< reco::TrackCollection > &tracks, edm::ESHandle< TransientTrackBuilder > &trackbuilder, unsigned int k)
Definition: TTHelpers.h:10
TemplatedVertexArbitrator< reco::TrackCollection, reco::Vertex > TrackVertexArbitrator
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1