CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions
PFSpecificAlgo Class Reference

#include <RecoMET/METAlgorithms/interface/PFSpecificAlgo.h>

Public Member Functions

 PFSpecificAlgo ()
 
SpecificPFMETData run (const edm::View< reco::Candidate > &pfCands)
 

Detailed Description

Description: Adds Particle Flow specific information to MET

Implementation: [Notes on implementation]

Definition at line 27 of file PFSpecificAlgo.h.

Constructor & Destructor Documentation

PFSpecificAlgo::PFSpecificAlgo ( )
inline

Definition at line 30 of file PFSpecificAlgo.h.

30 { }

Member Function Documentation

SpecificPFMETData PFSpecificAlgo::run ( const edm::View< reco::Candidate > &  pfCands)

Definition at line 15 of file PFSpecificAlgo.cc.

References edm::View< T >::begin(), SpecificPFMETData::ChargedEMFraction, SpecificPFMETData::ChargedHadFraction, alignCSCRings::e, edm::View< T >::end(), reco::LeafCandidate::energy(), SpecificPFMETData::MuonFraction, SpecificPFMETData::NeutralEMFraction, SpecificPFMETData::NeutralHadFraction, reco::PFCandidate::particleId(), funct::sin(), edm::View< T >::size(), timingPdfMaker::specific, reco::LeafCandidate::theta(), theta(), SpecificPFMETData::Type6Fraction, and SpecificPFMETData::Type7Fraction.

Referenced by cms::PFMETProducer::produce(), and cms::METProducer::produce_PFMET().

16 {
17  if(!pfCands.size()) return SpecificPFMETData();
18 
19  double NeutralEMEt = 0.0;
20  double NeutralHadEt = 0.0;
21  double ChargedEMEt = 0.0;
22  double ChargedHadEt = 0.0;
23  double MuonEt = 0.0;
24  double type6Et = 0.0;
25  double type7Et = 0.0;
26 
27  for( edm::View<reco::Candidate>::const_iterator iPfCand = pfCands.begin(); iPfCand != pfCands.end(); ++iPfCand)
28  {
29  const reco::PFCandidate* pfCand = dynamic_cast<const reco::PFCandidate*> (&(*iPfCand));
30  if (!pfCand) continue;
31  const double theta = pfCand->theta();
32  const double e = pfCand->energy();
33  const double et = e*sin(theta);
34 
35  if (pfCand->particleId() == 1) ChargedHadEt += et;
36  if (pfCand->particleId() == 2) ChargedEMEt += et;
37  if (pfCand->particleId() == 3) MuonEt += et;
38  if (pfCand->particleId() == 4) NeutralEMEt += et;
39  if (pfCand->particleId() == 5) NeutralHadEt += et;
40  if (pfCand->particleId() == 6) type6Et += et;
41  if (pfCand->particleId() == 7) type7Et += et;
42  }
43 
44  const double Et_total = NeutralEMEt + NeutralHadEt + ChargedEMEt + ChargedHadEt + MuonEt + type6Et + type7Et;
46  if (Et_total!=0.0)
47  {
48  specific.NeutralEMFraction = NeutralEMEt/Et_total;
49  specific.NeutralHadFraction = NeutralHadEt/Et_total;
50  specific.ChargedEMFraction = ChargedEMEt/Et_total;
51  specific.ChargedHadFraction = ChargedHadEt/Et_total;
52  specific.MuonFraction = MuonEt/Et_total;
53  specific.Type6Fraction = type6Et/Et_total;
54  specific.Type7Fraction = type7Et/Et_total;
55  }
56  return specific;
57 }
virtual double energy() const GCC11_FINAL
energy
dictionary specific
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
virtual double theta() const GCC11_FINAL
momentum polar angle
MET made from Particle Flow Candidates.
size_type size() const
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
const_iterator begin() const
const_iterator end() const
virtual ParticleType particleId() const
Definition: PFCandidate.h:355