Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-06-10 01:50:16

0001 #ifndef DataFormats_L1TParticleFlow_PFJet_h
0002 #define DataFormats_L1TParticleFlow_PFJet_h
0003 
0004 #include <vector>
0005 #include "DataFormats/L1Trigger/interface/L1Candidate.h"
0006 #include "DataFormats/Common/interface/Ptr.h"
0007 
0008 namespace l1t {
0009 
0010   class PFJet : public L1Candidate {
0011   public:
0012     /// constituent information. note that this is not going to be available in the hardware!
0013     typedef std::vector<edm::Ptr<l1t::L1Candidate>> Constituents;
0014 
0015     PFJet() {}
0016     PFJet(float pt, float eta, float phi, float mass = 0, int hwpt = 0, int hweta = 0, int hwphi = 0)
0017         : L1Candidate(PolarLorentzVector(pt, eta, phi, mass), hwpt, hweta, hwphi, /*hwQuality=*/0), rawPt_(pt) {}
0018 
0019     PFJet(const LorentzVector& p4, int hwpt = 0, int hweta = 0, int hwphi = 0)
0020         : L1Candidate(p4, hwpt, hweta, hwphi, /*hwQuality=*/0), rawPt_(p4.Pt()) {}
0021 
0022     // change the pt (but doesn't change the raw pt)
0023     void calibratePt(float newpt);
0024 
0025     // return the raw pT()
0026     float rawPt() const { return rawPt_; }
0027 
0028     /// constituent information. note that this is not going to be available in the hardware!
0029     const Constituents& constituents() const { return constituents_; }
0030     /// adds a candidate to this cluster; note that this only records the information, it's up to you to also set the 4-vector appropriately
0031     void addConstituent(const edm::Ptr<l1t::L1Candidate>& cand) { constituents_.emplace_back(cand); }
0032 
0033     // candidate interface
0034     size_t numberOfDaughters() const override { return constituents_.size(); }
0035     const reco::Candidate* daughter(size_type i) const override { return constituents_[i].get(); }
0036     using reco::LeafCandidate::daughter;  // avoid hiding the base
0037     edm::Ptr<l1t::L1Candidate> daughterPtr(size_type i) const { return constituents_[i]; }
0038 
0039     // Get and set the encodedJet_ bits. The Jet is encoded in 128 bits as a 2-element array of uint64_t
0040     std::array<uint64_t, 2> encodedJet() { return encodedJet_; }
0041     void setEncodedJet(std::array<uint64_t, 2> jet) { encodedJet_ = jet; }
0042 
0043   private:
0044     float rawPt_;
0045     Constituents constituents_;
0046     std::array<uint64_t, 2> encodedJet_ = {{0, 0}};
0047   };
0048 
0049   typedef std::vector<l1t::PFJet> PFJetCollection;
0050   typedef edm::Ref<l1t::PFJetCollection> PFJetRef;
0051   typedef std::vector<l1t::PFJetRef> PFJetVectorRef;
0052 }  // namespace l1t
0053 #endif