Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
RecoEgamma
ElectronIdentification
plugins
cuts
GsfEleHadronicOverEMLinearCut.cc
Go to the documentation of this file.
1
#include "
PhysicsTools/SelectorUtils/interface/CutApplicatorBase.h
"
2
#include "
DataFormats/EgammaCandidates/interface/GsfElectron.h
"
3
#include "
RecoEgamma/ElectronIdentification/interface/EBEECutValues.h
"
4
class
GsfEleHadronicOverEMLinearCut
:
public
CutApplicatorBase
{
5
public
:
6
GsfEleHadronicOverEMLinearCut
(
const
edm::ParameterSet
& params) :
7
CutApplicatorBase
(params),
8
slopeTerm_
(params,
"slopeTerm"
),
9
slopeStart_
(params,
"slopeStart"
),
10
constTerm_
(params,
"constTerm"
){}
11
12
result_type
operator()
(
const
reco::GsfElectronPtr
&)
const
override final;
13
14
CandidateType
candidateType
()
const
override final {
15
return
ELECTRON
;
16
}
17
18
private
:
19
EBEECutValues
slopeTerm_
;
20
EBEECutValues
slopeStart_
;
21
EBEECutValues
constTerm_
;
22
23
};
24
25
DEFINE_EDM_PLUGIN
(
CutApplicatorFactory
,
26
GsfEleHadronicOverEMLinearCut
,
27
"GsfEleHadronicOverEMLinearCut"
);
28
29
CutApplicatorBase::result_type
30
GsfEleHadronicOverEMLinearCut::
31
operator()
(
const
reco::GsfElectronPtr
& cand)
const
{
32
33
const
float
energy
= cand->superCluster()->energy();
34
const
float
cutValue = energy >
slopeStart_
(cand) ?
slopeTerm_
(cand)*(energy-
slopeStart_
(cand)) +
constTerm_
(cand) :
constTerm_
(cand);
35
36
return
cand->hadronicOverEm()*energy < cutValue;
37
}
GsfEleHadronicOverEMLinearCut
Definition:
GsfEleHadronicOverEMLinearCut.cc:4
CutApplicatorBase::ELECTRON
Definition:
CutApplicatorBase.h:48
GsfEleHadronicOverEMLinearCut::slopeTerm_
EBEECutValues slopeTerm_
Definition:
GsfEleHadronicOverEMLinearCut.cc:19
GsfEleHadronicOverEMLinearCut::GsfEleHadronicOverEMLinearCut
GsfEleHadronicOverEMLinearCut(const edm::ParameterSet ¶ms)
Definition:
GsfEleHadronicOverEMLinearCut.cc:6
EBEECutValues
Definition:
EBEECutValues.h:12
edmplugin::PluginFactory
Definition:
PluginFactory.h:31
GsfEleHadronicOverEMLinearCut::operator()
result_type operator()(const reco::GsfElectronPtr &) const overridefinal
Definition:
GsfEleHadronicOverEMLinearCut.cc:31
GsfEleHadronicOverEMLinearCut::candidateType
CandidateType candidateType() const overridefinal
Definition:
GsfEleHadronicOverEMLinearCut.cc:14
edm::Ptr< reco::GsfElectron >
GsfElectron.h
CutApplicatorBase
Definition:
CutApplicatorBase.h:45
EBEECutValues.h
GsfEleHadronicOverEMLinearCut::constTerm_
EBEECutValues constTerm_
Definition:
GsfEleHadronicOverEMLinearCut.cc:21
compareJSON.const
string const
Definition:
compareJSON.py:14
GsfEleHadronicOverEMLinearCut::slopeStart_
EBEECutValues slopeStart_
Definition:
GsfEleHadronicOverEMLinearCut.cc:20
CutApplicatorBase::CandidateType
CandidateType
Definition:
CutApplicatorBase.h:47
edm::ParameterSet
Definition:
ParameterSet.h:35
relval_parameters_module.energy
string energy
Definition:
relval_parameters_module.py:29
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition:
PluginFactory.h:101
CutApplicatorBase.h
Generated for CMSSW Reference Manual by
1.8.5