CMS 3D CMS Logo

VertexAssociatorByTracks.cc

Go to the documentation of this file.
00001 #include "DataFormats/Common/interface/EDProduct.h"
00002 #include "DataFormats/Common/interface/Ref.h"
00003 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00004 #include "DataFormats/VertexReco/interface/Vertex.h"
00005 
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00012 
00013 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00014 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
00015 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
00016 
00017 #include "SimTracker/TrackAssociation/interface/TrackAssociatorBase.h"
00018 #include "SimTracker/VertexAssociation/interface/VertexAssociatorByTracks.h"
00019 
00020 using namespace reco;
00021 using namespace edm;
00022 using namespace std;
00023 
00024 /* Constructor */
00025 VertexAssociatorByTracks::VertexAssociatorByTracks (const edm::ParameterSet& conf) :
00026   conf_(conf) {}
00027 
00028 
00029 /* Destructor */
00030 VertexAssociatorByTracks::~VertexAssociatorByTracks() {
00031   //do cleanup here
00032 }
00033 
00034 //
00035 //---member functions
00036 //
00037 
00038 VertexRecoToSimCollection VertexAssociatorByTracks::associateRecoToSim(
00039     edm::Handle<reco::VertexCollection>& vertexCollectionH,
00040     edm::Handle<TrackingVertexCollection>&  TVCollectionH,
00041     const edm::Event& event,
00042     reco::RecoToSimCollection& trackAssocResult) const {
00043 
00044 //  const double minHitFraction = theMinHitFraction;
00045 //  int nshared =0;
00046 //  float fraction=0;
00047 //  std::vector<unsigned int> SimTrackIds;
00048 //  std::vector<unsigned int> matchedIds;
00049 
00050   using reco::VertexRef;
00051 
00052   VertexRecoToSimCollection  outputCollection;
00053 
00054   const TrackingVertexCollection tVC = *(TVCollectionH.product());
00055   const   reco::VertexCollection  vC = *(vertexCollectionH.product());
00056 
00057 //  double minFraction = 0.01;
00058 //  double fraction = 0.8;
00059 
00060   std::map<TrackingVertexRef,int> tVCount;
00061 
00062   int iv = 0;
00063 
00064   // Loop over reco::Vertex
00065 
00066   for (reco::VertexCollection::const_iterator vertex = vC.begin();
00067        vertex != vC.end(); ++vertex, ++iv) {
00068     tVCount.clear();
00069     VertexRef rVertexR = VertexRef(vertexCollectionH,iv);
00070     double nRecoTracks = vertex->tracksSize();
00071 
00072     // Loop over daughter tracks of reco::Vertex
00073 
00074     for (reco::Vertex::trackRef_iterator recoDaughter = vertex->tracks_begin();
00075          recoDaughter != vertex->tracks_end(); ++recoDaughter) {
00076       RefToBase<reco::Track>  tr = recoDaughter->castTo<RefToBase<reco::Track> >();
00077       if (trackAssocResult[tr].size() > 0) {
00078         std::vector<std::pair<TrackingParticleRef, double> > tpV = trackAssocResult[tr];
00079 
00080         // Loop over TrackingParticles associated with reco::Track
00081 
00082         for (std::vector<std::pair<TrackingParticleRef, double> >::const_iterator match = tpV.begin();
00083              match != tpV.end(); ++match) {
00084           // ... and keep count of it's parent vertex
00085           TrackingParticleRef tp = match->first;
00086 //          double trackFraction = match->second;
00087           TrackingVertexRef   tv = tp->parentVertex();
00088           ++tVCount[tv]; // Count matches to this reco:Vertex for this TrackingVertex
00089         }
00090       }
00091     }
00092 
00093     // Loop over map, set score, add to outputCollection
00094 
00095     for (std::map<TrackingVertexRef,int>::const_iterator match = tVCount.begin();
00096          match != tVCount.end(); ++match) {
00097       TrackingVertexRef tV = match->first;
00098       double nMatches      = match->second;
00099       outputCollection.insert(rVertexR,std::make_pair(tV,nMatches/nRecoTracks));
00100     }
00101   } // Loop on reco::Vertex
00102 
00103   return outputCollection;
00104 }
00105 
00106 
00107 VertexSimToRecoCollection VertexAssociatorByTracks::associateSimToReco(
00108     edm::Handle<reco::VertexCollection>&   vertexCollectionH,
00109     edm::Handle<TrackingVertexCollection>& TVCollectionH,
00110     const edm::Event& e,
00111     reco::SimToRecoCollection& trackAssocResult) const {
00112 
00113   const TrackingVertexCollection tVC = *(TVCollectionH.product());
00114   const   reco::VertexCollection  vC = *(vertexCollectionH.product());
00115 
00116   VertexSimToRecoCollection  outputCollection; // return value
00117 
00118   // Loop over TrackingVertexes
00119 
00120   std::map<VertexRef,int> vCount;
00121   int iTV = 0;
00122   for (TrackingVertexCollection::const_iterator tV = tVC.begin();
00123        tV != tVC.end(); ++tV, ++iTV) {
00124     vCount.clear();
00125     TrackingVertexRef tVertexR = TrackingVertexRef(TVCollectionH,iTV);
00126     double nSimTracks = (tV->daughterTracks()).size();
00127 
00128     // Loop over daughter tracks of TrackingVertex
00129     for (TrackingVertex::tp_iterator simDaughter = tV->daughterTracks_begin();
00130          simDaughter != tV->daughterTracks_end(); ++simDaughter) {
00131       TrackingParticleRef tp = *simDaughter;
00132 
00133       SimToRecoCollection::const_iterator daughterPosition = trackAssocResult.find(*simDaughter);
00134       if (daughterPosition != trackAssocResult.end()) {
00135         std::vector<std::pair<RefToBase<reco::Track> , double> > recoTracks = trackAssocResult[*simDaughter];
00136 
00137        // Loop over reco::Tracks associated with TrackingParticle
00138         for (std::vector<std::pair<RefToBase<reco::Track> , double> >::const_iterator match = recoTracks.begin();
00139              match != recoTracks.end(); ++match) {
00140           // ... and keep count of it's parent vertex
00141 
00142           TrackBaseRef track(match->first);
00143 //          double   trackQuality = match->second;
00144 
00145           // Find vertex if any where this track comes from
00146           int iv = 0;
00147           for (reco::VertexCollection::const_iterator vertex = vC.begin();
00148               vertex != vC.end(); ++vertex,++iv) {
00149             VertexRef rVertexR = VertexRef(vertexCollectionH,iv);
00150             for (reco::Vertex::trackRef_iterator recoDaughter = vertex->tracks_begin();
00151                  recoDaughter != vertex->tracks_end(); ++recoDaughter) {
00152               if (*recoDaughter == track) {
00153                 ++vCount[rVertexR]; // Count matches to this TrackingVertex for this reco:Vertex
00154               }
00155             }
00156           }
00157         }
00158       }
00159     }
00160 
00161     // Loop over map, set score, add to outputCollection
00162     for (std::map<VertexRef,int>::const_iterator match = vCount.begin(); match != vCount.end(); ++match) {
00163       VertexRef v = match->first;
00164       double nMatches      = match->second;
00165       outputCollection.insert(tVertexR,std::make_pair(v,nMatches/nSimTracks));
00166     }
00167   }
00168 
00169   return outputCollection;
00170 }
00171 

Generated on Tue Jun 9 17:48:02 2009 for CMSSW by  doxygen 1.5.4