Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:33:07

0001 // -*- C++ -*-
0002 //
0003 // Package:    L1Trigger/L1CaloTrigger
0004 // Class:      Phase2L1CaloJetEmulator
0005 //
0006 /**\class Phase2L1CaloJetEmulator Phase2L1CaloJetEmulator.cc L1Trigger/L1CaloTrigger/plugins/Phase2L1CaloJetEmulator.cc
0007 
0008  Description: Producing GCT calo jets using GCT barrel, HGCal and HF towers, based on firmware logic.
0009 
0010  Implementation:
0011      Depends on producers for CaloTowerCollection, HGCalTowerBxCollection and HcalTrigPrimDigiCollection.
0012 */
0013 //
0014 // Original Author:  Pallabi Das
0015 //         Created:  Tue, 11 Apr 2023 11:27:33 GMT
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/StreamID.h"
0031 
0032 #include "DataFormats/L1TCalorimeterPhase2/interface/CaloCrystalCluster.h"
0033 #include "DataFormats/L1TCalorimeterPhase2/interface/CaloTower.h"
0034 #include "DataFormats/L1TCalorimeterPhase2/interface/CaloPFCluster.h"
0035 #include "DataFormats/L1TCalorimeterPhase2/interface/Phase2L1CaloJet.h"
0036 #include "DataFormats/L1Trigger/interface/EGamma.h"
0037 #include "DataFormats/L1THGCal/interface/HGCalTower.h"
0038 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
0039 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
0040 #include "CalibFormats/CaloTPG/interface/CaloTPGTranscoder.h"
0041 #include "CalibFormats/CaloTPG/interface/CaloTPGRecord.h"
0042 #include "L1Trigger/L1TCalorimeter/interface/CaloTools.h"
0043 
0044 #include "L1Trigger/L1CaloTrigger/interface/Phase2L1CaloJetEmulator.h"
0045 #include <ap_int.h>
0046 #include <cstdio>
0047 #include <fstream>
0048 #include <iomanip>
0049 #include <iostream>
0050 #include "TF1.h"
0051 
0052 //
0053 // class declaration
0054 //
0055 
0056 class Phase2L1CaloJetEmulator : public edm::stream::EDProducer<> {
0057 public:
0058   explicit Phase2L1CaloJetEmulator(const edm::ParameterSet&);
0059   ~Phase2L1CaloJetEmulator() override;
0060 
0061   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0062 
0063 private:
0064   void produce(edm::Event&, const edm::EventSetup&) override;
0065   float get_jet_pt_calibration(const float& jet_pt, const float& jet_eta) const;
0066   float get_tau_pt_calibration(const float& tau_pt, const float& tau_eta) const;
0067 
0068   // ----------member data ---------------------------
0069   edm::EDGetTokenT<l1tp2::CaloTowerCollection> caloTowerToken_;
0070   edm::EDGetTokenT<l1t::HGCalTowerBxCollection> hgcalTowerToken_;
0071   edm::EDGetTokenT<HcalTrigPrimDigiCollection> hfToken_;
0072   edm::ESGetToken<CaloTPGTranscoder, CaloTPGRecord> decoderTag_;
0073   std::vector<edm::ParameterSet> nHits_to_nvtx_params;
0074   std::vector<edm::ParameterSet> nvtx_to_PU_sub_params;
0075   std::map<std::string, TF1> nHits_to_nvtx_funcs;
0076   std::map<std::string, TF1> hgcalEM_nvtx_to_PU_sub_funcs;
0077   std::map<std::string, TF1> hgcalHad_nvtx_to_PU_sub_funcs;
0078   std::map<std::string, TF1> hf_nvtx_to_PU_sub_funcs;
0079   std::map<std::string, std::map<std::string, TF1>> all_nvtx_to_PU_sub_funcs;
0080 
0081   // For fetching jet pt calibrations
0082   std::vector<double> jetPtBins;
0083   std::vector<double> absEtaBinsBarrel;
0084   std::vector<double> jetCalibrationsBarrel;
0085   std::vector<double> absEtaBinsHGCal;
0086   std::vector<double> jetCalibrationsHGCal;
0087   std::vector<double> absEtaBinsHF;
0088   std::vector<double> jetCalibrationsHF;
0089 
0090   // For fetching tau pt calibrations
0091   std::vector<double> tauPtBins;
0092   std::vector<double> tauAbsEtaBinsBarrel;
0093   std::vector<double> tauCalibrationsBarrel;
0094   std::vector<double> tauAbsEtaBinsHGCal;
0095   std::vector<double> tauCalibrationsHGCal;
0096 
0097   // For storing jet calibrations
0098   std::vector<std::vector<double>> calibrationsBarrel;
0099   std::vector<std::vector<double>> calibrationsHGCal;
0100   std::vector<std::vector<double>> calibrationsHF;
0101 
0102   // For storing tau calibrations
0103   std::vector<std::vector<double>> tauPtCalibrationsBarrel;
0104   std::vector<std::vector<double>> tauPtCalibrationsHGCal;
0105 };
0106 
0107 //
0108 // constructors and destructor
0109 //
0110 Phase2L1CaloJetEmulator::Phase2L1CaloJetEmulator(const edm::ParameterSet& iConfig)
0111     : caloTowerToken_(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("gctFullTowers"))),
0112       hgcalTowerToken_(consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("hgcalTowers"))),
0113       hfToken_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("hcalDigis"))),
0114       decoderTag_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>(edm::ESInputTag("", ""))),
0115       nHits_to_nvtx_params(iConfig.getParameter<std::vector<edm::ParameterSet>>("nHits_to_nvtx_params")),
0116       nvtx_to_PU_sub_params(iConfig.getParameter<std::vector<edm::ParameterSet>>("nvtx_to_PU_sub_params")),
0117       jetPtBins(iConfig.getParameter<std::vector<double>>("jetPtBins")),
0118       absEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("absEtaBinsBarrel")),
0119       jetCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("jetCalibrationsBarrel")),
0120       absEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("absEtaBinsHGCal")),
0121       jetCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("jetCalibrationsHGCal")),
0122       absEtaBinsHF(iConfig.getParameter<std::vector<double>>("absEtaBinsHF")),
0123       jetCalibrationsHF(iConfig.getParameter<std::vector<double>>("jetCalibrationsHF")),
0124       tauPtBins(iConfig.getParameter<std::vector<double>>("tauPtBins")),
0125       tauAbsEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsBarrel")),
0126       tauCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("tauCalibrationsBarrel")),
0127       tauAbsEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsHGCal")),
0128       tauCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("tauCalibrationsHGCal")) {
0129   for (uint i = 0; i < nHits_to_nvtx_params.size(); i++) {
0130     edm::ParameterSet* pset = &nHits_to_nvtx_params.at(i);
0131     std::string calo = pset->getParameter<std::string>("fit");
0132     nHits_to_nvtx_funcs[calo.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
0133     nHits_to_nvtx_funcs[calo.c_str()].SetParameter(0, pset->getParameter<std::vector<double>>("nHits_params").at(0));
0134     nHits_to_nvtx_funcs[calo.c_str()].SetParameter(1, pset->getParameter<std::vector<double>>("nHits_params").at(1));
0135   }
0136   all_nvtx_to_PU_sub_funcs["hgcalEM"] = hgcalEM_nvtx_to_PU_sub_funcs;
0137   all_nvtx_to_PU_sub_funcs["hgcalHad"] = hgcalHad_nvtx_to_PU_sub_funcs;
0138   all_nvtx_to_PU_sub_funcs["hf"] = hf_nvtx_to_PU_sub_funcs;
0139   for (uint i = 0; i < nvtx_to_PU_sub_params.size(); i++) {
0140     edm::ParameterSet* pset = &nvtx_to_PU_sub_params.at(i);
0141     std::string calo = pset->getParameter<std::string>("calo");
0142     std::string iEta = pset->getParameter<std::string>("iEta");
0143     double p1 = pset->getParameter<std::vector<double>>("nvtx_params").at(0);
0144     double p2 = pset->getParameter<std::vector<double>>("nvtx_params").at(1);
0145 
0146     all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
0147     all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(0, p1);
0148     all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(1, p2);
0149   }
0150 
0151   // Fill the jet pt calibration 2D vector
0152   // Dimension 1 is AbsEta bin
0153   // Dimension 2 is jet pT bin which is filled with the actual callibration value
0154   // size()-1 b/c the inputs have lower and upper bounds
0155   // Do Barrel, then HGCal, then HF
0156   int index = 0;
0157   for (unsigned int abs_eta = 0; abs_eta < absEtaBinsBarrel.size() - 1; abs_eta++) {
0158     std::vector<double> pt_bin_calibs;
0159     for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
0160       pt_bin_calibs.push_back(jetCalibrationsBarrel.at(index));
0161       index++;
0162     }
0163     calibrationsBarrel.push_back(pt_bin_calibs);
0164   }
0165 
0166   index = 0;
0167   for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHGCal.size() - 1; abs_eta++) {
0168     std::vector<double> pt_bin_calibs;
0169     for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
0170       pt_bin_calibs.push_back(jetCalibrationsHGCal.at(index));
0171       index++;
0172     }
0173     calibrationsHGCal.push_back(pt_bin_calibs);
0174   }
0175 
0176   index = 0;
0177   for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHF.size() - 1; abs_eta++) {
0178     std::vector<double> pt_bin_calibs;
0179     for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
0180       pt_bin_calibs.push_back(jetCalibrationsHF.at(index));
0181       index++;
0182     }
0183     calibrationsHF.push_back(pt_bin_calibs);
0184   }
0185 
0186   // Fill the tau pt calibration 2D vector
0187   // Dimension 1 is AbsEta bin
0188   // Dimension 2 is tau pT bin which is filled with the actual calibration value
0189   // Do Barrel, then HGCal
0190   //
0191   // Note to future developers: be very concious of the order in which the calibrations are printed
0192   // out in tool which makse the cfg files.  You need to match that exactly when loading them and
0193   // using the calibrations below.
0194   index = 0;
0195   for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsBarrel.size() - 1; abs_eta++) {
0196     std::vector<double> pt_bin_calibs;
0197     for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
0198       pt_bin_calibs.push_back(tauCalibrationsBarrel.at(index));
0199       index++;
0200     }
0201     tauPtCalibrationsBarrel.push_back(pt_bin_calibs);
0202   }
0203 
0204   index = 0;
0205   for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsHGCal.size() - 1; abs_eta++) {
0206     std::vector<double> pt_bin_calibs;
0207     for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
0208       pt_bin_calibs.push_back(tauCalibrationsHGCal.at(index));
0209       index++;
0210     }
0211     tauPtCalibrationsHGCal.push_back(pt_bin_calibs);
0212   }
0213 
0214   produces<l1tp2::Phase2L1CaloJetCollection>("GCTJet");
0215 }
0216 
0217 Phase2L1CaloJetEmulator::~Phase2L1CaloJetEmulator() {}
0218 
0219 //
0220 // member functions
0221 //
0222 
0223 // ------------ method called to produce the data  ------------
0224 void Phase2L1CaloJetEmulator::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0225   using namespace edm;
0226   std::unique_ptr<l1tp2::Phase2L1CaloJetCollection> jetCands(make_unique<l1tp2::Phase2L1CaloJetCollection>());
0227 
0228   // Assign ETs to each eta-half of the barrel region (17x72 --> 18x72 to be able to make 3x3 super towers)
0229   edm::Handle<std::vector<l1tp2::CaloTower>> caloTowerCollection;
0230   if (!iEvent.getByToken(caloTowerToken_, caloTowerCollection))
0231     edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from caloTowerCollection!";
0232 
0233   iEvent.getByToken(caloTowerToken_, caloTowerCollection);
0234   float GCTintTowers[nBarrelEta][nBarrelPhi];
0235   float realEta[nBarrelEta][nBarrelPhi];
0236   float realPhi[nBarrelEta][nBarrelPhi];
0237   for (const l1tp2::CaloTower& i : *caloTowerCollection) {
0238     int ieta = i.towerIEta();
0239     int iphi = i.towerIPhi();
0240     if (i.ecalTowerEt() > 1.)
0241       GCTintTowers[ieta][iphi] = i.ecalTowerEt();  // suppress <= 1 GeV towers
0242     else
0243       GCTintTowers[ieta][iphi] = 0;
0244     realEta[ieta][iphi] = i.towerEta();
0245     realPhi[ieta][iphi] = i.towerPhi();
0246   }
0247 
0248   float temporary[nBarrelEta / 2][nBarrelPhi];
0249 
0250   // HGCal and HF info used for nvtx estimation
0251   edm::Handle<l1t::HGCalTowerBxCollection> hgcalTowerCollection;
0252   if (!iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection))
0253     edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from hgcalTowerCollection!";
0254   l1t::HGCalTowerBxCollection hgcalTowerColl;
0255   iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection);
0256   hgcalTowerColl = (*hgcalTowerCollection.product());
0257 
0258   edm::Handle<HcalTrigPrimDigiCollection> hfHandle;
0259   if (!iEvent.getByToken(hfToken_, hfHandle))
0260     edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get HcalTrigPrimDigi for HF!";
0261   iEvent.getByToken(hfToken_, hfHandle);
0262 
0263   int i_hgcalEM_hits_leq_threshold = 0;
0264   int i_hgcalHad_hits_leq_threshold = 0;
0265   int i_hf_hits_leq_threshold = 0;
0266   for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
0267     if (it->etEm() <= 1.75 && it->etEm() >= 1.25) {
0268       i_hgcalEM_hits_leq_threshold++;
0269     }
0270     if (it->etHad() <= 1.25 && it->etHad() >= 0.75) {
0271       i_hgcalHad_hits_leq_threshold++;
0272     }
0273   }
0274   const auto& decoder = iSetup.getData(decoderTag_);
0275   for (const auto& hit : *hfHandle.product()) {
0276     double et = decoder.hcaletValue(hit.id(), hit.t0());
0277     if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
0278       continue;
0279     if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
0280       continue;
0281     if (et <= 15.0 && et >= 10.0)
0282       i_hf_hits_leq_threshold++;
0283   }
0284 
0285   double hgcalEM_nvtx = nHits_to_nvtx_funcs["hgcalEM"].Eval(i_hgcalEM_hits_leq_threshold);
0286   if (hgcalEM_nvtx < 0)
0287     hgcalEM_nvtx = 0;
0288   double hgcalHad_nvtx = nHits_to_nvtx_funcs["hgcalHad"].Eval(i_hgcalHad_hits_leq_threshold);
0289   if (hgcalHad_nvtx < 0)
0290     hgcalHad_nvtx = 0;
0291   double hf_nvtx = nHits_to_nvtx_funcs["hf"].Eval(i_hf_hits_leq_threshold);
0292   if (hf_nvtx < 0)
0293     hf_nvtx = 0;
0294   double EstimatedNvtx = (hgcalEM_nvtx + hgcalHad_nvtx + hf_nvtx) / 3.;
0295 
0296   // Assign ETs to each eta-half of the endcap region (18x72)
0297   float hgcalTowers[nHgcalEta][nHgcalPhi];
0298   float hgcalEta[nHgcalEta][nHgcalPhi];
0299   float hgcalPhi[nHgcalEta][nHgcalPhi];
0300 
0301   for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
0302     for (int ieta = 0; ieta < nHgcalEta; ieta++) {
0303       hgcalTowers[ieta][iphi] = 0;
0304       if (ieta < nHgcalEta / 2)
0305         hgcalEta[ieta][iphi] = -3.045 + ieta * 0.087 + 0.0435;
0306       else
0307         hgcalEta[ieta][iphi] = 1.479 + (ieta - nHgcalEta / 2) * 0.087 + 0.0435;
0308       hgcalPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHgcalPhi) + (M_PI / nHgcalPhi);
0309     }
0310   }
0311 
0312   for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
0313     float eta = it->eta();
0314     int ieta;
0315     if (eta < 0)
0316       ieta = 19 - it->id().iEta();
0317     else
0318       ieta = 20 + it->id().iEta();
0319     if (eta > 1.479)
0320       ieta = ieta - 4;
0321     int iphi = it->id().iPhi();
0322 
0323     float hgcal_etEm = it->etEm();
0324     float hgcal_etHad = it->etHad();
0325     std::string etaKey = "";
0326     if (abs(eta) <= 1.8)
0327       etaKey = "er1p4to1p8";
0328     else if (abs(eta) <= 2.1 && abs(eta) > 1.8)
0329       etaKey = "er1p8to2p1";
0330     else if (abs(eta) <= 2.4 && abs(eta) > 2.1)
0331       etaKey = "er2p1to2p4";
0332     else if (abs(eta) <= 2.7 && abs(eta) > 2.4)
0333       etaKey = "er2p4to2p7";
0334     else if (abs(eta) <= 3.1 && abs(eta) > 2.7)
0335       etaKey = "er2p7to3p1";
0336     if (!etaKey.empty()) {
0337       hgcal_etEm = it->etEm() - all_nvtx_to_PU_sub_funcs["hgcalEM"][etaKey].Eval(EstimatedNvtx);
0338       hgcal_etHad = it->etHad() - all_nvtx_to_PU_sub_funcs["hgcalHad"][etaKey].Eval(EstimatedNvtx);
0339     }
0340 
0341     if (hgcal_etEm < 0)
0342       hgcal_etEm = 0;
0343     if (hgcal_etHad < 0)
0344       hgcal_etHad = 0;
0345     if ((it->etEm() + it->etHad() > 1.) && abs(eta) > 1.479)
0346       hgcalTowers[ieta][iphi] = hgcal_etEm + hgcal_etHad;  // suppress <= 1 GeV towers
0347   }
0348 
0349   float temporary_hgcal[nHgcalEta / 2][nHgcalPhi];
0350 
0351   // Assign ETs to each eta-half of the forward region (12x72)
0352   float hfTowers[2 * nHfEta][nHfPhi];  // split 12 -> 24
0353   float hfEta[nHfEta][nHfPhi];
0354   float hfPhi[nHfEta][nHfPhi];
0355   for (int iphi = 0; iphi < nHfPhi; iphi++) {
0356     for (int ieta = 0; ieta < nHfEta; ieta++) {
0357       hfTowers[2 * ieta][iphi] = 0;
0358       hfTowers[2 * ieta + 1][iphi] = 0;
0359       int temp;
0360       if (ieta < nHfEta / 2)
0361         temp = ieta - l1t::CaloTools::kHFEnd;
0362       else
0363         temp = ieta - nHfEta / 2 + l1t::CaloTools::kHFBegin + 1;
0364       hfEta[ieta][iphi] = l1t::CaloTools::towerEta(temp);
0365       hfPhi[ieta][iphi] = l1t::CaloTools::towerPhi(temp, iphi + 1);
0366     }
0367   }
0368 
0369   for (const auto& hit : *hfHandle.product()) {
0370     double et = decoder.hcaletValue(hit.id(), hit.t0());
0371     int ieta = 0;
0372     if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
0373       continue;
0374     if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
0375       continue;
0376     if (hit.id().ieta() <= -(l1t::CaloTools::kHFBegin + 1)) {
0377       ieta = hit.id().ieta() + l1t::CaloTools::kHFEnd;
0378     } else if (hit.id().ieta() >= (l1t::CaloTools::kHFBegin + 1)) {
0379       ieta = nHfEta / 2 + (hit.id().ieta() - (l1t::CaloTools::kHFBegin + 1));
0380     }
0381     int iphi = hit.id().iphi() - 1;  // HF phi runs between 1-72
0382     if (abs(hit.id().ieta()) <= 33 && abs(hit.id().ieta()) >= 29)
0383       et = et - all_nvtx_to_PU_sub_funcs["hf"]["er29to33"].Eval(EstimatedNvtx);
0384     if (abs(hit.id().ieta()) <= 37 && abs(hit.id().ieta()) >= 34)
0385       et = et - all_nvtx_to_PU_sub_funcs["hf"]["er34to37"].Eval(EstimatedNvtx);
0386     if (abs(hit.id().ieta()) <= 41 && abs(hit.id().ieta()) >= 38)
0387       et = et - all_nvtx_to_PU_sub_funcs["hf"]["er38to41"].Eval(EstimatedNvtx);
0388     if (et < 0)
0389       et = 0;
0390     // split tower energy
0391     hfTowers[2 * ieta][iphi] = et / 2;
0392     hfTowers[2 * ieta + 1][iphi] = et / 2;
0393     if ((ieta < 2 || ieta >= nHfEta - 2) && iphi % 4 == 2) {
0394       hfTowers[2 * ieta][iphi] = et / 8;
0395       hfTowers[2 * ieta + 1][iphi] = et / 8;
0396       hfTowers[2 * ieta][iphi + 1] = et / 8;
0397       hfTowers[2 * ieta + 1][iphi + 1] = et / 8;
0398       if (iphi + 2 == nHfPhi) {
0399         hfTowers[2 * ieta][0] = et / 8;
0400         hfTowers[2 * ieta + 1][0] = et / 8;
0401         hfTowers[2 * ieta][1] = et / 8;
0402         hfTowers[2 * ieta + 1][1] = et / 8;
0403       } else {
0404         hfTowers[2 * ieta][iphi + 2] = et / 8;
0405         hfTowers[2 * ieta + 1][iphi + 2] = et / 8;
0406         hfTowers[2 * ieta][iphi + 3] = et / 8;
0407         hfTowers[2 * ieta + 1][iphi + 3] = et / 8;
0408       }
0409     } else if ((ieta >= 2 && ieta < nHfEta - 2) && iphi % 2 == 0) {
0410       hfTowers[2 * ieta][iphi] = et / 4;
0411       hfTowers[2 * ieta + 1][iphi] = et / 4;
0412       hfTowers[2 * ieta][iphi + 1] = et / 4;
0413       hfTowers[2 * ieta + 1][iphi + 1] = et / 4;
0414     }
0415   }
0416 
0417   float temporary_hf[nHfEta][nHfPhi];
0418 
0419   // Begin creating jets
0420   // First create 3x3 super towers: 6x24 in barrel, endcap; 4x24 in forward
0421   // Then create up to 10 jets in each eta half of barrel, endcap, forward regions
0422 
0423   vector<l1tp2::Phase2L1CaloJet> halfBarrelJets, halfHgcalJets, halfHfJets;
0424   halfBarrelJets.clear();
0425   halfHgcalJets.clear();
0426   halfHfJets.clear();
0427   vector<l1tp2::Phase2L1CaloJet> allJets;
0428   allJets.clear();
0429 
0430   for (int k = 0; k < 2; k++) {
0431     halfBarrelJets.clear();
0432     halfHgcalJets.clear();
0433     halfHfJets.clear();
0434     gctobj::jetInfo jet[3 * nJets];
0435 
0436     // BARREL
0437     for (int iphi = 0; iphi < nBarrelPhi; iphi++) {
0438       for (int ieta = 0; ieta < nBarrelEta / 2; ieta++) {
0439         if (k == 0)
0440           temporary[ieta][iphi] = GCTintTowers[ieta][iphi];
0441         else
0442           temporary[ieta][iphi] = GCTintTowers[nBarrelEta / 2 + ieta][iphi];
0443       }
0444     }
0445 
0446     gctobj::GCTsupertower_t tempST[nSTEta][nSTPhi];
0447     gctobj::makeST(temporary, tempST);
0448     float TTseedThresholdBarrel = 5.;
0449 
0450     for (int i = 0; i < nJets; i++) {
0451       jet[i] = gctobj::getRegion(tempST, TTseedThresholdBarrel);
0452       l1tp2::Phase2L1CaloJet tempJet;
0453       int gctjeteta = jet[i].etaCenter;
0454       int gctjetphi = jet[i].phiCenter;
0455       tempJet.setJetIEta(gctjeteta + k * nBarrelEta / 2);
0456       tempJet.setJetIPhi(gctjetphi);
0457       float jeteta = realEta[gctjeteta + k * nBarrelEta / 2][gctjetphi];
0458       float jetphi = realPhi[gctjeteta + k * nBarrelEta / 2][gctjetphi];
0459       tempJet.setJetEta(jeteta);
0460       tempJet.setJetPhi(jetphi);
0461       tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
0462       tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
0463       tempJet.setTowerEt(jet[i].energyMax);
0464       int gcttowereta = jet[i].etaMax;
0465       int gcttowerphi = jet[i].phiMax;
0466       tempJet.setTowerIEta(gcttowereta + k * nBarrelEta / 2);
0467       tempJet.setTowerIPhi(gcttowerphi);
0468       float towereta = realEta[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
0469       float towerphi = realPhi[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
0470       tempJet.setTowerEta(towereta);
0471       tempJet.setTowerPhi(towerphi);
0472       reco::Candidate::PolarLorentzVector tempJetp4;
0473       tempJetp4.SetPt(tempJet.jetEt());
0474       tempJetp4.SetEta(tempJet.jetEta());
0475       tempJetp4.SetPhi(tempJet.jetPhi());
0476       tempJetp4.SetM(0.);
0477       tempJet.setP4(tempJetp4);
0478 
0479       if (jet[i].energy > 0.)
0480         halfBarrelJets.push_back(tempJet);
0481     }
0482 
0483     // ENDCAP
0484     for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
0485       for (int ieta = 0; ieta < nHgcalEta / 2; ieta++) {
0486         if (k == 0)
0487           temporary_hgcal[ieta][iphi] = hgcalTowers[ieta][iphi];
0488         else
0489           temporary_hgcal[ieta][iphi] = hgcalTowers[nHgcalEta / 2 + ieta][iphi];
0490       }
0491     }
0492 
0493     gctobj::GCTsupertower_t tempST_hgcal[nSTEta][nSTPhi];
0494     gctobj::makeST_hgcal(temporary_hgcal, tempST_hgcal);
0495     float TTseedThresholdEndcap = 3.;
0496     for (int i = nJets; i < 2 * nJets; i++) {
0497       jet[i] = gctobj::getRegion(tempST_hgcal, TTseedThresholdEndcap);
0498       l1tp2::Phase2L1CaloJet tempJet;
0499       int hgcaljeteta = jet[i].etaCenter;
0500       int hgcaljetphi = jet[i].phiCenter;
0501       tempJet.setJetIEta(hgcaljeteta + k * nHgcalEta / 2);
0502       tempJet.setJetIPhi(hgcaljetphi);
0503       float jeteta = hgcalEta[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
0504       float jetphi = hgcalPhi[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
0505       tempJet.setJetEta(jeteta);
0506       tempJet.setJetPhi(jetphi);
0507       tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
0508       tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
0509       tempJet.setTowerEt(jet[i].energyMax);
0510       int hgcaltowereta = jet[i].etaMax;
0511       int hgcaltowerphi = jet[i].phiMax;
0512       tempJet.setTowerIEta(hgcaltowereta + k * nHgcalEta / 2);
0513       tempJet.setTowerIPhi(hgcaltowerphi);
0514       float towereta = hgcalEta[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
0515       float towerphi = hgcalPhi[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
0516       tempJet.setTowerEta(towereta);
0517       tempJet.setTowerPhi(towerphi);
0518       reco::Candidate::PolarLorentzVector tempJetp4;
0519       tempJetp4.SetPt(tempJet.jetEt());
0520       tempJetp4.SetEta(tempJet.jetEta());
0521       tempJetp4.SetPhi(tempJet.jetPhi());
0522       tempJetp4.SetM(0.);
0523       tempJet.setP4(tempJetp4);
0524 
0525       if (jet[i].energy > 0.)
0526         halfHgcalJets.push_back(tempJet);
0527     }
0528 
0529     // HF
0530     for (int iphi = 0; iphi < nHfPhi; iphi++) {
0531       for (int ieta = 0; ieta < nHfEta; ieta++) {
0532         if (k == 0)
0533           temporary_hf[ieta][iphi] = hfTowers[ieta][iphi];
0534         else
0535           temporary_hf[ieta][iphi] = hfTowers[nHfEta + ieta][iphi];
0536       }
0537     }
0538 
0539     gctobj::GCTsupertower_t tempST_hf[nSTEta][nSTPhi];
0540     gctobj::makeST_hf(temporary_hf, tempST_hf);
0541     float TTseedThresholdHF = 3.;
0542     for (int i = 2 * nJets; i < 3 * nJets; i++) {
0543       jet[i] = gctobj::getRegion(tempST_hf, TTseedThresholdHF);
0544       l1tp2::Phase2L1CaloJet tempJet;
0545       int hfjeteta = jet[i].etaCenter / 2;  // 24 -> 12 towers
0546       int hfjetphi = jet[i].phiCenter;
0547       tempJet.setJetIEta(hfjeteta + k * nHfEta / 2);
0548       tempJet.setJetIPhi(hfjetphi);
0549       float jeteta = hfEta[hfjeteta + k * nHfEta / 2][hfjetphi];
0550       float jetphi = hfPhi[hfjeteta + k * nHfEta / 2][hfjetphi];
0551       tempJet.setJetEta(jeteta);
0552       tempJet.setJetPhi(jetphi);
0553       tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
0554       tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
0555       tempJet.setTowerEt(jet[i].energyMax);
0556       int hftowereta = jet[i].etaMax / 2;  // 24 -> 12 towers
0557       int hftowerphi = jet[i].phiMax;
0558       tempJet.setTowerIEta(hftowereta + k * nHfEta / 2);
0559       tempJet.setTowerIPhi(hftowerphi);
0560       float towereta = hfEta[hftowereta + k * nHfEta / 2][hftowerphi];
0561       float towerphi = hfPhi[hftowereta + k * nHfEta / 2][hftowerphi];
0562       tempJet.setTowerEta(towereta);
0563       tempJet.setTowerPhi(towerphi);
0564       reco::Candidate::PolarLorentzVector tempJetp4;
0565       tempJetp4.SetPt(tempJet.jetEt());
0566       tempJetp4.SetEta(tempJet.jetEta());
0567       tempJetp4.SetPhi(tempJet.jetPhi());
0568       tempJetp4.SetM(0.);
0569       tempJet.setP4(tempJetp4);
0570 
0571       if (jet[i].energy > 0.)
0572         halfHfJets.push_back(tempJet);
0573     }
0574 
0575     // Stitching:
0576     // if the jet eta is at the boundary: for HB it should be within 0-1 in -ve eta, 32-33 in +ve eta; for HE it should be within 0-1/16-17 in -ve eta, 34-35/18-19 in +ve eta; for HF it should be within 10-11 in -ve eta, 12-13 in +ve eta
0577     // then get the phi of that jet and check if there is a neighbouring jet with the same phi, then merge to the jet that has higher ET
0578     // in both eta/phi allow a maximum of one tower between jet centers for stitching
0579 
0580     for (size_t i = 0; i < halfHgcalJets.size(); i++) {
0581       if (halfHgcalJets.at(i).jetIEta() >= (nHgcalEta / 2 - 2) && halfHgcalJets.at(i).jetIEta() < (nHgcalEta / 2 + 2)) {
0582         float hgcal_ieta = k * nBarrelEta + halfHgcalJets.at(i).jetIEta();
0583         for (size_t j = 0; j < halfBarrelJets.size(); j++) {
0584           float barrel_ieta = nHgcalEta / 2 + halfBarrelJets.at(j).jetIEta();
0585           if (abs(barrel_ieta - hgcal_ieta) <= 2 &&
0586               abs(halfBarrelJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) <= 2) {
0587             float totalet = halfBarrelJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
0588             float totalTauEt = halfBarrelJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
0589             if (halfBarrelJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
0590               halfHgcalJets.at(i).setJetEt(0.);
0591               halfHgcalJets.at(i).setTauEt(0.);
0592               halfBarrelJets.at(j).setJetEt(totalet);
0593               halfBarrelJets.at(j).setTauEt(totalTauEt);
0594               reco::Candidate::PolarLorentzVector tempJetp4;
0595               tempJetp4.SetPt(totalet);
0596               tempJetp4.SetEta(halfBarrelJets.at(j).jetEta());
0597               tempJetp4.SetPhi(halfBarrelJets.at(j).jetPhi());
0598               tempJetp4.SetM(0.);
0599               halfBarrelJets.at(j).setP4(tempJetp4);
0600             } else {
0601               halfHgcalJets.at(i).setJetEt(totalet);
0602               halfHgcalJets.at(i).setTauEt(totalTauEt);
0603               halfBarrelJets.at(j).setJetEt(0.);
0604               halfBarrelJets.at(j).setTauEt(0.);
0605               reco::Candidate::PolarLorentzVector tempJetp4;
0606               tempJetp4.SetPt(totalet);
0607               tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
0608               tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
0609               tempJetp4.SetM(0.);
0610               halfHgcalJets.at(i).setP4(tempJetp4);
0611             }
0612           }
0613         }
0614       } else if (halfHgcalJets.at(i).jetIEta() < 2 || halfHgcalJets.at(i).jetIEta() >= (nHgcalEta - 2)) {
0615         float hgcal_ieta = k * nBarrelEta + nHfEta / 2 + halfHgcalJets.at(i).jetIEta();
0616         for (size_t j = 0; j < halfHfJets.size(); j++) {
0617           float hf_ieta = k * nBarrelEta + k * nHgcalEta + halfHfJets.at(j).jetIEta();
0618           if (abs(hgcal_ieta - hf_ieta) < 3 && abs(halfHfJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) < 3) {
0619             float totalet = halfHfJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
0620             float totalTauEt = halfHfJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
0621             if (halfHfJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
0622               halfHgcalJets.at(i).setJetEt(0.);
0623               halfHgcalJets.at(i).setTauEt(0.);
0624               halfHfJets.at(j).setJetEt(totalet);
0625               halfHfJets.at(j).setTauEt(totalTauEt);
0626               reco::Candidate::PolarLorentzVector tempJetp4;
0627               tempJetp4.SetPt(totalet);
0628               tempJetp4.SetEta(halfHfJets.at(j).jetEta());
0629               tempJetp4.SetPhi(halfHfJets.at(j).jetPhi());
0630               tempJetp4.SetM(0.);
0631               halfHfJets.at(j).setP4(tempJetp4);
0632             } else {
0633               halfHgcalJets.at(i).setJetEt(totalet);
0634               halfHgcalJets.at(i).setTauEt(totalTauEt);
0635               halfHfJets.at(j).setJetEt(0.);
0636               halfHfJets.at(j).setTauEt(0.);
0637               reco::Candidate::PolarLorentzVector tempJetp4;
0638               tempJetp4.SetPt(totalet);
0639               tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
0640               tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
0641               tempJetp4.SetM(0.);
0642               halfHgcalJets.at(i).setP4(tempJetp4);
0643             }
0644           }
0645         }
0646       }
0647     }
0648 
0649     // Write 6 leading jets from each eta half
0650 
0651     vector<l1tp2::Phase2L1CaloJet> halfAllJets;
0652     halfAllJets.clear();
0653 
0654     std::sort(halfBarrelJets.begin(), halfBarrelJets.end(), gctobj::compareByEt);
0655     for (size_t i = 0; i < halfBarrelJets.size(); i++) {
0656       if (halfBarrelJets.at(i).jetEt() > 0. && i < 6)
0657         halfAllJets.push_back(halfBarrelJets.at(i));
0658     }
0659 
0660     std::sort(halfHgcalJets.begin(), halfHgcalJets.end(), gctobj::compareByEt);
0661     for (size_t i = 0; i < halfHgcalJets.size(); i++) {
0662       if (halfHgcalJets.at(i).jetEt() > 0. && i < 6)
0663         halfAllJets.push_back(halfHgcalJets.at(i));
0664     }
0665 
0666     std::sort(halfHfJets.begin(), halfHfJets.end(), gctobj::compareByEt);
0667     for (size_t i = 0; i < halfHfJets.size(); i++) {
0668       if (halfHfJets.at(i).jetEt() > 0. && i < 6)
0669         halfAllJets.push_back(halfHfJets.at(i));
0670     }
0671 
0672     std::sort(halfAllJets.begin(), halfAllJets.end(), gctobj::compareByEt);
0673     for (size_t i = 0; i < halfAllJets.size(); i++) {
0674       if (halfAllJets.at(i).jetEt() > 0. && i < 6)
0675         allJets.push_back(halfAllJets.at(i));
0676     }
0677   }
0678 
0679   std::sort(allJets.begin(), allJets.end(), gctobj::compareByEt);
0680   for (size_t i = 0; i < allJets.size(); i++) {
0681     jetCands->push_back(allJets.at(i));
0682   }
0683 
0684   iEvent.put(std::move(jetCands), "GCTJet");
0685 }
0686 
0687 // Apply calibrations to HCAL energy based on Jet Eta, Jet pT
0688 float Phase2L1CaloJetEmulator::get_jet_pt_calibration(const float& jet_pt, const float& jet_eta) const {
0689   float abs_eta = std::abs(jet_eta);
0690   float tmp_jet_pt = jet_pt;
0691   if (tmp_jet_pt > 499)
0692     tmp_jet_pt = 499;
0693 
0694   // Different indices sizes in different calo regions.
0695   // Barrel...
0696   size_t eta_index = 0;
0697   size_t pt_index = 0;
0698   float calib = 1.0;
0699   if (abs_eta <= 1.5) {
0700     // Start loop checking 2nd value
0701     for (unsigned int i = 1; i < absEtaBinsBarrel.size(); i++) {
0702       if (abs_eta <= absEtaBinsBarrel.at(i))
0703         break;
0704       eta_index++;
0705     }
0706     // Start loop checking 2nd value
0707     for (unsigned int i = 1; i < jetPtBins.size(); i++) {
0708       if (tmp_jet_pt <= jetPtBins.at(i))
0709         break;
0710       pt_index++;
0711     }
0712     calib = calibrationsBarrel[eta_index][pt_index];
0713   }  // end Barrel
0714   else if (abs_eta <= 3.0)  // HGCal
0715   {
0716     // Start loop checking 2nd value
0717     for (unsigned int i = 1; i < absEtaBinsHGCal.size(); i++) {
0718       if (abs_eta <= absEtaBinsHGCal.at(i))
0719         break;
0720       eta_index++;
0721     }
0722     // Start loop checking 2nd value
0723     for (unsigned int i = 1; i < jetPtBins.size(); i++) {
0724       if (tmp_jet_pt <= jetPtBins.at(i))
0725         break;
0726       pt_index++;
0727     }
0728     calib = calibrationsHGCal[eta_index][pt_index];
0729   }  // end HGCal
0730   else  // HF
0731   {
0732     // Start loop checking 2nd value
0733     for (unsigned int i = 1; i < absEtaBinsHF.size(); i++) {
0734       if (abs_eta <= absEtaBinsHF.at(i))
0735         break;
0736       eta_index++;
0737     }
0738     // Start loop checking 2nd value
0739     for (unsigned int i = 1; i < jetPtBins.size(); i++) {
0740       if (tmp_jet_pt <= jetPtBins.at(i))
0741         break;
0742       pt_index++;
0743     }
0744     calib = calibrationsHF[eta_index][pt_index];
0745   }  // end HF
0746 
0747   return jet_pt * calib;
0748 }
0749 
0750 // Apply calibrations to tau pT based on L1EG info, EM Fraction, Tau Eta, Tau pT
0751 float Phase2L1CaloJetEmulator::get_tau_pt_calibration(const float& tau_pt, const float& tau_eta) const {
0752   float abs_eta = std::abs(tau_eta);
0753   float tmp_tau_pt = tau_pt;
0754   if (tmp_tau_pt > 199)
0755     tmp_tau_pt = 199;
0756 
0757   // Different indices sizes in different calo regions.
0758   // Barrel...
0759   size_t eta_index = 0;
0760   size_t pt_index = 0;
0761   float calib = 1.0;
0762   if (abs_eta <= 1.5) {
0763     // Start loop checking 2nd value
0764     for (unsigned int i = 1; i < tauAbsEtaBinsBarrel.size(); i++) {
0765       if (abs_eta <= tauAbsEtaBinsBarrel.at(i))
0766         break;
0767       eta_index++;
0768     }
0769     // Start loop checking 2nd value
0770     for (unsigned int i = 1; i < tauPtBins.size(); i++) {
0771       if (tmp_tau_pt <= tauPtBins.at(i))
0772         break;
0773       pt_index++;
0774     }
0775     calib = tauPtCalibrationsBarrel[eta_index][pt_index];
0776   }  // end Barrel
0777   else if (abs_eta <= 3.0)  // HGCal
0778   {
0779     // Start loop checking 2nd value
0780     for (unsigned int i = 1; i < tauAbsEtaBinsHGCal.size(); i++) {
0781       if (abs_eta <= tauAbsEtaBinsHGCal.at(i))
0782         break;
0783       eta_index++;
0784     }
0785     // Start loop checking 2nd value
0786     for (unsigned int i = 1; i < tauPtBins.size(); i++) {
0787       if (tmp_tau_pt <= tauPtBins.at(i))
0788         break;
0789       pt_index++;
0790     }
0791     calib = tauPtCalibrationsHGCal[eta_index][pt_index];
0792   }  // end HGCal
0793 
0794   return tau_pt * calib;
0795 }
0796 
0797 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0798 void Phase2L1CaloJetEmulator::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0799   edm::ParameterSetDescription desc;
0800   desc.add<edm::InputTag>("gctFullTowers", edm::InputTag("l1tPhase2L1CaloEGammaEmulator", "GCTFullTowers"));
0801   desc.add<edm::InputTag>("hgcalTowers", edm::InputTag("l1tHGCalTowerProducer", "HGCalTowerProcessor"));
0802   desc.add<edm::InputTag>("hcalDigis", edm::InputTag("simHcalTriggerPrimitiveDigis"));
0803 
0804   edm::ParameterSetDescription nHits_params_validator;
0805   nHits_params_validator.add<string>("fit", "type");
0806   nHits_params_validator.add<vector<double>>("nHits_params", {1., 1.});
0807   std::vector<edm::ParameterSet> nHits_params_default;
0808   edm::ParameterSet nHits_params1;
0809   nHits_params1.addParameter<string>("fit", "hgcalEM");
0810   nHits_params1.addParameter<vector<double>>("nHits_params", {157.522, 0.090});
0811   nHits_params_default.push_back(nHits_params1);
0812   edm::ParameterSet nHits_params2;
0813   nHits_params2.addParameter<string>("fit", "hgcalHad");
0814   nHits_params2.addParameter<vector<double>>("nHits_params", {159.295, 0.178});
0815   nHits_params_default.push_back(nHits_params2);
0816   edm::ParameterSet nHits_params3;
0817   nHits_params3.addParameter<string>("fit", "hf");
0818   nHits_params3.addParameter<vector<double>>("nHits_params", {165.706, 0.153});
0819   nHits_params_default.push_back(nHits_params3);
0820   desc.addVPSet("nHits_to_nvtx_params", nHits_params_validator, nHits_params_default);
0821 
0822   edm::ParameterSetDescription nvtx_params_validator;
0823   nvtx_params_validator.add<string>("calo", "type");
0824   nvtx_params_validator.add<string>("iEta", "etaregion");
0825   nvtx_params_validator.add<vector<double>>("nvtx_params", {1., 1.});
0826   std::vector<edm::ParameterSet> nvtx_params_default;
0827   edm::ParameterSet nvtx_params1;
0828   nvtx_params1.addParameter<string>("calo", "hgcalEM");
0829   nvtx_params1.addParameter<string>("iEta", "er1p4to1p8");
0830   nvtx_params1.addParameter<vector<double>>("nvtx_params", {-0.011772, 0.004142});
0831   nvtx_params_default.push_back(nvtx_params1);
0832   edm::ParameterSet nvtx_params2;
0833   nvtx_params2.addParameter<string>("calo", "hgcalEM");
0834   nvtx_params2.addParameter<string>("iEta", "er1p8to2p1");
0835   nvtx_params2.addParameter<vector<double>>("nvtx_params", {-0.015488, 0.005410});
0836   nvtx_params_default.push_back(nvtx_params2);
0837   edm::ParameterSet nvtx_params3;
0838   nvtx_params3.addParameter<string>("calo", "hgcalEM");
0839   nvtx_params3.addParameter<string>("iEta", "er2p1to2p4");
0840   nvtx_params3.addParameter<vector<double>>("nvtx_params", {-0.021150, 0.006078});
0841   nvtx_params_default.push_back(nvtx_params3);
0842   edm::ParameterSet nvtx_params4;
0843   nvtx_params4.addParameter<string>("calo", "hgcalEM");
0844   nvtx_params4.addParameter<string>("iEta", "er2p4to2p7");
0845   nvtx_params4.addParameter<vector<double>>("nvtx_params", {-0.015705, 0.005339});
0846   nvtx_params_default.push_back(nvtx_params4);
0847   edm::ParameterSet nvtx_params5;
0848   nvtx_params5.addParameter<string>("calo", "hgcalEM");
0849   nvtx_params5.addParameter<string>("iEta", "er2p7to3p1");
0850   nvtx_params5.addParameter<vector<double>>("nvtx_params", {-0.018492, 0.005620});
0851   nvtx_params_default.push_back(nvtx_params5);
0852   edm::ParameterSet nvtx_params6;
0853   nvtx_params6.addParameter<string>("calo", "hgcalHad");
0854   nvtx_params6.addParameter<string>("iEta", "er1p4to1p8");
0855   nvtx_params6.addParameter<vector<double>>("nvtx_params", {0.005675, 0.000615});
0856   nvtx_params_default.push_back(nvtx_params6);
0857   edm::ParameterSet nvtx_params7;
0858   nvtx_params7.addParameter<string>("calo", "hgcalHad");
0859   nvtx_params7.addParameter<string>("iEta", "er1p8to2p1");
0860   nvtx_params7.addParameter<vector<double>>("nvtx_params", {0.004560, 0.001099});
0861   nvtx_params_default.push_back(nvtx_params7);
0862   edm::ParameterSet nvtx_params8;
0863   nvtx_params8.addParameter<string>("calo", "hgcalHad");
0864   nvtx_params8.addParameter<string>("iEta", "er2p1to2p4");
0865   nvtx_params8.addParameter<vector<double>>("nvtx_params", {0.000036, 0.001608});
0866   nvtx_params_default.push_back(nvtx_params8);
0867   edm::ParameterSet nvtx_params9;
0868   nvtx_params9.addParameter<string>("calo", "hgcalHad");
0869   nvtx_params9.addParameter<string>("iEta", "er2p4to2p7");
0870   nvtx_params9.addParameter<vector<double>>("nvtx_params", {0.000869, 0.001754});
0871   nvtx_params_default.push_back(nvtx_params9);
0872   edm::ParameterSet nvtx_params10;
0873   nvtx_params10.addParameter<string>("calo", "hgcalHad");
0874   nvtx_params10.addParameter<string>("iEta", "er2p7to3p1");
0875   nvtx_params10.addParameter<vector<double>>("nvtx_params", {-0.006574, 0.003134});
0876   nvtx_params_default.push_back(nvtx_params10);
0877   edm::ParameterSet nvtx_params11;
0878   nvtx_params11.addParameter<string>("calo", "hf");
0879   nvtx_params11.addParameter<string>("iEta", "er29to33");
0880   nvtx_params11.addParameter<vector<double>>("nvtx_params", {-0.203291, 0.044096});
0881   nvtx_params_default.push_back(nvtx_params11);
0882   edm::ParameterSet nvtx_params12;
0883   nvtx_params12.addParameter<string>("calo", "hf");
0884   nvtx_params12.addParameter<string>("iEta", "er34to37");
0885   nvtx_params12.addParameter<vector<double>>("nvtx_params", {-0.210922, 0.045628});
0886   nvtx_params_default.push_back(nvtx_params12);
0887   edm::ParameterSet nvtx_params13;
0888   nvtx_params13.addParameter<string>("calo", "hf");
0889   nvtx_params13.addParameter<string>("iEta", "er38to41");
0890   nvtx_params13.addParameter<vector<double>>("nvtx_params", {-0.229562, 0.050560});
0891   nvtx_params_default.push_back(nvtx_params13);
0892   desc.addVPSet("nvtx_to_PU_sub_params", nvtx_params_validator, nvtx_params_default);
0893 
0894   desc.add<vector<double>>("jetPtBins", {0.0,   5.0,   7.5,   10.0,  12.5,  15.0,  17.5,  20.0,  22.5,  25.0,  27.5,
0895                                          30.0,  35.0,  40.0,  45.0,  50.0,  55.0,  60.0,  65.0,  70.0,  75.0,  80.0,
0896                                          85.0,  90.0,  95.0,  100.0, 110.0, 120.0, 130.0, 140.0, 150.0, 160.0, 170.0,
0897                                          180.0, 190.0, 200.0, 225.0, 250.0, 275.0, 300.0, 325.0, 400.0, 500.0});
0898   desc.add<vector<double>>("absEtaBinsBarrel", {0.00, 0.30, 0.60, 1.00, 1.50});
0899   desc.add<vector<double>>(
0900       "jetCalibrationsBarrel",
0901       {2.459, 2.320, 2.239, 2.166, 2.100, 2.040, 1.986, 1.937, 1.892, 1.852, 1.816, 1.768, 1.714, 1.670, 1.633, 1.603,
0902        1.578, 1.557, 1.540, 1.525, 1.513, 1.502, 1.493, 1.486, 1.479, 1.470, 1.460, 1.452, 1.445, 1.439, 1.433, 1.427,
0903        1.422, 1.417, 1.411, 1.403, 1.390, 1.377, 1.365, 1.352, 1.327, 1.284, 4.695, 3.320, 2.751, 2.361, 2.093, 1.908,
0904        1.781, 1.694, 1.633, 1.591, 1.562, 1.533, 1.511, 1.499, 1.492, 1.486, 1.482, 1.478, 1.474, 1.470, 1.467, 1.463,
0905        1.459, 1.456, 1.452, 1.447, 1.440, 1.433, 1.425, 1.418, 1.411, 1.404, 1.397, 1.390, 1.382, 1.370, 1.352, 1.334,
0906        1.316, 1.298, 1.262, 1.200, 5.100, 3.538, 2.892, 2.448, 2.143, 1.933, 1.789, 1.689, 1.620, 1.572, 1.539, 1.506,
0907        1.482, 1.469, 1.460, 1.455, 1.450, 1.446, 1.442, 1.438, 1.434, 1.431, 1.427, 1.423, 1.420, 1.414, 1.407, 1.400,
0908        1.392, 1.385, 1.378, 1.370, 1.363, 1.356, 1.348, 1.336, 1.317, 1.299, 1.281, 1.263, 1.226, 1.162, 3.850, 3.438,
0909        3.211, 3.017, 2.851, 2.708, 2.585, 2.479, 2.388, 2.310, 2.243, 2.159, 2.072, 2.006, 1.956, 1.917, 1.887, 1.863,
0910        1.844, 1.828, 1.814, 1.802, 1.791, 1.782, 1.773, 1.760, 1.744, 1.729, 1.714, 1.699, 1.685, 1.670, 1.656, 1.641,
0911        1.627, 1.602, 1.566, 1.530, 1.494, 1.458, 1.386, 1.260});
0912   desc.add<vector<double>>("absEtaBinsHGCal", {1.50, 1.90, 2.40, 3.00});
0913   desc.add<vector<double>>(
0914       "jetCalibrationsHGCal",
0915       {5.604,   4.578,  4.061,  3.647, 3.314, 3.047, 2.832, 2.660, 2.521, 2.410, 2.320, 2.216, 2.120, 2.056,
0916        2.013,   1.983,  1.961,  1.945, 1.932, 1.922, 1.913, 1.905, 1.898, 1.891, 1.884, 1.874, 1.861, 1.848,
0917        1.835,   1.822,  1.810,  1.797, 1.784, 1.771, 1.759, 1.736, 1.704, 1.673, 1.641, 1.609, 1.545, 1.434,
0918        4.385,   3.584,  3.177,  2.849, 2.584, 2.370, 2.197, 2.057, 1.944, 1.853, 1.780, 1.695, 1.616, 1.564,
0919        1.530,   1.507,  1.491,  1.480, 1.472, 1.466, 1.462, 1.459, 1.456, 1.453, 1.451, 1.447, 1.443, 1.439,
0920        1.435,   1.431,  1.427,  1.423, 1.419, 1.416, 1.412, 1.405, 1.395, 1.385, 1.376, 1.366, 1.346, 1.312,
0921        562.891, 68.647, 17.648, 5.241, 2.223, 1.490, 1.312, 1.270, 1.260, 1.259, 1.259, 1.260, 1.263, 1.265,
0922        1.267,   1.269,  1.271,  1.273, 1.275, 1.277, 1.279, 1.281, 1.283, 1.285, 1.287, 1.290, 1.295, 1.299,
0923        1.303,   1.307,  1.311,  1.315, 1.319, 1.323, 1.328, 1.335, 1.345, 1.355, 1.366, 1.376, 1.397, 1.433});
0924   desc.add<vector<double>>("absEtaBinsHF", {3.00, 3.60, 6.00});
0925   desc.add<vector<double>>(
0926       "jetCalibrationsHF",
0927       {8.169, 6.873, 6.155, 5.535, 5.001, 4.539, 4.141, 3.798, 3.501, 3.245, 3.024, 2.748, 2.463, 2.249,
0928        2.090, 1.971, 1.881, 1.814, 1.763, 1.725, 1.695, 1.673, 1.655, 1.642, 1.631, 1.618, 1.605, 1.596,
0929        1.588, 1.581, 1.575, 1.569, 1.563, 1.557, 1.551, 1.541, 1.527, 1.513, 1.498, 1.484, 1.456, 1.406,
0930        2.788, 2.534, 2.388, 2.258, 2.141, 2.037, 1.945, 1.862, 1.788, 1.722, 1.664, 1.587, 1.503, 1.436,
0931        1.382, 1.339, 1.305, 1.277, 1.255, 1.237, 1.223, 1.211, 1.201, 1.193, 1.186, 1.178, 1.170, 1.164,
0932        1.159, 1.154, 1.151, 1.147, 1.144, 1.141, 1.138, 1.133, 1.126, 1.118, 1.111, 1.104, 1.090, 1.064});
0933   desc.add<vector<double>>("tauPtBins", {0.0,  5.0,  7.5,  10.0, 12.5, 15.0, 20.0, 25.0,  30.0,  35.0,
0934                                          40.0, 45.0, 50.0, 55.0, 60.0, 70.0, 80.0, 100.0, 150.0, 200.0});
0935   desc.add<vector<double>>("tauAbsEtaBinsBarrel", {0.00, 0.30, 0.60, 1.00, 1.50});
0936   desc.add<vector<double>>("tauCalibrationsBarrel",
0937                            {1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.067,
0938                             1.067, 1.067, 1.067, 1.067, 1.067, 1.067, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106,
0939                             1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.106, 1.102,
0940                             1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102, 1.102,
0941                             1.102, 1.102, 1.102, 1.102, 1.102, 1.139, 1.139, 1.139, 1.139, 1.139, 1.139, 1.139, 1.139});
0942   desc.add<vector<double>>("tauAbsEtaBinsHGCal", {1.50, 1.90, 2.40, 3.00});
0943   desc.add<vector<double>>(
0944       "tauCalibrationsHGCal",
0945       {1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384, 1.384,
0946        1.384, 1.384, 1.384, 1.384, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473,
0947        1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.473, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133,
0948        1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133});
0949 
0950   descriptions.addWithDefaultLabel(desc);
0951 }
0952 
0953 //define this as a plug-in
0954 DEFINE_FWK_MODULE(Phase2L1CaloJetEmulator);