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