CMS 3D CMS Logo

PFMETAlgorithmMVA.cc
Go to the documentation of this file.
2 
5 
7 
8 #include <TFile.h>
9 
10 #include <iomanip>
11 
12 enum MVAType { kBaseline = 0 };
13 
14 const double Pi = std::cos(-1);
15 
17  if (input == "sumet")
18  return "particleFlow_SumET";
19  if (input == "npv")
20  return "nPV";
21  if (input == "pfu")
22  return "particleFlow_U";
23  if (input == "pfuphi")
24  return "particleFlow_UPhi";
25  if (input == "tksumet")
26  return "track_SumET";
27  if (input == "tku")
28  return "track_U";
29  if (input == "tkuphi")
30  return "track_UPhi";
31  if (input == "nopusumet")
32  return "noPileUp_SumET";
33  if (input == "nopuu")
34  return "noPileUp_U";
35  if (input == "nopuuphi")
36  return "noPileUp_UPhi";
37  if (input == "pusumet")
38  return "pileUp_SumET";
39  if (input == "pumet")
40  return "pileUp_MET";
41  if (input == "pumetphi")
42  return "pileUp_METPhi";
43  if (input == "pucsumet")
44  return "pileUpCorrected_SumET";
45  if (input == "pucu")
46  return "pileUpCorrected_U";
47  if (input == "pucuphi")
48  return "pileUpCorrected_UPhi";
49  if (input == "jetpt1")
50  return "jet1_pT";
51  if (input == "jeteta1")
52  return "jet1_eta";
53  if (input == "jetphi1")
54  return "jet1_Phi";
55  if (input == "jetpt2")
56  return "jet2_pT";
57  if (input == "jeteta2")
58  return "jet2_eta";
59  if (input == "jetphi2")
60  return "jet2_Phi";
61  if (input == "nalljet")
62  return "nJets";
63  if (input == "njet")
64  return "numJetsPtGt30";
65  if (input == "uphi_mva")
66  return "PhiCor_UPhi";
67  if (input == "uphix_mva")
68  return "PhiCor_UPhi";
69  if (input == "ux_mva")
70  return "RecoilCor_U";
71  return input;
72 }
73 
75  if (inputFileName.location() == edm::FileInPath::Unknown)
76  throw cms::Exception("PFMETAlgorithmMVA::loadMVA") << " Failed to find File = " << inputFileName << " !!\n";
77  std::unique_ptr<TFile> inputFile(new TFile(inputFileName.fullPath().data()));
78 
79  std::vector<std::string>* lVec = (std::vector<std::string>*)inputFile->Get("varlist");
80 
81  if (lVec == nullptr) {
82  throw cms::Exception("PFMETAlgorithmMVA::loadMVA")
83  << " Failed to load mva file : " << inputFileName.fullPath().data() << " is not a proper file !!\n";
84  }
85 
86  std::vector<std::string> variableNames;
87  for (unsigned int i = 0; i < lVec->size(); ++i) {
88  variableNames.push_back(updateVariableNames(lVec->at(i)));
89  }
90 
91  if (mvaName.find(mvaNameU_) != std::string::npos)
92  varForU_ = variableNames;
93  else if (mvaName.find(mvaNameDPhi_) != std::string::npos)
94  varForDPhi_ = variableNames;
95  else if (mvaName.find(mvaNameCovU1_) != std::string::npos)
96  varForCovU1_ = variableNames;
97  else if (mvaName.find(mvaNameCovU2_) != std::string::npos)
98  varForCovU2_ = variableNames;
99  else
100  throw cms::Exception("PFMETAlgorithmMVA::loadMVA")
101  << "MVA MET weight file tree names do not match specified inputs" << std::endl;
102 
103  const GBRForest* mva = (GBRForest*)inputFile->Get(mvaName.data());
104  if (!mva)
105  throw cms::Exception("PFMETAlgorithmMVA::loadMVA")
106  << " Failed to load MVA = " << mvaName.data() << " from file = " << inputFileName.fullPath().data() << " !!\n";
107 
108  return mva;
109 }
110 
112  : utils_(cfg),
113  mvaReaderU_(nullptr),
114  mvaReaderDPhi_(nullptr),
115  mvaReaderCovU1_(nullptr),
116  mvaReaderCovU2_(nullptr),
117  cfg_(cfg) {
119 
120  edm::ParameterSet cfgInputRecords = cfg_.getParameter<edm::ParameterSet>("inputRecords");
121  mvaNameU_ = cfgInputRecords.getParameter<std::string>("U");
122  mvaNameDPhi_ = cfgInputRecords.getParameter<std::string>("DPhi");
123  mvaNameCovU1_ = cfgInputRecords.getParameter<std::string>("CovU1");
124  mvaNameCovU2_ = cfgInputRecords.getParameter<std::string>("CovU2");
125 
126  loadMVAfromDB_ = cfg.getParameter<bool>("loadMVAfromDB");
127  if (loadMVAfromDB_) {
132  }
133 }
134 
136  if (!loadMVAfromDB_) {
137  delete mvaReaderU_;
138  delete mvaReaderDPhi_;
139  delete mvaReaderCovU1_;
140  delete mvaReaderCovU2_;
141  }
142 }
143 
144 //-------------------------------------------------------------------------------
146  if (loadMVAfromDB_) {
151  } else {
152  edm::ParameterSet cfgInputFileNames = cfg_.getParameter<edm::ParameterSet>("inputFileNames");
153 
154  edm::FileInPath inputFileNameU = cfgInputFileNames.getParameter<edm::FileInPath>("U");
155  mvaReaderU_ = loadMVAfromFile(inputFileNameU, mvaNameU_);
156  edm::FileInPath inputFileNameDPhi = cfgInputFileNames.getParameter<edm::FileInPath>("DPhi");
157  mvaReaderDPhi_ = loadMVAfromFile(inputFileNameDPhi, mvaNameDPhi_);
158  edm::FileInPath inputFileNameCovU1 = cfgInputFileNames.getParameter<edm::FileInPath>("CovU1");
159  mvaReaderCovU1_ = loadMVAfromFile(inputFileNameCovU1, mvaNameCovU1_);
160  edm::FileInPath inputFileNameCovU2 = cfgInputFileNames.getParameter<edm::FileInPath>("CovU2");
161  mvaReaderCovU2_ = loadMVAfromFile(inputFileNameCovU2, mvaNameCovU2_);
162  }
163 }
164 
165 //-------------------------------------------------------------------------------
166 void PFMETAlgorithmMVA::setInput(const std::vector<reco::PUSubMETCandInfo>& leptons,
167  const std::vector<reco::PUSubMETCandInfo>& jets,
168  const std::vector<reco::PUSubMETCandInfo>& pfCandidates,
169  const std::vector<reco::Vertex::Point>& vertices) {
171 
174 
177 
178  const std::vector<reco::PUSubMETCandInfo> jets_cleaned = utils_.getCleanedJets();
179 
185 
186  reco::Candidate::LorentzVector jet1P4 = utils_.leadJetP4(jets_cleaned);
187  reco::Candidate::LorentzVector jet2P4 = utils_.subleadJetP4(jets_cleaned);
188 
189  var_["particleFlow_U"] = pfRecoil_data.met;
190  var_["particleFlow_SumET"] = pfRecoil_data.sumet;
191  var_["particleFlow_UPhi"] = pfRecoil_data.phi;
192 
193  var_["track_SumET"] = chHSRecoil_data.sumet / var_["particleFlow_SumET"];
194  var_["track_U"] = chHSRecoil_data.met;
195  var_["track_UPhi"] = chHSRecoil_data.phi;
196 
197  var_["noPileUp_SumET"] = hsRecoil_data.sumet / var_["particleFlow_SumET"];
198  var_["noPileUp_U"] = hsRecoil_data.met;
199  var_["noPileUp_UPhi"] = hsRecoil_data.phi;
200 
201  var_["pileUp_SumET"] = puRecoil_data.sumet / var_["particleFlow_SumET"];
202  var_["pileUp_MET"] = puRecoil_data.met;
203  var_["pileUp_METPhi"] = puRecoil_data.phi;
204 
205  var_["pileUpCorrected_SumET"] = hsMinusNeutralPUMEt_data.sumet / var_["particleFlow_SumET"];
206  var_["pileUpCorrected_U"] = hsMinusNeutralPUMEt_data.met;
207  var_["pileUpCorrected_UPhi"] = hsMinusNeutralPUMEt_data.phi;
208 
209  var_["jet1_pT"] = jet1P4.pt();
210  var_["jet1_eta"] = jet1P4.eta();
211  var_["jet1_Phi"] = jet1P4.phi();
212  var_["jet2_pT"] = jet2P4.pt();
213  var_["jet2_eta"] = jet2P4.eta();
214  var_["jet2_Phi"] = jet2P4.phi();
215 
216  var_["numJetsPtGt30"] = utils_.numJetsAboveThreshold(jets_cleaned, 30.);
217  var_["nJets"] = jets_cleaned.size();
218  var_["nPV"] = vertices.size();
219 }
220 
221 //-------------------------------------------------------------------------------
222 std::unique_ptr<float[]> PFMETAlgorithmMVA::createFloatVector(std::vector<std::string> variableNames) {
223  std::unique_ptr<float[]> floatVector(new float[variableNames.size()]);
224  int i = 0;
225  for (const auto& variableName : variableNames) {
226  floatVector[i++] = var_[variableName];
227  }
228  return floatVector;
229 }
230 
231 //-------------------------------------------------------------------------------
233  // CV: MVAs needs to be evaluated in order { DPhi, U1, CovU1, CovU2 }
234  // as MVA for U1 (CovU1, CovU2) uses output of DPhi (DPhi and U1) MVA
236  var_["PhiCor_UPhi"] = var_["particleFlow_UPhi"] + mvaOutputDPhi_;
238  var_["RecoilCor_U"] = var_["particleFlow_U"] * mvaOutputU_;
239  var_["RecoilCor_UPhi"] = var_["PhiCor_UPhi"];
242 
243  // compute MET(Photon check)
244  if (hasPhotons_) {
245  //Fix events with unphysical properties
246  double sumLeptonPt = std::max(sqrt(sumLeptonPx_ * sumLeptonPx_ + sumLeptonPy_ * sumLeptonPy_), 1.);
247  if (var_["track_U"] / sumLeptonPt < 0.1 || var_["noPileUp_U"] / sumLeptonPt < 0.1) {
248  mvaOutputU_ = 1.;
249  mvaOutputDPhi_ = 0.;
250  }
251  }
252  computeMET();
253 }
254 //-------------------------------------------------------------------------------
255 
257  double U = var_["RecoilCor_U"];
258  double Phi = var_["PhiCor_UPhi"];
259  if (U < 0.)
260  Phi += Pi; //RF: No sign flip for U necessary in that case?
261  double cosPhi = std::cos(Phi);
262  double sinPhi = std::sin(Phi);
263  double metPx = U * cosPhi - sumLeptonPx_;
264  double metPy = U * sinPhi - sumLeptonPy_;
265  double metPt = sqrt(metPx * metPx + metPy * metPy);
266  mvaMEt_.SetCoordinates(metPx, metPy, 0., metPt);
267  // compute MET uncertainties in dirrections parallel and perpendicular to hadronic recoil
268  // (neglecting uncertainties on lepton momenta)
271  mvaMEtCov_(1, 0) = mvaMEtCov_(0, 1);
273 }
274 
275 //-------------------------------------------------------------------------------
276 const float PFMETAlgorithmMVA::GetResponse(const GBRForest* Reader, std::vector<std::string>& variableNames) {
277  std::unique_ptr<float[]> mvaInputVector = createFloatVector(variableNames);
278  double result = Reader->GetResponse(mvaInputVector.get());
279  return result;
280 }
281 
282 //-------------------------------------------------------------------------------
283 void PFMETAlgorithmMVA::print(std::ostream& stream) const {
284  stream << "<PFMETAlgorithmMVA::print>:" << std::endl;
285  for (const auto& entry : var_)
286  stream << entry.first << " = " << entry.second << std::endl;
287  stream << " covU1 = " << mvaOutputCovU1_ << ", covU2 = " << mvaOutputCovU2_ << std::endl;
288  stream << " sum(leptons): Pt = " << sqrt(sumLeptonPx_ * sumLeptonPx_ + sumLeptonPy_ * sumLeptonPy_) << ","
289  << " phi = " << atan2(sumLeptonPy_, sumLeptonPx_) << " "
290  << "(Px = " << sumLeptonPx_ << ", Py = " << sumLeptonPy_ << ")";
291  stream << " MVA output: U = " << mvaOutputU_ << ", dPhi = " << mvaOutputDPhi_ << ","
292  << " covU1 = " << mvaOutputCovU1_ << ", covU2 = " << mvaOutputCovU2_ << std::endl;
293  stream << std::endl;
294 }
std::vector< std::string > varForCovU2_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
unsigned numJetsAboveThreshold(const std::vector< reco::PUSubMETCandInfo > &, double)
MvaMEtUtilities utils_
const GBRForest * loadMVAfromFile(const edm::FileInPath &inputFileName, const std::string &mvaName)
std::vector< std::string > varForU_
double getLeptonsChSumMEX() const
const GBRForest * mvaReaderCovU1_
const float GetResponse(const GBRForest *Reader, std::vector< std::string > &variableNames)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const GBRForest * mvaReaderCovU2_
const std::string updateVariableNames(std::string input)
void print(std::ostream &) const
edm::ParameterSet cfg_
reco::Candidate::LorentzVector subleadJetP4(const std::vector< reco::PUSubMETCandInfo > &)
uint32_t T const *__restrict__ uint32_t const *__restrict__ int32_t int Histo::index_type cudaStream_t stream
std::string mvaNameCovU1_
void setInput(const std::vector< reco::PUSubMETCandInfo > &, const std::vector< reco::PUSubMETCandInfo > &, const std::vector< reco::PUSubMETCandInfo > &, const std::vector< reco::Vertex::Point > &)
static std::string const input
Definition: EdmProvDump.cc:50
const GBRForest * mvaReaderDPhi_
const GBRForest * mvaReaderU_
T sqrt(T t)
Definition: SSEVec.h:19
std::unique_ptr< float[]> createFloatVector(std::vector< std::string > variableNames)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double getLeptonsSumMEY() const
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
edm::ESGetToken< GBRForest, GBRWrapperRcd > mvaTokenCovU1_
void computeAllSums(const std::vector< reco::PUSubMETCandInfo > &jets, const std::vector< reco::PUSubMETCandInfo > &leptons, const std::vector< reco::PUSubMETCandInfo > &pfCandidates)
void initialize(const edm::EventSetup &)
reco::METCovMatrix mvaMEtCov_
std::vector< std::string > varForCovU1_
edm::ESGetToken< GBRForest, GBRWrapperRcd > mvaTokenCovU2_
const std::vector< reco::PUSubMETCandInfo > & getCleanedJets() const
reco::Candidate::LorentzVector leadJetP4(const std::vector< reco::PUSubMETCandInfo > &)
double getLeptonsSumMEX() const
const double Pi
CommonMETData computeRecoil(int metType)
PFMETAlgorithmMVA(const edm::ParameterSet &cfg, edm::ConsumesCollector iC)
std::vector< std::string > varForDPhi_
std::map< std::string, float > var_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
std::string mvaNameDPhi_
double getLeptonsChSumMEY() const
double GetResponse(const float *vector) const
Definition: GBRForest.h:48
reco::Candidate::LorentzVector mvaMEt_
edm::ESGetToken< GBRForest, GBRWrapperRcd > mvaTokenDPhi_
edm::ESGetToken< GBRForest, GBRWrapperRcd > mvaTokenU_
std::string mvaNameCovU2_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29