CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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)
 

Private Member Functions

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

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_
 
edm::EDGetTokenT
< reco::GenJetCollection
genJetsToken
 
edm::EDGetTokenT
< reco::GenMETCollection
genMetToken
 
edm::EDGetTokenT
< reco::GenParticleCollection
genParticlesToken
 
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
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- 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 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 57 of file GenToInputProducer.cc.

Constructor & Destructor Documentation

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

Definition at line 137 of file GenToInputProducer.cc.

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

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

Definition at line 176 of file GenToInputProducer.cc.

177  {
178  }

Member Function Documentation

void l1t::GenToInputProducer::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 619 of file GenToInputProducer.cc.

620 {
621 }
void l1t::GenToInputProducer::beginRun ( Run const &  iR,
EventSetup const &  iE 
)
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 630 of file GenToInputProducer.cc.

References LogDebug.

630  {
631 
632  LogDebug("l1t|Global") << "GenToInputProducer::beginRun function called...\n";
633 
634  counter_ = 0;
635 }
#define LogDebug(id)
unsigned int l1t::GenToInputProducer::convertEtaToHW ( double  ieta,
double  minEta,
double  maxEta,
int  steps,
unsigned int  bitMask 
)
private

Definition at line 654 of file GenToInputProducer.cc.

References benchmark_cfg::minEta.

654  {
655 
656  double binWidth = (maxEta - minEta)/steps;
657 
658  //if we are outside the limits, set error
659  if(ieta < minEta) return 99999;//ieta = minEta+binWidth/2.;
660  if(ieta > maxEta) return 99999;//ieta = maxEta-binWidth/2.;
661 
662  int binNum = (int)(ieta/binWidth);
663  if(ieta<0.) binNum--;
664 
665  unsigned int hwEta = binNum & bitMask;
666 
667 
668  return hwEta;
669 }
double maxEta
int l1t::GenToInputProducer::convertPhiToHW ( double  iphi,
int  steps 
)
private

Definition at line 644 of file GenToInputProducer.cc.

References M_PI.

644  {
645 
646  double phiMax = 2 * M_PI;
647  if( iphi < 0 ) iphi += 2*M_PI;
648  if( iphi > phiMax) iphi -= phiMax;
649 
650  int hwPhi = int( (iphi/phiMax)*steps + 0.00001 );
651  return hwPhi;
652 }
#define M_PI
int l1t::GenToInputProducer::convertPtToHW ( double  ipt,
int  maxPt,
double  step 
)
private

Definition at line 671 of file GenToInputProducer.cc.

671  {
672 
673  int hwPt = int( ipt/step + 0.0001 );
674  // if above max Pt, set to largest value
675  if( hwPt > maxPt ) hwPt = maxPt;
676 
677  return hwPt;
678 }
void l1t::GenToInputProducer::endJob ( void  )
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 625 of file GenToInputProducer.cc.

625  {
626 }
void l1t::GenToInputProducer::endRun ( Run const &  iR,
EventSetup const &  iE 
)
privatevirtual

Reimplemented from edm::EDProducer.

Definition at line 638 of file GenToInputProducer.cc.

638  {
639 
640 }
void l1t::GenToInputProducer::fillDescriptions ( ConfigurationDescriptions descriptions)
static

Definition at line 683 of file GenToInputProducer.cc.

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

683  {
684  //The following says we do not know what parameters are allowed so do no validation
685  // Please change this to state exactly what you do use, even if it is no parameters
687  desc.setUnknown();
688  descriptions.addDefault(desc);
689 }
void addDefault(ParameterSetDescription const &psetDescription)
void l1t::GenToInputProducer::produce ( Event iEvent,
EventSetup const &  iSetup 
)
privatevirtual

Implements edm::EDProducer.

Definition at line 188 of file GenToInputProducer.cc.

