src
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
using namespace
std
;
12
using namespace
edm
;
13
using namespace
reco
;
14
using namespace
math
;
15
using namespace
pf2pat
;
16
17
PFMETAlgo::PFMETAlgo(
const
edm::ParameterSet
& iConfig) {
18
verbose_ = iConfig.
getUntrackedParameter
<
bool
>(
"verbose"
,
false
);
19
20
hfCalibFactor_ = iConfig.
getParameter
<
double
>(
"hfCalibFactor"
);
21
}
22
23
reco::MET
PFMETAlgo::produce(
const
reco::PFCandidateCollection
&
pfCandidates
)
const
{
24
double
sumEx = 0;
25
double
sumEy = 0;
26
double
sumEt
= 0;
27
28
for
(
unsigned
i
= 0;
i
<
pfCandidates
.size();
i
++) {
29
const
reco::PFCandidate
&
cand
=
pfCandidates
[
i
];
30
31
double
E =
cand
.energy();
32
34
if
(
cand
.particleId() == PFCandidate::h_HF ||
cand
.particleId() == PFCandidate::egamma_HF)
35
E *= hfCalibFactor_;
36
37
double
phi =
cand
.phi();
38
double
cosphi =
cos
(phi);
39
double
sinphi =
sin
(phi);
40
41
double
theta
=
cand
.theta();
42
double
sintheta =
sin
(
theta
);
43
44
double
et
= E * sintheta;
45
double
ex =
et
* cosphi;
46
double
ey =
et
* sinphi;
47
48
sumEx += ex;
49
sumEy += ey;
50
sumEt
+=
et
;
51
}
52
53
double
Et =
sqrt
(sumEx * sumEx + sumEy * sumEy);
54
XYZTLorentzVector
missingEt(-sumEx, -sumEy, 0, Et);
55
56
if
(verbose_) {
57
cout
<<
"PFMETAlgo: mEx, mEy, mEt = "
<< missingEt.X() <<
", "
<< missingEt.Y() <<
", "
<< missingEt.T() << endl;
58
}
59
60
XYZPoint
vertex
;
// dummy vertex
61
return
MET
(
sumEt
, missingEt,
vertex
);
62
}
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition:
ParameterSet.h:303
EgHLTOffHistBins_cfi.et
et
Definition:
EgHLTOffHistBins_cfi.py:8
bphysicsOniaDQM_cfi.vertex
vertex
Definition:
bphysicsOniaDQM_cfi.py:7
mps_fire.i
i
Definition:
mps_fire.py:429
MET.h
funct::sin
Sin< T >::type sin(const T &t)
Definition:
Sin.h:22
EventSetup.h
std
Definition:
JetResolutionObject.h:76
METFwd.h
PFMETAlgo.h
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
reco::MET
Definition:
MET.h:41
mathSSE::sqrt
T sqrt(T t)
Definition:
SSEVec.h:19
objects.METAnalyzer.sumEt
sumEt
Definition:
METAnalyzer.py:97
funct::cos
Cos< T >::type cos(const T &t)
Definition:
Cos.h:22
LorentzVector.h
math
Definition:
choleskyInversion.h:8
zmumugammaAnalyzer_cfi.pfCandidates
pfCandidates
Definition:
zmumugammaAnalyzer_cfi.py:11
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
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition:
PFCandidate.h:41
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:46
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::ParameterSet
Definition:
ParameterSet.h:47
cand
Definition:
decayParser.h:32
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
PFCandidate.h
theta
Geom::Theta< T > theta() const
Definition:
Basic3DVectorLD.h:150
HLTTauDQMOffline_cfi.MET
MET
Definition:
HLTTauDQMOffline_cfi.py:64
XYZTLorentzVector
math::XYZTLorentzVector XYZTLorentzVector
Definition:
RawParticle.h:25
Generated for CMSSW Reference Manual by
1.8.14