Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 process = cms.Process("Gen")
0004 
0005 process.load("Configuration.StandardSequences.SimulationRandomNumberGeneratorSeeds_cff")
0006 
0007 process.source = cms.Source("LHESource",
0008     fileNames = cms.untracked.vstring('file:../../Pythia6Interface/test/ttbar_5flavours_xqcut20_10TeV.lhe')
0009     # fileNames = cms.untracked.vstring('file:/uscms_data/d2/yarba_j/lhe_for_tests/7TeV_Zbb_run45040_unweighted_events_qcut13_mgPostv2.lhe')
0010     # fileNames = cms.untracked.vstring('file:/uscms_data/d2/yarba_j/lhe_for_tests/7TeV_ttbarjets_run621_unweighted_events_qcut40_mgPost.lhe')
0011     # fileNames = cms.untracked.vstring('file:/storage/local/data1/condor/mrenna/lhe/7TeV_avjets_run50000_unweighted_events_qcut15_mgPost.lhe')
0012     # fileNames = cms.untracked.vstring('file:/storage/local/data1/condor/mrenna/lhe/7TeV_zvv_200_HT_inf_run114000_unweighted_events_qcut20_mgPostv2.lhe')
0013 )
0014 
0015 
0016 # process.load("Configuration.Generator.Hadronizer_MgmMatchTune4C_7TeV_madgraph_pythia8_cff")
0017 
0018 process.generator = cms.EDFilter("Pythia8HadronizerFilter",
0019     maxEventsToPrint = cms.untracked.int32(1),
0020     pythiaPylistVerbosity = cms.untracked.int32(1),
0021     filterEfficiency = cms.untracked.double(1.0),
0022     pythiaHepMCVerbosity = cms.untracked.bool(False),
0023     comEnergy = cms.double(7000.),
0024     jetMatching = cms.untracked.PSet(
0025        scheme = cms.string("MadgraphFastJet"),
0026        mode = cms.string("auto"),   # soup, or "inclusive"/"exclusive"
0027        #
0028        # ATTENTION PLEASE !
0029        # One can set some parameters to -1 to make the tool pock it up from LHE file.
0030        # However, -1 is ONLY possible if a givcen parameter is present in LHE file
0031        # - otherwise the code will throw. 
0032        # So the user should make sure what it is and what she/he wants to do.
0033        #
0034        MEMAIN_etaclmax = cms.double(5.),
0035        MEMAIN_qcut = cms.double(30.),       
0036        MEMAIN_minjets = cms.int32(-1),
0037        MEMAIN_maxjets = cms.int32(-1),
0038        MEMAIN_showerkt = cms.double(0),    # use 1=yes only for pt-ordered showers !
0039        MEMAIN_nqmatch = cms.int32(5),      # PID of the flavor until which the QCD radiation are kept in the matching procedure. 
0040                                            # If nqmatch=4, then all showered partons from b's are NOT taken into account.
0041                            # In many cases the D=5
0042        MEMAIN_excres = cms.string(""),
0043        outTree_flag = cms.int32(1)         # 1=yes, write out the tree for future sanity check
0044     ),    
0045     PythiaParameters = cms.PSet(
0046         pythia8_mg = cms.vstring(''), # this pset is for very initial testing
0047         # this pset below is actually used in large-scale (production-type) tests
0048     processParameters = cms.vstring(
0049             'Main:timesAllowErrors    = 10000', 
0050         'ParticleDecays:limitTau0 = on',
0051             'ParticleDecays:tauMax = 10',
0052         # '15:onMode = off', # tmp turn off tau decays, to process av sample (crash in Tau::decay in Py8)
0053         'Tune:ee 3',
0054         'Tune:pp 5',
0055     'ParticleDecays:sophisticatedTau = 0' ),
0056         parameterSets = cms.vstring('processParameters')
0057     )
0058 )
0059 
0060 process.load("FWCore.MessageLogger.MessageLogger_cfi")
0061 process.MessageLogger = cms.Service("MessageLogger",
0062     cerr = cms.untracked.PSet(
0063         enable = cms.untracked.bool(False)
0064     ),
0065     cout = cms.untracked.PSet(
0066         default = cms.untracked.PSet(
0067             limit = cms.untracked.int32(2)
0068         ),
0069         enable = cms.untracked.bool(True)
0070     )
0071 )
0072 
0073 process.RandomNumberGeneratorService = cms.Service("RandomNumberGeneratorService",
0074     generator = cms.PSet(
0075         initialSeed = cms.untracked.uint32(123456789),
0076     )
0077 )
0078 
0079 process.maxEvents = cms.untracked.PSet( input = cms.untracked.int32(10000) )
0080 
0081 process.GEN = cms.OutputModule("PoolOutputModule",
0082     fileName = cms.untracked.string('Py8Had_mgmatching.root')
0083 )
0084 
0085 process.p = cms.Path(process.generator)
0086 process.outpath = cms.EndPath(process.GEN)
0087 
0088 #process.schedule = cms.Schedule(process.p, process.outpath)
0089 process.schedule = cms.Schedule(process.p)
0090