CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFSpecificAlgo.cc
Go to the documentation of this file.
1 /*
2 class: PFSpecificAlgo.cc
3 description: MET made from Particle Flow candidates
4 authors: R. Remington (UF), R. Cavanaugh (UIC/Fermilab)
5  date: 10/27/08
6 */
7 
13 using namespace reco;
14 using namespace std;
15 
16 //--------------------------------------------------------------------------------------
17 // This algorithm adds Particle Flow specific global event information to the MET object
18 //--------------------------------------------------------------------------------------
19 
21 {
22  // Instantiate the container to hold the PF specific information
24  // Initialize the container
25  specific.NeutralEMFraction = 0.0;
26  specific.NeutralHadFraction = 0.0;
27  specific.ChargedEMFraction = 0.0;
28  specific.ChargedHadFraction = 0.0;
29  specific.MuonFraction = 0.0;
30  specific.Type6Fraction = 0.0;
31  specific.Type7Fraction = 0.0;
32 
33  if(!PFCandidates->size()) // if no Particle Flow candidates in the event
34  {
35  const LorentzVector p4( met.mex, met.mey, 0.0, met.met);
36  const Point vtx(0.0, 0.0, 0.0 );
37  PFMET specificPFMET( specific, met.sumet, p4, vtx);
38  return specificPFMET;
39  }
40 
41  double NeutralEMEt = 0.0;
42  double NeutralHadEt = 0.0;
43  double ChargedEMEt = 0.0;
44  double ChargedHadEt = 0.0;
45  double MuonEt = 0.0;
46  double type6Et = 0.0;
47  double type7Et = 0.0;
48 
49 
50  pfsignalgo_.useOriginalPtrs(PFCandidates.id());
51 
52  for( edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin() ; iParticle != (PFCandidates.product())->end() ; ++iParticle )
53  {
54  const Candidate* candidate = &(*iParticle);
55  if (candidate) {
56  //const PFCandidate* pfCandidate = static_cast<const PFCandidate*> (candidate);
57  const PFCandidate* pfCandidate = dynamic_cast<const PFCandidate*> (candidate);
58  if (pfCandidate)
59  {
60  //cout << pfCandidate->et() << " "
61  // << pfCandidate->hcalEnergy() << " "
62  // << pfCandidate->ecalEnergy() << endl;
63  //std::cout << "pfCandidate->particleId() = " << pfCandidate->particleId() << std::endl;
64  const double theta = iParticle->theta();
65  const double e = iParticle->energy();
66  const double et = e*sin(theta);
67  if(alsocalcsig){
68  reco::CandidatePtr dau(PFCandidates, iParticle - PFCandidates->begin());
69  if(dau.isNonnull () && dau.isAvailable()){
70  reco::PFCandidatePtr pf(dau.id(), pfCandidate, dau.key());
71  pfsignalgo_.addPFCandidate(pf);
72  }
73  //metSigInputVector.push_back(resolutions_.evalPF(pfCandidate));
74  }
75 
76  if (pfCandidate->particleId() == 1) ChargedHadEt += et;
77  if (pfCandidate->particleId() == 2) ChargedEMEt += et;
78  if (pfCandidate->particleId() == 3) MuonEt += et;
79  if (pfCandidate->particleId() == 4) NeutralEMEt += et;
80  if (pfCandidate->particleId() == 5) NeutralHadEt += et;
81  if (pfCandidate->particleId() == 6) type6Et += et;
82  if (pfCandidate->particleId() == 7) type7Et += et;
83  }
84  }
85  }
86 
87  const double Et_total=NeutralEMEt+NeutralHadEt+ChargedEMEt+ChargedHadEt+MuonEt+type6Et+type7Et;
88 
89  if (Et_total!=0.0)
90  {
91  specific.NeutralEMFraction = NeutralEMEt/Et_total;
92  specific.NeutralHadFraction = NeutralHadEt/Et_total;
93  specific.ChargedEMFraction = ChargedEMEt/Et_total;
94  specific.ChargedHadFraction = ChargedHadEt/Et_total;
95  specific.MuonFraction = MuonEt/Et_total;
96  specific.Type6Fraction = type6Et/Et_total;
97  specific.Type7Fraction = type7Et/Et_total;
98  }
99 
100  const LorentzVector p4(met.mex , met.mey, 0.0, met.met);
101  const Point vtx(0.0,0.0,0.0);
102  PFMET specificPFMET( specific, met.sumet, p4, vtx );
103 
104  specificPFMET.setSignificanceMatrix(pfsignalgo_.getSignifMatrix());
105 
106  return specificPFMET;
107 }
108 
110 {
111  alsocalcsig=true;
112  pfsignalgo_.setResolutions( &resolutions );
113  pfsignalgo_.addPFJets(jets);
114 }
dictionary specific
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
math::XYZPoint Point
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
double p4[4]
Definition: TauolaWrapper.h:92
vector< PseudoJet > jets
Structure containing data common to all types of MET.
Definition: CommonMETData.h:22
#define end
Definition: vmac.h:38
void setSignificanceMatrix(const TMatrixD &matrix)
Definition: MET.cc:185
math::XYZTLorentzVector LorentzVector
MET made from Particle Flow Candidates.
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:33
#define begin
Definition: vmac.h:31
void runSignificance(metsig::SignAlgoResolutions &resolutions, edm::Handle< edm::View< reco::PFJet > > jets)
reco::PFMET addInfo(edm::Handle< edm::View< reco::Candidate > > PFCandidates, CommonMETData met)
virtual ParticleType particleId() const
Definition: PFCandidate.h:324
tuple PFCandidates