CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
Phase2L1CaloJetEmulator Class Reference

#include <L1Trigger/L1CaloTrigger/plugins/Phase2L1CaloJetEmulator.cc>

Inheritance diagram for Phase2L1CaloJetEmulator:
edm::stream::EDProducer<>

Public Member Functions

 Phase2L1CaloJetEmulator (const edm::ParameterSet &)
 
 ~Phase2L1CaloJetEmulator () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

float get_jet_pt_calibration (const float &jet_pt, const float &jet_eta) const
 
float get_tau_pt_calibration (const float &tau_pt, const float &tau_eta) const
 
void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

std::vector< double > absEtaBinsBarrel
 
std::vector< double > absEtaBinsHF
 
std::vector< double > absEtaBinsHGCal
 
std::map< std::string, std::map< std::string, TF1 > > all_nvtx_to_PU_sub_funcs
 
std::vector< std::vector< double > > calibrationsBarrel
 
std::vector< std::vector< double > > calibrationsHF
 
std::vector< std::vector< double > > calibrationsHGCal
 
edm::EDGetTokenT< l1tp2::CaloTowerCollectioncaloTowerToken_
 
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecorddecoderTag_
 
std::map< std::string, TF1 > hf_nvtx_to_PU_sub_funcs
 
edm::EDGetTokenT< HcalTrigPrimDigiCollectionhfToken_
 
std::map< std::string, TF1 > hgcalEM_nvtx_to_PU_sub_funcs
 
std::map< std::string, TF1 > hgcalHad_nvtx_to_PU_sub_funcs
 
edm::EDGetTokenT< l1t::HGCalTowerBxCollectionhgcalTowerToken_
 
std::vector< double > jetCalibrationsBarrel
 
std::vector< double > jetCalibrationsHF
 
std::vector< double > jetCalibrationsHGCal
 
std::vector< double > jetPtBins
 
std::map< std::string, TF1 > nHits_to_nvtx_funcs
 
std::vector< edm::ParameterSetnHits_to_nvtx_params
 
std::vector< edm::ParameterSetnvtx_to_PU_sub_params
 
std::vector< double > tauAbsEtaBinsBarrel
 
std::vector< double > tauAbsEtaBinsHGCal
 
std::vector< double > tauCalibrationsBarrel
 
std::vector< double > tauCalibrationsHGCal
 
std::vector< double > tauPtBins
 
std::vector< std::vector< double > > tauPtCalibrationsBarrel
 
std::vector< std::vector< double > > tauPtCalibrationsHGCal
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Description: Producing GCT calo jets using GCT barrel, HGCal and HF towers, based on firmware logic.

Implementation: Depends on producers for CaloTowerCollection, HGCalTowerBxCollection and HcalTrigPrimDigiCollection.

Definition at line 56 of file Phase2L1CaloJetEmulator.cc.

Constructor & Destructor Documentation

◆ Phase2L1CaloJetEmulator()

Phase2L1CaloJetEmulator::Phase2L1CaloJetEmulator ( const edm::ParameterSet iConfig)
explicit

Definition at line 110 of file Phase2L1CaloJetEmulator.cc.

References absEtaBinsBarrel, absEtaBinsHF, absEtaBinsHGCal, all_nvtx_to_PU_sub_funcs, calibrationsBarrel, calibrationsHF, calibrationsHGCal, hf_nvtx_to_PU_sub_funcs, hgcalEM_nvtx_to_PU_sub_funcs, hgcalHad_nvtx_to_PU_sub_funcs, mps_fire::i, l1tPhase2CaloJetEmulator_cfi::iEta, jetCalibrationsBarrel, jetCalibrationsHF, jetCalibrationsHGCal, jetPtBins, nHits_to_nvtx_funcs, nHits_to_nvtx_params, nvtx_to_PU_sub_params, LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, muonDTDigis_cfi::pset, DiDispStaMuonMonitor_cfi::pt, AlCaHLTBitMon_QueryRunRegistry::string, tauAbsEtaBinsBarrel, tauAbsEtaBinsHGCal, tauCalibrationsBarrel, tauCalibrationsHGCal, tauPtBins, tauPtCalibrationsBarrel, tauPtCalibrationsHGCal, and parallelization::uint.

