SND@LHC Software
Loading...
Searching...
No Matches
pythia8_conf_utils.py
Go to the documentation of this file.
1from __future__ import print_function
2from __future__ import division
3import os
4import sys
5import re
6import six
7import numpy as np
8import scipy.interpolate
9import ROOT
10import shipunit as u
11
12def addHNLtoROOT(pid=9900015 ,m = 1.0, g=3.654203020370371E-21):
13 pdg = ROOT.TDatabasePDG.Instance()
14 pdg.AddParticle('N2','HNL', m, False, g, 0., 'N2', pid)
15
16def getbr_rpvsusy(h,histoname,mass,coupling):
17 if histoname in h:
18 normalized_br = h[histoname](mass)
19 br = normalized_br * coupling
20 else:
21 br = 0
22 return br
23
24def getmaxsumbrrpvsusy(h,histograms,mass,couplings):
25 #0 MeV< mass < 3.200 GeV
26 maxsumbr=0.0
27 sumbrs={}
28 for histoname in histograms:
29 item = histoname.split('_')
30 lepton = item[len(item)-1]
31 meson = item[0]
32 coupling=couplings[1]
33 try:
34 sumbrs[meson]+=getbr_rpvsusy(h,histoname,mass,coupling)
35 except:
36 sumbrs[meson]=getbr_rpvsusy(h,histoname,mass,coupling)
37 print(list(sumbrs.values()))
38 maxsumbr=max(sumbrs.values())
39 return maxsumbr
40
41def gettotalbrrpvsusy(h,histograms,mass,couplings):
42 totalbr=0.0
43 for histoname in histograms:
44 item = histoname.split('_')
45 coupling=couplings[1]
46 totalbr+=getbr_rpvsusy(h,histoname,mass,coupling)
47 return totalbr
48
49def make_particles_stable(P8gen, above_lifetime):
50 # FIXME: find time unit and add it to the docstring
51 """
52 Make the particles with a lifetime above the specified one stable, to allow
53 them to decay in Geant4 instead.
54 """
55 p8 = P8gen.getPythiaInstance()
56 n=1
57 while n!=0:
58 n = p8.particleData.nextId(n)
59 p = p8.particleData.particleDataEntryPtr(n)
60 if p.tau0() > above_lifetime:
61 command = "{}:mayDecay = false".format(n)
62 p8.readString(command)
63 print("Pythia8 configuration: Made {} stable for Pythia, should decay in Geant4".format(p.name()))
64
65def parse_histograms(filepath):
66 """
67 This function parses a file containing histograms of branching ratios.
68
69 It places them in a dictionary indexed by the decay string (e.g. 'd_K0_e'),
70 as a pair ([masses...], [branching ratios...]), where the mass is expressed
71 in GeV.
72 """
73 with open(filepath, 'r') as f:
74 lines = f.readlines()
75 # Define regular expressions matching (sub-)headers and data lines
76 th1f_exp = re.compile(r'^TH1F\|.+')
77 header_exp = re.compile(r'^TH1F\|(.+?)\|B(?:R|F)/U2(.+?)\|.+? mass \‍(GeV\‍)\|?')
78 subheader_exp = re.compile(r'^\s*?(\d+?),\s*(\d+?\.\d+?),\s*(\d+\.\d+)\s*$')
79 data_exp = re.compile(r'^\s*(\d+?)\s*,\s*(\d+\.\d+)\s*$')
80 # Locate beginning of each histogram
81 header_line_idx = [i for i in range(len(lines)) if th1f_exp.match(lines[i]) is not None]
82 # Iterate over histograms
83 histograms = {}
84 for offset in header_line_idx:
85 # Parse header
86 mh = header_exp.match(lines[offset])
87 if mh is None or len(mh.groups()) != 2:
88 raise ValueError("Malformed header encountered: {0}".format(lines[offset]))
89 decay_code = mh.group(1)
90 # Parse sub-header (min/max mass and number of points)
91 ms = subheader_exp.match(lines[offset+1])
92 if ms is None or len(ms.groups()) != 3:
93 raise ValueError("Malformed sub-header encountered: {0}".format(lines[offset+1]))
94 npoints = int(ms.group(1))
95 min_mass = float(ms.group(2))
96 max_mass = float(ms.group(3))
97 masses = np.linspace(min_mass, max_mass, npoints, endpoint=False)
98 branching_ratios = np.zeros(npoints)
99 # Now read the data lines (skipping the two header lines)
100 for line in lines[offset+2:offset+npoints+1]:
101 md = data_exp.match(line)
102 if md is None or len(md.groups()) != 2:
103 raise ValueError("Malformed data row encountered: {0}".format(line))
104 idx = int(md.group(1))
105 br = float(md.group(2))
106 branching_ratios[idx] = br
107 histograms[decay_code] = (masses, branching_ratios)
108 return histograms
109
110def make_interpolators(filepath, kind='linear'):
111 """
112 This function reads a file containing branching ratio histograms, and
113 returns a dictionary of interpolators of the branching ratios, indexed by
114 the decay string.
115 """
116 histogram_data = parse_histograms(filepath)
117 histograms = {}
118 for (hist_string, (masses, br)) in six.iteritems(histogram_data):
119 histograms[hist_string] = scipy.interpolate.interp1d(
120 masses, br, kind=kind, bounds_error=False, fill_value=0, assume_sorted=True)
121 return histograms
122
123def get_br(histograms, channel, mass, couplings):
124 """
125 Utility function used to reliably query the branching ratio for a given
126 channel at a given mass, taking into account the correct coupling.
127 """
128 hist = histograms[channel['decay']]
129 coupling = couplings[channel['coupling']]
130 normalized_br = hist(mass)
131 return normalized_br * coupling
132
133def add_particles(P8gen, particles, data):
134 """
135 Adds the corresponding particles to PYTHIA.
136
137 `particles` must be a list containing either the particles PDG IDs, or
138 their PYTHIA names. The commands needed to add the particles are queried
139 from `data`.
140
141 If the particle is not self-conjugate, the antiparticle is automatically
142 added by PYTHIA.
143 """
144 for particle_id in particles:
145 # Find particle in database (None: particle not found)
146 particle = next((p for p in data['particles']
147 if particle_id in [p['id'], p['name']]), None)
148 if particle is None:
149 raise ValueError("Could not find particle ID {0} in file {1}"
150 .format(particle, datafile))
151 # Add the particle
152 P8gen.SetParameters(particle['cmd'])
153
154def add_channel(P8gen, ch, histograms, mass, couplings, scale_factor):
155 "Add to PYTHIA a leptonic or semileptonic decay channel to HNL."
156 if 'idlepton' in ch:
157 br = get_br(histograms, ch, mass, couplings)
158 if br <= 0: # Ignore kinematically closed channels
159 return
160 if 'idhadron' in ch: # Semileptonic decay
161 P8gen.SetParameters('{}:addChannel 1 {:.17} 22 {} 9900015 {}'.format(ch['id'], br*scale_factor, ch['idlepton'], ch['idhadron']))
162 else: # Leptonic decay
163 P8gen.SetParameters('{}:addChannel 1 {:.17} 0 9900015 {}'.format(ch['id'], br*scale_factor, ch['idlepton']))
164 else: # Wrong decay
165 raise ValueError("Missing key 'idlepton' in channel {0}".format(ch))
166
167def add_tau_channel(P8gen, ch, histograms, mass, couplings, scale_factor):
168 "Add to PYTHIA a tau decay channel to HNL."
169 if 'idhadron' in ch:
170 br = get_br(histograms, ch, mass, couplings)
171 if br <= 0: # Ignore kinematically closed channels
172 return
173 if 'idlepton' in ch: # 3-body leptonic decay
174 P8gen.SetParameters('{}:addChannel 1 {:.16} 1531 9900015 {} {}'.format(ch['id'], br*scale_factor, ch['idlepton'], ch['idhadron']))
175 else: # 2-body semileptonic decay
176 P8gen.SetParameters('{}:addChannel 1 {:.16} 1521 9900015 {}'.format(ch['id'], br*scale_factor, ch['idhadron']))
177 else:
178 raise ValueError("Missing key 'idhadron' in channel {0}".format(ch))
179
180def fill_missing_channels(P8gen, max_total_br, decay_chains, epsilon=1e-6):
181 """
182 Add dummy channels for correct rejection sampling.
183
184 Even after rescaling the branching ratios, they do not sum up to unity
185 for most particles since we are ignoring SM processes.
186
187 This function adds a "filler" channel for each particle, in order to
188 preserve the ratios between different branching ratios.
189 """
190 top_level_particles = get_top_level_particles(decay_chains)
191 for particle in top_level_particles:
192 my_total_br = compute_total_br(particle, decay_chains)
193 remainder = 1 - my_total_br / max_total_br
194 assert(remainder > -epsilon)
195 assert(remainder < 1 + epsilon)
196 if remainder > epsilon:
197 add_dummy_channel(P8gen, particle, remainder)
198
199def add_dummy_channel(P8gen, particle, remainder):
200 """
201 Add a dummy channel to PYTHIA, with branching ratio equal to `remainder.`
202
203 The purpose of this function is to compensate for the absence of SM
204 channels, which are ignored when investigating rare processes. A dummy
205 decay channel is instead added to each particle in order to maintain the
206 correct ratios between the branching ratios of each particle to rare
207 processes. This is usually combined with a global reweighting of the
208 branching ratios.
209
210 In order to keep PYTHIA from complaining about charge conservation, a
211 suitable process is added which conserves electric charge.
212
213 All dummy channels can be identified by the presence of a photon among the
214 decay products.
215 """
216 pdg = P8gen.getPythiaInstance().particleData
217 charge = pdg.charge(particle)
218 if charge > 0:
219 P8gen.SetParameters('{}:addChannel 1 {:.16} 0 22 -11'.format(particle, remainder))
220 elif charge < 0:
221 P8gen.SetParameters('{}:addChannel 1 {:.16} 0 22 11'.format(particle, remainder))
222 else:
223 P8gen.SetParameters('{}:addChannel 1 {:.16} 0 22 22'.format(particle, remainder))
224
225def compute_max_total_br(decay_chains):
226 """
227 This function computes the maximum total branching ratio for all decay chains.
228
229 In order to make the event generation as efficient as possible when
230 studying very rare processes, it is necessary to rescale the branching
231 ratios, while enforcing the invariant that any total branching ratio must
232 remain lower that unity.
233
234 This is accomplished by computing, for each particle, the total branching
235 ratio to processes of interest, and then dividing all branching ratios by
236 the highest of those.
237
238 Note: the decay chain length must be at most 2.
239 """
240 # For each top-level charmed particle, sum BR over all its decay chains
241 top_level_particles = get_top_level_particles(decay_chains)
242 total_branching_ratios = [compute_total_br(particle, decay_chains)
243 for particle in top_level_particles]
244 # Find the maximum total branching ratio
245 return max(total_branching_ratios)
246
247def compute_total_br(particle, decay_chains):
248 """
249 Returns the total branching ratio to HNLs for a given particle.
250 """
251 return sum(np.prod(branching_ratios)
252 for (top, branching_ratios) in decay_chains
253 if top == particle)
254
255def get_top_level_particles(decay_chains):
256 """
257 Returns the set of particles which are at the top of a decay chain.
258 """
259 return set(top for (top, branching_ratios) in decay_chains)
260
261def exit_if_zero_br(max_total_br, selection, mass, particle='HNL'):
262 if max_total_br <= 0:
263 print("No phase space for {0} from {1} at this mass: {2}. Quitting."
264 .format(particle, selection, mass))
265 sys.exit()
266
267def print_scale_factor(scaling_factor):
268 "Prints the scale factor used to make event generation more efficient."
269 print("One simulated event per {0:.4g} meson decays".format(scaling_factor))
set(INCLUDE_DIRECTORIES ${SYSTEM_INCLUDE_DIRECTORIES} ${VMC_INCLUDE_DIRS} ${CMAKE_SOURCE_DIR}/shipdata ${CMAKE_SOURCE_DIR}/shipLHC ${CMAKE_SOURCE_DIR}/analysis/cuts ${CMAKE_SOURCE_DIR}/analysis/tools ${FMT_INCLUDE_DIR}) include_directories($
Definition CMakeLists.txt:1
compute_max_total_br(decay_chains)
add_particles(P8gen, particles, data)
make_interpolators(filepath, kind='linear')
fill_missing_channels(P8gen, max_total_br, decay_chains, epsilon=1e-6)
add_dummy_channel(P8gen, particle, remainder)
compute_total_br(particle, decay_chains)
gettotalbrrpvsusy(h, histograms, mass, couplings)
print_scale_factor(scaling_factor)
add_channel(P8gen, ch, histograms, mass, couplings, scale_factor)
getmaxsumbrrpvsusy(h, histograms, mass, couplings)
get_br(histograms, channel, mass, couplings)
get_top_level_particles(decay_chains)
make_particles_stable(P8gen, above_lifetime)
add_tau_channel(P8gen, ch, histograms, mass, couplings, scale_factor)
getbr_rpvsusy(h, histoname, mass, coupling)
exit_if_zero_br(max_total_br, selection, mass, particle='HNL')