#include <RecoMET/METAlgorithms/interface/PFSpecificAlgo.h>
Public Member Functions | |
reco::PFMET | addInfo (edm::Handle< edm::View< reco::Candidate > > PFCandidates, CommonMETData met) |
PFSpecificAlgo () | |
void | runSignificance (metsig::SignAlgoResolutions &resolutions, edm::Handle< edm::View< reco::PFJet > > jets) |
Private Types | |
typedef math::XYZTLorentzVector | LorentzVector |
typedef math::XYZPoint | Point |
Private Member Functions | |
void | initializeSpecificPFMETData (SpecificPFMETData &specific) |
TMatrixD | mkSignifMatrix (edm::Handle< edm::View< reco::Candidate > > &PFCandidates) |
SpecificPFMETData | mkSpecificPFMETData (edm::Handle< edm::View< reco::Candidate > > &PFCandidates) |
Private Attributes | |
bool | doSignificance |
metsig::SignPFSpecificAlgo | pfsignalgo_ |
Description: Adds Particle Flow specific information to MET
Implementation: [Notes on implementation]
Definition at line 40 of file PFSpecificAlgo.h.
typedef math::XYZTLorentzVector PFSpecificAlgo::LorentzVector [private] |
Definition at line 49 of file PFSpecificAlgo.h.
typedef math::XYZPoint PFSpecificAlgo::Point [private] |
Definition at line 50 of file PFSpecificAlgo.h.
PFSpecificAlgo::PFSpecificAlgo | ( | ) | [inline] |
Definition at line 43 of file PFSpecificAlgo.h.
: doSignificance(false) { }
reco::PFMET PFSpecificAlgo::addInfo | ( | edm::Handle< edm::View< reco::Candidate > > | PFCandidates, |
CommonMETData | met | ||
) |
Definition at line 20 of file PFSpecificAlgo.cc.
References CommonMETData::met, CommonMETData::mex, CommonMETData::mey, p4, benchmark_cfg::PFCandidates, pfMET_cfi::pfMET, reco::MET::setSignificanceMatrix(), timingPdfMaker::specific, and CommonMETData::sumet.
Referenced by cms::METProducer::produce_PFMET().
{ SpecificPFMETData specific = mkSpecificPFMETData(PFCandidates); const LorentzVector p4(met.mex , met.mey, 0.0, met.met); const Point vtx(0.0,0.0,0.0); PFMET pfMET(specific, met.sumet, p4, vtx ); if(doSignificance) pfMET.setSignificanceMatrix(mkSignifMatrix(PFCandidates)); return pfMET; }
void PFSpecificAlgo::initializeSpecificPFMETData | ( | SpecificPFMETData & | specific | ) | [private] |
Definition at line 42 of file PFSpecificAlgo.cc.
References SpecificPFMETData::ChargedEMFraction, SpecificPFMETData::ChargedHadFraction, SpecificPFMETData::MuonFraction, SpecificPFMETData::NeutralEMFraction, SpecificPFMETData::NeutralHadFraction, SpecificPFMETData::Type6Fraction, and SpecificPFMETData::Type7Fraction.
{ specific.NeutralEMFraction = 0.0; specific.NeutralHadFraction = 0.0; specific.ChargedEMFraction = 0.0; specific.ChargedHadFraction = 0.0; specific.MuonFraction = 0.0; specific.Type6Fraction = 0.0; specific.Type7Fraction = 0.0; }
TMatrixD PFSpecificAlgo::mkSignifMatrix | ( | edm::Handle< edm::View< reco::Candidate > > & | PFCandidates | ) | [private] |
Definition at line 107 of file PFSpecificAlgo.cc.
References begin, end, and benchmark_cfg::PFCandidates.
{ pfsignalgo_.useOriginalPtrs(PFCandidates.id()); for(edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin(); iParticle != (PFCandidates.product())->end(); ++iParticle ) { const PFCandidate* pfCandidate = dynamic_cast<const PFCandidate*> (&(*iParticle)); if (!pfCandidate) continue; reco::CandidatePtr dau(PFCandidates, iParticle - PFCandidates->begin()); if(dau.isNull()) continue; if(!dau.isAvailable()) continue; reco::PFCandidatePtr pf(dau.id(), pfCandidate, dau.key()); pfsignalgo_.addPFCandidate(pf); } return pfsignalgo_.getSignifMatrix(); }
SpecificPFMETData PFSpecificAlgo::mkSpecificPFMETData | ( | edm::Handle< edm::View< reco::Candidate > > & | PFCandidates | ) | [private] |
Definition at line 54 of file PFSpecificAlgo.cc.
References begin, SpecificPFMETData::ChargedEMFraction, SpecificPFMETData::ChargedHadFraction, alignCSCRings::e, end, reco::LeafCandidate::energy(), SpecificPFMETData::MuonFraction, SpecificPFMETData::NeutralEMFraction, SpecificPFMETData::NeutralHadFraction, reco::PFCandidate::particleId(), benchmark_cfg::PFCandidates, funct::sin(), timingPdfMaker::specific, theta(), reco::LeafCandidate::theta(), SpecificPFMETData::Type6Fraction, and SpecificPFMETData::Type7Fraction.
{ if(!PFCandidates->size()) { SpecificPFMETData specific; initializeSpecificPFMETData(specific); return specific; } double NeutralEMEt = 0.0; double NeutralHadEt = 0.0; double ChargedEMEt = 0.0; double ChargedHadEt = 0.0; double MuonEt = 0.0; double type6Et = 0.0; double type7Et = 0.0; for( edm::View<reco::Candidate>::const_iterator iParticle = (PFCandidates.product())->begin(); iParticle != (PFCandidates.product())->end(); ++iParticle ) { const PFCandidate* pfCandidate = dynamic_cast<const PFCandidate*> (&(*iParticle)); if (!pfCandidate) continue; const double theta = pfCandidate->theta(); const double e = pfCandidate->energy(); const double et = e*sin(theta); if (pfCandidate->particleId() == 1) ChargedHadEt += et; if (pfCandidate->particleId() == 2) ChargedEMEt += et; if (pfCandidate->particleId() == 3) MuonEt += et; if (pfCandidate->particleId() == 4) NeutralEMEt += et; if (pfCandidate->particleId() == 5) NeutralHadEt += et; if (pfCandidate->particleId() == 6) type6Et += et; if (pfCandidate->particleId() == 7) type7Et += et; } const double Et_total = NeutralEMEt + NeutralHadEt + ChargedEMEt + ChargedHadEt + MuonEt + type6Et + type7Et; SpecificPFMETData specific; initializeSpecificPFMETData(specific); if (Et_total!=0.0) { specific.NeutralEMFraction = NeutralEMEt/Et_total; specific.NeutralHadFraction = NeutralHadEt/Et_total; specific.ChargedEMFraction = ChargedEMEt/Et_total; specific.ChargedHadFraction = ChargedHadEt/Et_total; specific.MuonFraction = MuonEt/Et_total; specific.Type6Fraction = type6Et/Et_total; specific.Type7Fraction = type7Et/Et_total; } return specific; }
void PFSpecificAlgo::runSignificance | ( | metsig::SignAlgoResolutions & | resolutions, |
edm::Handle< edm::View< reco::PFJet > > | jets | ||
) |
Definition at line 34 of file PFSpecificAlgo.cc.
References analyzePatCleaning_cfg::jets.
Referenced by cms::METProducer::produce_PFMET().
{ doSignificance = true; pfsignalgo_.setResolutions( &resolutions ); pfsignalgo_.addPFJets(jets); }
bool PFSpecificAlgo::doSignificance [private] |
Definition at line 55 of file PFSpecificAlgo.h.
Definition at line 56 of file PFSpecificAlgo.h.