CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
l1t::GenToInputProducer Class Reference
Inheritance diagram for l1t::GenToInputProducer:

Public Member Functions

 GenToInputProducer (const ParameterSet &)
 
 ~GenToInputProducer () override
 

Static Public Member Functions

static void fillDescriptions (ConfigurationDescriptions &descriptions)
 

Private Member Functions

void beginJob () override
 
void beginRun (Run const &iR, EventSetup const &iE) override
 
int convertEtaToHW (double ieta, double minEta, double maxEta, int steps)
 
int convertPhiToHW (double iphi, int steps)
 
int convertPtToHW (double ipt, int maxPt, double step)
 
void endJob () override
 
void endRun (Run const &iR, EventSetup const &iE) override
 
void produce (Event &, EventSetup const &) override
 

Private Attributes

int bxFirst_
 
int bxLast_
 
std::vector< l1t::EGammaegammaVec_bx0
 
std::vector< l1t::EGammaegammaVec_bxm1
 
std::vector< l1t::EGammaegammaVec_bxm2
 
std::vector< l1t::EGammaegammaVec_bxp1
 
double egEtThreshold_
 
int emptyBxEvt_
 
int emptyBxTrailer_
 
std::vector< l1t::EtSumetsumVec_bx0
 
std::vector< l1t::EtSumetsumVec_bxm1
 
std::vector< l1t::EtSumetsumVec_bxm2
 
std::vector< l1t::EtSumetsumVec_bxp1
 
int eventCnt_
 
GlobalExtBlk extCond_bx0
 
GlobalExtBlk extCond_bxm1
 
GlobalExtBlk extCond_bxm2
 
GlobalExtBlk extCond_bxp1
 
edm::EDGetTokenT< reco::GenJetCollectiongenJetsToken
 
edm::EDGetTokenT< reco::GenMETCollectiongenMetToken
 
edm::EDGetTokenT< reco::GenParticleCollectiongenParticlesToken
 
TRandom3 * gRandom
 
double jetEtThreshold_
 
std::vector< l1t::JetjetVec_bx0
 
std::vector< l1t::JetjetVec_bxm1
 
std::vector< l1t::JetjetVec_bxm2
 
std::vector< l1t::JetjetVec_bxp1
 
unsigned long long m_paramsCacheId
 
int maxNumEGCands_
 
int maxNumJetCands_
 
int maxNumMuCands_
 
int maxNumMuShowerCands_
 
int maxNumTauCands_
 
double muEtThreshold_
 
std::vector< l1t::MuonShowermuonShowerVec_bx0
 
std::vector< l1t::MuonShowermuonShowerVec_bxm1
 
std::vector< l1t::MuonShowermuonShowerVec_bxm2
 
std::vector< l1t::MuonShowermuonShowerVec_bxp1
 
std::vector< l1t::MuonmuonVec_bx0
 
std::vector< l1t::MuonmuonVec_bxm1
 
std::vector< l1t::MuonmuonVec_bxm2
 
std::vector< l1t::MuonmuonVec_bxp1
 
TRandom3 * mushowerRandom
 
double tauEtThreshold_
 
std::vector< l1t::TautauVec_bx0
 
std::vector< l1t::TautauVec_bxm1
 
std::vector< l1t::TautauVec_bxm2
 
std::vector< l1t::TautauVec_bxp1
 

Detailed Description

Description: Create Input Collections for the GT from MC gen particles. Allows testing of emulation.

Author
: D. Puigh OSU

Modeled after FakeInputProducer.cc

features: R. Cavanaugh

Definition at line 64 of file GenToInputProducer.cc.

Constructor & Destructor Documentation

◆ GenToInputProducer()

l1t::GenToInputProducer::GenToInputProducer ( const ParameterSet iConfig)
explicit

Definition at line 155 of file GenToInputProducer.cc.

References edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

