CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Attributes
pf2pat::PFMETAlgo Class Reference

#include <PFMETAlgo.h>

Public Member Functions

 PFMETAlgo (const edm::ParameterSet &)
 
reco::MET produce (const reco::PFCandidateCollection &pfCandidates)
 
 ~PFMETAlgo ()
 

Private Attributes

double hfCalibFactor_
 HF calibration factor (in 31X applied by PFProducer) More...
 
bool verbose_
 verbose ? More...
 

Detailed Description

Definition at line 30 of file PFMETAlgo.h.

Constructor & Destructor Documentation

PFMETAlgo::PFMETAlgo ( const edm::ParameterSet iConfig)
explicit

Definition at line 17 of file PFMETAlgo.cc.

References edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

17  {
18  verbose_ = iConfig.getUntrackedParameter<bool>("verbose", false);
19 
20  hfCalibFactor_ = iConfig.getParameter<double>("hfCalibFactor");
21 }
T getUntrackedParameter(std::string const &, T const &) const
double hfCalibFactor_
HF calibration factor (in 31X applied by PFProducer)
Definition: PFMETAlgo.h:40
bool verbose_
verbose ?
Definition: PFMETAlgo.h:43
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
PFMETAlgo::~PFMETAlgo ( )

Definition at line 23 of file PFMETAlgo.cc.

23 {}

Member Function Documentation

reco::MET PFMETAlgo::produce ( const reco::PFCandidateCollection pfCandidates)

HF calibration factor (in 31X applied by PFProducer)

Definition at line 25 of file PFMETAlgo.cc.

References funct::cos(), gather_cfg::cout, reco::LeafCandidate::energy(), mps_fire::i, HLT_FULL_cff::MET, reco::PFCandidate::particleId(), reco::LeafCandidate::phi(), funct::sin(), mathSSE::sqrt(), theta(), and reco::LeafCandidate::theta().

25  {
26  double sumEx = 0;
27  double sumEy = 0;
28  double sumEt = 0;
29 
30  for (unsigned i = 0; i < pfCandidates.size(); i++) {
31  const reco::PFCandidate& cand = pfCandidates[i];
32 
33  double E = cand.energy();
34 
36  if (cand.particleId() == PFCandidate::h_HF || cand.particleId() == PFCandidate::egamma_HF)
37  E *= hfCalibFactor_;
38 
39  double phi = cand.phi();
40  double cosphi = cos(phi);
41  double sinphi = sin(phi);
42 
43  double theta = cand.theta();
44  double sintheta = sin(theta);
45 
46  double et = E * sintheta;
47  double ex = et * cosphi;
48  double ey = et * sinphi;
49 
50  sumEx += ex;
51  sumEy += ey;
52  sumEt += et;
53  }
54 
55  double Et = sqrt(sumEx * sumEx + sumEy * sumEy);
56  XYZTLorentzVector missingEt(-sumEx, -sumEy, 0, Et);
57 
58  if (verbose_) {
59  cout << "PFMETAlgo: mEx, mEy, mEt = " << missingEt.X() << ", " << missingEt.Y() << ", " << missingEt.T() << endl;
60  }
61 
62  XYZPoint vertex; // dummy vertex
63  return MET(sumEt, missingEt, vertex);
64 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
double hfCalibFactor_
HF calibration factor (in 31X applied by PFProducer)
Definition: PFMETAlgo.h:40
double theta() const final
momentum polar angle
T sqrt(T t)
Definition: SSEVec.h:19
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
bool verbose_
verbose ?
Definition: PFMETAlgo.h:43
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
tuple cout
Definition: gather_cfg.py:144
virtual ParticleType particleId() const
Definition: PFCandidate.h:392
double phi() const final
momentum azimuthal angle
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
double energy() const final
energy

Member Data Documentation

double pf2pat::PFMETAlgo::hfCalibFactor_
private

HF calibration factor (in 31X applied by PFProducer)

Definition at line 40 of file PFMETAlgo.h.

bool pf2pat::PFMETAlgo::verbose_
private

verbose ?

Definition at line 43 of file PFMETAlgo.h.