#include <METCorrectionAlgorithm.h>
Classes | |
struct | type2BinningEntryType |
Public Member Functions | |
CorrMETData | compMETCorrection (edm::Event &, const edm::EventSetup &) |
METCorrectionAlgorithm (const edm::ParameterSet &) | |
~METCorrectionAlgorithm () | |
Private Types | |
typedef std::vector < edm::InputTag > | vInputTag |
typedef std::vector< std::string > | vstring |
Private Attributes | |
bool | applyType0Corrections_ |
bool | applyType1Corrections_ |
bool | applyType2Corrections_ |
vInputTag | srcCHSSums_ |
vInputTag | srcType1Corrections_ |
double | type0Cuncl_ |
double | type0Rsoft_ |
std::vector < type2BinningEntryType * > | type2Binning_ |
Algorithm for o propagating jet energy corrections to MET (Type 1 MET corrections) o calibrating momentum of particles not within jets ("unclustered energy") and propagating those corrections to MET (Type 2 MET corrections)
Definition at line 33 of file METCorrectionAlgorithm.h.
typedef std::vector<edm::InputTag> METCorrectionAlgorithm::vInputTag [private] |
Definition at line 44 of file METCorrectionAlgorithm.h.
typedef std::vector<std::string> METCorrectionAlgorithm::vstring [private] |
Definition at line 55 of file METCorrectionAlgorithm.h.
METCorrectionAlgorithm::METCorrectionAlgorithm | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
Definition at line 9 of file METCorrectionAlgorithm.cc.
References applyType0Corrections_, applyType1Corrections_, applyType2Corrections_, AlCaHLTBitMon_QueryRunRegistry::data, Exception, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), srcCHSSums_, srcType1Corrections_, type0Cuncl_, type0Rsoft_, and type2Binning_.
{ applyType1Corrections_ = cfg.getParameter<bool>("applyType1Corrections"); if ( applyType1Corrections_ ) { srcType1Corrections_ = cfg.getParameter<vInputTag>("srcType1Corrections"); } applyType2Corrections_ = cfg.getParameter<bool>("applyType2Corrections"); if ( applyType2Corrections_ ) { vInputTag srcUnclEnergySums = cfg.getParameter<vInputTag>("srcUnclEnergySums"); 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, srcUnclEnergySums)); } } else { std::string type2CorrFormula = cfg.getParameter<std::string>("type2CorrFormula").data(); edm::ParameterSet type2CorrParameter = cfg.getParameter<edm::ParameterSet>("type2CorrParameter"); type2Binning_.push_back(new type2BinningEntryType(type2CorrFormula, type2CorrParameter, srcUnclEnergySums)); } } applyType0Corrections_ = cfg.exists("applyType0Corrections") ? cfg.getParameter<bool>("applyType0Corrections") : false; if ( applyType0Corrections_ ) { srcCHSSums_ = cfg.getParameter<vInputTag>("srcCHSSums"); type0Rsoft_ = cfg.getParameter<double>("type0Rsoft"); type0Cuncl_ = 1.0; if (applyType2Corrections_) { if (cfg.exists("type2Binning")) throw cms::Exception("Invalid Arg") << "Currently, applyType0Corrections and type2Binning cannot be used together!"; std::string type2CorrFormula = cfg.getParameter<std::string>("type2CorrFormula").data(); if (!(type2CorrFormula == "A")) throw cms::Exception("Invalid Arg") << "type2CorrFormula must be \"A\" if applyType0Corrections!"; edm::ParameterSet type2CorrParameter = cfg.getParameter<edm::ParameterSet>("type2CorrParameter"); type0Cuncl_ = type2CorrParameter.getParameter<double>("A"); } } }
METCorrectionAlgorithm::~METCorrectionAlgorithm | ( | ) |
Definition at line 49 of file METCorrectionAlgorithm.cc.
References type2Binning_.
{ for ( std::vector<type2BinningEntryType*>::const_iterator it = type2Binning_.begin(); it != type2Binning_.end(); ++it ) { delete (*it); } }
CorrMETData METCorrectionAlgorithm::compMETCorrection | ( | edm::Event & | evt, |
const edm::EventSetup & | es | ||
) |
Definition at line 57 of file METCorrectionAlgorithm.cc.
References applyType0Corrections_, applyType1Corrections_, applyType2Corrections_, edm::Event::getByLabel(), CorrMETData::mex, CorrMETData::mey, mathSSE::sqrt(), srcCHSSums_, srcType1Corrections_, CorrMETData::sumet, type0Cuncl_, type0Rsoft_, and type2Binning_.
Referenced by CorrectedMETProducerT< T >::produce().
{ CorrMETData metCorr; metCorr.mex = 0.; metCorr.mey = 0.; metCorr.sumet = 0.; if ( applyType0Corrections_ ) { //--- sum all Type 0 MET correction terms for ( vInputTag::const_iterator srcCHSSum = srcCHSSums_.begin(); srcCHSSum != srcCHSSums_.end(); ++srcCHSSum ) { edm::Handle<CorrMETData> chsSum; evt.getByLabel(*srcCHSSum, chsSum); metCorr.mex += type0Cuncl_*(1 - type0Rsoft_)*chsSum->mex; metCorr.mey += type0Cuncl_*(1 - type0Rsoft_)*chsSum->mey; metCorr.sumet += type0Cuncl_*(1 - type0Rsoft_)*chsSum->sumet; } } if ( applyType1Corrections_ ) { //--- sum all Type 1 MET correction terms for ( vInputTag::const_iterator srcType1Correction = srcType1Corrections_.begin(); srcType1Correction != srcType1Corrections_.end(); ++srcType1Correction ) { edm::Handle<CorrMETData> type1Correction; evt.getByLabel(*srcType1Correction, type1Correction); metCorr.mex += type1Correction->mex; metCorr.mey += type1Correction->mey; metCorr.sumet += type1Correction->sumet; } } if ( applyType2Corrections_ ) { //--- compute momentum sum of all "unclustered energy" in the event // // NOTE: calibration factors/formulas for Type 2 MET correction may depend on eta // (like the jet energy correction factors do) // for ( std::vector<type2BinningEntryType*>::const_iterator type2BinningEntry = type2Binning_.begin(); type2BinningEntry != type2Binning_.end(); ++type2BinningEntry ) { CorrMETData unclEnergySum; for ( vInputTag::const_iterator srcUnclEnergySum = (*type2BinningEntry)->srcUnclEnergySums_.begin(); srcUnclEnergySum != (*type2BinningEntry)->srcUnclEnergySums_.end(); ++srcUnclEnergySum ) { edm::Handle<CorrMETData> unclEnergySummand; evt.getByLabel(*srcUnclEnergySum, unclEnergySummand); unclEnergySum.mex += unclEnergySummand->mex; unclEnergySum.mey += unclEnergySummand->mey; unclEnergySum.sumet += unclEnergySummand->sumet; } //--- calibrate "unclustered energy" double unclEnergySumPt = sqrt(unclEnergySum.mex*unclEnergySum.mex + unclEnergySum.mey*unclEnergySum.mey); double unclEnergyScaleFactor = (*type2BinningEntry)->binCorrFormula_->Eval(unclEnergySumPt); //--- MET balances momentum of reconstructed particles, // hence correction to "unclustered energy" and corresponding Type 2 MET correction are of opposite sign metCorr.mex -= (unclEnergyScaleFactor - 1.)*unclEnergySum.mex; metCorr.mey -= (unclEnergyScaleFactor - 1.)*unclEnergySum.mey; metCorr.sumet += (unclEnergyScaleFactor - 1.)*unclEnergySum.sumet; } } return metCorr; }
bool METCorrectionAlgorithm::applyType0Corrections_ [private] |
Definition at line 48 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
bool METCorrectionAlgorithm::applyType1Corrections_ [private] |
Definition at line 49 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
bool METCorrectionAlgorithm::applyType2Corrections_ [private] |
Definition at line 50 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
vInputTag METCorrectionAlgorithm::srcCHSSums_ [private] |
Definition at line 45 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
Definition at line 46 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
double METCorrectionAlgorithm::type0Cuncl_ [private] |
Definition at line 53 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
double METCorrectionAlgorithm::type0Rsoft_ [private] |
Definition at line 52 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), and METCorrectionAlgorithm().
std::vector<type2BinningEntryType*> METCorrectionAlgorithm::type2Binning_ [private] |
Definition at line 123 of file METCorrectionAlgorithm.h.
Referenced by compMETCorrection(), METCorrectionAlgorithm(), and ~METCorrectionAlgorithm().