CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackGenAssociatorByChi2.h
Go to the documentation of this file.
1 #ifndef TrackGenAssociatorByChi2_h
2 #define TrackGenAssociatorByChi2_h
3 
27 
28 #include<map>
29 
30 //Note that the Association Map is filled with -ch2 and not chi2 because it is ordered using std::greater:
31 //the track with the lowest association chi2 will be the first in the output map.
32 
33 namespace reco{
40 }
41 
42 
44 
45  public:
46  typedef std::map<double, SimTrack> Chi2SimMap;
47  typedef std::pair< reco::Track, Chi2SimMap> RecoToSimPair;
48  typedef std::vector< RecoToSimPair > RecoToSimPairAssociation;
49 
52  chi2cut(conf.getParameter<double>("chi2cut")),
53  onlyDiagonal(conf.getParameter<bool>("onlyDiagonal")),
54  bsSrc(conf.getParameter<edm::InputTag>("beamSpot")) {
55  theMF=mF;
56  if (onlyDiagonal)
57  edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms = 0 ---- " << "\n";
58  else
59  edm::LogInfo("TrackAssociator") << " ---- Using Off Diagonal Covariance Terms != 0 ---- " << "\n";
60  }
61 
63  TrackGenAssociatorByChi2(const edm::ESHandle<MagneticField> mF, double chi2Cut, bool onlyDiag, const edm::InputTag& beamspotSrc){
64  chi2cut=chi2Cut;
65  onlyDiagonal=onlyDiag;
66  theMF=mF;
67  bsSrc = beamspotSrc;
68  }
69 
73  const edm::Event * event = 0,
74  const edm::EventSetup * setup = 0 ) const ;
78  const edm::Event * event = 0,
79  const edm::EventSetup * setup = 0 ) const ;
80 
84  const edm::Event * event = 0,
85  const edm::EventSetup * setup = 0) const {
87  for (unsigned int j=0; j<tCH->size();j++)
89 
91  for (unsigned int j=0; j<tPCH->size();j++)
93 
94  return associateRecoToGen(tc,tpc,event,setup);
95  }
96 
100  const edm::Event * event = 0,
101  const edm::EventSetup * setup = 0) const {
103  for (unsigned int j=0; j<tCH->size();j++)
105 
107  for (unsigned int j=0; j<tPCH->size();j++)
109 
110  return associateGenToReco(tc,tpc,event,setup);
111  }
112 
113 
114  private:
115 
117  double getChi2(const reco::TrackBase::ParameterVector& rParameters,
118  const reco::TrackBase::CovarianceMatrix& recoTrackCovMatrix,
119  const Basic3DVector<double>& momAtVtx,
120  const Basic3DVector<double>& vert,
121  int charge,
122  const reco::BeamSpot&) const;
123 
125  double associateRecoToSim(reco::TrackCollection::const_iterator,
126  TrackingParticleCollection::const_iterator,
127  const reco::BeamSpot&) const;
128 
129 
131  double chi2cut;
134 };
135 
136 #endif
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
reco::GenToRecoCollection associateGenToReco(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
Association Sim To Reco with Collections (Gen Particle version)
std::pair< reco::Track, Chi2SimMap > RecoToSimPair
virtual reco::GenToRecoCollection associateGenToReco(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< reco::GenParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
compare reco to sim the handle of reco::Track and GenParticle collections
ProductID id() const
Definition: HandleBase.cc:15
edm::AssociationMap< edm::OneToManyWithQualityGeneric< reco::GenParticleCollection, edm::View< reco::Track >, double > > GenToRecoCollection
edm::ESHandle< MagneticField > theMF
std::map< double, SimTrack > Chi2SimMap
math::Vector< dimension >::type ParameterVector
parameter vector
Definition: TrackBase.h:74
virtual reco::RecoToGenCollection associateRecoToGen(edm::Handle< edm::View< reco::Track > > &tCH, edm::Handle< reco::GenParticleCollection > &tPCH, const edm::Event *event=0, const edm::EventSetup *setup=0) const
compare reco to sim the handle of reco::Track and GenParticle collections
int j
Definition: DBlmapReader.cc:9
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
tuple conf
Definition: dbtoconf.py:185
reco::RecoToGenCollection associateRecoToGen(const edm::RefToBaseVector< reco::Track > &, const edm::RefVector< reco::GenParticleCollection > &, const edm::Event *event=0, const edm::EventSetup *setup=0) const
Association Sim To Reco with Collections (Gen Particle version)
TrackGenAssociatorByChi2(const edm::ESHandle< MagneticField > mF, const edm::ParameterSet &conf)
Constructor with PSet.
void push_back(const RefToBase< T > &)
double associateRecoToSim(reco::TrackCollection::const_iterator, TrackingParticleCollection::const_iterator, const reco::BeamSpot &) const
compare reco::TrackCollection and TrackingParticleCollection iterators: returns the chi2 ...
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
edm::AssociationMap< edm::OneToManyWithQualityGeneric< edm::View< reco::Track >, reco::GenParticleCollection, double > > RecoToGenCollection
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
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
TrackGenAssociatorByChi2(const edm::ESHandle< MagneticField > mF, double chi2Cut, bool onlyDiag, const edm::InputTag &beamspotSrc)
Constructor with magnetic field, double, bool and InputTag.
std::vector< RecoToSimPair > RecoToSimPairAssociation
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77