Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:03:35

0001 import FWCore.ParameterSet.Config as cms
0002 from Configuration.Generator.Pythia8CommonSettings_cfi import *
0003 from Configuration.Generator.Pythia8CUEP8M1Settings_cfi import *
0004 from GeneratorInterface.EvtGenInterface.EvtGenSetting_cff import *
0005 
0006 
0007 generator = cms.EDFilter("Pythia8GeneratorFilter",
0008                          pythiaPylistVerbosity = cms.untracked.int32(0),
0009                          filterEfficiency = cms.untracked.double(0.53),
0010                          pythiaHepMCVerbosity = cms.untracked.bool(False),
0011                          crossSection = cms.untracked.double(9090000.0),
0012                          comEnergy = cms.double(13000.0),
0013                          maxEventsToPrint = cms.untracked.int32(0),
0014                          ExternalDecays = cms.PSet(
0015         EvtGen130 = cms.untracked.PSet(
0016             decay_table = cms.string('GeneratorInterface/EvtGenInterface/data/DECAY_2010.DEC'),
0017             particle_property_file = cms.FileInPath('GeneratorInterface/EvtGenInterface/data/evt.pdl'),
0018             user_decay_file = cms.vstring('GeneratorInterface/ExternalDecays/data/LambdaB_Lambdamumu_ppi.dec'),
0019             list_forced_decays = cms.vstring('MyLambda_b0','Myanti-Lambda_b0','MyLambda0','Myanti-Lambda0'),
0020             operates_on_particles = cms.vint32()
0021             ),
0022         parameterSets = cms.vstring('EvtGen130')
0023         ),
0024                          PythiaParameters = cms.PSet(
0025         pythia8CommonSettingsBlock,
0026         pythia8CUEP8M1SettingsBlock,
0027         processParameters = cms.vstring(
0028             'HardQCD:gg2bbbar = on ',
0029             'HardQCD:qqbar2bbbar = on ',
0030             'HardQCD:hardbbbar = on',
0031             'PhaseSpace:pTHatMin = 20.',
0032             ),
0033         parameterSets = cms.vstring('pythia8CommonSettings',
0034                                     'pythia8CUEP8M1Settings',
0035                                     'processParameters',
0036                                     )
0037         )
0038                          )
0039 
0040 generator.PythiaParameters.processParameters.extend(EvtGenExtraParticles)
0041 
0042 configurationMetadata = cms.untracked.PSet(
0043     version = cms.untracked.string('$Revision: 1.1 $'),
0044     name = cms.untracked.string('$Source: Configuration/Generator/python/BuToKstarMuMu_forSTEAM_13TeV_TuneCUETP8M1_cfi.py $'),
0045     annotation = cms.untracked.string('Summer14: Pythia8+EvtGen130 generation of lambda_b --> lambda mumu ->ppimumu, 13TeV, Tune CUETP8M1')
0046     )
0047 
0048 ###########
0049 # Filters #
0050 ###########
0051 # Filter only pp events which produce a lambda_b:
0052 lambdabfilter = cms.EDFilter("PythiaFilter", ParticleID = cms.untracked.int32(5122))
0053 
0054 # Filter on final state muons
0055 mumugenfilter = cms.EDFilter("MCParticlePairFilter",
0056                              Status = cms.untracked.vint32(1, 1),
0057                              MinPt = cms.untracked.vdouble(2.8, 2.8),
0058                              MinP = cms.untracked.vdouble(2.8, 2.8),
0059                              MaxEta = cms.untracked.vdouble(2.3, 2.3),
0060                              MinEta = cms.untracked.vdouble(-2.3, -2.3),
0061                              ParticleID1 = cms.untracked.vint32(13,-13),
0062                              ParticleID2 = cms.untracked.vint32(13,-13)
0063                              )
0064 
0065 
0066 ProductionFilterSequence = cms.Sequence(generator*lambdabfilter*mumugenfilter)