FairShip
Loading...
Searching...
No Matches
pythia8_conf Namespace Reference

Functions

None configurerpvsusy (P8gen, mass, couplings, sfermionmass, benchmark, inclusive, bool deepCopy=False, bool debug=True)
 
None configure (P8gen, mass, production_couplings, decay_couplings, process_selection, bool deepCopy=False, bool debug=True)
 
None add_hnl (P8gen, mass, decay_couplings)
 
None setup_pythia_inclusive (P8gen)
 

Function Documentation

◆ add_hnl()

None pythia8_conf.add_hnl (   P8gen,
  mass,
  decay_couplings 
)

Definition at line 312 of file pythia8_conf.py.

312def add_hnl(P8gen, mass, decay_couplings) -> None:
313 "Adds the HNL to Pythia and ROOT"
314 hnl_instance = hnl.HNL(mass, decay_couplings, debug=True)
315 ctau = hnl_instance.computeNLifetime(system="FairShip") * u.c_light * u.cm
316 print(f"HNL ctau {ctau}")
317 P8gen.SetParameters(f"9900015:new = N2 N2 2 0 0 {mass:.12} 0.0 0.0 0.0 {ctau / u.mm:.12} 0 1 0 1 0")
318 P8gen.SetParameters("9900015:isResonance = false")
319 # Configuring decay modes...
321 P8gen, hnl_instance, conffile=os.path.expandvars("$FAIRSHIP/python/DecaySelection.conf"), verbose=False
322 )
323 # Finish HNL setup...
324 P8gen.SetParameters("9900015:mayDecay = on")
325 P8gen.SetHNLId(9900015)
326 # also add to PDG
327 gamma = u.hbarc / float(ctau) # 197.3269631e-16 / float(ctau) # hbar*c = 197 MeV*fm = 197e-16 GeV*cm
328 addHNLtoROOT(pid=9900015, m=mass, g=gamma)
329
330
Definition: hnl.py:700
def addHNLdecayChannels(P8Gen, hnl, conffile=os.path.expandvars("$FAIRSHIP/python/DecaySelection.conf"), verbose=True)

◆ configure()

None pythia8_conf.configure (   P8gen,
  mass,
  production_couplings,
  decay_couplings,
  process_selection,
bool   deepCopy = False,
bool   debug = True 
)
This function configures a HNLPythia8Generator instance for SHiP usage.

Definition at line 168 of file pythia8_conf.py.

