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 
18 #include "TMVA/Tools.h"
19 #include "TMVA/Reader.h"
20 #include "TMVA/Tools.h"
21 #include "TMVA/Reader.h"
22 
25 
26 // ----------------------------------------------------------------------------------------------------
28 public:
29  enum version_t { USER=-1, PHILv0=0 };
30 
32 
35 
37  float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, bool usePuppi);
38 
39  void set(const PileupJetIdentifier &);
40  float getMVAval(const std::vector<std::string> &, const std::unique_ptr<const GBRForest> &);
42  const std::string method() const { return cache_->tmvaMethod(); }
43 
44  std::string dumpVariables() const;
45 
46  typedef std::map<std::string,std::pair<float *,float> > variables_list_t;
47 
48  std::pair<int,int> getJetIdKey(float jetPt, float jetEta);
49  int computeCutIDflag(float betaStarClassic,float dR2Mean,float nvtx, float jetPt, float jetEta);
50  int computeIDflag (float mva, float jetPt, float jetEta);
51  int computeIDflag (float mva,int ptId,int etaId);
52 
54  const variables_list_t & getVariables() const { return variables_; };
55 
56  // In multithreaded mode, each PileupIdAlgo object will get duplicated
57  // on every stream. Some of the data it contains never changes after
58  // construction and can be shared by all streams. This nested class contains
59  // the data members that can be shared across streams. In particular
60  // the GBRForests take significant time to initialize and can be shared.
62  public:
63 
65 
66  std::unique_ptr<const GBRForest> const& reader() const { return reader_; }
67  std::vector<std::unique_ptr<const GBRForest>> const& etaReader() const { return etaReader_; }
68  bool cutBased() const { return cutBased_; }
69  bool etaBinnedWeights() const { return etaBinnedWeights_; }
70  bool runMvas() const { return runMvas_; }
71  int nEtaBins() const { return nEtaBins_; }
72  std::vector<double> const& jEtaMin() const { return jEtaMin_; }
73  std::vector<double> const& jEtaMax() const { return jEtaMax_; }
74  std::string const& label() const { return label_; }
75  std::string const& tmvaMethod() const { return tmvaMethod_; }
76  std::vector<std::string> const& tmvaVariables() const { return tmvaVariables_; }
77  std::vector<std::vector<std::string>> const& tmvaEtaVariables() const { return tmvaEtaVariables_; }
78 
79  typedef float array_t[3][4][4];
80  array_t const& mvacut() const { return mvacut_; }
81  array_t const& rmsCut() const { return rmsCut_; }
82  array_t const& betaStarCut() const { return betaStarCut_; }
83 
84  std::unique_ptr<const GBRForest> getMVA(std::vector<std::string> const& varList,
85  std::string const& tmvaWeights,
86  std::vector<std::string> const& tmvaSpectators);
87 
88  private:
89 
90  std::unique_ptr<const GBRForest> reader_;
91  std::vector<std::unique_ptr<const GBRForest>> etaReader_;
92  bool cutBased_;
94  bool runMvas_;
95  int nEtaBins_;
96  std::vector<double> jEtaMin_;
97  std::vector<double> jEtaMax_;
100  std::vector<std::string> tmvaVariables_;
101  std::vector<std::vector<std::string>> tmvaEtaVariables_;
102 
103  float mvacut_ [3][4][4]; //Keep the array fixed
104  float rmsCut_ [3][4][4]; //Keep the array fixed
105  float betaStarCut_[3][4][4]; //Keep the array fixed
106 
107  std::map<std::string,std::string> tmvaNames_;
108  };
109 
110 protected:
111 
112  void runMva();
113  void resetVariables();
114  void initVariables();
115 
117  variables_list_t variables_;
119 };
120 #endif
AlgoGBRForestsAndConstants const * cache_
AlgoGBRForestsAndConstants(edm::ParameterSet const &, bool runMvas)
std::vector< double > const & jEtaMin() const
std::map< std::string, std::string > tmvaNames_
std::vector< std::unique_ptr< const GBRForest > > etaReader_
Base class for all types of Jets.
Definition: Jet.h:20
std::string const & tmvaMethod() const
std::pair< int, int > getJetIdKey(float jetPt, float jetEta)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
const std::string method() const
int computeCutIDflag(float betaStarClassic, float dR2Mean, float nvtx, float jetPt, float jetEta)
std::vector< double > const & jEtaMax() const
const variables_list_t & getVariables() const
const PileupJetIdentifier::variables_list_t & getVariables() const { return variables_; }; ...
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, bool usePuppi)
PileupJetIdentifier internalId_
std::string dumpVariables() const
std::unique_ptr< const GBRForest > const & reader() const
std::vector< std::vector< std::string > > tmvaEtaVariables_
int computeIDflag(float mva, float jetPt, float jetEta)
def cache(function)
std::unique_ptr< const GBRForest > getMVA(std::vector< std::string > const &varList, std::string const &tmvaWeights, std::vector< std::string > const &tmvaSpectators)
PileupJetIdentifier computeMva()
std::vector< std::unique_ptr< const GBRForest > > const & etaReader() const
std::map< std::string, std::pair< float *, float > > variables_list_t
std::vector< std::vector< std::string > > const & tmvaEtaVariables() const
variables_list_t variables_
std::vector< std::string > const & tmvaVariables() const
float getMVAval(const std::vector< std::string > &, const std::unique_ptr< const GBRForest > &)
std::unique_ptr< const GBRForest > reader_
PileupJetIdAlgo(AlgoGBRForestsAndConstants const *cache)