CMS 3D CMS Logo

Classes | Public Member Functions | Private Member Functions | Private Attributes

PFJetMETcorrInputProducerT< T, Textractor > Class Template Reference

#include <PFJetMETcorrInputProducerT.h>

Inheritance diagram for PFJetMETcorrInputProducerT< T, Textractor >:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Classes

struct  type2BinningEntryType

Public Member Functions

 PFJetMETcorrInputProducerT (const edm::ParameterSet &cfg)
 ~PFJetMETcorrInputProducerT ()

Private Member Functions

void produce (edm::Event &evt, const edm::EventSetup &es)

Private Attributes

double jetCorrEtaMax_
Textractor jetCorrExtractor_
std::string jetCorrLabel_
std::string moduleLabel_
std::string offsetCorrLabel_
bool skipEM_
double skipEMfractionThreshold_
bool skipMuons_
StringCutObjectSelector
< reco::Muon > * 
skipMuonSelection_
edm::InputTag src_
double type1JetPtThreshold_
std::vector
< type2BinningEntryType * > 
type2Binning_

Detailed Description

template<typename T, typename Textractor>
class PFJetMETcorrInputProducerT< T, Textractor >

Produce Type 1 + 2 MET corrections corresponding to differences between raw PFJets and PFJets with jet energy corrections (JECs) applied

NOTE: class is templated to that it works with reco::PFJets as well as with pat::Jets of PF-type as input

Authors:
Michael Schmitt, Richard Cavanaugh, The University of Florida Florent Lacroix, University of Illinois at Chicago Christian Veelken, LLR
Version:
Revision:
1.8
Id:
PFJetMETcorrInputProducerT.h,v 1.8 2012/04/19 17:56:29 veelken Exp

Definition at line 66 of file PFJetMETcorrInputProducerT.h.


Constructor & Destructor Documentation

template<typename T , typename Textractor >
PFJetMETcorrInputProducerT< T, Textractor >::PFJetMETcorrInputProducerT ( const edm::ParameterSet cfg) [inline, explicit]

Definition at line 70 of file PFJetMETcorrInputProducerT.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, PFJetMETcorrInputProducerT< T, Textractor >::skipEM_, PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_, PFJetMETcorrInputProducerT< T, Textractor >::src_, PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_, and PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_.

    : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
      skipMuonSelection_(0)
  {
    src_ = cfg.getParameter<edm::InputTag>("src");
    
    offsetCorrLabel_ = ( cfg.exists("offsetCorrLabel") ) ?
      cfg.getParameter<std::string>("offsetCorrLabel") : "";
    jetCorrLabel_ = cfg.getParameter<std::string>("jetCorrLabel");
    
    jetCorrEtaMax_ = ( cfg.exists("jetCorrEtaMax") ) ?
      cfg.getParameter<double>("jetCorrEtaMax") : 9.9;
    
    type1JetPtThreshold_ = cfg.getParameter<double>("type1JetPtThreshold");
    
    skipEM_ = cfg.getParameter<bool>("skipEM");
    if ( skipEM_ ) {
      skipEMfractionThreshold_ = cfg.getParameter<double>("skipEMfractionThreshold");
    }
    
    skipMuons_ = cfg.getParameter<bool>("skipMuons");
    if ( skipMuons_ ) {
      std::string skipMuonSelection_string = cfg.getParameter<std::string>("skipMuonSelection");
      skipMuonSelection_ = new StringCutObjectSelector<reco::Muon>(skipMuonSelection_string);
    }

    if ( cfg.exists("type2Binning") ) {
      typedef std::vector<edm::ParameterSet> vParameterSet;
      vParameterSet cfgType2Binning = cfg.getParameter<vParameterSet>("type2Binning");
      for ( vParameterSet::const_iterator cfgType2BinningEntry = cfgType2Binning.begin();
            cfgType2BinningEntry != cfgType2Binning.end(); ++cfgType2BinningEntry ) {
        type2Binning_.push_back(new type2BinningEntryType(*cfgType2BinningEntry));
      }
    } else {
      type2Binning_.push_back(new type2BinningEntryType());
    }

    produces<CorrMETData>("type1");
    for ( typename std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
          type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) {   
      produces<CorrMETData>((*type2BinningEntry)->getInstanceLabel_full("type2"));
      produces<CorrMETData>((*type2BinningEntry)->getInstanceLabel_full("offset"));
    }
  }
