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.GeneratorSettingsSemantics import (
6 GeneratorSettingsLayer,
7 GeneratorSettingsPrecedence,
8)
9from GeneratorConfig.Sequences import EvgenSequence, EvgenSequenceFactory
10from Pythia8_i.Pythia8Tunes import (
11 a14_nnpdf23lo_tune_cmds,
12 a2_mstw2008lo_tune_cmds,
13)
14from AthenaCommon.SystemOfUnits import GeV
15
16# Get logger
17from AthenaCommon.Logging import logging
18log = logging.getLogger("Pythia8Config")
19
20
21def Pythia8CommandsCfg(flags, source, commands, precedence, name="Pythia8_i"):
22 """
23 Return a CA fragment that adds one command layer to Pythia8_i.
24 """
25 ca = ComponentAccumulator(EvgenSequenceFactory(EvgenSequence.Generator))
26 ca.addEventAlgo(
27 CompFactory.Pythia8_i(name, Commands=GeneratorSettingsLayer(
28 source=source,
29 values=tuple(commands or ()),
30 precedence=precedence,
31 report_context="Pythia8Cfg.Commands",
32 ))
33 )
34 return ca
35
36
37def Pythia8BaseCfg(flags, name="Pythia8_i", **kwargs):
38 """
39 The main Pythia8 configuration fragment that sets up the Pythia8_i algorithm
40 and returns a CA instance
41 """
42 # Baseline P8 settings
43 base_cmds = [
44 "Main:timesAllowErrors = 500",
45 "ParticleDecays:limitTau0 = on",
46 "ParticleDecays:tau0Max = 10.0"
47 ]
48
49 # Collision energy
50 if "CollisionEnergy" not in kwargs:
51 kwargs["CollisionEnergy"] = flags.Beam.Energy * 2 / GeV
52
53 # Load basic parameters
54 base_cmds.extend([
55 "6:m0 = 172.5",
56 "23:m0 = 91.1876",
57 "23:mWidth = 2.4952",
58 "24:m0 = 80.399",
59 "24:mWidth = 2.085",
60 "StandardModel:sin2thetaW = 0.23113",
61 "StandardModel:sin2thetaWbar = 0.23146",
62 ])
63
64 user_cmds = kwargs.pop("Commands", None)
65 kwargs["Commands"] = GeneratorSettingsLayer(
66 source="base_fragment_commands",
67 values=tuple(base_cmds),
68 precedence=GeneratorSettingsPrecedence.BASE,
69 report_context="Pythia8Cfg.Commands",
70 )
71
72 # Create CA object
73 ca = ComponentAccumulator(EvgenSequenceFactory(EvgenSequence.Generator))
74 ca.addEventAlgo(
75 CompFactory.Pythia8_i(name, **kwargs)
76 )
77
78 # Add the user commands
79 if user_cmds:
80 ca.merge(Pythia8CommandsCfg(
81 flags,
82 source="user_commands",
83 commands=user_cmds,
84 precedence=GeneratorSettingsPrecedence.USER,
85 name=name,
86 ))
87
88 # Announce generator to service
89 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
90 ca.merge(GeneratorInfoSvcCfg(flags, Generators=["Pythia8"]), sequenceName=EvgenSequence.Generator.value)
91
92 return ca
93
94
95def Pythia8EvtGenBaseCfg(flags, **kwargs):
96 """
97 Fragment for setting up EvtGen on top of Pythia 8
98 """
99 # Custom settings
100 auxfiles = ["inclusiveP8DsDPlus.pdt"]
101 # FHerwig has problems with omega b* (5334), so not present in the base EvtGen fragment.
102 whiteList = [-5334, 5334]
103
104 # Add custom EvtGenCfg
105 from EvtGen_i.EvtGenConfig import EvtGenCfg
106 ca = EvtGenCfg(
107 flags,
108 whiteList = whiteList,
109 auxfiles = auxfiles
110 )
111
112 return ca
113
114
116 """
117 Fragment for setting up A2 MSTW2008LO tune
118 """
119
120 # Remove any user command before calling base config
121 user_cmds = list(kwargs.pop("Commands", []))
122
123 # Get the base config
124 ca = Pythia8BaseCfg(flags, **kwargs)
125
126 # Get the tune commands and apply rapidity ordering
127 tune_cmds = a2_mstw2008lo_tune_cmds()
128 tune_cmds = ensureRapidityOrderMPI(tune_cmds)
129
130 # Add the tune commands
131 ca.merge(Pythia8CommandsCfg(
132 flags,
133 source="pythia_tune_A2_MSTW2008LO",
134 commands=tune_cmds,
135 precedence=GeneratorSettingsPrecedence.TUNE,
136 ))
137
138 # Now apply the user commands
139 if user_cmds:
140 ca.merge(Pythia8CommandsCfg(
141 flags,
142 source="job_options",
143 commands=user_cmds,
144 precedence=GeneratorSettingsPrecedence.USER,
145 ))
146
147 # Broadcast tune to service
148 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
149 ca.merge(GeneratorInfoSvcCfg(flags, Tune="A2 MSTW2008LO"), sequenceName=EvgenSequence.Generator.value)
150
151 # Call the base config
152 return ca
153
154
156 """
157 Fragment for setting up A14 tune with NNPDF23LO PDF
158 """
159
160 # Remove any user command before calling base config
161 user_cmds = list(kwargs.pop("Commands", []))
162
163 # Get the base config
164 ca = Pythia8BaseCfg(flags, **kwargs)
165
166 # Get the tune commands and apply rapidity ordering
167 tune_cmds = a14_nnpdf23lo_tune_cmds()
168 tune_cmds = ensureRapidityOrderMPI(tune_cmds)
169
170 # Add the tune commands
171 ca.merge(Pythia8CommandsCfg(
172 flags,
173 source="pythia_tune_A14_NNPDF23LO",
174 commands=tune_cmds,
175 precedence=GeneratorSettingsPrecedence.TUNE,
176 ))
177
178 # Now apply the user commands
179 if user_cmds:
180 ca.merge(Pythia8CommandsCfg(
181 flags,
182 source="job_options",
183 commands=user_cmds,
184 precedence=GeneratorSettingsPrecedence.USER,
185 ))
186
187 # Broadcast tune to service
188 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
189 ca.merge(GeneratorInfoSvcCfg(flags, Tune="A14 NNPDF23LO"), sequenceName=EvgenSequence.Generator.value)
190
191 # Call the base config
192 return ca
193
194
196 """
197 Config for Py8 tune A2 with MSTW2008LO tune
198 The default version of this includes EvtGen for standardised b fragmentation
199 This tune is generally only used for pile up samples at the start of run 2
200 for high pT physics at the start of run 2 the A14 tune is more appropriate.
201 There are also more recent soft QCD tunes, such as Monash,
202 but A2 was a conservative choice for initial 13 TeV pile up
203 """
204
205 # Add Pythia 8 to CA with correct tune settings
206 ca = Pythia8_A2_MSTW2008LO_Common_Cfg(flags, **kwargs)
207
208 # Add EvtGen
209 ca.merge(Pythia8EvtGenBaseCfg(flags, **kwargs))
210
211 return ca
212
213
215 """
216 Config for setting up Py8 with A14 tune
217 with EvtGen
218 """
219
220 # Add Pythia 8 to CA with correct tune settings
221 ca = Pythia8_A14_NNPDF23LO_Common_Cfg(flags, **kwargs)
222
223 # Add EvtGen
224 ca.merge(Pythia8EvtGenBaseCfg(flags, **kwargs))
225
226 return ca
227
228
230 """
231 A function that ensures rapidity ordering is set
232 """
233 # if MPI‐ordering already explicitly set, do nothing
234 if any("SpaceShower:rapidityOrderMPI" in c for c in cmds):
235 return cmds
236
237 # find the first rapidityOrder value
238 for c in cmds:
239 if "SpaceShower:rapidityOrder" in c and "MPI" not in c:
240 val = c.split("=", 1)[-1].strip()
241 cmds.append(f"SpaceShower:rapidityOrderMPI = {val}")
242 break
243 return cmds
244
245
246def Pythia8_MadGraph_Cfg(flags, ShowerCfg=Pythia8BaseCfg, **kwargs):
247 """
248 Modular fragment for MadGraph LHE input in Pythia8.
249 The Pythia8_i algorithm is configured through ShowerCfg (defaults to
250 Pythia8BaseCfg) so tune/EvtGen fragments can be injected without
251 instantiating Pythia8_i twice.
252 """
253 # Set LHE file name (override with LHEFile="myfile.lhe" if needed).
254 kwargs.setdefault("LHEFile", "events.lhe")
255
256 # Configure Pythia8 through the selected shower fragment.
257 ca = ShowerCfg(flags, **kwargs)
258
259 # Announce MadGraph to service
260 from GeneratorConfig.GeneratorInfoSvcConfig import GeneratorInfoSvcCfg
261 ca.merge(GeneratorInfoSvcCfg(flags, Generators=["MadGraph"]), sequenceName=EvgenSequence.Generator.value)
262
263 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)
Pythia8CommandsCfg(flags, source, commands, precedence, name="Pythia8_i")
Pythia8_A2_MSTW2008LO_Common_Cfg(flags, **kwargs)