111  : caloTowerToken_(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("gctFullTowers"))),
112  hgcalTowerToken_(consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("hgcalTowers"))),
113  hfToken_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("hcalDigis"))),
114  decoderTag_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>(edm::ESInputTag("", ""))),
115  nHits_to_nvtx_params(iConfig.getParameter<std::vector<edm::ParameterSet>>("nHits_to_nvtx_params")),
116  nvtx_to_PU_sub_params(iConfig.getParameter<std::vector<edm::ParameterSet>>("nvtx_to_PU_sub_params")),
117  jetPtBins(iConfig.getParameter<std::vector<double>>("jetPtBins")),
118  absEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("absEtaBinsBarrel")),
119  jetCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("jetCalibrationsBarrel")),
120  absEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("absEtaBinsHGCal")),
121  jetCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("jetCalibrationsHGCal")),
122  absEtaBinsHF(iConfig.getParameter<std::vector<double>>("absEtaBinsHF")),
123  jetCalibrationsHF(iConfig.getParameter<std::vector<double>>("jetCalibrationsHF")),
124  tauPtBins(iConfig.getParameter<std::vector<double>>("tauPtBins")),
125  tauAbsEtaBinsBarrel(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsBarrel")),
126  tauCalibrationsBarrel(iConfig.getParameter<std::vector<double>>("tauCalibrationsBarrel")),
127  tauAbsEtaBinsHGCal(iConfig.getParameter<std::vector<double>>("tauAbsEtaBinsHGCal")),
128  tauCalibrationsHGCal(iConfig.getParameter<std::vector<double>>("tauCalibrationsHGCal")) {
129  for (uint i = 0; i < nHits_to_nvtx_params.size(); i++) {
131  std::string calo = pset->getParameter<std::string>("fit");
132  nHits_to_nvtx_funcs[calo.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
133  nHits_to_nvtx_funcs[calo.c_str()].SetParameter(0, pset->getParameter<std::vector<double>>("nHits_params").at(0));
134  nHits_to_nvtx_funcs[calo.c_str()].SetParameter(1, pset->getParameter<std::vector<double>>("nHits_params").at(1));
135  }
139  for (uint i = 0; i < nvtx_to_PU_sub_params.size(); i++) {
141  std::string calo = pset->getParameter<std::string>("calo");
142  std::string iEta = pset->getParameter<std::string>("iEta");
143  double p1 = pset->getParameter<std::vector<double>>("nvtx_params").at(0);
144  double p2 = pset->getParameter<std::vector<double>>("nvtx_params").at(1);
145 
146  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()] = TF1(calo.c_str(), "[0] + [1] * x");
147  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(0, p1);
148  all_nvtx_to_PU_sub_funcs[calo.c_str()][iEta.c_str()].SetParameter(1, p2);
149  }
150 
151  // Fill the jet pt calibration 2D vector
152  // Dimension 1 is AbsEta bin
153  // Dimension 2 is jet pT bin which is filled with the actual callibration value
154  // size()-1 b/c the inputs have lower and upper bounds
155  // Do Barrel, then HGCal, then HF
156  int index = 0;
157  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsBarrel.size() - 1; abs_eta++) {
158  std::vector<double> pt_bin_calibs;
159  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
160  pt_bin_calibs.push_back(jetCalibrationsBarrel.at(index));
161  index++;
162  }
163  calibrationsBarrel.push_back(pt_bin_calibs);
164  }
165 
166  index = 0;
167  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHGCal.size() - 1; abs_eta++) {
168  std::vector<double> pt_bin_calibs;
169  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
170  pt_bin_calibs.push_back(jetCalibrationsHGCal.at(index));
171  index++;
172  }
173  calibrationsHGCal.push_back(pt_bin_calibs);
174  }
175 
176  index = 0;
177  for (unsigned int abs_eta = 0; abs_eta < absEtaBinsHF.size() - 1; abs_eta++) {
178  std::vector<double> pt_bin_calibs;
179  for (unsigned int pt = 0; pt < jetPtBins.size() - 1; pt++) {
180  pt_bin_calibs.push_back(jetCalibrationsHF.at(index));
181  index++;
182  }
183  calibrationsHF.push_back(pt_bin_calibs);
184  }
185 
186  // Fill the tau pt calibration 2D vector
187  // Dimension 1 is AbsEta bin
188  // Dimension 2 is tau pT bin which is filled with the actual calibration value
189  // Do Barrel, then HGCal
190  //
191  // Note to future developers: be very concious of the order in which the calibrations are printed
192  // out in tool which makse the cfg files. You need to match that exactly when loading them and
193  // using the calibrations below.
194  index = 0;
195  for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsBarrel.size() - 1; abs_eta++) {
196  std::vector<double> pt_bin_calibs;
197  for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
198  pt_bin_calibs.push_back(tauCalibrationsBarrel.at(index));
199  index++;
200  }
201  tauPtCalibrationsBarrel.push_back(pt_bin_calibs);
202  }
203 
204  index = 0;
205  for (unsigned int abs_eta = 0; abs_eta < tauAbsEtaBinsHGCal.size() - 1; abs_eta++) {
206  std::vector<double> pt_bin_calibs;
207  for (unsigned int pt = 0; pt < tauPtBins.size() - 1; pt++) {
208  pt_bin_calibs.push_back(tauCalibrationsHGCal.at(index));
209  index++;
210  }
211  tauPtCalibrationsHGCal.push_back(pt_bin_calibs);
212  }
213 
214  produces<l1tp2::Phase2L1CaloJetCollection>("GCTJet");
215 }
std::vector< double > jetPtBins
std::map< std::string, TF1 > hgcalEM_nvtx_to_PU_sub_funcs
std::map< std::string, TF1 > hgcalHad_nvtx_to_PU_sub_funcs
std::map< std::string, std::map< std::string, TF1 > > all_nvtx_to_PU_sub_funcs
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< double > tauPtBins
std::vector< double > jetCalibrationsBarrel
std::vector< std::vector< double > > calibrationsBarrel
std::map< std::string, TF1 > nHits_to_nvtx_funcs
std::vector< double > tauCalibrationsHGCal
std::vector< double > absEtaBinsHGCal
edm::EDGetTokenT< l1t::HGCalTowerBxCollection > hgcalTowerToken_
std::vector< edm::ParameterSet > nHits_to_nvtx_params
std::vector< std::vector< double > > calibrationsHGCal
std::vector< double > jetCalibrationsHF
std::vector< std::vector< double > > calibrationsHF
std::vector< double > tauAbsEtaBinsHGCal
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecord > decoderTag_
std::map< std::string, TF1 > hf_nvtx_to_PU_sub_funcs
std::vector< std::vector< double > > tauPtCalibrationsHGCal
edm::EDGetTokenT< HcalTrigPrimDigiCollection > hfToken_
std::vector< double > absEtaBinsBarrel
std::vector< double > tauCalibrationsBarrel
std::vector< double > tauAbsEtaBinsBarrel
std::vector< edm::ParameterSet > nvtx_to_PU_sub_params
std::vector< double > absEtaBinsHF
std::vector< std::vector< double > > tauPtCalibrationsBarrel
edm::EDGetTokenT< l1tp2::CaloTowerCollection > caloTowerToken_
Definition: Common.h:9
std::vector< double > jetCalibrationsHGCal

◆ ~Phase2L1CaloJetEmulator()

Phase2L1CaloJetEmulator::~Phase2L1CaloJetEmulator ( )
override

Definition at line 217 of file Phase2L1CaloJetEmulator.cc.

217 {}

Member Function Documentation

◆ fillDescriptions()

void Phase2L1CaloJetEmulator::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 768 of file Phase2L1CaloJetEmulator.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSet::addParameter(), edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

