10This module contains convenience functions to setup most commonly used physics
11generators correctly with their default settings. More information can be found
12in `BELLE2-NOTE-PH-2015-006`_
14.. _BELLE2-NOTE-PH-2015-006: https://docs.belle2.org/record/282
22def get_default_decayfile():
23 """Return the default DECAY.dec for Belle2"""
24 return b2.find_file(
"decfiles/dec/DECAY_BELLE2.DEC")
27def add_generator_preselection(
35 MaxChargedTheta=180.0,
42 stableParticles=False):
44 Adds generator preselection.
45 Should be added to the path after the generator.add_abc_generator but before simulation.add_simulation modules
46 It uses all particles from the event generator (i.e. primary, non-virtual, non-initial particles).
47 It checks if the required conditions are fulfilled.
48 If not, the events are given to the emptypath.
49 The main use case is a reduction of simulation time.
50 Note that you have to multiply the generated cross section by the retention fraction of the preselection.
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
70 generatorpreselection = path.add_module(
'GeneratorPreselection',
71 nChargedMin=nChargedMin,
72 nChargedMax=nChargedMax,
73 MinChargedP=MinChargedP,
74 MinChargedPt=MinChargedPt,
75 MinChargedTheta=MinChargedTheta,
76 MaxChargedTheta=MaxChargedTheta,
77 nPhotonMin=nPhotonMin,
78 nPhotonMax=nPhotonMax,
79 MinPhotonEnergy=MinPhotonEnergy,
80 MinPhotonTheta=MinPhotonTheta,
81 MaxPhotonTheta=MaxPhotonTheta,
82 applyInCMS=applyInCMS,
83 stableParticles=stableParticles
87 generatorpreselection.if_value(
'<11', emptypath)
90def add_aafh_generator(
101 Add the default two photon generator for four fermion final states
104 path (basf2.Path): path where the generator should be added
105 finalstate (str): either "e+e-e+e-", "e+e-mu+mu-", "e+e-tau+tau-", "mu+mu-mu+mu-" or "mu+mu-tau+tau-"
106 preselection (bool): if True, select events with at least one medium pt particle in the CDC acceptance
107 enableTauDecays (bool): if True, allow tau leptons to decay (using EvtGen)
108 minmass (float): minimum invariant mass
109 subweights (list(float)): list of four or eight values (first four are interpreted as WAP, rest as WBP)
110 which specify the relative weights for each of the four sub generators
111 maxsubweight (float): maximum expected subgenerator weight for rejection scheme
112 maxfinalweight (float): maximum expected final weight for rejection scheme
113 eventType (str) : event type information
116 if finalstate ==
'e+e-e+e-':
119 aafh_subgeneratorWeights = [1.0, 7.986e+01, 5.798e+04, 3.898e+05, 1.0, 1.664e+00, 2.812e+00, 7.321e-01]
121 aafh_subgeneratorWeights = subweights
122 if abs(minmass - 0.5) > 0.01
and not subweights:
123 b2.B2WARNING(
"add_aafh_generator: non default invariant mass cut without updated subweights requested!")
124 elif finalstate ==
'e+e-mu+mu-':
127 aafh_subgeneratorWeights = [1.000e+00, 1.520e+01, 3.106e+03, 6.374e+03, 1.000e+00, 1.778e+00, 6.075e+00, 6.512e+00]
129 aafh_subgeneratorWeights = subweights
130 if abs(minmass - 0.5) > 0.01
and not subweights:
131 b2.B2WARNING(
"add_aafh_generator: non default invariant mass cut without updated subweights requested!")
132 elif finalstate ==
'e+e-tau+tau-':
136 aafh_subgeneratorWeights = [1.000e+00, 2.214e+00, 1.202e+01, 1.536e+01, 1.000e+00, 1.664e+00, 1.680e+01, 6.934e+00]
138 aafh_subgeneratorWeights = subweights
141 f
"You requested a generator preselection for the final state {finalstate}: "
142 "please consider to remove it, since the cross section is small.")
143 elif finalstate ==
'mu+mu-mu+mu-':
148 aafh_subgeneratorWeights = [0.000e+00, 0.000e+00, 1.000e+00, 3.726e+00, 1.000e+00, 1.778e+00, 1.000e+00, 1.094e+00]
150 aafh_subgeneratorWeights = subweights
153 f
"You requested a generator preselection for the final state {finalstate}: "
154 "please consider to remove it, since the cross section is small.")
155 elif finalstate ==
'mu+mu-tau+tau-':
160 aafh_subgeneratorWeights = [0.000e+00, 0.000e+00, 1.000e+00, 1.715e+00, 1.000e+00, 1.778e+00, 1.000e+00, 6.257e-01]
162 aafh_subgeneratorWeights = subweights
165 f
"You requested a generator preselection for the final state {finalstate}: "
166 "please consider to remove it, since the cross section is small.")
167 elif finalstate ==
'tau+tau-tau+tau-':
168 b2.B2FATAL(f
"AAFH is not able to generate the {finalstate} final state. Please use KoralW instead.")
170 b2.B2FATAL(f
"add_aafh_generator final state not supported: {finalstate}")
172 aafh_maxSubgeneratorWeight = maxsubweight
173 aafh_maxFinalWeight = maxfinalweight
179 maxSubgeneratorWeight=aafh_maxSubgeneratorWeight,
180 maxFinalWeight=aafh_maxFinalWeight,
181 subgeneratorWeights=aafh_subgeneratorWeights,
182 suppressionLimits=[1e100] * 4,
188 generator_emptypath = b2.create_path()
189 add_generator_preselection(
191 emptypath=generator_emptypath,
194 MinChargedTheta=17.0,
195 MaxChargedTheta=150.0)
197 if 'tau+tau-' in finalstate:
199 path.add_module(
'EvtGenDecay')
201 b2.B2WARNING(
"The tau decays will not be generated.")
204def add_kkmc_generator(path, finalstate='', signalconfigfile='', useTauolaBelle=False, tauinputfile='', eventType=''):
206 Add the default muon pair and tau pair generator KKMC.
207 For tau decays, TauolaBelle and TauolaBelle2 are available.
208 Signal events can be produced setting a configuration file. Please notice that the configuration files for
209 TauolaBelle and TauolaBelle2 has a very different structure (see the examples below generators/examples).
212 path (basf2.Path): path where the generator should be added
213 finalstate(str): either "mu-mu+" or "tau-tau+"
214 signalconfigfile(str): File with configuration of the signal event to generate. It doesn't affect mu-mu+ decays.
215 useTauolaBelle(bool): If true, tau decay is driven by TauolaBelle. Otherwise TauolaBelle2 is used.
216 It doesn't affect mu-mu+ decays.
217 tauinputfile(str) : File to override KK2f_defaults. Only [sometimes] needed when tau decay is driven by TauolaBelle.
218 eventType (str) : event type information
222 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/tauola_bbb.input.dat')
225 kkmc_logfile =
'kkmc_tautau.txt'
228 kkmc_config = b2.find_file(
'data/generators/kkmc/KK2f_defaults.dat')
231 kkmc_tauconfigfile =
''
233 if finalstate ==
'tau+tau-':
234 b2.B2WARNING(
"add_kkmc_generator: please set finalstate as 'tau-tau+'. 'tau+tau-' will be deprecated in the future"
235 " for consistency in the configuration files.")
236 finalstate =
'tau-tau+'
237 if finalstate ==
'mu+mu-':
238 b2.B2WARNING(
"add_kkmc_generator: please set finalstate as 'mu-mu+'. 'mu+mu-' will be deprecated in the future for"
239 " consistency in the configuration files.")
240 finalstate =
'mu-mu+'
242 if finalstate ==
'tau-tau+':
244 b2.B2INFO(
"Generating tau pair events with TauolaBelle")
246 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/tau.input.dat')
247 kkmc_tauconfigfile = b2.find_file(
'data/generators/kkmc/tau_decaytable.dat')
249 if not signalconfigfile ==
'':
250 b2.B2INFO(f
"Using config file defined by user: {signalconfigfile}")
252 kkmc_tauconfigfile = b2.find_file(signalconfigfile)
254 kkmc_inputfile = b2.find_file(signalconfigfile)
256 if not tauinputfile ==
'':
257 kkmc_inputfile = b2.find_file(tauinputfile)
259 elif finalstate ==
'mu-mu+':
260 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/mu.input.dat')
261 kkmc_logfile =
'kkmc_mumu.txt'
264 b2.B2FATAL(
"add_kkmc_generator final state not supported: {}".format(finalstate))
269 tauinputFile=kkmc_inputfile,
270 KKdefaultFile=kkmc_config,
271 taudecaytableFile=kkmc_tauconfigfile,
272 kkmcoutputfilename=kkmc_logfile,
277def add_evtgen_generator(path, finalstate='', signaldecfile=None, coherentMixing=True, parentParticle='Upsilon(4S)
', eventType=''):
279 Add EvtGen for mixed and charged BB
282 path (basf2.Path): path where the generator should be added
283 finalstate (str): Either "charged" for generation of generic B+/B-, "mixed" for generic B0/anti-B0, or "signal" for
284 generation of user-defined signal mode
285 signaldecfile (str): path to decfile which defines the signal decay to be generated
286 (only needed if ``finalstate`` set to "signal")
287 coherentMixing: Either True or False. Switches on or off the coherent decay of the B0-B0bar pair.
288 It should always be True, unless you are generating Y(5,6S) -> BBar. In the latter case,
289 setting it False solves the internal limitation of Evtgen that allows to make a
290 coherent decay only starting from the Y(4S).
291 parentParticle (str): initial state (used only if it is not Upsilon(4S).
292 eventType (str) : event type information
294 evtgen_userdecfile = b2.find_file(
'data/generators/evtgen/charged.dec')
296 if parentParticle !=
'Upsilon(3S)' and parentParticle !=
'Upsilon(4S)'\
297 and parentParticle !=
'Upsilon(5S)' and parentParticle !=
'Upsilon(6S)':
298 b2.B2FATAL(
"add_evtgen_generator initial state not supported: {}".format(parentParticle))
300 if finalstate ==
'charged':
302 elif finalstate ==
'mixed':
303 evtgen_userdecfile = b2.find_file(
'data/generators/evtgen/mixed.dec')
304 elif finalstate ==
'signal':
305 evtgen_userdecfile = signaldecfile
307 b2.B2FATAL(
"add_evtgen_generator final state not supported: {}".format(finalstate))
309 if signaldecfile
and finalstate
in [
'charged',
'mixed']:
310 b2.B2WARNING(
"ignoring decfile: {}".format(signaldecfile))
313 if parentParticle ==
'Upsilon(3S)':
314 if finalstate !=
'signal':
315 b2.B2FATAL(
"add_evtgen_generator initial state {} is supported only with 'signal' final state".format(parentParticle))
317 coherentMixing =
False
318 b2.B2WARNING(
"add_evtgen_generator initial state {} has no BB mixing, now switching coherentMixing OFF"
319 .format(parentParticle))
321 if parentParticle ==
'Upsilon(5S)':
322 if finalstate !=
'signal':
323 b2.B2FATAL(
"add_evtgen_generator initial state {} is supported only with 'signal' final state".format(parentParticle))
325 coherentMixing =
False
327 "add_evtgen_generator initial state {} is supported only with false coherentMixing, now switching it OFF"
328 .format(parentParticle))
329 pdg.load(b2.find_file(
'decfiles/dec/Y5S.pdl'))
331 if parentParticle ==
'Upsilon(6S)':
332 if finalstate !=
'signal':
333 b2.B2FATAL(
"add_evtgen_generator initial state {} is supported only with 'signal' final state".format(parentParticle))
335 coherentMixing =
False
337 "add_evtgen_generator initial state {} is supported only with false coherentMixing, now switching it OFF"
338 .format(parentParticle))
339 pdg.load(b2.find_file(
'decfiles/dec/Y6S.pdl'))
343 userDECFile=evtgen_userdecfile,
344 CoherentMixing=coherentMixing,
345 ParentParticle=parentParticle,
350def add_continuum_generator(path, finalstate, userdecfile='', *, skip_on_failure=True, eventType=''):
352 Add the default continuum generators KKMC + PYTHIA including their default decfiles and PYTHIA settings
355 `add_inclusive_continuum_generator()` to add continuum generation with preselected particles
358 path (basf2.Path): path where the generator should be added
359 finalstate (str): uubar, ddbar, ssbar, ccbar
360 userdecfile (str): EvtGen decfile used for particle decays
361 skip_on_failure (bool): If True stop event processing right after
362 fragmentation fails. Otherwise continue normally
363 eventType (str) : event type information
367 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/uubar_nohadronization.input.dat')
370 kkmc_logfile =
'kkmc_uubar.txt'
373 pythia_config = b2.find_file(
'data/generators/modules/fragmentation/pythia_belle2.dat')
376 decay_user = b2.find_file(
'data/generators/modules/fragmentation/dec_belle2_qqbar.dec')
377 if userdecfile ==
'':
380 b2.B2INFO(
'Replacing default user decfile: {}'.format(userdecfile))
381 decay_user = userdecfile
384 kkmc_config = b2.find_file(
'data/generators/kkmc/KK2f_defaults.dat')
387 decay_file = get_default_decayfile()
389 if finalstate ==
'uubar':
391 elif finalstate ==
'ddbar':
392 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/ddbar_nohadronization.input.dat')
393 kkmc_logfile =
'kkmc_ddbar.txt'
394 elif finalstate ==
'ssbar':
395 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/ssbar_nohadronization.input.dat')
396 kkmc_logfile =
'kkmc_ssbar.txt'
397 elif finalstate ==
'ccbar':
398 kkmc_inputfile = b2.find_file(
'data/generators/kkmc/ccbar_nohadronization.input.dat')
399 pythia_config = b2.find_file(
'data/generators/modules/fragmentation/pythia_belle2_charm.dat')
400 kkmc_logfile =
'kkmc_ccbar.txt'
402 b2.B2FATAL(
"add_continuum_generator final state not supported: {}".format(finalstate))
407 tauinputFile=kkmc_inputfile,
408 KKdefaultFile=kkmc_config,
409 taudecaytableFile=
'',
410 kkmcoutputfilename=kkmc_logfile,
416 fragmentation = path.add_module(
418 ParameterFile=pythia_config,
422 UserDecFile=decay_user,
428 generator_emptypath = b2.create_path()
429 fragmentation.if_value(
'<1', generator_emptypath)
432def add_inclusive_continuum_generator(path, finalstate, particles, userdecfile='',
433 *, include_conjugates=True, max_iterations=100000, eventType=''):
435 Add continuum generation but require at least one of the given particles be
436 present in the event.
438 For example to only generate ccbar events which contain a "D*+" or an
439 electron one could would use
441 >>> add_inclusive_continuum_generator(path, "ccbar", ["D*+", 11])
443 If you are unsure how the particles are named in Belle II please have a look
444 at the ``b2help-particles`` executable or the `pdg` python module.
447 `add_continuum_generator()` to add continuum generation without preselection
450 path (basf2.Path): path to which the generator should be added
451 finalstate (str): uubar, ddbar, ssbar, ccbar
452 particles (list): A list of particle names or pdg codes. An event is
453 only accepted if at lease one of those particles appears in the event.
454 userdecfile (str): EvtGen decfile used for particle decays
455 include_conjugates (bool): If True (default) accept the event also if a
456 charge conjugate of the given particles is found
457 max_iterations (int): maximum tries per event to generate the requested
458 particle. If exceeded processing will be stopped with a
459 `FATAL <LogLevel.FATAL>` error so for rare particles one might need a
461 eventType (str) : event type information
463 loop_path = b2.create_path()
466 loop_path.add_module(
"PruneDataStore", keepMatchedEntries=
False, matchEntries=[
"MCParticles",
"EventExtraInfo"])
469 add_continuum_generator(loop_path, finalstate, userdecfile, skip_on_failure=
False, eventType=eventType)
471 loop_path.add_module(
"InclusiveParticleChecker", particles=particles, includeConjugates=include_conjugates)
473 path.do_while(loop_path, max_iterations=max_iterations)
476def add_bhwide_generator(path, minangle=0.5):
478 Add the high precision QED generator BHWIDE to the path. Settings are the default L1/HLT study settings
479 with a cross section of about 124000 nb (!)
482 path (basf2.Path): path where the generator should be added
483 minangle (float): minimum angle of the outgoing electron/positron in the CMS in degrees
486 if minangle < 0.0
or minangle > 180.0:
487 b2.B2FATAL(
"add_bhwide_generator minimum angle too small (<0.0) or too large (>180): {}".format(minangle))
489 bhwide = path.add_module(
"BHWideInput")
490 bhwide.param(
'ScatteringAngleRangePositron', [minangle, 180.0 - minangle])
491 bhwide.param(
'ScatteringAngleRangeElectron', [minangle, 180.0 - minangle])
492 bhwide.param(
'MaxAcollinearity', 180.0)
493 bhwide.param(
'MinEnergy', 0.10)
494 bhwide.param(
'VacuumPolarization',
'burkhardt')
495 bhwide.param(
'WeakCorrections',
True)
496 bhwide.param(
'WtMax', 3.0)
499def add_babayaganlo_generator(
505 generateInECLAcceptance=False,
508 Add the high precision QED generator BabaYaga@NLO to the path.
511 path (basf2.Path): path where the generator should be added.
512 finalstate (str): ee (e+e-), mm(mu+mu-) or gg (gammagamma).
513 minenergy (float): minimum particle (leptons for 'ee' or 'mm', photons for 'gg') energy in GeV.
514 minangle (float): angular range from minangle to 180-minangle for primary particles (in degrees).
515 fmax (float): maximum of differential cross section weight. This parameter should be set only by experts.
516 generateInECLAcceptance (bool): if True, the GeneratorPreselection module is used to select only events
517 with both the primary particles within the ECL acceptance.
518 eventType (str) : event type information
521 if finalstate !=
'ee' and finalstate !=
'gg' and finalstate !=
'mm':
522 b2.B2FATAL(f
'add_babayaganlo_generator final state not supported: {finalstate}')
524 if not (fmax == -1.0):
525 b2.B2WARNING(f
'The BabayagaNLOInput parameter "FMax" will be set to {fmax} instead to the default value (-1.0). '
526 'Please do not do this, unless you are extremely sure about this choice.')
531 FinalState=finalstate,
532 ScatteringAngleRange=[minangle, 180.0 - minangle],
537 if generateInECLAcceptance:
538 b2.B2INFO(f
'The final state {finalstate} is preselected requiring both primary particles within the ECL acceptance.')
539 emptypath = b2.Path()
540 add_generator_preselection(path=path,
544 if finalstate ==
'ee' or finalstate ==
'mm':
545 b2.set_module_parameters(path=path,
546 name=
'GeneratorPreselection',
548 MinChargedTheta=12.4,
549 MaxChargedTheta=155.1,)
551 elif finalstate ==
'gg':
552 b2.set_module_parameters(path=path,
553 name=
'GeneratorPreselection',
556 MaxPhotonTheta=155.1)
559def add_phokhara_generator(path, finalstate='', eventType=''):
561 Add the high precision QED generator PHOKHARA to the path. Almost full
562 acceptance settings for photons and hadrons/muons.
565 path (basf2.Path): path where the generator should be added
566 finalstate (str): One of the following final states: "mu+mu-", "pi+pi-", "pi+pi-pi0", "pi+pi-pi+pi-" (or "2(pi+pi-)"),
567 "pi+pi-pi0pi0" or ("pi+pi-2pi0"), "pi+pi-eta", "K+K-", "K0K0bar", "ppbar", "n0n0bar" or "Lambda0Lambda0bar"
568 eventType (str) : event type information
571 if finalstate ==
'mu+mu-':
576 MinInvMassHadrons=0.,
578 ).set_name(
'PHOKHARA_mumuISR')
580 elif finalstate ==
'pi+pi-':
585 MinInvMassHadrons=0.,
587 ).set_name(
'PHOKHARA_pipiISR')
589 elif finalstate ==
'pi+pi-pi0':
594 MinInvMassHadrons=0.,
596 ).set_name(
'PHOKHARA_pipipi0ISR')
598 elif finalstate ==
'pi+pi-pi+pi-' or finalstate ==
'2(pi+pi-)':
603 MinInvMassHadrons=0.,
605 ).set_name(
'PHOKHARA_pipipipiISR')
607 elif finalstate ==
'pi+pi-pi0pi0' or finalstate ==
'pi+pi-2pi0':
612 MinInvMassHadrons=0.,
614 ).set_name(
'PHOKHARA_pipipi0pi0ISR')
616 elif finalstate ==
'pi+pi-eta':
621 MinInvMassHadrons=0.,
623 ).set_name(
'PHOKHARA_pipietaISR')
625 elif finalstate ==
'K+K-':
630 MinInvMassHadrons=0.,
632 ).set_name(
'PHOKHARA_KKISR')
634 elif finalstate ==
'K0K0bar':
639 MinInvMassHadrons=0.,
641 ).set_name(
'PHOKHARA_K0K0barISR')
643 elif finalstate ==
'ppbar':
648 MinInvMassHadrons=0.,
650 ).set_name(
'PHOKHARA_ppbarISR')
652 elif finalstate ==
'n0n0bar':
657 MinInvMassHadrons=0.,
659 ).set_name(
'PHOKHARA_n0n0barISR')
661 elif finalstate ==
'Lambda0Lambda0bar':
663 "The Lambda pair production and decays are implemented at the leading order only,"
664 " but as the expected number of events is modest, the accuracy of the code is"
665 " sufficient for the description of this process.")
670 MinInvMassHadrons=0.,
672 ).set_name(
'PHOKHARA_Lambda0Lambda0barISR')
675 b2.B2FATAL(f
"add_phokhara_generator final state not supported: {finalstate}")
678def add_phokhara_evtgen_combination(
679 path, final_state_particles, user_decay_file,
680 beam_energy_spread=True, isr_events=False, min_inv_mass_vpho=0.0,
681 max_inv_mass_vpho=0.0, eventType=''):
683 Add combination of PHOKHARA and EvtGen to the path. Phokhara is
684 acting as ISR generator by generating e+ e- -> mu+ mu-, the muon pair is
685 then replaced by a virtual photon. Finally, the virtual photon is
689 path (basf2.Path): Path where the generator should be added.
690 final_state_particles (list): List of final-state particles of
691 the virtual-photon decay. It is necessary to define the correct
692 mass threshold in PHOKHARA. For example, for the process
693 e+ e- -> J/psi eta_c, the list should be ['J/psi', 'eta_c'];
694 it does not depend on subsequent J/psi or eta_c decays.
695 user_decay_file (str): Name of EvtGen user decay file. The initial
696 particle must be the virtual photon (vpho).
697 beam_energy_spread (bool): Whether beam-energy spread should be
699 isr_events (bool): If true, then PHOKHARA is used with ph0 (LO in
700 basf2) = 0, i.e. generation of processes with one photon.
701 min_inv_mass_hadrons(float): Minimum invariant mass of the virtual
702 photon. This parameter is used only if isr_events is true,
703 otherwise the minimum mass is computed from the masses of
704 the final-state particles.
705 max_inv_mass_hadrons(float): Maximum invariant mass of the virtual
706 photon. This parameter is used only if isr_events is true,
707 otherwise the maximum mass is not restricted.
708 eventType (str) : event type information
713 phokhara = path.add_module(
'PhokharaInput')
714 phokhara.param(
'eventType', eventType)
717 phokhara.param(
'FinalState', 0)
718 phokhara.param(
'ReplaceMuonsByVirtualPhoton',
True)
722 phokhara.param(
'BeamEnergySpread', beam_energy_spread)
725 phokhara.param(
'Epsilon', 0.0001)
728 phokhara.param(
'SearchMax', 5000)
731 phokhara.param(
'nMaxTrials', 25000)
735 phokhara.param(
'LO', 0)
737 phokhara.param(
'LO', 1)
738 phokhara.param(
'NLO', 1)
741 phokhara.param(
'QED', 0)
744 phokhara.param(
'IFSNLO', 0)
747 phokhara.param(
'Alpha', 1)
750 phokhara.param(
'NarrowRes', 0)
753 phokhara.param(
'ScatteringAngleRangePhoton', [0., 180.])
754 phokhara.param(
'ScatteringAngleRangeFinalStates', [0., 180.])
757 phokhara.param(
'MinInvMassHadronsGamma', 0.)
761 phokhara.param(
'MinInvMassHadrons',
762 min_inv_mass_vpho * min_inv_mass_vpho)
764 phokhara.param(
'MaxInvMassHadrons',
765 max_inv_mass_vpho * max_inv_mass_vpho)
770 for particle
in final_state_particles:
772 mass = mass + p.Mass()
773 phokhara.param(
'MinInvMassHadrons', mass * mass)
774 phokhara.param(
'ForceMinInvMassHadronsCut',
True)
777 phokhara.param(
'MaxInvMassHadrons', 200.0)
780 phokhara.param(
'MinEnergyGamma', 0.01)
783 evtgen_decay = path.add_module(
'EvtGenDecay')
784 evtgen_decay.param(
'UserDecFile', user_decay_file)
787def add_koralw_generator(path, finalstate='', enableTauDecays=True, eventType=''):
789 Add KoralW generator for radiative four fermion final states (only four leptons final states are currently supported).
792 path (basf2.Path): path where the generator should be added
793 finalstate (str): either 'e+e-e+e-', 'e+e-mu+mu-', 'e+e-tau+tau-', 'mu+mu-mu+mu-', 'mu+mu-tau+tau-' or 'tau+tau-tau+tau-'
794 enableTauDecays (bool): if True, allow tau leptons to decay (using EvtGen)
795 eventType (str) : event type information
799 if finalstate ==
'e+e-e+e-':
800 decayFile = b2.find_file(
'data/generators/koralw/KoralW_eeee.data')
801 elif finalstate ==
'e+e-mu+mu-':
802 decayFile = b2.find_file(
'data/generators/koralw/KoralW_eeMuMu.data')
803 elif finalstate ==
'e+e-tau+tau-':
804 decayFile = b2.find_file(
'data/generators/koralw/KoralW_eeTauTau.data')
805 elif finalstate ==
'mu+mu-mu+mu-':
806 decayFile = b2.find_file(
'data/generators/koralw/KoralW_MuMuMuMu.data')
807 elif finalstate ==
'mu+mu-tau+tau-':
808 decayFile = b2.find_file(
'data/generators/koralw/KoralW_MuMuTauTau.data')
809 elif finalstate ==
'tau+tau-tau+tau-':
810 decayFile = b2.find_file(
'data/generators/koralw/KoralW_TauTauTauTau.data')
812 b2.B2FATAL(f
'add_koralw_generator final state not supported: {finalstate}')
814 path.add_module(
'KoralWInput',
815 UserDataFile=decayFile,
818 if 'tau+tau-' in finalstate:
820 path.add_module(
'EvtGenDecay')
822 b2.B2WARNING(
'The tau decays will not be generated.')
825def add_cosmics_generator(path, components=None,
826 global_box_size=None, accept_box=None, keep_box=None,
827 geometry_xml_file='geometry/Beast2_phase2.xml',
828 cosmics_data_dir='data/generators/modules/cryinput/',
829 setup_file='generators/scripts/cry.setup',
830 data_taking_period='gcr2017', top_in_counter=False):
832 Add the cosmics generator CRY with the default parameters to the path.
835 Please remember to also change the reconstruction accordingly, if you
836 set "special" parameters here!
839 path (basf2.Path): path where the generator should be added
840 components (list(str)): list of geometry components to add in the
841 geometry module, or None for all components.
842 global_box_size (tuple(float, float, float)): sets global length, width
843 and height (in meters) in which to generate.
844 Default is ``[100, 100, 100]``
845 accept_box (tuple(float, float, float)): sets the size of the accept box in meter.
846 As a default it is set to ``[8.0, 8.0, 8.0]`` (the Belle II detector size).
847 keep_box (tuple(float, float, float)): sets the size of the keep box (keep box >= accept box).
848 geometry_xml_file (str): Name of the xml file to use for the geometry.
849 cosmics_data_dir (str): parameter CosmicDataDir for the cry module (absolute or relative to the basf2 repo).
850 setup_file (str): location of the cry.setup file (absolute or relative to the basf2 repo)
851 data_taking_period (str): The cosmics generation will be added using the
852 parameters, that where used in this period of data taking. The
853 periods can be found in ``cdc/cr/__init__.py``.
854 top_in_counter (bool): time of propagation from the hit point to the PMT in the trigger counter is subtracted
855 (assuming PMT is put at -z of the counter).
858 b2.B2FATAL(
'''The function "add_cosmics_generator()" is outdated and it is currently not working: please replace
860 add_cosmics_generator(path=path)
864 path.add_module('CRYInput')
866in your steering file (the module parameter "acceptance" has to be set, see the module documentation).''')
868 import cdc.cr
as cosmics_setup
870 if global_box_size
is None:
871 global_box_size = [100, 100, 100]
872 if accept_box
is None:
873 accept_box = [8, 8, 8]
877 cosmics_setup.set_cdc_cr_parameters(data_taking_period)
879 if cosmics_setup.cosmics_period ==
"201607":
880 b2.B2FATAL(
"The data taking period 201607 is very special (geometry setup, PMTs etc). This is not handled "
881 "by this script! Please ask the CDC group, if you want to simulate this.")
883 if 'Gearbox' not in path:
884 override = [(
"/Global/length", str(global_box_size[0]),
"m"),
885 (
"/Global/width", str(global_box_size[1]),
"m"),
886 (
"/Global/height", str(global_box_size[2]),
"m")]
888 if cosmics_setup.globalPhi:
889 override += [(
"/DetectorComponent[@name='CDC']//GlobalPhiRotation", str(cosmics_setup.globalPhi),
"deg")]
891 path.add_module(
'Gearbox', override=override, fileName=geometry_xml_file,)
894 if 'Geometry' not in path:
895 geometry = path.add_module(
'Geometry')
897 geometry.param(
'components', components)
899 cry = path.add_module(
'CRYInput')
902 cry.param(
'CosmicDataDir', b2.find_file(cosmics_data_dir))
905 cry.param(
'SetupFile', b2.find_file(setup_file))
908 cry.param(
'acceptLength', accept_box[0])
909 cry.param(
'acceptWidth', accept_box[1])
910 cry.param(
'acceptHeight', accept_box[2])
911 cry.param(
'maxTrials', 100000)
915 cry.param(
'keepLength', keep_box[0])
916 cry.param(
'keepWidth', keep_box[1])
917 cry.param(
'keepHeight', keep_box[2])
920 cry.param(
'kineticEnergyThreshold', 0.01)
923 if cosmics_setup.cosmics_period
not in [
"normal",
"gcr2017"]:
925 cosmics_selector = b2.register_module(
'CDCCosmicSelector',
926 lOfCounter=cosmics_setup.lengthOfCounter,
927 wOfCounter=cosmics_setup.widthOfCounter,
928 xOfCounter=cosmics_setup.triggerPos[0],
929 yOfCounter=cosmics_setup.triggerPos[1],
930 zOfCounter=cosmics_setup.triggerPos[2],
933 propSpeed=cosmics_setup.lightPropSpeed,
938 path.add_module(cosmics_selector)
940 empty_path = b2.create_path()
941 cosmics_selector.if_false(empty_path)
944def _get_treps_finalstates():
945 """Return the TREPS final-state bookkeeping table."""
946 tableFile = b2.find_file(
'generators/treps/data/treps_finalstates.json')
947 with open(tableFile)
as table:
948 return json.load(table)
951def add_treps_generator(path, finalstate='', useDiscreteAndSortedW=False, eventType=''):
953 Add TREPS generator to produce hadronic two-photon processes.
956 path (basf2.Path): path where the generator should be added
957 finalstate(str): "e+e-pi+pi-", "e+e-K+K-" or "e+e-ppbar"
958 useDiscreteAndSortedW(bool): if True, wListTableFile is used for discrete and sorted W. evtNumList must be set proper value.
959 eventType (str) : event type information
963 finalstateFiles = _get_treps_finalstates()[finalstate]
965 b2.B2FATAL(
"add_treps_generator final state not supported: {}".format(finalstate))
967 parameterFile = b2.find_file(finalstateFiles[
'parameter'])
968 differentialCrossSectionFile = b2.find_file(finalstateFiles[
'differential_cross_section'])
969 wListTableFile = b2.find_file(finalstateFiles[
'w_list'])
974 ParameterFile=parameterFile,
975 DifferentialCrossSectionFile=differentialCrossSectionFile,
976 WListTableFile=wListTableFile,
977 UseDiscreteAndSortedW=useDiscreteAndSortedW,
979 MaximalAbsCosTheta=1.01,
980 ApplyCosThetaCutCharged=
True,
981 MinimalTransverseMomentum=0,
982 ApplyTransverseMomentumCutCharged=
True,