Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:13:42

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 
0004 from Configuration.Eras.Era_Run2_2017_cff import Run2_2017
0005 process = cms.Process('GEN2',Run2_2017)
0006 
0007 # import of standard configurations
0008 process.load('Configuration.StandardSequences.Services_cff')
0009 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0010 process.load('FWCore.MessageService.MessageLogger_cfi')
0011 process.load('Configuration.EventContent.EventContent_cff')
0012 process.load('SimGeneral.MixingModule.mixNoPU_cfi')
0013 process.load('Configuration.StandardSequences.GeometryRecoDB_cff')
0014 process.load('Configuration.StandardSequences.MagneticField_cff')
0015 process.load('Configuration.StandardSequences.Generator_cff')
0016 process.load('IOMC.EventVertexGenerators.VtxSmearedRealistic25ns13TeVEarly2017Collision_cfi')
0017 process.load('GeneratorInterface.Core.genFilterSummary_cff')
0018 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0019 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0020 
0021 #----------
0022 # Input source
0023 #----------
0024 process.source = cms.Source("LHESource",
0025     # first few events from /eos/cms/store/lhe/5542/GJets_HT-200To400_8TeV-madgraph_10001.lhe
0026     fileNames = cms.untracked.vstring("file:GJets_HT-200To400_8TeV-madgraph.lhe")
0027 )
0028 
0029 process.options = cms.untracked.PSet()
0030 
0031 #----------
0032 # Output definition
0033 #----------
0034 process.output = cms.OutputModule("PoolOutputModule",
0035     fileName = cms.untracked.string('file:double-em-enrichment-test.root'),
0036 )
0037 
0038 process.output_step = cms.EndPath(process.output)
0039 
0040 #----------
0041 # PYTHIA hadronizer
0042 # with double EM enrichment filter
0043 #----------
0044 from Configuration.Generator.doubleEMEnrichingHepMCfilter_cfi import *
0045 from Configuration.Generator.Pythia8CommonSettings_cfi import *
0046 from Configuration.Generator.MCTunes2017.PythiaCP5Settings_cfi import *
0047 
0048 process.generator = cms.EDFilter("Pythia8HadronizerFilter",
0049 
0050     HepMCFilter = cms.PSet(
0051         filterName = cms.string('PythiaHepMCFilterGammaGamma'),
0052 
0053         # double EM enrichment filter parameters
0054         filterParameters = doubleEMenrichingHepMCfilterParams
0055     ),
0056                  
0057     PythiaParameters = cms.PSet(
0058         pythia8CommonSettingsBlock,
0059         pythia8CP5SettingsBlock,
0060    
0061         parameterSets = cms.vstring('pythia8CommonSettings', 
0062              'pythia8CP5Settings', 
0063              'processParameters'),
0064 
0065         # e.g. from /GJets_HT-40To100_TuneCP5_13TeV-madgraphMLM-pythia8 configuration
0066         processParameters = cms.vstring('JetMatching:setMad = off', 
0067             'JetMatching:scheme = 1', 
0068             'JetMatching:merge = on', 
0069             'JetMatching:jetAlgorithm = 2', 
0070             'JetMatching:etaJetMax = 5.', 
0071             'JetMatching:coneRadius = 1.', 
0072             'JetMatching:slowJetPower = 1', 
0073             'JetMatching:qCut = 19.', 
0074             'JetMatching:nQmatch = 5', 
0075             'JetMatching:nJetMax = 4', 
0076             'JetMatching:doShowerKt = off'),
0077         ),
0078 
0079         comEnergy = cms.double(13000.0),
0080         filterEfficiency = cms.untracked.double(1.0),
0081         maxEventsToPrint = cms.untracked.int32(0),
0082         
0083         # number of hadronization attempts per LHE event
0084         nAttempts = cms.uint32(20), 
0085         
0086         pythiaHepMCVerbosity = cms.untracked.bool(False),
0087         pythiaPylistVerbosity = cms.untracked.int32(1),
0088 )
0089 
0090 
0091 # Path and EndPath definitions
0092 process.generation_step = cms.Path(process.pgen)
0093 process.genfiltersummary_step = cms.EndPath(process.genFilterSummary)
0094 process.endjob_step = cms.EndPath(process.endOfProcess)
0095 
0096 # Schedule definition
0097 process.schedule = cms.Schedule(process.generation_step,
0098                                 process.genfiltersummary_step,
0099                                 process.endjob_step,
0100                                 process.output_step,
0101 )
0102 
0103 from PhysicsTools.PatAlgos.tools.helpers import associatePatAlgosToolsTask
0104 associatePatAlgosToolsTask(process)
0105 
0106 for path in process.paths:
0107     if path in ['lhe_step']: continue
0108     getattr(process,path)._seq = process.generator * getattr(process,path)._seq 
0109 
0110 process.maxEvents = cms.untracked.PSet( input = cms.untracked.int32(100))