CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalTriggerClusterIdentificationBDT.cc
Go to the documentation of this file.
1 #include <limits>
7 
9 public:
10  class Category {
11  public:
12  Category(float pt_min, float pt_max, float eta_min, float eta_max)
13  : pt_min_(pt_min), pt_max_(pt_max), eta_min_(eta_min), eta_max_(eta_max) {}
14  ~Category() {}
15  bool contains(float pt, float eta) const {
16  bool output = true;
17  if (pt < pt_min_ || pt >= pt_max_)
18  output = false;
19  if (std::abs(eta) < eta_min_ || std::abs(eta) >= eta_max_)
20  output = false;
21  return output;
22  }
23 
24  float pt_min() const { return pt_min_; }
25  float pt_max() const { return pt_max_; }
26  float eta_min() const { return eta_min_; }
27  float eta_max() const { return eta_max_; }
28 
29  private:
30  float pt_min_ = 0.;
32  float eta_min_ = 1.5;
33  float eta_max_ = 3.;
34  };
35 
36 public:
39  void initialize(const edm::ParameterSet& conf) final;
40  float value(const l1t::HGCalMulticluster& cluster) const final;
41  bool decision(const l1t::HGCalMulticluster& cluster) const final;
42 
43 private:
44  enum class ClusterVariable {
53  cl3d_szz,
57  };
58  std::vector<Category> categories_;
59  std::vector<std::unique_ptr<TMVAEvaluator>> bdts_;
60  std::vector<double> working_points_;
61  std::vector<std::string> input_variables_;
62  std::vector<ClusterVariable> input_variables_id_;
63 
65  int category(float pt, float eta) const;
66 };
67 
68 DEFINE_HGC_TPG_CLUSTER_ID(HGCalTriggerClusterIdentificationBDT, "HGCalTriggerClusterIdentificationBDT");
69 
72 
74  if (!bdts_.empty()) {
75  edm::LogWarning("HGCalTriggerClusterIdentificationBDT|Initialization") << "BDTs already initialized.";
76  return;
77  }
78  input_variables_ = conf.getParameter<std::vector<std::string>>("Inputs");
79  std::vector<std::string> bdt_files = conf.getParameter<std::vector<std::string>>("Weights");
80  std::vector<double> categories_etamin = conf.getParameter<std::vector<double>>("CategoriesEtaMin");
81  std::vector<double> categories_etamax = conf.getParameter<std::vector<double>>("CategoriesEtaMax");
82  std::vector<double> categories_ptmin = conf.getParameter<std::vector<double>>("CategoriesPtMin");
83  std::vector<double> categories_ptmax = conf.getParameter<std::vector<double>>("CategoriesPtMax");
84  working_points_ = conf.getParameter<std::vector<double>>("WorkingPoints");
85 
86  if (bdt_files.size() != categories_etamin.size() || categories_etamin.size() != categories_etamax.size() ||
87  categories_etamax.size() != categories_ptmin.size() || categories_ptmin.size() != categories_ptmax.size() ||
88  categories_ptmax.size() != working_points_.size()) {
89  throw cms::Exception("HGCalTriggerClusterIdentificationBDT|BadInitialization")
90  << "Inconsistent numbers of categories, BDT weight files and working points";
91  }
92  categories_.reserve(working_points_.size());
93  bdts_.reserve(working_points_.size());
94  for (unsigned cat = 0; cat < categories_etamin.size(); cat++) {
95  categories_.emplace_back(
96  categories_ptmin[cat], categories_ptmax[cat], categories_etamin[cat], categories_etamax[cat]);
97  }
98  std::vector<std::string> spectators = {};
99  for (const auto& file : bdt_files) {
100  bdts_.emplace_back(new TMVAEvaluator());
101  bdts_.back()->initialize("!Color:Silent:!Error",
102  "BDT::bdt",
105  spectators,
106  false,
107  false);
108  }
109 
110  // Transform input variable strings to enum values for later comparisons
111  input_variables_id_.reserve(input_variables_.size());
112  for (const auto& variable : input_variables_) {
113  if (variable == "cl3d_showerlength")
115  else if (variable == "cl3d_coreshowerlength")
117  else if (variable == "cl3d_firstlayer")
119  else if (variable == "cl3d_maxlayer")
121  else if (variable == "cl3d_seetot")
123  else if (variable == "cl3d_seemax")
125  else if (variable == "cl3d_spptot")
127  else if (variable == "cl3d_sppmax")
129  else if (variable == "cl3d_szz")
131  else if (variable == "cl3d_srrtot")
133  else if (variable == "cl3d_srrmax")
135  else if (variable == "cl3d_srrmean")
137  }
138 }
139 
141  std::map<std::string, float> inputs;
142  for (unsigned i = 0; i < input_variables_.size(); i++) {
144  }
145  float pt = cluster.pt();
146  float eta = cluster.eta();
147  int cat = category(pt, eta);
148  return (cat != -1 ? bdts_.at(cat)->evaluate(inputs) : -999.);
149 }
150 
152  float bdt_output = value(cluster);
153  float pt = cluster.pt();
154  float eta = cluster.eta();
155  int cat = category(pt, eta);
156  return (cat != -1 ? bdt_output > working_points_.at(cat) : true);
157 }
158 
160  for (unsigned cat = 0; cat < categories_.size(); cat++) {
161  if (categories_[cat].contains(pt, eta))
162  return static_cast<int>(cat);
163  }
164  return -1;
165 }
166 
168  const l1t::HGCalMulticluster& cluster) const {
169  switch (variable) {
171  return cluster.showerLength();
173  return cluster.coreShowerLength();
175  return cluster.firstLayer();
177  return cluster.maxLayer();
179  return cluster.sigmaEtaEtaTot();
181  return cluster.sigmaEtaEtaMax();
183  return cluster.sigmaPhiPhiTot();
185  return cluster.sigmaPhiPhiMax();
187  return cluster.sigmaZZ();
189  return cluster.sigmaRRTot();
191  return cluster.sigmaRRMax();
193  return cluster.sigmaRRMean();
194  default:
195  break;
196  }
197  return 0.;
198 }
int maxLayer() const
bool contains(EventRange const &lh, EventID const &rh)
Definition: EventRange.cc:37
double pt() const final
transverse momentum
float sigmaPhiPhiTot() const
float sigmaEtaEtaMax() const
float sigmaRRTot() const
float sigmaZZ() const
bool decision(const l1t::HGCalMulticluster &cluster) const final
int showerLength() const
def cat
Definition: eostools.py:401
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_HGC_TPG_CLUSTER_ID(type, name)
void initialize(const edm::ParameterSet &conf) final
Category(float pt_min, float pt_max, float eta_min, float eta_max)
float clusterVariable(ClusterVariable, const l1t::HGCalMulticluster &) const
float sigmaEtaEtaTot() const
int coreShowerLength() const
std::vector< std::unique_ptr< TMVAEvaluator > > bdts_
int firstLayer() const
float sigmaRRMean() const
float value(const l1t::HGCalMulticluster &cluster) const final
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float sigmaRRMax() const
Log< level::Warning, false > LogWarning
float sigmaPhiPhiMax() const
double eta() const final
momentum pseudorapidity