CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoMET/METAlgorithms/src/PFSpecificAlgo.cc

Go to the documentation of this file.
00001 /*
00002 class: PFSpecificAlgo.cc
00003 description:  MET made from Particle Flow candidates
00004 authors: R. Remington (UF), R. Cavanaugh (UIC/Fermilab)
00005   date: 10/27/08
00006 */
00007 
00008 #include "DataFormats/Math/interface/LorentzVector.h"
00009 #include "RecoMET/METAlgorithms/interface/PFSpecificAlgo.h"
00010 #include "RecoMET/METAlgorithms/interface/significanceAlgo.h"
00011 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00012 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00013 using namespace reco;
00014 using namespace std;
00015 
00016 //--------------------------------------------------------------------------------------
00017 // This algorithm adds Particle Flow specific global event information to the MET object
00018 //--------------------------------------------------------------------------------------
00019 
00020 reco::PFMET PFSpecificAlgo::addInfo(edm::Handle<edm::View<Candidate> > PFCandidates, CommonMETData met)
00021 {
00022   // Instantiate the container to hold the PF specific information
00023   SpecificPFMETData specific;
00024   // Initialize the container
00025   specific.NeutralEMFraction = 0.0;
00026   specific.NeutralHadFraction = 0.0;
00027   specific.ChargedEMFraction = 0.0;
00028   specific.ChargedHadFraction = 0.0;
00029   specific.MuonFraction = 0.0;
00030   specific.Type6Fraction = 0.0;
00031   specific.Type7Fraction = 0.0;
00032 
00033   if(!PFCandidates->size()) // if no Particle Flow candidates in the event
00034   {
00035     const LorentzVector p4( met.mex, met.mey, 0.0, met.met);
00036     const Point vtx(0.0, 0.0, 0.0 );
00037     PFMET specificPFMET( specific, met.sumet, p4, vtx);
00038     return specificPFMET;
00039   } 
00040 
00041   double NeutralEMEt = 0.0;
00042   double NeutralHadEt = 0.0;
00043   double ChargedEMEt = 0.0;
00044   double ChargedHadEt = 0.0;
00045   double MuonEt = 0.0;
00046   double type6Et = 0.0;
00047   double type7Et = 0.0;
00048   
00049   for( edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin() ; iParticle != (PFCandidates.product())->end() ; ++iParticle )
00050   {   
00051     const Candidate* candidate = &(*iParticle);
00052     if (candidate) {
00053       //const PFCandidate* pfCandidate = static_cast<const PFCandidate*> (candidate);
00054       const PFCandidate* pfCandidate = dynamic_cast<const PFCandidate*> (candidate);
00055       if (pfCandidate)
00056       {
00057         //cout << pfCandidate->et() << "     "
00058         //   << pfCandidate->hcalEnergy() << "    "
00059         //   << pfCandidate->ecalEnergy() << endl;
00060         //std::cout << "pfCandidate->particleId() = " << pfCandidate->particleId() << std::endl;
00061         const double theta = iParticle->theta();
00062         const double e     = iParticle->energy();
00063         const double et    = e*sin(theta);
00064         if(alsocalcsig){
00065             reco::CandidatePtr dau(PFCandidates, iParticle - PFCandidates->begin());
00066             if(dau.isNonnull () && dau.isAvailable()){
00067                 reco::PFCandidatePtr pf(dau.id(), pfCandidate, dau.key());
00068                 pfsignalgo_.addPFCandidate(pf);
00069             }
00070           //metSigInputVector.push_back(resolutions_.evalPF(pfCandidate));
00071         }
00072 
00073         if (pfCandidate->particleId() == 1) ChargedHadEt += et;
00074         if (pfCandidate->particleId() == 2) ChargedEMEt += et;
00075         if (pfCandidate->particleId() == 3) MuonEt += et;
00076         if (pfCandidate->particleId() == 4) NeutralEMEt += et;
00077         if (pfCandidate->particleId() == 5) NeutralHadEt += et;
00078         if (pfCandidate->particleId() == 6) type6Et += et;
00079         if (pfCandidate->particleId() == 7) type7Et += et;
00080       }
00081     } 
00082   }
00083 
00084   const double Et_total=NeutralEMEt+NeutralHadEt+ChargedEMEt+ChargedHadEt+MuonEt+type6Et+type7Et;
00085 
00086   if (Et_total!=0.0)
00087   {
00088     specific.NeutralEMFraction = NeutralEMEt/Et_total;
00089     specific.NeutralHadFraction = NeutralHadEt/Et_total;
00090     specific.ChargedEMFraction = ChargedEMEt/Et_total;
00091     specific.ChargedHadFraction = ChargedHadEt/Et_total;
00092     specific.MuonFraction = MuonEt/Et_total;
00093     specific.Type6Fraction = type6Et/Et_total;
00094     specific.Type7Fraction = type7Et/Et_total;
00095   }
00096   
00097   const LorentzVector p4(met.mex , met.mey, 0.0, met.met);
00098   const Point vtx(0.0,0.0,0.0);
00099   PFMET specificPFMET( specific, met.sumet, p4, vtx );
00100 
00101   specificPFMET.setSignificanceMatrix(pfsignalgo_.getSignifMatrix());
00102 
00103   return specificPFMET;
00104 }
00105 
00106 void PFSpecificAlgo::runSignificance(metsig::SignAlgoResolutions &resolutions, edm::Handle<edm::View<reco::PFJet> > jets)
00107 {
00108   alsocalcsig=true;
00109   pfsignalgo_.setResolutions( &resolutions );
00110   pfsignalgo_.addPFJets(jets);
00111 }