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 
59  pdesc.add<edm::InputTag>("beamSpot",edm::InputTag("offlineBeamSpot"));
60  pdesc.add<edm::InputTag>("primaryVertices",edm::InputTag("offlinePrimaryVertices"));
62  pdesc.add<edm::InputTag>("tracks",edm::InputTag("generalTracks"));
63  pdesc.add<edm::InputTag>("secondaryVertices",edm::InputTag("vertexMerger"));
65  pdesc.add<edm::InputTag>("tracks",edm::InputTag("particleFlow"));
66  pdesc.add<edm::InputTag>("secondaryVertices",edm::InputTag("candidateVertexMerger"));
67  } else {
68  pdesc.add<edm::InputTag>("tracks",edm::InputTag("generalTracks"));
69  pdesc.add<edm::InputTag>("secondaryVertices",edm::InputTag("vertexMerger"));
70  }
71  pdesc.add<double>("dLenFraction",0.333);
72  pdesc.add<double>("dRCut",0.4);
73  pdesc.add<double>("distCut",0.04);
74  pdesc.add<double>("sigCut",5.0);
75  pdesc.add<double>("fitterSigmacut",3.0);
76  pdesc.add<double>("fitterTini",256);
77  pdesc.add<double>("fitterRatio",0.25);
78  pdesc.add<int>("trackMinLayers",4);
79  pdesc.add<double>("trackMinPt",0.4);
80  pdesc.add<int>("trackMinPixels",1);
81  pdesc.add<double>("maxTimeSignificance",3.5);
83  cdesc.add("trackVertexArbitratorDefault",pdesc);
85  cdesc.add("candidateVertexArbitratorDefault",pdesc);
86  } else {
87  cdesc.addDefault(pdesc);
88  }
89  }
90 
91  virtual void produce(edm::Event &event, const edm::EventSetup &es) override ;
92 
93  private:
94  bool trackFilter(const reco::TrackRef &track) const;
95 
100  std::unique_ptr<TrackVertexArbitration<VTX> > theArbitrator;
101 };
102 
103 
104 template <class InputContainer, class VTX>
106 {
107  token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
108  token_secondaryVertex = consumes<Product>(params.getParameter<edm::InputTag>("secondaryVertices"));
109  token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
110  token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
111  produces<Product>();
112  theArbitrator.reset( new TrackVertexArbitration<VTX>(params) );
113 }
114 
115 template <class InputContainer, class VTX>
117 {
118  using namespace reco;
119 
121  event.getByToken(token_secondaryVertex, secondaryVertices);
122  Product theSecVertexColl = *(secondaryVertices.product());
123 
125  event.getByToken(token_primaryVertex, primaryVertices);
126 
127  auto recoVertices = std::make_unique<Product>();
128  if(primaryVertices->size()!=0){
129  const reco::Vertex &pv = (*primaryVertices)[0];
130 
132  event.getByToken(token_tracks, tracks);
133 
135  es.get<TransientTrackRecord>().get("TransientTrackBuilder",
136  trackBuilder);
137 
139  event.getByToken(token_beamSpot,beamSpot);
140 
141  std::vector<TransientTrack> selectedTracks;
142  for(typename InputContainer::const_iterator track = tracks->begin();
143  track != tracks->end(); ++track) {
144  selectedTracks.push_back(tthelpers::buildTT(tracks,trackBuilder,track - tracks->begin()));
145  }
146 
147 
148  // const edm::RefVector< TrackCollection > tracksForArbitration= selectedTracks;
149  Product theRecoVertices = theArbitrator->trackVertexArbitrator(beamSpot, pv, selectedTracks,
150  theSecVertexColl);
151 
152  for(unsigned int ivtx=0; ivtx < theRecoVertices.size(); ivtx++){
153  if ( !(nTracks(theRecoVertices[ivtx]) > 1) ) continue;
154  recoVertices->push_back(theRecoVertices[ivtx]);
155  }
156 
157  }
158  event.put(std::move(recoVertices));
159 
160 
161 
162 }
163 
164 
165 
166 #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
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:81
primaryVertices
Definition: jets_cff.py:130
unsigned int nTracks(float minWeight=0.5) const
Returns the number of tracks in the vertex with weight above minWeight.
Definition: Vertex.cc:169
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:62
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