Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-09 22:37:29

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/L1TParticleFlow/interface/PFCandidate.h"
0007 #include "DataFormats/L1TParticleFlow/interface/jets.h"
0008 #include "DataFormats/Common/interface/Ptr.h"
0009 
0010 namespace l1t {
0011 
0012   class PFJet : public L1Candidate {
0013   public:
0014     /// constituent information. note that this is not going to be available in the hardware!
0015     typedef std::vector<edm::Ptr<l1t::PFCandidate>> Constituents;
0016 
0017     PFJet() {}
0018     PFJet(float pt, float eta, float phi, float mass, int hwpt = 0, int hweta = 0, int hwphi = 0)
0019         : L1Candidate(PolarLorentzVector(pt, eta, phi, mass), hwpt, hweta, hwphi, /*hwQuality=*/0), rawPt_(pt) {}
0020 
0021     PFJet(const LorentzVector& p4, int hwpt = 0, int hweta = 0, int hwphi = 0)
0022         : L1Candidate(p4, hwpt, hweta, hwphi, /*hwQuality=*/0), rawPt_(p4.Pt()) {}
0023 
0024     // change the pt (but doesn't change the raw pt)
0025     void calibratePt(float newpt);
0026 
0027     // return the raw pT()
0028     float rawPt() const { return rawPt_; }
0029 
0030     /// constituent information. note that this is not going to be available in the hardware!
0031     const Constituents& constituents() const { return constituents_; }
0032     /// 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
0033     void addConstituent(const edm::Ptr<l1t::PFCandidate>& cand) { constituents_.emplace_back(cand); }
0034 
0035     // add jet tag prediction results
0036     void addTagScores(std::vector<float> scores, std::vector<l1ct::JetTagClass> classes, float ptcorrection) {
0037       tagScores_ = scores;
0038       tagClasses_ = classes;
0039       ptCorrection_ = ptcorrection;
0040     }
0041 
0042     std::vector<float> getTagScores() const { return tagScores_; }
0043     std::vector<l1ct::JetTagClass> getTagClasses() const { return tagClasses_; }
0044 
0045     float getTagScore(l1ct::JetTagClass tagClass) const {
0046       // get the tag score for a specific tagClass
0047       auto it = std::find(tagClasses_.begin(), tagClasses_.end(), tagClass);
0048       if (it != tagClasses_.end()) {
0049         return tagScores_[std::distance(tagClasses_.begin(), it)];
0050       }
0051       return -1.0f;  // Return an invalid/default score if tagClass is not found
0052     }
0053 
0054     float getTagScore(const std::string& tagClassStr) const {
0055       // get the tag score for a specific tagClass
0056       l1ct::JetTagClass tagClass(tagClassStr);
0057       return getTagScore(tagClass);
0058     }
0059 
0060     float getPtCorrection() const { return ptCorrection_; }
0061 
0062     // candidate interface
0063     size_t numberOfDaughters() const override { return constituents_.size(); }
0064     const reco::Candidate* daughter(size_type i) const override { return constituents_[i].get(); }
0065     using reco::LeafCandidate::daughter;  // avoid hiding the base
0066     edm::Ptr<l1t::PFCandidate> daughterPtr(size_type i) const { return constituents_[i]; }
0067 
0068     // Get and set the encodedJet_ bits. The Jet is encoded in 128 bits as a 2-element array of uint64_t
0069     // We store encodings both for Correlator internal usage and for Global Trigger
0070     enum class HWEncoding { CT, GT, GTWide };
0071     typedef std::array<uint64_t, 2> PackedJet;
0072     const PackedJet& encodedJet(const HWEncoding encoding = HWEncoding::GT) const {
0073       return encodedJet_[static_cast<int>(encoding)];
0074     }
0075     void setEncodedJet(const HWEncoding encoding, const PackedJet jet) {
0076       encodedJet_[static_cast<int>(encoding)] = jet;
0077     }
0078 
0079     // Accessors to HW objects with ap_* types from encoded words
0080     const PackedJet& getHWJetGT() const { return encodedJet(HWEncoding::GT); }
0081     const PackedJet& getHWJetGTWide() const { return encodedJet(HWEncoding::GTWide); }
0082     const PackedJet& getHWJetCT() const { return encodedJet(HWEncoding::CT); }
0083 
0084   private:
0085     float rawPt_;
0086     Constituents constituents_;
0087     std::vector<l1ct::JetTagClass> tagClasses_;
0088     std::vector<float> tagScores_;
0089     float ptCorrection_;
0090     std::array<PackedJet, 3> encodedJet_ = {{{{0, 0}}, {{0, 0}}, {{0, 0}}}};
0091   };
0092 
0093   typedef std::vector<l1t::PFJet> PFJetCollection;
0094   typedef edm::Ref<l1t::PFJetCollection> PFJetRef;
0095   typedef edm::RefVector<l1t::PFJetCollection> PFJetRefVector;
0096   typedef std::vector<l1t::PFJetRef> PFJetVectorRef;
0097 }  // namespace l1t
0098 #endif