CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
JetTracksAssociationDRVertexAssigned Class Reference

#include <JetTracksAssociationDRVertexAssigned.h>

Public Member Functions

 JetTracksAssociationDRVertexAssigned (double fDr)
 
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
 
 ~JetTracksAssociationDRVertexAssigned ()
 

Private Attributes

double mDeltaR2Threshold
 fidutial dR between track in the vertex and jet's reference direction More...
 

Detailed Description

Definition at line 14 of file JetTracksAssociationDRVertexAssigned.h.

Constructor & Destructor Documentation

JetTracksAssociationDRVertexAssigned::JetTracksAssociationDRVertexAssigned ( double  fDr)

Definition at line 15 of file JetTracksAssociationDRVertexAssigned.cc.

16 : mDeltaR2Threshold (fDr*fDr)
17 {}
double mDeltaR2Threshold
fidutial dR between track in the vertex and jet&#39;s reference direction
JetTracksAssociationDRVertexAssigned::~JetTracksAssociationDRVertexAssigned ( )
inline

Definition at line 17 of file JetTracksAssociationDRVertexAssigned.h.

17 {}

Member Function Documentation

void JetTracksAssociationDRVertexAssigned::produce ( reco::JetTracksAssociation::Container fAssociation,
const std::vector< edm::RefToBase< reco::Jet > > &  fJets,
const std::vector< reco::TrackRef > &  fTracks,
const reco::VertexCollection vertices 
) const

Definition at line 19 of file JetTracksAssociationDRVertexAssigned.cc.

References Geom::deltaR2(), reco::TrackBase::dz(), reco::TrackBase::eta(), reco::LeafCandidate::eta(), eta(), spr::find(), i, j, metsig::jet, reco::btau::jetEta, reco::btau::jetPhi, mDeltaR2Threshold, reco::TrackBase::p(), phi, reco::LeafCandidate::phi(), reco::TrackBase::phi(), reco::Jet::physicsEta(), edm::RefVector< C, T, F >::push_back(), findQualityFiles::rr, reco::JetTracksAssociation::setValue(), and lumiQTWidget::t.

Referenced by JetTracksAssociatorAtVertex::produce().

23 {
24  // cache tracks kinematics
25  std::vector <math::RhoEtaPhiVector> trackP3s;
26  std::map <int,double> trackvert;
27 
28  // std::cout<<" Number of vertices "<<vertices.size()<<std::endl;
29 
30  trackP3s.reserve (fTracks.size());
31  for (unsigned i = 0; i < fTracks.size(); ++i) {
32  const reco::Track* track = &*(fTracks[i]);
33  trackP3s.push_back (math::RhoEtaPhiVector (track->p(),track->eta(), track->phi()));
34 
35  // OK: Look for the tracks not associated with vertices
36 
37  const reco::TrackBaseRef ttr1(fTracks[i]);
38 
39 
40  int trackhasvert = -1;
41  for( reco::VertexCollection::const_iterator iv = vertices.begin(); iv != vertices.end(); iv++) {
42  std::vector<reco::TrackBaseRef>::const_iterator rr =
43  find((*iv).tracks_begin(),
44  (*iv).tracks_end(),
45  ttr1);
46  if( rr != (*iv).tracks_end() ) {
47  trackhasvert = 1;
48  trackvert[i] = (*iv).position().z();
49  // std::cout<<" Z "<<i<<" "<<trackhasvert<<" "<<(*iv).position().z()<<std::endl;
50  break;
51  }
52  } // all vertices
53  if(trackhasvert < 0) {
54  // Take impact parameter of the track as vertex position
55  math::XYZPoint ppt(0.,0.,0.);
56  trackvert[i] = track->dz(ppt);
57  // std::cout<<" Z "<<i<<" "<<trackhasvert<<" "<<track->dz(ppt)<<" "<<track->vz()<<" "<<track->vx()<<" "<<
58  // track->vy()<<" "<<track->pz()<<std::endl;
59  }
60  // OK
61  } // tracks
62 
63  for (unsigned j = 0; j < fJets.size(); ++j) {
64 
65  reco::TrackRefVector assoTracks;
66  const reco::Jet* jet = &*(fJets[j]);
67  double jetEta = jet->eta();
68  double jetPhi = jet->phi();
69  double neweta = 0;
70  for (unsigned t = 0; t < fTracks.size(); ++t) {
71 
72  std::map<int, double>::iterator cur = trackvert.find(t);
73  if(cur != trackvert.end()) {
74  neweta = jet->physicsEta((*cur).second,jetEta);
75  } else {
76  neweta = jetEta;
77  //std::cout<<" Lost track - not in map "<<std::endl;
78  }
79 
80  //std::cout<<" Old eta-new eta "<<(*cur).second<<" "<<jetEta<<" "<<neweta<<" Track "<<t<<" "<<trackP3s[t].eta()<<std::endl;
81 
82  double dR2 = deltaR2 (neweta, jetPhi, trackP3s[t].eta(), trackP3s[t].phi());
83  if (dR2 < mDeltaR2Threshold) assoTracks.push_back (fTracks[t]);
84  }
85 
86  reco::JetTracksAssociation::setValue (fAssociation, fJets[j], assoTracks);
87  }
88 
89 }
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
int i
Definition: DBlmapReader.cc:9
Base class for all types of Jets.
Definition: Jet.h:21
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
double mDeltaR2Threshold
fidutial dR between track in the vertex and jet&#39;s reference direction
T eta() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
RhoEtaPhiVectorD RhoEtaPhiVector
spatial vector with cylindrical internal representation using pseudorapidity
Definition: Vector3D.h:33
virtual double eta() const
momentum pseudorapidity
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
bool setValue(Container &, const reco::JetBaseRef &, reco::TrackRefVector)
associate jet with value. Returns false and associate nothing if jet is already associated ...
int j
Definition: DBlmapReader.cc:9
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
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:127
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
virtual double phi() const
momentum azimuthal angle
static float physicsEta(float fZVertex, float fDetectorEta)
static function to convert detector eta to physics eta
Definition: Jet.cc:319
Definition: DDAxes.h:10

Member Data Documentation

double JetTracksAssociationDRVertexAssigned::mDeltaR2Threshold
private

fidutial dR between track in the vertex and jet's reference direction

Definition at line 25 of file JetTracksAssociationDRVertexAssigned.h.

Referenced by produce().