CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TMVAEvaluator.cc
Go to the documentation of this file.
2 
5 
6 
8  mIsInitialized(false)
9 {
10 }
11 
12 
14 {
15 }
16 
17 
18 void TMVAEvaluator::initialize(const std::string & options, const std::string & method, const std::string & weightFile,
19  const std::vector<std::string> & variables, const std::vector<std::string> & spectators)
20 {
21  // initialize the TMVA reader
22  mReader.reset(new TMVA::Reader(options.c_str()));
23  mReader->SetVerbose(false);
24  mMethod = method;
25 
26  // add input variables
27  for(std::vector<std::string>::const_iterator it = variables.begin(); it!=variables.end(); ++it)
28  {
29  mVariables.insert( std::pair<std::string,float>(*it,0.) );
30  mReader->AddVariable(it->c_str(), &mVariables.at(*it));
31  }
32 
33  // add spectator variables
34  for(std::vector<std::string>::const_iterator it = spectators.begin(); it!=spectators.end(); ++it)
35  {
36  mSpectators.insert( std::pair<std::string,float>(*it,0.) );
37  mReader->AddSpectator(it->c_str(), &mSpectators.at(*it));
38  }
39 
40  // load the TMVA weights
41  reco::details::loadTMVAWeights(mReader.get(), mMethod.c_str(), weightFile.c_str());
42 
43  mIsInitialized = true;
44 }
45 
46 
47 float TMVAEvaluator::evaluate(const std::map<std::string,float> & inputs, const bool useSpectators)
48 {
49  if(!mIsInitialized)
50  {
51  edm::LogError("InitializationError") << "TMVAEvaluator not properly initialized.";
52  return -99.;
53  }
54 
55  if( useSpectators && inputs.size() < ( mVariables.size() + mSpectators.size() ) )
56  {
57  edm::LogError("MissingInputs") << "Too few inputs provided (" << inputs.size() << " provided but " << mVariables.size() << " input and " << mSpectators.size() << " spectator variables expected).";
58  return -99.;
59  }
60  else if( inputs.size() < mVariables.size() )
61  {
62  edm::LogError("MissingInputVariable(s)") << "Too few input variables provided (" << inputs.size() << " provided but " << mVariables.size() << " expected).";
63  return -99.;
64  }
65 
66  // set the input variable values
67  for(std::map<std::string,float>::iterator it = mVariables.begin(); it!=mVariables.end(); ++it)
68  {
69  if (inputs.count(it->first)>0)
70  it->second = inputs.at(it->first);
71  else
72  edm::LogError("MissingInputVariable") << "Input variable " << it->first << " is missing from the list of inputs. The returned discriminator value might not be sensible.";
73  }
74 
75  // if using spectator variables
76  if(useSpectators)
77  {
78  // set the spectator variable values
79  for(std::map<std::string,float>::iterator it = mSpectators.begin(); it!=mSpectators.end(); ++it)
80  {
81  if (inputs.count(it->first)>0)
82  it->second = inputs.at(it->first);
83  else
84  edm::LogError("MissingSpectatorVariable") << "Spectator variable " << it->first << " is missing from the list of inputs. The returned discriminator value might not be sensible.";
85  }
86  }
87 
88  // evaluate the MVA
89  float value = mReader->EvaluateMVA(mMethod.c_str());
90 
91  return value;
92 }
std::map< std::string, float > mVariables
Definition: TMVAEvaluator.h:28
float evaluate(const std::map< std::string, float > &inputs, const bool useSpectators=false)
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)
std::string mMethod
Definition: TMVAEvaluator.h:25
std::unique_ptr< TMVA::Reader > mReader
Definition: TMVAEvaluator.h:26
volatile std::atomic< bool > shutdown_flag false
void loadTMVAWeights(TMVA::Reader *reader, const std::string &method, const std::string &weightFile, bool verbose=false)
std::map< std::string, float > mSpectators
Definition: TMVAEvaluator.h:29