template<typename T , typename Textractor >
PFJetMETcorrInputProducerT< T, Textractor >::~PFJetMETcorrInputProducerT ( ) [inline]

Definition at line 114 of file PFJetMETcorrInputProducerT.h.

References PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_, and PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_.

  {
    delete skipMuonSelection_;

    for ( typename std::vector<type2BinningEntryType*>::const_iterator it = type2Binning_.begin();
          it != type2Binning_.end(); ++it ) {
      delete (*it);
    }
  }

Member Function Documentation

template<typename T , typename Textractor >
void PFJetMETcorrInputProducerT< T, Textractor >::produce ( edm::Event evt,
const edm::EventSetup es 
) [inline, private, virtual]

Implements edm::EDProducer.

Definition at line 126 of file PFJetMETcorrInputProducerT.h.

References edm::Event::getByLabel(), PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_, PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_, fwrapper::jets, PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_, edm::Event::put(), PFJetMETcorrInputProducerT< T, Textractor >::skipEM_, PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_, PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_, PFJetMETcorrInputProducerT< T, Textractor >::src_, PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_, and PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_.

  {
    std::auto_ptr<CorrMETData> type1Correction(new CorrMETData());
    for ( typename std::vector<type2BinningEntryType*>::iterator type2BinningEntry = type2Binning_.begin();
          type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) {
      (*type2BinningEntry)->binUnclEnergySum_   = CorrMETData();
      (*type2BinningEntry)->binOffsetEnergySum_ = CorrMETData();
    }

    typedef std::vector<T> JetCollection;
    edm::Handle<JetCollection> jets;
    evt.getByLabel(src_, jets);

    int numJets = jets->size();
    for ( int jetIndex = 0; jetIndex < numJets; ++jetIndex ) {
      const T& rawJet = jets->at(jetIndex);
      
      static PFJetMETcorrInputProducer_namespace::InputTypeCheckerT<T, Textractor> checkInputType;
      checkInputType(rawJet);
      
      double emEnergyFraction = rawJet.chargedEmEnergyFraction() + rawJet.neutralEmEnergyFraction();
      if ( skipEM_ && emEnergyFraction > skipEMfractionThreshold_ ) continue;
      
      static PFJetMETcorrInputProducer_namespace::RawJetExtractorT<T> rawJetExtractor;
      reco::Candidate::LorentzVector rawJetP4 = rawJetExtractor(rawJet);
      if ( skipMuons_ ) {
        std::vector<reco::PFCandidatePtr> cands = rawJet.getPFConstituents();
        for ( std::vector<reco::PFCandidatePtr>::const_iterator cand = cands.begin();
              cand != cands.end(); ++cand ) {
          if ( (*cand)->muonRef().isNonnull() && (*skipMuonSelection_)(*(*cand)->muonRef()) ) {
            reco::Candidate::LorentzVector muonP4 = (*cand)->p4();
            rawJetP4 -= muonP4;
          }
        }
      }

      reco::Candidate::LorentzVector corrJetP4 = jetCorrExtractor_(rawJet, jetCorrLabel_, &evt, &es, jetCorrEtaMax_, &rawJetP4);

      if ( corrJetP4.pt() > type1JetPtThreshold_ ) {
        
        reco::Candidate::LorentzVector rawJetP4offsetCorr = rawJetP4;
        if ( offsetCorrLabel_ != "" ) {
          rawJetP4offsetCorr = jetCorrExtractor_(rawJet, offsetCorrLabel_, &evt, &es, jetCorrEtaMax_, &rawJetP4);
          
          for ( typename std::vector<type2BinningEntryType*>::iterator type2BinningEntry = type2Binning_.begin();
                type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) {
            if ( !(*type2BinningEntry)->binSelection_ || (*(*type2BinningEntry)->binSelection_)(corrJetP4) ) {
              (*type2BinningEntry)->binOffsetEnergySum_.mex   += (rawJetP4.px() - rawJetP4offsetCorr.px());
              (*type2BinningEntry)->binOffsetEnergySum_.mey   += (rawJetP4.py() - rawJetP4offsetCorr.py());
              (*type2BinningEntry)->binOffsetEnergySum_.sumet += (rawJetP4.Et() - rawJetP4offsetCorr.Et());
            }
          }
        }

//--- MET balances momentum of reconstructed particles,
//    hence correction to jets and corresponding Type 1 MET correction are of opposite sign
        type1Correction->mex   -= (corrJetP4.px() - rawJetP4offsetCorr.px());
        type1Correction->mey   -= (corrJetP4.py() - rawJetP4offsetCorr.py());
        type1Correction->sumet += (corrJetP4.Et() - rawJetP4offsetCorr.Et());
      } else {
        for ( typename std::vector<type2BinningEntryType*>::iterator type2BinningEntry = type2Binning_.begin();
              type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) {
          if ( !(*type2BinningEntry)->binSelection_ || (*(*type2BinningEntry)->binSelection_)(corrJetP4) ) {
            (*type2BinningEntry)->binUnclEnergySum_.mex   += rawJetP4.px();
            (*type2BinningEntry)->binUnclEnergySum_.mey   += rawJetP4.py();
            (*type2BinningEntry)->binUnclEnergySum_.sumet += rawJetP4.Et();
          }
        }
      }
    }

//--- add 
//     o Type 1 MET correction                (difference corrected-uncorrected jet energy for jets of (corrected) Pt > 10 GeV)
//     o momentum sum of "unclustered energy" (jets of (corrected) Pt < 10 GeV)
//     o momentum sum of "offset energy"      (sum of energy attributed to pile-up/underlying event)
//    to the event
    evt.put(type1Correction, "type1");
    for ( typename std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin();
          type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) {
      evt.put(std::auto_ptr<CorrMETData>(new CorrMETData((*type2BinningEntry)->binUnclEnergySum_)), (*type2BinningEntry)->getInstanceLabel_full("type2"));
      evt.put(std::auto_ptr<CorrMETData>(new CorrMETData((*type2BinningEntry)->binOffsetEnergySum_)), (*type2BinningEntry)->getInstanceLabel_full("offset"));
    }
  }

