CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFMETAlgo.cc
Go to the documentation of this file.
2 
4 
8 
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 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
virtual double energy() const
energy
Definition: MET.h:32
T sqrt(T t)
Definition: SSEVec.h:28
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual double theta() const
momentum polar angle
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
tuple cout
Definition: gather_cfg.py:41
virtual ParticleType particleId() const
Definition: PFCandidate.h:314
virtual double phi() const
momentum azimuthal angle
Definition: DDAxes.h:10