CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
reco::tau::RecoTauMVAHelper Class Reference

#include <RecoTauMVAHelper.h>

Public Member Functions

const
PhysicsTools::Variable::ValueList
discriminants (const PFTauRef &tau) const
 
double operator() (const PFTauRef &tau) const
 
 RecoTauMVAHelper (const std::string &name, const std::string &eslabel, const edm::ParameterSet &pluginOptions)
 
void setEvent (const edm::Event &evt, const edm::EventSetup &es)
 
void train (const PFTauRef &tau, bool target, double weight=1.0) const
 
 ~RecoTauMVAHelper ()
 

Private Types

typedef boost::ptr_map
< PhysicsTools::AtomicId,
RecoTauDiscriminantPlugin
PluginMap
 

Private Member Functions

void fillValues (const reco::PFTauRef &tau) const
 
void loadDiscriminantPlugins (const PhysicsTools::Calibration::MVAComputer &computer)
 

Private Attributes

PhysicsTools::MVAComputerCache computer_
 
std::string eslabel_
 
std::string name_
 
edm::ParameterSet pluginOptions_
 
PluginMap plugins_
 
PhysicsTools::Variable::ValueList values_
 

Detailed Description

Definition at line 40 of file RecoTauMVAHelper.h.

Member Typedef Documentation

Definition at line 70 of file RecoTauMVAHelper.h.

Constructor & Destructor Documentation

reco::tau::RecoTauMVAHelper::RecoTauMVAHelper ( const std::string &  name,
const std::string &  eslabel,
const edm::ParameterSet pluginOptions 
)
explicit

Definition at line 16 of file RecoTauMVAHelper.cc.

reco::tau::RecoTauMVAHelper::~RecoTauMVAHelper ( )
inline

Definition at line 45 of file RecoTauMVAHelper.h.

45 {}

Member Function Documentation

const PhysicsTools::Variable::ValueList & reco::tau::RecoTauMVAHelper::discriminants ( const PFTauRef tau) const

Definition at line 108 of file RecoTauMVAHelper.cc.

References PhysicsTools::Variable::ValueList::clear(), fillValues(), and values_.

108  {
109  values_.clear();
110  fillValues(tau);
111  return values_;
112 }
void fillValues(const reco::PFTauRef &tau) const
PhysicsTools::Variable::ValueList values_
void reco::tau::RecoTauMVAHelper::fillValues ( const reco::PFTauRef tau) const
private

Definition at line 80 of file RecoTauMVAHelper.cc.

References PhysicsTools::Variable::ValueList::add(), error, instance, edm::detail::isnan(), pfTaus_cff::plugin, plugins_, metsig::tau, and values_.

Referenced by discriminants(), operator()(), and train().

80  {
81  // Loop over the relevant discriminators and the output
82  for (PluginMap::const_iterator plugin = plugins_.begin();
83  plugin != plugins_.end(); ++plugin) {
84  PhysicsTools::AtomicId id = plugin->first;
85  std::vector<double> pluginOutput = (plugin->second)->operator()(tau);
86  // Check for nans
87  for(size_t instance = 0; instance < pluginOutput.size(); ++instance) {
88  if (std::isnan(pluginOutput[instance])) {
89  std::ostringstream error;
90  error << "A nan was detected in"
91  << " the tau MVA variable " << id << " returning zero instead!"
92  << " The PFTau: " << *tau << std::endl;
93  tau->dump(error);
94  edm::LogError("CorruptedMVAInput") << error.str();
95  pluginOutput[instance] = 0.0;
96  }
97  }
98  //std::cout << "id: " << id << " first: " << pluginOutput[0] << std::endl;
99  // Build values and copy into values vector
100  std::for_each(pluginOutput.begin(), pluginOutput.end(),
102  boost::ref(values_), id, _1));
103  }
104 }
static PFTauRenderPlugin instance
Cheap generic unique keyword identifier class.
Definition: AtomicId.h:32
tuple plugin
Definition: pfTaus_cff.py:103
PhysicsTools::Variable::ValueList values_
bool isnan(float x)
Definition: math.h:13
void add(AtomicId id, double value)
Definition: Variable.h:107
void reco::tau::RecoTauMVAHelper::loadDiscriminantPlugins ( const PhysicsTools::Calibration::MVAComputer computer)
private

Definition at line 43 of file RecoTauMVAHelper.cc.

References edm::ParameterSet::addParameter(), SurfaceDeformationFactory::create(), reco::tau::discPluginName(), edm::ParameterSet::exists(), reco::get(), edm::ParameterSet::getParameter(), PhysicsTools::Calibration::MVAComputer::inputSet, AlCaHLTBitMon_ParallelJobs::options, pluginOptions_, and plugins_.

Referenced by setEvent().

