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

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

edm::EDGetTokenT< l1tp2::CaloTowerCollectioncaloTowerToken_
 
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecorddecoderTag_
 
edm::EDGetTokenT< HcalTrigPrimDigiCollectionhfToken_
 
edm::EDGetTokenT< l1t::HGCalTowerBxCollectionhgcalTowerToken_
 

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 55 of file Phase2L1CaloJetEmulator.cc.

Constructor & Destructor Documentation

◆ Phase2L1CaloJetEmulator()

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

Definition at line 75 of file Phase2L1CaloJetEmulator.cc.

76  : caloTowerToken_(consumes<l1tp2::CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("gctFullTowers"))),
77  hgcalTowerToken_(consumes<l1t::HGCalTowerBxCollection>(iConfig.getParameter<edm::InputTag>("hgcalTowers"))),
78  hfToken_(consumes<HcalTrigPrimDigiCollection>(iConfig.getParameter<edm::InputTag>("hcalDigis"))),
79  decoderTag_(esConsumes<CaloTPGTranscoder, CaloTPGRecord>(edm::ESInputTag("", ""))) {
80  produces<l1tp2::Phase2L1CaloJetCollection>("GCTJet");
81 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< l1t::HGCalTowerBxCollection > hgcalTowerToken_
edm::ESGetToken< CaloTPGTranscoder, CaloTPGRecord > decoderTag_
edm::EDGetTokenT< HcalTrigPrimDigiCollection > hfToken_
edm::EDGetTokenT< l1tp2::CaloTowerCollection > caloTowerToken_

◆ ~Phase2L1CaloJetEmulator()

Phase2L1CaloJetEmulator::~Phase2L1CaloJetEmulator ( )
override

Definition at line 83 of file Phase2L1CaloJetEmulator.cc.

83 {}

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 457 of file Phase2L1CaloJetEmulator.cc.

References edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

457  {
459  desc.add<edm::InputTag>("gctFullTowers", edm::InputTag("l1tPhase2L1CaloEGammaEmulator", "GCTFullTowers"));
460  desc.add<edm::InputTag>("hgcalTowers", edm::InputTag("l1tHGCalTowerProducer", "HGCalTowerProcessor"));
461  desc.add<edm::InputTag>("hcalDigis", edm::InputTag("simHcalTriggerPrimitiveDigis"));
462  descriptions.addWithDefaultLabel(desc);
463 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)

◆ produce()

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

Definition at line 90 of file Phase2L1CaloJetEmulator.cc.

References funct::abs(), reco::JetExtendedAssociation::allJets(), BXVector< T >::begin(), caloTowerToken_, gctobj::compareByEt(), decoderTag_, BXVector< T >::end(), hcalRecHitTable_cff::energy, EgHLTOffHistBins_cfi::et, PVValHelper::eta, edm::EventSetup::getData(), gctobj::getRegion(), hfToken_, l1tPhase2CaloJetEmulator_cfi::hgcalTowers, hgcalTowerToken_, mps_fire::i, hit::id, hcalRecHitTable_cff::ieta, iEvent, hcalRecHitTable_cff::iphi, 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, 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, groupFilesInBlocks::temp, and l1t::CaloTools::towerEta().

90  {
91  using namespace edm;
92  std::unique_ptr<l1tp2::Phase2L1CaloJetCollection> jetCands(make_unique<l1tp2::Phase2L1CaloJetCollection>());
93 
94  // Assign ETs to each eta-half of the barrel region (17x72 --> 18x72 to be able to make 3x3 super towers)
95  edm::Handle<std::vector<l1tp2::CaloTower>> caloTowerCollection;
96  if (!iEvent.getByToken(caloTowerToken_, caloTowerCollection))
97  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from caloTowerCollection!";
98 
99  iEvent.getByToken(caloTowerToken_, caloTowerCollection);
100  float GCTintTowers[nBarrelEta][nBarrelPhi];
101  float realEta[nBarrelEta][nBarrelPhi];
102  float realPhi[nBarrelEta][nBarrelPhi];
103  for (const l1tp2::CaloTower& i : *caloTowerCollection) {
104  int ieta = i.towerIEta();
105  int iphi = i.towerIPhi();
106  if (i.ecalTowerEt() > 1.)
107  GCTintTowers[ieta][iphi] = i.ecalTowerEt(); // suppress <= 1 GeV towers
108  else
109  GCTintTowers[ieta][iphi] = 0;
110  realEta[ieta][iphi] = i.towerEta();
111  realPhi[ieta][iphi] = i.towerPhi();
112  }
113 
114  float temporary[nBarrelEta / 2][nBarrelPhi];
115 
116  // Assign ETs to each eta-half of the endcap region (18x72)
117  edm::Handle<l1t::HGCalTowerBxCollection> hgcalTowerCollection;
118  if (!iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection))
119  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get towers from hgcalTowerCollection!";
120  l1t::HGCalTowerBxCollection hgcalTowerColl;
121  iEvent.getByToken(hgcalTowerToken_, hgcalTowerCollection);
122  hgcalTowerColl = (*hgcalTowerCollection.product());
124  float hgcalEta[nHgcalEta][nHgcalPhi];
125  float hgcalPhi[nHgcalEta][nHgcalPhi];
126 
127  for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
128  for (int ieta = 0; ieta < nHgcalEta; ieta++) {
129  hgcalTowers[ieta][iphi] = 0;
130  if (ieta < nHgcalEta / 2)
131  hgcalEta[ieta][iphi] = -3.045 + ieta * 0.087 + 0.0435;
132  else
133  hgcalEta[ieta][iphi] = 1.479 + (ieta - nHgcalEta / 2) * 0.087 + 0.0435;
134  hgcalPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHgcalPhi) + (M_PI / nHgcalPhi);
135  }
136  }
137 
138  for (auto it = hgcalTowerColl.begin(0); it != hgcalTowerColl.end(0); it++) {
139  float eta = it->eta();
140  int ieta;
141  if (eta < 0)
142  ieta = 19 - it->id().iEta();
143  else
144  ieta = 20 + it->id().iEta();
145  if (eta > 1.479)
146  ieta = ieta - 4;
147  int iphi = it->id().iPhi();
148  if ((it->etEm() + it->etHad() > 1.) && abs(eta) > 1.479)
149  hgcalTowers[ieta][iphi] = it->etEm() + it->etHad(); // suppress <= 1 GeV towers
150  }
151 
152  float temporary_hgcal[nHgcalEta / 2][nHgcalPhi];
153 
154  // Assign ETs to each eta-half of the forward region (12x72)
156  if (!iEvent.getByToken(hfToken_, hfHandle))
157  edm::LogError("Phase2L1CaloJetEmulator") << "Failed to get HcalTrigPrimDigi for HF!";
158  iEvent.getByToken(hfToken_, hfHandle);
159  float hfTowers[nHfEta][nHfPhi];
160  float hfEta[nHfEta][nHfPhi];
161  float hfPhi[nHfEta][nHfPhi];
162 
163  for (int iphi = 0; iphi < nHfPhi; iphi++) {
164  for (int ieta = 0; ieta < nHfEta; ieta++) {
165  hfTowers[ieta][iphi] = 0;
166  int temp;
167  if (ieta < nHfEta / 2)
169  else
172  hfPhi[ieta][iphi] = -M_PI + (iphi * 2 * M_PI / nHfPhi) + (M_PI / nHfPhi);
173  }
174  }
175 
176  const auto& decoder = iSetup.getData(decoderTag_);
177  for (const auto& hit : *hfHandle.product()) {
178  double et = decoder.hcaletValue(hit.id(), hit.t0());
179  int ieta = 0;
180  if (abs(hit.id().ieta()) < l1t::CaloTools::kHFBegin)
181  continue;
182  if (abs(hit.id().ieta()) > l1t::CaloTools::kHFEnd)
183  continue;
184  if (hit.id().ieta() <= -(l1t::CaloTools::kHFBegin + 1)) {
185  ieta = hit.id().ieta() + l1t::CaloTools::kHFEnd;
186  } else if (hit.id().ieta() >= (l1t::CaloTools::kHFBegin + 1)) {
187  ieta = nHfEta / 2 + (hit.id().ieta() - (l1t::CaloTools::kHFBegin + 1));
188  }
189  int iphi = 0;
190  if (hit.id().iphi() <= nHfPhi / 2)
191  iphi = hit.id().iphi() + (nHfPhi / 2 - 1);
192  else if (hit.id().iphi() > nHfPhi / 2)
193  iphi = hit.id().iphi() - (nHfPhi / 2 + 1);
194  if (et > 1.)
195  hfTowers[ieta][iphi] = et; // suppress <= 1 GeV towers
196  }
197 
198  float temporary_hf[nHfEta / 2][nHfPhi];
199 
200  // Begin creating jets
201  // First create 3x3 super towers: 6x24 in barrel, endcap; 4x24 in forward
202  // Then create up to 10 jets in each eta half of barrel, endcap, forward regions
203 
204  vector<l1tp2::Phase2L1CaloJet> halfBarrelJets, halfHgcalJets, halfHfJets;
205  halfBarrelJets.clear();
206  halfHgcalJets.clear();
207  halfHfJets.clear();
208  vector<l1tp2::Phase2L1CaloJet> allJets;
209  allJets.clear();
210 
211  for (int k = 0; k < 2; k++) {
212  halfBarrelJets.clear();
213  halfHgcalJets.clear();
214  halfHfJets.clear();
216 
217  // BARREL
218  for (int iphi = 0; iphi < nBarrelPhi; iphi++) {
219  for (int ieta = 0; ieta < nBarrelEta / 2; ieta++) {
220  if (k == 0)
221  temporary[ieta][iphi] = GCTintTowers[ieta][iphi];
222  else
223  temporary[ieta][iphi] = GCTintTowers[nBarrelEta / 2 + ieta][iphi];
224  }
225  }
226 
228  gctobj::makeST(temporary, tempST);
229 
230  for (int i = 0; i < nJets; i++) {
231  jet[i] = gctobj::getRegion(tempST);
232  l1tp2::Phase2L1CaloJet tempJet;
233  tempJet.setJetEt(jet[i].energy);
234  tempJet.setTauEt(jet[i].tauEt);
235  int gctjeteta = jet[i].etaCenter;
236  int gctjetphi = jet[i].phiCenter;
237  tempJet.setJetIEta(gctjeteta + k * nBarrelEta / 2);
238  tempJet.setJetIPhi(gctjetphi);
239  float jeteta = realEta[gctjeteta + k * nBarrelEta / 2][gctjetphi];
240  float jetphi = realPhi[gctjeteta + k * nBarrelEta / 2][gctjetphi];
241  tempJet.setJetEta(jeteta);
242  tempJet.setJetPhi(jetphi);
243  tempJet.setTowerEt(jet[i].energyMax);
244  int gcttowereta = jet[i].etaMax;
245  int gcttowerphi = jet[i].phiMax;
246  tempJet.setTowerIEta(gcttowereta + k * nBarrelEta / 2);
247  tempJet.setTowerIPhi(gcttowerphi);
248  float towereta = realEta[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
249  float towerphi = realPhi[gcttowereta + k * nBarrelEta / 2][gcttowerphi];
250  tempJet.setTowerEta(towereta);
251  tempJet.setTowerPhi(towerphi);
253  tempJetp4.SetPt(tempJet.jetEt());
254  tempJetp4.SetEta(tempJet.jetEta());
255  tempJetp4.SetPhi(tempJet.jetPhi());
256  tempJetp4.SetM(0.);
257  tempJet.setP4(tempJetp4);
258 
259  if (jet[i].energy > 0.)
260  halfBarrelJets.push_back(tempJet);
261  }
262 
263  // ENDCAP
264  for (int iphi = 0; iphi < nHgcalPhi; iphi++) {
265  for (int ieta = 0; ieta < nHgcalEta / 2; ieta++) {
266  if (k == 0)
267  temporary_hgcal[ieta][iphi] = hgcalTowers[ieta][iphi];
268  else
269  temporary_hgcal[ieta][iphi] = hgcalTowers[nHgcalEta / 2 + ieta][iphi];
270  }
271  }
272 
273  gctobj::GCTsupertower_t tempST_hgcal[nSTEta][nSTPhi];
274  gctobj::makeST_hgcal(temporary_hgcal, tempST_hgcal);
275  for (int i = nJets; i < 2 * nJets; i++) {
276  jet[i] = gctobj::getRegion(tempST_hgcal);
277  l1tp2::Phase2L1CaloJet tempJet;
278  tempJet.setJetEt(jet[i].energy);
279  tempJet.setTauEt(jet[i].tauEt);
280  int hgcaljeteta = jet[i].etaCenter;
281  int hgcaljetphi = jet[i].phiCenter;
282  tempJet.setJetIEta(hgcaljeteta + k * nHgcalEta / 2);
283  tempJet.setJetIPhi(hgcaljetphi);
284  float jeteta = hgcalEta[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
285  float jetphi = hgcalPhi[hgcaljeteta + k * nHgcalEta / 2][hgcaljetphi];
286  tempJet.setJetEta(jeteta);
287  tempJet.setJetPhi(jetphi);
288  tempJet.setTowerEt(jet[i].energyMax);
289  int hgcaltowereta = jet[i].etaMax;
290  int hgcaltowerphi = jet[i].phiMax;
291  tempJet.setTowerIEta(hgcaltowereta + k * nHgcalEta / 2);
292  tempJet.setTowerIPhi(hgcaltowerphi);
293  float towereta = hgcalEta[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
294  float towerphi = hgcalPhi[hgcaltowereta + k * nHgcalEta / 2][hgcaltowerphi];
295  tempJet.setTowerEta(towereta);
296  tempJet.setTowerPhi(towerphi);
298  tempJetp4.SetPt(tempJet.jetEt());
299  tempJetp4.SetEta(tempJet.jetEta());
300  tempJetp4.SetPhi(tempJet.jetPhi());
301  tempJetp4.SetM(0.);
302  tempJet.setP4(tempJetp4);
303 
304  if (jet[i].energy > 0.)
305  halfHgcalJets.push_back(tempJet);
306  }
307 
308  // HF
309  for (int iphi = 0; iphi < nHfPhi; iphi++) {
310  for (int ieta = 0; ieta < nHfEta / 2; ieta++) {
311  if (k == 0)
312  temporary_hf[ieta][iphi] = hfTowers[ieta][iphi];
313  else
314  temporary_hf[ieta][iphi] = hfTowers[nHfEta / 2 + ieta][iphi];
315  }
316  }
317 
319  gctobj::makeST_hf(temporary_hf, tempST_hf);
320  for (int i = 2 * nJets; i < 3 * nJets; i++) {
321  jet[i] = gctobj::getRegion(tempST_hf);
322  l1tp2::Phase2L1CaloJet tempJet;
323  tempJet.setJetEt(jet[i].energy);
324  tempJet.setTauEt(jet[i].tauEt);
325  int hfjeteta = jet[i].etaCenter;
326  int hfjetphi = jet[i].phiCenter;
327  tempJet.setJetIEta(hfjeteta + k * nHfEta / 2);
328  tempJet.setJetIPhi(hfjetphi);
329  float jeteta = hfEta[hfjeteta + k * nHfEta / 2][hfjetphi];
330  float jetphi = hfPhi[hfjeteta + k * nHfEta / 2][hfjetphi];
331  tempJet.setJetEta(jeteta);
332  tempJet.setJetPhi(jetphi);
333  tempJet.setTowerEt(jet[i].energyMax);
334  int hftowereta = jet[i].etaMax;
335  int hftowerphi = jet[i].phiMax;
336  tempJet.setTowerIEta(hftowereta + k * nHfEta / 2);
337  tempJet.setTowerIPhi(hftowerphi);
338  float towereta = hfEta[hftowereta + k * nHfEta / 2][hftowerphi];
339  float towerphi = hfPhi[hftowereta + k * nHfEta / 2][hftowerphi];
340  tempJet.setTowerEta(towereta);
341  tempJet.setTowerPhi(towerphi);
343  tempJetp4.SetPt(tempJet.jetEt());
344  tempJetp4.SetEta(tempJet.jetEta());
345  tempJetp4.SetPhi(tempJet.jetPhi());
346  tempJetp4.SetM(0.);
347  tempJet.setP4(tempJetp4);
348 
349  if (jet[i].energy > 0.)
350  halfHfJets.push_back(tempJet);
351  }
352 
353  // Stitching:
354  // 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
355  // 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
356  // in both eta/phi allow a maximum of one tower between jet centers for stitching
357 
358  for (size_t i = 0; i < halfHgcalJets.size(); i++) {
359  if (halfHgcalJets.at(i).jetIEta() >= (nHgcalEta / 2 - 2) && halfHgcalJets.at(i).jetIEta() < (nHgcalEta / 2 + 2)) {
360  float hgcal_ieta = k * nBarrelEta + halfHgcalJets.at(i).jetIEta();
361  for (size_t j = 0; j < halfBarrelJets.size(); j++) {
362  float barrel_ieta = nHgcalEta / 2 + halfBarrelJets.at(j).jetIEta();
363  if (abs(barrel_ieta - hgcal_ieta) <= 2 &&
364  abs(halfBarrelJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) <= 2) {
365  float totalet = halfBarrelJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
366  float totalTauEt = halfBarrelJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
367  if (halfBarrelJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
368  halfHgcalJets.at(i).setJetEt(0.);
369  halfHgcalJets.at(i).setTauEt(0.);
370  halfBarrelJets.at(j).setJetEt(totalet);
371  halfBarrelJets.at(j).setTauEt(totalTauEt);
373  tempJetp4.SetPt(totalet);
374  tempJetp4.SetEta(halfBarrelJets.at(j).jetEta());
375  tempJetp4.SetPhi(halfBarrelJets.at(j).jetPhi());
376  tempJetp4.SetM(0.);
377  halfBarrelJets.at(j).setP4(tempJetp4);
378  } else {
379  halfHgcalJets.at(i).setJetEt(totalet);
380  halfHgcalJets.at(i).setTauEt(totalTauEt);
381  halfBarrelJets.at(j).setJetEt(0.);
382  halfBarrelJets.at(j).setTauEt(0.);
384  tempJetp4.SetPt(totalet);
385  tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
386  tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
387  tempJetp4.SetM(0.);
388  halfHgcalJets.at(i).setP4(tempJetp4);
389  }
390  }
391  }
392  } else if (halfHgcalJets.at(i).jetIEta() < 2 || halfHgcalJets.at(i).jetIEta() >= (nHgcalEta - 2)) {
393  float hgcal_ieta = k * nBarrelEta + nHfEta / 2 + halfHgcalJets.at(i).jetIEta();
394  for (size_t j = 0; j < halfHfJets.size(); j++) {
395  float hf_ieta = k * nBarrelEta + k * nHgcalEta + halfHfJets.at(j).jetIEta();
396  if (abs(hgcal_ieta - hf_ieta) < 3 && abs(halfHfJets.at(j).jetIPhi() - halfHgcalJets.at(i).jetIPhi()) < 3) {
397  float totalet = halfHfJets.at(j).jetEt() + halfHgcalJets.at(i).jetEt();
398  float totalTauEt = halfHfJets.at(j).tauEt() + halfHgcalJets.at(i).tauEt();
399  if (halfHfJets.at(j).jetEt() > halfHgcalJets.at(i).jetEt()) {
400  halfHgcalJets.at(i).setJetEt(0.);
401  halfHgcalJets.at(i).setTauEt(0.);
402  halfHfJets.at(j).setJetEt(totalet);
403  halfHfJets.at(j).setTauEt(totalTauEt);
405  tempJetp4.SetPt(totalet);
406  tempJetp4.SetEta(halfHfJets.at(j).jetEta());
407  tempJetp4.SetPhi(halfHfJets.at(j).jetPhi());
408  tempJetp4.SetM(0.);
409  halfHfJets.at(j).setP4(tempJetp4);
410  } else {
411  halfHgcalJets.at(i).setJetEt(totalet);
412  halfHgcalJets.at(i).setTauEt(totalTauEt);
413  halfHfJets.at(j).setJetEt(0.);
414  halfHfJets.at(j).setTauEt(0.);
416  tempJetp4.SetPt(totalet);
417  tempJetp4.SetEta(halfHgcalJets.at(i).jetEta());
418  tempJetp4.SetPhi(halfHgcalJets.at(i).jetPhi());
419  tempJetp4.SetM(0.);
420  halfHgcalJets.at(i).setP4(tempJetp4);
421  }
422  }
423  }
424  }
425  }
426 
427  // Write up to 6 jets from each eta half of barrel, endcap, forward regions
428 
429  std::sort(halfBarrelJets.begin(), halfBarrelJets.end(), gctobj::compareByEt);
430  for (size_t i = 0; i < halfBarrelJets.size(); i++) {
431  if (halfBarrelJets.at(i).jetEt() > 0. && i < 6)
432  allJets.push_back(halfBarrelJets.at(i));
433  }
434 
435  std::sort(halfHgcalJets.begin(), halfHgcalJets.end(), gctobj::compareByEt);
436  for (size_t i = 0; i < halfHgcalJets.size(); i++) {
437  if (halfHgcalJets.at(i).jetEt() > 0. && i < 6)
438  allJets.push_back(halfHgcalJets.at(i));
439  }
440 
441  std::sort(halfHfJets.begin(), halfHfJets.end(), gctobj::compareByEt);
442  for (size_t i = 0; i < halfHfJets.size(); i++) {
443  if (halfHfJets.at(i).jetEt() > 0. && i < 6)
444  allJets.push_back(halfHfJets.at(i));
445  }
446  }
447 
448  std::sort(allJets.begin(), allJets.end(), gctobj::compareByEt);
449  for (size_t i = 0; i < allJets.size(); i++) {
450  jetCands->push_back(allJets.at(i));
451  }
452 
453  iEvent.put(std::move(jetCands), "GCTJet");
454 }
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
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
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
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)
void setJetEta(float jetEtaIn)
math::PtEtaPhiMLorentzVector PolarLorentzVector
Lorentz vector.
Definition: Candidate.h:38

Member Data Documentation

◆ caloTowerToken_

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

Definition at line 66 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ decoderTag_

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

Definition at line 69 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ hfToken_

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

Definition at line 68 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().

◆ hgcalTowerToken_

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

Definition at line 67 of file Phase2L1CaloJetEmulator.cc.

Referenced by produce().