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:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 GenToInputProducer (const ParameterSet &)
 
 ~GenToInputProducer () override
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducer () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 

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_
 
int counter_
 
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
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

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

Definition at line 60 of file GenToInputProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 147 of file GenToInputProducer.cc.

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

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

Definition at line 187 of file GenToInputProducer.cc.

188  {
189  }

Member Function Documentation

void l1t::GenToInputProducer::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 751 of file GenToInputProducer.cc.

752 {
753 }
void l1t::GenToInputProducer::beginRun ( Run const &  iR,
EventSetup const &  iE 
)
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 762 of file GenToInputProducer.cc.

References LogDebug.

762  {
763 
764  LogDebug("GtGenToInputProducer") << "GenToInputProducer::beginRun function called...\n";
765 
766  counter_ = 0;
767  srand( 0 );
768 
769  gRandom = new TRandom3();
770 }
#define LogDebug(id)
int l1t::GenToInputProducer::convertEtaToHW ( double  ieta,
double  minEta,
double  maxEta,
int  steps 
)
private

Definition at line 789 of file GenToInputProducer.cc.

References createfilelist::int, and trackingPOGFilters_cfi::minEta.

789  {
790 
791  double binWidth = (maxEta - minEta)/steps;
792 
793  //if we are outside the limits, set error
794  if(ieta < minEta) return 99999;//ieta = minEta+binWidth/2.;
795  if(ieta > maxEta) return 99999;//ieta = maxEta-binWidth/2.;
796 
797  int binNum = (int)(ieta/binWidth);
798  if(ieta<0.) binNum--;
799 
800 // unsigned int hwEta = binNum & bitMask;
801 // Remove masking for BXVectors...only assume in raw data
802 
803  return binNum;
804 }
double maxEta
int l1t::GenToInputProducer::convertPhiToHW ( double  iphi,
int  steps 
)
private

Definition at line 779 of file GenToInputProducer.cc.

References createfilelist::int, M_PI, and AlignmentTrackSelector_cfi::phiMax.

779  {
780 
781  double phiMax = 2 * M_PI;
782  if( iphi < 0 ) iphi += 2*M_PI;
783  if( iphi > phiMax) iphi -= phiMax;
784 
785  int hwPhi = int( (iphi/phiMax)*steps + 0.00001 );
786  return hwPhi;
787 }
#define M_PI
int l1t::GenToInputProducer::convertPtToHW ( double  ipt,
int  maxPt,
double  step 
)
private

Definition at line 806 of file GenToInputProducer.cc.

References createfilelist::int, and MuonErrorMatrixAnalyzer_cfi::maxPt.

806  {
807 
808  int hwPt = int( ipt/step + 0.0001 );
809  // if above max Pt, set to largest value
810  if( hwPt > maxPt ) hwPt = maxPt;
811 
812  return hwPt;
813 }
step
Definition: StallMonitor.cc:94
void l1t::GenToInputProducer::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 757 of file GenToInputProducer.cc.

757  {
758 }
void l1t::GenToInputProducer::endRun ( Run const &  iR,
EventSetup const &  iE 
)
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 773 of file GenToInputProducer.cc.

773  {
774 
775 }
void l1t::GenToInputProducer::fillDescriptions ( ConfigurationDescriptions descriptions)
static

Definition at line 818 of file GenToInputProducer.cc.

References edm::ConfigurationDescriptions::addDefault(), DEFINE_FWK_MODULE, and edm::ParameterSetDescription::setUnknown().

