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