00001 #ifndef TrackAssociatorByChi2_h
00002 #define TrackAssociatorByChi2_h
00003
00012 #include "SimTracker/TrackAssociation/interface/TrackAssociatorBase.h"
00013 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "MagneticField/Engine/interface/MagneticField.h"
00016 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "DataFormats/Math/interface/LorentzVector.h"
00022 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00023 #include "FWCore/Utilities/interface/InputTag.h"
00024 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00025 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00026
00027 #include<map>
00028
00029
00030
00031
00032 namespace reco{
00033 typedef edm::AssociationMap<edm::OneToManyWithQualityGeneric
00034 <reco::GenParticleCollection, edm::View<reco::Track>, double> >
00035 GenToRecoCollection;
00036 typedef edm::AssociationMap<edm::OneToManyWithQualityGeneric
00037 <edm::View<reco::Track>, reco::GenParticleCollection, double> >
00038 RecoToGenCollection;
00039 }
00040
00041
00042 class TrackAssociatorByChi2 : public TrackAssociatorBase {
00043
00044 public:
00045 typedef std::map<double, SimTrack> Chi2SimMap;
00046 typedef std::pair< reco::Track, Chi2SimMap> RecoToSimPair;
00047 typedef std::vector< RecoToSimPair > RecoToSimPairAssociation;
00048
00050 TrackAssociatorByChi2(const edm::ESHandle<MagneticField> mF, edm::ParameterSet conf):
00051 chi2cut(conf.getParameter<double>("chi2cut")),
00052 onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
00053 bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
00054 theMF=mF;
00055 if (onlyDiagonal)
00056 edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " << "\n";
00057 else
00058 edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " << "\n";
00059 }
00060
00062 TrackAssociatorByChi2(const edm::ESHandle<MagneticField> mF, double chi2Cut, bool onlyDiag, edm::InputTag beamspotSrc){
00063 chi2cut=chi2Cut;
00064 onlyDiagonal=onlyDiag;
00065 theMF=mF;
00066 bsSrc = beamspotSrc;
00067 }
00068
00070 ~TrackAssociatorByChi2(){}
00071
00073 double compareTracksParam(reco::TrackCollection::const_iterator,
00074 edm::SimTrackContainer::const_iterator,
00075 const math::XYZTLorentzVectorD,
00076 GlobalVector,
00077 reco::TrackBase::CovarianceMatrix,
00078 const reco::BeamSpot&) const;
00079
00081 RecoToSimPairAssociation compareTracksParam(const reco::TrackCollection&,
00082 const edm::SimTrackContainer&,
00083 const edm::SimVertexContainer&,
00084 const reco::BeamSpot&) const;
00085
00087 double getChi2(reco::TrackBase::ParameterVector& rParameters,
00088 reco::TrackBase::CovarianceMatrix& recoTrackCovMatrix,
00089 Basic3DVector<double>& momAtVtx,
00090 Basic3DVector<double>& vert,
00091 int& charge,
00092 const reco::BeamSpot&) const;
00093
00095 double associateRecoToSim(reco::TrackCollection::const_iterator,
00096 TrackingParticleCollection::const_iterator,
00097 const reco::BeamSpot&) const;
00098
00100 std::pair<bool,reco::TrackBase::ParameterVector> parametersAtClosestApproach(Basic3DVector<double>,
00101 Basic3DVector<double>,
00102 float,
00103 const reco::BeamSpot&) const;
00105 reco::RecoToSimCollection associateRecoToSim(const edm::RefToBaseVector<reco::Track>&,
00106 const edm::RefVector<TrackingParticleCollection>&,
00107 const edm::Event * event = 0,
00108 const edm::EventSetup * setup = 0 ) const ;
00110 reco::SimToRecoCollection associateSimToReco(const edm::RefToBaseVector<reco::Track>&,
00111 const edm::RefVector<TrackingParticleCollection>&,
00112 const edm::Event * event = 0,
00113 const edm::EventSetup * setup = 0 ) const ;
00114
00116 reco::RecoToSimCollection associateRecoToSim(edm::Handle<edm::View<reco::Track> >& tCH,
00117 edm::Handle<TrackingParticleCollection>& tPCH,
00118 const edm::Event * event = 0,
00119 const edm::EventSetup * setup = 0) const {
00120 return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event,setup);
00121 }
00122
00124 reco::SimToRecoCollection associateSimToReco(edm::Handle<edm::View<reco::Track> >& tCH,
00125 edm::Handle<TrackingParticleCollection>& tPCH,
00126 const edm::Event * event = 0,
00127 const edm::EventSetup * setup = 0) const {
00128 return TrackAssociatorBase::associateSimToReco(tCH,tPCH,event,setup);
00129 }
00130
00132 reco::RecoToGenCollection associateRecoToGen(const edm::RefToBaseVector<reco::Track>&,
00133 const edm::RefVector<reco::GenParticleCollection>&,
00134 const edm::Event * event = 0,
00135 const edm::EventSetup * setup = 0 ) const ;
00137 reco::GenToRecoCollection associateGenToReco(const edm::RefToBaseVector<reco::Track>&,
00138 const edm::RefVector<reco::GenParticleCollection>&,
00139 const edm::Event * event = 0,
00140 const edm::EventSetup * setup = 0 ) const ;
00141
00143 virtual reco::RecoToGenCollection associateRecoToGen(edm::Handle<edm::View<reco::Track> >& tCH,
00144 edm::Handle<reco::GenParticleCollection>& tPCH,
00145 const edm::Event * event = 0,
00146 const edm::EventSetup * setup = 0) const {
00147 edm::RefToBaseVector<reco::Track> tc(tCH);
00148 for (unsigned int j=0; j<tCH->size();j++)
00149 tc.push_back(edm::RefToBase<reco::Track>(tCH,j));
00150
00151 edm::RefVector<reco::GenParticleCollection> tpc(tPCH.id());
00152 for (unsigned int j=0; j<tPCH->size();j++)
00153 tpc.push_back(edm::Ref<reco::GenParticleCollection>(tPCH,j));
00154
00155 return associateRecoToGen(tc,tpc,event,setup);
00156 }
00157
00159 virtual reco::GenToRecoCollection associateGenToReco(edm::Handle<edm::View<reco::Track> >& tCH,
00160 edm::Handle<reco::GenParticleCollection>& tPCH,
00161 const edm::Event * event = 0,
00162 const edm::EventSetup * setup = 0) const {
00163 edm::RefToBaseVector<reco::Track> tc(tCH);
00164 for (unsigned int j=0; j<tCH->size();j++)
00165 tc.push_back(edm::RefToBase<reco::Track>(tCH,j));
00166
00167 edm::RefVector<reco::GenParticleCollection> tpc(tPCH.id());
00168 for (unsigned int j=0; j<tPCH->size();j++)
00169 tpc.push_back(edm::Ref<reco::GenParticleCollection>(tPCH,j));
00170
00171 return associateGenToReco(tc,tpc,event,setup);
00172 }
00173
00174
00175 private:
00176 edm::ESHandle<MagneticField> theMF;
00177 double chi2cut;
00178 bool onlyDiagonal;
00179 edm::InputTag bsSrc;
00180 };
00181
00182 #endif