File indexing completed on 2024-09-07 04:37:21
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036 #include "CondFormats/JetMETObjects/interface/FactorizedJetCorrector.h"
0037 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
0038 #include "DataFormats/Common/interface/ValueMap.h"
0039 #include "DataFormats/Common/interface/View.h"
0040 #include "DataFormats/JetReco/interface/CaloJet.h"
0041 #include "DataFormats/JetReco/interface/JPTJet.h"
0042 #include "DataFormats/JetReco/interface/Jet.h"
0043 #include "DataFormats/PatCandidates/interface/Jet.h"
0044 #include "DataFormats/PatCandidates/interface/JetCorrFactors.h"
0045 #include "DataFormats/VertexReco/interface/Vertex.h"
0046 #include "FWCore/Framework/interface/ESHandle.h"
0047 #include "FWCore/Framework/interface/ESWatcher.h"
0048 #include "FWCore/Framework/interface/Event.h"
0049 #include "FWCore/Framework/interface/EventSetup.h"
0050 #include "FWCore/Framework/interface/stream/EDProducer.h"
0051 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0052 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0053 #include "FWCore/ParameterSet/interface/ParameterDescription.h"
0054 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0055 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0056 #include "FWCore/Utilities/interface/InputTag.h"
0057 #include "JetMETCorrections/Objects/interface/JetCorrectionsRecord.h"
0058
0059 namespace pat {
0060
0061 class JetCorrFactorsProducer : public edm::stream::EDProducer<> {
0062 public:
0063
0064 typedef edm::ValueMap<pat::JetCorrFactors> JetCorrFactorsMap;
0065
0066 typedef std::map<JetCorrFactors::Flavor, std::vector<std::string> > FlavorCorrLevelMap;
0067
0068 public:
0069
0070 explicit JetCorrFactorsProducer(const edm::ParameterSet& cfg);
0071
0072 ~JetCorrFactorsProducer() override {}
0073
0074 void produce(edm::Event& event, const edm::EventSetup& setup) override;
0075
0076 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0077
0078 private:
0079
0080 bool flavorDependent() const { return (levels_.size() > 1); };
0081
0082 std::vector<JetCorrectorParameters> params(const JetCorrectorParametersCollection& parameters,
0083 const std::vector<std::string>& levels) const;
0084
0085
0086
0087 std::vector<std::string> expand(const std::vector<std::string>& levels, const JetCorrFactors::Flavor& flavor);
0088
0089 float evaluate(edm::View<reco::Jet>::const_iterator& jet, const JetCorrFactors::Flavor& flavor, int level);
0090
0091 int numberOf(const edm::Handle<std::vector<reco::Vertex> >& primaryVertices);
0092
0093
0094 private:
0095
0096 bool emf_;
0097
0098 edm::EDGetTokenT<edm::View<reco::Jet> > srcToken_;
0099
0100 std::string type_;
0101
0102 std::string label_;
0103
0104
0105 std::vector<std::string> extraJPTOffset_;
0106
0107 edm::InputTag primaryVertices_;
0108 edm::EDGetTokenT<std::vector<reco::Vertex> > primaryVerticesToken_;
0109
0110 edm::InputTag rho_;
0111 edm::EDGetTokenT<double> rhoToken_;
0112 const edm::ESGetToken<JetCorrectorParametersCollection, JetCorrectionsRecord> parametersToken_;
0113 edm::ESWatcher<JetCorrectionsRecord> parametersWatcher_;
0114
0115 bool useNPV_;
0116 bool useRho_;
0117
0118
0119
0120
0121
0122
0123
0124
0125
0126
0127 FlavorCorrLevelMap levels_;
0128
0129 std::map<JetCorrFactors::Flavor, std::unique_ptr<FactorizedJetCorrector> > correctors_;
0130
0131 std::unique_ptr<FactorizedJetCorrector> extraJPTOffsetCorrector_;
0132 };
0133
0134 inline int JetCorrFactorsProducer::numberOf(const edm::Handle<std::vector<reco::Vertex> >& primaryVertices) {
0135 int npv = 0;
0136 for (auto const& pv : *primaryVertices) {
0137 if (pv.ndof() >= 4)
0138 ++npv;
0139 }
0140 return npv;
0141 }
0142 }
0143
0144 using namespace pat;
0145
0146 JetCorrFactorsProducer::JetCorrFactorsProducer(const edm::ParameterSet& cfg)
0147 : emf_(cfg.getParameter<bool>("emf")),
0148 srcToken_(consumes(cfg.getParameter<edm::InputTag>("src"))),
0149 type_(cfg.getParameter<std::string>("flavorType")),
0150 label_(cfg.getParameter<std::string>("@module_label")),
0151 parametersToken_{esConsumes(edm::ESInputTag("", cfg.getParameter<std::string>("payload")))},
0152 useNPV_(cfg.getParameter<bool>("useNPV")),
0153 useRho_(cfg.getParameter<bool>("useRho")) {
0154 std::vector<std::string> levels = cfg.getParameter<std::vector<std::string> >("levels");
0155
0156
0157
0158
0159
0160
0161
0162
0163 if (std::find(levels.begin(), levels.end(), "L5Flavor") != levels.end() ||
0164 std::find(levels.begin(), levels.end(), "L7Parton") != levels.end()) {
0165 levels_[JetCorrFactors::GLUON] = expand(levels, JetCorrFactors::GLUON);
0166 levels_[JetCorrFactors::UDS] = expand(levels, JetCorrFactors::UDS);
0167 levels_[JetCorrFactors::CHARM] = expand(levels, JetCorrFactors::CHARM);
0168 levels_[JetCorrFactors::BOTTOM] = expand(levels, JetCorrFactors::BOTTOM);
0169 } else {
0170 levels_[JetCorrFactors::NONE] = levels;
0171 }
0172
0173
0174
0175
0176 if (std::find(levels.begin(), levels.end(), "L1JPTOffset") != levels.end()) {
0177 extraJPTOffset_.push_back(cfg.getParameter<std::string>("extraJPTOffset"));
0178 }
0179
0180
0181
0182
0183 if (useNPV_) {
0184 primaryVertices_ = cfg.getParameter<edm::InputTag>("primaryVertices");
0185 primaryVerticesToken_ = mayConsume<std::vector<reco::Vertex> >(primaryVertices_);
0186 }
0187
0188
0189
0190 if (useRho_) {
0191 if ((!extraJPTOffset_.empty() && extraJPTOffset_.front() == std::string("L1FastJet")) ||
0192 std::find(levels.begin(), levels.end(), "L1FastJet") != levels.end()) {
0193 rho_ = cfg.getParameter<edm::InputTag>("rho");
0194 rhoToken_ = mayConsume<double>(rho_);
0195 } else {
0196 edm::LogInfo message("Parameter rho not used");
0197 message << "Module is configured to use the parameter rho, but rho is only used \n"
0198 << "for L1FastJet corrections. The configuration of levels does not contain \n"
0199 << "L1FastJet corrections though, so rho will not be used by this module. \n";
0200 }
0201 }
0202 produces<JetCorrFactorsMap>();
0203 }
0204
0205 std::vector<std::string> JetCorrFactorsProducer::expand(const std::vector<std::string>& levels,
0206 const JetCorrFactors::Flavor& flavor) {
0207 std::vector<std::string> expand;
0208 for (std::vector<std::string>::const_iterator level = levels.begin(); level != levels.end(); ++level) {
0209 if ((*level) == "L5Flavor" || (*level) == "L7Parton") {
0210 if (flavor == JetCorrFactors::GLUON) {
0211 if (*level == "L7Parton" && type_ == "T") {
0212 edm::LogWarning message("L7Parton::GLUON not available");
0213 message << "Jet energy corrections requested for level: L7Parton and type: 'T'. \n"
0214 << "For this combination there is no GLUON correction available. The \n"
0215 << "correction for this flavor type will be taken from 'J'.";
0216 }
0217 expand.push_back(std::string(*level).append("_").append("g").append("J"));
0218 }
0219 if (flavor == JetCorrFactors::UDS)
0220 expand.push_back(std::string(*level).append("_").append("q").append(type_));
0221 if (flavor == JetCorrFactors::CHARM)
0222 expand.push_back(std::string(*level).append("_").append("c").append(type_));
0223 if (flavor == JetCorrFactors::BOTTOM)
0224 expand.push_back(std::string(*level).append("_").append("b").append(type_));
0225 } else {
0226 expand.push_back(*level);
0227 }
0228 }
0229 return expand;
0230 }
0231
0232 std::vector<JetCorrectorParameters> JetCorrFactorsProducer::params(const JetCorrectorParametersCollection& parameters,
0233 const std::vector<std::string>& levels) const {
0234 std::vector<JetCorrectorParameters> params;
0235 for (std::vector<std::string>::const_iterator level = levels.begin(); level != levels.end(); ++level) {
0236 const JetCorrectorParameters& ip = parameters[*level];
0237 params.push_back(ip);
0238 }
0239 return params;
0240 }
0241
0242 float JetCorrFactorsProducer::evaluate(edm::View<reco::Jet>::const_iterator& jet,
0243 const JetCorrFactors::Flavor& flavor,
0244 int level) {
0245 std::unique_ptr<FactorizedJetCorrector>& corrector = correctors_.find(flavor)->second;
0246
0247 const reco::JPTJet* jpt = dynamic_cast<reco::JPTJet const*>(&*jet);
0248 if (jpt) {
0249 TLorentzVector p4;
0250 p4.SetPtEtaPhiE(jpt->getCaloJetRef()->pt(),
0251 jpt->getCaloJetRef()->eta(),
0252 jpt->getCaloJetRef()->phi(),
0253 jpt->getCaloJetRef()->energy());
0254 if (extraJPTOffsetCorrector_) {
0255 extraJPTOffsetCorrector_->setJPTrawP4(p4);
0256 corrector->setJPTrawOff(extraJPTOffsetCorrector_->getSubCorrections()[0]);
0257 }
0258 corrector->setJPTrawP4(p4);
0259 }
0260
0261 const Jet* patjet = dynamic_cast<Jet const*>(&*jet);
0262 if (patjet) {
0263 corrector->setJetEta(patjet->correctedP4(0).eta());
0264 corrector->setJetPt(patjet->correctedP4(0).pt());
0265 corrector->setJetPhi(patjet->correctedP4(0).phi());
0266 corrector->setJetE(patjet->correctedP4(0).energy());
0267 } else {
0268 corrector->setJetEta(jet->eta());
0269 corrector->setJetPt(jet->pt());
0270 corrector->setJetPhi(jet->phi());
0271 corrector->setJetE(jet->energy());
0272 }
0273 if (emf_ && dynamic_cast<const reco::CaloJet*>(&*jet)) {
0274 corrector->setJetEMF(dynamic_cast<const reco::CaloJet*>(&*jet)->emEnergyFraction());
0275 }
0276 return corrector->getSubCorrections()[level];
0277 }
0278
0279 void JetCorrFactorsProducer::produce(edm::Event& event, const edm::EventSetup& setup) {
0280
0281 edm::Handle<edm::View<reco::Jet> > jets;
0282 event.getByToken(srcToken_, jets);
0283
0284
0285 edm::Handle<std::vector<reco::Vertex> > primaryVertices;
0286 if (!primaryVertices_.label().empty())
0287 event.getByToken(primaryVerticesToken_, primaryVertices);
0288
0289
0290 edm::Handle<double> rho;
0291 if (!rho_.label().empty())
0292 event.getByToken(rhoToken_, rho);
0293
0294 if (parametersWatcher_.check(setup)) {
0295
0296 auto const& parameters = setup.getData(parametersToken_);
0297
0298 for (auto const& flavor : levels_) {
0299 correctors_[flavor.first] = std::make_unique<FactorizedJetCorrector>(params(parameters, flavor.second));
0300 }
0301
0302 if (!extraJPTOffset_.empty()) {
0303 extraJPTOffsetCorrector_ = std::make_unique<FactorizedJetCorrector>(params(parameters, extraJPTOffset_));
0304 }
0305 }
0306
0307
0308 std::vector<JetCorrFactors> jcfs;
0309 for (edm::View<reco::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); ++jet) {
0310
0311
0312
0313
0314
0315
0316 std::vector<JetCorrFactors::CorrectionFactor> jec;
0317 jec.emplace_back("Uncorrected", std::vector<float>(1, 1));
0318
0319
0320
0321
0322
0323 FlavorCorrLevelMap::const_iterator corrLevel = levels_.begin();
0324 if (corrLevel == levels_.end()) {
0325 throw cms::Exception("No JECFactors")
0326 << "You request to create a jetCorrFactors object with no JEC Levels indicated. \n"
0327 << "This makes no sense, either you should correct this or drop the module from \n"
0328 << "the sequence.";
0329 }
0330 for (unsigned int idx = 0; idx < corrLevel->second.size(); ++idx) {
0331 std::vector<float> factors;
0332 if (corrLevel->second[idx].find("L5Flavor") != std::string::npos ||
0333 corrLevel->second[idx].find("L7Parton") != std::string::npos) {
0334 for (FlavorCorrLevelMap::const_iterator flavor = corrLevel; flavor != levels_.end(); ++flavor) {
0335 if (!primaryVertices_.label().empty()) {
0336
0337
0338 correctors_.find(flavor->first)->second->setNPV(numberOf(primaryVertices));
0339 }
0340 if (!rho_.label().empty()) {
0341
0342
0343 correctors_.find(flavor->first)->second->setRho(*rho);
0344 correctors_.find(flavor->first)->second->setJetA(jet->jetArea());
0345 }
0346 factors.push_back(evaluate(jet, flavor->first, idx));
0347 }
0348 } else {
0349 if (!primaryVertices_.label().empty()) {
0350
0351
0352 correctors_.find(corrLevel->first)->second->setNPV(numberOf(primaryVertices));
0353 }
0354 if (!rho_.label().empty()) {
0355
0356
0357 correctors_.find(corrLevel->first)->second->setRho(*rho);
0358 correctors_.find(corrLevel->first)->second->setJetA(jet->jetArea());
0359 }
0360 factors.push_back(evaluate(jet, corrLevel->first, idx));
0361 }
0362
0363
0364
0365
0366
0367
0368
0369
0370 jec.emplace_back((corrLevel->second[idx]).substr(0, (corrLevel->second[idx]).find('_')), factors);
0371 }
0372
0373
0374 jcfs.emplace_back(label_, jec);
0375 }
0376
0377 auto jetCorrsMap = std::make_unique<JetCorrFactorsMap>();
0378 JetCorrFactorsMap::Filler filler(*jetCorrsMap);
0379
0380 filler.insert(jets, jcfs.begin(), jcfs.end());
0381 filler.fill();
0382
0383 event.put(std::move(jetCorrsMap));
0384 }
0385
0386 void JetCorrFactorsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0387 edm::ParameterSetDescription iDesc;
0388 iDesc.add<bool>("emf", false);
0389 iDesc.add<std::string>("flavorType", "J");
0390 iDesc.add<edm::InputTag>("src", edm::InputTag("ak5CaloJets"));
0391 iDesc.add<std::string>("payload", "AK5Calo");
0392 iDesc.add<bool>("useNPV", true);
0393 iDesc.add<edm::InputTag>("primaryVertices", edm::InputTag("offlinePrimaryVertices"));
0394 iDesc.add<bool>("useRho", true);
0395 iDesc.add<edm::InputTag>("rho", edm::InputTag("fixedGridRhoFastjetAllCalo"));
0396 iDesc.add<std::string>("extraJPTOffset", "L1Offset");
0397
0398 iDesc.add<std::vector<std::string> >("levels",
0399 {
0400 "L1Offset",
0401 "L2Relative",
0402 "L3Absolute",
0403 "L2L3Residual",
0404 "L5Flavor",
0405 "L7Parton",
0406 });
0407 descriptions.add("JetCorrFactorsProducer", iDesc);
0408 }
0409
0410 #include "FWCore/Framework/interface/MakerMacros.h"
0411 DEFINE_FWK_MODULE(JetCorrFactorsProducer);