CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauMVAHelper.cc
Go to the documentation of this file.
3 
4 #include <boost/foreach.hpp>
5 #include <boost/bind.hpp>
6 #include <sstream>
7 
11 
14 
15 namespace reco { namespace tau {
16 
18  const std::string& eslabel,
19  const edm::ParameterSet& pluginOptions):
20  name_(name), eslabel_(eslabel), pluginOptions_(pluginOptions) {}
21 
23  const edm::EventSetup &es) {
24  // Update our MVA from the DB
27  if (eslabel_.size()) {
28  es.get<TauTagMVAComputerRcd>().get(eslabel_, handle);
29  } else {
30  es.get<TauTagMVAComputerRcd>().get(handle);
31  }
32  const MVAComputerContainer *container = handle.product();
33  // Load our MVA
34  bool reload = computer_.update(container, name_.c_str());
35  // If the MVA changed, update our list of discriminant plugins
36  if (reload && computer_.get())
37  loadDiscriminantPlugins(container->find(name_));
38  // Update the event info for all of our discriminators
39  BOOST_FOREACH(PluginMap::value_type plugin, plugins_) {
40  plugin.second->setup(evt, es);
41  }
42 }
43 
46  typedef std::vector<PhysicsTools::Calibration::Variable> VarList;
47  // List of input variables for this MVA.
48  const VarList &vars = comp.inputSet;
49  // Load the plugin for each of the Var if needed
50  BOOST_FOREACH(const VarList::value_type& var, vars) {
51  // Check to make sure it isn't a magic variable
52  if (std::strncmp(var.name.c_str(), "__", 2) != 0) {
53  // If we haven't added yet, build it.
54  PhysicsTools::AtomicId varId(var.name);
55  if (!plugins_.count(varId)) {
57  if (pluginOptions_.exists(var.name)) {
58  options = pluginOptions_.getParameter<edm::ParameterSet>(var.name);
59  };
60  // Make sure it has a name (required by base class)
61  if (!options.exists("name"))
62  options.addParameter("name", "MVA_" + var.name);
63  // Check if we want to specify the plugin name manually. This is
64  // required for things like the discriminant from discriminators, which
65  // take an InputTag. If we want to have more than one, we have to be
66  // able take the MVA name (like FlightPathSig) and map it to
67  // RecoTauDiscriminantFromDiscriminator[disc input tag = flight path sig]
68  // Otherwise we just keep our regular plugin mapping.
69  std::string pluginName = reco::tau::discPluginName(var.name);
70  if (options.exists("plugin")) {
71  pluginName = options.getParameter<std::string>("plugin");
72  }
73  plugins_.insert(
75  pluginName, options));
76  }
77  }
78  }
79 }
80 
82  // Loop over the relevant discriminators and the output
83  for (PluginMap::const_iterator plugin = plugins_.begin();
84  plugin != plugins_.end(); ++plugin) {
85  PhysicsTools::AtomicId id = plugin->first;
86  std::vector<double> pluginOutput = (plugin->second)->operator()(tau);
87  // Check for nans
88  for(size_t instance = 0; instance < pluginOutput.size(); ++instance) {
89  if (edm::isNotFinite(pluginOutput[instance])) {
90  std::ostringstream error;
91  error << "A nan was detected in"
92  << " the tau MVA variable " << id << " returning zero instead!"
93  << " The PFTau: " << *tau << std::endl;
94  tau->dump(error);
95  edm::LogError("CorruptedMVAInput") << error.str();
96  pluginOutput[instance] = 0.0;
97  }
98  }
99  //std::cout << "id: " << id << " first: " << pluginOutput[0] << std::endl;
100  // Build values and copy into values vector
101  std::for_each(pluginOutput.begin(), pluginOutput.end(),
103  boost::ref(values_), id, _1));
104  }
105 }
106 
107 // Get values
110  values_.clear();
111  fillValues(tau);
112  return values_;
113 }
114 
115 // Apply the MVA to a given tau
117  // Clear output
118  values_.clear();
119  // Build the values
120  fillValues(tau);
121  // Call the MVA
122  return computer_->eval(values_);
123 }
124 
126  double weight) const {
127  static const PhysicsTools::AtomicId kTargetId("__TARGET__");
128  static const PhysicsTools::AtomicId kWeightId("__WEIGHT__");
129  if (!computer_)
130  return;
131  values_.clear();
132  values_.add(kTargetId, target);
133  values_.add(kWeightId, weight);
134  // Build the discriminant values
135  fillValues(tau);
137 }
138 
139 }} // end namespace reco::tau
T getParameter(std::string const &) const
void setEvent(const edm::Event &evt, const edm::EventSetup &es)
double eval(Iterator_t first, Iterator_t last) const
evaluate variables given by a range of iterators given by first and last
auto_ptr< JetDefinition::Plugin > plugin
static PFTauRenderPlugin instance
void fillValues(const reco::PFTauRef &tau) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
Cheap generic unique keyword identifier class.
Definition: AtomicId.h:32
std::string discPluginName(const std::string &mvaName)
bool isNotFinite(T x)
Definition: isFinite.h:10
PhysicsTools::Variable::ValueList values_
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:145
tuple handle
Definition: patZpeak.py:22
double operator()(const PFTauRef &tau) const
edm::ParameterSet pluginOptions_
const PhysicsTools::Variable::ValueList & discriminants(const PFTauRef &tau) const
Container::value_type value_type
Helper class that can contain an list of identifier-value pairs.
Definition: Variable.h:82
RecoTauMVAHelper(const std::string &name, const std::string &eslabel, const edm::ParameterSet &pluginOptions)
bool update(const Calibration::MVAComputer *computer)
void train(const PFTauRef &tau, bool target, double weight=1.0) const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
std::vector< Variable > inputSet
Definition: MVAComputer.h:177
PhysicsTools::MVAComputerCache computer_
int weight
Definition: histoStyle.py:50
void loadDiscriminantPlugins(const PhysicsTools::Calibration::MVAComputer &computer)
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:56
void add(AtomicId id, double value)
Definition: Variable.h:107