![]() |
![]() |
Produces a ValueMap between JetCorrFactors and the to the originating reco jets. More...
#include <PhysicsTools/PatAlgos/interface/JetCorrFactorsProducer.h>
Public Types | |
typedef std::map < JetCorrFactors::Flavor, std::vector< std::string > > | FlavorCorrLevelMap |
map of correction levels to different flavors | |
typedef edm::ValueMap < pat::JetCorrFactors > | JetCorrFactorsMap |
value map for JetCorrFactors (to be written into the event) | |
Public Member Functions | |
JetCorrFactorsProducer (const edm::ParameterSet &cfg) | |
default constructor | |
virtual void | produce (edm::Event &event, const edm::EventSetup &setup) |
everything that needs to be done per event | |
~JetCorrFactorsProducer () | |
default destructor | |
Static Public Member Functions | |
static void | fillDescriptions (edm::ConfigurationDescriptions &descriptions) |
description of configuration file parameters | |
Private Member Functions | |
float | evaluate (edm::View< reco::Jet >::const_iterator &jet, boost::shared_ptr< FactorizedJetCorrector > &corrector, int level) |
evaluate jet correction factor up to a given level | |
std::vector< std::string > | expand (const std::vector< std::string > &levels, const JetCorrFactors::Flavor &flavor) |
bool | flavorDependent () const |
return true if the jec levels contain at least one flavor dependent correction level | |
int | numberOf (const edm::Handle< std::vector< reco::Vertex > > &primaryVertices) |
determines the number of valid primary vertices for the standard L1Offset correction of JetMET | |
std::vector < JetCorrectorParameters > | params (const JetCorrectorParametersCollection ¶meters, const std::vector< std::string > &levels) const |
return the jec parameters as input to the FactorizedJetCorrector for different flavors | |
std::string | payload () |
map jet algorithm to payload in DB | |
Private Attributes | |
bool | emf_ |
use electromagnetic fraction for jet energy corrections or not (will only have an effect for jets CaloJets) | |
std::string | label_ |
label of jec factors | |
FlavorCorrLevelMap | levels_ |
std::string | payload_ |
label of payload | |
edm::InputTag | primaryVertices_ |
label for L1Offset primaryVertex collection | |
edm::InputTag | rho_ |
label for L1FastJet energy density parameter rho | |
edm::InputTag | src_ |
input jet collection | |
std::string | type_ |
type of flavor dependent JEC factors (only 'J' and 'T' are allowed) | |
bool | useNPV_ |
use the NPV and rho with the JEC? (used for L1Offset/L1FastJet and L1FastJet, resp.) | |
bool | useRho_ |
Produces a ValueMap between JetCorrFactors and the to the originating reco jets.
The JetCorrFactorsProducer produces a set of correction factors, defined in the class pat::JetCorrFactors. This vector is linked to the originating reco jets through an edm::ValueMap. The initializing parameters of the module can be found in the recoLayer1/jetCorrFactors_cfi.py of the PatAlgos package. In the standard PAT workflow the module has to be run before the creation of the pat::Jet. The edm::ValueMap will then be embedded into the pat::Jet.
Jets corrected up to a given correction level can then be accessed via the pat::Jet member function correctedJet. For more details have a look into the class description of the pat::Jet.
ATTENTION: available options for flavor corrections are L5Flavor_gJ L7Parton_gJ gluon from dijets L5Flavor_qJ/_qT L7Parton_qJ/_qT quark from dijets/top L5Flavor_cJ/_cT L7Parton_cJ/_cT charm from dijets/top L5Flavor_bJ/_bT L7Parton_bJ/_bT beauty from dijets/top L7Parton_jJ/_tT mixture from dijets/top
where mixture refers to the flavor mixture as determined from the MC sample the flavor dependent corrections have been derived from. 'J' and 'T' stand for a typical dijet (ttbar) sample.
L1Offset corrections require the collection of _offlinePrimaryVertices_, which are supposed to be added as an additional optional parameter _primaryVertices_ in the jetCorrFactors_cfi.py file.
L1FastJet corrections, which are an alternative to the standard L1Offset correction as recommended by the JetMET PAG the energy density parameter _rho_ is supposed to be added as an additional optional parameter _rho_ in the jetCorrFactors_cfi.py file.
NOTE: the mixed mode (mc input mixture from dijets/ttbar) only exists for parton level corrections. jJ and tT are not covered in this implementation of the JetCorrFactorsProducer there are no gluon corrections available from the top sample on the L7Parton level.
Definition at line 59 of file JetCorrFactorsProducer.h.
typedef std::map<JetCorrFactors::Flavor, std::vector<std::string> > pat::JetCorrFactorsProducer::FlavorCorrLevelMap |
map of correction levels to different flavors
Definition at line 64 of file JetCorrFactorsProducer.h.
value map for JetCorrFactors (to be written into the event)
Definition at line 62 of file JetCorrFactorsProducer.h.
JetCorrFactorsProducer::JetCorrFactorsProducer | ( | const edm::ParameterSet & | cfg | ) | [explicit] |
default constructor
Definition at line 21 of file JetCorrFactorsProducer.cc.
References pat::JetCorrFactors::BOTTOM, pat::JetCorrFactors::CHARM, Exception, edm::ParameterSet::existsAs(), expand(), spr::find(), edm::ParameterSet::getParameter(), pat::JetCorrFactors::GLUON, levels_, argparse::message, NONE, primaryVertices_, rho_, pat::JetCorrFactors::UDS, useNPV_, and useRho_.
: emf_(cfg.getParameter<bool>( "emf" )), src_(cfg.getParameter<edm::InputTag>( "src" )), type_ (cfg.getParameter<std::string>("flavorType")), label_(cfg.getParameter<std::string>( "@module_label" )), payload_( cfg.getParameter<std::string>("payload") ), useNPV_(cfg.getParameter<bool>("useNPV")), useRho_(cfg.getParameter<bool>("useRho")) { std::vector<std::string> levels = cfg.getParameter<std::vector<std::string> >("levels"); // fill the std::map for levels_, which might be flavor dependent or not; // flavor dependency is determined from the fact whether the std::string // L5Flavor or L7Parton can be found in levels; if flavor dependent four // vectors of strings will be filled into the map corresponding to GLUON, // UDS, CHARM and BOTTOM (according to JetCorrFactors::Flavor), 'L5Flavor' // and 'L7Parton' will be expanded accordingly; if not levels_ is filled // with only one vector of strings according to NONE. This vector will be // equivalent to the original vector of strings. if(std::find(levels.begin(), levels.end(), "L5Flavor")!=levels.end() || std::find(levels.begin(), levels.end(), "L7Parton")!=levels.end()){ levels_[JetCorrFactors::GLUON ] = expand(levels, JetCorrFactors::GLUON ); levels_[JetCorrFactors::UDS ] = expand(levels, JetCorrFactors::UDS ); levels_[JetCorrFactors::CHARM ] = expand(levels, JetCorrFactors::CHARM ); levels_[JetCorrFactors::BOTTOM] = expand(levels, JetCorrFactors::BOTTOM); } else{ levels_[JetCorrFactors::NONE ] = levels; } // if the std::string L1Offset can be found in levels an additional para- // meter primaryVertices is needed, which should pass on the offline pri- // mary vertex collection. The size of this collection is needed for the // L1Offset correction. if(useNPV_){ if(cfg.existsAs<edm::InputTag>("primaryVertices")){ primaryVertices_=cfg.getParameter<edm::InputTag>("primaryVertices"); } else{ throw cms::Exception("No primaryVertices specified") << "The configured correction levels contain an L1Offset or L1FastJet correction, \n" << "which requires the number of offlinePrimaryVertices. Please specify this col- \n" << "lection as additional optional parameter primaryVertices in the jetCorrFactors\n" << "module. \n"; } } // if the std::string L1FastJet can be found in levels an additional // parameter rho is needed, which should pass on the energy density // parameter for the corresponding jet collection. if(useRho_){ if(std::find(levels.begin(), levels.end(), "L1FastJet")!=levels.end()){ if(cfg.existsAs<edm::InputTag>("rho")){ rho_=cfg.getParameter<edm::InputTag>("rho"); } else{ throw cms::Exception("No parameter rho specified") << "The configured correction levels contain a L1FastJet correction, which re- \n" << "quires the energy density parameter rho. Please specify this collection as \n" << "additional optional parameter rho in the jetCorrFactors module. \n"; } } else{ edm::LogWarning message( "Parameter rho not used" ); message << "Module is configured to use the parameter rho, but but rho is only used \n" << "for L1FastJet corrections at the moment. The configuration of levels \n" << "does not contain L1FastJet corrections though, so rho will not be used \n" << "throughout this module. \n"; } } produces<JetCorrFactorsMap>(); }
pat::JetCorrFactorsProducer::~JetCorrFactorsProducer | ( | ) | [inline] |
float JetCorrFactorsProducer::evaluate | ( | edm::View< reco::Jet >::const_iterator & | jet, |
boost::shared_ptr< FactorizedJetCorrector > & | corrector, | ||
int | level | ||
) | [private] |
evaluate jet correction factor up to a given level
Definition at line 129 of file JetCorrFactorsProducer.cc.
References emf_, and testEve_cfg::level.
Referenced by produce().
std::vector< std::string > JetCorrFactorsProducer::expand | ( | const std::vector< std::string > & | levels, |
const JetCorrFactors::Flavor & | flavor | ||
) | [private] |
return an expanded version of correction levels for different flavors; the result should be of type ['L2Relative', 'L3Absolute', 'L5FLavor_gJ', 'L7Parton_gJ']; L7Parton_gT will result in an empty string as this correction level is not available
Definition at line 92 of file JetCorrFactorsProducer.cc.
References pat::JetCorrFactors::BOTTOM, pat::JetCorrFactors::CHARM, pat::JetCorrFactors::GLUON, testEve_cfg::level, argparse::message, type_, and pat::JetCorrFactors::UDS.
Referenced by JetCorrFactorsProducer().
{ std::vector<std::string> expand; for(std::vector<std::string>::const_iterator level=levels.begin(); level!=levels.end(); ++level){ if((*level)=="L5Flavor" || (*level)=="L7Parton"){ if(flavor==JetCorrFactors::GLUON ){ if(*level=="L7Parton" && type_=="T"){ edm::LogWarning message( "L7Parton::GLUON not available" ); message << "Jet energy corrections requested for level: L7Parton and type: 'T'. \n" << "For this combination there is no GLUON correction available. The \n" << "correction for this flavor type will be taken from 'J'."; } expand.push_back(std::string(*level).append("_").append("g").append("J")); } if(flavor==JetCorrFactors::UDS ) expand.push_back(std::string(*level).append("_").append("q").append(type_)); if(flavor==JetCorrFactors::CHARM ) expand.push_back(std::string(*level).append("_").append("c").append(type_)); if(flavor==JetCorrFactors::BOTTOM) expand.push_back(std::string(*level).append("_").append("b").append(type_)); } else{ expand.push_back(*level); } } return expand; }
void JetCorrFactorsProducer::fillDescriptions | ( | edm::ConfigurationDescriptions & | descriptions | ) | [static] |
description of configuration file parameters
Reimplemented from edm::EDProducer.
Definition at line 255 of file JetCorrFactorsProducer.cc.
References edm::ParameterSetDescription::add(), and edm::ConfigurationDescriptions::add().
{ edm::ParameterSetDescription iDesc; iDesc.add<bool>("emf", false); iDesc.add<std::string>("flavorType", "J"); iDesc.add<edm::InputTag>("src", edm::InputTag("ak5CaloJets")); iDesc.add<std::string>("payload", "AK5Calo"); iDesc.add<bool>("useNPV", true); iDesc.add<edm::InputTag>("primaryVertices", edm::InputTag("offlinePrimaryVertices")); iDesc.add<bool>("useRho", false); iDesc.add<edm::InputTag>("rho", edm::InputTag("kt6PFJets", "rho")); std::vector<std::string> levels; levels.push_back(std::string("L1Offset" )); levels.push_back(std::string("L2Relative")); levels.push_back(std::string("L3Absolute")); levels.push_back(std::string("L5Flavor" )); levels.push_back(std::string("L7Parton" )); iDesc.add<std::vector<std::string> >("levels", levels); descriptions.add("JetCorrFactorsProducer", iDesc); }
bool pat::JetCorrFactorsProducer::flavorDependent | ( | ) | const [inline, private] |
return true if the jec levels contain at least one flavor dependent correction level
Definition at line 78 of file JetCorrFactorsProducer.h.
References levels_.
Referenced by produce().
{ return (levels_.size()>1); };
int pat::JetCorrFactorsProducer::numberOf | ( | const edm::Handle< std::vector< reco::Vertex > > & | primaryVertices | ) | [inline, private] |
determines the number of valid primary vertices for the standard L1Offset correction of JetMET
Definition at line 124 of file JetCorrFactorsProducer.h.
Referenced by produce().
{ int npv=0; for(std::vector<reco::Vertex>::const_iterator pv=primaryVertices->begin(); pv!=primaryVertices->end(); ++pv){ if(pv->ndof()>=4) ++npv; } return npv; }
std::vector< JetCorrectorParameters > JetCorrFactorsProducer::params | ( | const JetCorrectorParametersCollection & | parameters, |
const std::vector< std::string > & | levels | ||
) | const [private] |
return the jec parameters as input to the FactorizedJetCorrector for different flavors
Definition at line 118 of file JetCorrFactorsProducer.cc.
References testEve_cfg::level, and JetCorrectorParametersCollection::push_back().
Referenced by produce().
std::string JetCorrFactorsProducer::payload | ( | ) | [private] |
map jet algorithm to payload in DB
Definition at line 139 of file JetCorrFactorsProducer.cc.
References payload_.
Referenced by produce().
{ return payload_; }
void JetCorrFactorsProducer::produce | ( | edm::Event & | event, |
const edm::EventSetup & | setup | ||
) | [virtual] |
everything that needs to be done per event
Implements edm::EDProducer.
Definition at line 145 of file JetCorrFactorsProducer.cc.
References evaluate(), Exception, edm::helper::Filler< Map >::fill(), flavorDependent(), edm::EventSetup::get(), edm::Event::getByLabel(), edm::helper::Filler< Map >::insert(), metsig::jet, analyzePatCleaning_cfg::jets, edm::InputTag::label(), label_, levels_, numberOf(), Parameters::parameters, params(), payload(), primaryVertices_, rho, rho_, and src_.
{ // get jet collection from the event edm::Handle<edm::View<reco::Jet> > jets; event.getByLabel(src_, jets); // get primary vertices for L1Offset correction level if needed edm::Handle<std::vector<reco::Vertex> > primaryVertices; if(!primaryVertices_.label().empty()) event.getByLabel(primaryVertices_, primaryVertices); // get parameter rho for L1FastJet correction level if needed edm::Handle<double> rho; if(!rho_.label().empty()) event.getByLabel(rho_, rho); // retreive parameters from the DB edm::ESHandle<JetCorrectorParametersCollection> parameters; setup.get<JetCorrectionsRecord>().get(payload(), parameters); // initialize jet correctors std::map<JetCorrFactors::Flavor, boost::shared_ptr<FactorizedJetCorrector> > corrector; for(FlavorCorrLevelMap::const_iterator flavor=levels_.begin(); flavor!=levels_.end(); ++flavor){ corrector[flavor->first] = boost::shared_ptr<FactorizedJetCorrector>( new FactorizedJetCorrector(params(*parameters, flavor->second)) ); } // fill the jetCorrFactors std::vector<JetCorrFactors> jcfs; for(edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet!=jets->end(); ++jet){ // the JetCorrFactors::CorrectionFactor is a std::pair<std::string, std::vector<float> > // the string corresponds to the label of the correction level, the vector contains four // floats if flavor dependent and one float else. Per construction jet energy corrections // will be flavor independent up to the first flavor dependent correction and flavor de- // pendent afterwards. The first correction level is predefined with label 'Uncorrected'. // Per definition it is flavor independent. The correction factor is 1. std::vector<JetCorrFactors::CorrectionFactor> jec; jec.push_back(std::make_pair<std::string, std::vector<float> >(std::string("Uncorrected"), std::vector<float>(1, 1))); // pick the first element in the map (which could be the only one) and loop all jec // levels listed for that element. If this is not the only element all jec levels, which // are flavor independent will give the same correction factors until the first flavor // dependent correction level is reached. So the first element is still a good choice. FlavorCorrLevelMap::const_iterator corrLevel=levels_.begin(); if(corrLevel==levels_.end()){ throw cms::Exception("No JECFactors") << "You request to create a jetCorrFactors object with no JEC Levels indicated. \n" << "This makes no sense, either you should correct this or drop the module from \n" << "the sequence."; } for(unsigned int idx=0; idx<corrLevel->second.size(); ++idx){ bool flavorDependent=false; std::vector<float> factors; if(flavorDependent || corrLevel->second[idx].find("L5Flavor")!=std::string::npos || corrLevel->second[idx].find("L7Parton")!=std::string::npos){ flavorDependent=true; // after the first encounter all subsequent correction levels are flavor dependent for(FlavorCorrLevelMap::const_iterator flavor=corrLevel; flavor!=levels_.end(); ++flavor){ if(!primaryVertices_.label().empty()){ // if primaryVertices_ has a value the number of primary vertices needs to be // specified corrector.find(flavor->first)->second->setNPV(numberOf(primaryVertices)); } if(!rho_.label().empty()){ // if rho_ has a value the energy density parameter rho and the jet area need // to be specified corrector.find(flavor->first)->second->setRho(*rho); corrector.find(flavor->first)->second->setJetA(jet->jetArea()); } factors.push_back(evaluate(jet, corrector.find(flavor->first)->second, idx)); } } else{ if(!primaryVertices_.label().empty()){ // if primaryVertices_ has a value the number of primary vertices needs to be // specified corrector.find(corrLevel->first)->second->setNPV(numberOf(primaryVertices)); } if(!rho_.label().empty()){ // if rho_ has a value the energy density parameter rho and the jet area need // to be specified corrector.find(corrLevel->first)->second->setRho(*rho); corrector.find(corrLevel->first)->second->setJetA(jet->jetArea()); } factors.push_back(evaluate(jet, corrector.find(corrLevel->first)->second, idx)); } // push back the set of JetCorrFactors: the first entry corresponds to the label // of the correction level, which is taken from the first element in levels_. For // L5Flavor and L7Parton the part including the first '_' indicating the flavor // of the first element in levels_ is chopped of from the label to avoid confusion // of the correction levels. The second parameter corresponds to the set of jec // factors, which might be flavor dependent or not. In the default configuration // the CorrectionFactor will look like this: 'Uncorrected': 1 ; 'L2Relative': x ; // 'L3Absolute': x ; 'L5Flavor': v, x, y, z ; 'L7Parton': v, x, y, z jec.push_back(std::make_pair((corrLevel->second[idx]).substr(0, (corrLevel->second[idx]).find("_")), factors)); } // create the actual object with the scale factors we want the valuemap to refer to // label_ corresponds to the label of the module instance JetCorrFactors corrFactors(label_, jec); jcfs.push_back(corrFactors); } // build the value map std::auto_ptr<JetCorrFactorsMap> jetCorrsMap(new JetCorrFactorsMap()); JetCorrFactorsMap::Filler filler(*jetCorrsMap); // jets and jetCorrs have their indices aligned by construction filler.insert(jets, jcfs.begin(), jcfs.end()); filler.fill(); // do the actual filling // put our produced stuff in the event event.put(jetCorrsMap); }
bool pat::JetCorrFactorsProducer::emf_ [private] |
use electromagnetic fraction for jet energy corrections or not (will only have an effect for jets CaloJets)
Definition at line 94 of file JetCorrFactorsProducer.h.
Referenced by evaluate().
std::string pat::JetCorrFactorsProducer::label_ [private] |
label of jec factors
Definition at line 100 of file JetCorrFactorsProducer.h.
Referenced by produce().
jec levels for different flavors. In the default configuration this map would look like this: GLUON : 'L2Relative', 'L3Absolute', 'L5FLavor_jg', L7Parton_jg' UDS : 'L2Relative', 'L3Absolute', 'L5FLavor_jq', L7Parton_jq' CHARM : 'L2Relative', 'L3Absolute', 'L5FLavor_jc', L7Parton_jc' BOTTOM : 'L2Relative', 'L3Absolute', 'L5FLavor_jb', L7Parton_jb' or just like this: NONE : 'L2Relative', 'L3Absolute', 'L2L3Residual' per definition the vectors for all elements in this map should have the same size
Definition at line 120 of file JetCorrFactorsProducer.h.
Referenced by flavorDependent(), JetCorrFactorsProducer(), and produce().
std::string pat::JetCorrFactorsProducer::payload_ [private] |
label for L1Offset primaryVertex collection
Definition at line 104 of file JetCorrFactorsProducer.h.
Referenced by JetCorrFactorsProducer(), and produce().
label for L1FastJet energy density parameter rho
Definition at line 106 of file JetCorrFactorsProducer.h.
Referenced by JetCorrFactorsProducer(), and produce().
input jet collection
Definition at line 96 of file JetCorrFactorsProducer.h.
Referenced by produce().
std::string pat::JetCorrFactorsProducer::type_ [private] |
type of flavor dependent JEC factors (only 'J' and 'T' are allowed)
Definition at line 98 of file JetCorrFactorsProducer.h.
Referenced by expand().
bool pat::JetCorrFactorsProducer::useNPV_ [private] |
use the NPV and rho with the JEC? (used for L1Offset/L1FastJet and L1FastJet, resp.)
Definition at line 108 of file JetCorrFactorsProducer.h.
Referenced by JetCorrFactorsProducer().
bool pat::JetCorrFactorsProducer::useRho_ [private] |
Definition at line 109 of file JetCorrFactorsProducer.h.
Referenced by JetCorrFactorsProducer().