155  {
156  // register what you produce
157  produces<BXVector<l1t::EGamma>>();
158  produces<BXVector<l1t::Muon>>();
159  produces<BXVector<l1t::MuonShower>>();
160  produces<BXVector<l1t::Tau>>();
161  produces<BXVector<l1t::Jet>>();
162  produces<BXVector<l1t::EtSum>>();
163  produces<GlobalExtBlkBxCollection>();
164 
165  // Setup parameters
166  bxFirst_ = iConfig.getParameter<int>("bxFirst");
167  bxLast_ = iConfig.getParameter<int>("bxLast");
168 
169  maxNumMuCands_ = iConfig.getParameter<int>("maxMuCand");
170  maxNumMuShowerCands_ = iConfig.getParameter<int>("maxMuShowerCand");
171  maxNumJetCands_ = iConfig.getParameter<int>("maxJetCand");
172  maxNumEGCands_ = iConfig.getParameter<int>("maxEGCand");
173  maxNumTauCands_ = iConfig.getParameter<int>("maxTauCand");
174 
175  jetEtThreshold_ = iConfig.getParameter<double>("jetEtThreshold");
176  tauEtThreshold_ = iConfig.getParameter<double>("tauEtThreshold");
177  egEtThreshold_ = iConfig.getParameter<double>("egEtThreshold");
178  muEtThreshold_ = iConfig.getParameter<double>("muEtThreshold");
179 
180  emptyBxTrailer_ = iConfig.getParameter<int>("emptyBxTrailer");
181  emptyBxEvt_ = iConfig.getParameter<int>("emptyBxEvt");
182 
183  genParticlesToken = consumes<reco::GenParticleCollection>(std::string("genParticles"));
184  genJetsToken = consumes<reco::GenJetCollection>(std::string("ak4GenJets"));
185  genMetToken = consumes<reco::GenMETCollection>(std::string("genMetCalo"));
186 
187  // set cache id to zero, will be set at first beginRun:
188  m_paramsCacheId = 0;
189  eventCnt_ = 0;
190  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned long long m_paramsCacheId
edm::EDGetTokenT< reco::GenJetCollection > genJetsToken
edm::EDGetTokenT< reco::GenMETCollection > genMetToken
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken

◆ ~GenToInputProducer()

l1t::GenToInputProducer::~GenToInputProducer ( )
override

Definition at line 192 of file GenToInputProducer.cc.

192 {}

Member Function Documentation

◆ beginJob()

void l1t::GenToInputProducer::beginJob ( void  )
overrideprivate

Definition at line 832 of file GenToInputProducer.cc.

832 {}

◆ beginRun()

void l1t::GenToInputProducer::beginRun ( Run const &  iR,
EventSetup const &  iE 
)
overrideprivate

Definition at line 839 of file GenToInputProducer.cc.

References LogDebug.

839  {
840  LogDebug("GtGenToInputProducer") << "GenToInputProducer::beginRun function called...\n";
841 
842  srand(0);
843 
844  gRandom = new TRandom3();
845  mushowerRandom = new TRandom3();
846  }
#define LogDebug(id)

◆ convertEtaToHW()

int l1t::GenToInputProducer::convertEtaToHW ( double  ieta,
double  minEta,
double  maxEta,
int  steps 
)
private

Definition at line 863 of file GenToInputProducer.cc.

References hcalRecHitTable_cff::ieta, createfilelist::int, razorScouting_cff::maxEta, EgHLTOffEleSelection_cfi::minEta, and relval_machine::steps.

863  {
864  double binWidth = (maxEta - minEta) / steps;
865 
866  //if we are outside the limits, set error
867  if (ieta < minEta)
868  return 99999; //ieta = minEta+binWidth/2.;
869  if (ieta > maxEta)
870  return 99999; //ieta = maxEta-binWidth/2.;
871 
872  int binNum = (int)(ieta / binWidth);
873  if (ieta < 0.)
874  binNum--;
875 
876  // unsigned int hwEta = binNum & bitMask;
877  // Remove masking for BXVectors...only assume in raw data
878 
879  return binNum;
880  }

◆ convertPhiToHW()

int l1t::GenToInputProducer::convertPhiToHW ( double  iphi,
int  steps 
)
private

Definition at line 852 of file GenToInputProducer.cc.

References l1trig_cff::hwPhi, createfilelist::int, hcalRecHitTable_cff::iphi, M_PI, AlignmentTrackSelector_cfi::phiMax, and relval_machine::steps.

852  {
853  double phiMax = 2 * M_PI;
854  if (iphi < 0)
855  iphi += 2 * M_PI;
856  if (iphi > phiMax)
857  iphi -= phiMax;
858 
859  int hwPhi = int((iphi / phiMax) * steps + 0.00001);
860  return hwPhi;
861  }
#define M_PI

◆ convertPtToHW()

int l1t::GenToInputProducer::convertPtToHW ( double  ipt,
int  maxPt,
double  step 
)
private

Definition at line 882 of file GenToInputProducer.cc.

References l1trig_cff::hwPt, createfilelist::int, and PV_cfg::maxPt.

882  {
883  int hwPt = int(ipt / step + 0.0001);
884  // if above max Pt, set to largest value
885  if (hwPt > maxPt)
886  hwPt = maxPt;
887 
888  return hwPt;
889  }
maxPt
Definition: PV_cfg.py:224
step
Definition: StallMonitor.cc:83

◆ endJob()

void l1t::GenToInputProducer::endJob ( void  )
overrideprivate

Definition at line 835 of file GenToInputProducer.cc.

835 {}

◆ endRun()

void l1t::GenToInputProducer::endRun ( Run const &  iR,
EventSetup const &  iE 
)
overrideprivate

Definition at line 849 of file GenToInputProducer.cc.

849 {}

◆ fillDescriptions()

void l1t::GenToInputProducer::fillDescriptions ( ConfigurationDescriptions descriptions)
static

Definition at line 892 of file GenToInputProducer.cc.

References edm::ConfigurationDescriptions::addDefault(), and submitPVResolutionJobs::desc.

892  {
893  //The following says we do not know what parameters are allowed so do no validation
894  // Please change this to state exactly what you do use, even if it is no parameters
896  desc.setUnknown();
897  descriptions.addDefault(desc);
898  }
void addDefault(ParameterSetDescription const &psetDescription)

◆ produce()

void l1t::GenToInputProducer::produce ( Event iEvent,
EventSetup const &  iSetup 
)
overrideprivate

Definition at line 199 of file GenToInputProducer.cc.

References funct::abs(), nano_mu_local_reco_cff::bool, hcal_dqm_sourceclient-live_cfg::bxFirst, hcal_dqm_sourceclient-live_cfg::bxLast, ALCARECOTkAlJpsiMuMu_cff::charge, reco::Candidate::charge(), gather_cfg::cout, debug, HLT_2024v14_cff::dEta, PVValHelper::eta, reco::Candidate::eta(), l1trig_cff::etaAtVtx, l1tGTTFileWriter_cfi::etmiss, l1tCaloJetHTTProducer_cfi::genJets, HepMCValidationHelper::genMet(), AJJGenJetFilter_cfi::genParticles, l1tGTTFileWriter_cfi::htmiss, l1trig_cff::hwEtaAtVtx, l1trig_cff::hwPhiAtVtx, l1trig_cff::hwPtUnconstrained, mps_fire::i, iEvent, createfilelist::int, genparticles_cff::iso, objects.IsoTrackAnalyzer::isoSum, metsig::jet, PDWG_EXODelayedJetMET_cff::jets, isotrackApplyRegressor::k, L1Analysis::kAsymEt, L1Analysis::kAsymEtHF, L1Analysis::kAsymHt, L1Analysis::kAsymHtHF, L1Analysis::kCentrality, L1Analysis::kMinBiasHFM0, L1Analysis::kMinBiasHFM1, L1Analysis::kMinBiasHFP0, L1Analysis::kMinBiasHFP1, L1Analysis::kMissingEt, L1Analysis::kMissingEtHF, L1Analysis::kMissingHt, L1Analysis::kMissingHtHF, L1Analysis::kTotalEt, L1Analysis::kTotalEtEm, L1Analysis::kTotalHt, L1Analysis::kTowerCount, LogDebug, LogTrace, eostools::move(), amptDefaultParameters_cff::mu, DiMuonV_cfg::muons, ecaldqm::nTowers, EgammaValidation_cff::pdgId, reco::Candidate::pdgId(), reco::Candidate::phi(), l1trig_cff::phiAtVtx, DiDispStaMuonMonitor_cfi::pt, reco::Candidate::pt(), l1trig_cff::ptUnconstrained, GlobalExtBlk::reset(), GlobalExtBlk::setExternalDecision(), l1t::MuonShower::setMusOutOfTime0(), l1t::MuonShower::setMusOutOfTime1(), edm::shift, mps_update::status, reco::Candidate::status(), EcalPhiSymFlatTableProducers_cfi::sumEt, makeGlobalPositionRcd_cfg::tag, metsig::tau, and Tau3MuMonitor_cff::taus.

199  {
200  eventCnt_++;
201 
202  LogDebug("GtGenToInputProducer") << "GenToInputProducer::produce function called...\n";
203 
204  // Setup vectors
205  std::vector<l1t::Muon> muonVec;
206  std::vector<l1t::MuonShower> muonShowerVec;
207  std::vector<l1t::EGamma> egammaVec;
208  std::vector<l1t::Tau> tauVec;
209  std::vector<l1t::Jet> jetVec;
210  std::vector<l1t::EtSum> etsumVec;
211  GlobalExtBlk extCond_bx;
212 
213  // Set the range of BX....TO DO...move to Params or determine from param set.
214  int bxFirst = bxFirst_;
215  int bxLast = bxLast_;
216 
217  // Default values objects
218  double MaxLepPt_ = 255;
219  double MaxJetPt_ = 1023;
220  double MaxEt_ = 2047;
221 
222  double MaxCaloEta_ = 5.0;
223  double MaxMuonEta_ = 2.45;
224 
225  double PhiStepCalo_ = 144;
226  double PhiStepMuon_ = 576;
227 
228  // eta scale
229  double EtaStepCalo_ = 230;
230  double EtaStepMuon_ = 450;
231 
232  // Et scale (in GeV)
233  double PtStep_ = 0.5;
234 
235  //outputs
236  std::unique_ptr<l1t::EGammaBxCollection> egammas(new l1t::EGammaBxCollection(0, bxFirst, bxLast));
237  std::unique_ptr<l1t::MuonBxCollection> muons(new l1t::MuonBxCollection(0, bxFirst, bxLast));
238  std::unique_ptr<l1t::MuonShowerBxCollection> muonShowers(new l1t::MuonShowerBxCollection(0, bxFirst, bxLast));
239  std::unique_ptr<l1t::TauBxCollection> taus(new l1t::TauBxCollection(0, bxFirst, bxLast));
240  std::unique_ptr<l1t::JetBxCollection> jets(new l1t::JetBxCollection(0, bxFirst, bxLast));
241  std::unique_ptr<l1t::EtSumBxCollection> etsums(new l1t::EtSumBxCollection(0, bxFirst, bxLast));
242  std::unique_ptr<GlobalExtBlkBxCollection> extCond(new GlobalExtBlkBxCollection(0, bxFirst, bxLast));
243 
244  std::vector<int> mu_cands_index;
245  std::vector<int> eg_cands_index;
246  std::vector<int> tau_cands_index;
248  // Make sure that you can get genParticles
249  if (iEvent.getByToken(genParticlesToken, genParticles)) {
250  for (size_t k = 0; k < genParticles->size(); k++) {
251  const reco::Candidate& mcParticle = (*genParticles)[k];
252 
253  int status = mcParticle.status();
254  int pdgId = mcParticle.pdgId();
255  double pt = mcParticle.pt();
256 
257  // Only use status 1 particles (Tau's need to be allowed through..take status 2 taus)
258  if (status != 1 && !(abs(pdgId) == 15 && status == 2))
259  continue;
260 
261  int absId = abs(pdgId);
262 
263  if (absId == 11 && pt >= egEtThreshold_)
264  eg_cands_index.push_back(k);
265  else if (absId == 13 && pt >= muEtThreshold_)
266  mu_cands_index.push_back(k);
267  else if (absId == 15 && pt >= tauEtThreshold_)
268  tau_cands_index.push_back(k);
269  }
270  } else {
271  LogTrace("GtGenToInputProducer") << ">>> GenParticles collection not found!" << std::endl;
272  }
273 
274  // Muon Shower Collection
275 
276  bool mus0 = (bool)mushowerRandom->Integer(2); // should be [0,1] = 1 bit;
277  bool mus1 = (bool)mushowerRandom->Integer(2);
278  bool mus2 = (bool)mushowerRandom->Integer(2);
279  bool mus2Loose = (bool)mushowerRandom->Integer(2);
280  bool musoot0 = (bool)mushowerRandom->Integer(2);
281  bool musoot1 = (bool)mushowerRandom->Integer(2);
282  bool musoot2Loose = (bool)mushowerRandom->Integer(2);
283  //fill a vector of MuonShower objs with only one obj per BX, not one obj per muon obj
284  cout << "GenToInputProducer MuonShower = (MUS0, MUS1, MUS2, MUSOOT0, MUSOOT1) = (" << mus0 << "," << mus1 << ","
285  << mus2 << "," << musoot0 << "," << musoot1 << ")" << endl;
286  l1t::MuonShower muShower(mus0, musoot0, mus2Loose, musoot2Loose, mus1, musoot1, mus2);
287  muShower.setMusOutOfTime0(musoot0);
288  muShower.setMusOutOfTime1(musoot1);
289  muonShowerVec.push_back(muShower);
290 
291  // Muon Collection
292  int numMuCands = int(mu_cands_index.size());
293  Int_t idxMu[numMuCands];
294  double muPtSorted[numMuCands];
295  for (int iMu = 0; iMu < numMuCands; iMu++)
296  muPtSorted[iMu] = genParticles->at(mu_cands_index[iMu]).pt();
297 
298  TMath::Sort(numMuCands, muPtSorted, idxMu);
299  for (int iMu = 0; iMu < numMuCands; iMu++) {
300  if (iMu >= maxNumMuCands_)
301  continue;
302 
303  const reco::Candidate& mcParticle = (*genParticles)[mu_cands_index[idxMu[iMu]]];
304 
305  int pt = convertPtToHW(mcParticle.pt(), MaxLepPt_, PtStep_);
306  int eta = convertEtaToHW(mcParticle.eta(), -MaxMuonEta_, MaxMuonEta_, EtaStepMuon_);
307  int phi = convertPhiToHW(mcParticle.phi(), PhiStepMuon_);
308  int qual = gRandom->Integer(16); //4;
309  int iso = gRandom->Integer(4) % 2; //1;
310  int charge = (mcParticle.charge() < 0) ? 1 : 0;
311  int chargeValid = 1;
312  int tfMuIdx = 0;
313  int tag = 1;
314  bool debug = false;
315  int isoSum = 0;
316  int dPhi = 0;
317  int dEta = 0;
318  int rank = 0;
319  int hwEtaAtVtx = eta;
320  int hwPhiAtVtx = phi;
321  double etaAtVtx = 0.0;
322  double phiAtVtx = 0.0;
323  int hwPtUnconstrained =
324  convertPtToHW(mcParticle.pt(), MaxLepPt_, PtStep_) / 2; // word is 8 bits wide so divide 9 bit word by 2
325  double ptUnconstrained = 0.0;
326  int dXY = gRandom->Integer(4); // should be [0,3] = 2 bits
327 
328  // Eta outside of acceptance
329  if (eta >= 9999)
330  continue;
331 
332  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>* p4 =
333  new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>();
334 
335  l1t::Muon mu(*p4,
336  pt,
337  eta,
338  phi,
339  qual,
340  charge,
341  chargeValid,
342  iso,
343  tfMuIdx,
344  tag,
345  debug,
346  isoSum,
347  dPhi,
348  dEta,
349  rank,
350  hwEtaAtVtx,
351  hwPhiAtVtx,
352  etaAtVtx,
353  phiAtVtx,
356  dXY); // modified to conform to latest Muon.h interface
357  muonVec.push_back(mu);
358  }
359 
360  // EG Collection
361  int numEgCands = int(eg_cands_index.size());
362  Int_t idxEg[numEgCands];
363  double egPtSorted[numEgCands];
364  for (int iEg = 0; iEg < numEgCands; iEg++)
365  egPtSorted[iEg] = genParticles->at(eg_cands_index[iEg]).pt();
366 
367  TMath::Sort(numEgCands, egPtSorted, idxEg);
368  for (int iEg = 0; iEg < numEgCands; iEg++) {
369  if (iEg >= maxNumEGCands_)
370  continue;
371 
372  const reco::Candidate& mcParticle = (*genParticles)[eg_cands_index[idxEg[iEg]]];
373 
374  int pt = convertPtToHW(mcParticle.pt(), MaxLepPt_, PtStep_);
375  int eta = convertEtaToHW(mcParticle.eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_);
376  int phi = convertPhiToHW(mcParticle.phi(), PhiStepCalo_);
377  int qual = gRandom->Integer(2); // modified for LLP Jets
378  int iso = gRandom->Integer(4) % 2;
379 
380  // Eta outside of acceptance
381  if (eta >= 9999)
382  continue;
383 
384  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>* p4 =
385  new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>();
386 
387  l1t::EGamma eg(*p4, pt, eta, phi, qual, iso);
388  egammaVec.push_back(eg);
389  }
390 
391  // Tau Collection
392  int numTauCands = int(tau_cands_index.size());
393  Int_t idxTau[numTauCands];
394  double tauPtSorted[numTauCands];
395  for (int iTau = 0; iTau < numTauCands; iTau++)
396  tauPtSorted[iTau] = genParticles->at(tau_cands_index[iTau]).pt();
397 
398  TMath::Sort(numTauCands, tauPtSorted, idxTau);
399  for (int iTau = 0; iTau < numTauCands; iTau++) {
400  if (iTau >= maxNumTauCands_)
401  continue;
402 
403  const reco::Candidate& mcParticle = (*genParticles)[tau_cands_index[idxTau[iTau]]];
404 
405  int pt = convertPtToHW(mcParticle.pt(), MaxLepPt_, PtStep_);
406  int eta = convertEtaToHW(mcParticle.eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_);
407  int phi = convertPhiToHW(mcParticle.phi(), PhiStepCalo_);
408  int qual = gRandom->Integer(2); // modified for LLP Jets
409  int iso = gRandom->Integer(4) % 2;
410 
411  // Eta outside of acceptance
412  if (eta >= 9999)
413  continue;
414 
415  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>* p4 =
416  new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>();
417 
418  l1t::Tau tau(*p4, pt, eta, phi, qual, iso);
419  tauVec.push_back(tau);
420  }
421 
422  // Temporary hack to increase number of EGs and taus
423  int maxOtherEGs = 4;
424  int maxOtherTaus = 8;
425  int numCurrentEGs = int(egammaVec.size());
426  int numCurrentTaus = int(tauVec.size());
427 
428  int numExtraEGs = 0, numExtraTaus = 0;
429  // end hack
430 
431  // Use to sum the energy of the objects in the event for ETT and HTT
432  // sum all jets
433  double sumEt = 0;
434 
435  int nJet = 0;
437  // Make sure that you can get genJets
438  if (iEvent.getByToken(genJetsToken, genJets)) { // Jet Collection
439  for (reco::GenJetCollection::const_iterator genJet = genJets->begin(); genJet != genJets->end(); ++genJet) {
440  //Keep running sum of total Et
441  sumEt += genJet->et();
442 
443  // Apply pt and eta cut?
444  if (genJet->pt() < jetEtThreshold_)
445  continue;
446 
447  //
448  if (nJet >= maxNumJetCands_)
449  continue;
450  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>* p4 =
451  new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>();
452 
453  int pt = convertPtToHW(genJet->et(), MaxJetPt_, PtStep_);
454  int eta = convertEtaToHW(genJet->eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_);
455  int phi = convertPhiToHW(genJet->phi(), PhiStepCalo_);
456 
457  // Eta outside of acceptance
458  if (eta >= 9999)
459  continue;
460 
461  int qual = gRandom->Integer(2); // modified for LLP Jets
462 
463  l1t::Jet jet(*p4, pt, eta, phi, qual);
464  jetVec.push_back(jet);
465 
466  nJet++;
467 
468  // Temporary hack to increase number of EGs and taus
469  if ((numExtraEGs + numCurrentEGs) < maxNumEGCands_ && numExtraEGs < maxOtherEGs) {
470  numExtraEGs++;
471 
472  int EGpt = convertPtToHW(genJet->et(), MaxLepPt_, PtStep_);
473  int EGeta = convertEtaToHW(genJet->eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_);
474  int EGphi = convertPhiToHW(genJet->phi(), PhiStepCalo_);
475 
476  int EGqual = gRandom->Integer(2); // modified for LLP Jets
477  int EGiso = gRandom->Integer(4) % 2;
478 
479  l1t::EGamma eg(*p4, EGpt, EGeta, EGphi, EGqual, EGiso);
480  egammaVec.push_back(eg);
481  }
482 
483  if ((numExtraTaus + numCurrentTaus) < maxNumTauCands_ && numExtraTaus < maxOtherTaus) {
484  numExtraTaus++;
485 
486  int Taupt = convertPtToHW(genJet->et(), MaxLepPt_, PtStep_);
487  int Taueta = convertEtaToHW(genJet->eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_);
488  int Tauphi = convertPhiToHW(genJet->phi(), PhiStepCalo_);
489  int Tauqual = gRandom->Integer(2); // modified for LLP Jets
490  int Tauiso = gRandom->Integer(4) % 2;
491 
492  l1t::Tau tau(*p4, Taupt, Taueta, Tauphi, Tauqual, Tauiso);
493  tauVec.push_back(tau);
494  }
495  // end hack
496  }
497  } else {
498  LogTrace("GtGenToInputProducer") << ">>> GenJets collection not found!" << std::endl;
499  }
500 
501  // Put the total Et into EtSums (Make HTT slightly smaller to tell them apart....not supposed to be realistic)
502  int pt = convertPtToHW(sumEt, 2047, PtStep_);
503  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>* p4 =
504  new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double>>();
505  l1t::EtSum etTotal(*p4, l1t::EtSum::EtSumType::kTotalEt, pt, 0, 0, 0);
506 
507  // Scale down ETTem as an estimate
508  pt = convertPtToHW(sumEt * 0.6, 2047, PtStep_);
509  l1t::EtSum etEmTotal(*p4, l1t::EtSum::EtSumType::kTotalEtEm, pt, 0, 0, 0);
510 
511  //ccla Generate uniform distribution of tower counts
512  int nTowers = 4095 * gRandom->Rndm();
513  l1t::EtSum towerCounts(*p4, l1t::EtSum::EtSumType::kTowerCount, nTowers, 0, 0, 0);
514 
515  //ccla Generate uniform distributions of AsymEt, AsymHt, AsymEtHF, AsymHtF
516  int nAsymEt = 255 * gRandom->Rndm();
517  l1t::EtSum AsymEt(*p4, l1t::EtSum::EtSumType::kAsymEt, nAsymEt, 0, 0, 0);
518  int nAsymHt = 255 * gRandom->Rndm();
519  l1t::EtSum AsymHt(*p4, l1t::EtSum::EtSumType::kAsymHt, nAsymHt, 0, 0, 0);
520  int nAsymEtHF = 255 * gRandom->Rndm();
521  l1t::EtSum AsymEtHF(*p4, l1t::EtSum::EtSumType::kAsymEtHF, nAsymEtHF, 0, 0, 0);
522  int nAsymHtHF = 255 * gRandom->Rndm();
523  l1t::EtSum AsymHtHF(*p4, l1t::EtSum::EtSumType::kAsymHtHF, nAsymHtHF, 0, 0, 0);
524 
525  pt = convertPtToHW(sumEt * 0.9, 2047, PtStep_);
526  l1t::EtSum htTotal(*p4, l1t::EtSum::EtSumType::kTotalHt, pt, 0, 0, 0);
527 
528  // Add EtSums for testing the MinBias Trigger (use some random numbers)
529  int hfP0val = gRandom->Poisson(4.);
530  if (hfP0val > 15)
531  hfP0val = 15;
532  l1t::EtSum hfP0(*p4, l1t::EtSum::EtSumType::kMinBiasHFP0, hfP0val, 0, 0, 0);
533 
534  int hfM0val = gRandom->Poisson(4.);
535  if (hfM0val > 15)
536  hfM0val = 15;
537  l1t::EtSum hfM0(*p4, l1t::EtSum::EtSumType::kMinBiasHFM0, hfM0val, 0, 0, 0);
538 
539  int hfP1val = gRandom->Poisson(4.);
540  if (hfP1val > 15)
541  hfP1val = 15;
542  l1t::EtSum hfP1(*p4, l1t::EtSum::EtSumType::kMinBiasHFP1, hfP1val, 0, 0, 0);
543 
544  int hfM1val = gRandom->Poisson(4.);
545  if (hfM1val > 15)
546  hfM1val = 15;
547  l1t::EtSum hfM1(*p4, l1t::EtSum::EtSumType::kMinBiasHFM1, hfM1val, 0, 0, 0);
548 
549  // Do same for Centrality
550  int cent30val(0), cent74val(0);
551  int centa = gRandom->Poisson(2.);
552  int centb = gRandom->Poisson(2.);
553  if (centa >= centb) {
554  cent30val = centa;
555  cent74val = centb;
556  } else {
557  cent30val = centb;
558  cent74val = centa;
559  }
560 
561  if (cent30val > 15)
562  cent30val = 15;
563  if (cent74val > 15)
564  cent74val = 15;
565 
566  int shift = 4;
567  int centralval = 0;
568  centralval |= cent30val & 0xF;
569  centralval |= (cent74val & 0xF) << shift;
570 
571  l1t::EtSum centrality(*p4, l1t::EtSum::EtSumType::kCentrality, centralval, 0, 0, 0);
572 
573  int mpt = 0;
574  int mphi = 0;
575  int mptHf = 0;
576  int mphiHf = 0;
577  int mhpt = 0;
578  int mhphi = 0;
579  int mhptHf = 0;
580  int mhphiHf = 0;
581 
583  // Make sure that you can get genMET
584  if (iEvent.getByToken(genMetToken, genMet)) {
585  mpt = convertPtToHW(genMet->front().pt(), MaxEt_, PtStep_);
586  mphi = convertPhiToHW(genMet->front().phi(), PhiStepCalo_);
587 
588  // Make Missing Et with HF slightly largeer and rotated (These are all fake inputs anyway...not supposed to be realistic)
589  mptHf = convertPtToHW(genMet->front().pt() * 1.1, MaxEt_, PtStep_);
590  mphiHf = convertPhiToHW(genMet->front().phi() + 3.14 / 7., PhiStepCalo_);
591 
592  // Make Missing Ht slightly smaller and rotated (These are all fake inputs anyway...not supposed to be realistic)
593  mhpt = convertPtToHW(genMet->front().pt() * 0.9, MaxEt_, PtStep_);
594  mhphi = convertPhiToHW(genMet->front().phi() + 3.14 / 5., PhiStepCalo_);
595 
596  // Ditto with Hissing Ht with HF
597  mhptHf = convertPtToHW(genMet->front().pt() * 0.95, MaxEt_, PtStep_);
598  mhphiHf = convertPhiToHW(genMet->front().phi() + 3.14 / 6., PhiStepCalo_);
599  } else {
600  LogTrace("GtGenToInputProducer") << ">>> GenMet collection not found!" << std::endl;
601  }
602 
603  // Missing Et and missing htt
604  l1t::EtSum etmiss(*p4, l1t::EtSum::EtSumType::kMissingEt, mpt, 0, mphi, 0);
605  l1t::EtSum etmissHF(*p4, l1t::EtSum::EtSumType::kMissingEtHF, mptHf, 0, mphiHf, 0);
606  l1t::EtSum htmiss(*p4, l1t::EtSum::EtSumType::kMissingHt, mhpt, 0, mhphi, 0);
607  l1t::EtSum htmissHF(*p4, l1t::EtSum::EtSumType::kMissingHtHF, mhptHf, 0, mhphiHf, 0);
608 
609  // Fill the EtSums in the Correct order
610  etsumVec.push_back(etTotal);
611  etsumVec.push_back(etEmTotal);
612  etsumVec.push_back(hfP0); // Frame0
613 
614  etsumVec.push_back(htTotal);
615  etsumVec.push_back(towerCounts);
616  etsumVec.push_back(hfM0); //Frame1
617 
618  etsumVec.push_back(etmiss);
619  etsumVec.push_back(AsymEt);
620  etsumVec.push_back(hfP1); //Frame2
621 
622  etsumVec.push_back(htmiss);
623  etsumVec.push_back(AsymHt);
624  etsumVec.push_back(hfM1); //Frame3
625 
626  etsumVec.push_back(etmissHF);
627  etsumVec.push_back(AsymEtHF); // Frame4
628 
629  etsumVec.push_back(htmissHF);
630  etsumVec.push_back(AsymHtHF);
631  etsumVec.push_back(centrality); // Frame5
632 
633  // Fill in some external conditions for testing
634  if ((iEvent.id().event()) % 2 == 0) {
635  for (int i = 0; i < 255; i = i + 2)
636  extCond_bx.setExternalDecision(i, true);
637  } else {
638  for (int i = 1; i < 255; i = i + 2)
639  extCond_bx.setExternalDecision(i, true);
640  }
641 
642  // Insert all the bx into the L1 Collections
643  //printf("Event %i EmptyBxEvt %i emptyBxTrailer %i diff %i \n",eventCnt_,emptyBxEvt_,emptyBxTrailer_,(emptyBxEvt_ - eventCnt_));
644 
645  // Fill Muons
646  for (int iMu = 0; iMu < int(muonVec_bxm2.size()); iMu++) {
647  muons->push_back(-2, muonVec_bxm2[iMu]);
648  }
649  for (int iMu = 0; iMu < int(muonVec_bxm1.size()); iMu++) {
650  muons->push_back(-1, muonVec_bxm1[iMu]);
651  }
652  for (int iMu = 0; iMu < int(muonVec_bx0.size()); iMu++) {
653  muons->push_back(0, muonVec_bx0[iMu]);
654  }
655  for (int iMu = 0; iMu < int(muonVec_bxp1.size()); iMu++) {
656  muons->push_back(1, muonVec_bxp1[iMu]);
657  }
659  for (int iMu = 0; iMu < int(muonVec.size()); iMu++) {
660  muons->push_back(2, muonVec[iMu]);
661  }
662  } else {
663  // this event is part of empty trailer...clear out data
664  muonVec.clear();
665  }
666 
667  // Fill MuonShowers
668  for (int iMuShower = 0; iMuShower < int(muonShowerVec_bxm2.size()); iMuShower++) {
669  muonShowers->push_back(-2, muonShowerVec_bxm2[iMuShower]);
670  }
671  for (int iMuShower = 0; iMuShower < int(muonShowerVec_bxm1.size()); iMuShower++) {
672  muonShowers->push_back(-1, muonShowerVec_bxm1[iMuShower]);
673  }
674  for (int iMuShower = 0; iMuShower < int(muonShowerVec_bx0.size()); iMuShower++) {
675  muonShowers->push_back(0, muonShowerVec_bx0[iMuShower]);
676  }
677  for (int iMuShower = 0; iMuShower < int(muonShowerVec_bxp1.size()); iMuShower++) {
678  muonShowers->push_back(1, muonShowerVec_bxp1[iMuShower]);
679  }
681  for (int iMuShower = 0; iMuShower < int(muonShowerVec.size()); iMuShower++) {
682  muonShowers->push_back(2, muonShowerVec[iMuShower]);
683  }
684  } else {
685  // this event is part of empty trailer...clear out data
686  muonShowerVec.clear();
687  }
688 
689  // Fill Egammas
690  for (int iEG = 0; iEG < int(egammaVec_bxm2.size()); iEG++) {
691  egammas->push_back(-2, egammaVec_bxm2[iEG]);
692  }
693  for (int iEG = 0; iEG < int(egammaVec_bxm1.size()); iEG++) {
694  egammas->push_back(-1, egammaVec_bxm1[iEG]);
695  }
696  for (int iEG = 0; iEG < int(egammaVec_bx0.size()); iEG++) {
697  egammas->push_back(0, egammaVec_bx0[iEG]);
698  }
699  for (int iEG = 0; iEG < int(egammaVec_bxp1.size()); iEG++) {
700  egammas->push_back(1, egammaVec_bxp1[iEG]);
701  }
703  for (int iEG = 0; iEG < int(egammaVec.size()); iEG++) {
704  egammas->push_back(2, egammaVec[iEG]);
705  }
706  } else {
707  // this event is part of empty trailer...clear out data
708  egammaVec.clear();
709  }
710 
711  // Fill Taus
712  for (int iTau = 0; iTau < int(tauVec_bxm2.size()); iTau++) {
713  taus->push_back(-2, tauVec_bxm2[iTau]);
714  }
715  for (int iTau = 0; iTau < int(tauVec_bxm1.size()); iTau++) {
716  taus->push_back(-1, tauVec_bxm1[iTau]);
717  }
718  for (int iTau = 0; iTau < int(tauVec_bx0.size()); iTau++) {
719  taus->push_back(0, tauVec_bx0[iTau]);
720  }
721  for (int iTau = 0; iTau < int(tauVec_bxp1.size()); iTau++) {
722  taus->push_back(1, tauVec_bxp1[iTau]);
723  }
725  for (int iTau = 0; iTau < int(tauVec.size()); iTau++) {
726  taus->push_back(2, tauVec[iTau]);
727  }
728  } else {
729  // this event is part of empty trailer...clear out data
730  tauVec.clear();
731  }
732 
733  // Fill Jets
734  for (int iJet = 0; iJet < int(jetVec_bxm2.size()); iJet++) {
735  jets->push_back(-2, jetVec_bxm2[iJet]);
736  }
737  for (int iJet = 0; iJet < int(jetVec_bxm1.size()); iJet++) {
738  jets->push_back(-1, jetVec_bxm1[iJet]);
739  }
740  for (int iJet = 0; iJet < int(jetVec_bx0.size()); iJet++) {
741  jets->push_back(0, jetVec_bx0[iJet]);
742  }
743  for (int iJet = 0; iJet < int(jetVec_bxp1.size()); iJet++) {
744  jets->push_back(1, jetVec_bxp1[iJet]);
745  }
747  for (int iJet = 0; iJet < int(jetVec.size()); iJet++) {
748  jets->push_back(2, jetVec[iJet]);
749  }
750  } else {
751  // this event is part of empty trailer...clear out data
752  jetVec.clear();
753  }
754 
755  // Fill Etsums
756  for (int iETsum = 0; iETsum < int(etsumVec_bxm2.size()); iETsum++) {
757  etsums->push_back(-2, etsumVec_bxm2[iETsum]);
758  }
759  for (int iETsum = 0; iETsum < int(etsumVec_bxm1.size()); iETsum++) {
760  etsums->push_back(-1, etsumVec_bxm1[iETsum]);
761  }
762  for (int iETsum = 0; iETsum < int(etsumVec_bx0.size()); iETsum++) {
763  etsums->push_back(0, etsumVec_bx0[iETsum]);
764  }
765  for (int iETsum = 0; iETsum < int(etsumVec_bxp1.size()); iETsum++) {
766  etsums->push_back(1, etsumVec_bxp1[iETsum]);
767  }
769  for (int iETsum = 0; iETsum < int(etsumVec.size()); iETsum++) {
770  etsums->push_back(2, etsumVec[iETsum]);
771  }
772  } else {
773  // this event is part of empty trailer...clear out data
774  etsumVec.clear();
775  }
776 
777  // Fill Externals
778  extCond->push_back(-2, extCond_bxm2);
779  extCond->push_back(-1, extCond_bxm1);
780  extCond->push_back(0, extCond_bx0);
781  extCond->push_back(1, extCond_bxp1);
783  extCond->push_back(2, extCond_bx);
784  } else {
785  // this event is part of the empty trailer...clear out data
786  extCond_bx.reset();
787  }
788 
789  iEvent.put(std::move(egammas));
790  iEvent.put(std::move(muons));
791  iEvent.put(std::move(muonShowers));
792  iEvent.put(std::move(taus));
793  iEvent.put(std::move(jets));
794  iEvent.put(std::move(etsums));
795  iEvent.put(std::move(extCond));
796 
797  // Now shift the bx data by one to prepare for next event.
805 
813 
821 
822  muonVec_bxp1 = muonVec;
823  muonShowerVec_bxp1 = muonShowerVec;
824  egammaVec_bxp1 = egammaVec;
825  tauVec_bxp1 = tauVec;
826  jetVec_bxp1 = jetVec;
827  etsumVec_bxp1 = etsumVec;
828  extCond_bxp1 = extCond_bx;
829  }
BXVector< GlobalExtBlk > GlobalExtBlkBxCollection
Definition: GlobalExtBlk.h:29
std::vector< l1t::EtSum > etsumVec_bxm1
void reset()
reset the content of a GlobalExtBlk
Definition: GlobalExtBlk.cc:58
std::vector< l1t::MuonShower > muonShowerVec_bxm1
std::vector< l1t::Muon > muonVec_bxm1
static constexpr int nTowers
std::vector< l1t::EtSum > etsumVec_bx0
virtual double pt() const =0
transverse momentum
std::vector< l1t::Muon > muonVec_bxp1
virtual int status() const =0
status word
Definition: Tau.h:20
std::vector< l1t::Tau > tauVec_bx0
std::vector< l1t::Tau > tauVec_bxp1
std::vector< l1t::Tau > tauVec_bxm2
std::vector< l1t::Muon > muonVec_bxm2
std::vector< l1t::EtSum > etsumVec_bxp1
muons
the two sets of parameters below are mutually exclusive, depending if RECO or ALCARECO is used the us...
Definition: DiMuonV_cfg.py:214
#define LogTrace(id)
std::vector< l1t::EGamma > egammaVec_bxm1
Definition: Jet.h:20
int convertPtToHW(double ipt, int maxPt, double step)
int iEvent
Definition: GenABIO.cc:224
std::vector< l1t::EGamma > egammaVec_bxm2
std::vector< l1t::EGamma > egammaVec_bxp1
edm::EDGetTokenT< reco::GenJetCollection > genJetsToken
std::vector< l1t::Muon > muonVec_bx0
std::vector< l1t::Tau > tauVec_bxm1
std::vector< l1t::Jet > jetVec_bxm1
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual int charge() const =0
electric charge
std::vector< l1t::Jet > jetVec_bxp1
std::vector< l1t::EGamma > egammaVec_bx0
void setExternalDecision(unsigned int bit, bool val)
Set decision bits.
Definition: GlobalExtBlk.cc:40
int convertPhiToHW(double iphi, int steps)
Definition: Muon.h:21
virtual int pdgId() const =0
PDG identifier.
#define debug
Definition: HDRShower.cc:19
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
int convertEtaToHW(double ieta, double minEta, double maxEta, int steps)
edm::EDGetTokenT< reco::GenMETCollection > genMetToken
std::vector< l1t::MuonShower > muonShowerVec_bxm2
static unsigned int const shift
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken
isoSum
===> compute the isolation and find the most isolated track
std::vector< l1t::Jet > jetVec_bxm2
std::vector< l1t::EtSum > etsumVec_bxm2
def move(src, dest)
Definition: eostools.py:511
std::vector< l1t::Jet > jetVec_bx0
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity
std::vector< l1t::MuonShower > muonShowerVec_bxp1
#define LogDebug(id)
std::vector< l1t::MuonShower > muonShowerVec_bx0

