CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Static Public Attributes | Private Attributes
ElectronDNNEstimator Class Reference

#include <ElectronDNNEstimator.h>

Public Member Functions

 ElectronDNNEstimator (const egammaTools::DNNConfiguration &, const bool useEBModelInGap)
 
std::vector< std::vector< float > > evaluate (const reco::GsfElectronCollection &ele, const std::vector< tensorflow::Session * > &sessions) const
 
std::map< std::string, float > getInputsVars (const reco::GsfElectron &ele) const
 
std::vector
< tensorflow::Session * > 
getSessions () const
 

Static Public Attributes

static const std::vector
< std::string > 
dnnAvaibleInputs
 
static constexpr float ecalBarrelMaxEtaNoGap = 1.485
 
static constexpr float ecalBarrelMaxEtaWithGap = 1.566
 
static constexpr float ptThreshold = 10.
 

Private Attributes

const egammaTools::EgammaDNNHelper dnnHelper_
 
const bool useEBModelInGap_
 

Detailed Description

Definition at line 12 of file ElectronDNNEstimator.h.

Constructor & Destructor Documentation

ElectronDNNEstimator::ElectronDNNEstimator ( const egammaTools::DNNConfiguration cfg,
const bool  useEBModelInGap 
)

Definition at line 29 of file ElectronDNNEstimator.cc.

30  : dnnHelper_(cfg,
31  std::bind(electronModelSelector,
32  _1,
37  useEBModelInGap_(useEBModelInGap) {}
uint electronModelSelector(const std::map< std::string, float > &vars, float ptThr, float etaThr)
static const std::vector< std::string > dnnAvaibleInputs
static constexpr float ptThreshold
static constexpr float ecalBarrelMaxEtaWithGap
const egammaTools::EgammaDNNHelper dnnHelper_
static constexpr float ecalBarrelMaxEtaNoGap

Member Function Documentation

std::vector< std::vector< float > > ElectronDNNEstimator::evaluate ( const reco::GsfElectronCollection ele,
const std::vector< tensorflow::Session * > &  sessions 
) const

Definition at line 144 of file ElectronDNNEstimator.cc.

References dnnHelper_, egammaTools::EgammaDNNHelper::evaluate(), getInputsVars(), and PixelMapPlotter::inputs.

145  {
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
std::map< std::string, float > getInputsVars(const reco::GsfElectron &ele) const
const egammaTools::EgammaDNNHelper dnnHelper_
std::map< std::string, float > ElectronDNNEstimator::getInputsVars ( const reco::GsfElectron ele) const

Definition at line 87 of file ElectronDNNEstimator.cc.

References funct::abs(), reco::GsfElectron::closestCtfTrackNLayers(), reco::GsfElectron::closestCtfTrackNormChi2(), reco::GsfElectron::convVtxFitProb(), reco::GsfElectron::deltaEtaSeedClusterTrackAtCalo(), reco::GsfElectron::deltaEtaSuperClusterTrackAtVtx(), reco::GsfElectron::deltaPhiSuperClusterTrackAtVtx(), reco::GsfElectron::dr03EcalRecHitSumEt(), reco::GsfElectron::dr03HcalTowerSumEt(), reco::GsfElectron::dr03TkSumPt(), reco::GsfElectron::ecalEnergy(), reco::GsfElectron::ecalEnergyError(), reco::GsfElectron::ecalPFClusterIso(), reco::GsfElectron::eEleClusterOverPout(), reco::GsfElectron::eSuperClusterOverP(), reco::LeafCandidate::eta(), reco::GsfElectron::fbrem(), reco::GsfElectron::full5x5_e1x5(), reco::GsfElectron::full5x5_e2x5Max(), reco::GsfElectron::full5x5_e5x5(), reco::GsfElectron::full5x5_hcalOverEcal(), reco::GsfElectron::full5x5_hcalOverEcalValid(), reco::GsfElectron::full5x5_r9(), reco::GsfElectron::full5x5_sigmaIetaIeta(), reco::GsfElectron::full5x5_sigmaIphiIphi(), reco::GsfElectron::gsfTrack(), reco::GsfElectron::hadronicOverEm(), reco::GsfElectron::hcalOverEcalValid(), reco::GsfElectron::hcalPFClusterIso(), reco::HitPattern::MISSING_INNER_HITS, reco::GsfElectron::numberOfBrems(), reco::LeafCandidate::pt(), reco::GsfElectron::superCluster(), reco::GsfElectron::superClusterFbrem(), reco::GsfElectron::trackMomentumAtVtx(), reco::GsfElectron::trackMomentumError(), and L1TEGammaDiff_cfi::variables.

Referenced by evaluate().

87  {
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 }
float trackMomentumError() const
Definition: GsfElectron.h:884
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
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
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
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
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
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
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
double eta() const final
momentum pseudorapidity
std::vector< tensorflow::Session * > ElectronDNNEstimator::getSessions ( ) const

Definition at line 39 of file ElectronDNNEstimator.cc.

References dnnHelper_, and egammaTools::EgammaDNNHelper::getSessions().

39 { return dnnHelper_.getSessions(); };
std::vector< tensorflow::Session * > getSessions() const
const egammaTools::EgammaDNNHelper dnnHelper_

Member Data Documentation

const std::vector< std::string > ElectronDNNEstimator::dnnAvaibleInputs
static

Definition at line 29 of file ElectronDNNEstimator.h.

const egammaTools::EgammaDNNHelper ElectronDNNEstimator::dnnHelper_
private

Definition at line 36 of file ElectronDNNEstimator.h.

Referenced by evaluate(), and getSessions().

constexpr float ElectronDNNEstimator::ecalBarrelMaxEtaNoGap = 1.485
static

Definition at line 33 of file ElectronDNNEstimator.h.

constexpr float ElectronDNNEstimator::ecalBarrelMaxEtaWithGap = 1.566
static

Definition at line 32 of file ElectronDNNEstimator.h.

constexpr float ElectronDNNEstimator::ptThreshold = 10.
static

Definition at line 31 of file ElectronDNNEstimator.h.

const bool ElectronDNNEstimator::useEBModelInGap_
private

Definition at line 38 of file ElectronDNNEstimator.h.