Go to the documentation of this file.00001
00002
00003 #include "CommonTools/ParticleFlow/plugins/Type1PFMET.h"
00004
00005 #include "DataFormats/METReco/interface/MET.h"
00006 #include "DataFormats/METReco/interface/CaloMET.h"
00007 #include "DataFormats/JetReco/interface/PFJet.h"
00008
00009
00010
00011 using namespace reco;
00012
00013
00014 Type1PFMET::Type1PFMET( const edm::ParameterSet& iConfig )
00015 {
00016 inputUncorMetLabel = iConfig.getParameter<std::string>("inputUncorMetLabel");
00017 inputUncorJetsTag = iConfig.getParameter<edm::InputTag>("inputUncorJetsTag");
00018 correctorLabel = iConfig.getParameter<std::string>("corrector");
00019 jetPTthreshold = iConfig.getParameter<double>("jetPTthreshold");
00020 jetEMfracLimit = iConfig.getParameter<double>("jetEMfracLimit");
00021 jetMufracLimit = iConfig.getParameter<double>("jetMufracLimit");
00022 produces<METCollection>();
00023 }
00024 Type1PFMET::Type1PFMET() {}
00025
00026
00027 Type1PFMET::~Type1PFMET() {}
00028
00029
00030 void Type1PFMET::produce( edm::Event& iEvent, const edm::EventSetup& iSetup )
00031 {
00032 using namespace edm;
00033 Handle<PFJetCollection> inputUncorJets;
00034 iEvent.getByLabel( inputUncorJetsTag, inputUncorJets );
00035 const JetCorrector* corrector = JetCorrector::getJetCorrector (correctorLabel, iSetup);
00036 Handle<METCollection> inputUncorMet;
00037 iEvent.getByLabel( inputUncorMetLabel, inputUncorMet );
00038 std::auto_ptr<METCollection> output( new METCollection() );
00039 run( *(inputUncorMet.product()), *corrector, *(inputUncorJets.product()),
00040 jetPTthreshold, jetEMfracLimit, jetMufracLimit,
00041 &*output );
00042 iEvent.put( output );
00043 }
00044
00045 void Type1PFMET::run(const METCollection& uncorMET,
00046 const JetCorrector& corrector,
00047 const PFJetCollection& uncorJet,
00048 double jetPTthreshold,
00049 double jetEMfracLimit,
00050 double jetMufracLimit,
00051 METCollection* corMET)
00052 {
00053 if (!corMET) {
00054 std::cerr << "Type1METAlgo_run-> undefined output MET collection. Stop. " << std::endl;
00055 return;
00056 }
00057
00058 double DeltaPx = 0.0;
00059 double DeltaPy = 0.0;
00060 double DeltaSumET = 0.0;
00061
00062
00063
00064 for( PFJetCollection::const_iterator jet = uncorJet.begin(); jet != uncorJet.end(); ++jet) {
00065 if( jet->pt() > jetPTthreshold ) {
00066 double emEFrac =
00067 jet->chargedEmEnergyFraction() + jet->neutralEmEnergyFraction();
00068 double muEFrac = jet->chargedMuEnergyFraction();
00069 if( emEFrac < jetEMfracLimit
00070 && muEFrac < jetMufracLimit ) {
00071 double corr = corrector.correction (*jet) - 1.;
00072 DeltaPx += jet->px() * corr;
00073 DeltaPy += jet->py() * corr;
00074 DeltaSumET += jet->et() * corr;
00075 }
00076 }
00077 }
00078
00079 CorrMETData delta;
00080 delta.mex = - DeltaPx;
00081 delta.mey = - DeltaPy;
00082 delta.sumet = DeltaSumET;
00083
00084 const MET* u = &(uncorMET.front());
00085 double corrMetPx = u->px()+delta.mex;
00086 double corrMetPy = u->py()+delta.mey;
00087 MET::LorentzVector correctedMET4vector( corrMetPx, corrMetPy, 0.,
00088 sqrt (corrMetPx*corrMetPx + corrMetPy*corrMetPy)
00089 );
00090
00091 std::vector<CorrMETData> corrections = u->mEtCorr();
00092 corrections.push_back( delta );
00093
00094 MET result = MET(u->sumEt()+delta.sumet,
00095 corrections,
00096 correctedMET4vector,
00097 u->vertex() );
00098 corMET->push_back(result);
00099
00100 return;
00101 }
00102
00103
00104
00105