CMS 3D CMS Logo

TMVAEvaluator.cc
Go to the documentation of this file.
4 
8 
9 TMVAEvaluator::TMVAEvaluator() : mIsInitialized(false), mUsingGBRForest(false), mUseAdaBoost(false) {}
10 
12  const std::string& method,
13  const std::string& weightFile,
14  const std::vector<std::string>& variables,
15  const std::vector<std::string>& spectators,
16  bool useGBRForest,
17  bool useAdaBoost) {
18  // initialize the TMVA reader
19  mReader.reset(new TMVA::Reader(options.c_str()));
20  mReader->SetVerbose(false);
21  mMethod = method;
22 
23  // add input variables
24  for (std::vector<std::string>::const_iterator it = variables.begin(); it != variables.end(); ++it) {
25  mVariables.insert(std::make_pair(*it, std::make_pair(it - variables.begin(), 0.)));
26  mReader->AddVariable(it->c_str(), &(mVariables.at(*it).second));
27  }
28 
29  // add spectator variables
30  for (std::vector<std::string>::const_iterator it = spectators.begin(); it != spectators.end(); ++it) {
31  mSpectators.insert(std::make_pair(*it, std::make_pair(it - spectators.begin(), 0.)));
32  mReader->AddSpectator(it->c_str(), &(mSpectators.at(*it).second));
33  }
34 
35  // load the TMVA weights
37 
38  if (useGBRForest) {
39  mGBRForest = createGBRForest(weightFile);
40 
41  // now can free some memory
42  mReader.reset(nullptr);
43 
44  mUsingGBRForest = true;
46  }
47 
48  mIsInitialized = true;
49 }
50 
52  const std::vector<std::string>& variables,
53  const std::vector<std::string>& spectators,
54  bool useAdaBoost) {
55  // add input variables
56  for (std::vector<std::string>::const_iterator it = variables.begin(); it != variables.end(); ++it)
57  mVariables.insert(std::make_pair(*it, std::make_pair(it - variables.begin(), 0.)));
58 
59  // add spectator variables
60  for (std::vector<std::string>::const_iterator it = spectators.begin(); it != spectators.end(); ++it)
61  mSpectators.insert(std::make_pair(*it, std::make_pair(it - spectators.begin(), 0.)));
62 
63  // do not take ownership if getting GBRForest from an external source
64  mGBRForest = std::shared_ptr<const GBRForest>(gbrForest, [](const GBRForest*) {});
65 
66  mIsInitialized = true;
67  mUsingGBRForest = true;
69 }
70 
72  const std::string& label,
73  const std::vector<std::string>& variables,
74  const std::vector<std::string>& spectators,
75  bool useAdaBoost) {
76  edm::ESHandle<GBRForest> gbrForestHandle;
77 
78  iSetup.get<GBRWrapperRcd>().get(label.c_str(), gbrForestHandle);
79 
81 }
82 
83 float TMVAEvaluator::evaluateTMVA(const std::map<std::string, float>& inputs, bool useSpectators) const {
84  // default value
85  float value = -99.;
86 
87  // TMVA::Reader is not thread safe
88  std::lock_guard<std::mutex> lock(m_mutex);
89 
90  // set the input variable values
91  for (auto it = mVariables.begin(); it != mVariables.end(); ++it) {
92  if (inputs.count(it->first) > 0)
93  it->second.second = inputs.at(it->first);
94  else
95  edm::LogError("MissingInputVariable")
96  << "Input variable " << it->first
97  << " is missing from the list of inputs. The returned discriminator value might not be sensible.";
98  }
99 
100  // if using spectator variables
101  if (useSpectators) {
102  // set the spectator variable values
103  for (auto it = mSpectators.begin(); it != mSpectators.end(); ++it) {
104  if (inputs.count(it->first) > 0)
105  it->second.second = inputs.at(it->first);
106  else
107  edm::LogError("MissingSpectatorVariable")
108  << "Spectator variable " << it->first
109  << " is missing from the list of inputs. The returned discriminator value might not be sensible.";
110  }
111  }
112 
113  // evaluate the MVA
114  value = mReader->EvaluateMVA(mMethod.c_str());
115 
116  return value;
117 }
118 
119 float TMVAEvaluator::evaluateGBRForest(const std::map<std::string, float>& inputs) const {
120  // default value
121  float value = -99.;
122 
123  std::unique_ptr<float[]> vars(new float[mVariables.size()]); // allocate n floats
124 
125  // set the input variable values
126  for (auto it = mVariables.begin(); it != mVariables.end(); ++it) {
127  if (inputs.count(it->first) > 0)
128  vars[it->second.first] = inputs.at(it->first);
129  else
130  edm::LogError("MissingInputVariable")
131  << "Input variable " << it->first
132  << " is missing from the list of inputs. The returned discriminator value might not be sensible.";
133  }
134 
135  // evaluate the MVA
136  if (mUseAdaBoost)
137  value = mGBRForest->GetAdaBoostClassifier(vars.get());
138  else
139  value = mGBRForest->GetGradBoostClassifier(vars.get());
140 
141  return value;
142 }
143 
144 float TMVAEvaluator::evaluate(const std::map<std::string, float>& inputs, bool useSpectators) const {
145  // default value
146  float value = -99.;
147 
148  if (!mIsInitialized) {
149  edm::LogError("InitializationError") << "TMVAEvaluator not properly initialized.";
150  return value;
151  }
152 
153  if (useSpectators && inputs.size() < (mVariables.size() + mSpectators.size())) {
154  edm::LogError("MissingInputs") << "Too few inputs provided (" << inputs.size() << " provided but "
155  << mVariables.size() << " input and " << mSpectators.size()
156  << " spectator variables expected).";
157  return value;
158  } else if (inputs.size() < mVariables.size()) {
159  edm::LogError("MissingInputVariable(s)") << "Too few input variables provided (" << inputs.size()
160  << " provided but " << mVariables.size() << " expected).";
161  return value;
162  }
163 
164  if (mUsingGBRForest) {
165  if (useSpectators)
166  edm::LogWarning("UnsupportedFunctionality")
167  << "Use of spectator variables with GBRForest is not supported. Spectator variables will be ignored.";
168  value = evaluateGBRForest(inputs);
169  } else
170  value = evaluateTMVA(inputs, useSpectators);
171 
172  return value;
173 }
void initializeGBRForest(const GBRForest *gbrForest, const std::vector< std::string > &variables, const std::vector< std::string > &spectators, bool useAdaBoost=false)
std::map< std::string, std::pair< size_t, float > > mVariables
Definition: TMVAEvaluator.h:53
bool mUsingGBRForest
Definition: TMVAEvaluator.h:45
void initialize(const std::string &options, const std::string &method, const std::string &weightFile, const std::vector< std::string > &variables, const std::vector< std::string > &spectators, bool useGBRForest=false, bool useAdaBoost=false)
char const * label
float evaluateTMVA(const std::map< std::string, float > &inputs, bool useSpectators) const
std::shared_ptr< const GBRForest > mGBRForest
Definition: TMVAEvaluator.h:51
Definition: value.py:1
std::string mMethod
Definition: TMVAEvaluator.h:48
std::map< std::string, std::pair< size_t, float > > mSpectators
Definition: TMVAEvaluator.h:54
std::mutex m_mutex
Definition: TMVAEvaluator.h:49
TMVA::IMethod * loadTMVAWeights(TMVA::Reader *reader, const std::string &method, const std::string &weightFile, bool verbose=false)
std::unique_ptr< TMVA::Reader > mReader
Definition: TMVAEvaluator.h:50
float evaluate(const std::map< std::string, float > &inputs, bool useSpectators=false) const
T get() const
Definition: EventSetup.h:71
vars
Definition: DeepTauId.cc:77
float evaluateGBRForest(const std::map< std::string, float > &inputs) const
T const * product() const
Definition: ESHandle.h:86
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)