Belle II Software prerelease-11-00-00a
generators.py
1
8
9'''
10This module contains convenience functions to setup most commonly used physics
11generators correctly with their default settings. More information can be found
12in `BELLE2-NOTE-PH-2015-006`_
13
14.. _BELLE2-NOTE-PH-2015-006: https://docs.belle2.org/record/282
15'''
16
17import basf2 as b2
18import json
19import pdg
20
21
22def get_default_decayfile():
23 """Return the default DECAY.dec for Belle2"""
24 return b2.find_file("decfiles/dec/DECAY_BELLE2.DEC")
25
26
27def add_generator_preselection(
28 path,
29 emptypath,
30 nChargedMin=0,
31 nChargedMax=999,
32 MinChargedP=-1.0,
33 MinChargedPt=-1.0,
34 MinChargedTheta=0.0,
35 MaxChargedTheta=180.0,
36 nPhotonMin=0,
37 nPhotonMax=999,
38 MinPhotonEnergy=-1,
39 MinPhotonTheta=0.0,
40 MaxPhotonTheta=180.0,
41 applyInCMS=False,
42 stableParticles=False):
43 """
44 Adds generator preselection.
45 Should be added to the path after the generator.add_abc_generator but before simulation.add_simulation modules
46 It uses all particles from the event generator (i.e. primary, non-virtual, non-initial particles).
47 It checks if the required conditions are fulfilled.
48 If not, the events are given to the emptypath.
49 The main use case is a reduction of simulation time.
50 Note that you have to multiply the generated cross section by the retention fraction of the preselection.
51
52 Parameters:
53 path (basf2.Path): path where the generator should be added
54 emptypath (basf2.Path): path where the skipped events are given to
55 nChargedMin (int): minimum number of charged particles
56 nChargedMax (int): maximum number of charged particles
57 MinChargedP (float): minimum charged momentum [GeV]
58 MinChargedPt (float): minimum charged transverse momentum (pt) [GeV]
59 MinChargedTheta (float): minimum polar angle of charged particle [deg]
60 MaxChargedTheta (float): maximum polar angle of charged particle [deg]
61 nPhotonMin (int): minimum number of photons
62 nPhotonMax (int): maximum number of photons
63 MinPhotonEnergy (float): minimum photon energy [GeV]
64 MinPhotonTheta (float): minimum polar angle of photon [deg]
65 MaxPhotonTheta (float): maximum polar angle of photon [deg]
66 applyInCMS (bool): if true apply the P,Pt,theta, and energy cuts in the center of mass frame
67 stableParticles (bool): if true apply the selection criteria for stable particles in the generator
68 """
69
70 generatorpreselection = path.add_module('GeneratorPreselection',
71 nChargedMin=nChargedMin,
72 nChargedMax=nChargedMax,
73 MinChargedP=MinChargedP,
74 MinChargedPt=MinChargedPt,
75 MinChargedTheta=MinChargedTheta,
76 MaxChargedTheta=MaxChargedTheta,
77 nPhotonMin=nPhotonMin,
78 nPhotonMax=nPhotonMax,
79 MinPhotonEnergy=MinPhotonEnergy,
80 MinPhotonTheta=MinPhotonTheta,
81 MaxPhotonTheta=MaxPhotonTheta,
82 applyInCMS=applyInCMS,
83 stableParticles=stableParticles
84 )
85
86 # empty path for unwanted events
87 generatorpreselection.if_value('<11', emptypath)
88
89
90def add_aafh_generator(
91 path,
92 finalstate='',
93 preselection=False,
94 enableTauDecays=True,
95 minmass=0.5,
96 subweights=[],
97 maxsubweight=1,
98 maxfinalweight=3.0,
99 eventType=''):
100 """
101 Add the default two photon generator for four fermion final states
102
103 Parameters:
104 path (basf2.Path): path where the generator should be added
105 finalstate (str): either "e+e-e+e-", "e+e-mu+mu-", "e+e-tau+tau-", "mu+mu-mu+mu-" or "mu+mu-tau+tau-"
106 preselection (bool): if True, select events with at least one medium pt particle in the CDC acceptance
107 enableTauDecays (bool): if True, allow tau leptons to decay (using EvtGen)
108 minmass (float): minimum invariant mass
109 subweights (list(float)): list of four or eight values (first four are interpreted as WAP, rest as WBP)
110 which specify the relative weights for each of the four sub generators
111 maxsubweight (float): maximum expected subgenerator weight for rejection scheme
112 maxfinalweight (float): maximum expected final weight for rejection scheme
113 eventType (str) : event type information
114 """
115
116 if finalstate == 'e+e-e+e-':
117 aafh_mode = 5
118 if not subweights: # default subweights are for minmass=0.5
119 aafh_subgeneratorWeights = [1.0, 7.986e+01, 5.798e+04, 3.898e+05, 1.0, 1.664e+00, 2.812e+00, 7.321e-01]
120 else:
121 aafh_subgeneratorWeights = subweights
122 if abs(minmass - 0.5) > 0.01 and not subweights:
123 b2.B2WARNING("add_aafh_generator: non default invariant mass cut without updated subweights requested!")
124 elif finalstate == 'e+e-mu+mu-':
125 aafh_mode = 3
126 if not subweights: # default subweights are for minmass=0.5
127 aafh_subgeneratorWeights = [1.000e+00, 1.520e+01, 3.106e+03, 6.374e+03, 1.000e+00, 1.778e+00, 6.075e+00, 6.512e+00]
128 else:
129 aafh_subgeneratorWeights = subweights
130 if abs(minmass - 0.5) > 0.01 and not subweights:
131 b2.B2WARNING("add_aafh_generator: non default invariant mass cut without updated subweights requested!")
132 elif finalstate == 'e+e-tau+tau-':
133 aafh_mode = 4
134 minmass = 0
135 if not subweights:
136 aafh_subgeneratorWeights = [1.000e+00, 2.214e+00, 1.202e+01, 1.536e+01, 1.000e+00, 1.664e+00, 1.680e+01, 6.934e+00]
137 else:
138 aafh_subgeneratorWeights = subweights
139 if preselection:
140 b2.B2WARNING(
141 f"You requested a generator preselection for the final state {finalstate}: "
142 "please consider to remove it, since the cross section is small.")
143 elif finalstate == 'mu+mu-mu+mu-':
144 aafh_mode = 2
145 minmass = 0
146 maxsubweight = 3
147 if not subweights:
148 aafh_subgeneratorWeights = [0.000e+00, 0.000e+00, 1.000e+00, 3.726e+00, 1.000e+00, 1.778e+00, 1.000e+00, 1.094e+00]
149 else:
150 aafh_subgeneratorWeights = subweights
151 if preselection:
152 b2.B2WARNING(
153 f"You requested a generator preselection for the final state {finalstate}: "
154 "please consider to remove it, since the cross section is small.")
155 elif finalstate == 'mu+mu-tau+tau-':
156 aafh_mode = 1
157 minmass = 0
158 maxsubweight = 3
159 if not subweights:
160 aafh_subgeneratorWeights = [0.000e+00, 0.000e+00, 1.000e+00, 1.715e+00, 1.000e+00, 1.778e+00, 1.000e+00, 6.257e-01]
161 else:
162 aafh_subgeneratorWeights = subweights
163 if preselection:
164 b2.B2WARNING(
165 f"You requested a generator preselection for the final state {finalstate}: "
166 "please consider to remove it, since the cross section is small.")
167 elif finalstate == 'tau+tau-tau+tau-':
168 b2.B2FATAL(f"AAFH is not able to generate the {finalstate} final state. Please use KoralW instead.")
169 else:
170 b2.B2FATAL(f"add_aafh_generator final state not supported: {finalstate}")
171
172 aafh_maxSubgeneratorWeight = maxsubweight
173 aafh_maxFinalWeight = maxfinalweight
174
175 path.add_module(
176 'AafhInput',
177 mode=aafh_mode,
178 rejection=2,
179 maxSubgeneratorWeight=aafh_maxSubgeneratorWeight,
180 maxFinalWeight=aafh_maxFinalWeight,
181 subgeneratorWeights=aafh_subgeneratorWeights,
182 suppressionLimits=[1e100] * 4,
183 minMass=minmass,
184 eventType=eventType,
185 )
186
187 if preselection:
188 generator_emptypath = b2.create_path()
189 add_generator_preselection(
190 path=path,
191 emptypath=generator_emptypath,
192 nChargedMin=1,
193 MinChargedPt=0.1,
194 MinChargedTheta=17.0,
195 MaxChargedTheta=150.0)
196
197 if 'tau+tau-' in finalstate:
198 if enableTauDecays:
199 path.add_module('EvtGenDecay')
200 else:
201 b2.B2WARNING("The tau decays will not be generated.")
202
203
204def add_kkmc_generator(path, finalstate='', signalconfigfile='', useTauolaBelle=False, tauinputfile='', eventType=''):
205 """
206 Add the default muon pair and tau pair generator KKMC.
207 For tau decays, TauolaBelle and TauolaBelle2 are available.
208 Signal events can be produced setting a configuration file. Please notice that the configuration files for
209 TauolaBelle and TauolaBelle2 has a very different structure (see the examples below generators/examples).
210
211 Parameters:
212 path (basf2.Path): path where the generator should be added
213 finalstate(str): either "mu-mu+" or "tau-tau+"
214 signalconfigfile(str): File with configuration of the signal event to generate. It doesn't affect mu-mu+ decays.
215 useTauolaBelle(bool): If true, tau decay is driven by TauolaBelle. Otherwise TauolaBelle2 is used.
216 It doesn't affect mu-mu+ decays.
217 tauinputfile(str) : File to override KK2f_defaults. Only [sometimes] needed when tau decay is driven by TauolaBelle.
218 eventType (str) : event type information
219 """
220
221
222 kkmc_inputfile = b2.find_file('data/generators/kkmc/tauola_bbb.input.dat')
223
224
225 kkmc_logfile = 'kkmc_tautau.txt'
226
227
228 kkmc_config = b2.find_file('data/generators/kkmc/KK2f_defaults.dat')
229
230
231 kkmc_tauconfigfile = ''
232
233 if finalstate == 'tau+tau-':
234 b2.B2WARNING("add_kkmc_generator: please set finalstate as 'tau-tau+'. 'tau+tau-' will be deprecated in the future"
235 " for consistency in the configuration files.")
236 finalstate = 'tau-tau+'
237 if finalstate == 'mu+mu-':
238 b2.B2WARNING("add_kkmc_generator: please set finalstate as 'mu-mu+'. 'mu+mu-' will be deprecated in the future for"
239 " consistency in the configuration files.")
240 finalstate = 'mu-mu+'
241
242 if finalstate == 'tau-tau+':
243 if useTauolaBelle:
244 b2.B2INFO("Generating tau pair events with TauolaBelle")
245
246 kkmc_inputfile = b2.find_file('data/generators/kkmc/tau.input.dat')
247 kkmc_tauconfigfile = b2.find_file('data/generators/kkmc/tau_decaytable.dat')
248
249 if not signalconfigfile == '':
250 b2.B2INFO(f"Using config file defined by user: {signalconfigfile}")
251 if useTauolaBelle:
252 kkmc_tauconfigfile = b2.find_file(signalconfigfile)
253 else:
254 kkmc_inputfile = b2.find_file(signalconfigfile)
255
256 if not tauinputfile == '':
257 kkmc_inputfile = b2.find_file(tauinputfile)
258
259 elif finalstate == 'mu-mu+':
260 kkmc_inputfile = b2.find_file('data/generators/kkmc/mu.input.dat')
261 kkmc_logfile = 'kkmc_mumu.txt'
262
263 else:
264 b2.B2FATAL("add_kkmc_generator final state not supported: {}".format(finalstate))
265
266 # use KKMC to generate lepton pairs
267 path.add_module(
268 'KKGenInput',
269 tauinputFile=kkmc_inputfile,
270 KKdefaultFile=kkmc_config,
271 taudecaytableFile=kkmc_tauconfigfile,
272 kkmcoutputfilename=kkmc_logfile,
273 eventType=eventType,
274 )
275
276
277def add_evtgen_generator(path, finalstate='', signaldecfile=None, coherentMixing=True, parentParticle='Upsilon(4S)', eventType=''):
278 """
279 Add EvtGen for mixed and charged BB
280
281 Parameters:
282 path (basf2.Path): path where the generator should be added
283 finalstate (str): Either "charged" for generation of generic B+/B-, "mixed" for generic B0/anti-B0, or "signal" for
284 generation of user-defined signal mode
285 signaldecfile (str): path to decfile which defines the signal decay to be generated
286 (only needed if ``finalstate`` set to "signal")
287 coherentMixing: Either True or False. Switches on or off the coherent decay of the B0-B0bar pair.
288 It should always be True, unless you are generating Y(5,6S) -> BBar. In the latter case,
289 setting it False solves the internal limitation of Evtgen that allows to make a
290 coherent decay only starting from the Y(4S).
291 parentParticle (str): initial state (used only if it is not Upsilon(4S).
292 eventType (str) : event type information
293 """
294 evtgen_userdecfile = b2.find_file('data/generators/evtgen/charged.dec')
295
296 if parentParticle != 'Upsilon(3S)' and parentParticle != 'Upsilon(4S)'\
297 and parentParticle != 'Upsilon(5S)' and parentParticle != 'Upsilon(6S)':
298 b2.B2FATAL("add_evtgen_generator initial state not supported: {}".format(parentParticle))
299
300 if finalstate == 'charged':
301 pass
302 elif finalstate == 'mixed':
303 evtgen_userdecfile = b2.find_file('data/generators/evtgen/mixed.dec')
304 elif finalstate == 'signal':
305 evtgen_userdecfile = signaldecfile
306 else:
307 b2.B2FATAL("add_evtgen_generator final state not supported: {}".format(finalstate))
308
309 if signaldecfile and finalstate in ['charged', 'mixed']:
310 b2.B2WARNING("ignoring decfile: {}".format(signaldecfile))
311
312 # use EvtGen
313 if parentParticle == 'Upsilon(3S)':
314 if finalstate != 'signal':
315 b2.B2FATAL("add_evtgen_generator initial state {} is supported only with 'signal' final state".format(parentParticle))
316 if coherentMixing:
317 coherentMixing = False
318 b2.B2WARNING("add_evtgen_generator initial state {} has no BB mixing, now switching coherentMixing OFF"
319 .format(parentParticle))
320
321 if parentParticle == 'Upsilon(5S)':
322 if finalstate != 'signal':
323 b2.B2FATAL("add_evtgen_generator initial state {} is supported only with 'signal' final state".format(parentParticle))
324 if coherentMixing:
325 coherentMixing = False
326 b2.B2WARNING(
327 "add_evtgen_generator initial state {} is supported only with false coherentMixing, now switching it OFF"
328 .format(parentParticle))
329 pdg.load(b2.find_file('decfiles/dec/Y5S.pdl'))
330
331 if parentParticle == 'Upsilon(6S)':
332 if finalstate != 'signal':
333 b2.B2FATAL("add_evtgen_generator initial state {} is supported only with 'signal' final state".format(parentParticle))
334 if coherentMixing:
335 coherentMixing = False
336 b2.B2WARNING(
337 "add_evtgen_generator initial state {} is supported only with false coherentMixing, now switching it OFF"
338 .format(parentParticle))
339 pdg.load(b2.find_file('decfiles/dec/Y6S.pdl'))
340
341 path.add_module(
342 'EvtGenInput',
343 userDECFile=evtgen_userdecfile,
344 CoherentMixing=coherentMixing,
345 ParentParticle=parentParticle,
346 eventType=eventType,
347 )
348
349
350def add_continuum_generator(path, finalstate, userdecfile='', *, skip_on_failure=True, eventType=''):
351 """
352 Add the default continuum generators KKMC + PYTHIA including their default decfiles and PYTHIA settings
353
354 See Also:
355 `add_inclusive_continuum_generator()` to add continuum generation with preselected particles
356
357 Parameters:
358 path (basf2.Path): path where the generator should be added
359 finalstate (str): uubar, ddbar, ssbar, ccbar
360 userdecfile (str): EvtGen decfile used for particle decays
361 skip_on_failure (bool): If True stop event processing right after
362 fragmentation fails. Otherwise continue normally
363 eventType (str) : event type information
364 """
365
366
367 kkmc_inputfile = b2.find_file('data/generators/kkmc/uubar_nohadronization.input.dat')
368
369
370 kkmc_logfile = 'kkmc_uubar.txt'
371
372
373 pythia_config = b2.find_file('data/generators/modules/fragmentation/pythia_belle2.dat')
374
375
376 decay_user = b2.find_file('data/generators/modules/fragmentation/dec_belle2_qqbar.dec')
377 if userdecfile == '':
378 pass
379 else:
380 b2.B2INFO('Replacing default user decfile: {}'.format(userdecfile))
381 decay_user = userdecfile
382
383
384 kkmc_config = b2.find_file('data/generators/kkmc/KK2f_defaults.dat')
385
386
387 decay_file = get_default_decayfile()
388
389 if finalstate == 'uubar':
390 pass
391 elif finalstate == 'ddbar':
392 kkmc_inputfile = b2.find_file('data/generators/kkmc/ddbar_nohadronization.input.dat')
393 kkmc_logfile = 'kkmc_ddbar.txt'
394 elif finalstate == 'ssbar':
395 kkmc_inputfile = b2.find_file('data/generators/kkmc/ssbar_nohadronization.input.dat')
396 kkmc_logfile = 'kkmc_ssbar.txt'
397 elif finalstate == 'ccbar':
398 kkmc_inputfile = b2.find_file('data/generators/kkmc/ccbar_nohadronization.input.dat')
399 pythia_config = b2.find_file('data/generators/modules/fragmentation/pythia_belle2_charm.dat')
400 kkmc_logfile = 'kkmc_ccbar.txt'
401 else:
402 b2.B2FATAL("add_continuum_generator final state not supported: {}".format(finalstate))
403
404 # use KKMC to generate qqbar events (no fragmentation at this stage)
405 path.add_module(
406 'KKGenInput',
407 tauinputFile=kkmc_inputfile,
408 KKdefaultFile=kkmc_config,
409 taudecaytableFile='',
410 kkmcoutputfilename=kkmc_logfile,
411 eventType=eventType,
412 )
413
414 # add the fragmentation module to fragment the generated quarks into hadrons
415 # using PYTHIA8
416 fragmentation = path.add_module(
417 'Fragmentation',
418 ParameterFile=pythia_config,
419 ListPYTHIAEvent=0,
420 UseEvtGen=1,
421 DecFile=decay_file,
422 UserDecFile=decay_user,
423 )
424
425 if skip_on_failure:
426 # branch to an empty path if PYTHIA failed, this will change the number of events
427 # but the file meta data will contain the total number of generated events
428 generator_emptypath = b2.create_path()
429 fragmentation.if_value('<1', generator_emptypath)
430
431
432def add_inclusive_continuum_generator(path, finalstate, particles, userdecfile='',
433 *, include_conjugates=True, max_iterations=100000, eventType=''):
434 """
435 Add continuum generation but require at least one of the given particles be
436 present in the event.
437
438 For example to only generate ccbar events which contain a "D*+" or an
439 electron one could would use
440
441 >>> add_inclusive_continuum_generator(path, "ccbar", ["D*+", 11])
442
443 If you are unsure how the particles are named in Belle II please have a look
444 at the ``b2help-particles`` executable or the `pdg` python module.
445
446 See Also:
447 `add_continuum_generator()` to add continuum generation without preselection
448
449 Parameters:
450 path (basf2.Path): path to which the generator should be added
451 finalstate (str): uubar, ddbar, ssbar, ccbar
452 particles (list): A list of particle names or pdg codes. An event is
453 only accepted if at lease one of those particles appears in the event.
454 userdecfile (str): EvtGen decfile used for particle decays
455 include_conjugates (bool): If True (default) accept the event also if a
456 charge conjugate of the given particles is found
457 max_iterations (int): maximum tries per event to generate the requested
458 particle. If exceeded processing will be stopped with a
459 `FATAL <LogLevel.FATAL>` error so for rare particles one might need a
460 larger number.
461 eventType (str) : event type information
462 """
463 loop_path = b2.create_path()
464 # we might run this more than once so make sure we remove any particles
465 # before generating new ones
466 loop_path.add_module("PruneDataStore", keepMatchedEntries=False, matchEntries=["MCParticles", "EventExtraInfo"])
467 # add the generator but make sure it doesn't stop processing on
468 # fragmentation failure as is this currently not supported by do_while
469 add_continuum_generator(loop_path, finalstate, userdecfile, skip_on_failure=False, eventType=eventType)
470 # check for the particles we want
471 loop_path.add_module("InclusiveParticleChecker", particles=particles, includeConjugates=include_conjugates)
472 # Done, add this to the path and iterate it until we found our particle
473 path.do_while(loop_path, max_iterations=max_iterations)
474
475
476def add_bhwide_generator(path, minangle=0.5):
477 """
478 Add the high precision QED generator BHWIDE to the path. Settings are the default L1/HLT study settings
479 with a cross section of about 124000 nb (!)
480
481 Parameters:
482 path (basf2.Path): path where the generator should be added
483 minangle (float): minimum angle of the outgoing electron/positron in the CMS in degrees
484 """
485
486 if minangle < 0.0 or minangle > 180.0:
487 b2.B2FATAL("add_bhwide_generator minimum angle too small (<0.0) or too large (>180): {}".format(minangle))
488
489 bhwide = path.add_module("BHWideInput")
490 bhwide.param('ScatteringAngleRangePositron', [minangle, 180.0 - minangle])
491 bhwide.param('ScatteringAngleRangeElectron', [minangle, 180.0 - minangle])
492 bhwide.param('MaxAcollinearity', 180.0)
493 bhwide.param('MinEnergy', 0.10)
494 bhwide.param('VacuumPolarization', 'burkhardt')
495 bhwide.param('WeakCorrections', True)
496 bhwide.param('WtMax', 3.0)
497
498
499def add_babayaganlo_generator(
500 path,
501 finalstate='',
502 minenergy=0.01,
503 minangle=10.0,
504 fmax=-1.0,
505 generateInECLAcceptance=False,
506 eventType=''):
507 '''
508 Add the high precision QED generator BabaYaga@NLO to the path.
509
510 Parameters:
511 path (basf2.Path): path where the generator should be added.
512 finalstate (str): ee (e+e-), mm(mu+mu-) or gg (gammagamma).
513 minenergy (float): minimum particle (leptons for 'ee' or 'mm', photons for 'gg') energy in GeV.
514 minangle (float): angular range from minangle to 180-minangle for primary particles (in degrees).
515 fmax (float): maximum of differential cross section weight. This parameter should be set only by experts.
516 generateInECLAcceptance (bool): if True, the GeneratorPreselection module is used to select only events
517 with both the primary particles within the ECL acceptance.
518 eventType (str) : event type information
519 '''
520
521 if finalstate != 'ee' and finalstate != 'gg' and finalstate != 'mm':
522 b2.B2FATAL(f'add_babayaganlo_generator final state not supported: {finalstate}')
523
524 if not (fmax == -1.0):
525 b2.B2WARNING(f'The BabayagaNLOInput parameter "FMax" will be set to {fmax} instead to the default value (-1.0). '
526 'Please do not do this, unless you are extremely sure about this choice.')
527
528 path.add_module(
529 'BabayagaNLOInput',
530 eventType=eventType,
531 FinalState=finalstate,
532 ScatteringAngleRange=[minangle, 180.0 - minangle],
533 MinEnergy=minenergy,
534 FMax=fmax
535 )
536
537 if generateInECLAcceptance:
538 b2.B2INFO(f'The final state {finalstate} is preselected requiring both primary particles within the ECL acceptance.')
539 emptypath = b2.Path()
540 add_generator_preselection(path=path,
541 emptypath=emptypath,
542 applyInCMS=False)
543
544 if finalstate == 'ee' or finalstate == 'mm':
545 b2.set_module_parameters(path=path,
546 name='GeneratorPreselection',
547 nChargedMin=2,
548 MinChargedTheta=12.4,
549 MaxChargedTheta=155.1,)
550
551 elif finalstate == 'gg':
552 b2.set_module_parameters(path=path,
553 name='GeneratorPreselection',
554 nPhotonMin=2,
555 MinPhotonTheta=12.4,
556 MaxPhotonTheta=155.1)
557
558
559def add_phokhara_generator(path, finalstate='', eventType=''):
560 """
561 Add the high precision QED generator PHOKHARA to the path. Almost full
562 acceptance settings for photons and hadrons/muons.
563
564 Parameters:
565 path (basf2.Path): path where the generator should be added
566 finalstate (str): One of the following final states: "mu+mu-", "pi+pi-", "pi+pi-pi0", "pi+pi-pi+pi-" (or "2(pi+pi-)"),
567 "pi+pi-pi0pi0" or ("pi+pi-2pi0"), "pi+pi-eta", "K+K-", "K0K0bar", "ppbar", "n0n0bar" or "Lambda0Lambda0bar"
568 eventType (str) : event type information
569 """
570
571 if finalstate == 'mu+mu-':
572 path.add_module(
573 'PhokharaInput',
574 FinalState=0, # mu+mu-
575 LO=0, NLO=1, QED=0, # use full two loop corrections
576 MinInvMassHadrons=0.,
577 eventType=eventType,
578 ).set_name('PHOKHARA_mumuISR')
579
580 elif finalstate == 'pi+pi-':
581 path.add_module(
582 'PhokharaInput',
583 FinalState=1, # pi+pi-
584 LO=0, NLO=1, QED=0, # use full two loop corrections
585 MinInvMassHadrons=0.,
586 eventType=eventType,
587 ).set_name('PHOKHARA_pipiISR')
588
589 elif finalstate == 'pi+pi-pi0':
590 path.add_module(
591 'PhokharaInput',
592 FinalState=8, # pi+pi-pi0
593 LO=0, NLO=1, QED=0, # use full two loop corrections
594 MinInvMassHadrons=0.,
595 eventType=eventType,
596 ).set_name('PHOKHARA_pipipi0ISR')
597
598 elif finalstate == 'pi+pi-pi+pi-' or finalstate == '2(pi+pi-)':
599 path.add_module(
600 'PhokharaInput',
601 FinalState=3, # 2(pi+pi-)
602 LO=0, NLO=1, QED=0, # use full two loop corrections
603 MinInvMassHadrons=0.,
604 eventType=eventType,
605 ).set_name('PHOKHARA_pipipipiISR')
606
607 elif finalstate == 'pi+pi-pi0pi0' or finalstate == 'pi+pi-2pi0':
608 path.add_module(
609 'PhokharaInput',
610 FinalState=2, # pi+pi-2pi0
611 LO=0, NLO=1, QED=0, # use full two loop corrections
612 MinInvMassHadrons=0.,
613 eventType=eventType,
614 ).set_name('PHOKHARA_pipipi0pi0ISR')
615
616 elif finalstate == 'pi+pi-eta':
617 path.add_module(
618 'PhokharaInput',
619 FinalState=10, # pi+pi-eta
620 LO=0, NLO=1, QED=0, # use full two loop corrections
621 MinInvMassHadrons=0.,
622 eventType=eventType,
623 ).set_name('PHOKHARA_pipietaISR')
624
625 elif finalstate == 'K+K-':
626 path.add_module(
627 'PhokharaInput',
628 FinalState=6, # K+K-
629 LO=0, NLO=1, QED=0, # use full two loop corrections
630 MinInvMassHadrons=0.,
631 eventType=eventType,
632 ).set_name('PHOKHARA_KKISR')
633
634 elif finalstate == 'K0K0bar':
635 path.add_module(
636 'PhokharaInput',
637 FinalState=7, # K0K0bar
638 LO=0, NLO=1, QED=0, # use full two loop corrections
639 MinInvMassHadrons=0.,
640 eventType=eventType,
641 ).set_name('PHOKHARA_K0K0barISR')
642
643 elif finalstate == 'ppbar':
644 path.add_module(
645 'PhokharaInput',
646 FinalState=4, # ppbar
647 LO=0, NLO=1, QED=0, # use full two loop corrections
648 MinInvMassHadrons=0.,
649 eventType=eventType,
650 ).set_name('PHOKHARA_ppbarISR')
651
652 elif finalstate == 'n0n0bar':
653 path.add_module(
654 'PhokharaInput',
655 FinalState=5, # n0n0bar
656 LO=0, NLO=1, QED=0, # use full two loop corrections
657 MinInvMassHadrons=0.,
658 eventType=eventType,
659 ).set_name('PHOKHARA_n0n0barISR')
660
661 elif finalstate == 'Lambda0Lambda0bar':
662 b2.B2WARNING(
663 "The Lambda pair production and decays are implemented at the leading order only,"
664 " but as the expected number of events is modest, the accuracy of the code is"
665 " sufficient for the description of this process.")
666 path.add_module(
667 'PhokharaInput',
668 FinalState=9, # Lambda0Lambda0bar
669 LO=0, NLO=0, QED=0, # use only one loop corrections
670 MinInvMassHadrons=0.,
671 eventType=eventType,
672 ).set_name('PHOKHARA_Lambda0Lambda0barISR')
673
674 else:
675 b2.B2FATAL(f"add_phokhara_generator final state not supported: {finalstate}")
676
677
678def add_phokhara_evtgen_combination(
679 path, final_state_particles, user_decay_file,
680 beam_energy_spread=True, isr_events=False, min_inv_mass_vpho=0.0,
681 max_inv_mass_vpho=0.0, eventType=''):
682 """
683 Add combination of PHOKHARA and EvtGen to the path. Phokhara is
684 acting as ISR generator by generating e+ e- -> mu+ mu-, the muon pair is
685 then replaced by a virtual photon. Finally, the virtual photon is
686 decayed by EvtGen.
687
688 Parameters:
689 path (basf2.Path): Path where the generator should be added.
690 final_state_particles (list): List of final-state particles of
691 the virtual-photon decay. It is necessary to define the correct
692 mass threshold in PHOKHARA. For example, for the process
693 e+ e- -> J/psi eta_c, the list should be ['J/psi', 'eta_c'];
694 it does not depend on subsequent J/psi or eta_c decays.
695 user_decay_file (str): Name of EvtGen user decay file. The initial
696 particle must be the virtual photon (vpho).
697 beam_energy_spread (bool): Whether beam-energy spread should be
698 simulated
699 isr_events (bool): If true, then PHOKHARA is used with ph0 (LO in
700 basf2) = 0, i.e. generation of processes with one photon.
701 min_inv_mass_hadrons(float): Minimum invariant mass of the virtual
702 photon. This parameter is used only if isr_events is true,
703 otherwise the minimum mass is computed from the masses of
704 the final-state particles.
705 max_inv_mass_hadrons(float): Maximum invariant mass of the virtual
706 photon. This parameter is used only if isr_events is true,
707 otherwise the maximum mass is not restricted.
708 eventType (str) : event type information
709 """
710
711 import pdg
712
713 phokhara = path.add_module('PhokharaInput')
714 phokhara.param('eventType', eventType)
715
716 # Generate muons and replace them by a virtual photon.
717 phokhara.param('FinalState', 0)
718 phokhara.param('ReplaceMuonsByVirtualPhoton', True)
719
720 # Simulate beam-energy spread. This performs initialization for every
721 # event, and, thus, very slow, but usually necessary except for testing.
722 phokhara.param('BeamEnergySpread', beam_energy_spread)
723
724 # Soft photon cutoff.
725 phokhara.param('Epsilon', 0.0001)
726
727 # Maximum search.
728 phokhara.param('SearchMax', 5000)
729
730 # Number of generation attempts for event.
731 phokhara.param('nMaxTrials', 25000)
732
733 # Use NNLO.
734 if isr_events:
735 phokhara.param('LO', 0)
736 else:
737 phokhara.param('LO', 1)
738 phokhara.param('NLO', 1)
739
740 # Use ISR only.
741 phokhara.param('QED', 0)
742
743 # No interference.
744 phokhara.param('IFSNLO', 0)
745
746 # Vacuum polarization by Fred Jegerlehner.
747 phokhara.param('Alpha', 1)
748
749 # Do not include narrow resonances.
750 phokhara.param('NarrowRes', 0)
751
752 # Angular ragnes.
753 phokhara.param('ScatteringAngleRangePhoton', [0., 180.])
754 phokhara.param('ScatteringAngleRangeFinalStates', [0., 180.])
755
756 # Minimal invariant mass of the muons and tagged photon combination.
757 phokhara.param('MinInvMassHadronsGamma', 0.)
758
759 if isr_events:
760 # Minimum squared invariant mass of muons (final state).
761 phokhara.param('MinInvMassHadrons',
762 min_inv_mass_vpho * min_inv_mass_vpho)
763 # Maximum squared invariant mass of muons (final state).
764 phokhara.param('MaxInvMassHadrons',
765 max_inv_mass_vpho * max_inv_mass_vpho)
766 else:
767 # Minimum squared invariant mass of muons (final state).
768 # Force application of this cut.
769 mass = 0
770 for particle in final_state_particles:
771 p = pdg.get(particle)
772 mass = mass + p.Mass()
773 phokhara.param('MinInvMassHadrons', mass * mass)
774 phokhara.param('ForceMinInvMassHadronsCut', True)
775 # Maximum squared invariant mass of muons (final state).
776 # It is set to a large value, resulting in unrestricted mass.
777 phokhara.param('MaxInvMassHadrons', 200.0)
778
779 # Minimal photon energy.
780 phokhara.param('MinEnergyGamma', 0.01)
781
782 # EvtGen.
783 evtgen_decay = path.add_module('EvtGenDecay')
784 evtgen_decay.param('UserDecFile', user_decay_file)
785
786
787def add_koralw_generator(path, finalstate='', enableTauDecays=True, eventType=''):
788 """
789 Add KoralW generator for radiative four fermion final states (only four leptons final states are currently supported).
790
791 Parameters:
792 path (basf2.Path): path where the generator should be added
793 finalstate (str): either 'e+e-e+e-', 'e+e-mu+mu-', 'e+e-tau+tau-', 'mu+mu-mu+mu-', 'mu+mu-tau+tau-' or 'tau+tau-tau+tau-'
794 enableTauDecays (bool): if True, allow tau leptons to decay (using EvtGen)
795 eventType (str) : event type information
796 """
797
798 decayFile = ''
799 if finalstate == 'e+e-e+e-':
800 decayFile = b2.find_file('data/generators/koralw/KoralW_eeee.data')
801 elif finalstate == 'e+e-mu+mu-':
802 decayFile = b2.find_file('data/generators/koralw/KoralW_eeMuMu.data')
803 elif finalstate == 'e+e-tau+tau-':
804 decayFile = b2.find_file('data/generators/koralw/KoralW_eeTauTau.data')
805 elif finalstate == 'mu+mu-mu+mu-':
806 decayFile = b2.find_file('data/generators/koralw/KoralW_MuMuMuMu.data')
807 elif finalstate == 'mu+mu-tau+tau-':
808 decayFile = b2.find_file('data/generators/koralw/KoralW_MuMuTauTau.data')
809 elif finalstate == 'tau+tau-tau+tau-':
810 decayFile = b2.find_file('data/generators/koralw/KoralW_TauTauTauTau.data')
811 else:
812 b2.B2FATAL(f'add_koralw_generator final state not supported: {finalstate}')
813
814 path.add_module('KoralWInput',
815 UserDataFile=decayFile,
816 eventType=eventType)
817
818 if 'tau+tau-' in finalstate:
819 if enableTauDecays:
820 path.add_module('EvtGenDecay')
821 else:
822 b2.B2WARNING('The tau decays will not be generated.')
823
824
825def add_cosmics_generator(path, components=None,
826 global_box_size=None, accept_box=None, keep_box=None,
827 geometry_xml_file='geometry/Beast2_phase2.xml',
828 cosmics_data_dir='data/generators/modules/cryinput/',
829 setup_file='generators/scripts/cry.setup',
830 data_taking_period='gcr2017', top_in_counter=False):
831 """
832 Add the cosmics generator CRY with the default parameters to the path.
833
834 Warning:
835 Please remember to also change the reconstruction accordingly, if you
836 set "special" parameters here!
837
838 Parameters:
839 path (basf2.Path): path where the generator should be added
840 components (list(str)): list of geometry components to add in the
841 geometry module, or None for all components.
842 global_box_size (tuple(float, float, float)): sets global length, width
843 and height (in meters) in which to generate.
844 Default is ``[100, 100, 100]``
845 accept_box (tuple(float, float, float)): sets the size of the accept box in meter.
846 As a default it is set to ``[8.0, 8.0, 8.0]`` (the Belle II detector size).
847 keep_box (tuple(float, float, float)): sets the size of the keep box (keep box >= accept box).
848 geometry_xml_file (str): Name of the xml file to use for the geometry.
849 cosmics_data_dir (str): parameter CosmicDataDir for the cry module (absolute or relative to the basf2 repo).
850 setup_file (str): location of the cry.setup file (absolute or relative to the basf2 repo)
851 data_taking_period (str): The cosmics generation will be added using the
852 parameters, that where used in this period of data taking. The
853 periods can be found in ``cdc/cr/__init__.py``.
854 top_in_counter (bool): time of propagation from the hit point to the PMT in the trigger counter is subtracted
855 (assuming PMT is put at -z of the counter).
856 """
857
858 b2.B2FATAL('''The function "add_cosmics_generator()" is outdated and it is currently not working: please replace
859
860 add_cosmics_generator(path=path)
861
862with
863
864 path.add_module('CRYInput')
865
866in your steering file (the module parameter "acceptance" has to be set, see the module documentation).''')
867
868 import cdc.cr as cosmics_setup
869
870 if global_box_size is None:
871 global_box_size = [100, 100, 100]
872 if accept_box is None:
873 accept_box = [8, 8, 8]
874 if keep_box is None:
875 keep_box = [8, 8, 8]
876
877 cosmics_setup.set_cdc_cr_parameters(data_taking_period)
878
879 if cosmics_setup.cosmics_period == "201607":
880 b2.B2FATAL("The data taking period 201607 is very special (geometry setup, PMTs etc). This is not handled "
881 "by this script! Please ask the CDC group, if you want to simulate this.")
882
883 if 'Gearbox' not in path:
884 override = [("/Global/length", str(global_box_size[0]), "m"),
885 ("/Global/width", str(global_box_size[1]), "m"),
886 ("/Global/height", str(global_box_size[2]), "m")]
887
888 if cosmics_setup.globalPhi:
889 override += [("/DetectorComponent[@name='CDC']//GlobalPhiRotation", str(cosmics_setup.globalPhi), "deg")]
890
891 path.add_module('Gearbox', override=override, fileName=geometry_xml_file,)
892
893 # detector geometry
894 if 'Geometry' not in path:
895 geometry = path.add_module('Geometry')
896 if components:
897 geometry.param('components', components)
898
899 cry = path.add_module('CRYInput')
900
901 # cosmic data input
902 cry.param('CosmicDataDir', b2.find_file(cosmics_data_dir))
903
904 # user input file
905 cry.param('SetupFile', b2.find_file(setup_file))
906
907 # acceptance half-lengths - at least one particle has to enter that box to use that event
908 cry.param('acceptLength', accept_box[0])
909 cry.param('acceptWidth', accept_box[1])
910 cry.param('acceptHeight', accept_box[2])
911 cry.param('maxTrials', 100000)
912
913 # keep half-lengths - all particles that do not enter the box are removed (keep box >= accept box)
914 # default was 6.0
915 cry.param('keepLength', keep_box[0])
916 cry.param('keepWidth', keep_box[1])
917 cry.param('keepHeight', keep_box[2])
918
919 # minimal kinetic energy - all particles below that energy are ignored
920 cry.param('kineticEnergyThreshold', 0.01)
921
922 # TODO: I still do not fully understand, when the cosmics selector is needed and when not
923 if cosmics_setup.cosmics_period not in ["normal", "gcr2017"]:
924 # Selector module.
925 cosmics_selector = b2.register_module('CDCCosmicSelector',
926 lOfCounter=cosmics_setup.lengthOfCounter,
927 wOfCounter=cosmics_setup.widthOfCounter,
928 xOfCounter=cosmics_setup.triggerPos[0],
929 yOfCounter=cosmics_setup.triggerPos[1],
930 zOfCounter=cosmics_setup.triggerPos[2],
931 phiOfCounter=0.,
932 TOP=top_in_counter,
933 propSpeed=cosmics_setup.lightPropSpeed,
934 TOF=1,
935 cryGenerator=True
936 )
937
938 path.add_module(cosmics_selector)
939
940 empty_path = b2.create_path()
941 cosmics_selector.if_false(empty_path)
942
943
944def _get_treps_finalstates():
945 """Return the TREPS final-state bookkeeping table."""
946 tableFile = b2.find_file('generators/treps/data/treps_finalstates.json')
947 with open(tableFile) as table:
948 return json.load(table)
949
950
951def add_treps_generator(path, finalstate='', useDiscreteAndSortedW=False, eventType=''):
952 """
953 Add TREPS generator to produce hadronic two-photon processes.
954
955 Parameters:
956 path (basf2.Path): path where the generator should be added
957 finalstate(str): "e+e-pi+pi-", "e+e-K+K-" or "e+e-ppbar"
958 useDiscreteAndSortedW(bool): if True, wListTableFile is used for discrete and sorted W. evtNumList must be set proper value.
959 eventType (str) : event type information
960 """
961
962 try:
963 finalstateFiles = _get_treps_finalstates()[finalstate]
964 except KeyError:
965 b2.B2FATAL("add_treps_generator final state not supported: {}".format(finalstate))
966
967 parameterFile = b2.find_file(finalstateFiles['parameter'])
968 differentialCrossSectionFile = b2.find_file(finalstateFiles['differential_cross_section'])
969 wListTableFile = b2.find_file(finalstateFiles['w_list'])
970
971 # use TREPS to generate two-photon events.
972 path.add_module(
973 'TrepsInput',
974 ParameterFile=parameterFile,
975 DifferentialCrossSectionFile=differentialCrossSectionFile,
976 WListTableFile=wListTableFile,
977 UseDiscreteAndSortedW=useDiscreteAndSortedW,
978 MaximalQ2=1.0,
979 MaximalAbsCosTheta=1.01,
980 ApplyCosThetaCutCharged=True,
981 MinimalTransverseMomentum=0,
982 ApplyTransverseMomentumCutCharged=True,
983 eventType=eventType,
984 )
get(name)
Definition pdg.py:48
load(filename)
Definition pdg.py:122