1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
|
import FWCore.ParameterSet.Config as cms
from Configuration.Generator.Pythia8CommonSettings_cfi import *
from Configuration.Generator.MCTunes2017.PythiaCP5Settings_cfi import *
from GeneratorInterface.EvtGenInterface.EvtGenSetting_cff import *
generator = cms.EDFilter("Pythia8GeneratorFilter",
pythiaHepMCVerbosity = cms.untracked.bool(False),
maxEventsToPrint = cms.untracked.int32(0),
pythiaPylistVerbosity = cms.untracked.int32(0),
filterEfficiency = cms.untracked.double(1.38e-3),
crossSection = cms.untracked.double(540000000.),
comEnergy = cms.double(14000.0),
ExternalDecays = cms.PSet(
EvtGen130 = cms.untracked.PSet(
decay_table = cms.string('GeneratorInterface/EvtGenInterface/data/DECAY_2010.DEC'),
particle_property_file = cms.FileInPath('GeneratorInterface/EvtGenInterface/data/evt.pdl'),
user_decay_embedded= cms.vstring(
'#',
'Alias MyBs B_s0',
'Alias Myanti-Bs anti-B_s0',
'ChargeConj Myanti-Bs MyBs',
'Alias MyJpsi J/psi',
'ChargeConj MyJpsi MyJpsi',
'#',
'Decay MyBs',
'1.000 MyJpsi gamma SVP_HELAMP 1.0 0.0 1.0 0.0;',
'Enddecay',
'CDecay Myanti-Bs',
'#',
'Decay MyJpsi',
'1.000 mu+ mu- PHOTOS VLL;',
'Enddecay',
'End'
),
list_forced_decays = cms.vstring('MyBs','Myanti-Bs'),
operates_on_particles = cms.vint32(),
),
parameterSets = cms.vstring('EvtGen130')
),
PythiaParameters = cms.PSet(pythia8CommonSettingsBlock,
pythia8CP5SettingsBlock,
processParameters = cms.vstring(
"SoftQCD:nonDiffractive = on",
'PTFilter:filter = on', # this turn on the filter
'PTFilter:quarkToFilter = 5', # PDG id of q quark
'PTFilter:scaleToFilter = 1.0'
),
parameterSets = cms.vstring('pythia8CommonSettings',
'pythia8CP5Settings',
'processParameters',
)
)
)
generator.PythiaParameters.processParameters.extend(EvtGenExtraParticles)
configurationMetadata = cms.untracked.PSet(
version = cms.untracked.string('$Revision: 1.0 $'),
name = cms.untracked.string('$Source: /Configuration/Generator/python/BsJpsiGamma_PhaseII_cfi.py $'),
annotation = cms.untracked.string('PhaseII: Pythia8+EvtGen130 generation of Bs --> Jpsi gamma, 14TeV, Tune CP5')
)
bfilter = cms.EDFilter(
"PythiaFilter",
MaxEta = cms.untracked.double(9999.),
MinEta = cms.untracked.double(-9999.),
ParticleID = cms.untracked.int32(531)
)
decayfilter = cms.EDFilter(
"PythiaDauVFilter",
verbose = cms.untracked.int32(1),
NumberDaughters = cms.untracked.int32(2),
ParticleID = cms.untracked.int32(531),
DaughterIDs = cms.untracked.vint32(443, 22), ## jpsi gamma
MinPt = cms.untracked.vdouble(3.5, 4.0),
MinEta = cms.untracked.vdouble(-3.5, -2.5),
MaxEta = cms.untracked.vdouble( 3.5, 2.5)
)
jpsifilter = cms.EDFilter(
"PythiaDauVFilter",
MotherID = cms.untracked.int32(531),
ParticleID = cms.untracked.int32(443),
NumberDaughters = cms.untracked.int32(2),
DaughterIDs = cms.untracked.vint32(13, -13),
MinPt = cms.untracked.vdouble(2.5, 2.5),
MinEta = cms.untracked.vdouble(-2.9, -2.9),
MaxEta = cms.untracked.vdouble(2.9, 2.9),
verbose = cms.untracked.int32(1)
)
ProductionFilterSequence = cms.Sequence(generator*decayfilter*jpsifilter)
|