ATLAS Offline Software
Loading...
Searching...
No Matches
Pythia8Config.py
Go to the documentation of this file.
1# Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
2
3from AthenaConfiguration.ComponentAccumulator import ComponentAccumulator
4from AthenaConfiguration.ComponentFactory import CompFactory
5from GeneratorConfig.Sequences import EvgenSequence, EvgenSequenceFactory
6from AthenaCommon.SystemOfUnits import GeV
7
8# Get logger
9from AthenaCommon.Logging import logging
10log = logging.getLogger("Pythia8Config")
11
12
13def Pythia8BaseCfg(flags, name="Pythia8_i", **kwargs):
14 """
15 The main Pythia8 configuration fragment that sets up the Pythia8_i algorithm
16 and returns a CA instance
17 """
18 # Default commands
19 commands = kwargs.get("Commands", [])
20
21 # Baseline P8 settings
22 base_cmds = [
23 "Main:timesAllowErrors = 500",
24 "ParticleDecays:limitTau0 = on",
25 "ParticleDecays:tau0Max = 10.0"
26 ]
27 commands.extend(base_cmds)
28
29 # Collision energy
30 if "CollisionEnergy" not in kwargs:
31 kwargs["CollisionEnergy"] = flags.Beam.Energy * 2 / GeV
32
33 # Extended settings
34 if flags.Generator.PDGparams:
35
36 from EvgenProdTools.offline_dict import parameters
37
38
39 particle_params = parameters.get("particles")
40 if particle_params:
41 for pdg_str, vals in particle_params.items():
42
43 pdg = int(pdg_str)
44 if 6 <= pdg < 26:
45 commands.append(f"{pdg}:m0 = {vals['mass']}")
46 commands.append(f"{pdg}:mWidth = {vals['width']}")
47 else:
48 log.warning("Could not retrieve standard ATLAS particle parameters")
49
50
51 ew_params = parameters.get("EW_parameters")
52 if ew_params:
53
54 for key, val in ew_params.items():
55 if key[1] in ('sin2thetaW', 'sin2thetaWbar'):
56 commands.append(f"StandardModel:{key[1]} = {val}")
57 else:
58 log.warning("Could not retrieve standard ATLAS EW parameters")
59 else:
60
61 commands.extend([
62 "6:m0 = 172.5",
63 "23:m0 = 91.1876",
64 "23:mWidth = 2.4952",
65 "24:m0 = 80.399",
66 "24:mWidth = 2.085",
67 "StandardModel:sin2thetaW = 0.23113",
68 "StandardModel:sin2thetaWbar = 0.23146",
69 ])
70
71 # Remove duplicate commands and update kwargs
72 commands = list(dict.fromkeys(commands))
73 kwargs["Commands"] = commands
74
75 # Create CA object
76 ca = ComponentAccumulator(EvgenSequenceFactory(EvgenSequence.Generator))
77 ca.addEventAlgo(
78 CompFactory.Pythia8_i("Pythia8_i", **kwargs)
79 )
80
81 # Announce generator to service
82 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
83 ca.merge(GeneratorInfoSvcCfg(flags, Generators=["Pythia8"]), sequenceName=EvgenSequence.Generator.value)
84
85 return ca
86
87
88def Pythia8EvtGenBaseCfg(flags, **kwargs):
89 """
90 Fragment for setting up EvtGen on top of Pythia 8
91 """
92 # Custom settings
93 auxfiles = ["inclusiveP8DsDPlus.pdt"]
94 # FHerwig has problems with omega b* (5334), so not present in the base EvtGen fragment.
95 whiteList = [-5334, 5334]
96
97 # Add custom EvtGenCfg
98 from EvtGen_i.EvtGenConfig import EvtGenCfg
99 ca = EvtGenCfg(
100 flags,
101 whiteList = whiteList,
102 auxfiles = auxfiles
103 )
104
105 return ca
106
107
109 """
110 Fragment for setting up A2 MSTW2008LO tune
111 """
112 # Defaults
113 cmds = kwargs.get("Commands", [])
114
115 # Tune parameters
116 cmds.extend([
117 "Tune:pp = 5",
118 "MultipartonInteractions:bProfile = 4",
119 "MultipartonInteractions:a1 = 0.03",
120 "MultipartonInteractions:pT0Ref = 1.90",
121 "MultipartonInteractions:ecmPow = 0.30",
122 "SpaceShower:rapidityOrder = 0",
123 "PDF:pSet = LHAPDF6:MSTW2008lo68cl",
124 "ColourReconnection:range = 2.28"
125 ])
126
127 # Now call rapidity ordering
128 cmds = ensureRapidityOrderMPI(cmds)
129
130 # Update kwargs
131 kwargs["Commands"] = list(dict.fromkeys(cmds))
132
133 # Now get the base config
134 ca = Pythia8BaseCfg(flags, **kwargs)
135
136 # Broadcast tune to service
137 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
138 ca.merge(GeneratorInfoSvcCfg(flags, Tune="A2 MSTW2008LO"), sequenceName=EvgenSequence.Generator.value)
139
140 # Call the base config
141 return ca
142
143
145 """
146 Fragment for setting up A14 tune with NNPDF23LO PDF
147 """
148 # Defaults
149 cmds = kwargs.get("Commands", [])
150
151 # Tune parameters
152 cmds.extend([
153 "Tune:ee = 7",
154 "Tune:pp = 14",
155 "SpaceShower:rapidityOrder = on",
156 "SigmaProcess:alphaSvalue = 0.140",
157 "SpaceShower:pT0Ref = 1.56",
158 "SpaceShower:pTmaxFudge = 0.91",
159 "SpaceShower:pTdampFudge = 1.05",
160 "SpaceShower:alphaSvalue = 0.127",
161 "TimeShower:alphaSvalue = 0.127",
162 "BeamRemnants:primordialKThard = 1.88",
163 "MultipartonInteractions:pT0Ref = 2.09",
164 "MultipartonInteractions:alphaSvalue = 0.126",
165 "PDF:pSet=LHAPDF6:NNPDF23_lo_as_0130_qed",
166 "ColourReconnection:range = 1.71"
167 ])
168
169 # Now call rapidity ordering
170 cmds = ensureRapidityOrderMPI(cmds)
171
172 # Update kwargs
173 kwargs["Commands"] = list(dict.fromkeys(cmds))
174
175 # Now get the base config
176 ca = Pythia8BaseCfg(flags, **kwargs)
177
178 # Broadcast tune to service
179 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
180 ca.merge(GeneratorInfoSvcCfg(flags, Tune="A14 NNPDF23LO"), sequenceName=EvgenSequence.Generator.value)
181
182 # Call the base config
183 return ca
184
185
187 """
188 Config for Py8 tune A2 with MSTW2008LO tune
189 The default version of this includes EvtGen for standardised b fragmentation
190 This tune is generally only used for pile up samples at the start of run 2
191 for high pT physics at the start of run 2 the A14 tune is more appropriate.
192 There are also more recent soft QCD tunes, such as Monash,
193 but A2 was a conservative choice for initial 13 TeV pile up
194 """
195
196 # Add Pythia 8 to CA with correct tune settings
197 ca = Pythia8_A2_MSTW2008LO_Common_Cfg(flags, **kwargs)
198
199 # Add EvtGen
200 ca.merge(Pythia8EvtGenBaseCfg(flags, **kwargs))
201
202 return ca
203
204
206 """
207 Config for setting up Py8 with A14 tune
208 with EvtGen
209 """
210
211 # Add Pythia 8 to CA with correct tune settings
212 ca = Pythia8_A14_NNPDF23LO_Common_Cfg(flags, **kwargs)
213
214 # Add EvtGen
215 ca.merge(Pythia8EvtGenBaseCfg(flags, **kwargs))
216
217 return ca
218
219
221 """
222 A function that ensures rapidity ordering is set
223 """
224 # if MPI‐ordering already explicitly set, do nothing
225 if any("SpaceShower:rapidityOrderMPI" in c for c in cmds):
226 return cmds
227
228 # find the first rapidityOrder value
229 for c in cmds:
230 if "SpaceShower:rapidityOrder" in c and "MPI" not in c:
231 val = c.split("=", 1)[-1].strip()
232 cmds.append(f"SpaceShower:rapidityOrderMPI = {val}")
233 break
234 return cmds
235
236
237def Pythia8_MadGraph_Cfg(flags, ShowerCfg=Pythia8BaseCfg, **kwargs):
238 """
239 Modular fragment for MadGraph LHE input in Pythia8.
240 The Pythia8_i algorithm is configured through ShowerCfg (defaults to
241 Pythia8BaseCfg) so tune/EvtGen fragments can be injected without
242 instantiating Pythia8_i twice.
243 """
244 # Set LHE file name (override with LHEFile="myfile.lhe" if needed).
245 kwargs.setdefault("LHEFile", "events.lhe")
246
247 # Configure Pythia8 through the selected shower fragment.
248 ca = ShowerCfg(flags, **kwargs)
249
250 # Announce MadGraph to service
251 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
252 ca.merge(GeneratorInfoSvcCfg(flags, Generators=["MadGraph"]), sequenceName=EvgenSequence.Generator.value)
253
254 return ca
Pythia8_A2_MSTW2008LO_EvtGen_Common_Cfg(flags, **kwargs)
Pythia8_MadGraph_Cfg(flags, ShowerCfg=Pythia8BaseCfg, **kwargs)
Pythia8_A14_NNPDF23LO_Common_Cfg(flags, **kwargs)
Pythia8_A14_NNPDF23LO_EvtGen_Common_Cfg(flags, **kwargs)
Pythia8BaseCfg(flags, name="Pythia8_i", **kwargs)
Pythia8EvtGenBaseCfg(flags, **kwargs)
Pythia8_A2_MSTW2008LO_Common_Cfg(flags, **kwargs)