References funct::abs(), runGlobalFakeInputProducer::bxFirst, runGlobalFakeInputProducer::bxLast, DeDxDiscriminatorTools::charge(), reco::Candidate::charge(), reco::Candidate::eta(), eta(), EtaStepCalo_, EtaStepMuon_, HepMCValidationHelper::genMet(), genParticleCandidates2GenParticles_cfi::genParticles, edm::Event::getByToken(), metsig::jet, fwrapper::jets, gen::k, LogDebug, LogTrace, MaxCaloEta_, MaxEt_, MaxJetPt_, MaxLepPt_, MaxMuonEta_, RPCpg::mu, patZpeak::muons, p4, benchmark_cfg::pdgId, reco::Candidate::pdgId(), phi, reco::Candidate::phi(), PhiStepCalo_, PhiStepMuon_, RecoTauCleanerPlugins::pt, reco::Candidate::pt(), PtStep_, edm::Event::put(), reco::Candidate::status(), ntuplemaker::status, GlobalPosition_Frontier_DevDB_cff::tag, and metsig::tau.

189 {
190 
191  eventCnt_++;
192 
193  LogDebug("l1t|Global") << "GenToInputProducer::produce function called...\n";
194 
195  // Setup vectors
196  std::vector<l1t::Muon> muonVec;
197  std::vector<l1t::EGamma> egammaVec;
198  std::vector<l1t::Tau> tauVec;
199  std::vector<l1t::Jet> jetVec;
200  std::vector<l1t::EtSum> etsumVec;
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 
207  // Default values objects
208  double MaxLepPt_ = 255;
209  double MaxJetPt_ = 1023;
210  double MaxEt_ = 2047;
211 
212  double MaxCaloEta_ = 5.0;
213  double MaxMuonEta_ = 2.45;
214 
215  double PhiStepCalo_ = 144;
216  double PhiStepMuon_ = 576;
217 
218  // eta scale
219  double EtaStepCalo_ = 230;
220  double EtaStepMuon_ = 450;
221 
222  // Et scale (in GeV)
223  double PtStep_ = 0.5;
224 
225 
226  //outputs
227  std::auto_ptr<l1t::EGammaBxCollection> egammas (new l1t::EGammaBxCollection(0, bxFirst, bxLast));
228  std::auto_ptr<l1t::MuonBxCollection> muons (new l1t::MuonBxCollection(0, bxFirst, bxLast));
229  std::auto_ptr<l1t::TauBxCollection> taus (new l1t::TauBxCollection(0, bxFirst, bxLast));
230  std::auto_ptr<l1t::JetBxCollection> jets (new l1t::JetBxCollection(0, bxFirst, bxLast));
231  std::auto_ptr<l1t::EtSumBxCollection> etsums (new l1t::EtSumBxCollection(0, bxFirst, bxLast));
232 
233  std::vector<int> mu_cands_index;
234  std::vector<int> eg_cands_index;
235  std::vector<int> tau_cands_index;
237  // Make sure that you can get genParticles
238  if( iEvent.getByToken(genParticlesToken, genParticles) ){
239 
240  for( size_t k = 0; k < genParticles->size(); k++ ){
241  const reco::Candidate & mcParticle = (*genParticles)[k];
242 
243  int status = mcParticle.status();
244  int pdgId = mcParticle.pdgId();
245  double pt = mcParticle.pt();
246 
247  // Only use status 1 particles (Tau's need to be allowed through..take status 2 taus)
248  if( status!=1 && !(abs(pdgId)==15 && status==2) ) continue;
249 
250  int absId = abs(pdgId);
251 
252  if( absId==11 && pt>=egEtThreshold_ ) eg_cands_index.push_back(k);
253  else if( absId==13 && pt>=muEtThreshold_ ) mu_cands_index.push_back(k);
254  else if( absId==15 && pt>=tauEtThreshold_ ) tau_cands_index.push_back(k);
255  }
256  }
257  else {
258  LogTrace("l1t|Global") << ">>> GenParticles collection not found!" << std::endl;
259  }
260 
261 
262 
263  // Muon Collection
264  int numMuCands = int( mu_cands_index.size() );
265  Int_t idxMu[numMuCands];
266  double muPtSorted[numMuCands];
267  for( int iMu=0; iMu<numMuCands; iMu++ ) 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 
272  if( iMu>=maxNumMuCands_ ) 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_, 0x1ff );
278  int phi = convertPhiToHW( mcParticle.phi(), PhiStepMuon_ );
279  int qual = 4;
280  int iso = 1;
281  int charge = ( mcParticle.charge()>0 ) ? 1 : 0;
282  int chargeValid = 1;
283  int mip = 1;
284  int tag = 1;
285 
286  // Eta outside of acceptance
287  if( eta>=9999 ) continue;
288 
289  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > *p4 = new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> >();
290 
291  l1t::Muon mu(*p4, pt, eta, phi, qual, charge, chargeValid, iso, mip, tag);
292  muonVec.push_back(mu);
293  }
294 
295 
296  // EG Collection
297  int numEgCands = int( eg_cands_index.size() );
298  Int_t idxEg[numEgCands];
299  double egPtSorted[numEgCands];
300  for( int iEg=0; iEg<numEgCands; iEg++ ) egPtSorted[iEg] = genParticles->at(eg_cands_index[iEg]).pt();
301 
302  TMath::Sort(numEgCands,egPtSorted,idxEg);
303  for( int iEg=0; iEg<numEgCands; iEg++ ){
304 
305  if( iEg>=maxNumEGCands_ ) continue;
306 
307  const reco::Candidate & mcParticle = (*genParticles)[eg_cands_index[idxEg[iEg]]];
308 
309  int pt = convertPtToHW( mcParticle.pt(), MaxLepPt_, PtStep_ );
310  int eta = convertEtaToHW( mcParticle.eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_ , 0xff);
311  int phi = convertPhiToHW( mcParticle.phi(), PhiStepCalo_ );
312  int qual = 1;
313  int iso = 1;
314 
315  // Eta outside of acceptance
316  if( eta>=9999 ) continue;
317 
318  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > *p4 = new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> >();
319 
320  l1t::EGamma eg(*p4, pt, eta, phi, qual, iso);
321  egammaVec.push_back(eg);
322  }
323 
324 
325 
326  // Tau Collection
327  int numTauCands = int( tau_cands_index.size() );
328  Int_t idxTau[numTauCands];
329  double tauPtSorted[numTauCands];
330  for( int iTau=0; iTau<numTauCands; iTau++ ) tauPtSorted[iTau] = genParticles->at(tau_cands_index[iTau]).pt();
331 
332  TMath::Sort(numTauCands,tauPtSorted,idxTau);
333  for( int iTau=0; iTau<numTauCands; iTau++ ){
334 
335  if( iTau>=maxNumTauCands_ ) continue;
336 
337  const reco::Candidate & mcParticle = (*genParticles)[tau_cands_index[idxTau[iTau]]];
338 
339  int pt = convertPtToHW( mcParticle.pt(), MaxLepPt_, PtStep_ );
340  int eta = convertEtaToHW( mcParticle.eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_ , 0xff);
341  int phi = convertPhiToHW( mcParticle.phi(), PhiStepCalo_ );
342  int qual = 1;
343  int iso = 1;
344 
345  // Eta outside of acceptance
346  if( eta>=9999 ) continue;
347 
348  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > *p4 = new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> >();
349 
350  l1t::Tau tau(*p4, pt, eta, phi, qual, iso);
351  tauVec.push_back(tau);
352  }
353 
354 
355  // Temporary hack to increase number of EGs and taus
356  int maxOtherEGs = 4;
357  int maxOtherTaus = 8;
358  int numCurrentEGs = int( egammaVec.size() );
359  int numCurrentTaus = int( tauVec.size() );
360 
361  int numExtraEGs=0, numExtraTaus=0;
362  // end hack
363 
364  // Use to sum the energy of the objects in the event for ETT and HTT
365  // sum all jets
366  double sumEt = 0;
367 
368  int nJet = 0;
370  // Make sure that you can get genJets
371  if( iEvent.getByToken(genJetsToken, genJets) ){ // Jet Collection
372  for(reco::GenJetCollection::const_iterator genJet = genJets->begin(); genJet!=genJets->end(); ++genJet ){
373 
374  //Keep running sum of total Et
375  sumEt += genJet->et();
376 
377  // Apply pt and eta cut?
378  if( genJet->pt()<jetEtThreshold_ ) continue;
379 
380  //
381  if( nJet>=maxNumJetCands_ ) continue;
382  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > *p4 = new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> >();
383 
384  int pt = convertPtToHW( genJet->et(), MaxJetPt_, PtStep_ );
385  int eta = convertEtaToHW( genJet->eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_ , 0xff);
386  int phi = convertPhiToHW( genJet->phi(), PhiStepCalo_ );
387 
388  // Eta outside of acceptance
389  if( eta>=9999 ) continue;
390 
391  int qual = 0;
392 
393  l1t::Jet jet(*p4, pt, eta, phi, qual);
394  jetVec.push_back(jet);
395 
396  nJet++;
397 
398  // Temporary hack to increase number of EGs and taus
399  if( (numExtraEGs+numCurrentEGs)<maxNumEGCands_ && numExtraEGs<maxOtherEGs ){
400  numExtraEGs++;
401 
402  int EGpt = convertPtToHW( genJet->et(), MaxLepPt_, PtStep_ );
403  int EGeta = convertEtaToHW( genJet->eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_ , 0xff);
404  int EGphi = convertPhiToHW( genJet->phi(), PhiStepCalo_ );
405 
406  int EGqual = 1;
407  int EGiso = 1;
408 
409  l1t::EGamma eg(*p4, EGpt, EGeta, EGphi, EGqual, EGiso);
410  egammaVec.push_back(eg);
411  }
412 
413  if( (numExtraTaus+numCurrentTaus)<maxNumTauCands_ && numExtraTaus<maxOtherTaus ){
414  numExtraTaus++;
415 
416  int Taupt = convertPtToHW( genJet->et(), MaxLepPt_, PtStep_ );
417  int Taueta = convertEtaToHW( genJet->eta(), -MaxCaloEta_, MaxCaloEta_, EtaStepCalo_ , 0xff);
418  int Tauphi = convertPhiToHW( genJet->phi(), PhiStepCalo_ );
419  int Tauqual = 1;
420  int Tauiso = 1;
421 
422  l1t::Tau tau(*p4, Taupt, Taueta, Tauphi, Tauqual, Tauiso);
423  tauVec.push_back(tau);
424  }
425  // end hack
426  }
427  }
428  else {
429  LogTrace("l1t|Global") << ">>> GenJets collection not found!" << std::endl;
430  }
431 
432 
434  // Make sure that you can get genMET
435  if( iEvent.getByToken(genMetToken, genMet) ){
436  int pt = convertPtToHW( genMet->front().pt(), MaxEt_, PtStep_ );
437  int phi = convertPhiToHW( genMet->front().phi(), PhiStepCalo_ );
438 
439  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > *p4 = new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> >();
440 
441  // Missing Et
442  l1t::EtSum etmiss(*p4, l1t::EtSum::EtSumType::kMissingEt,pt, 0,phi, 0);
443  etsumVec.push_back(etmiss);
444 
445  // Make Missing Ht slightly smaller and rotated (These are all fake inputs anyway...not supposed to be realistic)
446  pt = convertPtToHW( genMet->front().pt()*0.9, MaxEt_, PtStep_ );
447  phi = convertPhiToHW( genMet->front().phi()+ 3.14/5., PhiStepCalo_ );
448 
449  l1t::EtSum htmiss(*p4, l1t::EtSum::EtSumType::kMissingHt,pt, 0,phi, 0);
450  etsumVec.push_back(htmiss);
451 
452 
453  }
454  else {
455  LogTrace("l1t|Global") << ">>> GenMet collection not found!" << std::endl;
456  }
457 
458 
459 // Put the total Et into EtSums (Make HTT slightly smaller to tell them apart....not supposed to be realistic)
460  int pt = convertPtToHW( sumEt, 2047, PtStep_ );
461  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> > *p4 = new ROOT::Math::LorentzVector<ROOT::Math::PxPyPzE4D<double> >();
462  l1t::EtSum etTotal(*p4, l1t::EtSum::EtSumType::kTotalEt,pt, 0, 0, 0);
463  etsumVec.push_back(etTotal);
464 
465  pt = convertPtToHW( sumEt*0.9, 2047, PtStep_ );
466  l1t::EtSum htTotal(*p4, l1t::EtSum::EtSumType::kTotalHt,pt, 0, 0, 0);
467  etsumVec.push_back(htTotal);
468 
469 
470  // Insert all the bx into the L1 Collections
471  printf("Event %i EmptyBxEvt %i emptyBxTrailer %i diff %i \n",eventCnt_,emptyBxEvt_,emptyBxTrailer_,(emptyBxEvt_ - eventCnt_));
472 
473  // Fill Muons
474  for( int iMu=0; iMu<int(muonVec_bxm2.size()); iMu++ ){
475  muons->push_back(-2, muonVec_bxm2[iMu]);
476  }
477  for( int iMu=0; iMu<int(muonVec_bxm1.size()); iMu++ ){
478  muons->push_back(-1, muonVec_bxm1[iMu]);
479  }
480  for( int iMu=0; iMu<int(muonVec_bx0.size()); iMu++ ){
481  muons->push_back(0, muonVec_bx0[iMu]);
482  }
483  for( int iMu=0; iMu<int(muonVec_bxp1.size()); iMu++ ){
484  muons->push_back(1, muonVec_bxp1[iMu]);
485  }
487  for( int iMu=0; iMu<int(muonVec.size()); iMu++ ){
488  muons->push_back(2, muonVec[iMu]);
489  }
490  } else {
491  // this event is part of empty trailer...clear out data
492  muonVec.clear();
493  }
494 
495  // Fill Egammas
496  for( int iEG=0; iEG<int(egammaVec_bxm2.size()); iEG++ ){
497  egammas->push_back(-2, egammaVec_bxm2[iEG]);
498  }
499  for( int iEG=0; iEG<int(egammaVec_bxm1.size()); iEG++ ){
500  egammas->push_back(-1, egammaVec_bxm1[iEG]);
501  }
502  for( int iEG=0; iEG<int(egammaVec_bx0.size()); iEG++ ){
503  egammas->push_back(0, egammaVec_bx0[iEG]);
504  }
505  for( int iEG=0; iEG<int(egammaVec_bxp1.size()); iEG++ ){
506  egammas->push_back(1, egammaVec_bxp1[iEG]);
507  }
509  for( int iEG=0; iEG<int(egammaVec.size()); iEG++ ){
510  egammas->push_back(2, egammaVec[iEG]);
511  }
512  } else {
513  // this event is part of empty trailer...clear out data
514  egammaVec.clear();
515  }
516 
517  // Fill Taus
518  for( int iTau=0; iTau<int(tauVec_bxm2.size()); iTau++ ){
519  taus->push_back(-2, tauVec_bxm2[iTau]);
520  }
521  for( int iTau=0; iTau<int(tauVec_bxm1.size()); iTau++ ){
522  taus->push_back(-1, tauVec_bxm1[iTau]);
523  }
524  for( int iTau=0; iTau<int(tauVec_bx0.size()); iTau++ ){
525  taus->push_back(0, tauVec_bx0[iTau]);
526  }
527  for( int iTau=0; iTau<int(tauVec_bxp1.size()); iTau++ ){
528  taus->push_back(1, tauVec_bxp1[iTau]);
529  }
531  for( int iTau=0; iTau<int(tauVec.size()); iTau++ ){
532  taus->push_back(2, tauVec[iTau]);
533  }
534  } else {
535  // this event is part of empty trailer...clear out data
536  tauVec.clear();
537  }
538 
539  // Fill Jets
540  for( int iJet=0; iJet<int(jetVec_bxm2.size()); iJet++ ){
541  jets->push_back(-2, jetVec_bxm2[iJet]);
542  }
543  for( int iJet=0; iJet<int(jetVec_bxm1.size()); iJet++ ){
544  jets->push_back(-1, jetVec_bxm1[iJet]);
545  }
546  for( int iJet=0; iJet<int(jetVec_bx0.size()); iJet++ ){
547  jets->push_back(0, jetVec_bx0[iJet]);
548  }
549  for( int iJet=0; iJet<int(jetVec_bxp1.size()); iJet++ ){
550  jets->push_back(1, jetVec_bxp1[iJet]);
551  }
553  for( int iJet=0; iJet<int(jetVec.size()); iJet++ ){
554  jets->push_back(2, jetVec[iJet]);
555  }
556  } else {
557  // this event is part of empty trailer...clear out data
558  jetVec.clear();
559  }
560 
561  // Fill Etsums
562  for( int iETsum=0; iETsum<int(etsumVec_bxm2.size()); iETsum++ ){
563  etsums->push_back(-2, etsumVec_bxm2[iETsum]);
564  }
565  for( int iETsum=0; iETsum<int(etsumVec_bxm1.size()); iETsum++ ){
566  etsums->push_back(-1, etsumVec_bxm1[iETsum]);
567  }
568  for( int iETsum=0; iETsum<int(etsumVec_bx0.size()); iETsum++ ){
569  etsums->push_back(0, etsumVec_bx0[iETsum]);
570  }
571  for( int iETsum=0; iETsum<int(etsumVec_bxp1.size()); iETsum++ ){
572  etsums->push_back(1, etsumVec_bxp1[iETsum]);
573  }
575  for( int iETsum=0; iETsum<int(etsumVec.size()); iETsum++ ){
576  etsums->push_back(2, etsumVec[iETsum]);
577  }
578  } else {
579  // this event is part of empty trailer...clear out data
580  etsumVec.clear();
581  }
582 
583 
584  iEvent.put(egammas);
585  iEvent.put(muons);
586  iEvent.put(taus);
587  iEvent.put(jets);
588  iEvent.put(etsums);
589 
590  // Now shift the bx data by one to prepare for next event.
596 
602 
608 
609  muonVec_bxp1 = muonVec;
610  egammaVec_bxp1 = egammaVec;
611  tauVec_bxp1 = tauVec;
612  jetVec_bxp1 = jetVec;
613  etsumVec_bxp1 = etsumVec;
614 
615 }
#define LogDebug(id)
double MaxEt_
std::vector< l1t::EtSum > etsumVec_bxm1
std::vector< l1t::Muon > muonVec_bxm1
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::vector< l1t::EtSum > etsumVec_bx0
virtual float eta() const =0
momentum pseudorapidity
std::vector< l1t::Muon > muonVec_bxp1
virtual int status() const =0
status word
Definition: Tau.h:13
std::vector< l1t::Tau > tauVec_bx0
std::vector< l1t::Tau > tauVec_bxp1
std::vector< l1t::Tau > tauVec_bxm2
std::vector< l1t::Muon > muonVec_bxm2
virtual float phi() const =0
momentum azimuthal angle
std::vector< l1t::EtSum > etsumVec_bxp1
T eta() const
double charge(const std::vector< uint8_t > &Ampls)
double MaxMuonEta_
std::vector< l1t::EGamma > egammaVec_bxm1
Definition: Jet.h:13
double MaxJetPt_
int convertPtToHW(double ipt, int maxPt, double step)
virtual float pt() const =0
transverse momentum
double PhiStepCalo_
double EtaStepMuon_
std::vector< l1t::EGamma > egammaVec_bxm2
std::vector< l1t::EGamma > egammaVec_bxp1
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::EDGetTokenT< reco::GenJetCollection > genJetsToken
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
double PtStep_
double MaxCaloEta_
double PhiStepMuon_
virtual int charge() const =0
electric charge
const int mu
Definition: Constants.h:22
std::vector< l1t::Jet > jetVec_bxp1
std::vector< l1t::EGamma > egammaVec_bx0
#define LogTrace(id)
int convertPhiToHW(double iphi, int steps)
Definition: Muon.h:12
int k[5][pyjets_maxn]
virtual int pdgId() const =0
PDG identifier.
double MaxLepPt_
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
edm::EDGetTokenT< reco::GenMETCollection > genMetToken
tuple muons
Definition: patZpeak.py:38
unsigned int convertEtaToHW(double ieta, double minEta, double maxEta, int steps, unsigned int bitMask)
edm::EDGetTokenT< reco::GenParticleCollection > genParticlesToken
tuple status
Definition: ntuplemaker.py:245
std::vector< l1t::Jet > jetVec_bxm2
std::vector< l1t::EtSum > etsumVec_bxm2
std::vector< l1t::Jet > jetVec_bx0
Definition: DDAxes.h:10
double EtaStepCalo_

