Main Page
Namespaces
Classes
Package Documentation
CommonTools
ParticleFlow
src
PFMETAlgo.cc
Go to the documentation of this file.
1
#include "
CommonTools/ParticleFlow/interface/PFMETAlgo.h
"
2
3
#include "
DataFormats/ParticleFlowCandidate/interface/PFCandidate.h
"
4
5
#include "
DataFormats/METReco/interface/MET.h
"
6
#include "
DataFormats/METReco/interface/METFwd.h
"
7
#include "
DataFormats/Math/interface/LorentzVector.h
"
8
9
#include "
FWCore/Framework/interface/EventSetup.h
"
10
11
12
using namespace
std
;
13
using namespace
edm
;
14
using namespace
reco
;
15
using namespace
math
;
16
using namespace
pf2pat
;
17
18
PFMETAlgo::PFMETAlgo(
const
edm::ParameterSet
& iConfig) {
19
20
21
verbose_ =
22
iConfig.
getUntrackedParameter
<
bool
>(
"verbose"
,
false
);
23
24
hfCalibFactor_ =
25
iConfig.
getParameter
<
double
>(
"hfCalibFactor"
);
26
27
}
28
29
30
31
PFMETAlgo::~PFMETAlgo() { }
32
33
reco::MET
PFMETAlgo::produce(
const
reco::PFCandidateCollection
&
pfCandidates
) {
34
35
double
sumEx = 0;
36
double
sumEy = 0;
37
double
sumEt
= 0;
38
39
for
(
unsigned
i
=0;
i
<pfCandidates.size();
i
++ ) {
40
41
const
reco::PFCandidate
&
cand
= pfCandidates[
i
];
42
43
double
E = cand.
energy
();
44
46
if
( cand.
particleId
()==PFCandidate::h_HF ||
47
cand.
particleId
()==PFCandidate::egamma_HF )
48
E *= hfCalibFactor_;
49
50
double
phi = cand.
phi
();
51
double
cosphi =
cos
(phi);
52
double
sinphi =
sin
(phi);
53
54
double
theta
= cand.
theta
();
55
double
sintheta =
sin
(theta);
56
57
double
et
= E*sintheta;
58
double
ex = et*cosphi;
59
double
ey = et*sinphi;
60
61
sumEx += ex;
62
sumEy += ey;
63
sumEt +=
et
;
64
}
65
66
double
Et =
sqrt
( sumEx*sumEx + sumEy*sumEy);
67
XYZTLorentzVector
missingEt( -sumEx, -sumEy, 0, Et);
68
69
if
(verbose_) {
70
cout
<<
"PFMETAlgo: mEx, mEy, mEt = "
71
<< missingEt.X() <<
", "
72
<< missingEt.Y() <<
", "
73
<< missingEt.T() <<endl;
74
}
75
76
XYZPoint
vertex;
// dummy vertex
77
return
MET
(sumEt, missingEt, vertex);
78
79
80
}
81
82
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
mps_fire.i
i
Definition:
mps_fire.py:330
reco::LeafCandidate::theta
double theta() const final
momentum polar angle
Definition:
LeafCandidate.h:135
MET.h
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
EventSetup.h
theta
Geom::Theta< T > theta() const
Definition:
Basic3DVectorLD.h:179
std
Definition:
JetResolutionObject.h:80
METFwd.h
PFMETAlgo.h
slimmedMuons_cfi.pfCandidates
pfCandidates
Definition:
slimmedMuons_cfi.py:6
PFCandidate.h
reco::MET
Definition:
MET.h:42
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:18
objects.METAnalyzer.sumEt
sumEt
Definition:
METAnalyzer.py:97
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
reco::LeafCandidate::energy
double energy() const final
energy
Definition:
LeafCandidate.h:110
LorentzVector.h
math
Definition:
Error.h:16
nanoDQM_cfi.MET
MET
Definition:
nanoDQM_cfi.py:315
reco::PFCandidateCollection
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
Definition:
PFCandidateFwd.h:12
pf2pat
Definition:
ElectronIDPFCandidateSelectorDefinition.h:22
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition:
Point3D.h:12
stringResolutionProvider_cfi.et
et
define resolution functions of each parameter
Definition:
stringResolutionProvider_cfi.py:13
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition:
PFCandidate.h:40
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:44
edm
HLT enums.
Definition:
AlignableModifier.h:17
edm::ParameterSet
Definition:
ParameterSet.h:36
cand
Definition:
decayParser.h:34
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
reco::PFCandidate::particleId
virtual ParticleType particleId() const
Definition:
PFCandidate.h:374
reco::LeafCandidate::phi
double phi() const final
momentum azimuthal angle
Definition:
LeafCandidate.h:133
XYZTLorentzVector
math::XYZTLorentzVector XYZTLorentzVector
Definition:
RawParticle.h:27
Generated for CMSSW Reference Manual by
1.8.11