CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackAssociatorByChi2.h
Go to the documentation of this file.
1 #ifndef TrackAssociatorByChi2_h
2 #define TrackAssociatorByChi2_h
3 
22 
23 #include<map>
24 
25 //Note that the Association Map is filled with -ch2 and not chi2 because it is ordered using std::greater:
26 //the track with the lowest association chi2 will be the first in the output map.
27 
29 
30  public:
31  typedef std::map<double, SimTrack> Chi2SimMap;
32  typedef std::pair< reco::Track, Chi2SimMap> RecoToSimPair;
33  typedef std::vector< RecoToSimPair > RecoToSimPairAssociation;
34 
37  chi2cut(conf.getParameter<double>("chi2cut")),
38  onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
39  bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
40  theMF=mF;
41  if (onlyDiagonal)
42  edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " << "\n";
43  else
44  edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " << "\n";
45  }
46 
48  TrackAssociatorByChi2(const edm::ESHandle<MagneticField> mF, double chi2Cut, bool onlyDiag, const edm::InputTag& beamspotSrc){
49  chi2cut=chi2Cut;
50  onlyDiagonal=onlyDiag;
51  theMF=mF;
52  bsSrc = beamspotSrc;
53  }
54 
57 
58 
60  virtual
63  const edm::Event * event = 0,
64  const edm::EventSetup * setup = 0 ) const override;
66  virtual
69  const edm::Event * event = 0,
70  const edm::EventSetup * setup = 0 ) const override;
71 
73  virtual
76  const edm::Event * event = 0,
77  const edm::EventSetup * setup = 0) const override {
79  }
80 
82  virtual
85  const edm::Event * event = 0,
86  const edm::EventSetup * setup = 0) const override {
88  }
89 
90 
91  private:
92 
94  std::pair<bool,reco::TrackBase::ParameterVector> parametersAtClosestApproach(const Basic3DVector<double>&,// vertex
95  const Basic3DVector<double>&,// momAtVtx
96  float,// charge
97  const reco::BeamSpot&) const;//beam spot
98 
100  double compareTracksParam(reco::TrackCollection::const_iterator,
101  edm::SimTrackContainer::const_iterator,
102  const math::XYZTLorentzVectorD&,
103  const GlobalVector&,
105  const reco::BeamSpot&) const;
106 
109  const edm::SimTrackContainer&,
111  const reco::BeamSpot&) const;
112 
114  double getChi2(const reco::TrackBase::ParameterVector& rParameters,
115  const reco::TrackBase::CovarianceMatrix& recoTrackCovMatrix,
116  const Basic3DVector<double>& momAtVtx,
117  const Basic3DVector<double>& vert,
118  int charge,
119  const reco::BeamSpot&) const;
120 
122  double associateRecoToSim(reco::TrackCollection::const_iterator,
123  TrackingParticleCollection::const_iterator,
124  const reco::BeamSpot&) const;
125 
126 
128  double chi2cut;
131 };
132 
133 #endif
std::pair< bool, reco::TrackBase::ParameterVector > parametersAtClosestApproach(const Basic3DVector< double > &, const Basic3DVector< double > &, float, const reco::BeamSpot &) const
propagate the track parameters of TrackinParticle from production vertex to the point of closest appr...
virtual reco::SimToRecoCollection associateSimToReco(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const override
Association Sim To Reco with Collections.
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
TrackAssociatorByChi2(const edm::ESHandle< MagneticField > mF, const edm::ParameterSet &conf)
Constructor with PSet.
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
double compareTracksParam(reco::TrackCollection::const_iterator, edm::SimTrackContainer::const_iterator, const math::XYZTLorentzVectorD &, const GlobalVector &, const reco::TrackBase::CovarianceMatrix &, const reco::BeamSpot &) const
compare reco::TrackCollection and edm::SimTrackContainer iterators: returns the chi2 ...
std::pair< reco::Track, Chi2SimMap > RecoToSimPair
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:73
TrackAssociatorByChi2(const edm::ESHandle< MagneticField > mF, double chi2Cut, bool onlyDiag, const edm::InputTag &beamspotSrc)
Constructor with magnetic field, double, bool and InputTag.
~TrackAssociatorByChi2()
Destructor.
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
edm::ESHandle< MagneticField > theMF
tuple conf
Definition: dbtoconf.py:185
virtual reco::SimToRecoCollection associateSimToReco(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event, const edm::EventSetup *setup) const
compare reco to sim the handle of reco::Track and TrackingParticle collections
virtual reco::RecoToSimCollection associateRecoToSim(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< TrackingParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const override
Association Reco To Sim with Collections.
std::vector< SimVertex > SimVertexContainer
double getChi2(const reco::TrackBase::ParameterVector &rParameters, const reco::TrackBase::CovarianceMatrix &recoTrackCovMatrix, const Basic3DVector< double > &momAtVtx, const Basic3DVector< double > &vert, int charge, const reco::BeamSpot &) const
basic method where chi2 is computed
std::map< double, SimTrack > Chi2SimMap
std::vector< RecoToSimPair > RecoToSimPairAssociation
virtual reco::RecoToSimCollection associateRecoToSim(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
virtual reco::SimToRecoCollection associateSimToReco(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< TrackingParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const override
compare reco to sim the handle of reco::Track and TrackingParticle collections
std::vector< SimTrack > SimTrackContainer
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:76