Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 from Configuration.Generator.Pythia8CommonSettings_cfi import *
0004 from Configuration.Generator.Pythia8CUEP8M1Settings_cfi import *
0005 from GeneratorInterface.EvtGenInterface.EvtGenSetting_cff import *
0006 
0007 generator = cms.EDFilter("Pythia8GeneratorFilter",
0008                          pythiaHepMCVerbosity = cms.untracked.bool(False),
0009                          maxEventsToPrint = cms.untracked.int32(0),
0010                          pythiaPylistVerbosity = cms.untracked.int32(0),
0011                          filterEfficiency = cms.untracked.double(1.38e-3),
0012                          crossSection = cms.untracked.double(540000000.),
0013                          comEnergy = cms.double(13000.0),
0014 
0015                          ExternalDecays = cms.PSet(
0016         EvtGen130 = cms.untracked.PSet(
0017             decay_table = cms.string('GeneratorInterface/EvtGenInterface/data/DECAY_2010.DEC'),
0018             particle_property_file = cms.FileInPath('GeneratorInterface/EvtGenInterface/data/evt_Bmm.pdl'),
0019             user_decay_file = cms.vstring('GeneratorInterface/EvtGenInterface/data/Bu_JpsiK.dec'),
0020             list_forced_decays = cms.vstring('MyB+',
0021                                              'MyB-'),
0022        #     operates_on_particles = cms.vint32(),
0023             ),
0024         parameterSets = cms.vstring('EvtGen130')
0025         ),
0026                          
0027                          PythiaParameters = cms.PSet(pythia8CommonSettingsBlock,
0028                                                      pythia8CUEP8M1SettingsBlock,
0029                                                      processParameters = cms.vstring("SoftQCD:nonDiffractive = on"),
0030                                                      parameterSets = cms.vstring('pythia8CommonSettings',
0031                                                                                  'pythia8CUEP8M1Settings',
0032                                                                                  'processParameters',
0033                                                                                  )
0034                                                      )
0035                          )
0036 
0037 generator.PythiaParameters.processParameters.extend(EvtGenExtraParticles)
0038 
0039 configurationMetadata = cms.untracked.PSet(
0040     version = cms.untracked.string('$Revision: 1.1 $'),
0041     name = cms.untracked.string('$Source: Configuration/Generator/python/PYTHIA8_BuJpsiK_EtaPtFilter_CUEP8M1_13TeV_cff.py $'),
0042     annotation = cms.untracked.string('Spring 2015: Pythia8+EvtGen130 generation of B+ --> J/psi K+, 13TeV, Tune CUETP8M1')
0043     )
0044 
0045 bfilter = cms.EDFilter(
0046     "PythiaFilter",
0047     MaxEta = cms.untracked.double(9999.),
0048     MinEta = cms.untracked.double(-9999.),
0049     ParticleID = cms.untracked.int32(521)
0050     )
0051 
0052 jpsifilter = cms.EDFilter(
0053     "PythiaDauVFilter",
0054     verbose         = cms.untracked.int32(1), 
0055     NumberDaughters = cms.untracked.int32(2), 
0056     MotherID        = cms.untracked.int32(521),  
0057     ParticleID      = cms.untracked.int32(443),  
0058     DaughterIDs     = cms.untracked.vint32(13, -13),
0059     MinPt           = cms.untracked.vdouble(2.5, 2.5), 
0060     MinEta          = cms.untracked.vdouble(-2.5, -2.5), 
0061     MaxEta          = cms.untracked.vdouble( 2.5,  2.5)
0062     )
0063 
0064 kfilter = cms.EDFilter(
0065     "PythiaDauVFilter",
0066     verbose         = cms.untracked.int32(1), 
0067     NumberDaughters = cms.untracked.int32(2), 
0068     MotherID        = cms.untracked.int32(0),  
0069     ParticleID      = cms.untracked.int32(521),  
0070     DaughterIDs     = cms.untracked.vint32(443, 321),
0071     MinPt           = cms.untracked.vdouble(-99., 0.4), 
0072     MinEta          = cms.untracked.vdouble(-9999., -2.5), 
0073     MaxEta          = cms.untracked.vdouble(9999.,   2.5)
0074     )
0075 
0076 ProductionFilterSequence = cms.Sequence(generator*bfilter*jpsifilter*kfilter)
0077