Member Data Documentation

◆ bxFirst_

int l1t::GenToInputProducer::bxFirst_
private

Definition at line 92 of file GenToInputProducer.cc.

◆ bxLast_

int l1t::GenToInputProducer::bxLast_
private

Definition at line 93 of file GenToInputProducer.cc.

◆ egammaVec_bx0

std::vector<l1t::EGamma> l1t::GenToInputProducer::egammaVec_bx0
private

Definition at line 128 of file GenToInputProducer.cc.

◆ egammaVec_bxm1

std::vector<l1t::EGamma> l1t::GenToInputProducer::egammaVec_bxm1
private

Definition at line 127 of file GenToInputProducer.cc.

◆ egammaVec_bxm2

std::vector<l1t::EGamma> l1t::GenToInputProducer::egammaVec_bxm2
private

Definition at line 126 of file GenToInputProducer.cc.

◆ egammaVec_bxp1

std::vector<l1t::EGamma> l1t::GenToInputProducer::egammaVec_bxp1
private

Definition at line 129 of file GenToInputProducer.cc.

◆ egEtThreshold_

double l1t::GenToInputProducer::egEtThreshold_
private

Definition at line 103 of file GenToInputProducer.cc.

◆ emptyBxEvt_

int l1t::GenToInputProducer::emptyBxEvt_
private

