CMS 3D CMS Logo

IsolatedPixelTrackCandidate.h
Go to the documentation of this file.
1 #ifndef HcalIsolatedTrack_IsolatedPixelTrackCandidate_h
2 #define HcalIsolatedTrack_IsolatedPixelTrackCandidate_h
3 
14 
16 
17 #include <vector>
18 #include <map>
19 #include <utility>
20 
21 namespace reco {
22 
24 
25  public:
26 
27  // default constructor
29  enIn_=-1;
30  enOut_=-1;
31  nhitIn_=-1;
32  nhitOut_=-1;
33  maxPtPxl_=-1;
34  sumPtPxl_=-1;
35  etaEcal_=0;
36  phiEcal_=0;
37  etaPhiEcal_=false;
38  }
41  enIn_=-1;
42  enOut_=-1;
43  nhitIn_=-1;
44  nhitOut_=-1;
45  maxPtPxl_=-1;
46  sumPtPxl_=-1;
47  etaEcal_=0;
48  phiEcal_=0;
49  etaPhiEcal_=false;
50  }
52  IsolatedPixelTrackCandidate(const reco::TrackRef & tr, const l1extra::L1JetParticleRef & tauRef, double max, double sum):
53  RecoCandidate( 0, LorentzVector((tr.get()->px()),(tr.get())->py(),(tr.get())->pz(),(tr.get())->p()) ),
54  track_(tr), l1tauJet_(tauRef), maxPtPxl_(max), sumPtPxl_(sum) {
55  enIn_=-1;
56  enOut_=-1;
57  nhitIn_=-1;
58  nhitOut_=-1;
59  etaEcal_=0;
60  phiEcal_=0;
61  etaPhiEcal_=false;
62  }
63  // constructor from a track using l1t
64  IsolatedPixelTrackCandidate(const reco::TrackRef & tr, const l1t::TauRef & tauRef, double max, double sum):
65  RecoCandidate( 0, LorentzVector((tr.get()->px()),(tr.get())->py(),(tr.get())->pz(),(tr.get())->p()) ),
66  track_(tr), l1ttauJet_(tauRef), maxPtPxl_(max), sumPtPxl_(sum) {
67  enIn_=-1;
68  enOut_=-1;
69  nhitIn_=-1;
70  nhitOut_=-1;
71  etaEcal_=0;
72  phiEcal_=0;
73  etaPhiEcal_=false;
74  }
75 
77  IsolatedPixelTrackCandidate(const l1extra::L1JetParticleRef & tauRef, double enIn, double enOut, int nhitIn, int nhitOut):
78  RecoCandidate( 0, LorentzVector(tauRef->px(),tauRef->py(),tauRef->pz(),tauRef->p()) ),
79  l1tauJet_(tauRef), enIn_(enIn), enOut_(enOut), nhitIn_(nhitIn), nhitOut_(nhitOut) {
80  maxPtPxl_=-1;
81  sumPtPxl_=-1;
82  etaEcal_=0;
83  phiEcal_=0;
84  etaPhiEcal_=false;
85  }
87  IsolatedPixelTrackCandidate(const l1t::TauRef & tauRef, double enIn, double enOut, int nhitIn, int nhitOut):
88  RecoCandidate( 0, LorentzVector(tauRef->px(),tauRef->py(),tauRef->pz(),tauRef->p()) ),
89  l1ttauJet_(tauRef), enIn_(enIn), enOut_(enOut), nhitIn_(nhitIn), nhitOut_(nhitOut) {
90  maxPtPxl_=-1;
91  sumPtPxl_=-1;
92  etaEcal_=0;
93  phiEcal_=0;
94  etaPhiEcal_=false;
95  }
98 
100  ~IsolatedPixelTrackCandidate() override;
101 
103  IsolatedPixelTrackCandidate * clone() const override;
104 
106  reco::TrackRef track() const override;
107  void setTrack( const reco::TrackRef & tr ) { track_ = tr; }
108 
110  double maxPtPxl() const {return maxPtPxl_;}
111  void setMaxPtPxl(double mptpxl) {maxPtPxl_=mptpxl;}
112 
114  double sumPtPxl() const {return sumPtPxl_;}
115  void setSumPtPxl(double sumptpxl) {sumPtPxl_=sumptpxl;}
116 
118  virtual l1extra::L1JetParticleRef l1tau() const;
119  void setL1TauJet( const l1extra::L1JetParticleRef & tauRef ) { l1tauJet_ = tauRef; }
120 
122  virtual l1t::TauRef l1ttau() const;
123  void setL1TTauJet( const l1t::TauRef & tauRef ) { l1ttauJet_ = tauRef; }
124 
126  double energyIn() const {return enIn_; }
127  void setEnergyIn(double a) {enIn_=a;}
128 
130  double energyOut() const {return enOut_;}
131  void setEnergyOut(double a) {enOut_=a;}
132 
134  int nHitIn() const {return nhitIn_;}
135  void setNHitIn(int a) {nhitIn_=a;}
136 
138  int nHitOut() const {return nhitOut_;}
139  void setNHitOut(int a) {nhitOut_=a;}
140 
142  std::pair<int,int> towerIndex() const;
143 
145  void setEtaPhiEcal(double eta, double phi) {
147  }
148  std::pair<double,double> etaPhiEcal() const {
149  return ((etaPhiEcal_) ? std::pair<double,double>(etaEcal_,phiEcal_) : std::pair<double,double>(0,0));
150  }
151  bool etaPhiEcalValid() const {return etaPhiEcal_;}
152 
153  private:
155  bool overlap( const Candidate & ) const override;
163  double maxPtPxl_;
165  double sumPtPxl_;
167  double enIn_;
169  double enOut_;
171  int nhitIn_;
173  int nhitOut_;
177  };
178 
179 
180 }
181 
182 #endif
double enOut_
energy in outer cone around L1 tau jet
IsolatedPixelTrackCandidate(const l1t::TauRef &tauRef, double enIn, double enOut, int nhitIn, int nhitOut)
constructor from tau jet using l1t
double eta() const final
momentum pseudorapidity
reco::TrackRef track_
reference to a Track
double energyIn() const
ECAL energy in the inner cone around tau jet.
l1t::TauRef l1ttauJet_
reference to a S2 L1 tau jet
double px() const final
x coordinate of momentum vector
IsolatedPixelTrackCandidate * clone() const override
returns a clone of the candidate
IsolatedPixelTrackCandidate(const l1extra::L1JetParticleRef &tauRef, double enIn, double enOut, int nhitIn, int nhitOut)
constructor from tau jet
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.
double pz() const final
z coordinate of momentum vector
int nHitOut() const
number of ECAL hits in the outer cone around tau jet
void setL1TauJet(const l1extra::L1JetParticleRef &tauRef)
void setL1TTauJet(const l1t::TauRef &tauRef)
int nHitIn() const
number of ECAL hits in the inner cone around tau jet
T get() const
get a component
IsolatedPixelTrackCandidate(const LorentzVector &v)
constructor from LorentzVector
void setEtaPhiEcal(double eta, double phi)
eta, phi at ECAL surface
double maxPtPxl_
highest Pt of other pixel tracks in the cone around the candidate
int nhitIn_
number of hits in inner cone
bool overlap(const Candidate &) const override
check overlap with another candidate
virtual l1t::TauRef l1ttau() const
get reference to L1 tau jet from lt1
virtual l1extra::L1JetParticleRef l1tau() const
get reference to L1 tau jet
std::pair< double, double > etaPhiEcal() const
double p() const final
magnitude of momentum vector
int nhitOut_
number of hits in inner cone
double energyOut() const
ECAL energy in the outer cone around tau jet.
reco::TrackRef track() const override
refrence to a Track
double py() const final
y coordinate of momentum vector
double enIn_
energy in inner cone around L1 tau jet
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
IsolatedPixelTrackCandidate(const reco::TrackRef &tr, const l1t::TauRef &tauRef, double max, double sum)
double maxPtPxl() const
highest Pt of other pixel tracks in the cone around the candidate
fixed size matrix
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
double phi() const final
momentum azimuthal angle
void setTrack(const reco::TrackRef &tr)
std::pair< int, int > towerIndex() const
get index of tower which track is hitting