File indexing completed on 2024-09-07 04:35:30
0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002
0003 #include "CondCore/Utilities/interface/PayloadInspectorModule.h"
0004 #include "CondCore/Utilities/interface/PayloadInspector.h"
0005 #include "CondCore/CondDB/interface/Time.h"
0006
0007
0008 #include "CondFormats/JetMETObjects/interface/JetCorrectorParameters.h"
0009 #include "CondFormats/JetMETObjects/interface/JetCorrectorParametersHelper.h"
0010 #include "CondFormats/JetMETObjects/interface/JetCorrectionUncertainty.h"
0011 #include "CondFormats/JetMETObjects/interface/FactorizedJetCorrector.h"
0012 #include "CondFormats/JetMETObjects/interface/SimpleJetCorrectionUncertainty.h"
0013 #include "CondFormats/JetMETObjects/interface/Utilities.h"
0014
0015 #include "CommonTools/Utils/interface/FormulaEvaluator.h"
0016
0017 #include <memory>
0018 #include <sstream>
0019 #include <fstream>
0020 #include <iostream>
0021
0022
0023 #include "TH2F.h"
0024 #include "TLegend.h"
0025 #include "TCanvas.h"
0026 #include "TLine.h"
0027 #include "TStyle.h"
0028 #include "TLatex.h"
0029 #include "TPave.h"
0030 #include "TPaveStats.h"
0031
0032 #define MIN_ETA -5.05
0033 #define MAX_ETA 5.05
0034 #define NBIN_ETA 51
0035 #define MIN_PT -5.
0036 #define MAX_PT 3005.
0037 #define NBIN_PT 301
0038
0039 namespace {
0040
0041 using namespace cond::payloadInspector;
0042
0043 enum levelt {
0044 L1Offset = 0,
0045 L1JPTOffset = 7,
0046 L1FastJet = 10,
0047 L2Relative = 1,
0048 L3Absolute = 2,
0049 L2L3Residual = 8,
0050 L4EMF = 3,
0051 L5Flavor = 4,
0052 L6UE = 5,
0053 L7Parton = 6,
0054 Uncertainty = 9,
0055 UncertaintyAbsolute = 11,
0056 UncertaintyHighPtExtra = 12,
0057 UncertaintySinglePionECAL = 13,
0058 UncertaintySinglePionHCAL = 27,
0059 UncertaintyFlavor = 14,
0060 UncertaintyTime = 15,
0061 UncertaintyRelativeJEREC1 = 16,
0062 UncertaintyRelativeJEREC2 = 17,
0063 UncertaintyRelativeJERHF = 18,
0064 UncertaintyRelativePtEC1 = 28,
0065 UncertaintyRelativePtEC2 = 29,
0066 UncertaintyRelativePtHF = 30,
0067 UncertaintyRelativeStatEC2 = 19,
0068 UncertaintyRelativeStatHF = 20,
0069 UncertaintyRelativeFSR = 21,
0070 UncertaintyRelativeSample = 31,
0071 UncertaintyPileUpDataMC = 22,
0072 UncertaintyPileUpOOT = 23,
0073 UncertaintyPileUpPtBB = 24,
0074 UncertaintyPileUpPtEC = 32,
0075 UncertaintyPileUpPtHF = 33,
0076 UncertaintyPileUpBias = 25,
0077 UncertaintyPileUpJetRate = 26,
0078 L1RC = 34,
0079 L1Residual = 35,
0080 UncertaintyAux3 = 36,
0081 UncertaintyAux4 = 37,
0082 N_LEVELS = 38
0083 };
0084
0085 const std::vector<std::string> labels_ = {
0086 "L1Offset",
0087 "L2Relative",
0088 "L3Absolute",
0089 "L4EMF",
0090 "L5Flavor",
0091 "L6UE",
0092 "L7Parton",
0093 "L1JPTOffset",
0094 "L2L3Residual",
0095 "Uncertainty",
0096 "L1FastJet",
0097 "UncertaintyAbsolute",
0098 "UncertaintyHighPtExtra",
0099 "UncertaintySinglePionECAL",
0100 "UncertaintyFlavor",
0101 "UncertaintyTime",
0102 "UncertaintyRelativeJEREC1",
0103 "UncertaintyRelativeJEREC2",
0104 "UncertaintyRelativeJERHF",
0105 "UncertaintyRelativeStatEC2",
0106 "UncertaintyRelativeStatHF",
0107 "UncertaintyRelativeFSR",
0108 "UncertaintyPileUpDataMC",
0109 "UncertaintyPileUpOOT",
0110 "UncertaintyPileUpPtBB",
0111 "UncertaintyPileUpBias",
0112 "UncertaintyPileUpJetRate",
0113 "UncertaintySinglePionHCAL",
0114 "UncertaintyRelativePtEC1",
0115 "UncertaintyRelativePtEC2",
0116 "UncertaintyRelativePtHF",
0117 "UncertaintyRelativeSample",
0118 "UncertaintyPileUpPtEC",
0119 "UncertaintyPileUpPtHF",
0120 "L1RC",
0121 "L1Residual",
0122 "UncertaintyAux3",
0123 "UncertaintyAux4",
0124 };
0125
0126 bool fill_eta_hist(JetCorrectorParameters const& JCParam,
0127 TH1D* hist,
0128 const std::map<std::string, std::string>& paramValues) {
0129 if (!(JCParam.isValid())) {
0130 edm::LogWarning("JEC_PI") << "JetCorrectorParameter is not valid.";
0131 return false;
0132 }
0133
0134 std::vector<double> vars;
0135 std::vector<float> bins;
0136 std::vector<double> params;
0137
0138 double par_JetPt = 100.;
0139 double par_JetEta = 0.;
0140 double par_JetA = 0.5;
0141 double par_Rho = 20.;
0142 double par_JetPhi = 0.;
0143 double par_JetE = 150.;
0144
0145
0146 auto ip = paramValues.find("Jet_Pt");
0147 if (ip != paramValues.end()) {
0148 par_JetPt = std::stod(ip->second);
0149 }
0150 ip = paramValues.find("Jet_Rho");
0151 if (ip != paramValues.end()) {
0152 par_Rho = std::stod(ip->second);
0153 }
0154
0155 int ir = -1;
0156
0157 for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
0158 par_JetEta = (idx + 0.5) * (MAX_ETA - MIN_ETA) / NBIN_ETA + MIN_ETA;
0159
0160 if (JCParam.definitions().formula().compare("1") == 0) {
0161 hist->SetBinContent(idx + 1, 1.);
0162 continue;
0163 } else if (JCParam.definitions().level().compare("Uncertainty") == 0 &&
0164 JCParam.definitions().formula().compare("\"\"") == 0) {
0165 JetCorrectionUncertainty JCU(JCParam);
0166 JCU.setJetEta(par_JetEta);
0167 JCU.setJetPt(par_JetPt);
0168 JCU.setJetPhi(par_JetPhi);
0169 JCU.setJetE(par_JetE);
0170
0171 float unc = JCU.getUncertainty(true);
0172 hist->SetBinContent(idx + 1, unc);
0173
0174 continue;
0175 }
0176
0177 reco::FormulaEvaluator formula(JCParam.definitions().formula());
0178
0179 vars.clear();
0180 bins.clear();
0181 params.clear();
0182
0183 for (size_t i = 0; i < JCParam.definitions().nBinVar(); i++) {
0184
0185 if (JCParam.definitions().binVar(i).compare("JetPt") == 0)
0186 bins.push_back(par_JetPt);
0187 if (JCParam.definitions().binVar(i).compare("JetEta") == 0)
0188 bins.push_back(par_JetEta);
0189 if (JCParam.definitions().binVar(i).compare("JetA") == 0)
0190 bins.push_back(par_JetA);
0191 if (JCParam.definitions().binVar(i).compare("Rho") == 0)
0192 bins.push_back(par_Rho);
0193 }
0194
0195 for (size_t i = 0; i < JCParam.definitions().nParVar(); i++) {
0196
0197 if (JCParam.definitions().parVar(i).compare("JetPt") == 0)
0198 vars.push_back(par_JetPt);
0199 if (JCParam.definitions().parVar(i).compare("JetEta") == 0)
0200 vars.push_back(par_JetEta);
0201 if (JCParam.definitions().parVar(i).compare("JetA") == 0)
0202 vars.push_back(par_JetA);
0203 if (JCParam.definitions().parVar(i).compare("Rho") == 0)
0204 vars.push_back(par_Rho);
0205 }
0206
0207 ir = JCParam.binIndex(bins);
0208
0209 if (ir < 0 || ir > (int)JCParam.size())
0210 continue;
0211
0212
0213 for (size_t i = 2 * vars.size(); i < JCParam.record(ir).nParameters(); i++) {
0214 double par = JCParam.record(ir).parameter(i);
0215 params.push_back(par);
0216 }
0217
0218 double jec = formula.evaluate(vars, params);
0219 hist->SetBinContent(idx + 1, jec);
0220
0221 }
0222 return true;
0223
0224 }
0225
0226
0227
0228
0229
0230
0231
0232
0233
0234 template <levelt ii>
0235 class JetCorrectorVsEta : public cond::payloadInspector::Histogram1D<JetCorrectorParametersCollection, SINGLE_IOV> {
0236 public:
0237 JetCorrectorVsEta()
0238 : cond::payloadInspector::Histogram1D<JetCorrectorParametersCollection, SINGLE_IOV>(
0239 "Jet Corrector", "#eta", NBIN_ETA, MIN_ETA, MAX_ETA, "Corrector") {
0240 cond::payloadInspector::PlotBase::addInputParam("Jet_Pt");
0241 cond::payloadInspector::PlotBase::addInputParam("Jet_Rho");
0242 }
0243
0244 bool fill() override {
0245 double par_JetPt = 100.;
0246 double par_Rho = 20.;
0247
0248
0249 auto paramValues = cond::payloadInspector::PlotBase::inputParamValues();
0250 auto ip = paramValues.find("Jet_Pt");
0251 if (ip != paramValues.end()) {
0252 par_JetPt = std::stod(ip->second);
0253 edm::LogPrint("JEC_PI") << "Jet Pt: " << par_JetPt;
0254 }
0255 ip = paramValues.find("Jet_Rho");
0256 if (ip != paramValues.end()) {
0257 par_Rho = std::stod(ip->second);
0258 edm::LogPrint("JEC_PI") << "Rho: " << par_Rho;
0259 }
0260
0261 TH1D* jec_hist = new TH1D("JEC vs. #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0262
0263 auto tag = PlotBase::getTag<0>();
0264 for (auto const& iov : tag.iovs) {
0265 std::shared_ptr<JetCorrectorParametersCollection> payload = Base::fetchPayload(std::get<1>(iov));
0266 if (payload.get()) {
0267
0268 std::vector<key_t> keys;
0269
0270 (*payload).validKeys(keys);
0271
0272 if (std::find(keys.begin(), keys.end(), ii) == keys.end()) {
0273 edm::LogWarning("JEC_PI") << "Jet corrector level " << (*payload).findLabel(ii) << " is not available.";
0274 return false;
0275 }
0276 auto JCParams = (*payload)[ii];
0277
0278 edm::LogInfo("JEC_PI") << "Jet corrector level " << (*payload).findLabel(ii) << " as "
0279 << JCParams.definitions().level() << " has " << JCParams.definitions().nParVar()
0280 << " parameters with " << JCParams.definitions().nBinVar() << " bin(s).";
0281 edm::LogInfo("JEC_PI") << "JCParam size: " << JCParams.size();
0282 edm::LogInfo("JEC_PI") << "JCParam def: " << JCParams.definitions().formula();
0283
0284 fill_eta_hist(JCParams, jec_hist, paramValues);
0285
0286 for (size_t idx = 0; idx <= NBIN_ETA; idx++) {
0287 fillWithBinAndValue(idx + 1, jec_hist->GetBinContent(idx + 1));
0288 }
0289
0290 return true;
0291 }
0292 return false;
0293 }
0294 return false;
0295 }
0296 };
0297
0298 class JetCorrectorVsEtaSummary
0299 : public cond::payloadInspector::PlotImage<JetCorrectorParametersCollection, SINGLE_IOV> {
0300 public:
0301 JetCorrectorVsEtaSummary()
0302 : cond::payloadInspector::PlotImage<JetCorrectorParametersCollection, SINGLE_IOV>("Jet Correction Summary") {
0303 cond::payloadInspector::PlotBase::addInputParam("Jet_Pt");
0304 cond::payloadInspector::PlotBase::addInputParam("Jet_Rho");
0305 }
0306
0307 bool fill() override {
0308 double par_JetPt = 100.;
0309 double par_JetA = 0.5;
0310 double par_Rho = 40.;
0311
0312
0313 auto paramValues = cond::payloadInspector::PlotBase::inputParamValues();
0314 auto ip = paramValues.find("Jet_Pt");
0315 if (ip != paramValues.end()) {
0316 par_JetPt = std::stod(ip->second);
0317 }
0318 ip = paramValues.find("Jet_Rho");
0319 if (ip != paramValues.end()) {
0320 par_Rho = std::stod(ip->second);
0321 }
0322
0323 TH1D* jec_l1fj = new TH1D("JEC L1FastJet vs. #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0324 TH1D* jec_l2rel = new TH1D("JEC L2Relative vs. #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0325 TH1D* jec_l2l3 = new TH1D("JEC L2L3Residual vs. #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0326 TH1D* jec_l1rc = new TH1D("JEC L1RC vs. #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0327 TH1D* jec_uncert = new TH1D("JEC Uncertainty vs. #eta", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0328
0329 TLegend* leg_eta = new TLegend(0.50, 0.73, 0.935, 0.90);
0330 TLegend* leg_eta2 = new TLegend(0.50, 0.83, 0.935, 0.90);
0331
0332 leg_eta->SetBorderSize(0);
0333 leg_eta->SetLineStyle(0);
0334 leg_eta->SetFillStyle(0);
0335 leg_eta->SetTextFont(42);
0336
0337 leg_eta2->SetBorderSize(0);
0338 leg_eta2->SetLineStyle(0);
0339 leg_eta2->SetFillStyle(0);
0340 leg_eta2->SetTextFont(42);
0341
0342 auto tag = PlotBase::getTag<0>();
0343 auto iov = tag.iovs.front();
0344
0345 std::shared_ptr<JetCorrectorParametersCollection> payload = fetchPayload(std::get<1>(iov));
0346 unsigned int run = std::get<0>(iov);
0347 std::string tagname = tag.name;
0348 std::stringstream ss_tagname(tag.name);
0349 std::string stmp;
0350
0351 std::string tag_ver;
0352
0353 getline(ss_tagname, stmp, '_');
0354 getline(ss_tagname, stmp);
0355 tag_ver = stmp;
0356
0357 if (payload.get()) {
0358 auto JCParam_L1FJ = (*payload)[L1FastJet];
0359 auto JCParam_L2 = (*payload)[L2Relative];
0360 auto JCParam_L2L3 = (*payload)[L2L3Residual];
0361 auto JCParam_L1RC = (*payload)[L1RC];
0362 auto JCParam_Unc = (*payload)[Uncertainty];
0363
0364 auto paramValues = cond::payloadInspector::PlotBase::inputParamValues();
0365
0366 fill_eta_hist(JCParam_L1FJ, jec_l1fj, paramValues);
0367 fill_eta_hist(JCParam_L2, jec_l2rel, paramValues);
0368 fill_eta_hist(JCParam_L2L3, jec_l2l3, paramValues);
0369 fill_eta_hist(JCParam_L1RC, jec_l1rc, paramValues);
0370 fill_eta_hist(JCParam_Unc, jec_uncert, paramValues);
0371
0372 gStyle->SetOptStat(0);
0373 gStyle->SetLabelFont(42, "XYZ");
0374 gStyle->SetLabelSize(0.05, "XYZ");
0375 gStyle->SetFrameLineWidth(3);
0376
0377 std::string title = Form("Summary Run %i", run);
0378 TCanvas canvas("Jet Energy Correction", title.c_str(), 800, 1200);
0379 canvas.Divide(1, 2);
0380
0381 canvas.cd(1);
0382 jec_l1fj->SetTitle(tag_ver.c_str());
0383 jec_l1fj->SetXTitle("#eta");
0384 jec_l1fj->SetMaximum(1.6);
0385 jec_l1fj->SetMinimum(0.0);
0386 jec_l1fj->SetLineWidth(3);
0387 jec_l1fj->Draw("][");
0388
0389 jec_l2rel->SetLineColor(2);
0390 jec_l2rel->SetLineWidth(3);
0391 jec_l2rel->Draw("][same");
0392
0393 jec_l2l3->SetLineColor(8);
0394 jec_l2l3->SetLineWidth(3);
0395 jec_l2l3->Draw("][same");
0396
0397 jec_l1rc->SetLineColor(9);
0398 jec_l1rc->SetLineWidth(3);
0399 jec_l1rc->Draw("][same");
0400
0401 leg_eta->AddEntry(jec_l1fj, (*payload).findLabel(L1FastJet).c_str(), "l");
0402 leg_eta->AddEntry(jec_l2rel, (*payload).findLabel(L2Relative).c_str(), "l");
0403 leg_eta->AddEntry(jec_l2l3, (*payload).findLabel(L2L3Residual).c_str(), "l");
0404 leg_eta->AddEntry(jec_l1rc, (*payload).findLabel(L1RC).c_str(), "l");
0405 leg_eta->AddEntry((TObject*)nullptr, Form("JetPt=%.2f; JetA=%.2f; Rho=%.2f", par_JetPt, par_JetA, par_Rho), "");
0406 leg_eta->Draw();
0407
0408 canvas.cd(2);
0409 jec_uncert->SetTitle(tag_ver.c_str());
0410 jec_uncert->SetXTitle("#eta");
0411 jec_uncert->SetMaximum(0.1);
0412 jec_uncert->SetMinimum(0.0);
0413 jec_uncert->SetLineColor(6);
0414 jec_uncert->SetLineWidth(3);
0415 jec_uncert->Draw("][");
0416
0417 leg_eta2->AddEntry(jec_uncert, (*payload).findLabel(Uncertainty).c_str(), "l");
0418 leg_eta2->AddEntry(
0419 (TObject*)nullptr, Form("JetPt=%.2f; JetA=%.2f; Rho=%.2f", par_JetPt, par_JetA, par_Rho), "");
0420 leg_eta2->Draw();
0421
0422 canvas.SaveAs(m_imageFileName.c_str());
0423
0424 return true;
0425 } else
0426 return false;
0427 }
0428
0429 };
0430
0431 class JetCorrectorVsEtaCompare
0432 : public cond::payloadInspector::PlotImage<JetCorrectorParametersCollection, SINGLE_IOV, 2> {
0433 public:
0434 JetCorrectorVsEtaCompare()
0435 : cond::payloadInspector::PlotImage<JetCorrectorParametersCollection, SINGLE_IOV, 2>(
0436 "Jet Correction Compare Two Tags") {
0437 cond::payloadInspector::PlotBase::addInputParam("Corr_Level");
0438 cond::payloadInspector::PlotBase::addInputParam("Jet_Pt");
0439 cond::payloadInspector::PlotBase::addInputParam("Jet_Rho");
0440 }
0441
0442 bool fill() override {
0443 double par_JetPt = 100.;
0444 double par_JetA = 0.5;
0445 double par_Rho = 40.;
0446 key_t par_Level = L1FastJet;
0447
0448
0449 auto paramValues = cond::payloadInspector::PlotBase::inputParamValues();
0450 auto ip = paramValues.find("Jet_Pt");
0451 if (ip != paramValues.end()) {
0452 par_JetPt = std::stod(ip->second);
0453 }
0454 ip = paramValues.find("Jet_Rho");
0455 if (ip != paramValues.end()) {
0456 par_Rho = std::stod(ip->second);
0457 }
0458
0459 TH1D* jec_one = new TH1D("JEC vs. #eta one", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0460 TH1D* jec_two = new TH1D("JEC vs. #eta two", "", NBIN_ETA, MIN_ETA, MAX_ETA);
0461
0462 TLegend* leg_eta = new TLegend(0.50, 0.73, 0.935, 0.90);
0463
0464 leg_eta->SetBorderSize(0);
0465 leg_eta->SetLineStyle(0);
0466 leg_eta->SetFillStyle(0);
0467 leg_eta->SetTextFont(42);
0468
0469 auto tag1 = PlotBase::getTag<0>();
0470 auto iov1 = tag1.iovs.front();
0471
0472 auto tag2 = PlotBase::getTag<1>();
0473 auto iov2 = tag2.iovs.front();
0474
0475 std::shared_ptr<JetCorrectorParametersCollection> payload1 = fetchPayload(std::get<1>(iov1));
0476 std::string tagname1 = tag1.name;
0477
0478 std::shared_ptr<JetCorrectorParametersCollection> payload2 = fetchPayload(std::get<1>(iov2));
0479 std::string tagname2 = tag2.name;
0480
0481 ip = paramValues.find("Corr_Level");
0482 if (ip != paramValues.end()) {
0483
0484 if (ip->second.length() < 3) {
0485 int ii = std::stoi(ip->second);
0486 par_Level = static_cast<key_t>(ii);
0487 } else {
0488
0489 std::vector<std::string>::const_iterator found3 = find(labels_.begin(), labels_.end(), ip->second);
0490 if (found3 != labels_.end()) {
0491 par_Level = static_cast<key_t>(found3 - labels_.begin());
0492 }
0493 }
0494 }
0495
0496 std::string stmp;
0497 std::stringstream ss_tagname(tag1.name);
0498
0499 std::string tag_ver1;
0500 std::string tag_ver2;
0501
0502 getline(ss_tagname, stmp, '_');
0503 getline(ss_tagname, stmp);
0504 tag_ver1 = stmp;
0505
0506 std::stringstream ss_tagname2(tag2.name);
0507
0508 getline(ss_tagname2, stmp, '_');
0509 getline(ss_tagname2, stmp);
0510 tag_ver2 = stmp;
0511
0512 if (payload1.get() && payload2.get()) {
0513
0514 std::vector<key_t> keys;
0515
0516 (*payload1).validKeys(keys);
0517
0518 if (std::find(keys.begin(), keys.end(), par_Level) == keys.end()) {
0519 edm::LogWarning("JEC_PI") << "Jet corrector level " << (*payload1).findLabel(par_Level)
0520 << " is not available for tag one.";
0521 return false;
0522 }
0523
0524 keys.clear();
0525 (*payload2).validKeys(keys);
0526 if (std::find(keys.begin(), keys.end(), par_Level) == keys.end()) {
0527 edm::LogWarning("JEC_PI") << "Jet corrector level " << (*payload2).findLabel(par_Level)
0528 << " is not available for twg two.";
0529 return false;
0530 }
0531
0532 auto JCParam_one = (*payload1)[par_Level];
0533 auto JCParam_two = (*payload2)[par_Level];
0534
0535 fill_eta_hist(JCParam_one, jec_one, paramValues);
0536 fill_eta_hist(JCParam_two, jec_two, paramValues);
0537
0538 gStyle->SetOptStat(0);
0539 gStyle->SetLabelFont(42, "XYZ");
0540 gStyle->SetLabelSize(0.05, "XYZ");
0541 gStyle->SetFrameLineWidth(3);
0542
0543 std::string title = Form("Comparison between %s and %s", tag1.name.c_str(), tag2.name.c_str());
0544 TCanvas canvas("Jet Energy Correction", title.c_str(), 800, 600);
0545
0546 canvas.cd();
0547 jec_one->SetTitle(("JetCorrector comparison for " + (*payload1).findLabel(par_Level)).c_str());
0548 jec_one->SetXTitle("#eta");
0549 jec_one->SetMaximum(1.6);
0550 if (par_Level == Uncertainty) {
0551 jec_one->SetMaximum(0.1);
0552 }
0553 jec_one->SetMinimum(0.0);
0554 jec_one->SetLineWidth(3);
0555 jec_one->Draw("][");
0556
0557 jec_two->SetLineColor(2);
0558 jec_two->SetLineWidth(3);
0559 jec_two->SetLineStyle(2);
0560 jec_two->Draw("][same");
0561
0562 leg_eta->AddEntry(jec_one, tag_ver1.c_str(), "l");
0563 leg_eta->AddEntry(jec_two, tag_ver2.c_str(), "l");
0564 leg_eta->AddEntry((TObject*)nullptr, Form("JetPt=%.2f; JetA=%.2f; Rho=%.2f", par_JetPt, par_JetA, par_Rho), "");
0565 leg_eta->Draw();
0566
0567 canvas.SaveAs(m_imageFileName.c_str());
0568
0569 return true;
0570 } else
0571 return false;
0572 }
0573
0574 };
0575
0576 typedef JetCorrectorVsEta<L1Offset> JetCorrectorVsEtaL1Offset;
0577 typedef JetCorrectorVsEta<L1FastJet> JetCorrectorVsEtaL1FastJet;
0578 typedef JetCorrectorVsEta<L2Relative> JetCorrectorVsEtaL2Relative;
0579 typedef JetCorrectorVsEta<L3Absolute> JetCorrectorVsEtaL3Absolute;
0580 typedef JetCorrectorVsEta<L2L3Residual> JetCorrectorVsEtaL2L3Residual;
0581 typedef JetCorrectorVsEta<Uncertainty> JetCorrectorVsEtaUncertainty;
0582 typedef JetCorrectorVsEta<L1RC> JetCorrectorVsEtaL1RC;
0583
0584
0585 PAYLOAD_INSPECTOR_MODULE(JetCorrectorParametersCollection) {
0586 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaL1Offset);
0587 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaL1FastJet);
0588 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaL2Relative);
0589 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaL3Absolute);
0590 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaL2L3Residual);
0591 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaUncertainty);
0592 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaL1RC);
0593
0594 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaSummary);
0595 PAYLOAD_INSPECTOR_CLASS(JetCorrectorVsEtaCompare);
0596 }
0597
0598 }