Definition at line 108 of file GenToInputProducer.cc.

◆ emptyBxTrailer_

int l1t::GenToInputProducer::emptyBxTrailer_
private

Definition at line 107 of file GenToInputProducer.cc.

◆ etsumVec_bx0

std::vector<l1t::EtSum> l1t::GenToInputProducer::etsumVec_bx0
private

Definition at line 143 of file GenToInputProducer.cc.

◆ etsumVec_bxm1

std::vector<l1t::EtSum> l1t::GenToInputProducer::etsumVec_bxm1
private

Definition at line 142 of file GenToInputProducer.cc.

◆ etsumVec_bxm2

std::vector<l1t::EtSum> l1t::GenToInputProducer::etsumVec_bxm2
private

Definition at line 141 of file GenToInputProducer.cc.

◆ etsumVec_bxp1

std::vector<l1t::EtSum> l1t::GenToInputProducer::etsumVec_bxp1
private

Definition at line 144 of file GenToInputProducer.cc.

◆ eventCnt_

int l1t::GenToInputProducer::eventCnt_
private

Definition at line 109 of file GenToInputProducer.cc.

◆ extCond_bx0

GlobalExtBlk l1t::GenToInputProducer::extCond_bx0
private

Definition at line 148 of file GenToInputProducer.cc.

