Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:02

0001 #include "RecoVertex/AdaptiveVertexFinder/interface/VertexMerging.h"
0002 
0003 VertexMerging::VertexMerging(const edm::ParameterSet &params)
0004     : maxFraction(params.getParameter<double>("maxFraction")),
0005       minSignificance(params.getParameter<double>("minSignificance")) {}
0006 
0007 static double computeSharedTracks(const reco::Vertex &pv, const reco::Vertex &sv) {
0008   std::set<reco::TrackRef> pvTracks;
0009   for (std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin(); iter != pv.tracks_end(); iter++) {
0010     if (pv.trackWeight(*iter) >= 0.5)
0011       pvTracks.insert(iter->castTo<reco::TrackRef>());
0012   }
0013 
0014   unsigned int count = 0, total = 0;
0015   for (std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin(); iter != sv.tracks_end(); iter++) {
0016     if (sv.trackWeight(*iter) >= 0.5) {
0017       total++;
0018       count += pvTracks.count(iter->castTo<reco::TrackRef>());
0019     }
0020   }
0021 
0022   return (double)count / (double)total;
0023 }
0024 
0025 reco::VertexCollection VertexMerging::mergeVertex(reco::VertexCollection &secondaryVertices) {
0026   VertexDistance3D dist;
0027   reco::VertexCollection recoVertices;
0028   for (std::vector<reco::Vertex>::const_iterator sv = secondaryVertices.begin(); sv != secondaryVertices.end(); ++sv) {
0029     recoVertices.push_back(*sv);
0030   }
0031   for (std::vector<reco::Vertex>::iterator sv = recoVertices.begin(); sv != recoVertices.end(); ++sv) {
0032     bool shared = false;
0033     for (std::vector<reco::Vertex>::iterator sv2 = recoVertices.begin(); sv2 != recoVertices.end(); ++sv2) {
0034       double fr = computeSharedTracks(*sv2, *sv);
0035       //        std::cout << sv2-recoVertices->begin() << " vs " << sv-recoVertices->begin() << " : " << fr << " "  <<  computeSharedTracks(*sv, *sv2) << " sig " << dist.distance(*sv,*sv2).significance() << std::endl;
0036       //      std::cout << (fr > maxFraction) << " && " << (dist.distance(*sv,*sv2).significance() < 2)  <<  " && " <<  (sv-sv2!=0)  << " && " <<  (fr >= computeSharedTracks(*sv2, *sv))  << std::endl;
0037       if (fr > maxFraction && dist.distance(*sv, *sv2).significance() < minSignificance && sv - sv2 != 0 &&
0038           fr >= computeSharedTracks(*sv, *sv2)) {
0039         shared = true;
0040         // std::cout << "shared " << sv-recoVertices->begin() << " and "  << sv2-recoVertices->begin() << " fractions: " << fr << " , "  << computeSharedTracks(*sv2, *sv) << " sig: " <<  dist.distance(*sv,*sv2).significance() <<  std::endl;
0041       }
0042     }
0043     if (shared) {
0044       sv = recoVertices.erase(sv) - 1;
0045     }
0046     //    std::cout << "it = " <<  sv-recoVertices->begin() << " new size is: " << recoVertices->size() <<   std::endl;
0047   }
0048 
0049   return recoVertices;
0050 }