CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoVertex/AdaptiveVertexFinder/plugins/VertexMerger.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include <set>
00003 
00004 #include "FWCore/Framework/interface/EDProducer.h"
00005 #include "FWCore/Framework/interface/Event.h"
00006 #include "FWCore/Framework/interface/MakerMacros.h"
00007 #include "FWCore/Utilities/interface/InputTag.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "DataFormats/TrackReco/interface/Track.h"
00012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00013 #include "DataFormats/VertexReco/interface/Vertex.h"
00014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00015 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
00016 
00017 class VertexMerger : public edm::EDProducer {
00018     public:
00019         VertexMerger(const edm::ParameterSet &params);
00020 
00021         virtual void produce(edm::Event &event, const edm::EventSetup &es);
00022 
00023     private:
00024         bool trackFilter(const reco::TrackRef &track) const;
00025 
00026 //      edm::InputTag                           primaryVertexCollection;
00027         edm::InputTag                           secondaryVertexCollection;
00028         double                                  maxFraction;
00029         double                                  minSignificance;
00030 };
00031 
00032 VertexMerger::VertexMerger(const edm::ParameterSet &params) :
00033 //      primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
00034         secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
00035         maxFraction(params.getParameter<double>("maxFraction")),
00036         minSignificance(params.getParameter<double>("minSignificance"))
00037 {
00038         produces<reco::VertexCollection>();
00039 }
00040 
00041 static double computeSharedTracks(const reco::Vertex &pv,
00042                                         const reco::Vertex &sv)
00043 {
00044         std::set<reco::TrackRef> pvTracks;
00045         for(std::vector<reco::TrackBaseRef>::const_iterator iter = pv.tracks_begin();
00046             iter != pv.tracks_end(); iter++) {
00047                 if (pv.trackWeight(*iter) >= 0.5)
00048                         pvTracks.insert(iter->castTo<reco::TrackRef>());
00049         }
00050 
00051         unsigned int count = 0, total = 0;
00052         for(std::vector<reco::TrackBaseRef>::const_iterator iter = sv.tracks_begin();
00053             iter != sv.tracks_end(); iter++) {
00054                 if (sv.trackWeight(*iter) >= 0.5) {
00055                         total++;
00056                         count += pvTracks.count(iter->castTo<reco::TrackRef>());
00057                 }
00058         }
00059 
00060         return (double)count / (double)total;
00061 }
00062 
00063 void VertexMerger::produce(edm::Event &event, const edm::EventSetup &es)
00064 {
00065         using namespace reco;
00066 
00067         edm::Handle<VertexCollection> secondaryVertices;
00068         event.getByLabel(secondaryVertexCollection, secondaryVertices);
00069 
00070         VertexDistance3D dist;
00071         std::auto_ptr<VertexCollection> recoVertices(new VertexCollection);
00072         for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
00073             sv != secondaryVertices->end(); ++sv) {
00074           recoVertices->push_back(*sv);
00075         }
00076        for(std::vector<reco::Vertex>::iterator sv = recoVertices->begin();
00077             sv != recoVertices->end(); ++sv) {
00078 
00079         bool shared=false;
00080        for(std::vector<reco::Vertex>::iterator sv2 = recoVertices->begin();
00081             sv2 != recoVertices->end(); ++sv2) {
00082                   double fr=computeSharedTracks(*sv2, *sv);
00083         //        std::cout << sv2-recoVertices->begin() << " vs " << sv-recoVertices->begin() << " : " << fr << " "  <<  computeSharedTracks(*sv, *sv2) << " sig " << dist.distance(*sv,*sv2).significance() << std::endl;
00084           //      std::cout << (fr > maxFraction) << " && " << (dist.distance(*sv,*sv2).significance() < 2)  <<  " && " <<  (sv-sv2!=0)  << " && " <<  (fr >= computeSharedTracks(*sv2, *sv))  << std::endl;
00085                 if (fr > maxFraction && dist.distance(*sv,*sv2).significance() < minSignificance && sv-sv2!=0 
00086                     && fr >= computeSharedTracks(*sv, *sv2) )
00087                   {
00088                       shared=true; 
00089                      // std::cout << "shared " << sv-recoVertices->begin() << " and "  << sv2-recoVertices->begin() << " fractions: " << fr << " , "  << computeSharedTracks(*sv2, *sv) << " sig: " <<  dist.distance(*sv,*sv2).significance() <<  std::endl;
00090          
00091                   }
00092                  
00093 
00094         }
00095         if(shared) { sv=recoVertices->erase(sv)-1; }
00096      //    std::cout << "it = " <<  sv-recoVertices->begin() << " new size is: " << recoVertices->size() <<   std::endl;
00097        }
00098 
00099         event.put(recoVertices);
00100 
00101 
00102 
00103 }
00104 
00105 DEFINE_FWK_MODULE(VertexMerger);