◆ extCond_bxm1

GlobalExtBlk l1t::GenToInputProducer::extCond_bxm1
private

Definition at line 147 of file GenToInputProducer.cc.

◆ extCond_bxm2

GlobalExtBlk l1t::GenToInputProducer::extCond_bxm2
private

Definition at line 146 of file GenToInputProducer.cc.

◆ extCond_bxp1

GlobalExtBlk l1t::GenToInputProducer::extCond_bxp1
private

Definition at line 149 of file GenToInputProducer.cc.

◆ genJetsToken

edm::EDGetTokenT<reco::GenJetCollection> l1t::GenToInputProducer::genJetsToken
private

Definition at line 113 of file GenToInputProducer.cc.

◆ genMetToken

edm::EDGetTokenT<reco::GenMETCollection> l1t::GenToInputProducer::genMetToken
private

Definition at line 114 of file GenToInputProducer.cc.

◆ genParticlesToken

edm::EDGetTokenT<reco::GenParticleCollection> l1t::GenToInputProducer::genParticlesToken
private

Definition at line 112 of file GenToInputProducer.cc.

◆ gRandom

TRandom3* l1t::GenToInputProducer::gRandom
private

Definition at line 88 of file GenToInputProducer.cc.

◆ jetEtThreshold_

double l1t::GenToInputProducer::jetEtThreshold_
private