768  {
770  desc.add<edm::InputTag>("gctFullTowers", edm::InputTag("l1tPhase2L1CaloEGammaEmulator", "GCTFullTowers"));
771  desc.add<edm::InputTag>("hgcalTowers", edm::InputTag("l1tHGCalTowerProducer", "HGCalTowerProcessor"));
772  desc.add<edm::InputTag>("hcalDigis", edm::InputTag("simHcalTriggerPrimitiveDigis"));
773 
774  edm::ParameterSetDescription nHits_params_validator;
775  nHits_params_validator.add<string>("fit", "type");
776  nHits_params_validator.add<vector<double>>("nHits_params", {1., 1.});
777  std::vector<edm::ParameterSet> nHits_params_default;
778  edm::ParameterSet nHits_params1;
779  nHits_params1.addParameter<string>("fit", "hgcalEM");
780  nHits_params1.addParameter<vector<double>>("nHits_params", {157.522, 0.090});
781  nHits_params_default.push_back(nHits_params1);
782  edm::ParameterSet nHits_params2;
783  nHits_params2.addParameter<string>("fit", "hgcalHad");
784  nHits_params2.addParameter<vector<double>>("nHits_params", {159.295, 0.178});
785  nHits_params_default.push_back(nHits_params2);
786  edm::ParameterSet nHits_params3;
787  nHits_params3.addParameter<string>("fit", "hf");
788  nHits_params3.addParameter<vector<double>>("nHits_params", {165.706, 0.153});
789  nHits_params_default.push_back(nHits_params3);
790  desc.addVPSet("nHits_to_nvtx_params", nHits_params_validator, nHits_params_default);
791 
792  edm::ParameterSetDescription nvtx_params_validator;
793  nvtx_params_validator.add<string>("calo", "type");
794  nvtx_params_validator.add<string>("iEta", "etaregion");
795  nvtx_params_validator.add<vector<double>>("nvtx_params", {1., 1.});
796  std::vector<edm::ParameterSet> nvtx_params_default;
797  edm::ParameterSet nvtx_params1;
798  nvtx_params1.addParameter<string>("calo", "hgcalEM");
799  nvtx_params1.addParameter<string>("iEta", "er1p4to1p8");
800  nvtx_params1.addParameter<vector<double>>("nvtx_params", {-0.011772, 0.004142});
801  nvtx_params_default.push_back(nvtx_params1);
802  edm::ParameterSet nvtx_params2;
803  nvtx_params2.addParameter<string>("calo", "hgcalEM");
804  nvtx_params2.addParameter<string>("iEta", "er1p8to2p1");
805  nvtx_params2.addParameter<vector<double>>("nvtx_params", {-0.015488, 0.005410});
806  nvtx_params_default.push_back(nvtx_params2);
807  edm::ParameterSet nvtx_params3;
808  nvtx_params3.addParameter<string>("calo", "hgcalEM");
809  nvtx_params3.addParameter<string>("iEta", "er2p1to2p4");
810  nvtx_params3.addParameter<vector<double>>("nvtx_params", {-0.021150, 0.006078});
811  nvtx_params_default.push_back(nvtx_params3);
812  edm::ParameterSet nvtx_params4;
813  nvtx_params4.addParameter<string>("calo", "hgcalEM");
814  nvtx_params4.addParameter<string>("iEta", "er2p4to2p7");
815  nvtx_params4.addParameter<vector<double>>("nvtx_params", {-0.015705, 0.005339});
816  nvtx_params_default.push_back(nvtx_params4);
817  edm::ParameterSet nvtx_params5;
818  nvtx_params5.addParameter<string>("calo", "hgcalEM");
819  nvtx_params5.addParameter<string>("iEta", "er2p7to3p1");
820  nvtx_params5.addParameter<vector<double>>("nvtx_params", {-0.018492, 0.005620});
821  nvtx_params_default.push_back(nvtx_params5);
822  edm::ParameterSet nvtx_params6;
823  nvtx_params6.addParameter<string>("calo", "hgcalHad");
824  nvtx_params6.addParameter<string>("iEta", "er1p4to1p8");
825  nvtx_params6.addParameter<vector<double>>("nvtx_params", {0.005675, 0.000615});
826  nvtx_params_default.push_back(nvtx_params6);
827  edm::ParameterSet nvtx_params7;
828  nvtx_params7.addParameter<string>("calo", "hgcalHad");
829  nvtx_params7.addParameter<string>("iEta", "er1p8to2p1");
830  nvtx_params7.addParameter<vector<double>>("nvtx_params", {0.004560, 0.001099});
831  nvtx_params_default.push_back(nvtx_params7);
832  edm::ParameterSet nvtx_params8;
833  nvtx_params8.addParameter<string>("calo", "hgcalHad");
834  nvtx_params8.addParameter<string>("iEta", "er2p1to2p4");
835  nvtx_params8.addParameter<vector<double>>("nvtx_params", {0.000036, 0.001608});
836  nvtx_params_default.push_back(nvtx_params8);
837  edm::ParameterSet nvtx_params9;
838  nvtx_params9.addParameter<string>("calo", "hgcalHad");
839  nvtx_params9.addParameter<string>("iEta", "er2p4to2p7");
840  nvtx_params9.addParameter<vector<double>>("nvtx_params", {0.000869, 0.001754});
841  nvtx_params_default.push_back(nvtx_params9);
842  edm::ParameterSet nvtx_params10;
843  nvtx_params10.addParameter<string>("calo", "hgcalHad");
844  nvtx_params10.addParameter<string>("iEta", "er2p7to3p1");
845  nvtx_params10.addParameter<vector<double>>("nvtx_params", {-0.006574, 0.003134});
846  nvtx_params_default.push_back(nvtx_params10);
847  edm::ParameterSet nvtx_params11;
848  nvtx_params11.addParameter<string>("calo", "hf");
849  nvtx_params11.addParameter<string>("iEta", "er29to33");
850  nvtx_params11.addParameter<vector<double>>("nvtx_params", {-0.203291, 0.044096});
851  nvtx_params_default.push_back(nvtx_params11);
852  edm::ParameterSet nvtx_params12;
853  nvtx_params12.addParameter<string>("calo", "hf");
854  nvtx_params12.addParameter<string>("iEta", "er34to37");
855  nvtx_params12.addParameter<vector<double>>("nvtx_params", {-0.210922, 0.045628});
856  nvtx_params_default.push_back(nvtx_params12);
857  edm::ParameterSet nvtx_params13;
858  nvtx_params13.addParameter<string>("calo", "hf");
859  nvtx_params13.addParameter<string>("iEta", "er38to41");
860  nvtx_params13.addParameter<vector<double>>("nvtx_params", {-0.229562, 0.050560});
861  nvtx_params_default.push_back(nvtx_params13);
862  desc.addVPSet("nvtx_to_PU_sub_params", nvtx_params_validator, nvtx_params_default);
863 
864  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,
865  30.0, 35.0, 40.0, 45.0, 50.0, 55.0, 60.0, 65.0, 70.0, 75.0, 80.0,
866  85.0, 90.0, 95.0, 100.0, 110.0, 120.0, 130.0, 140.0, 150.0, 160.0, 170.0,
867  180.0, 190.0, 200.0, 225.0, 250.0, 275.0, 300.0, 325.0, 400.0, 500.0});
868  desc.add<vector<double>>("absEtaBinsBarrel", {0.00, 0.30, 0.60, 1.00, 1.50});
869  desc.add<vector<double>>(
870  "jetCalibrationsBarrel",
871  {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,
872  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,
873  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,
874  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,
875  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,
876  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,
877  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,
878  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,
879  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,
880  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,
881  1.627, 1.602, 1.566, 1.530, 1.494, 1.458, 1.386, 1.260});
882  desc.add<vector<double>>("absEtaBinsHGCal", {1.50, 1.90, 2.40, 3.00});
883  desc.add<vector<double>>(
884  "jetCalibrationsHGCal",
885  {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,
886  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,
887  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,
888  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,
889  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,
890  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,
891  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,
892  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,
893  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});
894  desc.add<vector<double>>("absEtaBinsHF", {3.00, 3.60, 6.00});
895  desc.add<vector<double>>(
896  "jetCalibrationsHF",
897  {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,
898  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,
899  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,
900  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,
901  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,
902  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});
903  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,
904  40.0, 45.0, 50.0, 55.0, 60.0, 70.0, 80.0, 100.0, 150.0, 200.0});
905  desc.add<vector<double>>("tauAbsEtaBinsBarrel", {0.00, 0.30, 0.60, 1.00, 1.50});
906  desc.add<vector<double>>("tauCalibrationsBarrel",
907  {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,
908  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,
909  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,
910  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,
911  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});
912  desc.add<vector<double>>("tauAbsEtaBinsHGCal", {1.50, 1.90, 2.40, 3.00});
913  desc.add<vector<double>>(
914  "tauCalibrationsHGCal",
915  {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,
916  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,
917  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,
918  1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133, 1.133});
919 
920  descriptions.addWithDefaultLabel(desc);
921 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:136
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ get_jet_pt_calibration()

float Phase2L1CaloJetEmulator::get_jet_pt_calibration ( const float &  jet_pt,
const float &  jet_eta 
) const
private

Definition at line 656 of file Phase2L1CaloJetEmulator.cc.

References funct::abs(), absEtaBinsBarrel, absEtaBinsHF, absEtaBinsHGCal, l1tEGammaCrystalsEmulatorProducer_cfi::calib, calibrationsBarrel, calibrationsHF, calibrationsHGCal, mps_fire::i, and jetPtBins.

Referenced by produce().