Member Data Documentation

template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::jetCorrEtaMax_ [private]
template<typename T , typename Textractor >
Textractor PFJetMETcorrInputProducerT< T, Textractor >::jetCorrExtractor_ [private]
template<typename T , typename Textractor >
std::string PFJetMETcorrInputProducerT< T, Textractor >::jetCorrLabel_ [private]
template<typename T , typename Textractor >
std::string PFJetMETcorrInputProducerT< T, Textractor >::moduleLabel_ [private]

Definition at line 210 of file PFJetMETcorrInputProducerT.h.

template<typename T , typename Textractor >
std::string PFJetMETcorrInputProducerT< T, Textractor >::offsetCorrLabel_ [private]
template<typename T , typename Textractor >
bool PFJetMETcorrInputProducerT< T, Textractor >::skipEM_ [private]
template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::skipEMfractionThreshold_ [private]
template<typename T , typename Textractor >
bool PFJetMETcorrInputProducerT< T, Textractor >::skipMuons_ [private]
template<typename T , typename Textractor >
StringCutObjectSelector<reco::Muon>* PFJetMETcorrInputProducerT< T, Textractor >::skipMuonSelection_ [private]
template<typename T , typename Textractor >
edm::InputTag PFJetMETcorrInputProducerT< T, Textractor >::src_ [private]
template<typename T , typename Textractor >
double PFJetMETcorrInputProducerT< T, Textractor >::type1JetPtThreshold_ [private]
template<typename T , typename Textractor >
std::vector<type2BinningEntryType*> PFJetMETcorrInputProducerT< T, Textractor >::type2Binning_ [private]