CMS 3D CMS Logo

MVAValueMapProducer.h
Go to the documentation of this file.
1 #ifndef __RecoEgamma_EgammaTools_MVAValueMapProducer_H__
2 #define __RecoEgamma_EgammaTools_MVAValueMapProducer_H__
3 
15 
16 #include <memory>
17 #include <vector>
18 #include <cmath>
19 
20 template <class ParticleType>
22 
23  public:
24 
26  ~MVAValueMapProducer() override {}
27 
28  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
29 
30  private:
31 
32  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
33 
34  // for AOD and MiniAOD case
36 
37  // MVA estimators
38  const std::vector<std::unique_ptr<AnyMVAEstimatorRun2Base>> mvaEstimators_;
39 
40  // Value map names
41  const std::vector<std::string> mvaValueMapNames_;
42  const std::vector<std::string> mvaRawValueMapNames_;
43  const std::vector<std::string> mvaCategoriesMapNames_;
44 
45  // To get the auxiliary MVA variables
47 
48 };
49 
50 namespace {
51 
52  auto getMVAEstimators(const edm::VParameterSet& vConfig)
53  {
54  std::vector<std::unique_ptr<AnyMVAEstimatorRun2Base>> mvaEstimators;
55 
56  // Loop over the list of MVA configurations passed here from python and
57  // construct all requested MVA estimators.
58  for( auto &imva : vConfig )
59  {
60 
61  // The factory below constructs the MVA of the appropriate type based
62  // on the "mvaName" which is the name of the derived MVA class (plugin)
63  if( !imva.empty() ) {
64 
65  mvaEstimators.emplace_back(AnyMVAEstimatorRun2Factory::get()->create(
66  imva.getParameter<std::string>("mvaName"), imva));
67 
68  } else
69  throw cms::Exception(" MVA configuration not found: ")
70  << " failed to find proper configuration for one of the MVAs in the main python script " << std::endl;
71  }
72 
73  return mvaEstimators;
74  }
75 
76  std::vector<std::string> getValueMapNames(const edm::VParameterSet& vConfig, std::string && suffix)
77  {
78  std::vector<std::string> names;
79  for( auto &imva : vConfig )
80  {
81  names.push_back(imva.getParameter<std::string>("mvaName") + imva.getParameter<std::string>("mvaTag") + suffix);
82  }
83 
84  return names;
85  }
86 }
87 
88 template <class ParticleType>
90  : src_(consumesCollector(), iConfig, "src", "srcMiniAOD")
91  , mvaEstimators_(getMVAEstimators(iConfig.getParameterSetVector("mvaConfigurations")))
92  , mvaValueMapNames_ (getValueMapNames(iConfig.getParameterSetVector("mvaConfigurations"), "Values" ))
93  , mvaRawValueMapNames_ (getValueMapNames(iConfig.getParameterSetVector("mvaConfigurations"), "RawValues" ))
94  , mvaCategoriesMapNames_(getValueMapNames(iConfig.getParameterSetVector("mvaConfigurations"), "Categories"))
96 {
97  for( auto const& name : mvaValueMapNames_ ) produces<edm::ValueMap<float>>(name);
98  for( auto const& name : mvaRawValueMapNames_ ) produces<edm::ValueMap<float>>(name);
99  for( auto const& name : mvaCategoriesMapNames_ ) produces<edm::ValueMap<int >>(name);
100 }
101 
102 template <class ParticleType>
104 {
105 
107 
108  // Loop over MVA estimators
109  for( unsigned iEstimator = 0; iEstimator < mvaEstimators_.size(); iEstimator++ ){
110 
111  std::vector<float> mvaValues;
112  std::vector<float> mvaRawValues;
113  std::vector<int> mvaCategories;
114 
115  // Loop over particles
116  for (size_t i = 0; i < src->size(); ++i)
117  {
118  auto iCand = src->ptrAt(i);
119 
120  std::vector<float> auxVariables = variableHelper_.getAuxVariables(iCand, iEvent);
121 
122  int cat = -1; // Passed by reference to the mvaValue function to store the category
123  const float response = mvaEstimators_[iEstimator]->mvaValue( &(*iCand), auxVariables, cat );
124  mvaRawValues.push_back( response ); // The MVA score
125  mvaValues.push_back( 2.0/(1.0+exp(-2.0*response))-1 ); // MVA output between -1 and 1
126  mvaCategories.push_back( cat );
127  } // end loop over particles
128 
129  writeValueMap(iEvent, src, mvaValues , mvaValueMapNames_ [iEstimator] );
130  writeValueMap(iEvent, src, mvaRawValues , mvaRawValueMapNames_ [iEstimator] );
131  writeValueMap(iEvent, src, mvaCategories, mvaCategoriesMapNames_[iEstimator] );
132 
133  } // end loop over estimators
134 
135 }
136 
137 template <class ParticleType>
139  //The following says we do not know what parameters are allowed so do no validation
140  // Please change this to state exactly what you do use, even if it is no parameters
142  desc.setUnknown();
143  descriptions.addDefault(desc);
144 }
145 
146 #endif
const std::vector< std::string > mvaRawValueMapNames_
def create(alignables, pedeDump, additionalData, outputFile, config)
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
const std::vector< std::unique_ptr< AnyMVAEstimatorRun2Base > > mvaEstimators_
const std::string names[nVars_]
void writeValueMap(edm::Event &iEvent, const edm::Handle< HandleType > &handle, const std::vector< ValueType > &values, const std::string &label)
Definition: Utils.h:13
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
def cat(path)
Definition: eostools.py:401
void addDefault(ParameterSetDescription const &psetDescription)
const MultiTokenT< edm::View< ParticleType > > src_
const std::vector< std::string > mvaCategoriesMapNames_
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const std::vector< std::string > mvaValueMapNames_
edm::Handle< T > getValidHandle(const edm::Event &iEvent) const
Definition: MultiToken.h:95
MVAValueMapProducer(const edm::ParameterSet &)
const MVAVariableHelper< ParticleType > variableHelper_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
T get(const Candidate &c)
Definition: component.h:55