Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:36

0001 /**
0002  *  @file     L1TPhase2CorrelatorOffline.cc
0003  *  @authors  Dylan Rankin (MIT)
0004  *  @date     20/10/2020
0005  *  @version  0.1
0006  *
0007  */
0008 
0009 #include "DQMOffline/L1Trigger/interface/L1TPhase2CorrelatorOffline.h"
0010 
0011 #include "TLorentzVector.h"
0012 #include "TGraph.h"
0013 
0014 #include <iomanip>
0015 #include <cstdio>
0016 #include <string>
0017 #include <sstream>
0018 #include <cmath>
0019 #include <algorithm>
0020 
0021 using namespace reco;
0022 using namespace trigger;
0023 using namespace edm;
0024 using namespace std;
0025 
0026 const std::map<std::string, unsigned int> L1TPhase2CorrelatorOffline::PlotConfigNames = {
0027     {"resVsPt", PlotConfig::resVsPt},
0028     {"resVsEta", PlotConfig::resVsEta},
0029     {"ptDist", PlotConfig::ptDist},
0030     {"etaDist", PlotConfig::etaDist}};
0031 
0032 //
0033 // -------------------------------------- Constructor --------------------------------------------
0034 //
0035 L1TPhase2CorrelatorOffline::L1TPhase2CorrelatorOffline(const edm::ParameterSet& ps)
0036     : genJetToken_(consumes<std::vector<reco::GenJet>>(ps.getUntrackedParameter<edm::InputTag>("genJetsInputTag"))),
0037       genParticleToken_(
0038           consumes<std::vector<reco::GenParticle>>(ps.getUntrackedParameter<edm::InputTag>("genParticlesInputTag"))),
0039       BFieldTag_{esConsumes<MagneticField, IdealMagneticFieldRecord, edm::Transition::BeginRun>()},
0040       objs_(ps.getParameter<edm::ParameterSet>("objects")),
0041       isParticleGun_(ps.getParameter<bool>("isParticleGun")),
0042       histFolder_(ps.getParameter<std::string>("histFolder")),
0043       respresolFolder_(histFolder_ + "/respresol_raw"),
0044       histDefinitions_(dqmoffline::l1t::readHistDefinitions(ps.getParameterSet("histDefinitions"), PlotConfigNames)),
0045       h_L1PF_pt_(),
0046       h_L1PF_eta_(),
0047       h_L1Puppi_pt_(),
0048       h_L1Puppi_eta_(),
0049       h_L1PF_pt_mu_(),
0050       h_L1PF_eta_mu_(),
0051       h_L1Puppi_pt_mu_(),
0052       h_L1Puppi_eta_mu_(),
0053       h_L1PF_pt_el_(),
0054       h_L1PF_eta_el_(),
0055       h_L1Puppi_pt_el_(),
0056       h_L1Puppi_eta_el_(),
0057       h_L1PF_pt_pho_(),
0058       h_L1PF_eta_pho_(),
0059       h_L1Puppi_pt_pho_(),
0060       h_L1Puppi_eta_pho_(),
0061       h_L1PF_pt_ch_(),
0062       h_L1PF_eta_ch_(),
0063       h_L1Puppi_pt_ch_(),
0064       h_L1Puppi_eta_ch_(),
0065       h_L1PF_pt_nh_(),
0066       h_L1PF_eta_nh_(),
0067       h_L1Puppi_pt_nh_(),
0068       h_L1Puppi_eta_nh_(),
0069       h_L1PF_part_ptratio_0p2_vs_pt_barrel_(),
0070       h_L1PF_part_ptratio_0p2_vs_pt_endcap_(),
0071       h_L1PF_part_ptratio_0p2_vs_pt_ecnotk_(),
0072       h_L1PF_part_ptratio_0p2_vs_pt_hf_(),
0073       h_L1PF_part_ptratio_0p2_vs_eta_(),
0074       h_L1Puppi_part_ptratio_0p2_vs_pt_barrel_(),
0075       h_L1Puppi_part_ptratio_0p2_vs_pt_endcap_(),
0076       h_L1Puppi_part_ptratio_0p2_vs_pt_ecnotk_(),
0077       h_L1Puppi_part_ptratio_0p2_vs_pt_hf_(),
0078       h_L1Puppi_part_ptratio_0p2_vs_eta_(),
0079       h_L1PF_jet_ptratio_vs_pt_barrel_(),
0080       h_L1PF_jet_ptratio_vs_pt_endcap_(),
0081       h_L1PF_jet_ptratio_vs_pt_ecnotk_(),
0082       h_L1PF_jet_ptratio_vs_pt_hf_(),
0083       h_L1PF_jet_ptratio_vs_eta_(),
0084       h_L1Puppi_jet_ptratio_vs_pt_barrel_(),
0085       h_L1Puppi_jet_ptratio_vs_pt_endcap_(),
0086       h_L1Puppi_jet_ptratio_vs_pt_ecnotk_(),
0087       h_L1Puppi_jet_ptratio_vs_pt_hf_(),
0088       h_L1Puppi_jet_ptratio_vs_eta_() {
0089   edm::LogInfo("L1TPhase2CorrelatorOffline") << "Constructor "
0090                                              << "L1TPhase2CorrelatorOffline::L1TPhase2CorrelatorOffline " << std::endl;
0091 
0092   auto reconames = objs_.getParameterNamesForType<std::vector<edm::InputTag>>();
0093   for (const std::string& name : reconames) {
0094     reco_.emplace_back(L1TPhase2CorrelatorOffline::MultiCollection(objs_, name, consumesCollector()), RecoVars());
0095   }
0096   for (auto& obj : objs_.getParameter<std::vector<edm::InputTag>>("L1PF")) {
0097     phase2PFToken_.push_back(consumes<std::vector<l1t::PFCandidate>>(obj));
0098   }
0099   for (auto& obj : objs_.getParameter<std::vector<edm::InputTag>>("L1Puppi")) {
0100     phase2PuppiToken_.push_back(consumes<std::vector<l1t::PFCandidate>>(obj));
0101   }
0102 }
0103 
0104 //
0105 // -- Destructor
0106 //
0107 L1TPhase2CorrelatorOffline::~L1TPhase2CorrelatorOffline() {
0108   edm::LogInfo("L1TPhase2CorrelatorOffline")
0109       << "Destructor L1TPhase2CorrelatorOffline::~L1TPhase2CorrelatorOffline " << std::endl;
0110 }
0111 
0112 //
0113 // -------------------------------------- beginRun --------------------------------------------
0114 //
0115 void L1TPhase2CorrelatorOffline::dqmBeginRun(const edm::Run& run, const edm::EventSetup& iSetup) {
0116   edm::LogInfo("L1TPhase2CorrelatorOffline") << "L1TPhase2CorrelatorOffline::beginRun" << std::endl;
0117 
0118   bZ_ = iSetup.getData(BFieldTag_).inTesla(GlobalPoint(0, 0, 0)).z();
0119 }
0120 
0121 //
0122 // -------------------------------------- bookHistos --------------------------------------------
0123 //
0124 void L1TPhase2CorrelatorOffline::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) {
0125   edm::LogInfo("L1TPhase2CorrelatorOffline") << "L1TPhase2CorrelatorOffline::bookHistograms" << std::endl;
0126 
0127   // book at beginRun
0128   bookPhase2CorrelatorHistos(ibooker);
0129 }
0130 
0131 //
0132 // -------------------------------------- Analyze --------------------------------------------
0133 //
0134 void L1TPhase2CorrelatorOffline::analyze(edm::Event const& e, edm::EventSetup const& eSetup) {
0135   edm::Handle<std::vector<reco::GenJet>> genjets;
0136   edm::Handle<std::vector<reco::GenParticle>> genparticles;
0137   e.getByToken(genJetToken_, genjets);
0138   e.getByToken(genParticleToken_, genparticles);
0139 
0140   std::vector<const reco::GenParticle*> prompts, taus;
0141   for (const reco::GenParticle& gen : *genparticles) {
0142     if (isParticleGun_) {
0143       if (gen.statusFlags().isPrompt() == 1)
0144         prompts.push_back(&gen);
0145       continue;
0146     }
0147     if ((gen.isPromptFinalState() || gen.isDirectPromptTauDecayProductFinalState()) &&
0148         (std::abs(gen.pdgId()) == 11 || std::abs(gen.pdgId()) == 13) && gen.pt() > 5) {
0149       prompts.push_back(&gen);
0150     } else if (gen.isPromptFinalState() && std::abs(gen.pdgId()) == 22 && gen.pt() > 10) {
0151       prompts.push_back(&gen);
0152     } else if (abs(gen.pdgId()) == 15 && gen.isPromptDecayed()) {
0153       taus.push_back(&gen);
0154     }
0155   }
0156 
0157   for (auto& recopair : reco_) {
0158     recopair.first.get(e);
0159   }
0160 
0161   for (const reco::GenJet& j1 : *genjets) {
0162     bool ok = true;
0163     const reco::Candidate* match = nullptr;
0164     for (const reco::GenParticle* gp : prompts) {
0165       if (::deltaR2(*gp, j1) < 0.16f) {
0166         if (match != nullptr) {
0167           ok = false;
0168           break;
0169         } else {
0170           match = gp;
0171         }
0172       }
0173     }
0174     if (!ok)
0175       continue;
0176     if (!match) {
0177       // look for a tau
0178       for (const reco::GenParticle* gp : taus) {
0179         if (::deltaR2(*gp, j1) < 0.16f) {
0180           if (match != nullptr) {
0181             ok = false;
0182             break;
0183           } else {
0184             match = gp;
0185           }
0186         }
0187       }
0188       if (!ok)
0189         continue;
0190       if (match != nullptr && match->numberOfDaughters() == 2 &&
0191           std::abs(match->daughter(0)->pdgId()) + std::abs(match->daughter(1)->pdgId()) == 211 + 16) {
0192         // one-prong tau, consider it a pion
0193         match = (std::abs(match->daughter(0)->pdgId()) == 211 ? match->daughter(0) : match->daughter(1));
0194       }
0195     }
0196     if (match != nullptr) {
0197       if (std::abs(match->pdgId()) == 15) {
0198         reco::Particle::LorentzVector pvis;
0199         for (unsigned int i = 0, n = match->numberOfDaughters(); i < n; ++i) {
0200           const reco::Candidate* dau = match->daughter(i);
0201           if (std::abs(dau->pdgId()) == 12 || std::abs(dau->pdgId()) == 14 || std::abs(dau->pdgId()) == 16) {
0202             continue;
0203           }
0204           pvis += dau->p4();
0205         }
0206         mc_.pt = pvis.Pt();
0207         mc_.eta = pvis.Eta();
0208         mc_.phi = pvis.Phi();
0209       } else {
0210         mc_.fillP4(*match);
0211         mc_.fillPropagated(*match, bZ_);
0212       }
0213       mc_.id = std::abs(match->pdgId());
0214       mc_.iso04 = j1.pt() / mc_.pt - 1;
0215       mc_.iso02 = 0;
0216       for (const auto& dptr : j1.daughterPtrVector()) {
0217         if (::deltaR2(*dptr, *match) < 0.04f) {
0218           mc_.iso02 += dptr->pt();
0219         }
0220       }
0221       mc_.iso02 = mc_.iso02 / mc_.pt - 1;
0222     } else {
0223       if (j1.pt() < 20)
0224         continue;
0225       mc_.fillP4(j1);
0226       mc_.id = 0;
0227       mc_.iso02 = 0;
0228       mc_.iso04 = 0;
0229     }
0230     mc_.iso08 = mc_.iso04;
0231     for (const reco::GenJet& j2 : *genjets) {
0232       if (&j2 == &j1)
0233         continue;
0234       if (::deltaR2(j1, j2) < 0.64f)
0235         mc_.iso08 += j2.pt() / mc_.pt;
0236     }
0237     for (auto& recopair : reco_) {
0238       recopair.second.fill(recopair.first.objects(),
0239                            recopair.first.prop() ? mc_.caloeta : mc_.eta,
0240                            recopair.first.prop() ? mc_.calophi : mc_.phi);
0241 
0242       if (abs(mc_.id) != 0 && fabs(mc_.iso04) < 0.05) {
0243         if (recopair.first.name() == "L1PF") {
0244           if (fabs(mc_.eta) < 1.5) {
0245             h_L1PF_part_ptratio_0p2_vs_pt_barrel_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0246           } else if (fabs(mc_.eta) < 2.5) {
0247             h_L1PF_part_ptratio_0p2_vs_pt_endcap_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0248           } else if (fabs(mc_.eta) < 3.) {
0249             h_L1PF_part_ptratio_0p2_vs_pt_ecnotk_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0250           } else if (fabs(mc_.eta) < 5.) {
0251             h_L1PF_part_ptratio_0p2_vs_pt_hf_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0252           }
0253           h_L1PF_part_ptratio_0p2_vs_eta_->Fill(mc_.eta, recopair.second.pt02 / mc_.pt);
0254         }
0255         if (recopair.first.name() == "L1Puppi") {
0256           if (fabs(mc_.eta) < 1.5) {
0257             h_L1Puppi_part_ptratio_0p2_vs_pt_barrel_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0258           } else if (fabs(mc_.eta) < 2.5) {
0259             h_L1Puppi_part_ptratio_0p2_vs_pt_endcap_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0260           } else if (fabs(mc_.eta) < 3.) {
0261             h_L1Puppi_part_ptratio_0p2_vs_pt_ecnotk_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0262           } else if (fabs(mc_.eta) < 5.) {
0263             h_L1Puppi_part_ptratio_0p2_vs_pt_hf_->Fill(mc_.pt, recopair.second.pt02 / mc_.pt);
0264           }
0265           h_L1Puppi_part_ptratio_0p2_vs_eta_->Fill(mc_.eta, recopair.second.pt02 / mc_.pt);
0266         }
0267       }
0268       if (abs(mc_.id) == 0) {
0269         if (recopair.first.name() == "L1PF") {
0270           if (fabs(mc_.eta) < 1.5) {
0271             h_L1PF_jet_ptratio_vs_pt_barrel_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0272           } else if (fabs(mc_.eta) < 2.5) {
0273             h_L1PF_jet_ptratio_vs_pt_endcap_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0274           } else if (fabs(mc_.eta) < 3.) {
0275             h_L1PF_jet_ptratio_vs_pt_ecnotk_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0276           } else if (fabs(mc_.eta) < 5.) {
0277             h_L1PF_jet_ptratio_vs_pt_hf_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0278           }
0279           h_L1PF_jet_ptratio_vs_eta_->Fill(mc_.eta, recopair.second.pt / mc_.pt);
0280         }
0281         if (recopair.first.name() == "L1Puppi") {
0282           if (fabs(mc_.eta) < 1.5) {
0283             h_L1Puppi_jet_ptratio_vs_pt_barrel_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0284           } else if (fabs(mc_.eta) < 2.5) {
0285             h_L1Puppi_jet_ptratio_vs_pt_endcap_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0286           } else if (fabs(mc_.eta) < 3.) {
0287             h_L1Puppi_jet_ptratio_vs_pt_ecnotk_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0288           } else if (fabs(mc_.eta) < 5.) {
0289             h_L1Puppi_jet_ptratio_vs_pt_hf_->Fill(mc_.pt, recopair.second.pt / mc_.pt);
0290           }
0291           h_L1Puppi_jet_ptratio_vs_eta_->Fill(mc_.eta, recopair.second.pt / mc_.pt);
0292         }
0293       }
0294     }
0295   }
0296   for (auto& recopair : reco_) {
0297     recopair.first.clear();
0298     recopair.second.clear();
0299   }
0300 
0301   for (auto& pfToken : phase2PFToken_) {
0302     edm::Handle<std::vector<l1t::PFCandidate>> l1pfs;
0303     e.getByToken(pfToken, l1pfs);
0304 
0305     for (const auto& pfc : *l1pfs) {
0306       h_L1PF_pt_->Fill(pfc.pt());
0307       h_L1PF_eta_->Fill(pfc.eta());
0308       if (abs(pfc.pdgId()) == 13) {
0309         h_L1PF_pt_mu_->Fill(pfc.pt());
0310         h_L1PF_eta_mu_->Fill(pfc.eta());
0311       } else if (abs(pfc.pdgId()) == 11) {
0312         h_L1PF_pt_el_->Fill(pfc.pt());
0313         h_L1PF_eta_el_->Fill(pfc.eta());
0314       } else if (abs(pfc.pdgId()) == 22) {
0315         h_L1PF_pt_pho_->Fill(pfc.pt());
0316         h_L1PF_eta_pho_->Fill(pfc.eta());
0317       } else if (abs(pfc.pdgId()) == 211) {
0318         h_L1PF_pt_ch_->Fill(pfc.pt());
0319         h_L1PF_eta_ch_->Fill(pfc.eta());
0320       } else if (abs(pfc.pdgId()) == 130) {
0321         h_L1PF_pt_nh_->Fill(pfc.pt());
0322         h_L1PF_eta_nh_->Fill(pfc.eta());
0323       }
0324     }  // loop over L1 PF
0325   }
0326   for (auto& pupToken : phase2PuppiToken_) {
0327     edm::Handle<std::vector<l1t::PFCandidate>> l1pups;
0328     e.getByToken(pupToken, l1pups);
0329     for (const auto& pupc : *l1pups) {
0330       h_L1Puppi_pt_->Fill(pupc.pt());
0331       h_L1Puppi_eta_->Fill(pupc.eta());
0332       if (abs(pupc.pdgId()) == 13) {
0333         h_L1Puppi_pt_mu_->Fill(pupc.pt());
0334         h_L1Puppi_eta_mu_->Fill(pupc.eta());
0335       } else if (abs(pupc.pdgId()) == 11) {
0336         h_L1Puppi_pt_el_->Fill(pupc.pt());
0337         h_L1Puppi_eta_el_->Fill(pupc.eta());
0338       } else if (abs(pupc.pdgId()) == 22) {
0339         h_L1Puppi_pt_pho_->Fill(pupc.pt());
0340         h_L1Puppi_eta_pho_->Fill(pupc.eta());
0341       } else if (abs(pupc.pdgId()) == 211) {
0342         h_L1Puppi_pt_ch_->Fill(pupc.pt());
0343         h_L1Puppi_eta_ch_->Fill(pupc.eta());
0344       } else if (abs(pupc.pdgId()) == 130) {
0345         h_L1Puppi_pt_nh_->Fill(pupc.pt());
0346         h_L1Puppi_eta_nh_->Fill(pupc.eta());
0347       }
0348     }  // loop over L1 Puppi
0349   }
0350 }
0351 
0352 //
0353 // -------------------------------------- endRun --------------------------------------------
0354 //
0355 //
0356 // -------------------------------------- book histograms --------------------------------------------
0357 //
0358 void L1TPhase2CorrelatorOffline::bookPhase2CorrelatorHistos(DQMStore::IBooker& ibooker) {
0359   ibooker.cd();
0360   ibooker.setCurrentFolder(histFolder_);
0361   ibooker.setScope(MonitorElementData::Scope::RUN);
0362 
0363   dqmoffline::l1t::HistDefinition resVsPtDef = histDefinitions_[PlotConfig::resVsPt];
0364   dqmoffline::l1t::HistDefinition resVsEtaDef = histDefinitions_[PlotConfig::resVsEta];
0365   dqmoffline::l1t::HistDefinition ptDistDef = histDefinitions_[PlotConfig::ptDist];
0366   dqmoffline::l1t::HistDefinition etaDistDef = histDefinitions_[PlotConfig::etaDist];
0367 
0368   int ptratio_nbins = 300;
0369   float ptratio_lo = 0.;
0370   float ptratio_hi = 3.;
0371 
0372   h_L1PF_pt_ = ibooker.book1D("PF_pt", "L1 PF p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0373   h_L1PF_eta_ = ibooker.book1D("PF_eta", "L1 PF #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0374   h_L1Puppi_pt_ = ibooker.book1D("Puppi_pt", "L1 PUPPI p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0375   h_L1Puppi_eta_ = ibooker.book1D("Puppi_eta", "L1 PUPPI #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0376 
0377   h_L1PF_pt_mu_ = ibooker.book1D("PF_pt_mu", "L1 PF Muon p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0378   h_L1PF_eta_mu_ = ibooker.book1D("PF_eta_mu", "L1 PF Muon #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0379   h_L1Puppi_pt_mu_ =
0380       ibooker.book1D("Puppi_pt_mu", "L1 PUPPI Muon p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0381   h_L1Puppi_eta_mu_ =
0382       ibooker.book1D("Puppi_eta_mu", "L1 PUPPI Muon #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0383 
0384   h_L1PF_pt_el_ = ibooker.book1D("PF_pt_el", "L1 PF Electron p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0385   h_L1PF_eta_el_ =
0386       ibooker.book1D("PF_eta_el", "L1 PF Electron #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0387   h_L1Puppi_pt_el_ =
0388       ibooker.book1D("Puppi_pt_el", "L1 PUPPI Electron p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0389   h_L1Puppi_eta_el_ =
0390       ibooker.book1D("Puppi_eta_el", "L1 PUPPI Electron #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0391 
0392   h_L1PF_pt_pho_ = ibooker.book1D("PF_pt_pho", "L1 PF Photon p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0393   h_L1PF_eta_pho_ =
0394       ibooker.book1D("PF_eta_pho", "L1 PF Photon #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0395   h_L1Puppi_pt_pho_ =
0396       ibooker.book1D("Puppi_pt_pho", "L1 PUPPI Photon p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0397   h_L1Puppi_eta_pho_ =
0398       ibooker.book1D("Puppi_eta_pho", "L1 PUPPI Photon #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0399 
0400   h_L1PF_pt_ch_ =
0401       ibooker.book1D("PF_pt_ch", "L1 PF Charged Hadron p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0402   h_L1PF_eta_ch_ =
0403       ibooker.book1D("PF_eta_ch", "L1 PF Charged Hadron #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0404   h_L1Puppi_pt_ch_ =
0405       ibooker.book1D("Puppi_pt_ch", "L1 PUPPI Charged Hadron p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0406   h_L1Puppi_eta_ch_ = ibooker.book1D(
0407       "Puppi_eta_ch", "L1 PUPPI Charged Hadron #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0408 
0409   h_L1PF_pt_nh_ =
0410       ibooker.book1D("PF_pt_nh", "L1 PF Neutral Hadron p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0411   h_L1PF_eta_nh_ =
0412       ibooker.book1D("PF_eta_nh", "L1 PF Neutral Hadron #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0413   h_L1Puppi_pt_nh_ =
0414       ibooker.book1D("Puppi_pt_nh", "L1 PUPPI Neutral Hadron p_{T}", ptDistDef.nbinsX, ptDistDef.xmin, ptDistDef.xmax);
0415   h_L1Puppi_eta_nh_ = ibooker.book1D(
0416       "Puppi_eta_nh", "L1 PUPPI Neutral Hadron #eta", etaDistDef.nbinsX, etaDistDef.xmin, etaDistDef.xmax);
0417 
0418   ibooker.setCurrentFolder(respresolFolder_);
0419 
0420   h_L1PF_part_ptratio_0p2_vs_pt_barrel_ = ibooker.book2D("L1PFParticlePtRatio0p2VsPtBarrel",
0421                                                          "L1 PF Particle L1/Gen (#Delta R < 0.2) vs p_{T}, Barrel",
0422                                                          resVsPtDef.nbinsX,
0423                                                          resVsPtDef.xmin,
0424                                                          resVsPtDef.xmax,
0425                                                          ptratio_nbins,
0426                                                          ptratio_lo,
0427                                                          ptratio_hi);
0428 
0429   h_L1PF_part_ptratio_0p2_vs_pt_endcap_ = ibooker.book2D("L1PFParticlePtRatio0p2VsPtEndcap",
0430                                                          "L1 PF Particle L1/Gen (#Delta R < 0.2) vs p_{T}, Endcap",
0431                                                          resVsPtDef.nbinsX,
0432                                                          resVsPtDef.xmin,
0433                                                          resVsPtDef.xmax,
0434                                                          ptratio_nbins,
0435                                                          ptratio_lo,
0436                                                          ptratio_hi);
0437 
0438   h_L1PF_part_ptratio_0p2_vs_pt_ecnotk_ =
0439       ibooker.book2D("L1PFParticlePtRatio0p2VsPtEndcapNoTk",
0440                      "L1 PF Particle L1/Gen (#Delta R < 0.2) vs p_{T}, Endcap No Tk",
0441                      resVsPtDef.nbinsX,
0442                      resVsPtDef.xmin,
0443                      resVsPtDef.xmax,
0444                      ptratio_nbins,
0445                      ptratio_lo,
0446                      ptratio_hi);
0447 
0448   h_L1PF_part_ptratio_0p2_vs_pt_hf_ = ibooker.book2D("L1PFParticlePtRatio0p2VsPtHF",
0449                                                      "L1 PF Particle L1/Gen (#Delta R < 0.2) vs p_{T}, HF",
0450                                                      resVsPtDef.nbinsX,
0451                                                      resVsPtDef.xmin,
0452                                                      resVsPtDef.xmax,
0453                                                      ptratio_nbins,
0454                                                      ptratio_lo,
0455                                                      ptratio_hi);
0456 
0457   h_L1PF_part_ptratio_0p2_vs_eta_ = ibooker.book2D("L1PFParticlePtRatio0p2VsEta",
0458                                                    "L1 PF Particle L1/Gen (#Delta R < 0.2) vs #eta",
0459                                                    resVsEtaDef.nbinsX,
0460                                                    resVsEtaDef.xmin,
0461                                                    resVsEtaDef.xmax,
0462                                                    ptratio_nbins,
0463                                                    ptratio_lo,
0464                                                    ptratio_hi);
0465 
0466   h_L1Puppi_part_ptratio_0p2_vs_pt_barrel_ =
0467       ibooker.book2D("L1PUPPIParticlePtRatio0p2VsPtBarrel",
0468                      "L1 PUPPI Particle L1/Gen (#Delta R < 0.2) vs p_{T}, Barrel",
0469                      resVsPtDef.nbinsX,
0470                      resVsPtDef.xmin,
0471                      resVsPtDef.xmax,
0472                      ptratio_nbins,
0473                      ptratio_lo,
0474                      ptratio_hi);
0475 
0476   h_L1Puppi_part_ptratio_0p2_vs_pt_endcap_ =
0477       ibooker.book2D("L1PUPPIParticlePtRatio0p2VsPtEndcap",
0478                      "L1 PUPPI Particle L1/Gen (#Delta R < 0.2) vs p_{T}, Endcap",
0479                      resVsPtDef.nbinsX,
0480                      resVsPtDef.xmin,
0481                      resVsPtDef.xmax,
0482                      ptratio_nbins,
0483                      ptratio_lo,
0484                      ptratio_hi);
0485 
0486   h_L1Puppi_part_ptratio_0p2_vs_pt_ecnotk_ =
0487       ibooker.book2D("L1PUPPIParticlePtRatio0p2VsPtEndcapNoTk",
0488                      "L1 PUPPI Particle L1/Gen (#Delta R < 0.2) vs p_{T}, Endcap No Tk",
0489                      resVsPtDef.nbinsX,
0490                      resVsPtDef.xmin,
0491                      resVsPtDef.xmax,
0492                      ptratio_nbins,
0493                      ptratio_lo,
0494                      ptratio_hi);
0495 
0496   h_L1Puppi_part_ptratio_0p2_vs_pt_hf_ = ibooker.book2D("L1PUPPIParticlePtRatio0p2VsPtHF",
0497                                                         "L1 PUPPI Particle L1/Gen (#Delta R < 0.2) vs p_{T}, HF",
0498                                                         resVsPtDef.nbinsX,
0499                                                         resVsPtDef.xmin,
0500                                                         resVsPtDef.xmax,
0501                                                         ptratio_nbins,
0502                                                         ptratio_lo,
0503                                                         ptratio_hi);
0504 
0505   h_L1Puppi_part_ptratio_0p2_vs_eta_ = ibooker.book2D("L1PUPPIParticlePtRatio0p2VsEta",
0506                                                       "L1 PUPPI Particle L1/Gen (#Delta R < 0.2) vs #eta",
0507                                                       resVsEtaDef.nbinsX,
0508                                                       resVsEtaDef.xmin,
0509                                                       resVsEtaDef.xmax,
0510                                                       ptratio_nbins,
0511                                                       ptratio_lo,
0512                                                       ptratio_hi);
0513 
0514   h_L1PF_jet_ptratio_vs_pt_barrel_ = ibooker.book2D("L1PFJetPtRatioVsPtBarrel",
0515                                                     "L1 PF Jet L1/Gen vs p_{T}, Barrel",
0516                                                     resVsPtDef.nbinsX,
0517                                                     resVsPtDef.xmin,
0518                                                     resVsPtDef.xmax,
0519                                                     ptratio_nbins,
0520                                                     ptratio_lo,
0521                                                     ptratio_hi);
0522 
0523   h_L1PF_jet_ptratio_vs_pt_endcap_ = ibooker.book2D("L1PFJetPtRatioVsPtEndcap",
0524                                                     "L1 PF Jet L1/Gen vs p_{T}, Endcap",
0525                                                     resVsPtDef.nbinsX,
0526                                                     resVsPtDef.xmin,
0527                                                     resVsPtDef.xmax,
0528                                                     ptratio_nbins,
0529                                                     ptratio_lo,
0530                                                     ptratio_hi);
0531 
0532   h_L1PF_jet_ptratio_vs_pt_ecnotk_ = ibooker.book2D("L1PFJetPtRatioVsPtEndcapNoTk",
0533                                                     "L1 PF Jet L1/Gen vs p_{T}, Endcap No Tk",
0534                                                     resVsPtDef.nbinsX,
0535                                                     resVsPtDef.xmin,
0536                                                     resVsPtDef.xmax,
0537                                                     ptratio_nbins,
0538                                                     ptratio_lo,
0539                                                     ptratio_hi);
0540 
0541   h_L1PF_jet_ptratio_vs_pt_hf_ = ibooker.book2D("L1PFJetPtRatioVsPtHF",
0542                                                 "L1 PF Jet L1/Gen vs p_{T}, HF",
0543                                                 resVsPtDef.nbinsX,
0544                                                 resVsPtDef.xmin,
0545                                                 resVsPtDef.xmax,
0546                                                 ptratio_nbins,
0547                                                 ptratio_lo,
0548                                                 ptratio_hi);
0549 
0550   h_L1PF_jet_ptratio_vs_eta_ = ibooker.book2D("L1PFJetPtRatioVsEta",
0551                                               "L1 PF Jet L1/Gen vs #eta",
0552                                               resVsEtaDef.nbinsX,
0553                                               resVsEtaDef.xmin,
0554                                               resVsEtaDef.xmax,
0555                                               ptratio_nbins,
0556                                               ptratio_lo,
0557                                               ptratio_hi);
0558 
0559   h_L1Puppi_jet_ptratio_vs_pt_barrel_ = ibooker.book2D("L1PUPPIJetPtRatioVsPtBarrel",
0560                                                        "L1 PUPPI Jet L1/Gen vs p_{T}, Barrel",
0561                                                        resVsPtDef.nbinsX,
0562                                                        resVsPtDef.xmin,
0563                                                        resVsPtDef.xmax,
0564                                                        ptratio_nbins,
0565                                                        ptratio_lo,
0566                                                        ptratio_hi);
0567 
0568   h_L1Puppi_jet_ptratio_vs_pt_endcap_ = ibooker.book2D("L1PUPPIJetPtRatioVsPtEndcap",
0569                                                        "L1 PUPPI Jet L1/Gen vs p_{T}, Endcap",
0570                                                        resVsPtDef.nbinsX,
0571                                                        resVsPtDef.xmin,
0572                                                        resVsPtDef.xmax,
0573                                                        ptratio_nbins,
0574                                                        ptratio_lo,
0575                                                        ptratio_hi);
0576 
0577   h_L1Puppi_jet_ptratio_vs_pt_ecnotk_ = ibooker.book2D("L1PUPPIJetPtRatioVsPtEndcapNoTk",
0578                                                        "L1 PUPPI Jet L1/Gen vs p_{T}, EndcapNoTk",
0579                                                        resVsPtDef.nbinsX,
0580                                                        resVsPtDef.xmin,
0581                                                        resVsPtDef.xmax,
0582                                                        ptratio_nbins,
0583                                                        ptratio_lo,
0584                                                        ptratio_hi);
0585 
0586   h_L1Puppi_jet_ptratio_vs_pt_hf_ = ibooker.book2D("L1PUPPIJetPtRatioVsPtHF",
0587                                                    "L1 PUPPI Jet L1/Gen vs p_{T}, HF",
0588                                                    resVsPtDef.nbinsX,
0589                                                    resVsPtDef.xmin,
0590                                                    resVsPtDef.xmax,
0591                                                    ptratio_nbins,
0592                                                    ptratio_lo,
0593                                                    ptratio_hi);
0594 
0595   h_L1Puppi_jet_ptratio_vs_eta_ = ibooker.book2D("L1PUPPIJetPtRatioVsEta",
0596                                                  "L1 PUPPI Jet L1/Gen vs #eta",
0597                                                  resVsEtaDef.nbinsX,
0598                                                  resVsEtaDef.xmin,
0599                                                  resVsEtaDef.xmax,
0600                                                  ptratio_nbins,
0601                                                  ptratio_lo,
0602                                                  ptratio_hi);
0603 
0604   ibooker.setCurrentFolder(histFolder_);
0605   //Response
0606 
0607   h_L1PF_part_response_0p2_pt_barrel_ = ibooker.book1D("L1PFParticleResponse0p2VsPtBarrel",
0608                                                        "L1 PF Particle Response (#Delta R < 0.2) vs p_{T}, Barrel",
0609                                                        resVsPtDef.nbinsX,
0610                                                        resVsPtDef.xmin,
0611                                                        resVsPtDef.xmax);
0612 
0613   h_L1PF_part_response_0p2_pt_endcap_ = ibooker.book1D("L1PFParticleResponse0p2VsPtEndcap",
0614                                                        "L1 PF Particle Response (#Delta R < 0.2) vs p_{T}, Endcap",
0615                                                        resVsPtDef.nbinsX,
0616                                                        resVsPtDef.xmin,
0617                                                        resVsPtDef.xmax);
0618 
0619   h_L1PF_part_response_0p2_pt_ecnotk_ =
0620       ibooker.book1D("L1PFParticleResponse0p2VsPtEndcapNoTk",
0621                      "L1 PF Particle Response (#Delta R < 0.2) vs p_{T}, Endcap No Tk",
0622                      resVsPtDef.nbinsX,
0623                      resVsPtDef.xmin,
0624                      resVsPtDef.xmax);
0625 
0626   h_L1PF_part_response_0p2_pt_hf_ = ibooker.book1D("L1PFParticleResponse0p2VsPtHF",
0627                                                    "L1 PF Particle Response (#Delta R < 0.2) vs p_{T}, HF",
0628                                                    resVsPtDef.nbinsX,
0629                                                    resVsPtDef.xmin,
0630                                                    resVsPtDef.xmax);
0631 
0632   h_L1PF_part_response_0p2_eta_ = ibooker.book1D("L1PFParticleResponse0p2VsEta",
0633                                                  "L1 PF Particle Response (#Delta R < 0.2) vs #eta",
0634                                                  resVsEtaDef.nbinsX,
0635                                                  resVsEtaDef.xmin,
0636                                                  resVsEtaDef.xmax);
0637 
0638   h_L1Puppi_part_response_0p2_pt_barrel_ =
0639       ibooker.book1D("L1PUPPIParticleResponse0p2VsPtBarrel",
0640                      "L1 PUPPI Particle Response (#Delta R < 0.2) vs p_{T}, Barrel",
0641                      resVsPtDef.nbinsX,
0642                      resVsPtDef.xmin,
0643                      resVsPtDef.xmax);
0644 
0645   h_L1Puppi_part_response_0p2_pt_endcap_ =
0646       ibooker.book1D("L1PUPPIParticleResponse0p2VsPtEndcap",
0647                      "L1 PUPPI Particle Response (#Delta R < 0.2) vs p_{T}, Endcap",
0648                      resVsPtDef.nbinsX,
0649                      resVsPtDef.xmin,
0650                      resVsPtDef.xmax);
0651 
0652   h_L1Puppi_part_response_0p2_pt_ecnotk_ =
0653       ibooker.book1D("L1PUPPIParticleResponse0p2VsPtEndcapNoTk",
0654                      "L1 PUPPI Particle Response (#Delta R < 0.2) vs p_{T}, Endcap No Tk",
0655                      resVsPtDef.nbinsX,
0656                      resVsPtDef.xmin,
0657                      resVsPtDef.xmax);
0658 
0659   h_L1Puppi_part_response_0p2_pt_hf_ = ibooker.book1D("L1PUPPIParticleResponse0p2VsPtHF",
0660                                                       "L1 PUPPI Particle Response (#Delta R < 0.2) vs p_{T}, HF",
0661                                                       resVsPtDef.nbinsX,
0662                                                       resVsPtDef.xmin,
0663                                                       resVsPtDef.xmax);
0664 
0665   h_L1Puppi_part_response_0p2_eta_ = ibooker.book1D("L1PUPPIParticleResponse0p2VsEta",
0666                                                     "L1 PUPPI Particle Response (#Delta R < 0.2) vs #eta",
0667                                                     resVsEtaDef.nbinsX,
0668                                                     resVsEtaDef.xmin,
0669                                                     resVsEtaDef.xmax);
0670 
0671   h_L1PF_jet_response_pt_barrel_ = ibooker.book1D("L1PFJetResponseVsPtBarrel",
0672                                                   "L1 PF Jet Response vs p_{T}, Barrel",
0673                                                   resVsPtDef.nbinsX,
0674                                                   resVsPtDef.xmin,
0675                                                   resVsPtDef.xmax);
0676 
0677   h_L1PF_jet_response_pt_endcap_ = ibooker.book1D("L1PFJetResponseVsPtEndcap",
0678                                                   "L1 PF Jet Response vs p_{T}, Endcap",
0679                                                   resVsPtDef.nbinsX,
0680                                                   resVsPtDef.xmin,
0681                                                   resVsPtDef.xmax);
0682 
0683   h_L1PF_jet_response_pt_ecnotk_ = ibooker.book1D("L1PFJetResponseVsPtEndcapNoTk",
0684                                                   "L1 PF Jet Response vs p_{T}, Endcap No Tk",
0685                                                   resVsPtDef.nbinsX,
0686                                                   resVsPtDef.xmin,
0687                                                   resVsPtDef.xmax);
0688 
0689   h_L1PF_jet_response_pt_hf_ = ibooker.book1D(
0690       "L1PFJetResponseVsPtHF", "L1 PF Jet Response vs p_{T}, HF", resVsPtDef.nbinsX, resVsPtDef.xmin, resVsPtDef.xmax);
0691 
0692   h_L1PF_jet_response_eta_ = ibooker.book1D(
0693       "L1PFJetResponseVsEta", "L1 PF Jet Response vs #eta", resVsEtaDef.nbinsX, resVsEtaDef.xmin, resVsEtaDef.xmax);
0694 
0695   h_L1Puppi_jet_response_pt_barrel_ = ibooker.book1D("L1PUPPIJetResponseVsPtBarrel",
0696                                                      "L1 PUPPI Jet Response vs p_{T}, Barrel",
0697                                                      resVsPtDef.nbinsX,
0698                                                      resVsPtDef.xmin,
0699                                                      resVsPtDef.xmax);
0700 
0701   h_L1Puppi_jet_response_pt_endcap_ = ibooker.book1D("L1PUPPIJetResponseVsPtEndcap",
0702                                                      "L1 PUPPI Jet Response vs p_{T}, Endcap",
0703                                                      resVsPtDef.nbinsX,
0704                                                      resVsPtDef.xmin,
0705                                                      resVsPtDef.xmax);
0706 
0707   h_L1Puppi_jet_response_pt_ecnotk_ = ibooker.book1D("L1PUPPIJetResponseVsPtEndcapNoTk",
0708                                                      "L1 PUPPI Jet Response vs p_{T}, EndcapNoTk",
0709                                                      resVsPtDef.nbinsX,
0710                                                      resVsPtDef.xmin,
0711                                                      resVsPtDef.xmax);
0712 
0713   h_L1Puppi_jet_response_pt_hf_ = ibooker.book1D("L1PUPPIJetResponseVsPtHF",
0714                                                  "L1 PUPPI Jet Response vs p_{T}, HF",
0715                                                  resVsPtDef.nbinsX,
0716                                                  resVsPtDef.xmin,
0717                                                  resVsPtDef.xmax);
0718 
0719   h_L1Puppi_jet_response_eta_ = ibooker.book1D("L1PUPPIJetResponseVsEta",
0720                                                "L1 PUPPI Jet Response vs #eta",
0721                                                resVsEtaDef.nbinsX,
0722                                                resVsEtaDef.xmin,
0723                                                resVsEtaDef.xmax);
0724 
0725   //Resolution
0726 
0727   h_L1PF_part_resolution_0p2_pt_barrel_ = ibooker.book1D("L1PFParticleResolution0p2VsPtBarrel",
0728                                                          "L1 PF Particle Resolution (#Delta R < 0.2) vs p_{T}, Barrel",
0729                                                          resVsPtDef.nbinsX,
0730                                                          resVsPtDef.xmin,
0731                                                          resVsPtDef.xmax);
0732 
0733   h_L1PF_part_resolution_0p2_pt_endcap_ = ibooker.book1D("L1PFParticleResolution0p2VsPtEndcap",
0734                                                          "L1 PF Particle Resolution (#Delta R < 0.2) vs p_{T}, Endcap",
0735                                                          resVsPtDef.nbinsX,
0736                                                          resVsPtDef.xmin,
0737                                                          resVsPtDef.xmax);
0738 
0739   h_L1PF_part_resolution_0p2_pt_ecnotk_ =
0740       ibooker.book1D("L1PFParticleResolution0p2VsPtEndcapNoTk",
0741                      "L1 PF Particle Resolution (#Delta R < 0.2) vs p_{T}, Endcap No Tk",
0742                      resVsPtDef.nbinsX,
0743                      resVsPtDef.xmin,
0744                      resVsPtDef.xmax);
0745 
0746   h_L1PF_part_resolution_0p2_pt_hf_ = ibooker.book1D("L1PFParticleResolution0p2VsPtHF",
0747                                                      "L1 PF Particle Resolution (#Delta R < 0.2) vs p_{T}, HF",
0748                                                      resVsPtDef.nbinsX,
0749                                                      resVsPtDef.xmin,
0750                                                      resVsPtDef.xmax);
0751 
0752   h_L1Puppi_part_resolution_0p2_pt_barrel_ =
0753       ibooker.book1D("L1PUPPIParticleResolution0p2VsPtBarrel",
0754                      "L1 PUPPI Particle Resolution (#Delta R < 0.2) vs p_{T}, Barrel",
0755                      resVsPtDef.nbinsX,
0756                      resVsPtDef.xmin,
0757                      resVsPtDef.xmax);
0758 
0759   h_L1Puppi_part_resolution_0p2_pt_endcap_ =
0760       ibooker.book1D("L1PUPPIParticleResolution0p2VsPtEndcap",
0761                      "L1 PUPPI Particle Resolution (#Delta R < 0.2) vs p_{T}, Endcap",
0762                      resVsPtDef.nbinsX,
0763                      resVsPtDef.xmin,
0764                      resVsPtDef.xmax);
0765 
0766   h_L1Puppi_part_resolution_0p2_pt_ecnotk_ =
0767       ibooker.book1D("L1PUPPIParticleResolution0p2VsPtEndcapNoTk",
0768                      "L1 PUPPI Particle Resolution (#Delta R < 0.2) vs p_{T}, Endcap No Tk",
0769                      resVsPtDef.nbinsX,
0770                      resVsPtDef.xmin,
0771                      resVsPtDef.xmax);
0772 
0773   h_L1Puppi_part_resolution_0p2_pt_hf_ = ibooker.book1D("L1PUPPIParticleResolution0p2VsPtHF",
0774                                                         "L1 PUPPI Particle Resolution (#Delta R < 0.2) vs p_{T}, HF",
0775                                                         resVsPtDef.nbinsX,
0776                                                         resVsPtDef.xmin,
0777                                                         resVsPtDef.xmax);
0778 
0779   h_L1PF_jet_resolution_pt_barrel_ = ibooker.book1D("L1PFJetResolutionVsPtBarrel",
0780                                                     "L1 PF Jet Resolution vs p_{T}, Barrel",
0781                                                     resVsPtDef.nbinsX,
0782                                                     resVsPtDef.xmin,
0783                                                     resVsPtDef.xmax);
0784 
0785   h_L1PF_jet_resolution_pt_endcap_ = ibooker.book1D("L1PFJetResolutionVsPtEndcap",
0786                                                     "L1 PF Jet Resolution vs p_{T}, Endcap",
0787                                                     resVsPtDef.nbinsX,
0788                                                     resVsPtDef.xmin,
0789                                                     resVsPtDef.xmax);
0790 
0791   h_L1PF_jet_resolution_pt_ecnotk_ = ibooker.book1D("L1PFJetResolutionVsPtEndcapNoTk",
0792                                                     "L1 PF Jet Resolution vs p_{T}, Endcap No Tk",
0793                                                     resVsPtDef.nbinsX,
0794                                                     resVsPtDef.xmin,
0795                                                     resVsPtDef.xmax);
0796 
0797   h_L1PF_jet_resolution_pt_hf_ = ibooker.book1D("L1PFJetResolutionVsPtHF",
0798                                                 "L1 PF Jet Resolution vs p_{T}, HF",
0799                                                 resVsPtDef.nbinsX,
0800                                                 resVsPtDef.xmin,
0801                                                 resVsPtDef.xmax);
0802 
0803   h_L1Puppi_jet_resolution_pt_barrel_ = ibooker.book1D("L1PUPPIJetResolutionVsPtBarrel",
0804                                                        "L1 PUPPI Jet Resolution vs p_{T}, Barrel",
0805                                                        resVsPtDef.nbinsX,
0806                                                        resVsPtDef.xmin,
0807                                                        resVsPtDef.xmax);
0808 
0809   h_L1Puppi_jet_resolution_pt_endcap_ = ibooker.book1D("L1PUPPIJetResolutionVsPtEndcap",
0810                                                        "L1 PUPPI Jet Resolution vs p_{T}, Endcap",
0811                                                        resVsPtDef.nbinsX,
0812                                                        resVsPtDef.xmin,
0813                                                        resVsPtDef.xmax);
0814 
0815   h_L1Puppi_jet_resolution_pt_ecnotk_ = ibooker.book1D("L1PUPPIJetResolutionVsPtEndcapNoTk",
0816                                                        "L1 PUPPI Jet Resolution vs p_{T}, EndcapNoTk",
0817                                                        resVsPtDef.nbinsX,
0818                                                        resVsPtDef.xmin,
0819                                                        resVsPtDef.xmax);
0820 
0821   h_L1Puppi_jet_resolution_pt_hf_ = ibooker.book1D("L1PUPPIJetResolutionVsPtHF",
0822                                                    "L1 PUPPI Jet Resolution vs p_{T}, HF",
0823                                                    resVsPtDef.nbinsX,
0824                                                    resVsPtDef.xmin,
0825                                                    resVsPtDef.xmax);
0826 
0827   ibooker.cd();
0828 
0829   return;
0830 }
0831 
0832 //
0833 // -------------------------------------- endRun --------------------------------------------
0834 //
0835 void L1TPhase2CorrelatorOffline::dqmEndRun(const edm::Run& run, const edm::EventSetup& iSetup) {
0836   computeResponseResolution();
0837 }
0838 
0839 void L1TPhase2CorrelatorOffline::computeResponseResolution() {
0840   std::vector<MonitorElement*> monElementstoComputeIn = {h_L1PF_part_ptratio_0p2_vs_pt_barrel_,
0841                                                          h_L1PF_part_ptratio_0p2_vs_pt_endcap_,
0842                                                          h_L1PF_part_ptratio_0p2_vs_pt_ecnotk_,
0843                                                          h_L1PF_part_ptratio_0p2_vs_pt_hf_,
0844                                                          h_L1PF_part_ptratio_0p2_vs_eta_,
0845                                                          h_L1Puppi_part_ptratio_0p2_vs_pt_barrel_,
0846                                                          h_L1Puppi_part_ptratio_0p2_vs_pt_endcap_,
0847                                                          h_L1Puppi_part_ptratio_0p2_vs_pt_ecnotk_,
0848                                                          h_L1Puppi_part_ptratio_0p2_vs_pt_hf_,
0849                                                          h_L1Puppi_part_ptratio_0p2_vs_eta_,
0850                                                          h_L1PF_jet_ptratio_vs_pt_barrel_,
0851                                                          h_L1PF_jet_ptratio_vs_pt_endcap_,
0852                                                          h_L1PF_jet_ptratio_vs_pt_ecnotk_,
0853                                                          h_L1PF_jet_ptratio_vs_pt_hf_,
0854                                                          h_L1PF_jet_ptratio_vs_eta_,
0855                                                          h_L1Puppi_jet_ptratio_vs_pt_barrel_,
0856                                                          h_L1Puppi_jet_ptratio_vs_pt_endcap_,
0857                                                          h_L1Puppi_jet_ptratio_vs_pt_ecnotk_,
0858                                                          h_L1Puppi_jet_ptratio_vs_pt_hf_,
0859                                                          h_L1Puppi_jet_ptratio_vs_eta_};
0860   std::vector<MonitorElement*> monElementstoComputeResp = {h_L1PF_part_response_0p2_pt_barrel_,
0861                                                            h_L1PF_part_response_0p2_pt_endcap_,
0862                                                            h_L1PF_part_response_0p2_pt_ecnotk_,
0863                                                            h_L1PF_part_response_0p2_pt_hf_,
0864                                                            h_L1PF_part_response_0p2_eta_,
0865                                                            h_L1Puppi_part_response_0p2_pt_barrel_,
0866                                                            h_L1Puppi_part_response_0p2_pt_endcap_,
0867                                                            h_L1Puppi_part_response_0p2_pt_ecnotk_,
0868                                                            h_L1Puppi_part_response_0p2_pt_hf_,
0869                                                            h_L1Puppi_part_response_0p2_eta_,
0870                                                            h_L1PF_jet_response_pt_barrel_,
0871                                                            h_L1PF_jet_response_pt_endcap_,
0872                                                            h_L1PF_jet_response_pt_ecnotk_,
0873                                                            h_L1PF_jet_response_pt_hf_,
0874                                                            h_L1PF_jet_response_eta_,
0875                                                            h_L1Puppi_jet_response_pt_barrel_,
0876                                                            h_L1Puppi_jet_response_pt_endcap_,
0877                                                            h_L1Puppi_jet_response_pt_ecnotk_,
0878                                                            h_L1Puppi_jet_response_pt_hf_,
0879                                                            h_L1Puppi_jet_response_eta_};
0880   std::vector<MonitorElement*> monElementstoComputeResol = {h_L1PF_part_resolution_0p2_pt_barrel_,
0881                                                             h_L1PF_part_resolution_0p2_pt_endcap_,
0882                                                             h_L1PF_part_resolution_0p2_pt_ecnotk_,
0883                                                             h_L1PF_part_resolution_0p2_pt_hf_,
0884                                                             nullptr,
0885                                                             h_L1Puppi_part_resolution_0p2_pt_barrel_,
0886                                                             h_L1Puppi_part_resolution_0p2_pt_endcap_,
0887                                                             h_L1Puppi_part_resolution_0p2_pt_ecnotk_,
0888                                                             h_L1Puppi_part_resolution_0p2_pt_hf_,
0889                                                             nullptr,
0890                                                             h_L1PF_jet_resolution_pt_barrel_,
0891                                                             h_L1PF_jet_resolution_pt_endcap_,
0892                                                             h_L1PF_jet_resolution_pt_ecnotk_,
0893                                                             h_L1PF_jet_resolution_pt_hf_,
0894                                                             nullptr,
0895                                                             h_L1Puppi_jet_resolution_pt_barrel_,
0896                                                             h_L1Puppi_jet_resolution_pt_endcap_,
0897                                                             h_L1Puppi_jet_resolution_pt_ecnotk_,
0898                                                             h_L1Puppi_jet_resolution_pt_hf_,
0899                                                             nullptr};
0900 
0901   for (unsigned int i = 0; i < monElementstoComputeIn.size(); i++) {
0902     if (monElementstoComputeIn[i] != nullptr && monElementstoComputeResp[i] != nullptr) {
0903       medianResponseCorrResolution(
0904           monElementstoComputeIn[i], monElementstoComputeResp[i], monElementstoComputeResol[i]);
0905     }
0906   }
0907 }
0908 
0909 std::vector<float> L1TPhase2CorrelatorOffline::getQuantile(float quant, TH2F* hist) {
0910   std::vector<float> quantiles(hist->GetNbinsX(), 1.);
0911   for (int ix = 1; ix < hist->GetNbinsX() + 1; ix++) {
0912     float thresh = quant * (hist->Integral(ix, ix, 0, -1));
0913     if (hist->Integral(ix, ix, 0, -1) == 0.) {
0914     } else if (quant <= 0. || thresh < hist->GetBinContent(ix, 0)) {
0915       quantiles[ix - 1] = hist->GetYaxis()->GetBinLowEdge(1);
0916     } else if (quant >= 1. || thresh >= hist->Integral(ix, ix, 0, hist->GetNbinsY())) {
0917       quantiles[ix - 1] = hist->GetYaxis()->GetBinUpEdge(hist->GetNbinsY());
0918     } else {
0919       float sum = hist->GetBinContent(ix, 0);
0920       for (int iy = 1; iy < hist->GetNbinsY() + 1; iy++) {
0921         float add = hist->GetBinContent(ix, iy);
0922         if (sum + add >= thresh) {
0923           quantiles[ix - 1] =
0924               hist->GetYaxis()->GetBinLowEdge(iy) + hist->GetYaxis()->GetBinWidth(iy) * ((thresh - sum) / add);
0925           break;
0926         }
0927         sum += add;
0928       }
0929     }
0930   }
0931   return quantiles;
0932 }
0933 
0934 void L1TPhase2CorrelatorOffline::medianResponseCorrResolution(MonitorElement* in2D,
0935                                                               MonitorElement* response,
0936                                                               MonitorElement* resolution) {
0937   auto hbase = in2D->getTH2F();
0938   auto hresp = response->getTH1F();
0939   if (hbase != nullptr && hresp != nullptr) {
0940     if (hbase->GetNbinsX() == hresp->GetNbinsX()) {
0941       auto med = getQuantile(0.5, hbase);
0942       TGraph* ptrecgen = new TGraph(hbase->GetNbinsX());
0943       for (int ib = 1; ib < hbase->GetNbinsX() + 1; ib++) {
0944         float corr = med[ib - 1];
0945         float xval = hbase->GetXaxis()->GetBinCenter(ib);
0946         ptrecgen->SetPoint(ib - 1, xval * corr, xval);
0947         hresp->SetBinContent(ib, corr);
0948       }
0949       if (resolution != nullptr) {
0950         auto hresol = resolution->getTH1F();
0951         if (hresol != nullptr) {
0952           if (hbase->GetNbinsX() == hresol->GetNbinsX()) {
0953             ptrecgen->Sort();
0954             TH2F* ch = new TH2F(*hbase);
0955             ch->Reset("ICE");
0956             for (int ibx = 1; ibx < ch->GetNbinsX() + 1; ibx++) {
0957               float xval = hbase->GetXaxis()->GetBinCenter(ibx);
0958               for (int iby = 1; iby < ch->GetNbinsY() + 1; iby++) {
0959                 float yval = hbase->GetYaxis()->GetBinCenter(iby);
0960                 float newyval = ptrecgen->Eval(yval * xval) / xval;
0961                 int ycb = ch->FindBin(xval, newyval);
0962                 ch->SetBinContent(ycb, ch->GetBinContent(ycb) + hbase->GetBinContent(ibx, iby));
0963               }
0964             }
0965             delete ptrecgen;
0966             auto qc = getQuantile(0.5, ch);
0967             auto qhi = getQuantile(0.84, ch);
0968             auto qlo = getQuantile(0.16, ch);
0969             delete ch;
0970             for (int ibx = 1; ibx < hbase->GetNbinsX() + 1; ibx++) {
0971               hresol->SetBinContent(ibx, qc[ibx - 1] > 0.2 ? (qhi[ibx - 1] - qlo[ibx - 1]) / 2. : 0.);
0972             }
0973           }
0974         }
0975       } else {
0976         delete ptrecgen;
0977       }
0978     }
0979   }
0980 }
0981 
0982 // define this as a plug-in
0983 DEFINE_FWK_MODULE(L1TPhase2CorrelatorOffline);