CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PhotonMVAEstimator.cc
Go to the documentation of this file.
14 
16 public:
17  // Constructor and destructor
19  ~PhotonMVAEstimator() override{};
20 
21  // Calculation of the MVA value
22  float mvaValue(const reco::Candidate* candPtr, std::vector<float> const& auxVars, int& iCategory) const override;
23 
24  int findCategory(const reco::Candidate* candPtr) const override;
25 
26  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
27 
28 private:
29  int findCategory(reco::Photon const& photon) const;
30 
31  // The number of categories and number of variables per category
33  std::vector<ThreadSafeFunctor<StringCutObjectSelector<reco::Photon>>> categoryFunctions_;
34  std::vector<int> nVariables_;
35 
36  // Data members
37  std::vector<std::unique_ptr<const GBRForest>> gbrForests_;
38 
39  // There might be different variables for each category, so the variables
40  // names vector is itself a vector of length nCategories
41  std::vector<std::vector<int>> variables_;
42 
43  // The variable manager which stores how to obtain the variables
45 
46  // Other objects needed by the MVA
47  std::unique_ptr<EffectiveAreas> effectiveAreas_;
48  std::vector<double> phoIsoPtScalingCoeff_;
49  double phoIsoCutoff_;
50 };
51 
54  mvaVarMngr_(conf.getParameter<std::string>("variableDefinition"), MVAVariableHelper::indexMap()) {
55  //
56  // Construct the MVA estimators
57  //
58  if (getTag() == "Run2Spring16NonTrigV1") {
60  std::make_unique<EffectiveAreas>((conf.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath());
61  phoIsoPtScalingCoeff_ = conf.getParameter<std::vector<double>>("phoIsoPtScalingCoeff");
62  phoIsoCutoff_ = conf.getParameter<double>("phoIsoCutoff");
63  }
64 
65  const auto weightFileNames = conf.getParameter<std::vector<std::string>>("weightFileNames");
66  const auto categoryCutStrings = conf.getParameter<std::vector<std::string>>("categoryCuts");
67 
68  if ((int)(categoryCutStrings.size()) != getNCategories())
69  throw cms::Exception("MVA config failure: ")
70  << "wrong number of category cuts in PhotonMVAEstimator" << getTag() << std::endl;
71 
72  for (auto const& cut : categoryCutStrings)
73  categoryFunctions_.emplace_back(cut);
74 
75  // Initialize GBRForests
76  if (static_cast<int>(weightFileNames.size()) != getNCategories())
77  throw cms::Exception("MVA config failure: ")
78  << "wrong number of weightfiles in PhotonMVAEstimator" << getTag() << std::endl;
79 
80  gbrForests_.clear();
81  // Create a TMVA reader object for each category
82  for (int i = 0; i < getNCategories(); i++) {
83  std::vector<int> variablesInCategory;
84 
85  std::vector<std::string> variableNamesInCategory;
86  gbrForests_.push_back(createGBRForest(weightFileNames[i], variableNamesInCategory));
87 
88  nVariables_.push_back(variableNamesInCategory.size());
89 
90  variables_.push_back(variablesInCategory);
91 
92  for (int j = 0; j < nVariables_[i]; ++j) {
93  int index = mvaVarMngr_.getVarIndex(variableNamesInCategory[j]);
94  if (index == -1) {
95  throw cms::Exception("MVA config failure: ")
96  << "Concerning PhotonMVAEstimator" << getTag() << std::endl
97  << "Variable " << variableNamesInCategory[j] << " not found in variable definition file!" << std::endl;
98  }
99  variables_[i].push_back(index);
100  }
101  }
102 }
103 
105  std::vector<float> const& auxVars,
106  int& iCategory) const {
107  const reco::Photon* phoPtr = dynamic_cast<const reco::Photon*>(candPtr);
108  if (phoPtr == nullptr) {
109  throw cms::Exception("MVA failure: ")
110  << " given particle is expected to be reco::Photon or pat::Photon," << std::endl
111  << " but appears to be neither" << std::endl;
112  }
113 
114  iCategory = findCategory(phoPtr);
115 
116  std::vector<float> vars;
117 
118  vars.reserve(nVariables_[iCategory]);
119  for (int i = 0; i < nVariables_[iCategory]; ++i) {
120  vars.push_back(mvaVarMngr_.getValue(variables_[iCategory][i], *phoPtr, auxVars));
121  }
122 
123  // Special case for Spring16!
124  if (getTag() == "Run2Spring16NonTrigV1" and iCategory == 1) { // Endcap category
125  // Raw value for EB only, because of loss of transparency in EE
126  // for endcap MVA only in 2016
127  double eA = effectiveAreas_->getEffectiveArea(std::abs(phoPtr->superCluster()->eta()));
128  double phoIsoCorr = vars[10] - eA * (double)vars[9] - phoIsoPtScalingCoeff_.at(1) * phoPtr->pt();
129  vars[10] = std::max(phoIsoCorr, phoIsoCutoff_);
130  }
131 
132  if (isDebug()) {
133  std::cout << " *** Inside PhotonMVAEstimator" << getTag() << std::endl;
134  std::cout << " category " << iCategory << std::endl;
135  for (int i = 0; i < nVariables_[iCategory]; ++i) {
136  std::cout << " " << mvaVarMngr_.getName(variables_[iCategory][i]) << " " << vars[i] << std::endl;
137  }
138  }
139 
140  const float response = gbrForests_.at(iCategory)->GetResponse(vars.data());
141 
142  if (isDebug()) {
143  std::cout << " ### MVA " << response << std::endl << std::endl;
144  }
145 
146  return response;
147 }
148 
150  const reco::Photon* phoPtr = dynamic_cast<const reco::Photon*>(candPtr);
151  if (phoPtr == nullptr) {
152  throw cms::Exception("MVA failure: ")
153  << " given particle is expected to be reco::Photon or pat::Photon," << std::endl
154  << " but appears to be neither" << std::endl;
155  }
156 
157  return findCategory(*phoPtr);
158 }
159 
161  for (int i = 0; i < getNCategories(); ++i) {
162  if (categoryFunctions_[i](photon))
163  return i;
164  }
165 
166  edm::LogWarning("MVA warning") << "category not defined for particle with pt " << photon.pt() << " GeV, eta "
167  << photon.superCluster()->eta() << " in PhotonMVAEstimator" << getTag();
168 
169  return -1;
170 }
171 
std::vector< std::unique_ptr< const GBRForest > > gbrForests_
const std::string & getName(int index) const
double pt() const final
transverse momentum
int findCategory(const reco::Candidate *candPtr) const override
PhotonMVAEstimator(const edm::ParameterSet &conf)
MVAVariableManager< reco::Photon > mvaVarMngr_
float getValue(int index, const ParticleType &particle, const std::vector< float > &auxVariables) const
const std::string & getTag() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
reco::SuperClusterRef superCluster() const override
Ref to SuperCluster.
std::vector< ThreadSafeFunctor< StringCutObjectSelector< reco::Photon > > > categoryFunctions_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > nVariables_
std::vector< double > phoIsoPtScalingCoeff_
int getVarIndex(const std::string &name)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
float mvaValue(const reco::Candidate *candPtr, std::vector< float > const &auxVars, int &iCategory) const override
std::unique_ptr< EffectiveAreas > effectiveAreas_
tuple cout
Definition: gather_cfg.py:144
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::vector< int > > variables_
~PhotonMVAEstimator() override
vars
Definition: DeepTauId.cc:164
Log< level::Warning, false > LogWarning
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)