656  {
657  float abs_eta = std::abs(jet_eta);
658  float tmp_jet_pt = jet_pt;
659  if (tmp_jet_pt > 499)
660  tmp_jet_pt = 499;
661 
662  // Different indices sizes in different calo regions.
663  // Barrel...
664  size_t eta_index = 0;
665  size_t pt_index = 0;
666  float calib = 1.0;
667  if (abs_eta <= 1.5) {
668  // Start loop checking 2nd value
669  for (unsigned int i = 1; i < absEtaBinsBarrel.size(); i++) {
670  if (abs_eta <= absEtaBinsBarrel.at(i))
671  break;
672  eta_index++;
673  }
674  // Start loop checking 2nd value
675  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
676  if (tmp_jet_pt <= jetPtBins.at(i))
677  break;
678  pt_index++;
679  }
680  calib = calibrationsBarrel[eta_index][pt_index];
681  } // end Barrel
682  else if (abs_eta <= 3.0) // HGCal
683  {
684  // Start loop checking 2nd value
685  for (unsigned int i = 1; i < absEtaBinsHGCal.size(); i++) {
686  if (abs_eta <= absEtaBinsHGCal.at(i))
687  break;
688  eta_index++;
689  }
690  // Start loop checking 2nd value
691  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
692  if (tmp_jet_pt <= jetPtBins.at(i))
693  break;
694  pt_index++;
695  }
696  calib = calibrationsHGCal[eta_index][pt_index];
697  } // end HGCal
698  else // HF
699  {
700  // Start loop checking 2nd value
701  for (unsigned int i = 1; i < absEtaBinsHF.size(); i++) {
702  if (abs_eta <= absEtaBinsHF.at(i))
703  break;
704  eta_index++;
705  }
706  // Start loop checking 2nd value
707  for (unsigned int i = 1; i < jetPtBins.size(); i++) {
708  if (tmp_jet_pt <= jetPtBins.at(i))
709  break;
710  pt_index++;
711  }
712  calib = calibrationsHF[eta_index][pt_index];
713  } // end HF
714 
715  return jet_pt * calib;
716 }
std::vector< double > jetPtBins
std::vector< std::vector< double > > calibrationsBarrel
std::vector< double > absEtaBinsHGCal
std::vector< std::vector< double > > calibrationsHGCal
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::vector< double > > calibrationsHF
std::vector< double > absEtaBinsBarrel
std::vector< double > absEtaBinsHF

◆ get_tau_pt_calibration()

float Phase2L1CaloJetEmulator::get_tau_pt_calibration ( const float &  tau_pt,
const float &  tau_eta 
) const
private

Definition at line 719 of file Phase2L1CaloJetEmulator.cc.

References funct::abs(), l1tEGammaCrystalsEmulatorProducer_cfi::calib, mps_fire::i, runTauDisplay::tau_eta, runTauDisplay::tau_pt, tauAbsEtaBinsBarrel, tauAbsEtaBinsHGCal, tauPtBins, tauPtCalibrationsBarrel, and tauPtCalibrationsHGCal.

Referenced by produce().

719  {
720  float abs_eta = std::abs(tau_eta);
721  float tmp_tau_pt = tau_pt;
722  if (tmp_tau_pt > 199)
723  tmp_tau_pt = 199;
724 
725  // Different indices sizes in different calo regions.
726  // Barrel...
727  size_t eta_index = 0;
728  size_t pt_index = 0;
729  float calib = 1.0;
730  if (abs_eta <= 1.5) {
731  // Start loop checking 2nd value
732  for (unsigned int i = 1; i < tauAbsEtaBinsBarrel.size(); i++) {
733  if (abs_eta <= tauAbsEtaBinsBarrel.at(i))
734  break;
735  eta_index++;
736  }
737  // Start loop checking 2nd value
738  for (unsigned int i = 1; i < tauPtBins.size(); i++) {
739  if (tmp_tau_pt <= tauPtBins.at(i))
740  break;
741  pt_index++;
742  }
743  calib = tauPtCalibrationsBarrel[eta_index][pt_index];
744  } // end Barrel
745  else if (abs_eta <= 3.0) // HGCal
746  {
747  // Start loop checking 2nd value
748  for (unsigned int i = 1; i < tauAbsEtaBinsHGCal.size(); i++) {
749  if (abs_eta <= tauAbsEtaBinsHGCal.at(i))
750  break;
751  eta_index++;
752  }
753  // Start loop checking 2nd value
754  for (unsigned int i = 1; i < tauPtBins.size(); i++) {
755  if (tmp_tau_pt <= tauPtBins.at(i))
756  break;
757  pt_index++;
758  }
759  calib = tauPtCalibrationsHGCal[eta_index][pt_index];
760  } // end HGCal
761  else
762  return calib;
763 
764  return tau_pt * calib;
765 }
std::vector< double > tauPtBins
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > tauAbsEtaBinsHGCal
std::vector< std::vector< double > > tauPtCalibrationsHGCal
std::vector< double > tauAbsEtaBinsBarrel
std::vector< std::vector< double > > tauPtCalibrationsBarrel

◆ produce()

void Phase2L1CaloJetEmulator::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 224 of file Phase2L1CaloJetEmulator.cc.

References funct::abs(), all_nvtx_to_PU_sub_funcs, reco::JetExtendedAssociation::allJets(), BXVector< T >::begin(), caloTowerToken_, gctobj::compareByEt(), decoderTag_, BXVector< T >::end(), hcalRecHitTable_cff::energy, EgHLTOffHistBins_cfi::et, PVValHelper::eta, get_jet_pt_calibration(), get_tau_pt_calibration(), edm::EventSetup::getData(), gctobj::getRegion(), hfToken_, l1tPhase2CaloJetEmulator_cfi::hgcalTowers, hgcalTowerToken_, mps_fire::i, hit::id, hcalRecHitTable_cff::ieta, iEvent, hcalRecHitTable_cff::iphi, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, dqmiolumiharvest::j, metsig::jet, l1tp2::Phase2L1CaloJet::jetEt(), l1tp2::Phase2L1CaloJet::jetEta(), l1tp2::Phase2L1CaloJet::jetPhi(), dqmdumpme::k, l1t::CaloTools::kHFBegin, l1t::CaloTools::kHFEnd, M_PI, gctobj::makeST(), gctobj::makeST_hf(), gctobj::makeST_hgcal(), eostools::move(), nBarrelEta, nBarrelPhi, nHfEta, nHfPhi, nHgcalEta, nHgcalPhi, nHits_to_nvtx_funcs, nJets, nSTEta, nSTPhi, edm::Handle< T >::product(), l1tp2::Phase2L1CaloJet::setJetEt(), l1tp2::Phase2L1CaloJet::setJetEta(), l1tp2::Phase2L1CaloJet::setJetIEta(), l1tp2::Phase2L1CaloJet::setJetIPhi(), l1tp2::Phase2L1CaloJet::setJetPhi(), reco::LeafCandidate::setP4(), l1tp2::Phase2L1CaloJet::setTauEt(), l1tp2::Phase2L1CaloJet::setTowerEt(), l1tp2::Phase2L1CaloJet::setTowerEta(), l1tp2::Phase2L1CaloJet::setTowerIEta(), l1tp2::Phase2L1CaloJet::setTowerIPhi(), l1tp2::Phase2L1CaloJet::setTowerPhi(), jetUpdater_cfi::sort, AlCaHLTBitMon_QueryRunRegistry::string, groupFilesInBlocks::temp, and l1t::CaloTools::towerEta().