170) -> None:
171 """
172 This function configures a HNLPythia8Generator instance for SHiP usage.
173 """
174
175 if process_selection is True: # For backward compatibility
176 process_selection = "inclusive"
177
178 # Wrap the Pythia8 object into a class logging all of its method calls
179 _exit_stack = contextlib.ExitStack()
180 if debug:
181 pythia_log = _exit_stack.enter_context(open("pythia8_conf.txt", "w")) # noqa: SIM115
182 P8gen = MethodLogger(P8gen, sink=pythia_log)
183
184 fairship_root = os.environ["FAIRSHIP"]
185 histograms = make_interpolators(fairship_root + "/shipgen/branchingratios.dat")
186 P8gen.UseRandom3() # TRandom1 or TRandom3 ?
187 P8gen.SetMom(400) # beam momentum in GeV
188 if deepCopy:
189 P8gen.UseDeepCopy()
190 ROOT.TDatabasePDG.Instance()
191 P8gen.SetParameters("Next:numberCount = 0")
192 # let strange particle decay in Geant4
193 make_particles_stable(P8gen, above_lifetime=1)
194
195 # Load particle & decay data
196 # ==========================
197
198 datafile = fairship_root + "/python/hnl_production.yaml"
199 with open(datafile) as f:
200 data = yaml.load(f, Loader=yaml.FullLoader)
201 all_channels = data["channels"]
202
203 # Inclusive
204 # =========
205
206 if process_selection == "inclusive":
207 setup_pythia_inclusive(P8gen)
208
209 # Charm decays only (with secondary production from tau)
210 # ======================================================
211
212 if process_selection == "c":
213 selection = data["selections"]["c"]
214 for cmd in selection["parameters"]:
215 P8gen.SetParameters(cmd)
216 add_hnl(P8gen, mass, decay_couplings)
217
218 # Add new charmed particles
219 # -------------------------
220
221 # Select all charmed particles (+ tau lepton)
222 c_particles = selection["particles"]
223 tau_id = 15 # tau- Monte-Carlo ID
224 add_particles(P8gen, c_particles + [tau_id], data)
225
226 # Add HNL production channels from charmed particles
227 # --------------------------------------------------
228
229 # Find charm and tau decays to HNLs
230 c_channels = [ch for ch in all_channels if ch["id"] in c_particles]
231 tau_channels = [ch for ch in all_channels if ch["id"] == tau_id]
232 # Standard model process: tau+ production from D_s+ decay
233 ds_id = 431 # D_s+ Monte-Carlo ID
234 ds_tau_br = 0.0548 # SM branching ratio Br(D_s+ -> tau+ nu_tau) (source: PDG 2018)
235
236 # Compute the branching ratio scaling factor, taking into account
237 # secondary production from tau
238 # Decay chains are encoded as follows:
239 # [(top level id A, [br A -> B, br B -> C, ...]), ...]
240
241 # Most charm particles directly decay to HNLs
242 primary_decays = [(ch["id"], [get_br(histograms, ch, mass, production_couplings)]) for ch in c_channels]
243 # The D_s+ can indirectly produce a HNL by first producing a tau+
244 secondary_decays = [
245 (ds_id, [ds_tau_br, get_br(histograms, ch, mass, production_couplings)]) for ch in tau_channels
246 ]
247 all_decays = primary_decays + secondary_decays
248
249 # Compute maximum total branching ratio (to rescale all BRs)
250 max_total_br = compute_max_total_br(all_decays)
251 exit_if_zero_br(max_total_br, process_selection, mass)
252 print_scale_factor(1 / max_total_br)
253
254 # Add charm decays
255 for ch in c_channels:
256 add_channel(P8gen, ch, histograms, mass, production_couplings, 1 / max_total_br)
257 # Add tau production from D_s+
258 # We can freely rescale Br(Ds -> tau) and Br(tau -> N X...) as long as
259 # Br(Ds -> tau -> N X...) remains the same.
260 # Here, we set Br(tau -> N) to unity to make event generation more efficient.
261 # The implicit assumption here is that we will disregard the tau during the analysis.
262 total_tau_br = sum(branching_ratios[1] for (_, branching_ratios) in secondary_decays)
263 assert ds_tau_br * total_tau_br <= max_total_br + 1e-12
264 P8gen.SetParameters(
265 f"431:addChannel 1 {ds_tau_br * total_tau_br / max_total_br:.12} 0 -15 16"
266 )
267 # Add secondary HNL production from tau
268 for ch in tau_channels:
269 # Rescale branching ratios only if some are non-zero. Otherwise leave them at zero.
270 add_tau_channel(P8gen, ch, histograms, mass, production_couplings, 1 / (total_tau_br or 1))
271
272 # Add dummy channels in place of SM processes
273 fill_missing_channels(P8gen, max_total_br, all_decays)
274
275 # List channels to confirm that Pythia has been properly set up
276 P8gen.List(9900015)
277
278 # B/Bc decays only
279 # ================
280
281 if process_selection in ["b", "bc"]:
282 selection = data["selections"][process_selection]
283 for cmd in selection["parameters"]:
284 P8gen.SetParameters(cmd)
285 add_hnl(P8gen, mass, decay_couplings)
286
287 # Add particles
288 particles = selection["particles"]
289 add_particles(P8gen, particles, data)
290
291 # Find all decay channels
292 channels = [ch for ch in all_channels if ch["id"] in particles]
293 decays = [(ch["id"], [get_br(histograms, ch, mass, production_couplings)]) for ch in channels]
294
295 # Compute scaling factor
296 max_total_br = compute_max_total_br(decays)
297 exit_if_zero_br(max_total_br, process_selection, mass)
298 print_scale_factor(1 / max_total_br)
299
300 # Add beauty decays
301 for ch in channels:
302 add_channel(P8gen, ch, histograms, mass, production_couplings, 1 / max_total_br)
303
304 # Add dummy channels in place of SM processes
305 fill_missing_channels(P8gen, max_total_br, decays)
306
307 P8gen.List(9900015)
308
309 _exit_stack.close()
310
311
int open(const char *, int)
Opens a file descriptor.

◆ configurerpvsusy()

