CMS 3D CMS Logo

PFJet.h
Go to the documentation of this file.
1 #ifndef DataFormats_L1TParticleFlow_PFJet_h
2 #define DataFormats_L1TParticleFlow_PFJet_h
3 
4 #include <vector>
8 
9 namespace l1t {
10 
11  class PFJet : public L1Candidate {
12  public:
14  typedef std::vector<edm::Ptr<l1t::PFCandidate>> Constituents;
15 
16  PFJet() {}
17  PFJet(float pt, float eta, float phi, float mass = 0, int hwpt = 0, int hweta = 0, int hwphi = 0)
18  : L1Candidate(PolarLorentzVector(pt, eta, phi, mass), hwpt, hweta, hwphi, /*hwQuality=*/0), rawPt_(pt) {}
19 
20  PFJet(const LorentzVector& p4, int hwpt = 0, int hweta = 0, int hwphi = 0)
21  : L1Candidate(p4, hwpt, hweta, hwphi, /*hwQuality=*/0), rawPt_(p4.Pt()) {}
22 
23  // change the pt (but doesn't change the raw pt)
24  void calibratePt(float newpt);
25 
26  // return the raw pT()
27  float rawPt() const { return rawPt_; }
28 
30  const Constituents& constituents() const { return constituents_; }
33 
34  // candidate interface
35  size_t numberOfDaughters() const override { return constituents_.size(); }
36  const reco::Candidate* daughter(size_type i) const override { return constituents_[i].get(); }
37  using reco::LeafCandidate::daughter; // avoid hiding the base
39 
40  // Get and set the encodedJet_ bits. The Jet is encoded in 128 bits as a 2-element array of uint64_t
41  // We store encodings both for Correlator internal usage and for Global Trigger
42  enum class HWEncoding { CT, GT };
43  typedef std::array<uint64_t, 2> PackedJet;
44  const PackedJet& encodedJet(const HWEncoding encoding = HWEncoding::GT) const {
45  return encodedJet_[static_cast<int>(encoding)];
46  }
47  void setEncodedJet(const HWEncoding encoding, const PackedJet jet) {
48  encodedJet_[static_cast<int>(encoding)] = jet;
49  }
50 
51  // Accessors to HW objects with ap_* types from encoded words
52  const PackedJet& getHWJetGT() const { return encodedJet(HWEncoding::GT); }
53  const PackedJet& getHWJetCT() const { return encodedJet(HWEncoding::CT); }
54 
55  private:
56  float rawPt_;
58  std::array<PackedJet, 2> encodedJet_ = {{{{0, 0}}, {{0, 0}}}};
59  };
60 
61  typedef std::vector<l1t::PFJet> PFJetCollection;
63  typedef std::vector<l1t::PFJetRef> PFJetVectorRef;
64 } // namespace l1t
65 #endif
double pt() const final
transverse momentum
float rawPt_
Definition: PFJet.h:56
size_t size_type
Definition: Candidate.h:29
edm::Ptr< l1t::PFCandidate > daughterPtr(size_type i) const
Definition: PFJet.h:38
PFJet(const LorentzVector &p4, int hwpt=0, int hweta=0, int hwphi=0)
Definition: PFJet.h:20
const PackedJet & getHWJetCT() const
Definition: PFJet.h:53
const PackedJet & encodedJet(const HWEncoding encoding=HWEncoding::GT) const
Definition: PFJet.h:44
const reco::Candidate * daughter(size_type i) const override
return daughter at a given position (throws an exception)
Definition: PFJet.h:36
void setEncodedJet(const HWEncoding encoding, const PackedJet jet)
Definition: PFJet.h:47
std::array< PackedJet, 2 > encodedJet_
Definition: PFJet.h:58
delete x;
Definition: CaloConfig.h:22
Constituents constituents_
Definition: PFJet.h:57
const LorentzVector & p4() const final
four-momentum Lorentz vector
std::vector< l1t::PFJet > PFJetCollection
Definition: PFJet.h:61
const Candidate * daughter(size_type) const override
return daughter at a given position (throws an exception)
float rawPt() const
Definition: PFJet.h:27
const PackedJet & getHWJetGT() const
Definition: PFJet.h:52
PFJet(float pt, float eta, float phi, float mass=0, int hwpt=0, int hweta=0, int hwphi=0)
Definition: PFJet.h:17
void calibratePt(float newpt)
Definition: PFJet.cc:3
edm::Ref< l1t::PFJetCollection > PFJetRef
Definition: PFJet.h:62
PFJet()
Definition: PFJet.h:16
std::array< uint64_t, 2 > PackedJet
Definition: PFJet.h:43
const Constituents & constituents() const
constituent information. note that this is not going to be available in the hardware! ...
Definition: PFJet.h:30
std::vector< edm::Ptr< l1t::PFCandidate > > Constituents
constituent information. note that this is not going to be available in the hardware! ...
Definition: PFJet.h:14
void addConstituent(const edm::Ptr< l1t::PFCandidate > &cand)
adds a candidate to this cluster; note that this only records the information, it&#39;s up to you to also...
Definition: PFJet.h:32
std::vector< l1t::PFJetRef > PFJetVectorRef
Definition: PFJet.h:63
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
double mass() const final
mass
double phi() const final
momentum azimuthal angle
size_t numberOfDaughters() const override
number of daughters
Definition: PFJet.h:35
HWEncoding
Definition: PFJet.h:42
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38
double eta() const final
momentum pseudorapidity