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  public:
25  // default constructor
27  enIn_ = -1;
28  enOut_ = -1;
29  nhitIn_ = -1;
30  nhitOut_ = -1;
31  maxPtPxl_ = -1;
32  sumPtPxl_ = -1;
33  etaEcal_ = 0;
34  phiEcal_ = 0;
35  etaPhiEcal_ = false;
36  }
39  enIn_ = -1;
40  enOut_ = -1;
41  nhitIn_ = -1;
42  nhitOut_ = -1;
43  maxPtPxl_ = -1;
44  sumPtPxl_ = -1;
45  etaEcal_ = 0;
46  phiEcal_ = 0;
47  etaPhiEcal_ = false;
48  }
51  const l1extra::L1JetParticleRef& tauRef,
52  double max,
53  double sum)
54  : RecoCandidate(0, LorentzVector((tr.get()->px()), (tr.get())->py(), (tr.get())->pz(), (tr.get())->p())),
55  track_(tr),
56  l1tauJet_(tauRef),
57  maxPtPxl_(max),
58  sumPtPxl_(sum) {
59  enIn_ = -1;
60  enOut_ = -1;
61  nhitIn_ = -1;
62  nhitOut_ = -1;
63  etaEcal_ = 0;
64  phiEcal_ = 0;
65  etaPhiEcal_ = false;
66  }
67  // constructor from a track using l1t
68  IsolatedPixelTrackCandidate(const reco::TrackRef& tr, const l1t::TauRef& tauRef, double max, double sum)
69  : RecoCandidate(0, LorentzVector((tr.get()->px()), (tr.get())->py(), (tr.get())->pz(), (tr.get())->p())),
70  track_(tr),
71  l1ttauJet_(tauRef),
72  maxPtPxl_(max),
73  sumPtPxl_(sum) {
74  enIn_ = -1;
75  enOut_ = -1;
76  nhitIn_ = -1;
77  nhitOut_ = -1;
78  etaEcal_ = 0;
79  phiEcal_ = 0;
80  etaPhiEcal_ = false;
81  }
82 
85  const l1extra::L1JetParticleRef& tauRef, double enIn, double enOut, int nhitIn, int nhitOut)
86  : RecoCandidate(0, LorentzVector(tauRef->px(), tauRef->py(), tauRef->pz(), tauRef->p())),
87  l1tauJet_(tauRef),
88  enIn_(enIn),
89  enOut_(enOut),
90  nhitIn_(nhitIn),
91  nhitOut_(nhitOut) {
92  maxPtPxl_ = -1;
93  sumPtPxl_ = -1;
94  etaEcal_ = 0;
95  phiEcal_ = 0;
96  etaPhiEcal_ = false;
97  }
99  IsolatedPixelTrackCandidate(const l1t::TauRef& tauRef, double enIn, double enOut, int nhitIn, int nhitOut)
100  : RecoCandidate(0, LorentzVector(tauRef->px(), tauRef->py(), tauRef->pz(), tauRef->p())),
101  l1ttauJet_(tauRef),
102  enIn_(enIn),
103  enOut_(enOut),
104  nhitIn_(nhitIn),
105  nhitOut_(nhitOut) {
106  maxPtPxl_ = -1;
107  sumPtPxl_ = -1;
108  etaEcal_ = 0;
109  phiEcal_ = 0;
110  etaPhiEcal_ = false;
111  }
114 
116  ~IsolatedPixelTrackCandidate() override;
117 
119  IsolatedPixelTrackCandidate* clone() const override;
120 
122  reco::TrackRef track() const override;
123  void setTrack(const reco::TrackRef& tr) { track_ = tr; }
124 
126  double maxPtPxl() const { return maxPtPxl_; }
127  void setMaxPtPxl(double mptpxl) { maxPtPxl_ = mptpxl; }
128 
130  double sumPtPxl() const { return sumPtPxl_; }
131  void setSumPtPxl(double sumptpxl) { sumPtPxl_ = sumptpxl; }
132 
134  virtual l1extra::L1JetParticleRef l1tau() const;
135  void setL1TauJet(const l1extra::L1JetParticleRef& tauRef) { l1tauJet_ = tauRef; }
136 
138  virtual l1t::TauRef l1ttau() const;
139  void setL1TTauJet(const l1t::TauRef& tauRef) { l1ttauJet_ = tauRef; }
140 
142  double energyIn() const { return enIn_; }
143  void setEnergyIn(double a) { enIn_ = a; }
144 
146  double energyOut() const { return enOut_; }
147  void setEnergyOut(double a) { enOut_ = a; }
148 
150  int nHitIn() const { return nhitIn_; }
151  void setNHitIn(int a) { nhitIn_ = a; }
152 
154  int nHitOut() const { return nhitOut_; }
155  void setNHitOut(int a) { nhitOut_ = a; }
156 
158  std::pair<int, int> towerIndex() const;
159 
161  void setEtaPhiEcal(double eta, double phi) {
162  etaEcal_ = eta;
163  phiEcal_ = phi;
164  etaPhiEcal_ = true;
165  }
166  std::pair<double, double> etaPhiEcal() const {
167  return ((etaPhiEcal_) ? std::pair<double, double>(etaEcal_, phiEcal_) : std::pair<double, double>(0, 0));
168  }
169  bool etaPhiEcalValid() const { return etaPhiEcal_; }
170 
171  private:
173  bool overlap(const Candidate&) const override;
181  double maxPtPxl_;
183  double sumPtPxl_;
185  double enIn_;
187  double enOut_;
189  int nhitIn_;
191  int nhitOut_;
195  };
196 
197 } // namespace reco
198 
199 #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 pz() const final
z coordinate of momentum vector
reco::TrackRef track_
reference to a Track
virtual l1extra::L1JetParticleRef l1tau() const
get reference to L1 tau jet
l1t::TauRef l1ttauJet_
reference to a S2 L1 tau jet
IsolatedPixelTrackCandidate(const l1extra::L1JetParticleRef &tauRef, double enIn, double enOut, int nhitIn, int nhitOut)
constructor from tau jet
T get() const
get a component
virtual l1t::TauRef l1ttau() const
get reference to L1 tau jet from lt1
double maxPtPxl() const
highest Pt 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 px() const final
x coordinate of momentum vector
IsolatedPixelTrackCandidate * clone() const override
returns a clone of the candidate
double p() const final
magnitude of momentum vector
void setL1TauJet(const l1extra::L1JetParticleRef &tauRef)
void setL1TTauJet(const l1t::TauRef &tauRef)
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
double py() const final
y coordinate of momentum vector
int nHitOut() const
number of ECAL hits in the outer cone around tau jet
int nhitOut_
number of hits in inner cone
std::pair< int, int > towerIndex() const
get index of tower which track is hitting
double enIn_
energy in inner cone around L1 tau jet
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
IsolatedPixelTrackCandidate(const reco::TrackRef &tr, const l1t::TauRef &tauRef, double max, double sum)
double energyIn() const
ECAL energy in the inner cone around tau jet.
reco::TrackRef track() const override
refrence to a Track
fixed size matrix
double a
Definition: hdecay.h:119
l1extra::L1JetParticleRef l1tauJet_
reference to a L1 tau jet
int nHitIn() const
number of ECAL hits in the inner cone around tau jet
double phi() const final
momentum azimuthal angle
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< double, double > etaPhiEcal() const
double energyOut() const
ECAL energy in the outer cone around tau jet.
bool overlap(const Candidate &) const override
check overlap with another candidate
double sumPtPxl() const
Pt sum of other pixel tracks in the cone around the candidate.
double eta() const final
momentum pseudorapidity