CMS 3D CMS Logo

PhotonDNNEstimator.cc
Go to the documentation of this file.
4 
5 #include <iostream>
6 #include <fstream>
7 #include <memory>
8 
9 using namespace std::placeholders;
10 
11 inline uint photonModelSelector(const std::map<std::string, float>& vars, float etaThr) {
12  /*
13  Selection of the model to be applied on the photon based on eta limit
14  */
15  const auto absEta = std::abs(vars.at("eta"));
16  if (absEta <= etaThr) {
17  return 0;
18  } else {
19  return 1;
20  }
21 }
22 
24  : dnnHelper_(cfg,
26  _1,
27  (useEBModelInGap) ? PhotonDNNEstimator::ecalBarrelMaxEtaWithGap
28  : PhotonDNNEstimator::ecalBarrelMaxEtaNoGap),
29  PhotonDNNEstimator::dnnAvaibleInputs),
30  useEBModelInGap_(useEBModelInGap) {}
31 
32 std::vector<tensorflow::Session*> PhotonDNNEstimator::getSessions() const { return dnnHelper_.getSessions(); };
33 
34 const std::vector<std::string> PhotonDNNEstimator::dnnAvaibleInputs = {{"pt",
35  "eta",
36  "hadTowOverEm",
37  "trkSumPtHollowConeDR03",
38  "EcalRecHit",
39  "SigmaIetaIeta",
40  "SigmaIetaIetaFull5x5",
41  "SigmaIEtaIPhiFull5x5",
42  "EcalPFClusterIso",
43  "HcalPFClusterIso",
44  "HasPixelSeed",
45  "R9Full5x5",
46  "hcalTower"}};
47 
48 std::map<std::string, float> PhotonDNNEstimator::getInputsVars(const reco::Photon& photon) const {
49  // Prepare a map with all the defined variables
50  std::map<std::string, float> variables;
51  variables["pt"] = photon.pt();
52  variables["eta"] = photon.eta();
53  variables["hadTowOverEm"] = photon.hadTowOverEmValid() ? photon.hadTowOverEm() : 0;
54  variables["trkSumPtHollowConeDR03"] = photon.trkSumPtHollowConeDR03();
55  variables["EcalRecHit"] = photon.ecalRecHitSumEtConeDR03();
56  variables["SigmaIetaIeta"] = photon.sigmaIetaIeta();
57  variables["SigmaIetaIetaFull5x5"] = photon.full5x5_sigmaIetaIeta();
58  variables["SigmaIEtaIPhiFull5x5"] = photon.full5x5_showerShapeVariables().sigmaIetaIphi;
59  variables["EcalPFClusterIso"] = photon.ecalPFClusterIso();
60  variables["HcalPFClusterIso"] = photon.hcalPFClusterIso();
61  variables["HasPixelSeed"] = (Int_t)photon.hasPixelSeed();
62  variables["R9Full5x5"] = photon.full5x5_r9();
63  variables["hcalTower"] = photon.hcalTowerSumEtConeDR03();
64  variables["R9Full5x5"] = photon.full5x5_r9();
65  // Define more variables here and use them directly in the model config!
66  return variables;
67 }
68 
69 std::vector<std::pair<uint, std::vector<float>>> PhotonDNNEstimator::evaluate(
70  const reco::PhotonCollection& photons, const std::vector<tensorflow::Session*>& sessions) const {
71  // Collect the map of variables for each candidate and call the dnnHelper
72  // Scaling, model selection and running is performed in the helper
73  std::vector<std::map<std::string, float>> inputs;
74  for (const auto& photon : photons) {
75  inputs.push_back(getInputsVars(photon));
76  }
77  return dnnHelper_.evaluate(inputs, sessions);
78 }
uint photonModelSelector(const std::map< std::string, float > &vars, float etaThr)
std::vector< tensorflow::Session * > getSessions() const
std::map< std::string, float > getInputsVars(const reco::Photon &ele) const
std::vector< std::pair< uint, std::vector< float > > > evaluate(const reco::PhotonCollection &ele, const std::vector< tensorflow::Session *> &sessions) const
std::vector< tensorflow::Session * > getSessions() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const egammaTools::EgammaDNNHelper dnnHelper_
std::vector< Photon > PhotonCollection
collectin of Photon objects
Definition: PhotonFwd.h:9
static const std::vector< std::string > dnnAvaibleInputs
std::vector< std::pair< uint, std::vector< float > > > evaluate(const std::vector< std::map< std::string, float >> &candidates, const std::vector< tensorflow::Session *> &sessions) const
vars
Definition: DeepTauId.cc:166
PhotonDNNEstimator(const egammaTools::DNNConfiguration &, const bool useEBModelInGap)