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 maxNumTauCands_
 
double muEtThreshold_
 
std::vector< l1t::MuonmuonVec_bx0
 
std::vector< l1t::MuonmuonVec_bxm1
 
std::vector< l1t::MuonmuonVec_bxm2
 
std::vector< l1t::MuonmuonVec_bxp1
 
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 63 of file GenToInputProducer.cc.

Constructor & Destructor Documentation

◆ GenToInputProducer()

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

Definition at line 147 of file GenToInputProducer.cc.

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

147  {
148  // register what you produce
149  produces<BXVector<l1t::EGamma>>();
150  produces<BXVector<l1t::Muon>>();
151  produces<BXVector<l1t::Tau>>();
152  produces<BXVector<l1t::Jet>>();
153  produces<BXVector<l1t::EtSum>>();
154  produces<GlobalExtBlkBxCollection>();
155 
156  // Setup parameters
157  bxFirst_ = iConfig.getParameter<int>("bxFirst");
158  bxLast_ = iConfig.getParameter<int>("bxLast");
159 
160  maxNumMuCands_ = iConfig.getParameter<int>("maxMuCand");
161  maxNumJetCands_ = iConfig.getParameter<int>("maxJetCand");
162  maxNumEGCands_ = iConfig.getParameter<int>("maxEGCand");
163  maxNumTauCands_ = iConfig.getParameter<int>("maxTauCand");
164 
165  jetEtThreshold_ = iConfig.getParameter<double>("jetEtThreshold");
166  tauEtThreshold_ = iConfig.getParameter<double>("tauEtThreshold");
167  egEtThreshold_ = iConfig.getParameter<double>("egEtThreshold");
168  muEtThreshold_ = iConfig.getParameter<double>("muEtThreshold");
169 
170  emptyBxTrailer_ = iConfig.getParameter<int>("emptyBxTrailer");
171  emptyBxEvt_ = iConfig.getParameter<int>("emptyBxEvt");
172 
173  genParticlesToken = consumes<reco::GenParticleCollection>(std::string("genParticles"));
174  genJetsToken = consumes<reco::GenJetCollection>(std::string("ak4GenJets"));
175  genMetToken = consumes<reco::GenMETCollection>(std::string("genMetCalo"));
176 
177  // set cache id to zero, will be set at first beginRun:
178  m_paramsCacheId = 0;
179  eventCnt_ = 0;
180  }
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 182 of file GenToInputProducer.cc.

182 {}

Member Function Documentation

◆ beginJob()

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

Definition at line 776 of file GenToInputProducer.cc.

776 {}

◆ beginRun()

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

Definition at line 783 of file GenToInputProducer.cc.

References LogDebug.

783  {
784  LogDebug("GtGenToInputProducer") << "GenToInputProducer::beginRun function called...\n";
785 
786  srand(0);
787 
788  gRandom = new TRandom3();
789  }
#define LogDebug(id)

◆ convertEtaToHW()

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

Definition at line 806 of file GenToInputProducer.cc.

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

806  {
807  double binWidth = (maxEta - minEta) / steps;
808 
809  //if we are outside the limits, set error
810  if (ieta < minEta)
811  return 99999; //ieta = minEta+binWidth/2.;
812  if (ieta > maxEta)
813  return 99999; //ieta = maxEta-binWidth/2.;
814 
815  int binNum = (int)(ieta / binWidth);
816  if (ieta < 0.)
817  binNum--;
818 
819  // unsigned int hwEta = binNum & bitMask;
820  // Remove masking for BXVectors...only assume in raw data
821 
822  return binNum;
823  }

◆ convertPhiToHW()

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

Definition at line 795 of file GenToInputProducer.cc.

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

795  {
796  double phiMax = 2 * M_PI;
797  if (iphi < 0)
798  iphi += 2 * M_PI;
799  if (iphi > phiMax)
800  iphi -= phiMax;
801 
802  int hwPhi = int((iphi / phiMax) * steps + 0.00001);
803  return hwPhi;
804  }
#define M_PI

◆ convertPtToHW()

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

Definition at line 825 of file GenToInputProducer.cc.

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

825  {
826  int hwPt = int(ipt / step + 0.0001);
827  // if above max Pt, set to largest value
828  if (hwPt > maxPt)
829  hwPt = maxPt;
830 
831  return hwPt;
832  }
maxPt
Definition: PV_cfg.py:224
step
Definition: StallMonitor.cc:83

◆ endJob()

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

Definition at line 779 of file GenToInputProducer.cc.

779 {}

◆ endRun()

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

Definition at line 792 of file GenToInputProducer.cc.

792 {}

◆ fillDescriptions()

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

Definition at line 835 of file GenToInputProducer.cc.

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

835  {
836  //The following says we do not know what parameters are allowed so do no validation
837  // Please change this to state exactly what you do use, even if it is no parameters
839  desc.setUnknown();
840  descriptions.addDefault(desc);
841  }
void addDefault(ParameterSetDescription const &psetDescription)

◆ produce()

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

Definition at line 189 of file GenToInputProducer.cc.

References funct::abs(), simCaloStage2Layer1Digis_cfi::bxFirst, simCaloStage2Layer1Digis_cfi::bxLast, ALCARECOTkAlJpsiMuMu_cff::charge, reco::Candidate::charge(), debug, HLT_2024v10_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, objects.IsoTrackAnalyzer::isoSum, metsig::jet, PDWG_EXODelayedJetMET_cff::jets, dqmdumpme::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, l1tGTMenu_lepSeeds_cff::qual, GlobalExtBlk::reset(), GlobalExtBlk::setExternalDecision(), edm::shift, mps_update::status, reco::Candidate::status(), objects.METAnalyzer::sumEt, makeGlobalPositionRcd_cfg::tag, metsig::tau, and Tau3MuMonitor_cff::taus.

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

