CMS 3D CMS Logo

PhotonMVAEstimator.cc
Go to the documentation of this file.
2 
5  name_(conf.getParameter<std::string>("mvaName")),
6  tag_(conf.getParameter<std::string>("mvaTag")),
7  methodName_("BDTG method"),
8  mvaVarMngr_(conf.getParameter<std::string>("variableDefinition")),
9  ebeeSplit_ (conf.getParameter<double> ("ebeeSplit")),
10  debug_(conf.getUntrackedParameter<bool>("debug", false))
11 {
12  //
13  // Construct the MVA estimators
14  //
15  if (tag_ == "Run2Spring16NonTrigV1") {
16  effectiveAreas_ = std::make_unique<EffectiveAreas>((conf.getParameter<edm::FileInPath>("effAreasConfigFile")).fullPath());
17  phoIsoPtScalingCoeff_ = conf.getParameter<std::vector<double >>("phoIsoPtScalingCoeff");
18  phoIsoCutoff_ = conf.getParameter<double>("phoIsoCutoff");
19  }
20 
21  const std::vector <std::string> weightFileNames
22  = conf.getParameter<std::vector<std::string> >("weightFileNames");
23 
24  // Initialize GBRForests
25  if( (int)(weightFileNames.size()) != nCategories_ )
26  throw cms::Exception("MVA config failure: ")
27  << "wrong number of weightfiles in " << name_ << tag_ << std::endl;
28 
29  gbrForests_.clear();
30  // Create a TMVA reader object for each category
31  for(int i=0; i<nCategories_; i++){
32 
33  std::vector<std::string> variableNamesInCategory;
34  std::vector<int> variablesInCategory;
35 
36  // Use unique_ptr so that all readers are properly cleaned up
37  // when the vector clear() is called in the destructor
38 
39  gbrForests_.push_back( GBRForestTools::createGBRForest( weightFileNames[i], variableNamesInCategory ) );
40 
41  nVariables_.push_back(variableNamesInCategory.size());
42 
43  variables_.push_back(variablesInCategory);
44 
45  for (int j=0; j<nVariables_[i];++j) {
46  int index = mvaVarMngr_.getVarIndex(variableNamesInCategory[j]);
47  if(index == -1) {
48  throw cms::Exception("MVA config failure: ")
49  << "Concerning " << name_ << tag_ << std::endl
50  << "Variable " << variableNamesInCategory[j]
51  << " not found in variable definition file!" << std::endl;
52  }
53  variables_[i].push_back(index);
54 
55  }
56  }
57 }
58 
61 }
62 
65 
66  const int iCategory = findCategory( candPtr );
67 
68  const edm::Ptr<reco::Photon> phoPtr{ candPtr };
69  if( phoPtr.get() == nullptr) {
70  throw cms::Exception("MVA failure: ")
71  << " given particle is expected to be reco::Photon or pat::Photon," << std::endl
72  << " but appears to be neither" << std::endl;
73  }
74 
75  std::vector<float> vars;
76 
77  for (int i = 0; i < nVariables_[iCategory]; ++i) {
78  vars.push_back(mvaVarMngr_.getValue(variables_[iCategory][i], phoPtr, iEvent));
79  }
80 
81  // Special case for Spring16!
82  if (tag_ == "Run2Spring16NonTrigV1" and iCategory == CAT_EE) {
83  // Raw value for EB only, because of loss of transparency in EE
84  // for endcap MVA only in 2016
85  double eA = effectiveAreas_->getEffectiveArea( std::abs(phoPtr->superCluster()->eta()) );
86  double phoIsoCorr = vars[10] - eA*(double)vars[9] - phoIsoPtScalingCoeff_.at(1) * phoPtr->pt();
87  vars[10] = TMath::Max( phoIsoCorr, phoIsoCutoff_);
88  }
89 
90  if(debug_) {
91  std::cout << " *** Inside " << name_ << tag_ << std::endl;
92  std::cout << " category " << iCategory << std::endl;
93  for (int i = 0; i < nVariables_[iCategory]; ++i) {
94  std::cout << " " << mvaVarMngr_.getName(variables_[iCategory][i]) << " " << vars[i] << std::endl;
95  }
96  }
97 
98  const float response = gbrForests_.at(iCategory)->GetResponse(vars.data());
99 
100  if(debug_) {
101  std::cout << " ### MVA " << response << std::endl << std::endl;
102  }
103 
104  return response;
105 }
106 
108 
109  // Try to cast the particle into a reco particle.
110  // This should work for both reco and pat.
111  auto pho = dynamic_cast<reco::Photon const*>(candPtr.get());
112 
113  if( pho == nullptr) {
114  throw cms::Exception("MVA failure: ")
115  << " given particle is expected to be reco::Photon or pat::Photon," << std::endl
116  << " but appears to be neither" << std::endl;
117  }
118 
119  float eta = pho->superCluster()->eta();
120 
121  //
122  // Determine the category
123  //
124  if ( std::abs(eta) < ebeeSplit_)
125  return CAT_EB;
126  else
127  return CAT_EE;
128 }
129 
131  // All tokens for event content needed by this MVA
132  // Tags from the variable helper
133  for (auto &tag : mvaVarMngr_.getHelperInputTags()) {
134  cc.consumes<edm::ValueMap<float>>(tag);
135  }
136  for (auto &tag : mvaVarMngr_.getGlobalInputTags()) {
137  cc.consumes<double>(tag);
138  }
139 }
140 
143  "PhotonMVAEstimator");
const std::string getName(int index) const
T getParameter(std::string const &) const
void setConsumes(edm::ConsumesCollector &&) const override
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:159
PhotonMVAEstimator(const edm::ParameterSet &conf)
MVAVariableManager< reco::Photon > mvaVarMngr_
std::vector< edm::InputTag > getHelperInputTags() const
const std::string tag_
float mvaValue(const edm::Ptr< reco::Candidate > &candPtr, const edm::EventBase &) const override
int iEvent
Definition: GenABIO.cc:230
int findCategory(const edm::Ptr< reco::Candidate > &candPtr) const override
static std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightFile)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int getVarIndex(std::string &name)
std::vector< std::unique_ptr< const GBRForest > > gbrForests_
T Max(T a, T b)
Definition: MathUtil.h:44
std::vector< int > nVariables_
std::vector< double > phoIsoPtScalingCoeff_
float getValue(int index, const edm::Ptr< ParticleType > &ptclPtr, const edm::EventBase &iEvent) const
const std::string name_
std::unique_ptr< EffectiveAreas > effectiveAreas_
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::vector< int > > variables_
~PhotonMVAEstimator() override
std::vector< edm::InputTag > getGlobalInputTags() const