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 from Configuration.Generator.Pythia8CommonSettings_cfi import *
0003 from Configuration.Generator.MCTunes2017.PythiaCP5Settings_cfi import *
0004 
0005 generator = cms.EDFilter("Pythia8GeneratorFilter",
0006                          pythiaPylistVerbosity = cms.untracked.int32(0),
0007                          pythiaHepMCVerbosity = cms.untracked.bool(False),
0008                          comEnergy = cms.double(14000.0),
0009                          ##crossSection = cms.untracked.double(54000000000), # Given by PYTHIA after running
0010                          ##filterEfficiency = cms.untracked.double(0.004), # Given by PYTHIA after running
0011                          maxEventsToPrint = cms.untracked.int32(0),
0012                          ExternalDecays = cms.PSet(
0013         EvtGen130 = cms.untracked.PSet(
0014             decay_table = cms.string('GeneratorInterface/EvtGenInterface/data/DECAY_2014_NOLONGLIFE.DEC'),
0015             particle_property_file = cms.FileInPath('GeneratorInterface/EvtGenInterface/data/evt_2014.pdl'),
0016             user_decay_file = cms.vstring('GeneratorInterface/ExternalDecays/data/Bd_Kstarmumu_Kpi.dec'),
0017             list_forced_decays = cms.vstring('MyB0','Myanti-B0'),
0018             convertPythiaCodes = cms.untracked.bool(False),
0019             operates_on_particles = cms.vint32()
0020             ),
0021         parameterSets = cms.vstring('EvtGen130')
0022         ),
0023                          PythiaParameters = cms.PSet(
0024         pythia8CommonSettingsBlock,
0025         pythia8CP5SettingsBlock,
0026         ## check this (need extra parameters?)
0027         processParameters = cms.vstring('SoftQCD:nonDiffractive = on',
0028                                         'PTFilter:filter = on', # this turn on the filter
0029                                         'PTFilter:quarkToFilter = 5', # PDG id of q quark (can be any other)
0030                                         'PTFilter:scaleToFilter = 1.0'),
0031         parameterSets = cms.vstring('pythia8CommonSettings',
0032                                     'pythia8CP5Settings',
0033                                     'processParameters',
0034         )
0035         )
0036 )
0037 
0038 
0039 ###########
0040 # Filters #
0041 ###########
0042 
0043 bfilter = cms.EDFilter(
0044     "PythiaFilter", 
0045     MaxEta = cms.untracked.double(9999.),
0046     MinEta = cms.untracked.double(-9999.),
0047     ParticleID = cms.untracked.int32(511)  ## Bd
0048     )
0049 
0050 decayfilter = cms.EDFilter(
0051     "PythiaDauVFilter",
0052     verbose         = cms.untracked.int32(1),
0053     NumberDaughters = cms.untracked.int32(3),
0054     ParticleID      = cms.untracked.int32(511),
0055     DaughterIDs     = cms.untracked.vint32(-13, 13, 313),  ## mu+, mu-, K*^0(892)
0056     MinPt           = cms.untracked.vdouble(2.5, 2.5, -1.),
0057     MinEta          = cms.untracked.vdouble(-2.5, -2.5, -9999.),
0058     MaxEta          = cms.untracked.vdouble( 2.5,  2.5,  9999.)
0059     )
0060 
0061 kstarfilter = cms.EDFilter(
0062     "PythiaDauVFilter",
0063     verbose         = cms.untracked.int32(1), 
0064     NumberDaughters = cms.untracked.int32(2), 
0065     MotherID        = cms.untracked.int32(511),  ## Bd
0066     ParticleID      = cms.untracked.int32(313),  ## K*^0(892)
0067     DaughterIDs     = cms.untracked.vint32(321, -211), ## K+, pi-
0068     MinPt           = cms.untracked.vdouble(0.4, 0.4), 
0069     MinEta          = cms.untracked.vdouble(-2.5, -2.5), 
0070     MaxEta          = cms.untracked.vdouble( 2.5,  2.5)
0071     )
0072 
0073 
0074 ProductionFilterSequence = cms.Sequence(generator*bfilter*decayfilter*kstarfilter)