Definition at line 101 of file GenToInputProducer.cc.

◆ jetVec_bx0

std::vector<l1t::Jet> l1t::GenToInputProducer::jetVec_bx0
private

Definition at line 138 of file GenToInputProducer.cc.

◆ jetVec_bxm1

std::vector<l1t::Jet> l1t::GenToInputProducer::jetVec_bxm1
private

Definition at line 137 of file GenToInputProducer.cc.

◆ jetVec_bxm2

std::vector<l1t::Jet> l1t::GenToInputProducer::jetVec_bxm2
private

Definition at line 136 of file GenToInputProducer.cc.

◆ jetVec_bxp1

std::vector<l1t::Jet> l1t::GenToInputProducer::jetVec_bxp1
private

Definition at line 139 of file GenToInputProducer.cc.

◆ m_paramsCacheId

unsigned long long l1t::GenToInputProducer::m_paramsCacheId
private

Definition at line 83 of file GenToInputProducer.cc.

◆ maxNumEGCands_

int l1t::GenToInputProducer::maxNumEGCands_
private

Definition at line 98 of file GenToInputProducer.cc.

◆ maxNumJetCands_

int l1t::GenToInputProducer::maxNumJetCands_
private

Definition at line 97 of file GenToInputProducer.cc.

◆ maxNumMuCands_

