CMS 3D CMS Logo

PFCluster.h
Go to the documentation of this file.
1 #ifndef DataFormats_L1TParticleFlow_PFCluster_h
2 #define DataFormats_L1TParticleFlow_PFCluster_h
3 
4 #include <vector>
7 
8 namespace l1t {
9 
10  class PFCluster : public L1Candidate {
11  public:
13  typedef std::pair<edm::Ptr<l1t::L1Candidate>, float> ConstituentAndFraction;
14  typedef std::vector<ConstituentAndFraction> ConstituentsAndFractions;
15 
16  PFCluster() {}
17  PFCluster(float pt,
18  float eta,
19  float phi,
20  float hOverE = 0,
21  bool isEM = false,
22  float ptError = 0,
23  int hwpt = 0,
24  int hweta = 0,
25  int hwphi = 0)
26  : L1Candidate(PolarLorentzVector(pt, eta, phi, 0), hwpt, hweta, hwphi, /*hwQuality=*/isEM ? 1 : 0),
27  hOverE_(hOverE),
28  ptError_(ptError) {
29  setPdgId(isEM ? 22 : 130); // photon : non-photon(K0)
30  }
32  const LorentzVector& p4, float hOverE, bool isEM, float ptError = 0, int hwpt = 0, int hweta = 0, int hwphi = 0)
33  : L1Candidate(p4, hwpt, hweta, hwphi, /*hwQuality=*/isEM ? 1 : 0), hOverE_(hOverE), ptError_(ptError) {
34  setPdgId(isEM ? 22 : 130); // photon : non-photon(K0)
35  }
36 
37  float hOverE() const { return hOverE_; }
38  void setHOverE(float hOverE) { hOverE_ = hOverE; }
39 
40  float emEt() const {
41  if (hOverE_ == -1)
42  return 0;
43  return pt() / (1 + hOverE_);
44  }
45 
46  // change the pt. H/E also recalculated to keep emEt constant
47  void calibratePt(float newpt, float preserveEmEt = true);
48 
53  constituents_.emplace_back(cand, fraction);
54  }
55 
56  float ptError() const { return ptError_; }
57  void setPtError(float ptError) { ptError_ = ptError; }
58 
59  bool isEM() const { return hwQual(); }
60  void setIsEM(bool isEM) { setHwQual(isEM); }
61  unsigned int hwEmID() const { return hwQual(); }
62 
63  float egVsPionMVAOut() const { return egVsPionMVAOut_; }
65 
66  float egVsPUMVAOut() const { return egVsPUMVAOut_; }
68 
69  private:
72  };
73 
74  typedef std::vector<l1t::PFCluster> PFClusterCollection;
76 } // namespace l1t
77 #endif
std::vector< ConstituentAndFraction > ConstituentsAndFractions
Definition: PFCluster.h:14
float hOverE() const
Definition: PFCluster.h:37
double pt() const final
transverse momentum
void setEgVsPionMVAOut(float egVsPionMVAOut)
Definition: PFCluster.h:64
float ptError_
Definition: PFCluster.h:70
void setHOverE(float hOverE)
Definition: PFCluster.h:38
void setHwQual(int qual)
Definition: L1Candidate.h:31
delete x;
Definition: CaloConfig.h:22
float egVsPionMVAOut_
Definition: PFCluster.h:70
void setEgVsPUMVAOut(float egVsPUMVAOut)
Definition: PFCluster.h:67
ConstituentsAndFractions constituents_
Definition: PFCluster.h:71
void calibratePt(float newpt, float preserveEmEt=true)
Definition: PFCluster.cc:3
const LorentzVector & p4() const final
four-momentum Lorentz vector
int hwQual() const
Definition: L1Candidate.h:38
float ptError() const
Definition: PFCluster.h:56
float egVsPUMVAOut() const
Definition: PFCluster.h:66
std::pair< edm::Ptr< l1t::L1Candidate >, float > ConstituentAndFraction
constituent information. note that this is not going to be available in the hardware! ...
Definition: PFCluster.h:13
PFCluster(const LorentzVector &p4, float hOverE, bool isEM, float ptError=0, int hwpt=0, int hweta=0, int hwphi=0)
Definition: PFCluster.h:31
PFCluster(float pt, float eta, float phi, float hOverE=0, bool isEM=false, float ptError=0, int hwpt=0, int hweta=0, int hwphi=0)
Definition: PFCluster.h:17
std::vector< l1t::PFCluster > PFClusterCollection
Definition: PFCluster.h:74
const ConstituentsAndFractions & constituentsAndFractions() const
constituent information. note that this is not going to be available in the hardware! ...
Definition: PFCluster.h:50
float egVsPionMVAOut() const
Definition: PFCluster.h:63
bool isEM() const
Definition: PFCluster.h:59
void setIsEM(bool isEM)
Definition: PFCluster.h:60
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
float egVsPUMVAOut_
Definition: PFCluster.h:70
void setPtError(float ptError)
Definition: PFCluster.h:57
unsigned int hwEmID() const
Definition: PFCluster.h:61
void addConstituent(const edm::Ptr< l1t::L1Candidate > &cand, float fraction=1.0)
adds a candidate to this cluster; note that this only records the information, it&#39;s up to you to also...
Definition: PFCluster.h:52
double phi() const final
momentum azimuthal angle
void setPdgId(int pdgId) final
float emEt() const
Definition: PFCluster.h:40
edm::Ref< l1t::PFClusterCollection > PFClusterRef
Definition: PFCluster.h:75
float hOverE_
Definition: PFCluster.h:70
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38
double eta() const final
momentum pseudorapidity