Member Data Documentation

int l1t::GenToInputProducer::bxFirst_
private

Definition at line 82 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::bxLast_
private

Definition at line 83 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::counter_
private

Definition at line 105 of file GenToInputProducer.cc.

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

Definition at line 114 of file GenToInputProducer.cc.

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

Definition at line 113 of file GenToInputProducer.cc.

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

Definition at line 112 of file GenToInputProducer.cc.

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

Definition at line 115 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::egEtThreshold_
private

Definition at line 92 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::emptyBxEvt_
private

Definition at line 97 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::emptyBxTrailer_
private

Definition at line 96 of file GenToInputProducer.cc.

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

Definition at line 129 of file GenToInputProducer.cc.

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

Definition at line 128 of file GenToInputProducer.cc.

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

Definition at line 127 of file GenToInputProducer.cc.

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

Definition at line 130 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::eventCnt_
private

Definition at line 98 of file GenToInputProducer.cc.

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

Definition at line 102 of file GenToInputProducer.cc.

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

Definition at line 103 of file GenToInputProducer.cc.

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

Definition at line 101 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::jetEtThreshold_
private

Definition at line 90 of file GenToInputProducer.cc.

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

Definition at line 124 of file GenToInputProducer.cc.

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

Definition at line 123 of file GenToInputProducer.cc.

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

Definition at line 122 of file GenToInputProducer.cc.

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

Definition at line 125 of file GenToInputProducer.cc.

unsigned long long l1t::GenToInputProducer::m_paramsCacheId
private

Definition at line 76 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumEGCands_
private

Definition at line 87 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumJetCands_
private

Definition at line 86 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumMuCands_
private

Definition at line 85 of file GenToInputProducer.cc.

int l1t::GenToInputProducer::maxNumTauCands_
private

Definition at line 88 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::muEtThreshold_
private

Definition at line 93 of file GenToInputProducer.cc.

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

Definition at line 109 of file GenToInputProducer.cc.

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

Definition at line 108 of file GenToInputProducer.cc.

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

Definition at line 107 of file GenToInputProducer.cc.

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

Definition at line 110 of file GenToInputProducer.cc.

double l1t::GenToInputProducer::tauEtThreshold_
private

Definition at line 91 of file GenToInputProducer.cc.

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

Definition at line 119 of file GenToInputProducer.cc.

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

Definition at line 118 of file GenToInputProducer.cc.

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

Definition at line 117 of file GenToInputProducer.cc.

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

Definition at line 120 of file GenToInputProducer.cc.