CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MVAValueMapProducer.h
Go to the documentation of this file.
3 
6 
8 
11 
13 
14 #include <memory>
15 #include <vector>
16 
17 template <class ParticleType>
19 
20  public:
21 
22  explicit MVAValueMapProducer(const edm::ParameterSet&);
24 
25  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
26 
27  private:
28 
29  virtual void produce(edm::Event&, const edm::EventSetup&) override;
30 
31  template<typename T>
34  const std::vector<T> & values,
35  const std::string & label) const ;
36 
37  // for AOD case
39 
40  // for miniAOD case
42 
43  // MVA estimator
44  std::vector<std::unique_ptr<AnyMVAEstimatorRun2Base>> mvaEstimators_;
45 
46  // Value map names
47  std::vector <std::string> mvaValueMapNames_;
48  std::vector <std::string> mvaCategoriesMapNames_;
49 
50 };
51 
52 template <class ParticleType>
54 {
55 
56  //
57  // Declare consummables, handle both AOD and miniAOD case
58  //
59  src_ = mayConsume<edm::View<ParticleType> >(iConfig.getParameter<edm::InputTag>("src"));
60  srcMiniAOD_ = mayConsume<edm::View<ParticleType> >(iConfig.getParameter<edm::InputTag>("srcMiniAOD"));
61 
62  // Loop over the list of MVA configurations passed here from python and
63  // construct all requested MVA esimtators.
64  const std::vector<edm::ParameterSet>& mvaEstimatorConfigs
65  = iConfig.getParameterSetVector("mvaConfigurations");
66  for( auto &imva : mvaEstimatorConfigs ){
67 
68  std::unique_ptr<AnyMVAEstimatorRun2Base> thisEstimator;
69  thisEstimator.reset(NULL);
70  if( !imva.empty() ) {
71  const std::string& pName = imva.getParameter<std::string>("mvaName");
72  // The factory below constructs the MVA of the appropriate type based
73  // on the "mvaName" which is the name of the derived MVA class (plugin)
75  // Declare all event content, such as ValueMaps produced upstream or other,
76  // original event data pieces, that is needed (if any is implemented in the specific
77  // MVA classes)
78  //edm::ConsumesCollector &cc = consumesCollector();
79  estimator->setConsumes( consumesCollector() );
80 
81  thisEstimator.reset(estimator);
82 
83  } else
84  throw cms::Exception(" MVA configuration not found: ")
85  << " failed to find proper configuration for one of the MVAs in the main python script " << std::endl;
86 
87  // The unique pointer control is passed to the vector in the line below.
88  // Don't use thisEstimator pointer beyond the next line.
89  mvaEstimators_.emplace_back( thisEstimator.release() );
90 
91  //
92  // Compose and save the names of the value maps to be produced
93  //
94  const auto& currentEstimator = mvaEstimators_.back();
95  std::string thisValueMapName = currentEstimator->getName() + "Values";
96  std::string thisCategoriesMapName = currentEstimator->getName() + "Categories";
97  mvaValueMapNames_.push_back( thisValueMapName );
98  mvaCategoriesMapNames_.push_back( thisCategoriesMapName );
99 
100  // Declare the maps to the framework
101  produces<edm::ValueMap<float> >(thisValueMapName);
102  produces<edm::ValueMap<int> >(thisCategoriesMapName);
103 
104  }
105 
106 
107 }
108 
109 template <class ParticleType>
111 }
112 
113 template <class ParticleType>
115 
116  using namespace edm;
117 
119 
120  // Retrieve the collection of particles from the event.
121  // If we fail to retrieve the collection with the standard AOD
122  // name, we next look for the one with the stndard miniAOD name.
123  iEvent.getByToken(src_, src);
124  if( !src.isValid() ){
125  iEvent.getByToken(srcMiniAOD_,src);
126  if( !src.isValid() )
127  throw cms::Exception(" Collection not found: ")
128  << " failed to find a standard AOD or miniAOD particle collection " << std::endl;
129  }
130 
131 
132  // Loop over MVA estimators
133  for( unsigned iEstimator = 0; iEstimator < mvaEstimators_.size(); iEstimator++ ){
134 
135  // Set up all event content, such as ValueMaps produced upstream or other,
136  // original event data pieces, that is needed (if any is implemented in the specific
137  // MVA classes)
138  mvaEstimators_[iEstimator]->getEventContent( iEvent );
139 
140  std::vector<float> mvaValues;
141  std::vector<int> mvaCategories;
142 
143  // Loop over particles
144  for (size_t i = 0; i < src->size(); ++i){
145  auto iCand = src->ptrAt(i);
146 
147  mvaValues.push_back( mvaEstimators_[iEstimator]->mvaValue( iCand ) );
148  mvaCategories.push_back( mvaEstimators_[iEstimator]->findCategory( iCand ) );
149  } // end loop over particles
150 
151  writeValueMap(iEvent, src, mvaValues, mvaValueMapNames_[iEstimator] );
152  writeValueMap(iEvent, src, mvaCategories, mvaCategoriesMapNames_[iEstimator] );
153 
154  } // end loop over estimators
155 
156 
157 }
158 
159 template<class ParticleType> template<typename T>
162  const std::vector<T> & values,
163  const std::string & label) const
164 {
165  using namespace edm;
166  using namespace std;
167  auto_ptr<ValueMap<T> > valMap(new ValueMap<T>());
168  typename edm::ValueMap<T>::Filler filler(*valMap);
169  filler.insert(handle, values.begin(), values.end());
170  filler.fill();
171  iEvent.put(valMap, label);
172 }
173 
174 template <class ParticleType>
176  //The following says we do not know what parameters are allowed so do no validation
177  // Please change this to state exactly what you do use, even if it is no parameters
179  desc.setUnknown();
180  descriptions.addDefault(desc);
181 }
182 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
VParameterSet const & getParameterSetVector(std::string const &name) const
std::vector< std::string > mvaValueMapNames_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
#define NULL
Definition: scimark2.h:8
std::vector< std::unique_ptr< AnyMVAEstimatorRun2Base > > mvaEstimators_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
virtual void produce(edm::Event &, const edm::EventSetup &) override
tuple handle
Definition: patZpeak.py:22
void writeValueMap(edm::Event &iEvent, const edm::Handle< edm::View< ParticleType > > &handle, const std::vector< T > &values, const std::string &label) const
std::vector< std::string > mvaCategoriesMapNames_
MVAValueMapProducer(const edm::ParameterSet &)
edm::EDGetToken srcMiniAOD_
T get(const Candidate &c)
Definition: component.h:55