Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-03 04:17:53

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 process = cms.Process("TQAF")
0004 
0005 process.load("FWCore.MessageLogger.MessageLogger_cfi")
0006 process.load("Configuration.StandardSequences.SimulationRandomNumberGeneratorSeeds_cff")
0007 
0008 process.source = cms.Source("EmptySource")
0009 
0010 process.generator = cms.EDFilter("Pythia8GeneratorFilter",
0011     maxEventsToPrint = cms.untracked.int32(1),
0012     pythiaPylistVerbosity = cms.untracked.int32(1),
0013     filterEfficiency = cms.untracked.double(1.0),
0014     pythiaHepMCVerbosity = cms.untracked.bool(True),
0015     comEnergy = cms.double(13000.),
0016     PythiaParameters = cms.PSet(
0017         processParameters = cms.vstring('Top:all = on'),
0018         parameterSets = cms.vstring('processParameters')
0019     )
0020 )
0021 
0022 process.genParticles = cms.EDProducer("GenParticleProducer",
0023     abortOnUnknownPDGCode = cms.untracked.bool(False),
0024     saveBarCodes = cms.untracked.bool(True),
0025     src = cms.InputTag("generator:unsmeared")
0026 )
0027 
0028 
0029 ## define maximal number of events to loop over
0030 process.maxEvents = cms.untracked.PSet(
0031     input = cms.untracked.int32(100)
0032 )
0033 ## configure process options
0034 process.options = cms.untracked.PSet(
0035     allowUnscheduled = cms.untracked.bool(True),
0036     wantSummary      = cms.untracked.bool(True)
0037 )
0038 
0039 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0040 process.load("GeneratorInterface.RivetInterface.genParticles2HepMC_cfi")
0041 process.load("GeneratorInterface.RivetInterface.particleLevel_cfi")
0042 process.particleLevel.src = cms.InputTag("genParticles2HepMC:unsmeared")
0043 
0044 process.path = cms.Path(process.generator*process.genParticles*process.genParticles2HepMC*process.particleLevel)
0045 
0046 process.out = cms.OutputModule("PoolOutputModule",
0047     fileName = cms.untracked.string("particleLevel.root"),
0048     outputCommands = cms.untracked.vstring(
0049         "drop *",
0050         "keep *_genParticles_*_*",
0051         "keep *_particleLevel_*_*",
0052     ),
0053 )
0054 process.outPath = cms.EndPath(process.out)