CMS 3D CMS Logo

PileupJetIdAlgo.h
Go to the documentation of this file.
1 //--------------------------------------------------------------------------------------------------
2 //
3 // PileupJetIdAlgo
4 //
5 // Author: P. Musella, P. Harris
6 //--------------------------------------------------------------------------------------------------
7 
8 #ifndef RecoJets_JetProducers_plugins_PileupJetIdAlgo_h
9 #define RecoJets_JetProducers_plugins_PileupJetIdAlgo_h
10 
17 
20 
21 // ----------------------------------------------------------------------------------------------------
23 public:
24  enum version_t { USER = -1, PHILv0 = 0 };
25 
26  class AlgoGBRForestsAndConstants;
27 
28  PileupJetIdAlgo(AlgoGBRForestsAndConstants const* cache);
30 
32  const reco::Jet* jet, float jec, const reco::Vertex*, const reco::VertexCollection&, double rho, bool usePuppi);
33 
34  void set(const PileupJetIdentifier&);
35  float getMVAval(const std::vector<std::string>&, const std::unique_ptr<const GBRForest>&);
37  const std::string method() const { return cache_->tmvaMethod(); }
38 
39  std::string dumpVariables() const;
40 
41  typedef std::map<std::string, std::pair<float*, float>> variables_list_t;
42 
43  std::pair<int, int> getJetIdKey(float jetPt, float jetEta);
44  int computeCutIDflag(float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta);
45  int computeIDflag(float mva, float jetPt, float jetEta);
46  int computeIDflag(float mva, int ptId, int etaId);
47 
49  const variables_list_t& getVariables() const { return variables_; };
50 
51  // In multithreaded mode, each PileupIdAlgo object will get duplicated
52  // on every stream. Some of the data it contains never changes after
53  // construction and can be shared by all streams. This nested class contains
54  // the data members that can be shared across streams. In particular
55  // the GBRForests take significant time to initialize and can be shared.
57  public:
59 
60  std::unique_ptr<const GBRForest> const& reader() const { return reader_; }
61  std::vector<std::unique_ptr<const GBRForest>> const& etaReader() const { return etaReader_; }
62  bool cutBased() const { return cutBased_; }
63  bool etaBinnedWeights() const { return etaBinnedWeights_; }
64  bool runMvas() const { return runMvas_; }
65  int nEtaBins() const { return nEtaBins_; }
66  std::vector<double> const& jEtaMin() const { return jEtaMin_; }
67  std::vector<double> const& jEtaMax() const { return jEtaMax_; }
68  std::string const& label() const { return label_; }
69  std::string const& tmvaMethod() const { return tmvaMethod_; }
70  std::vector<std::string> const& tmvaVariables() const { return tmvaVariables_; }
71  std::vector<std::vector<std::string>> const& tmvaEtaVariables() const { return tmvaEtaVariables_; }
72 
73  typedef float array_t[3][5][4];
74  array_t const& mvacut() const { return mvacut_; }
75  array_t const& rmsCut() const { return rmsCut_; }
76  array_t const& betaStarCut() const { return betaStarCut_; }
77 
78  private:
79  std::unique_ptr<const GBRForest> reader_;
80  std::vector<std::unique_ptr<const GBRForest>> etaReader_;
81  bool cutBased_;
83  bool runMvas_;
84  int nEtaBins_;
85  std::vector<double> jEtaMin_;
86  std::vector<double> jEtaMax_;
89  std::vector<std::string> tmvaVariables_;
90  std::vector<std::vector<std::string>> tmvaEtaVariables_;
91 
92  float mvacut_[3][5][4]; //Keep the array fixed
93  float rmsCut_[3][5][4]; //Keep the array fixed
94  float betaStarCut_[3][5][4]; //Keep the array fixed
95 
96  std::map<std::string, std::string> tmvaNames_;
97  };
98 
99 protected:
100  void runMva();
101  void resetVariables();
102  void initVariables();
103 
107 };
108 #endif
const std::string method() const
std::map< std::string, std::pair< float *, float > > variables_list_t
std::vector< double > const & jEtaMin() const
AlgoGBRForestsAndConstants const * cache_
AlgoGBRForestsAndConstants(edm::ParameterSet const &, bool runMvas)
std::vector< std::unique_ptr< const GBRForest > > etaReader_
Base class for all types of Jets.
Definition: Jet.h:20
std::string dumpVariables() const
std::pair< int, int > getJetIdKey(float jetPt, float jetEta)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
int computeCutIDflag(float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta)
std::vector< std::vector< std::string > > const & tmvaEtaVariables() const
std::vector< std::string > tmvaVariables_
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, bool usePuppi)
PileupJetIdentifier internalId_
std::map< std::string, std::string > tmvaNames_
std::vector< std::string > const & tmvaVariables() const
std::vector< std::vector< std::string > > tmvaEtaVariables_
int computeIDflag(float mva, float jetPt, float jetEta)
std::vector< std::unique_ptr< const GBRForest > > const & etaReader() const
std::string const & tmvaMethod() const
const variables_list_t & getVariables() const
const PileupJetIdentifier::variables_list_t & getVariables() const { return variables_; }; ...
PileupJetIdentifier computeMva()
std::unique_ptr< const GBRForest > const & reader() const
def cache(function)
Definition: utilities.py:3
variables_list_t variables_
std::vector< double > const & jEtaMax() const
float getMVAval(const std::vector< std::string > &, const std::unique_ptr< const GBRForest > &)
std::unique_ptr< const GBRForest > reader_
PileupJetIdAlgo(AlgoGBRForestsAndConstants const *cache)