44  {
45  typedef std::vector<PhysicsTools::Calibration::Variable> VarList;
46  // List of input variables for this MVA.
47  const VarList &vars = comp.inputSet;
48  // Load the plugin for each of the Var if needed
49  BOOST_FOREACH(const VarList::value_type& var, vars) {
50  // Check to make sure it isn't a magic variable
51  if (std::strncmp(var.name.c_str(), "__", 2) != 0) {
52  // If we haven't added yet, build it.
53  PhysicsTools::AtomicId varId(var.name);
54  if (!plugins_.count(varId)) {
56  if (pluginOptions_.exists(var.name)) {
57  options = pluginOptions_.getParameter<edm::ParameterSet>(var.name);
58  };
59  // Make sure it has a name (required by base class)
60  if (!options.exists("name"))
61  options.addParameter("name", "MVA_" + var.name);
62  // Check if we want to specify the plugin name manually. This is
63  // required for things like the discriminant from discriminators, which
64  // take an InputTag. If we want to have more than one, we have to be
65  // able take the MVA name (like FlightPathSig) and map it to
66  // RecoTauDiscriminantFromDiscriminator[disc input tag = flight path sig]
67  // Otherwise we just keep our regular plugin mapping.
68  std::string pluginName = reco::tau::discPluginName(var.name);
69  if (options.exists("plugin")) {
70  pluginName = options.getParameter<std::string>("plugin");
71  }
72  plugins_.insert(
74  pluginName, options));
75  }
76  }
77  }
78 }
T getParameter(std::string const &) 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)
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:145
edm::ParameterSet pluginOptions_
Container::value_type value_type
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:56
double reco::tau::RecoTauMVAHelper::operator() ( const PFTauRef tau) const

Definition at line 115 of file RecoTauMVAHelper.cc.

References PhysicsTools::Variable::ValueList::clear(), computer_, PhysicsTools::MVAComputer::eval(), fillValues(), and values_.

115  {
116  // Clear output
117  values_.clear();
118  // Build the values
119  fillValues(tau);
120  // Call the MVA
121  return computer_->eval(values_);
122 }
double eval(Iterator_t first, Iterator_t last) const
evaluate variables given by a range of iterators given by first and last
void fillValues(const reco::PFTauRef &tau) const
PhysicsTools::Variable::ValueList values_
PhysicsTools::MVAComputerCache computer_
void reco::tau::RecoTauMVAHelper::setEvent ( const edm::Event evt,
const edm::EventSetup es 
)

Definition at line 21 of file RecoTauMVAHelper.cc.

References computer_, eslabel_, edm::EventSetup::get(), PhysicsTools::MVAComputerCache::get(), patZpeak::handle, loadDiscriminantPlugins(), name_, pfTaus_cff::plugin, plugins_, edm::ESHandle< class >::product(), and PhysicsTools::MVAComputerCache::update().

Referenced by RecoTauMVATrainer::analyze().

22  {
23  // Update our MVA from the DB
26  if (eslabel_.size()) {
27  es.get<TauTagMVAComputerRcd>().get(eslabel_, handle);
28  } else {
29  es.get<TauTagMVAComputerRcd>().get(handle);
30  }
31  const MVAComputerContainer *container = handle.product();
32  // Load our MVA
33  bool reload = computer_.update(container, name_.c_str());
34  // If the MVA changed, update our list of discriminant plugins
35  if (reload && computer_.get())
36  loadDiscriminantPlugins(container->find(name_));
37  // Update the event info for all of our discriminators
38  BOOST_FOREACH(PluginMap::value_type plugin, plugins_) {
39  plugin.second->setup(evt, es);
40  }
41 }
tuple plugin
Definition: pfTaus_cff.py:103
tuple handle
Definition: patZpeak.py:22
Container::value_type value_type
bool update(const Calibration::MVAComputer *computer)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
PhysicsTools::MVAComputerCache computer_
void loadDiscriminantPlugins(const PhysicsTools::Calibration::MVAComputer &computer)
void reco::tau::RecoTauMVAHelper::train ( const PFTauRef tau,
bool  target,
double  weight = 1.0 
) const

Definition at line 124 of file RecoTauMVAHelper.cc.

References PhysicsTools::Variable::ValueList::add(), PhysicsTools::Variable::ValueList::clear(), computer_, PhysicsTools::MVAComputer::eval(), fillValues(), and values_.

125  {
126  static const PhysicsTools::AtomicId kTargetId("__TARGET__");
127  static const PhysicsTools::AtomicId kWeightId("__WEIGHT__");
128  if (!computer_)
129  return;
130  values_.clear();
131  values_.add(kTargetId, target);
132  values_.add(kWeightId, weight);
133  // Build the discriminant values
134  fillValues(tau);
136 }
double eval(Iterator_t first, Iterator_t last) const
evaluate variables given by a range of iterators given by first and last
void fillValues(const reco::PFTauRef &tau) const
Cheap generic unique keyword identifier class.
Definition: AtomicId.h:32
PhysicsTools::Variable::ValueList values_
PhysicsTools::MVAComputerCache computer_
void add(AtomicId id, double value)
Definition: Variable.h:107

Member Data Documentation

PhysicsTools::MVAComputerCache reco::tau::RecoTauMVAHelper::computer_
private

Definition at line 67 of file RecoTauMVAHelper.h.

Referenced by operator()(), setEvent(), and train().

std::string reco::tau::RecoTauMVAHelper::eslabel_
private

Definition at line 63 of file RecoTauMVAHelper.h.

Referenced by setEvent().

std::string reco::tau::RecoTauMVAHelper::name_
private

Definition at line 61 of file RecoTauMVAHelper.h.

Referenced by setEvent().

edm::ParameterSet reco::tau::RecoTauMVAHelper::pluginOptions_
private

Definition at line 65 of file RecoTauMVAHelper.h.

Referenced by loadDiscriminantPlugins().

PluginMap reco::tau::RecoTauMVAHelper::plugins_
private

Definition at line 71 of file RecoTauMVAHelper.h.

Referenced by fillValues(), loadDiscriminantPlugins(), and setEvent().

PhysicsTools::Variable::ValueList reco::tau::RecoTauMVAHelper::values_
mutableprivate

Definition at line 77 of file RecoTauMVAHelper.h.

Referenced by discriminants(), fillValues(), operator()(), and train().