Member Data Documentation

◆ bxFirst_

int l1t::GenToInputProducer::bxFirst_
private

Definition at line 90 of file GenToInputProducer.cc.

◆ bxLast_

int l1t::GenToInputProducer::bxLast_
private

Definition at line 91 of file GenToInputProducer.cc.

◆ egammaVec_bx0

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

Definition at line 120 of file GenToInputProducer.cc.

◆ egammaVec_bxm1

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

Definition at line 119 of file GenToInputProducer.cc.

◆ egammaVec_bxm2

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

Definition at line 118 of file GenToInputProducer.cc.

◆ egammaVec_bxp1

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

Definition at line 121 of file GenToInputProducer.cc.

◆ egEtThreshold_

double l1t::GenToInputProducer::egEtThreshold_
private

Definition at line 100 of file GenToInputProducer.cc.

◆ emptyBxEvt_

int l1t::GenToInputProducer::emptyBxEvt_
private

Definition at line 105 of file GenToInputProducer.cc.

◆ emptyBxTrailer_

int l1t::GenToInputProducer::emptyBxTrailer_
private

Definition at line 104 of file GenToInputProducer.cc.

◆ etsumVec_bx0

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

Definition at line 135 of file GenToInputProducer.cc.

◆ etsumVec_bxm1

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

Definition at line 134 of file GenToInputProducer.cc.

◆ etsumVec_bxm2

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

Definition at line 133 of file GenToInputProducer.cc.

◆ etsumVec_bxp1

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

Definition at line 136 of file GenToInputProducer.cc.

◆ eventCnt_

int l1t::GenToInputProducer::eventCnt_
private

Definition at line 106 of file GenToInputProducer.cc.

◆ extCond_bx0

GlobalExtBlk l1t::GenToInputProducer::extCond_bx0
private

Definition at line 140 of file GenToInputProducer.cc.

◆ extCond_bxm1

GlobalExtBlk l1t::GenToInputProducer::extCond_bxm1
private

Definition at line 139 of file GenToInputProducer.cc.

◆ extCond_bxm2

GlobalExtBlk l1t::GenToInputProducer::extCond_bxm2
private

Definition at line 138 of file GenToInputProducer.cc.

◆ extCond_bxp1

GlobalExtBlk l1t::GenToInputProducer::extCond_bxp1
private

Definition at line 141 of file GenToInputProducer.cc.

◆ genJetsToken

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

Definition at line 110 of file GenToInputProducer.cc.

◆ genMetToken

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

Definition at line 111 of file GenToInputProducer.cc.

◆ genParticlesToken

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

Definition at line 109 of file GenToInputProducer.cc.

◆ gRandom

TRandom3* l1t::GenToInputProducer::gRandom
private

Definition at line 87 of file GenToInputProducer.cc.

◆ jetEtThreshold_

double l1t::GenToInputProducer::jetEtThreshold_
private

Definition at line 98 of file GenToInputProducer.cc.

◆ jetVec_bx0

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

Definition at line 130 of file GenToInputProducer.cc.

◆ jetVec_bxm1

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

Definition at line 129 of file GenToInputProducer.cc.

◆ jetVec_bxm2

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

Definition at line 128 of file GenToInputProducer.cc.

◆ jetVec_bxp1

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

Definition at line 131 of file GenToInputProducer.cc.

◆ m_paramsCacheId

unsigned long long l1t::GenToInputProducer::m_paramsCacheId
private

Definition at line 82 of file GenToInputProducer.cc.

◆ maxNumEGCands_

int l1t::GenToInputProducer::maxNumEGCands_
private

Definition at line 95 of file GenToInputProducer.cc.

◆ maxNumJetCands_

int l1t::GenToInputProducer::maxNumJetCands_
private

Definition at line 94 of file GenToInputProducer.cc.

◆ maxNumMuCands_

int l1t::GenToInputProducer::maxNumMuCands_
private

Definition at line 93 of file GenToInputProducer.cc.

◆ maxNumTauCands_

int l1t::GenToInputProducer::maxNumTauCands_
private

Definition at line 96 of file GenToInputProducer.cc.

◆ muEtThreshold_

double l1t::GenToInputProducer::muEtThreshold_
private

Definition at line 101 of file GenToInputProducer.cc.

◆ muonVec_bx0

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

Definition at line 115 of file GenToInputProducer.cc.

◆ muonVec_bxm1

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

Definition at line 114 of file GenToInputProducer.cc.

◆ muonVec_bxm2

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

Definition at line 113 of file GenToInputProducer.cc.

◆ muonVec_bxp1

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

Definition at line 116 of file GenToInputProducer.cc.

◆ tauEtThreshold_

double l1t::GenToInputProducer::tauEtThreshold_
private

Definition at line 99 of file GenToInputProducer.cc.

◆ tauVec_bx0

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

Definition at line 125 of file GenToInputProducer.cc.

◆ tauVec_bxm1

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

Definition at line 124 of file GenToInputProducer.cc.

◆ tauVec_bxm2

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

Definition at line 123 of file GenToInputProducer.cc.

◆ tauVec_bxp1

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

Definition at line 126 of file GenToInputProducer.cc.