CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/SimTracker/VertexAssociation/src/VertexAssociatorByTracks.cc

Go to the documentation of this file.
00001 #include "DataFormats/Common/interface/Ref.h"
00002 #include "DataFormats/Common/interface/RefToBase.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 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00021 
00022 
00023 /* Constructor */
00024 VertexAssociatorByTracks::VertexAssociatorByTracks (const edm::ParameterSet & config) : config_(config)
00025 {
00026     R2SMatchedSimRatio_ = config.getParameter<double>("R2SMatchedSimRatio");
00027     R2SMatchedRecoRatio_ = config.getParameter<double>("R2SMatchedRecoRatio");
00028     S2RMatchedSimRatio_ = config.getParameter<double>("S2RMatchedSimRatio");
00029     S2RMatchedRecoRatio_ = config.getParameter<double>("S2RMatchedRecoRatio");
00030 
00031     std::string trackQualityType = config.getParameter<std::string>("trackQuality");
00032     trackQuality_ = reco::TrackBase::qualityByName(trackQualityType);
00033 
00034     edm::ParameterSet param = config.getParameter<edm::ParameterSet>("trackingParticleSelector");
00035 
00036     selector_ = TrackingParticleSelector(
00037                     param.getParameter<double>("ptMinTP"),
00038                     param.getParameter<double>("minRapidityTP"),
00039                     param.getParameter<double>("maxRapidityTP"),
00040                     param.getParameter<double>("tipTP"),
00041                     param.getParameter<double>("lipTP"),
00042                     param.getParameter<int>("minHitTP"),
00043                     param.getParameter<bool>("signalOnlyTP"),
00044                     param.getParameter<bool>("chargedOnlyTP"),
00045                     param.getParameter<bool>("stableOnlyTP"),
00046                     param.getParameter<std::vector<int> >("pdgIdTP")
00047                 );
00048 }
00049 
00050 
00051 /* Destructor */
00052 VertexAssociatorByTracks::~VertexAssociatorByTracks()
00053 {
00054     //do cleanup here
00055 }
00056 
00057 
00058 reco::VertexRecoToSimCollection VertexAssociatorByTracks::associateRecoToSim(
00059     edm::Handle<edm::View<reco::Vertex> > & recoVertexes,
00060     edm::Handle<TrackingVertexCollection> & trackingVertexes,
00061     const edm::Event& event,
00062     reco::RecoToSimCollection & associator
00063 ) const
00064 {
00065     reco::VertexRecoToSimCollection  outputCollection;
00066 
00067     std::map<TrackingVertexRef,std::pair<double, std::size_t> > matches;
00068 
00069     std::cout << "reco::VertexCollection size = " << recoVertexes->size();
00070     std::cout << " ; TrackingVertexCollection size = " << trackingVertexes->size() << std::endl << std::endl;
00071 
00072     // Loop over RecoVertex
00073     for (std::size_t recoIndex = 0; recoIndex < recoVertexes->size(); ++recoIndex)
00074     {
00075         reco::VertexBaseRef recoVertex(recoVertexes, recoIndex); 
00076 
00077         matches.clear();
00078 
00079         double recoDaughterWeight = 0.;
00080 
00081         // Loop over daughter tracks of RecoVertex
00082         for (
00083             reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
00084             recoDaughter != recoVertex->tracks_end();
00085             ++recoDaughter
00086         )
00087         {
00088             // Check the quality of the RecoDoughters
00089             if ( !(*recoDaughter)->quality(trackQuality_) ) continue;
00090 
00091             // Check for association for the given RecoDaughter
00092             if ( associator.numberOfAssociations(*recoDaughter) > 0 )
00093             {
00094                 std::vector<std::pair<TrackingParticleRef,double> > associations = associator[*recoDaughter];
00095 
00096                 // Loop over TrackingParticles associated with RecoDaughter
00097                 for (
00098                     std::vector<std::pair<TrackingParticleRef,double> >::const_iterator association = associations.begin();
00099                     association != associations.end();
00100                     ++association
00101                 )
00102                 {
00103                     // Get a reference to parent vertex of TrackingParticle associated to the RecoDaughter
00104                     TrackingVertexRef trackingVertex = association->first->parentVertex();
00105                     // Store matched RecoDaughter to the trackingVertex
00106                     matches[trackingVertex].first += recoVertex->trackWeight(*recoDaughter);
00107                     matches[trackingVertex].second++;
00108                 }
00109             }
00110 
00111             recoDaughterWeight += recoVertex->trackWeight(*recoDaughter);
00112         } //loop over daughter tracks of RecoVertex
00113 
00114         std::size_t assoIndex = 0;
00115 
00116         // Loop over map between TrackingVertexes and matched RecoDaugther
00117         for (
00118             std::map<TrackingVertexRef,std::pair<double, std::size_t> >::const_iterator match = matches.begin();
00119             match != matches.end();
00120             ++match
00121         )
00122         {
00123             // Getting the TrackingVertex information
00124             TrackingVertexRef trackingVertex = match->first;
00125             double matchedDaughterWeight = match->second.first;
00126             std::size_t matchedDaughterCounter = match->second.second;
00127 
00128             // Count for only reconstructible SimDaughters           
00129             std::size_t simDaughterCounter = 0;
00130 
00131             for (
00132                 TrackingVertex::tp_iterator simDaughter = trackingVertex->daughterTracks_begin();
00133                 simDaughter != trackingVertex->daughterTracks_end();
00134                 ++simDaughter
00135             )
00136                 if ( selector_(**simDaughter) ) simDaughterCounter++;
00137 
00138             // Sanity condition in case that reconstructable condition is too tight
00139             if ( simDaughterCounter < matchedDaughterCounter )
00140                 simDaughterCounter = matchedDaughterCounter;
00141 
00142             // Condition over S2RMatchedSimRatio
00143             if ( (double)matchedDaughterCounter/simDaughterCounter < R2SMatchedSimRatio_ ) continue;
00144 
00145             double quality = (double)matchedDaughterWeight/recoDaughterWeight;
00146 
00147             // Condition over R2SMatchedRecoRatio
00148             if (quality < R2SMatchedRecoRatio_) continue;
00149 
00150             outputCollection.insert(recoVertex, std::make_pair(trackingVertex, quality));
00151 
00152             std::cout << "R2S: INDEX " << assoIndex;
00153             std::cout << " ; SimDaughterCounter = " << simDaughterCounter;
00154             std::cout << " ; RecoDaughterWeight = " << recoDaughterWeight;
00155             std::cout << " ; MatchedDaughterCounter = " << matchedDaughterCounter;
00156             std::cout << " ; MatchedDaughterWeight = " << matchedDaughterWeight;
00157             std::cout << " ; quality = " << quality << std::endl;
00158 
00159             assoIndex++;
00160         }
00161     } // Loop on RecoVertex
00162 
00163     std::cout << std::endl;
00164     std::cout << "RecoToSim OUTPUT COLLECTION: outputCollection.size() = " << outputCollection.size() << std::endl << std::endl;
00165 
00166     return outputCollection;
00167 }
00168 
00169 
00170 
00171 reco::VertexSimToRecoCollection VertexAssociatorByTracks::associateSimToReco(
00172     edm::Handle<edm::View<reco::Vertex> > & recoVertexes,
00173     edm::Handle<TrackingVertexCollection> & trackingVertexes,
00174     const edm::Event& event,
00175     reco::SimToRecoCollection & associator
00176 ) const
00177 {
00178     reco::VertexSimToRecoCollection  outputCollection;
00179 
00180     // Loop over TrackingVertexes
00181     std::map<std::size_t,std::pair<double, std::size_t> > matches;
00182 
00183     // Loop over TrackingVertexes
00184     for (std::size_t simIndex = 0; simIndex < trackingVertexes->size(); ++simIndex)
00185     {
00186         TrackingVertexRef trackingVertex(trackingVertexes, simIndex);
00187 
00188         matches.clear();
00189 
00190         std::size_t simDaughterCounter = 0;
00191 
00192         // Loop over daughter tracks of TrackingVertex
00193         for (
00194             TrackingVertex::tp_iterator simDaughter = trackingVertex->daughterTracks_begin();
00195             simDaughter != trackingVertex->daughterTracks_end();
00196             ++simDaughter
00197         )
00198         {
00199             // Select only reconstructible SimDaughters
00200             if ( !selector_(**simDaughter) ) continue;
00201 
00202             // Check for association for the given RecoDaughter
00203             if ( associator.numberOfAssociations(*simDaughter) > 0 )
00204             {
00205                 std::vector<std::pair<reco::TrackBaseRef, double> > associations = associator[*simDaughter];
00206 
00207                 // Loop over RecoTracks associated with TrackingParticle
00208                 for (
00209                     std::vector<std::pair<reco::TrackBaseRef,double> >::const_iterator association = associations.begin();
00210                     association != associations.end();
00211                     ++association
00212                 )
00213                 {
00214                     reco::TrackBaseRef recoTrack = association->first;
00215 
00216                     for (std::size_t recoIndex = 0; recoIndex < recoVertexes->size(); ++recoIndex)
00217                     {
00218                         reco::VertexBaseRef recoVertex(recoVertexes, recoIndex);
00219 
00220                         for (
00221                             reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
00222                             recoDaughter != recoVertex->tracks_end();
00223                             ++recoDaughter
00224                         )
00225                         {
00226                             // Store matched RecoDaughter to the RecoVertex
00227                             if (
00228                                 recoDaughter->id() == recoTrack.id() &&
00229                                 recoDaughter->key() == recoTrack.key()
00230                             )
00231                             {
00232                                 matches[recoIndex].first += recoVertex->trackWeight(*recoDaughter);
00233                                 matches[recoIndex].second++; 
00234                             }
00235                         }
00236                     }
00237                 } // loop a recotracks
00238             }
00239 
00240             simDaughterCounter++;
00241         } // loop a simDaughter
00242 
00243         std::size_t assoIndex = 0;
00244 
00245         // Loop over map, set score, add to outputCollection
00246         for (
00247             std::map<std::size_t,std::pair<double,std::size_t> >::const_iterator match = matches.begin();
00248             match != matches.end();
00249             ++match
00250         )
00251         {
00252             // Getting the TrackingVertex information
00253             reco::VertexBaseRef recoVertex(recoVertexes, match->first);
00254             double matchedDaughterWeight = match->second.first;
00255             std::size_t matchedDaughterCounter = match->second.second;
00256 
00257             double recoDaughterWeight = 0.;
00258 
00259             // Weighted count of those tracks with a given quality of the RecoDoughters
00260             for (
00261                 reco::Vertex::trackRef_iterator recoDaughter = recoVertex->tracks_begin();
00262                 recoDaughter != recoVertex->tracks_end();
00263                 ++recoDaughter
00264             )
00265                 if ( (*recoDaughter)->quality(trackQuality_) ) 
00266                     recoDaughterWeight += recoVertex->trackWeight(*recoDaughter);
00267 
00268             // Sanity condition in case that track quality condition is too tight
00269             if ( recoDaughterWeight < matchedDaughterWeight )
00270                 recoDaughterWeight = matchedDaughterWeight;
00271 
00272             // Condition over S2RMatchedRecoRatio
00273             if ( matchedDaughterWeight/recoDaughterWeight < S2RMatchedRecoRatio_ ) continue;
00274 
00275             double quality = (double)matchedDaughterCounter/simDaughterCounter;
00276 
00277             // Condition over S2RMatchedSimRatio
00278             if (quality < S2RMatchedSimRatio_) continue;
00279 
00280             outputCollection.insert(trackingVertex, std::make_pair(recoVertex, quality));
00281 
00282             std::cout << "R2S: INDEX " << assoIndex;
00283             std::cout << " ; SimDaughterCounter = " << simDaughterCounter;
00284             std::cout << " ; RecoDaughterWeight = " << recoDaughterWeight;
00285             std::cout << " ; MatchedDaughterCounter = " << matchedDaughterCounter;
00286             std::cout << " ; MatchedDaughterWeight = " << matchedDaughterWeight;
00287             std::cout << " ; quality = " << quality << std::endl;
00288 
00289             assoIndex++;
00290         }
00291     } // loop over TrackingVertexes
00292 
00293     std::cout << "SimToReco OUTPUT COLLECTION: outputCollection.size() = " << outputCollection.size() << std::endl << std::endl;
00294 
00295     return outputCollection;
00296 }
00297