CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/RecoVertex/AdaptiveVertexFinder/src/VertexMerging.cc

Go to the documentation of this file.
00001 #include "RecoVertex/AdaptiveVertexFinder/interface/VertexMerging.h"
00002 
00003 
00004 VertexMerging::VertexMerging(const edm::ParameterSet &params) :
00005         maxFraction(params.getParameter<double>("maxFraction")),
00006         minSignificance(params.getParameter<double>("minSignificance"))
00007 {
00008         
00009 }
00010 
00011 static double computeSharedTracks(const reco::Vertex &pv,
00012                                         const reco::Vertex &sv)
00013 {
00014         std::set<reco::TrackRef> pvTracks;
00015         for(std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin();
00016             iter != pv.tracks_end(); iter++) {
00017                 if (pv.trackWeight(*iter) >= 0.5)
00018                         pvTracks.insert(iter->castTo<reco::TrackRef>());
00019         }
00020 
00021         unsigned int count = 0, total = 0;
00022         for(std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin();
00023             iter != sv.tracks_end(); iter++) {
00024                 if (sv.trackWeight(*iter) >= 0.5) {
00025                         total++;
00026                         count += pvTracks.count(iter->castTo<reco::TrackRef>());
00027                 }
00028         }
00029 
00030         return (double)count / (double)total;
00031 }
00032 
00033 
00034 
00035 
00036 
00037 reco::VertexCollection VertexMerging::mergeVertex(
00038         reco::VertexCollection & secondaryVertices){
00039 
00040 
00041 
00042         
00043         VertexDistance3D dist;
00044         reco::VertexCollection recoVertices;
00045         for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices.begin();
00046             sv != secondaryVertices.end(); ++sv) {
00047           recoVertices.push_back(*sv);
00048         }
00049        for(std::vector<reco::Vertex>::iterator sv = recoVertices.begin();
00050             sv != recoVertices.end(); ++sv) {
00051 
00052         bool shared=false;
00053        for(std::vector<reco::Vertex>::iterator sv2 = recoVertices.begin();
00054             sv2 != recoVertices.end(); ++sv2) {
00055                   double fr=computeSharedTracks(*sv2, *sv);
00056         //        std::cout << sv2-recoVertices->begin() << " vs " << sv-recoVertices->begin() << " : " << fr << " "  <<  computeSharedTracks(*sv, *sv2) << " sig " << dist.distance(*sv,*sv2).significance() << std::endl;
00057           //      std::cout << (fr > maxFraction) << " && " << (dist.distance(*sv,*sv2).significance() < 2)  <<  " && " <<  (sv-sv2!=0)  << " && " <<  (fr >= computeSharedTracks(*sv2, *sv))  << std::endl;
00058                 if (fr > maxFraction && dist.distance(*sv,*sv2).significance() < minSignificance && sv-sv2!=0 
00059                     && fr >= computeSharedTracks(*sv, *sv2) )
00060                   {
00061                       shared=true; 
00062                      // std::cout << "shared " << sv-recoVertices->begin() << " and "  << sv2-recoVertices->begin() << " fractions: " << fr << " , "  << computeSharedTracks(*sv2, *sv) << " sig: " <<  dist.distance(*sv,*sv2).significance() <<  std::endl;
00063          
00064                   }
00065                  
00066 
00067         }
00068         if(shared) { sv=recoVertices.erase(sv)-1; }
00069      //    std::cout << "it = " <<  sv-recoVertices->begin() << " new size is: " << recoVertices->size() <<   std::endl;
00070        }
00071 
00072         return recoVertices;
00073 
00074 
00075 }
00076