CMS 3D CMS Logo

StGenEventReco.cc
Go to the documentation of this file.
2 
4  : srcToken_(consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("src"))),
5  initToken_(consumes<reco::GenParticleCollection>(cfg.getParameter<edm::InputTag>("init"))) {
6  produces<StGenEvent>();
7 }
8 
10 
13  evt.getByToken(srcToken_, parts);
14 
16  evt.getByToken(initToken_, inits);
17 
18  //add TopDecayTree
20 
21  //add InitialStatePartons
22  reco::GenParticleRefProd initParts(inits);
23 
24  //add genEvt to the output stream
25  StGenEvent* genEvt = new StGenEvent(cands, initParts);
26  std::unique_ptr<StGenEvent> gen(genEvt);
27  evt.put(std::move(gen));
28 }
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
edm::EDGetTokenT< reco::GenParticleCollection > initToken_
edm::EDGetTokenT< reco::GenParticleCollection > srcToken_
~StGenEventReco() override
Class derived from the TopGenEvent for single-top events.
Definition: StGenEvent.h:17
def gen(fragment, howMuch)
Production test section ####.
fixed size matrix
HLT enums.
StGenEventReco(const edm::ParameterSet &)
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &, const edm::EventSetup &) override