224  {
225  using namespace edm;
226  std::unique_ptr<l1tp2::Phase2L1CaloJetCollection> jetCands(make_unique<l1tp2::Phase2L1CaloJetCollection>());
227 
228  // Assign ETs to each eta-half of the barrel region (17x72 --> 18x72 to be able to make 3x3 super towers)
229  edm::Handle<std::vector<l1tp2::CaloTower>> caloTowerCollection;
230  if (!iEvent.getByToken(caloTowerToken_, caloTowerCollection))
231  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from caloTowerCollection!";
232 
233  iEvent.getByToken(caloTowerToken_, caloTowerCollection);
234  float GCTintTowers[nBarrelEta][nBarrelPhi];
235  float realEta[nBarrelEta][nBarrelPhi];
236  float realPhi[nBarrelEta][nBarrelPhi];
237  for (const l1tp2::CaloTower& i : *caloTowerCollection) {
238  int ieta = i.towerIEta();
239  int iphi = i.towerIPhi();
240  if (i.ecalTowerEt() > 1.)
241  GCTintTowers[ieta][iphi] = i.ecalTowerEt(); // suppress <= 1 GeV towers
242  else
243  GCTintTowers[ieta][iphi] = 0;
244  realEta[ieta][iphi] = i.towerEta();
245  realPhi[ieta][iphi] = i.towerPhi();
246  }
247 
248  float temporary[nBarrelEta / 2][nBarrelPhi];
249 
250  // HGCal and HF info used for nvtx estimation
251  edm::Handle<l1t::HGCalTowerBxCollection> hgcalTowerCollection;
252  if (!iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection))
253  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from hgcalTowerCollection!";
254  l1t::HGCalTowerBxCollection hgcalTowerColl;
255  iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection);
256  hgcalTowerColl = (*hgcalTowerCollection.product());
257 
259  if (!iEvent.getByToken(hfToken_, hfHandle))
260  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get HcalTrigPrimDigi for HF!";
261  iEvent.getByToken(hfToken_, hfHandle);
262 
263  int i_hgcalEM_hits_leq_threshold = 0;
264  int i_hgcalHad_hits_leq_threshold = 0;
265  int i_hf_hits_leq_threshold = 0;
266  for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
267  if (it->etEm() <= 1.75 && it->etEm() >= 1.25) {
268  i_hgcalEM_hits_leq_threshold++;
269  }
270  if (it->etHad() <= 1.25 && it->etHad() >= 0.75) {
271  i_hgcalHad_hits_leq_threshold++;
272  }
273  }
274  const auto& decoder = iSetup.getData(decoderTag_);
275  for (const auto& hit : *hfHandle.product()) {
276  double et = decoder.hcaletValue(hit.id(), hit.t0());
277  if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
278  continue;
279  if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
280  continue;
281  if (et <= 15.0 && et >= 10.0)
282  i_hf_hits_leq_threshold++;
283  }
284 
285  double hgcalEM_nvtx = nHits_to_nvtx_funcs["hgcalEM"].Eval(i_hgcalEM_hits_leq_threshold);
286  if (hgcalEM_nvtx < 0)
287  hgcalEM_nvtx = 0;
288  double hgcalHad_nvtx = nHits_to_nvtx_funcs["hgcalHad"].Eval(i_hgcalHad_hits_leq_threshold);
289  if (hgcalHad_nvtx < 0)
290  hgcalHad_nvtx = 0;
291  double hf_nvtx = nHits_to_nvtx_funcs["hf"].Eval(i_hf_hits_leq_threshold);
292  if (hf_nvtx < 0)
293  hf_nvtx = 0;
294  double EstimatedNvtx = (hgcalEM_nvtx + hgcalHad_nvtx + hf_nvtx) / 3.;
295 
296  // Assign ETs to each eta-half of the endcap region (18x72)
298  float hgcalEta[nHgcalEta][nHgcalPhi];
299  float hgcalPhi[nHgcalEta][nHgcalPhi];
300 
301  for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
302  for (int ieta = 0; ieta < nHgcalEta; ieta++) {
303  hgcalTowers[ieta][iphi] = 0;
304  if (ieta < nHgcalEta / 2)
305  hgcalEta[ieta][iphi] = -3.045 + ieta * 0.087 + 0.0435;
306  else
307  hgcalEta[ieta][iphi] = 1.479 + (ieta - nHgcalEta / 2) * 0.087 + 0.0435;
308  hgcalPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHgcalPhi) + (M_PI / nHgcalPhi);
309  }
310  }
311 
312  for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
313  float eta = it->eta();
314  int ieta;
315  if (eta < 0)
316  ieta = 19 - it->id().iEta();
317  else
318  ieta = 20 + it->id().iEta();
319  if (eta > 1.479)
320  ieta = ieta - 4;
321  int iphi = it->id().iPhi();
322 
323  float hgcal_etEm = it->etEm();
324  float hgcal_etHad = it->etHad();
325  std::string etaKey = "";
326  if (abs(eta) <= 1.8)
327  etaKey = "er1p4to1p8";
328  else if (abs(eta) <= 2.1 && abs(eta) > 1.8)
329  etaKey = "er1p8to2p1";
330  else if (abs(eta) <= 2.4 && abs(eta) > 2.1)
331  etaKey = "er2p1to2p4";
332  else if (abs(eta) <= 2.7 && abs(eta) > 2.4)
333  etaKey = "er2p4to2p7";
334  else if (abs(eta) <= 3.1 && abs(eta) > 2.7)
335  etaKey = "er2p7to3p1";
336  if (!etaKey.empty()) {
337  hgcal_etEm = it->etEm() - all_nvtx_to_PU_sub_funcs["hgcalEM"][etaKey].Eval(EstimatedNvtx);
338  hgcal_etHad = it->etHad() - all_nvtx_to_PU_sub_funcs["hgcalHad"][etaKey].Eval(EstimatedNvtx);
339  }
340 
341  if (hgcal_etEm < 0)
342  hgcal_etEm = 0;
343  if (hgcal_etHad < 0)
344  hgcal_etHad = 0;
345  if ((it->etEm() + it->etHad() > 1.) && abs(eta) > 1.479)
346  hgcalTowers[ieta][iphi] = hgcal_etEm + hgcal_etHad; // suppress <= 1 GeV towers
347  }
348 
349  float temporary_hgcal[nHgcalEta / 2][nHgcalPhi];
350 
351  // Assign ETs to each eta-half of the forward region (12x72)
352  float hfTowers[nHfEta][nHfPhi];
353  float hfEta[nHfEta][nHfPhi];
354  float hfPhi[nHfEta][nHfPhi];
355  for (int iphi = 0; iphi < nHfPhi; iphi++) {
356  for (int ieta = 0; ieta < nHfEta; ieta++) {
357  hfTowers[ieta][iphi] = 0;
358  int temp;
359  if (ieta < nHfEta / 2)
361  else
364  hfPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHfPhi) + (M_PI / nHfPhi);
365  }
366  }
367 
368  for (const auto& hit : *hfHandle.product()) {
369  double et = decoder.hcaletValue(hit.id(), hit.t0());
370  int ieta = 0;
371  if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
372  continue;
373  if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
374  continue;
375  if (hit.id().ieta() <= -(l1t::CaloTools::kHFBegin + 1)) {
376  ieta = hit.id().ieta() + l1t::CaloTools::kHFEnd;
377  } else if (hit.id().ieta() >= (l1t::CaloTools::kHFBegin + 1)) {
378  ieta = nHfEta / 2 + (hit.id().ieta() - (l1t::CaloTools::kHFBegin + 1));
379  }
380  int iphi = 0;
381  if (hit.id().iphi() <= nHfPhi / 2)
382  iphi = hit.id().iphi() + (nHfPhi / 2 - 1);
383  else if (hit.id().iphi() > nHfPhi / 2)
384  iphi = hit.id().iphi() - (nHfPhi / 2 + 1);
385  if (abs(hit.id().ieta()) <= 33 && abs(hit.id().ieta()) >= 29)
386  et = et - all_nvtx_to_PU_sub_funcs["hf"]["er29to33"].Eval(EstimatedNvtx);
387  if (abs(hit.id().ieta()) <= 37 && abs(hit.id().ieta()) >= 34)
388  et = et - all_nvtx_to_PU_sub_funcs["hf"]["er34to37"].Eval(EstimatedNvtx);
389  if (abs(hit.id().ieta()) <= 41 && abs(hit.id().ieta()) >= 38)
390  et = et - all_nvtx_to_PU_sub_funcs["hf"]["er38to41"].Eval(EstimatedNvtx);
391  if (et < 0)
392  et = 0;
393  if (et > 1.)
394  hfTowers[ieta][iphi] = et; // suppress <= 1 GeV towers
395  }
396 
397  float temporary_hf[nHfEta / 2][nHfPhi];
398 
399  // Begin creating jets
400  // First create 3x3 super towers: 6x24 in barrel, endcap; 4x24 in forward
401  // Then create up to 10 jets in each eta half of barrel, endcap, forward regions
402 
403  vector<l1tp2::Phase2L1CaloJet> halfBarrelJets, halfHgcalJets, halfHfJets;
404  halfBarrelJets.clear();
405  halfHgcalJets.clear();
406  halfHfJets.clear();
407  vector<l1tp2::Phase2L1CaloJet> allJets;
408  allJets.clear();
409 
410  for (int k = 0; k < 2; k++) {
411  halfBarrelJets.clear();
412  halfHgcalJets.clear();
413  halfHfJets.clear();
415 
416  // BARREL
417  for (int iphi = 0; iphi < nBarrelPhi; iphi++) {
418  for (int ieta = 0; ieta < nBarrelEta / 2; ieta++) {
419  if (k == 0)
420  temporary[ieta][iphi] = GCTintTowers[ieta][iphi];
421  else
422  temporary[ieta][iphi] = GCTintTowers[nBarrelEta / 2 + ieta][iphi];
423  }
424  }
425 
427  gctobj::makeST(temporary, tempST);
428 
429  for (int i = 0; i < nJets; i++) {
430  jet[i] = gctobj::getRegion(tempST);
431  l1tp2::Phase2L1CaloJet tempJet;
432  int gctjeteta = jet[i].etaCenter;
433  int gctjetphi = jet[i].phiCenter;
434  tempJet.setJetIEta(gctjeteta + k * nBarrelEta / 2);
435  tempJet.setJetIPhi(gctjetphi);
436  float jeteta = realEta[gctjeteta + k * nBarrelEta / 2][gctjetphi];
437  float jetphi = realPhi[gctjeteta + k * nBarrelEta / 2][gctjetphi];
438  tempJet.setJetEta(jeteta);
439  tempJet.setJetPhi(jetphi);
440  tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
441  tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
442  tempJet.setTowerEt(jet[i].energyMax);
443  int gcttowereta = jet[i].etaMax;
444  int gcttowerphi = jet[i].phiMax;
445  tempJet.setTowerIEta(gcttowereta + k * nBarrelEta / 2);
446  tempJet.setTowerIPhi(gcttowerphi);
447  float towereta = realEta[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
448  float towerphi = realPhi[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
449  tempJet.setTowerEta(towereta);
450  tempJet.setTowerPhi(towerphi);
452  tempJetp4.SetPt(tempJet.jetEt());
453  tempJetp4.SetEta(tempJet.jetEta());
454  tempJetp4.SetPhi(tempJet.jetPhi());
455  tempJetp4.SetM(0.);
456  tempJet.setP4(tempJetp4);
457 
458  if (jet[i].energy > 0.)
459  halfBarrelJets.push_back(tempJet);
460  }
461 
462  // ENDCAP
463  for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
464  for (int ieta = 0; ieta < nHgcalEta / 2; ieta++) {
465  if (k == 0)
466  temporary_hgcal[ieta][iphi] = hgcalTowers[ieta][iphi];
467  else
468  temporary_hgcal[ieta][iphi] = hgcalTowers[nHgcalEta / 2 + ieta][iphi];
469  }
470  }
471 
472  gctobj::GCTsupertower_t tempST_hgcal[nSTEta][nSTPhi];
473  gctobj::makeST_hgcal(temporary_hgcal, tempST_hgcal);
474  for (int i = nJets; i < 2 * nJets; i++) {
475  jet[i] = gctobj::getRegion(tempST_hgcal);
476  l1tp2::Phase2L1CaloJet tempJet;
477  int hgcaljeteta = jet[i].etaCenter;
478  int hgcaljetphi = jet[i].phiCenter;
479  tempJet.setJetIEta(hgcaljeteta + k * nHgcalEta / 2);
480  tempJet.setJetIPhi(hgcaljetphi);
481  float jeteta = hgcalEta[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
482  float jetphi = hgcalPhi[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
483  tempJet.setJetEta(jeteta);
484  tempJet.setJetPhi(jetphi);
485  tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
486  tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
487  tempJet.setTowerEt(jet[i].energyMax);
488  int hgcaltowereta = jet[i].etaMax;
489  int hgcaltowerphi = jet[i].phiMax;
490  tempJet.setTowerIEta(hgcaltowereta + k * nHgcalEta / 2);
491  tempJet.setTowerIPhi(hgcaltowerphi);
492  float towereta = hgcalEta[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
493  float towerphi = hgcalPhi[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
494  tempJet.setTowerEta(towereta);
495  tempJet.setTowerPhi(towerphi);
497  tempJetp4.SetPt(tempJet.jetEt());
498  tempJetp4.SetEta(tempJet.jetEta());
499  tempJetp4.SetPhi(tempJet.jetPhi());
500  tempJetp4.SetM(0.);
501  tempJet.setP4(tempJetp4);
502 
503  if (jet[i].energy > 0.)
504  halfHgcalJets.push_back(tempJet);
505  }
506 
507  // HF
508  for (int iphi = 0; iphi < nHfPhi; iphi++) {
509  for (int ieta = 0; ieta < nHfEta / 2; ieta++) {
510  if (k == 0)
511  temporary_hf[ieta][iphi] = hfTowers[ieta][iphi];
512  else
513  temporary_hf[ieta][iphi] = hfTowers[nHfEta / 2 + ieta][iphi];
514  }
515  }
516 
518  gctobj::makeST_hf(temporary_hf, tempST_hf);
519  for (int i = 2 * nJets; i < 3 * nJets; i++) {
520  jet[i] = gctobj::getRegion(tempST_hf);
521  l1tp2::Phase2L1CaloJet tempJet;
522  int hfjeteta = jet[i].etaCenter;
523  int hfjetphi = jet[i].phiCenter;
524  tempJet.setJetIEta(hfjeteta + k * nHfEta / 2);
525  tempJet.setJetIPhi(hfjetphi);
526  float jeteta = hfEta[hfjeteta + k * nHfEta / 2][hfjetphi];
527  float jetphi = hfPhi[hfjeteta + k * nHfEta / 2][hfjetphi];
528  tempJet.setJetEta(jeteta);
529  tempJet.setJetPhi(jetphi);
530  tempJet.setJetEt(get_jet_pt_calibration(jet[i].energy, jeteta));
531  tempJet.setTauEt(get_tau_pt_calibration(jet[i].tauEt, jeteta));
532  tempJet.setTowerEt(jet[i].energyMax);
533  int hftowereta = jet[i].etaMax;
534  int hftowerphi = jet[i].phiMax;
535  tempJet.setTowerIEta(hftowereta + k * nHfEta / 2);
536  tempJet.setTowerIPhi(hftowerphi);
537  float towereta = hfEta[hftowereta + k * nHfEta / 2][hftowerphi];
538  float towerphi = hfPhi[hftowereta + k * nHfEta / 2][hftowerphi];
539  tempJet.setTowerEta(towereta);
540  tempJet.setTowerPhi(towerphi);
542  tempJetp4.SetPt(tempJet.jetEt());
543  tempJetp4.SetEta(tempJet.jetEta());
544  tempJetp4.SetPhi(tempJet.jetPhi());
545  tempJetp4.SetM(0.);
546  tempJet.setP4(tempJetp4);
547 
548  if (jet[i].energy > 0.)
549  halfHfJets.push_back(tempJet);
550  }
551 
552  // Stitching:
553  // 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
554  // 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
555  // in both eta/phi allow a maximum of one tower between jet centers for stitching
556 
557  for (size_t i = 0; i < halfHgcalJets.size(); i++) {
558  if (halfHgcalJets.at(i).jetIEta() >= (nHgcalEta / 2 - 2) && halfHgcalJets.at(i).jetIEta() < (nHgcalEta / 2 + 2)) {
559  float hgcal_ieta = k * nBarrelEta + halfHgcalJets.at(i).jetIEta();
560  for (size_t j = 0; j < halfBarrelJets.size(); j++) {
561  float barrel_ieta = nHgcalEta / 2 + halfBarrelJets.at(j).jetIEta();
562  if (abs(barrel_ieta - hgcal_ieta) <= 2 &&
563  abs(halfBarrelJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) <= 2) {
564  float totalet = halfBarrelJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
565  float totalTauEt = halfBarrelJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
566  if (halfBarrelJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
567  halfHgcalJets.at(i).setJetEt(0.);
568  halfHgcalJets.at(i).setTauEt(0.);
569  halfBarrelJets.at(j).setJetEt(totalet);
570  halfBarrelJets.at(j).setTauEt(totalTauEt);
572  tempJetp4.SetPt(totalet);
573  tempJetp4.SetEta(halfBarrelJets.at(j).jetEta());
574  tempJetp4.SetPhi(halfBarrelJets.at(j).jetPhi());
575  tempJetp4.SetM(0.);
576  halfBarrelJets.at(j).setP4(tempJetp4);
577  } else {
578  halfHgcalJets.at(i).setJetEt(totalet);
579  halfHgcalJets.at(i).setTauEt(totalTauEt);
580  halfBarrelJets.at(j).setJetEt(0.);
581  halfBarrelJets.at(j).setTauEt(0.);
583  tempJetp4.SetPt(totalet);
584  tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
585  tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
586  tempJetp4.SetM(0.);
587  halfHgcalJets.at(i).setP4(tempJetp4);
588  }
589  }
590  }
591  } else if (halfHgcalJets.at(i).jetIEta() < 2 || halfHgcalJets.at(i).jetIEta() >= (nHgcalEta - 2)) {
592  float hgcal_ieta = k * nBarrelEta + nHfEta / 2 + halfHgcalJets.at(i).jetIEta();
593  for (size_t j = 0; j < halfHfJets.size(); j++) {
594  float hf_ieta = k * nBarrelEta + k * nHgcalEta + halfHfJets.at(j).jetIEta();
595  if (abs(hgcal_ieta - hf_ieta) < 3 && abs(halfHfJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) < 3) {
596  float totalet = halfHfJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
597  float totalTauEt = halfHfJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
598  if (halfHfJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
599  halfHgcalJets.at(i).setJetEt(0.);
600  halfHgcalJets.at(i).setTauEt(0.);
601  halfHfJets.at(j).setJetEt(totalet);
602  halfHfJets.at(j).setTauEt(totalTauEt);
604  tempJetp4.SetPt(totalet);
605  tempJetp4.SetEta(halfHfJets.at(j).jetEta());
606  tempJetp4.SetPhi(halfHfJets.at(j).jetPhi());
607  tempJetp4.SetM(0.);
608  halfHfJets.at(j).setP4(tempJetp4);
609  } else {
610  halfHgcalJets.at(i).setJetEt(totalet);
611  halfHgcalJets.at(i).setTauEt(totalTauEt);
612  halfHfJets.at(j).setJetEt(0.);
613  halfHfJets.at(j).setTauEt(0.);
615  tempJetp4.SetPt(totalet);
616  tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
617  tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
618  tempJetp4.SetM(0.);
619  halfHgcalJets.at(i).setP4(tempJetp4);
620  }
621  }
622  }
623  }
624  }
625 
626  // Write up to 6 jets from each eta half of barrel, endcap, forward regions
627 
628  std::sort(halfBarrelJets.begin(), halfBarrelJets.end(), gctobj::compareByEt);
629  for (size_t i = 0; i < halfBarrelJets.size(); i++) {
630  if (halfBarrelJets.at(i).jetEt() > 0. && i < 6)
631  allJets.push_back(halfBarrelJets.at(i));
632  }
633 
634  std::sort(halfHgcalJets.begin(), halfHgcalJets.end(), gctobj::compareByEt);
635  for (size_t i = 0; i < halfHgcalJets.size(); i++) {
636  if (halfHgcalJets.at(i).jetEt() > 0. && i < 6)
637  allJets.push_back(halfHgcalJets.at(i));
638  }
639 
640  std::sort(halfHfJets.begin(), halfHfJets.end(), gctobj::compareByEt);
641  for (size_t i = 0; i < halfHfJets.size(); i++) {
642  if (halfHfJets.at(i).jetEt() > 0. && i < 6)
643  allJets.push_back(halfHfJets.at(i));
644  }
645  }
646 
647  std::sort(allJets.begin(), allJets.end(), gctobj::compareByEt);
648  for (size_t i = 0; i < allJets.size(); i++) {
649  jetCands->push_back(allJets.at(i));
650  }
651 
652  iEvent.put(std::move(jetCands), "GCTJet");
653 }
static float towerEta(int ieta)
Definition: CaloTools.cc:201
std::vector< reco::JetBaseRef > allJets(const Container &)
fill list of all jets associated with values. Return # of jets in the list
std::map< std::string, std::map< std::string, TF1 > > all_nvtx_to_PU_sub_funcs
static constexpr int nBarrelPhi
void setTowerEt(float towerEtIn)
void setTowerIEta(int towerIEtaIn)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static constexpr int nSTEta
T const * product() const
Definition: Handle.h:70
std::map< std::string, TF1 > nHits_to_nvtx_funcs
static constexpr int nHgcalEta
static constexpr int nJets
static const int kHFBegin
Definition: CaloTools.h:39
edm::EDGetTokenT< l1t::HGCalTowerBxCollection > hgcalTowerToken_
const_iterator begin(int bx) const
static constexpr int nHfPhi
float get_tau_pt_calibration(const float &tau_pt, const float &tau_eta) const
int iEvent
Definition: GenABIO.cc:224
void makeST_hf(const float hfTowers[nHfEta/2][nHfPhi], GCTsupertower_t supertower_return[nSTEta][nSTPhi])
static constexpr int nSTPhi
static const int kHFEnd
Definition: CaloTools.h:40
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void setTauEt(float tauEtIn)
void setTowerPhi(float towerPhiIn)
#define M_PI
void setJetIPhi(int jetIPhiIn)
unsigned int id
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecord > decoderTag_
void makeST(const float GCTintTowers[nBarrelEta/2][nBarrelPhi], GCTsupertower_t supertower_return[nSTEta][nSTPhi])
jetInfo getRegion(GCTsupertower_t temp[nSTEta][nSTPhi])
static constexpr int nHgcalPhi
void setTowerEta(float towerEtaIn)
const_iterator end(int bx) const
HLT enums.
void makeST_hgcal(const float hgcalTowers[nHgcalEta/2][nHgcalPhi], GCTsupertower_t supertower_return[nSTEta][nSTPhi])
edm::EDGetTokenT< HcalTrigPrimDigiCollection > hfToken_
static constexpr int nBarrelEta
static constexpr int nHfEta
void setTowerIPhi(int towerIPhiIn)
void setJetEt(float jetEtIn)
void setJetIEta(int jetIEtaIn)
edm::EDGetTokenT< l1tp2::CaloTowerCollection > caloTowerToken_
bool compareByEt(l1tp2::Phase2L1CaloJet i, l1tp2::Phase2L1CaloJet j)
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
void setJetPhi(float jetPhiIn)
float get_jet_pt_calibration(const float &jet_pt, const float &jet_eta) const
void setJetEta(float jetEtaIn)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38

Member Data Documentation

◆ absEtaBinsBarrel

std::vector<double> Phase2L1CaloJetEmulator::absEtaBinsBarrel
private

Definition at line 83 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ absEtaBinsHF

std::vector<double> Phase2L1CaloJetEmulator::absEtaBinsHF
private

Definition at line 87 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ absEtaBinsHGCal

std::vector<double> Phase2L1CaloJetEmulator::absEtaBinsHGCal
private

Definition at line 85 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ all_nvtx_to_PU_sub_funcs

std::map<std::string, std::map<std::string, TF1> > Phase2L1CaloJetEmulator::all_nvtx_to_PU_sub_funcs
private

Definition at line 79 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator(), and produce().

◆ calibrationsBarrel

std::vector<std::vector<double> > Phase2L1CaloJetEmulator::calibrationsBarrel
private

Definition at line 98 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ calibrationsHF

std::vector<std::vector<double> > Phase2L1CaloJetEmulator::calibrationsHF
private

Definition at line 100 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ calibrationsHGCal

std::vector<std::vector<double> > Phase2L1CaloJetEmulator::calibrationsHGCal
private

Definition at line 99 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ caloTowerToken_

edm::EDGetTokenT<l1tp2::CaloTowerCollection> Phase2L1CaloJetEmulator::caloTowerToken_
private

Definition at line 69 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ decoderTag_

edm::ESGetToken<CaloTPGTranscoder, CaloTPGRecord> Phase2L1CaloJetEmulator::decoderTag_
private

Definition at line 72 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ hf_nvtx_to_PU_sub_funcs

std::map<std::string, TF1> Phase2L1CaloJetEmulator::hf_nvtx_to_PU_sub_funcs
private

Definition at line 78 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ hfToken_

edm::EDGetTokenT<HcalTrigPrimDigiCollection> Phase2L1CaloJetEmulator::hfToken_
private

Definition at line 71 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ hgcalEM_nvtx_to_PU_sub_funcs

std::map<std::string, TF1> Phase2L1CaloJetEmulator::hgcalEM_nvtx_to_PU_sub_funcs
private

Definition at line 76 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ hgcalHad_nvtx_to_PU_sub_funcs

std::map<std::string, TF1> Phase2L1CaloJetEmulator::hgcalHad_nvtx_to_PU_sub_funcs
private

Definition at line 77 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ hgcalTowerToken_

edm::EDGetTokenT<l1t::HGCalTowerBxCollection> Phase2L1CaloJetEmulator::hgcalTowerToken_
private

Definition at line 70 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ jetCalibrationsBarrel

std::vector<double> Phase2L1CaloJetEmulator::jetCalibrationsBarrel
private

Definition at line 84 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ jetCalibrationsHF

std::vector<double> Phase2L1CaloJetEmulator::jetCalibrationsHF
private

Definition at line 88 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ jetCalibrationsHGCal

std::vector<double> Phase2L1CaloJetEmulator::jetCalibrationsHGCal
private

Definition at line 86 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ jetPtBins

std::vector<double> Phase2L1CaloJetEmulator::jetPtBins
private

Definition at line 82 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_jet_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ nHits_to_nvtx_funcs

std::map<std::string, TF1> Phase2L1CaloJetEmulator::nHits_to_nvtx_funcs
private

Definition at line 75 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator(), and produce().

◆ nHits_to_nvtx_params

std::vector<edm::ParameterSet> Phase2L1CaloJetEmulator::nHits_to_nvtx_params
private

Definition at line 73 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ nvtx_to_PU_sub_params

std::vector<edm::ParameterSet> Phase2L1CaloJetEmulator::nvtx_to_PU_sub_params
private

Definition at line 74 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ tauAbsEtaBinsBarrel

std::vector<double> Phase2L1CaloJetEmulator::tauAbsEtaBinsBarrel
private

Definition at line 92 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_tau_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ tauAbsEtaBinsHGCal

std::vector<double> Phase2L1CaloJetEmulator::tauAbsEtaBinsHGCal
private

Definition at line 94 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_tau_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ tauCalibrationsBarrel

std::vector<double> Phase2L1CaloJetEmulator::tauCalibrationsBarrel
private

Definition at line 93 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ tauCalibrationsHGCal

std::vector<double> Phase2L1CaloJetEmulator::tauCalibrationsHGCal
private

Definition at line 95 of file Phase2L1CaloJetEmulator.cc.

Referenced by Phase2L1CaloJetEmulator().

◆ tauPtBins

std::vector<double> Phase2L1CaloJetEmulator::tauPtBins
private

Definition at line 91 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_tau_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ tauPtCalibrationsBarrel

std::vector<std::vector<double> > Phase2L1CaloJetEmulator::tauPtCalibrationsBarrel
private

Definition at line 103 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_tau_pt_calibration(), and Phase2L1CaloJetEmulator().

◆ tauPtCalibrationsHGCal

std::vector<std::vector<double> > Phase2L1CaloJetEmulator::tauPtCalibrationsHGCal
private

Definition at line 104 of file Phase2L1CaloJetEmulator.cc.

Referenced by get_tau_pt_calibration(), and Phase2L1CaloJetEmulator().