None pythia8_conf.configurerpvsusy (   P8gen,
  mass,
  couplings,
  sfermionmass,
  benchmark,
  inclusive,
bool   deepCopy = False,
bool   debug = True 
)

Definition at line 32 of file pythia8_conf.py.

34) -> None:
35 # configure pythia8 for Ship usage
36 _exit_stack = contextlib.ExitStack()
37 if debug:
38 pythia_log = _exit_stack.enter_context(open("pythia8_conf.txt", "w")) # noqa: SIM115
39 P8gen = MethodLogger(P8gen, sink=pythia_log)
40 h = make_interpolators(os.path.expandvars(f"$FAIRSHIP/shipgen/branchingratiosrpvsusybench{benchmark}.dat"))
41 P8gen.UseRandom3()
42 P8gen.SetMom(400) # beam momentum in GeV
43 if deepCopy:
44 P8gen.UseDeepCopy()
45 ROOT.TDatabasePDG.Instance()
46 # let strange particle decay in Geant4
47 make_particles_stable(P8gen, above_lifetime=1)
48
49 if inclusive == "True":
50 setup_pythia_inclusive(P8gen)
51
52 # generate RPV neutralino from inclusive charm hadrons
53 if inclusive == "c":
54 P8gen.SetParameters("HardQCD::hardccbar = on")
55 # add RPVSUSY
56 rpvsusy_instance = rpvsusy.RPVSUSY(mass, couplings, sfermionmass, benchmark, debug=True)
57 ctau = rpvsusy_instance.computeNLifetime(system="FairShip") * u.c_light * u.cm
58 print("RPVSUSY ctau ", ctau)
59 P8gen.SetParameters(f"9900015:new = N2 N2 2 0 0 {mass:.12} 0.0 0.0 0.0 {ctau / u.mm:.12} 0 1 0 1 0")
60 P8gen.SetParameters("9900015:isResonance = false")
61 P8gen.SetParameters("Next:numberCount = 0")
62 # Configuring decay modes...
63 rpvsusy_instance.AddChannelsToPythia(P8gen)
64
65 # Finish HNL setup...
66 P8gen.SetParameters("9900015:mayDecay = on")
67 P8gen.SetHNLId(9900015)
68 # also add to PDG
69 gamma = u.hbarc / float(ctau) # 197.3269631e-16 / float(ctau) # hbar*c = 197 MeV*fm = 197e-16 GeV*cm
70 addHNLtoROOT(pid=9900015, m=mass, g=gamma)
71 # 12 14 16 neutrinos replace with N2
72 charmhistograms = ["d_mu", "ds_mu"]
73 # no tau decay here to consider
74 maxsumBR = getmaxsumbrrpvsusy(h, charmhistograms, mass, couplings)
75 exit_if_zero_br(maxsumBR, inclusive, mass, particle="RPV neutralino")
76 gettotalbrrpvsusy(h, charmhistograms, mass, couplings)
77
78 # overwrite D_s+ decays
79 P8gen.SetParameters(
80 "431:new D_s+ D_s- 1 3 0 1.96849"
81 " 0.00000 0.00000 0.00000 1.49900e-01 0 1 0 1 0"
82 )
83 sumBR = 0.0
84 if getbr_rpvsusy(h, "ds_mu", mass, couplings[1]) > 0.0:
85 P8gen.SetParameters(
86 "431:addChannel 1 {:.12} 0 -13 9900015".format(
87 getbr_rpvsusy(h, "ds_mu", mass, couplings[1]) / maxsumBR
88 )
89 )
90 sumBR += float(getbr_rpvsusy(h, "ds_mu", mass, couplings[1]) / maxsumBR)
91 if sumBR < 1.0 and sumBR > 0.0:
92 P8gen.SetParameters(f"431:addChannel 1 {1.0 - sumBR:.12} 0 22 -11")
93
94 # overwrite D+ decays
95 P8gen.SetParameters(
96 "411:new D+ D- 1 3 0 1.86962 0.00000 0.00000 0.00000 3.11800e-01 0 1 0 1 0"
97 )
98 sumBR = 0.0
99 if getbr_rpvsusy(h, "d_mu", mass, couplings[1]) > 0.0:
100 P8gen.SetParameters(
101 "411:addChannel 1 {:.12} 0 -13 9900015".format(
102 getbr_rpvsusy(h, "d_mu", mass, couplings[1]) / maxsumBR
103 )
104 )
105 sumBR += float(getbr_rpvsusy(h, "d_mu", mass, couplings[1]) / maxsumBR)
106 if sumBR < 1.0 and sumBR > 0.0:
107 P8gen.SetParameters(f"411:addChannel 1 {1.0 - sumBR:.12} 0 22 -11")
108
109 P8gen.List(9900015)
110
111 if inclusive == "b":
112 P8gen.SetParameters("HardQCD::hardbbbar = on")
113 # add RPVSUSY
114 rpvsusy_instance = rpvsusy.RPVSUSY(mass, couplings, sfermionmass, benchmark, debug=True)
115 ctau = rpvsusy_instance.computeNLifetime(system="FairShip") * u.c_light * u.cm
116 P8gen.SetParameters(f"9900015:new = N2 N2 2 0 0 {mass:.12} 0.0 0.0 0.0 {ctau / u.mm:.12} 0 1 0 1 0")
117 P8gen.SetParameters("9900015:isResonance = false")
118 # Configuring decay modes...
119 rpvsusy_instance.AddChannelsToPythia(P8gen)
120 # Finish HNL setup...
121 P8gen.SetParameters("9900015:mayDecay = on")
122 P8gen.SetHNLId(9900015)
123 # also add to PDG
124 gamma = u.hbarc / float(ctau) # 197.3269631e-16 / float(ctau) # hbar*c = 197 MeV*fm = 197e-16 GeV*cm
125 addHNLtoROOT(pid=9900015, m=mass, g=gamma)
126 # 12 14 16 neutrinos replace with N2
127 beautyhistograms = ["b_mu", "b_tau", "b0_nu_mu", "b0_nu_tau"]
128 maxsumBR = getmaxsumbrrpvsusy(h, beautyhistograms, mass, couplings)
129 exit_if_zero_br(maxsumBR, inclusive, mass, particle="RPV neutralino")
130 gettotalbrrpvsusy(h, beautyhistograms, mass, couplings)
131
132 # overwrite B+ decays
133 P8gen.SetParameters(
134 "521:new B+ B- 1 3 0 5.27925"
135 " 0.00000 0.00000 0.00000 4.91100e-01 0 1 0 1 0"
136 )
137 sumBR = 0.0
138 if getbr_rpvsusy(h, "b_tau", mass, couplings[1]) > 0.0:
139 P8gen.SetParameters(
140 "521:addChannel 1 {:.12} 0 9900015 -15".format(
141 getbr_rpvsusy(h, "b_tau", mass, couplings[1]) / maxsumBR
142 )
143 )
144 sumBR += float(getbr_rpvsusy(h, "b_tau", mass, couplings[1]) / maxsumBR)
145 if sumBR < 1.0 and sumBR > 0.0:
146 P8gen.SetParameters(f"521:addChannel 1 {1.0 - sumBR:.12} 0 22 22")
147
148 # overwrite B0 decays
149 P8gen.SetParameters(
150 "511:new B0 Bbar0 1 0 0 5.27958"
151 " 0.00000 0.00000 0.00000 4.58700e-01 0 1 0 1 0"
152 )
153 sumBR = 0.0
154 if getbr_rpvsusy(h, "b0_nu_tau", mass, couplings[1]) > 0.0:
155 P8gen.SetParameters(
156 "511:addChannel 1 {:.12} 22 9900015 16".format(
157 getbr_rpvsusy(h, "b0_nu_tau", mass, couplings[1]) / maxsumBR
158 )
159 )
160 if sumBR < 1.0 and sumBR > 0.0:
161 P8gen.SetParameters(f"511:addChannel 1 {1.0 - sumBR:.12} 0 22 22")
162
163 P8gen.List(9900015)
164
165 _exit_stack.close()
166
167

◆ setup_pythia_inclusive()

None pythia8_conf.setup_pythia_inclusive (   P8gen)

Definition at line 331 of file pythia8_conf.py.

331def setup_pythia_inclusive(P8gen) -> None:
332 P8gen.SetParameters("SoftQCD:inelastic = on")
333 P8gen.SetParameters("PhotonCollision:gmgm2mumu = on")
334 P8gen.SetParameters("PromptPhoton:all = on")
335 P8gen.SetParameters("WeakBosonExchange:all = on")