CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
IsolatedPixelTrackCandidate.h
Go to the documentation of this file.
1 #ifndef HcalIsolatedTrack_IsolatedPixelTrackCandidate_h
2 #define HcalIsolatedTrack_IsolatedPixelTrackCandidate_h
3 
12 
14 
15 #include <vector>
16 #include <map>
17 #include <utility>
18 
19 namespace reco {
20 
22 
23  public:
24 
29  {
30  enIn_=-1;
31  enOut_=-1;
32  nhitIn_=-1;
33  nhitOut_=-1;
34  maxPtPxl_=-1;
35  sumPtPxl_=-1;
36  }
38  IsolatedPixelTrackCandidate(const reco::TrackRef & tr, const l1extra::L1JetParticleRef & tauRef, double max, double sum):
39  RecoCandidate( 0, LorentzVector((tr.get()->px()),(tr.get())->py(),(tr.get())->pz(),(tr.get())->p()) ),
40  track_(tr), l1tauJet_(tauRef), maxPtPxl_(max), sumPtPxl_(sum)
41  {
42  enIn_=-1;
43  enOut_=-1;
44  nhitIn_=-1;
45  nhitOut_=-1;
46  }
47 
49  IsolatedPixelTrackCandidate(const l1extra::L1JetParticleRef & tauRef, double enIn, double enOut, int nhitIn, int nhitOut):
50  RecoCandidate( 0, LorentzVector(tauRef->px(),tauRef->py(),tauRef->pz(),tauRef->p()) ),
51  l1tauJet_(tauRef), enIn_(enIn), enOut_(enOut), nhitIn_(nhitIn), nhitOut_(nhitOut)
52  {
53  maxPtPxl_=-1;
54  sumPtPxl_=-1;
55  }
56 
57 
58 
59 
62 
64  virtual IsolatedPixelTrackCandidate * clone() const;
65 
67  virtual reco::TrackRef track() const;
68  void setTrack( const reco::TrackRef & tr ) { track_ = tr; }
69 
71  double maxPtPxl() const {return maxPtPxl_;}
72  void SetMaxPtPxl(double mptpxl) {maxPtPxl_=mptpxl;}
73 
75  double sumPtPxl() const {return sumPtPxl_;}
76  void SetSumPtPxl(double sumptpxl) {sumPtPxl_=sumptpxl;}
77 
79  virtual l1extra::L1JetParticleRef l1tau() const;
80  void setL1TauJet( const l1extra::L1JetParticleRef & tauRef ) { l1tauJet_ = tauRef; }
81 
83  double energyIn() const {return enIn_; }
84  void SetEnergyIn(double a) {enIn_=a;}
85 
87  double energyOut() const {return enOut_;}
88  void SetEnergyOut(double a) {enOut_=a;}
89 
91  int nHitIn() const {return nhitIn_;}
92  void SetNHitIn(int a) {nhitIn_=a;}
93 
95  int nHitOut() const {return nhitOut_;}
96  void SetNHitOut(int a) {nhitOut_=a;}
97 
99  std::pair<int,int> towerIndex() const;
100 
101  private:
103  virtual bool overlap( const Candidate & ) const;
109  double maxPtPxl_;
111  double sumPtPxl_;
113  double enIn_;
115  double enOut_;
117  int nhitIn_;
119  int nhitOut_;
120 
121  };
122 
123 
124 }
125 
126 #endif
double enOut_
energy in outer cone around L1 tau jet
virtual reco::TrackRef track() const
refrence to a Track
reco::TrackRef track_
reference to a Track
virtual double p() const GCC11_FINAL
magnitude of momentum vector
double energyIn() const
ECAL energy in the inner cone around tau jet.
IsolatedPixelTrackCandidate(const l1extra::L1JetParticleRef &tauRef, double enIn, double enOut, int nhitIn, int nhitOut)
constructor from tau jet
virtual double pz() const GCC11_FINAL
z coordinate of momentum vector
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
double sumPtPxl() const
Pt sum of other pixel tracks in the cone around the candidate.
double sumPtPxl_
Pt sum of other pixel tracks in the cone around the candidate.
const T & max(const T &a, const T &b)
virtual IsolatedPixelTrackCandidate * clone() const
returns a clone of the candidate
int nHitOut() const
number of ECAL hits in the outer cone around tau jet
void setL1TauJet(const l1extra::L1JetParticleRef &tauRef)
int nHitIn() const
number of ECAL hits in the inner cone around tau jet
T get() const
get a component
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
IsolatedPixelTrackCandidate(const LorentzVector &v)
constructor from LorentzVector
double maxPtPxl_
highest Pt of other pixel tracks in the cone around the candidate
int nhitIn_
number of hits in inner cone
virtual bool overlap(const Candidate &) const
check overlap with another candidate
virtual l1extra::L1JetParticleRef l1tau() const
get reference to L1 tau jet
int nhitOut_
number of hits in inner cone
double energyOut() const
ECAL energy in the outer cone around tau jet.
double enIn_
energy in inner cone around L1 tau jet
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
double maxPtPxl() const
highest Pt of other pixel tracks in the cone around the candidate
double a
Definition: hdecay.h:121
l1extra::L1JetParticleRef l1tauJet_
reference to a L1 tau jet
IsolatedPixelTrackCandidate(const reco::TrackRef &tr, const l1extra::L1JetParticleRef &tauRef, double max, double sum)
constructor from a track
void setTrack(const reco::TrackRef &tr)
std::pair< int, int > towerIndex() const
get index of tower which track is hitting