Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-20 01:53:33

0001 #ifndef L1Trigger_Phase2L1ParticleFlow_corrector_h
0002 #define L1Trigger_Phase2L1ParticleFlow_corrector_h
0003 #include <TGraph.h>
0004 #include <TH1.h>
0005 #include <string>
0006 #include <vector>
0007 
0008 class TDirectory;
0009 
0010 namespace l1t {
0011   class PFCluster;
0012 }
0013 
0014 namespace l1tpf {
0015   class corrector {
0016   public:
0017     enum class EmulationMode { None, Correction, CorrectedPt };
0018 
0019     corrector()
0020         : is2d_(false),
0021           neta_(0),
0022           nemf_(0),
0023           emfMax_(-1),
0024           emulate_(false),
0025           debug_(false),
0026           emulationMode_(l1tpf::corrector::EmulationMode::CorrectedPt) {}
0027     corrector(const std::string &iFile,
0028               float emfMax = -1,
0029               bool debug = false,
0030               bool emulate = false,
0031               l1tpf::corrector::EmulationMode emulationMode = l1tpf::corrector::EmulationMode::CorrectedPt);
0032     corrector(const std::string &iFile,
0033               const std::string &directory,
0034               float emfMax = -1,
0035               bool debug = false,
0036               bool emulate = false,
0037               l1tpf::corrector::EmulationMode emulationMode = l1tpf::corrector::EmulationMode::CorrectedPt);
0038     corrector(TDirectory *src,
0039               float emfMax = -1,
0040               bool debug = false,
0041               bool emulate = false,
0042               l1tpf::corrector::EmulationMode emulationMode = l1tpf::corrector::EmulationMode::CorrectedPt);
0043     // create an empty corrector (you'll need to fill the graphs later)
0044     corrector(const TH1 *index, float emfMax = -1);
0045     ~corrector();
0046 
0047     // no copy, but can move
0048     corrector(const corrector &corr) = delete;
0049     corrector &operator=(const corrector &corr) = delete;
0050     corrector(corrector &&corr);
0051     corrector &operator=(corrector &&corr);
0052 
0053     float correctedPt(float et, float emEt, float eta) const;
0054     float correctedPt(float et, float eta) const { return correctedPt(et, 0, eta); }
0055     void correctPt(l1t::PFCluster &cluster, float preserveEmEt = true) const;
0056 
0057     bool valid() const { return (index_.get() != nullptr); }
0058 
0059     // set the graph (note: it is cloned, and the corrector owns the clone)
0060     void setGraph(const TGraph &graph, int ieta, int iemf = 0);
0061 
0062     bool is2d() const { return is2d_; }
0063     unsigned int neta() const { return neta_; }
0064     unsigned int nemf() const { return nemf_; }
0065     // access the index histogram
0066     const TH1 &getIndex() const { return *index_; }
0067     // access the graphs (owned by the corrector, may be null)
0068     TGraph *getGraph(int ieta, int iemf = 0) { return corrections_[ieta * nemf_ + iemf]; }
0069     const TGraph *getGraph(int ieta, int iemf = 0) const { return corrections_[ieta * nemf_ + iemf]; }
0070 
0071     // store the corrector
0072     void writeToFile(const std::string &filename, const std::string &directory) const;
0073     // store the corrector
0074     void writeToFile(TDirectory *dest) const;
0075 
0076   private:
0077     std::unique_ptr<TH1> index_;
0078     std::vector<TGraph *> corrections_;
0079     std::vector<TH1 *> correctionsEmulated_;
0080     bool is2d_;
0081     unsigned int neta_, nemf_;
0082     float emfMax_;
0083     bool emulate_;
0084     bool debug_;
0085     l1tpf::corrector::EmulationMode emulationMode_;
0086 
0087     void init_(const std::string &iFile, const std::string &directory, bool debug, bool emulate);
0088     void init_(TDirectory *src, bool debug);
0089     void initGraphs_(TDirectory *src, bool debug);
0090     void initEmulation_(TDirectory *src, bool debug);
0091   };
0092 }  // namespace l1tpf
0093 #endif