CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ElectronDNNEstimator.cc
Go to the documentation of this file.
5 
6 #include <iostream>
7 #include <fstream>
8 #include <memory>
9 
10 using namespace std::placeholders;
11 
12 inline uint electronModelSelector(const std::map<std::string, float>& vars, float ptThr, float etaThr) {
13  /*
14  Selection of the model to be applied on the electron based on pt/eta cuts or whatever selection
15  */
16  const auto pt = vars.at("pt");
17  const auto absEta = std::abs(vars.at("eta"));
18  if (pt < ptThr)
19  return 0;
20  else {
21  if (absEta <= etaThr) {
22  return 1;
23  } else {
24  return 2;
25  }
26  }
27 }
28 
30  : dnnHelper_(cfg,
31  std::bind(electronModelSelector,
32  _1,
33  ElectronDNNEstimator::ptThreshold,
34  (useEBModelInGap) ? ElectronDNNEstimator::ecalBarrelMaxEtaWithGap
35  : ElectronDNNEstimator::ecalBarrelMaxEtaNoGap),
36  ElectronDNNEstimator::dnnAvaibleInputs),
37  useEBModelInGap_(useEBModelInGap) {}
38 
39 std::vector<tensorflow::Session*> ElectronDNNEstimator::getSessions() const { return dnnHelper_.getSessions(); };
40 
41 const std::vector<std::string> ElectronDNNEstimator::dnnAvaibleInputs = {
42  {"pt",
43  "eta",
44  "fbrem",
45  "abs(deltaEtaSuperClusterTrackAtVtx)",
46  "abs(deltaPhiSuperClusterTrackAtVtx)",
47  "full5x5_sigmaIetaIeta",
48  "full5x5_hcalOverEcal",
49  "eSuperClusterOverP",
50  "full5x5_e1x5",
51  "eEleClusterOverPout",
52  "closestCtfTrackNormChi2",
53  "closestCtfTrackNLayers",
54  "gsfTrack.missing_inner_hits",
55  "dr03TkSumPt",
56  "dr03EcalRecHitSumEt",
57  "dr03HcalTowerSumEt",
58  "gsfTrack.normalizedChi2",
59  "superCluster.eta",
60  "ecalPFClusterIso",
61  "hcalPFClusterIso",
62  "numberOfBrems",
63  "abs(deltaEtaSeedClusterTrackAtCalo)",
64  "hadronicOverEm",
65  "full5x5_e2x5Max",
66  "full5x5_e5x5",
67  "full5x5_sigmaIphiIphi",
68  "1_minus_full5x5_e1x5_ratio_full5x5_e5x5",
69  "full5x5_e1x5_ratio_full5x5_e5x5",
70  "full5x5_r9",
71  "gsfTrack.trackerLayersWithMeasurement",
72  "superCluster.energy",
73  "superCluster.rawEnergy",
74  "superClusterFbrem",
75  "1_ratio_ecalEnergy_minus_1_ratio_trackMomentumAtVtx.R",
76  "superCluster.preshowerEnergy_ratio_superCluster.rawEnergy",
77  "convVtxFitProb",
78  "superCluster.clustersSize",
79  "ecalEnergyError_ratio_ecalEnergy",
80  "superClusterFbrem_plus_superCluster.energy",
81  "superClusterFbrem_plus_superCluster.rawEnergy",
82  "trackMomentumError",
83  "trackMomentumError_ratio_pt",
84  "full5x5_e5x5_ratio_superCluster.rawEnergy",
85  "full5x5_e5x5_ratio_superCluster.energy"}};
86 
87 std::map<std::string, float> ElectronDNNEstimator::getInputsVars(const reco::GsfElectron& ele) const {
88  // Prepare a map with all the defined variables
89  std::map<std::string, float> variables;
90  variables["pt"] = ele.pt();
91  variables["eta"] = ele.eta();
92  variables["fbrem"] = ele.fbrem();
93  variables["abs(deltaEtaSuperClusterTrackAtVtx)"] = std::abs(ele.deltaEtaSuperClusterTrackAtVtx());
94  variables["abs(deltaPhiSuperClusterTrackAtVtx)"] = std::abs(ele.deltaPhiSuperClusterTrackAtVtx());
95  variables["full5x5_sigmaIetaIeta"] = ele.full5x5_sigmaIetaIeta();
96  variables["full5x5_hcalOverEcal"] = ele.full5x5_hcalOverEcalValid() ? ele.full5x5_hcalOverEcal() : 0;
97  variables["eSuperClusterOverP"] = ele.eSuperClusterOverP();
98  variables["full5x5_e1x5"] = ele.full5x5_e1x5();
99  variables["eEleClusterOverPout"] = ele.eEleClusterOverPout();
100  variables["closestCtfTrackNormChi2"] = ele.closestCtfTrackNormChi2();
101  variables["closestCtfTrackNLayers"] = ele.closestCtfTrackNLayers();
102  variables["gsfTrack.missing_inner_hits"] =
103  ele.gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
104  variables["dr03TkSumPt"] = ele.dr03TkSumPt();
105  variables["dr03EcalRecHitSumEt"] = ele.dr03EcalRecHitSumEt();
106  variables["dr03HcalTowerSumEt"] = ele.dr03HcalTowerSumEt();
107  variables["gsfTrack.normalizedChi2"] = ele.gsfTrack()->normalizedChi2();
108  variables["superCluster.eta"] = ele.superCluster()->eta();
109  variables["ecalPFClusterIso"] = ele.ecalPFClusterIso();
110  variables["hcalPFClusterIso"] = ele.hcalPFClusterIso();
111  variables["numberOfBrems"] = ele.numberOfBrems();
112  variables["abs(deltaEtaSeedClusterTrackAtCalo)"] = std::abs(ele.deltaEtaSeedClusterTrackAtCalo());
113  variables["hadronicOverEm"] = ele.hcalOverEcalValid() ? ele.hadronicOverEm() : 0;
114  variables["full5x5_e2x5Max"] = ele.full5x5_e2x5Max();
115  variables["full5x5_e5x5"] = ele.full5x5_e5x5();
116 
117  variables["full5x5_sigmaIphiIphi"] = ele.full5x5_sigmaIphiIphi();
118  variables["1_minus_full5x5_e1x5_ratio_full5x5_e5x5"] = 1 - ele.full5x5_e1x5() / ele.full5x5_e5x5();
119  variables["full5x5_e1x5_ratio_full5x5_e5x5"] = ele.full5x5_e1x5() / ele.full5x5_e5x5();
120  variables["full5x5_r9"] = ele.full5x5_r9();
121  variables["gsfTrack.trackerLayersWithMeasurement"] = ele.gsfTrack()->hitPattern().trackerLayersWithMeasurement();
122  variables["superCluster.energy"] = ele.superCluster()->energy();
123  variables["superCluster.rawEnergy"] = ele.superCluster()->rawEnergy();
124  variables["superClusterFbrem"] = ele.superClusterFbrem();
125  variables["1_ratio_ecalEnergy_minus_1_ratio_trackMomentumAtVtx.R"] =
126  1 / ele.ecalEnergy() - 1 / ele.trackMomentumAtVtx().R();
127  variables["superCluster.preshowerEnergy_ratio_superCluster.rawEnergy"] =
128  ele.superCluster()->preshowerEnergy() / ele.superCluster()->rawEnergy();
129  variables["convVtxFitProb"] = ele.convVtxFitProb();
130  variables["superCluster.clustersSize"] = ele.superCluster()->clustersSize();
131  variables["ecalEnergyError_ratio_ecalEnergy"] = ele.ecalEnergyError() / ele.ecalEnergy();
132  variables["superClusterFbrem_plus_superCluster.energy"] = ele.superClusterFbrem() + ele.superCluster()->energy();
133  variables["superClusterFbrem_plus_superCluster.rawEnergy"] =
134  ele.superClusterFbrem() + ele.superCluster()->rawEnergy();
135  variables["trackMomentumError"] = ele.trackMomentumError();
136  variables["trackMomentumError_ratio_pt"] = ele.trackMomentumError() / ele.pt();
137  variables["full5x5_e5x5_ratio_superCluster.rawEnergy"] = ele.full5x5_e5x5() / ele.superCluster()->rawEnergy();
138  variables["full5x5_e5x5_ratio_superCluster.energy"] = ele.full5x5_e5x5() / ele.superCluster()->energy();
139 
140  // Define more variables here and use them directly in the model config!
141  return variables;
142 }
143 
144 std::vector<std::vector<float>> ElectronDNNEstimator::evaluate(
145  const reco::GsfElectronCollection& electrons, const std::vector<tensorflow::Session*>& sessions) const {
146  // Collect the map of variables for each candidate and call the dnnHelper
147  // Scaling, model selection and running is performed in the helper
148  std::vector<std::map<std::string, float>> inputs;
149  for (const auto& ele : electrons) {
150  inputs.push_back(getInputsVars(ele));
151  }
152  return dnnHelper_.evaluate(inputs, sessions);
153 }
std::vector< std::vector< float > > evaluate(const std::vector< std::map< std::string, float >> &candidates, const std::vector< tensorflow::Session * > &sessions) const
uint electronModelSelector(const std::map< std::string, float > &vars, float ptThr, float etaThr)
ElectronDNNEstimator(const egammaTools::DNNConfiguration &, const bool useEBModelInGap)
float trackMomentumError() const
Definition: GsfElectron.h:884
tuple cfg
Definition: looper.py:296
double pt() const final
transverse momentum
float eSuperClusterOverP() const
Definition: GsfElectron.h:221
float full5x5_e5x5() const
Definition: GsfElectron.h:475
float full5x5_e1x5() const
Definition: GsfElectron.h:473
math::XYZVectorF trackMomentumAtVtx() const
Definition: GsfElectron.h:268
std::map< std::string, float > getInputsVars(const reco::GsfElectron &ele) const
float dr03HcalTowerSumEt(int depth=0) const
Definition: GsfElectron.h:576
float full5x5_sigmaIphiIphi() const
Definition: GsfElectron.h:472
float fbrem() const
Definition: GsfElectron.h:809
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
float superClusterFbrem() const
Definition: GsfElectron.h:803
int closestCtfTrackNLayers() const
Definition: GsfElectron.h:168
float full5x5_sigmaIetaIeta() const
Definition: GsfElectron.h:471
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:156
static const std::vector< std::string > dnnAvaibleInputs
float deltaEtaSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:225
float hadronicOverEm() const
Definition: GsfElectron.h:500
float closestCtfTrackNormChi2() const
Definition: GsfElectron.h:165
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:228
int numberOfBrems() const
Definition: GsfElectron.h:808
std::vector< tensorflow::Session * > getSessions() const
float dr03TkSumPt() const
Definition: GsfElectron.h:557
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float eEleClusterOverPout() const
Definition: GsfElectron.h:224
float hcalPFClusterIso() const
Definition: GsfElectron.h:732
float full5x5_hcalOverEcal(int depth=0) const
Definition: GsfElectron.h:477
float ecalEnergyError() const
Definition: GsfElectron.h:897
std::vector< tensorflow::Session * > getSessions() const
float full5x5_e2x5Max() const
Definition: GsfElectron.h:474
float ecalPFClusterIso() const
Definition: GsfElectron.h:731
float ecalEnergy() const
Definition: GsfElectron.h:896
float dr03EcalRecHitSumEt() const
Definition: GsfElectron.h:559
bool full5x5_hcalOverEcalValid() const
Definition: GsfElectron.h:479
float full5x5_r9() const
Definition: GsfElectron.h:476
float deltaEtaSeedClusterTrackAtCalo() const
Definition: GsfElectron.h:226
float convVtxFitProb() const
Definition: GsfElectron.h:648
bool hcalOverEcalValid() const
Definition: GsfElectron.h:462
const egammaTools::EgammaDNNHelper dnnHelper_
vars
Definition: DeepTauId.cc:164
std::vector< std::vector< float > > evaluate(const reco::GsfElectronCollection &ele, const std::vector< tensorflow::Session * > &sessions) const
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
double eta() const final
momentum pseudorapidity