CMS 3D CMS Logo

PFTau.h
Go to the documentation of this file.
1 #ifndef DataFormats_L1TParticleFlow_PFTau_h
2 #define DataFormats_L1TParticleFlow_PFTau_h
3 
4 #include <algorithm>
5 #include <vector>
9 
10 namespace l1t {
11 
12  static constexpr float PFTAU_NN_OFFSET = 0.1;
13  static constexpr float PFTAU_NN_SLOPE = 0.2;
14  static constexpr float PFTAU_NN_OVERALL_SCALE = 1. / 20.1;
15 
16  static constexpr float PFTAU_NN_LOOSE_CUT = 0.28;
17  static constexpr float PFTAU_NN_TIGHT_CUT = 0.25;
18 
19  static constexpr float PFTAU_PF_LOOSE_CUT = 10.0;
20  static constexpr float PFTAU_PF_TIGHT_CUT = 5.0;
21 
22  static constexpr float PTSCALING_MASSCUT = 40.0;
23 
24  static constexpr double PFTAU_NN_PT_CUTOFF = 100.0;
25 
26  class PFTau : public L1Candidate {
27  public:
28  PFTau() {}
29  enum { unidentified = 0, oneprong = 1, oneprongpi0 = 2, threeprong = 3 };
31  float iVector[80],
32  float iso = -1,
33  float fulliso = -1,
34  int id = 0,
35  int hwpt = 0,
36  int hweta = 0,
37  int hwphi = 0)
38  : PFTau(PolarLorentzVector(p), iVector, iso, id, hwpt, hweta, hwphi) {}
39  PFTau(const PolarLorentzVector& p,
40  float iVector[80],
41  float iso = -1,
42  float fulliso = -1,
43  int id = 0,
44  int hwpt = 0,
45  int hweta = 0,
46  int hwphi = 0);
47  float chargedIso() const { return iso_; }
48  float fullIso() const { return fullIso_; }
49  int id() const { return id_; }
50 
51  void setZ0(float z0) { setVertex(reco::Particle::Point(0, 0, z0)); }
52  void setDxy(float dxy) { dxy_ = dxy; }
53 
54  float z0() const { return vz(); }
55  float dxy() const { return dxy_; }
56  const float* NNValues() const { return NNValues_; }
57 
58  bool passMass() const { return (mass() < 2 + pt() / PTSCALING_MASSCUT); }
59  bool passLooseNN() const { return iso_ > PFTAU_NN_LOOSE_CUT; }
60  bool passLooseNNMass() const {
61  if (!passMass())
62  return false;
63  return passLooseNN();
64  }
65  bool passLoosePF() const { return fullIso_ < PFTAU_PF_LOOSE_CUT; }
66  bool passTightNN() const {
69  }
70  bool passTightNNMass() const {
71  if (!passMass())
72  return false;
73  return passTightNN();
74  }
75  bool passTightPF() const { return fullIso_ < PFTAU_PF_TIGHT_CUT; }
76 
77  //Tau encoding for GT
79  l1gt::PackedTau encodedTau() const { return encodedTau_; } //Can be unpacked using l1gt::Tau::unpack()
80 
81  //Return the l1gt Tau object from the encoded objects
83 
84  private:
85  float NNValues_[80]; // Values for each of the 80 NN inputs
86  float iso_;
87  float fullIso_;
88  int id_;
89  float dxy_;
91  };
92 
93  typedef std::vector<l1t::PFTau> PFTauCollection;
94 
97  typedef std::vector<l1t::PFTauRef> PFTauVectorRef;
98 } // namespace l1t
99 #endif
void set_encodedTau(l1gt::PackedTau encodedTau)
Definition: PFTau.h:78
static Tau unpack(const PackedTau &src)
Definition: gt_datatypes.h:232
std::vector< l1t::PFTau > PFTauCollection
Definition: PFTau.h:93
float fullIso() const
Definition: PFTau.h:48
l1gt::PackedTau encodedTau() const
Definition: PFTau.h:79
std::vector< l1t::PFTauRef > PFTauVectorRef
Definition: PFTau.h:97
double pt() const final
transverse momentum
double vz() const override
z coordinate of vertex position
void setDxy(float dxy)
Definition: PFTau.h:52
static constexpr float PFTAU_NN_OVERALL_SCALE
Definition: PFTau.h:14
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: LeafCandidate.h:25
static constexpr float PTSCALING_MASSCUT
Definition: PFTau.h:22
bool passLooseNNMass() const
Definition: PFTau.h:60
float dxy() const
Definition: PFTau.h:55
delete x;
Definition: CaloConfig.h:22
static constexpr double PFTAU_NN_PT_CUTOFF
Definition: PFTau.h:24
void setVertex(const Point &vertex) override
set vertex
bool passMass() const
Definition: PFTau.h:58
std::array< uint64_t, 2 > PackedTau
Definition: gt_datatypes.h:38
bool passLoosePF() const
Definition: PFTau.h:65
edm::RefVector< l1t::PFTauCollection > PFTauRefVector
Definition: PFTau.h:96
double p() const final
magnitude of momentum vector
int id_
Definition: PFTau.h:88
static constexpr float PFTAU_NN_OFFSET
Definition: PFTau.h:12
bool passLooseNN() const
Definition: PFTau.h:59
bool passTightNN() const
Definition: PFTau.h:66
l1gt::Tau getHWTauGT() const
Definition: PFTau.h:82
int id() const
Definition: PFTau.h:49
float fullIso_
Definition: PFTau.h:87
math::XYZPoint Point
point in the space
Definition: Particle.h:25
float dxy_
Definition: PFTau.h:89
static constexpr float PFTAU_NN_SLOPE
Definition: PFTau.h:13
PFTau(const LorentzVector &p, float iVector[80], float iso=-1, float fulliso=-1, int id=0, int hwpt=0, int hweta=0, int hwphi=0)
Definition: PFTau.h:30
edm::Ref< l1t::PFTauCollection > PFTauRef
Definition: PFTau.h:95
float NNValues_[80]
Definition: PFTau.h:85
static constexpr float PFTAU_PF_LOOSE_CUT
Definition: PFTau.h:19
bool passTightNNMass() const
Definition: PFTau.h:70
bool passTightPF() const
Definition: PFTau.h:75
static constexpr float PFTAU_NN_TIGHT_CUT
Definition: PFTau.h:17
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
float iso_
Definition: PFTau.h:86
double mass() const final
mass
static constexpr float PFTAU_NN_LOOSE_CUT
Definition: PFTau.h:16
void setZ0(float z0)
Definition: PFTau.h:51
float chargedIso() const
Definition: PFTau.h:47
const float * NNValues() const
Definition: PFTau.h:56
static constexpr float PFTAU_PF_TIGHT_CUT
Definition: PFTau.h:20
PFTau()
Definition: PFTau.h:28
float z0() const
Definition: PFTau.h:54
l1gt::PackedTau encodedTau_
Definition: PFTau.h:90
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38