CMS 3D CMS Logo

PFMETAlgo.cc
Go to the documentation of this file.
2 
4 
8 
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++) {
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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T getUntrackedParameter(std::string const &, T const &) const
Definition: MET.h:41
T sqrt(T t)
Definition: SSEVec.h:23
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
fixed size matrix
HLT enums.
Geom::Theta< T > theta() const
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25