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