CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VertexMerger.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <set>
3 
9 
16 
17 class VertexMerger : public edm::EDProducer {
18  public:
19  VertexMerger(const edm::ParameterSet &params);
20 
21  virtual void produce(edm::Event &event, const edm::EventSetup &es);
22 
23  private:
24  bool trackFilter(const reco::TrackRef &track) const;
25 
26 // edm::InputTag primaryVertexCollection;
28  double maxFraction;
30 };
31 
33 // primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
34  secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
35  maxFraction(params.getParameter<double>("maxFraction")),
36  minSignificance(params.getParameter<double>("minSignificance"))
37 {
38  produces<reco::VertexCollection>();
39 }
40 
41 static double computeSharedTracks(const reco::Vertex &pv,
42  const reco::Vertex &sv)
43 {
44  std::set<reco::TrackRef> pvTracks;
45  for(std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin();
46  iter != pv.tracks_end(); iter++) {
47  if (pv.trackWeight(*iter) >= 0.5)
48  pvTracks.insert(iter->castTo<reco::TrackRef>());
49  }
50 
51  unsigned int count = 0, total = 0;
52  for(std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin();
53  iter != sv.tracks_end(); iter++) {
54  if (sv.trackWeight(*iter) >= 0.5) {
55  total++;
56  count += pvTracks.count(iter->castTo<reco::TrackRef>());
57  }
58  }
59 
60  return (double)count / (double)total;
61 }
62 
64 {
65  using namespace reco;
66 
67  edm::Handle<VertexCollection> secondaryVertices;
68  event.getByLabel(secondaryVertexCollection, secondaryVertices);
69 
70  VertexDistance3D dist;
71  std::auto_ptr<VertexCollection> recoVertices(new VertexCollection);
72  for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
73  sv != secondaryVertices->end(); ++sv) {
74  recoVertices->push_back(*sv);
75  }
76  for(std::vector<reco::Vertex>::iterator sv = recoVertices->begin();
77  sv != recoVertices->end(); ++sv) {
78 
79  bool shared=false;
80  for(std::vector<reco::Vertex>::iterator sv2 = recoVertices->begin();
81  sv2 != recoVertices->end(); ++sv2) {
82  double fr=computeSharedTracks(*sv2, *sv);
83  // std::cout << sv2-recoVertices->begin() << " vs " << sv-recoVertices->begin() << " : " << fr << " " << computeSharedTracks(*sv, *sv2) << " sig " << dist.distance(*sv,*sv2).significance() << std::endl;
84  // std::cout << (fr > maxFraction) << " && " << (dist.distance(*sv,*sv2).significance() < 2) << " && " << (sv-sv2!=0) << " && " << (fr >= computeSharedTracks(*sv2, *sv)) << std::endl;
85  if (fr > maxFraction && dist.distance(*sv,*sv2).significance() < minSignificance && sv-sv2!=0
86  && fr >= computeSharedTracks(*sv, *sv2) )
87  {
88  shared=true;
89  // std::cout << "shared " << sv-recoVertices->begin() << " and " << sv2-recoVertices->begin() << " fractions: " << fr << " , " << computeSharedTracks(*sv2, *sv) << " sig: " << dist.distance(*sv,*sv2).significance() << std::endl;
90 
91  }
92 
93 
94  }
95  if(shared) { sv=recoVertices->erase(sv)-1; }
96  // std::cout << "it = " << sv-recoVertices->begin() << " new size is: " << recoVertices->size() << std::endl;
97  }
98 
99  event.put(recoVertices);
100 
101 
102 
103 }
104 
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:45
VertexMerger(const edm::ParameterSet &params)
Definition: VertexMerger.cc:32
edm::InputTag secondaryVertexCollection
Definition: VertexMerger.cc:27
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool trackFilter(const reco::TrackRef &track) const
static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv)
Definition: VertexMerger.cc:41
double maxFraction
Definition: VertexMerger.cc:28
virtual void produce(edm::Event &event, const edm::EventSetup &es)
Definition: VertexMerger.cc:63
float trackWeight(const TrackBaseRef &r) const
returns the weight with which a Track has contributed to the vertex-fit.
double minSignificance
Definition: VertexMerger.cc:29
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:40