CMS 3D CMS Logo

ZElectronsSelector.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 // Class: ZElectronsSelector
3 //
4 // Original Author: Silvia Taroni
5 // Created: Wed, 29 Nov 2017 18:23:54 GMT
6 //
7 //
8 
10 
11 // system include files
12 
13 #include <algorithm>
14 #include <iostream>
15 #include <memory>
16 #include <string>
17 #include <vector>
18 
19 // user include files
22 
25 
27 
33 //
38 
39 using namespace std;
40 using namespace reco;
41 namespace edm {
42  class EventSetup;
43 }
44 
46 public:
48  bool operator()(const reco::GsfElectron&) const;
49  void newEvent(const edm::Event&, const edm::EventSetup&);
50  const float getEffectiveArea(float eta) const;
51  void printEffectiveAreas() const;
52 
56 
57  std::vector<double> absEtaMin_; // low limit of the eta range
58  std::vector<double> absEtaMax_; // upper limit of the eta range
59  std::vector<double> effectiveAreaValues_; // effective area for this eta range
60 
62 
63  vector<int> missHits;
64  vector<double> sigmaIEtaIEtaCut;
65  vector<double> dEtaInSeedCut;
66  vector<double> dPhiInCut;
67  vector<double> hOverECut;
68  vector<double> relCombIso;
69  vector<double> EInvMinusPInv;
70 };
71 
73  printf(" eta_min eta_max effective area\n");
74  uint nEtaBins = absEtaMin_.size();
75  for (uint iEta = 0; iEta < nEtaBins; iEta++) {
76  printf(" %8.4f %8.4f %8.5f\n", absEtaMin_[iEta], absEtaMax_[iEta], effectiveAreaValues_[iEta]);
77  }
78 }
79 const float ZElectronsSelector::getEffectiveArea(float eta) const {
80  float effArea = 0;
81  uint nEtaBins = absEtaMin_.size();
82  for (uint iEta = 0; iEta < nEtaBins; iEta++) {
83  if (std::abs(eta) >= absEtaMin_[iEta] && std::abs(eta) < absEtaMax_[iEta]) {
84  effArea = effectiveAreaValues_[iEta];
85  break;
86  }
87  }
88 
89  return effArea;
90 }
91 
93  : theRhoToken(iC.consumes<double>(cfg.getParameter<edm::InputTag>("rho"))) {
94  absEtaMin_ = cfg.getParameter<std::vector<double> >("absEtaMin");
95  absEtaMax_ = cfg.getParameter<std::vector<double> >("absEtaMax");
96  effectiveAreaValues_ = cfg.getParameter<std::vector<double> >("effectiveAreaValues");
97  //printEffectiveAreas();
98  eleIDWP = cfg.getParameter<edm::ParameterSet>("eleID");
99 
100  missHits = eleIDWP.getParameter<std::vector<int> >("missingHitsCut");
101  sigmaIEtaIEtaCut = eleIDWP.getParameter<std::vector<double> >("full5x5_sigmaIEtaIEtaCut");
102  dEtaInSeedCut = eleIDWP.getParameter<std::vector<double> >("dEtaInSeedCut");
103  dPhiInCut = eleIDWP.getParameter<std::vector<double> >("dPhiInCut");
104  hOverECut = eleIDWP.getParameter<std::vector<double> >("hOverECut");
105  relCombIso = eleIDWP.getParameter<std::vector<double> >("relCombIsolationWithEACut");
106  EInvMinusPInv = eleIDWP.getParameter<std::vector<double> >("EInverseMinusPInverseCut");
107 }
108 
111 }
112 
114  float pt_e = el.pt();
115  unsigned int ind = 0;
116  auto etrack = el.gsfTrack();
117  float abseta = fabs((el.superCluster().get())->position().eta());
118 
119  if (el.isEB()) {
120  if (abseta > 1.479)
121  return false; // check if it is really needed
122  }
123  if (el.isEE()) {
124  ind = 1;
125  if (abseta < 1.479)
126  return false; // check if it is really needed
127  if (abseta >= 2.5)
128  return false; // check if it is really needed
129  }
130 
131  if (etrack->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > missHits[ind])
132  return false;
133  if (el.full5x5_sigmaIetaIeta() > sigmaIEtaIEtaCut[ind])
134  return false;
135  if (fabs(el.deltaPhiSuperClusterTrackAtVtx()) > dPhiInCut[ind])
136  return false;
137  if (fabs(el.deltaEtaSeedClusterTrackAtVtx()) > dEtaInSeedCut[ind])
138  return false;
139  if (el.hadronicOverEm() > hOverECut[ind])
140  return false;
141  const float eA = getEffectiveArea(abseta);
142  const float rho = _rhoHandle.isValid() ? (float)(*_rhoHandle.product()) : 0;
144  std::max(float(0.0),
146  relCombIso[ind] * pt_e)
147  return false;
148  const float ecal_energy_inverse = 1.0 / el.ecalEnergy();
149  const float eSCoverP = el.eSuperClusterOverP();
150  if (std::abs(1.0 - eSCoverP) * ecal_energy_inverse > EInvMinusPInv[ind])
151  return false;
152 
153  return true;
154 }
155 
157 
159 
const PflowIsolationVariables & pfIsolationVariables() const
Definition: GsfElectron.h:650
T getParameter(std::string const &) const
GsfTrackRef gsfTrack() const override
reference to a GsfTrack
Definition: GsfElectron.h:156
float eSuperClusterOverP() const
Definition: GsfElectron.h:221
vector< double > EInvMinusPInv
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::Handle< double > _rhoHandle
std::vector< double > absEtaMax_
vector< double > dEtaInSeedCut
std::vector< double > effectiveAreaValues_
double pt() const final
transverse momentum
bool ev
bool operator()(const reco::GsfElectron &) const
constexpr int nEtaBins
Definition: Common.h:11
bool isEE() const
Definition: GsfElectron.h:329
bool isEB() const
Definition: GsfElectron.h:328
float full5x5_sigmaIetaIeta() const
Definition: GsfElectron.h:435
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float hadronicOverEm() const
Definition: GsfElectron.h:468
edm::EDGetTokenT< double > theRhoToken
float sumPhotonEt
sum pt of PF photons // old float photonIso ;
Definition: GsfElectron.h:602
float deltaPhiSuperClusterTrackAtVtx() const
Definition: GsfElectron.h:228
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define EVENTSETUP_STD_INIT(SELECTOR)
float sumNeutralHadronEt
sum pt of neutral hadrons // old float neutralHadronIso ;
Definition: GsfElectron.h:601
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
bool isValid() const
Definition: HandleBase.h:70
float deltaEtaSeedClusterTrackAtVtx() const
Definition: GsfElectron.h:231
vector< double > relCombIso
SingleObjectSelector< edm::View< reco::GsfElectron >, ZElectronsSelector > ZElectronsSelectorAndSkim
T const * product() const
Definition: Handle.h:69
void newEvent(const edm::Event &, const edm::EventSetup &)
float ecalEnergy() const
Definition: GsfElectron.h:812
void printEffectiveAreas() const
edm::ParameterSet eleIDWP
fixed size matrix
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:289
SuperClusterRef superCluster() const override
reference to a SuperCluster
Definition: GsfElectron.h:155
vector< double > dPhiInCut
float sumChargedHadronPt
sum-pt of charged Hadron // old float chargedHadronIso ;
Definition: GsfElectron.h:600
std::vector< double > absEtaMin_
edm::EDGetTokenT< reco::GsfElectronCollection > theGsfEToken
const float getEffectiveArea(float eta) const
vector< double > hOverECut
ZElectronsSelector(const edm::ParameterSet &, edm::ConsumesCollector &iC)
vector< double > sigmaIEtaIEtaCut