int l1t::GenToInputProducer::maxNumMuCands_
private

Definition at line 95 of file GenToInputProducer.cc.

◆ maxNumMuShowerCands_

int l1t::GenToInputProducer::maxNumMuShowerCands_
private

Definition at line 96 of file GenToInputProducer.cc.

◆ maxNumTauCands_

int l1t::GenToInputProducer::maxNumTauCands_
private

Definition at line 99 of file GenToInputProducer.cc.

◆ muEtThreshold_

double l1t::GenToInputProducer::muEtThreshold_
private

Definition at line 104 of file GenToInputProducer.cc.

◆ muonShowerVec_bx0

std::vector<l1t::MuonShower> l1t::GenToInputProducer::muonShowerVec_bx0
private

Definition at line 123 of file GenToInputProducer.cc.

◆ muonShowerVec_bxm1

std::vector<l1t::MuonShower> l1t::GenToInputProducer::muonShowerVec_bxm1
private

Definition at line 122 of file GenToInputProducer.cc.

◆ muonShowerVec_bxm2

std::vector<l1t::MuonShower> l1t::GenToInputProducer::muonShowerVec_bxm2
private

Definition at line 121 of file GenToInputProducer.cc.

◆ muonShowerVec_bxp1

std::vector<l1t::MuonShower> l1t::GenToInputProducer::muonShowerVec_bxp1
private

