CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/SimTracker/TrackAssociation/interface/TrackAssociatorByChi2.h

Go to the documentation of this file.
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 
00025 #include<map>
00026 
00027 //Note that the Association Map is filled with -ch2 and not chi2 because it is ordered using std::greater:
00028 //the track with the lowest association chi2 will be the first in the output map.
00029 
00030 class TrackAssociatorByChi2 : public TrackAssociatorBase {
00031 
00032  public:
00033   typedef std::map<double,  SimTrack> Chi2SimMap;
00034   typedef std::pair< reco::Track, Chi2SimMap> RecoToSimPair;
00035   typedef std::vector< RecoToSimPair > RecoToSimPairAssociation;
00036 
00038   TrackAssociatorByChi2(const edm::ESHandle<MagneticField> mF, edm::ParameterSet conf):
00039     chi2cut(conf.getParameter<double>("chi2cut")),
00040     onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
00041     bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
00042     theMF=mF;  
00043     if (onlyDiagonal)
00044       edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " <<  "\n";
00045     else 
00046       edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " <<  "\n";
00047   }
00048 
00050   TrackAssociatorByChi2(const edm::ESHandle<MagneticField> mF, double chi2Cut, bool onlyDiag, edm::InputTag beamspotSrc){
00051     chi2cut=chi2Cut;
00052     onlyDiagonal=onlyDiag;
00053     theMF=mF;  
00054     bsSrc = beamspotSrc;
00055   }
00056 
00058   ~TrackAssociatorByChi2(){}
00059 
00061   double compareTracksParam(reco::TrackCollection::const_iterator, 
00062                             edm::SimTrackContainer::const_iterator, 
00063                             const math::XYZTLorentzVectorD, 
00064                             GlobalVector,
00065                             reco::TrackBase::CovarianceMatrix,
00066                             const reco::BeamSpot&) const;
00067 
00069   RecoToSimPairAssociation compareTracksParam(const reco::TrackCollection&, 
00070                                               const edm::SimTrackContainer&, 
00071                                               const edm::SimVertexContainer&,
00072                                               const reco::BeamSpot&) const;
00073 
00075   double associateRecoToSim(reco::TrackCollection::const_iterator,
00076                             TrackingParticleCollection::const_iterator,
00077                             const reco::BeamSpot&) const;
00078 
00080   std::pair<bool,reco::TrackBase::ParameterVector> parametersAtClosestApproach(Basic3DVector<double>,// vertex
00081                                                                                Basic3DVector<double>,// momAtVtx
00082                                                                                float,// charge
00083                                                                                const reco::BeamSpot&) const;//beam spot
00085   reco::RecoToSimCollection associateRecoToSim(const edm::RefToBaseVector<reco::Track>&,
00086                                                const edm::RefVector<TrackingParticleCollection>&,
00087                                                const edm::Event * event = 0,
00088                                                const edm::EventSetup * setup = 0 ) const ;
00090   reco::SimToRecoCollection associateSimToReco(const edm::RefToBaseVector<reco::Track>&,
00091                                                const edm::RefVector<TrackingParticleCollection>&,
00092                                                const edm::Event * event = 0,
00093                                                const edm::EventSetup * setup = 0 ) const ;
00094   
00096   reco::RecoToSimCollection associateRecoToSim(edm::Handle<edm::View<reco::Track> >& tCH, 
00097                                                edm::Handle<TrackingParticleCollection>& tPCH, 
00098                                                const edm::Event * event = 0,
00099                                                const edm::EventSetup * setup = 0) const {
00100     return TrackAssociatorBase::associateRecoToSim(tCH,tPCH,event,setup);
00101   }
00102   
00104   reco::SimToRecoCollection associateSimToReco(edm::Handle<edm::View<reco::Track> >& tCH, 
00105                                                edm::Handle<TrackingParticleCollection>& tPCH,
00106                                                const edm::Event * event = 0,
00107                                                const edm::EventSetup * setup = 0) const {
00108     return TrackAssociatorBase::associateSimToReco(tCH,tPCH,event,setup);
00109   }  
00110  private:
00111   edm::ESHandle<MagneticField> theMF;
00112   double chi2cut;
00113   bool onlyDiagonal;
00114   edm::InputTag bsSrc;
00115 };
00116 
00117 #endif