818  {
819  //The following says we do not know what parameters are allowed so do no validation
820  // Please change this to state exactly what you do use, even if it is no parameters
822  desc.setUnknown();
823  descriptions.addDefault(desc);
824 }
void addDefault(ParameterSetDescription const &psetDescription)
void l1t::GenToInputProducer::produce ( Event ,
EventSetup const &   
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 199 of file GenToInputProducer.cc.

References funct::abs(), ALCARECOTkAlJpsiMuMu_cff::charge, reco::Candidate::charge(), debug, particleFlow_cfi::dEta, particleFlow_cfi::dPhi, PVValHelper::eta, reco::Candidate::eta(), edm::EventID::event(), ttbarCategorization_cff::genJets, HepMCValidationHelper::genMet(), GenHFHadronMatcher_cfi::genParticles, edm::Event::getByToken(), mps_fire::i, edm::EventBase::id(), createfilelist::int, objects.IsoTrackAnalyzer::isoSum, metsig::jet, fwrapper::jets, gen::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(), RPCpg::mu, extraflags_cff::muons, ecaldqm::nTowers, p4, common_cff::pdgId, reco::Candidate::pdgId(), reco::Candidate::phi(), EnergyCorrector::pt, reco::Candidate::pt(), edm::Event::put(), GlobalExtBlk::reset(), GlobalExtBlk::setExternalDecision(), edm::shift, mps_update::status, reco::Candidate::status(), objects.METAnalyzer::sumEt, GlobalPosition_Frontier_DevDB_cff::tag, metsig::tau, and nano_cff::taus.

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

Member Data Documentation

int l1t::GenToInputProducer::bxFirst_
private

Definition at line 87 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::bxLast_
private

Definition at line 88 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::counter_
private

Definition at line 110 of file GenToInputProducer.cc.

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

Definition at line 119 of file GenToInputProducer.cc.

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

Definition at line 118 of file GenToInputProducer.cc.

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

Definition at line 117 of file GenToInputProducer.cc.

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

Definition at line 120 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::egEtThreshold_
private

Definition at line 97 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::emptyBxEvt_
private

Definition at line 102 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::emptyBxTrailer_
private

Definition at line 101 of file GenToInputProducer.cc.

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

Definition at line 134 of file GenToInputProducer.cc.

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

Definition at line 133 of file GenToInputProducer.cc.

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

Definition at line 132 of file GenToInputProducer.cc.

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

Definition at line 135 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::eventCnt_
private

Definition at line 103 of file GenToInputProducer.cc.

GlobalExtBlk l1t::GenToInputProducer::extCond_bx0
private

Definition at line 139 of file GenToInputProducer.cc.

GlobalExtBlk l1t::GenToInputProducer::extCond_bxm1
private

Definition at line 138 of file GenToInputProducer.cc.

GlobalExtBlk l1t::GenToInputProducer::extCond_bxm2
private

Definition at line 137 of file GenToInputProducer.cc.

GlobalExtBlk l1t::GenToInputProducer::extCond_bxp1
private

Definition at line 140 of file GenToInputProducer.cc.

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

Definition at line 107 of file GenToInputProducer.cc.

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

Definition at line 108 of file GenToInputProducer.cc.

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

Definition at line 106 of file GenToInputProducer.cc.

TRandom3* l1t::GenToInputProducer::gRandom
private

Definition at line 84 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::jetEtThreshold_
private

Definition at line 95 of file GenToInputProducer.cc.

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

Definition at line 129 of file GenToInputProducer.cc.

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

Definition at line 128 of file GenToInputProducer.cc.

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

Definition at line 127 of file GenToInputProducer.cc.

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

Definition at line 130 of file GenToInputProducer.cc.

unsigned long long l1t::GenToInputProducer::m_paramsCacheId
private

Definition at line 79 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumEGCands_
private

Definition at line 92 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumJetCands_
private

Definition at line 91 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumMuCands_
private

Definition at line 90 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumTauCands_
private

Definition at line 93 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::muEtThreshold_
private

Definition at line 98 of file GenToInputProducer.cc.

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

Definition at line 114 of file GenToInputProducer.cc.

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

Definition at line 113 of file GenToInputProducer.cc.

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

Definition at line 112 of file GenToInputProducer.cc.

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

Definition at line 115 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::tauEtThreshold_
private

Definition at line 96 of file GenToInputProducer.cc.

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

Definition at line 124 of file GenToInputProducer.cc.

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

Definition at line 123 of file GenToInputProducer.cc.

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

Definition at line 122 of file GenToInputProducer.cc.

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

Definition at line 125 of file GenToInputProducer.cc.