CMS 3D CMS Logo

JetTracksAssociationDRVertexAssigned.cc
Go to the documentation of this file.
1 // Associate jets with tracks by simple "dR" criteria
2 // Fedor Ratnikov (UMd), Aug. 28, 2007
3 
5 
8 
13 
15 : mDeltaR2Threshold (fDr*fDr)
16 {}
17 
19  const std::vector <edm::RefToBase<reco::Jet> >& fJets,
20  const std::vector <reco::TrackRef>& fTracks,
21  const reco::VertexCollection& vertices) const
22 {
23  // cache tracks kinematics
24  std::vector <math::RhoEtaPhiVector> trackP3s;
25  std::map <int,double> trackvert;
26 
27  // std::cout<<" Number of vertices "<<vertices.size()<<std::endl;
28 
29  trackP3s.reserve (fTracks.size());
30  for (unsigned i = 0; i < fTracks.size(); ++i) {
31  const reco::Track* track = &*(fTracks[i]);
32  trackP3s.push_back (math::RhoEtaPhiVector (track->p(),track->eta(), track->phi()));
33 
34  // OK: Look for the tracks not associated with vertices
35 
36  const reco::TrackBaseRef ttr1(fTracks[i]);
37 
38 
39  int trackhasvert = -1;
40  for( reco::VertexCollection::const_iterator iv = vertices.begin(); iv != vertices.end(); iv++) {
41  std::vector<reco::TrackBaseRef>::const_iterator rr =
42  find((*iv).tracks_begin(),
43  (*iv).tracks_end(),
44  ttr1);
45  if( rr != (*iv).tracks_end() ) {
46  trackhasvert = 1;
47  trackvert[i] = (*iv).position().z();
48  // std::cout<<" Z "<<i<<" "<<trackhasvert<<" "<<(*iv).position().z()<<std::endl;
49  break;
50  }
51  } // all vertices
52  if(trackhasvert < 0) {
53  // Take impact parameter of the track as vertex position
54  math::XYZPoint ppt(0.,0.,0.);
55  trackvert[i] = track->dz(ppt);
56  // std::cout<<" Z "<<i<<" "<<trackhasvert<<" "<<track->dz(ppt)<<" "<<track->vz()<<" "<<track->vx()<<" "<<
57  // track->vy()<<" "<<track->pz()<<std::endl;
58  }
59  // OK
60  } // tracks
61 
62  for (unsigned j = 0; j < fJets.size(); ++j) {
63 
64  reco::TrackRefVector assoTracks;
65  const reco::Jet* jet = &*(fJets[j]);
66  double jetEta = jet->eta();
67  double jetPhi = jet->phi();
68  double neweta = 0;
69  for (unsigned t = 0; t < fTracks.size(); ++t) {
70 
71  std::map<int, double>::iterator cur = trackvert.find(t);
72  if(cur != trackvert.end()) {
73  neweta = jet->physicsEta((*cur).second,jetEta);
74  } else {
75  neweta = jetEta;
76  //std::cout<<" Lost track - not in map "<<std::endl;
77  }
78 
79  //std::cout<<" Old eta-new eta "<<(*cur).second<<" "<<jetEta<<" "<<neweta<<" Track "<<t<<" "<<trackP3s[t].eta()<<std::endl;
80 
81  double dR2 = deltaR2 (neweta, jetPhi, trackP3s[t].eta(), trackP3s[t].phi());
82  if (dR2 < mDeltaR2Threshold) assoTracks.push_back (fTracks[t]);
83  }
84 
85  reco::JetTracksAssociation::setValue (fAssociation, fJets[j], assoTracks);
86  }
87 
88 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:654
double eta() const final
momentum pseudorapidity
Base class for all types of Jets.
Definition: Jet.h:20
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:684
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
double mDeltaR2Threshold
fidutial dR between track in the vertex and jet&#39;s reference direction
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
RhoEtaPhiVectorD RhoEtaPhiVector
spatial vector with cylindrical internal representation using pseudorapidity
Definition: Vector3D.h:32
void produce(reco::JetTracksAssociation::Container *fAssociation, const std::vector< edm::RefToBase< reco::Jet > > &fJets, const std::vector< reco::TrackRef > &fTracks, const reco::VertexCollection &vertices) const
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:690
bool setValue(Container &, const reco::JetBaseRef &, reco::TrackRefVector)
associate jet with value. Returns false and associate nothing if jet is already associated ...
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:648
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:69
double phi() const final
momentum azimuthal angle
static float physicsEta(float fZVertex, float fDetectorEta)
static function to convert detector eta to physics eta