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