Definition at line 124 of file GenToInputProducer.cc.

◆ muonVec_bx0

std::vector<l1t::Muon> l1t::GenToInputProducer::muonVec_bx0
private

Definition at line 118 of file GenToInputProducer.cc.

◆ muonVec_bxm1

std::vector<l1t::Muon> l1t::GenToInputProducer::muonVec_bxm1
private

Definition at line 117 of file GenToInputProducer.cc.

◆ muonVec_bxm2

std::vector<l1t::Muon> l1t::GenToInputProducer::muonVec_bxm2
private

Definition at line 116 of file GenToInputProducer.cc.

◆ muonVec_bxp1

std::vector<l1t::Muon> l1t::GenToInputProducer::muonVec_bxp1
private

Definition at line 119 of file GenToInputProducer.cc.

◆ mushowerRandom

TRandom3* l1t::GenToInputProducer::mushowerRandom
private

Definition at line 89 of file GenToInputProducer.cc.

◆ tauEtThreshold_

double l1t::GenToInputProducer::tauEtThreshold_
private

Definition at line 102 of file GenToInputProducer.cc.

◆ tauVec_bx0

std::vector<l1t::Tau> l1t::GenToInputProducer::tauVec_bx0
private

Definition at line 133 of file GenToInputProducer.cc.

◆ tauVec_bxm1

std::vector<l1t::Tau> l1t::GenToInputProducer::tauVec_bxm1
private

Definition at line 132 of file GenToInputProducer.cc.

◆ tauVec_bxm2

std::vector<l1t::Tau> l1t::GenToInputProducer::tauVec_bxm2
private

Definition at line 131 of file GenToInputProducer.cc.

◆ tauVec_bxp1

std::vector<l1t::Tau> l1t::GenToInputProducer::tauVec_bxp1
private

Definition at line 134 of file GenToInputProducer.cc.