11 #include <simulation/modules/fullsim/FullSimModule.h>
12 #include <simulation/kernel/RunManager.h>
13 #include <simulation/kernel/DetectorConstruction.h>
15 #include <simulation/physicslist/Belle2PhysicsList.h>
16 #include <simulation/kernel/ExtPhysicsConstructor.h>
17 #include <simulation/kernel/MagneticField.h>
18 #include <simulation/kernel/PrimaryGeneratorAction.h>
19 #include <simulation/kernel/EventAction.h>
20 #include <simulation/kernel/TrackingAction.h>
21 #include <simulation/kernel/SteppingAction.h>
22 #include <simulation/kernel/StackingAction.h>
24 #include <mdst/dataobjects/MCParticle.h>
25 #include <framework/datastore/StoreObjPtr.h>
26 #include <framework/datastore/StoreArray.h>
27 #include <framework/dataobjects/EventMetaData.h>
28 #include <framework/gearbox/Unit.h>
30 #include <CLHEP/Units/SystemOfUnits.h>
32 #include <simulation/monopoles/G4MonopolePhysics.h>
33 #include <simulation/longlivedneutral/G4LongLivedNeutralPhysics.h>
35 #include <G4TransportationManager.hh>
36 #include <G4Transportation.hh>
37 #include <G4PhysListFactory.hh>
38 #include <G4ProcessVector.hh>
39 #include <G4OpticalPhysics.hh>
40 #include <G4ParticleDefinition.hh>
41 #include <G4ParticleTable.hh>
42 #include <G4EventManager.hh>
43 #include <G4RunManager.hh>
44 #include <G4UImanager.hh>
45 #include <G4VisExecutive.hh>
46 #include <G4StepLimiter.hh>
47 #include <G4EmParameters.hh>
48 #include <G4HadronicProcessStore.hh>
49 #include <G4InuclParticleNames.hh>
51 #include <G4Mag_UsualEqRhs.hh>
52 #include <G4NystromRK4.hh>
53 #include <G4HelixExplicitEuler.hh>
54 #include <G4HelixSimpleRunge.hh>
55 #include <G4CachedMagneticField.hh>
59 using namespace Belle2::Simulation;
60 using namespace Belle2::Monopoles;
61 using namespace G4InuclParticleNames;
75 setDescription(
"Performs the full Geant4 detector simulation. Requires a valid geometry in memory.");
76 setPropertyFlags(c_ParallelProcessingCertified | c_TerminateInAllProcesses);
79 addParam(
"InputMCParticleCollection", m_mcParticleInputColName,
"The name of the input MCParticle collection.",
string(
""));
80 addParam(
"ThresholdImportantEnergy", m_thresholdImportantEnergy,
81 "[GeV] A particle which got 'stuck' and has less than this energy will be killed after 'ThresholdTrials' trials.", 0.250);
82 addParam(
"ThresholdTrials", m_thresholdTrials,
83 "Geant4 will try 'ThresholdTrials' times to move a particle which got 'stuck' and has an energy less than 'ThresholdImportantEnergy'.",
85 addParam(
"RunEventVerbosity", m_runEventVerbosity,
"Geant4 run/event verbosity: 0=silent; 1=info level; 2=debug level", 0);
86 addParam(
"TrackingVerbosity", m_trackingVerbosity,
87 "Tracking verbosity: 0=Silent; 1=Min info per step; 2=sec particles; 3=pre/post step info; 4=like 3 but more info; 5=proposed step length info.",
89 addParam(
"HadronProcessVerbosity", m_hadronProcessVerbosity,
"Hadron Process verbosity: 0=Silent; 1=info level; 2=debug level", 0);
90 addParam(
"EmProcessVerbosity", m_emProcessVerbosity,
"Em Process verbosity: 0=Silent; 1=info level; 2=debug level", 0);
91 addParam(
"PhysicsList", m_physicsList,
"The name of the physics list which is used for the simulation.",
string(
"Belle2"));
92 addParam(
"StandardEM", m_standardEM,
"If true, replaces fast EM physics with standard EM physics.",
false);
93 addParam(
"RegisterOptics", m_optics,
"If true, G4OpticalPhysics is registered in Geant4 PhysicsList.",
true);
94 addParam(
"UseHighPrecisionNeutrons", m_HPneutrons,
"If true, high precision neutron models used below 20 MeV.",
false);
95 addParam(
"RegisterMonopoles", m_monopoles,
"If set to true, G4MonopolePhysics is registered in Geant4 PhysicsList.",
false);
96 addParam(
"MonopoleMagCharge", m_monopoleMagneticCharge,
"The value of monopole magnetic charge in units of e+.", 1.0);
97 addParam(
"ProductionCut", m_productionCut,
98 "[cm] Apply continuous energy loss to primary particle which has no longer enough energy to produce secondaries which travel at least the specified productionCut distance.",
100 addParam(
"PXDProductionCut", m_pxdProductionCut,
"[cm] Secondary production threshold in PXD envelope.", 0.0);
101 addParam(
"SVDProductionCut", m_svdProductionCut,
"[cm] Secondary production threshold in SVD envelope.", 0.0);
102 addParam(
"CDCProductionCut", m_cdcProductionCut,
"[cm] Secondary production threshold in CDC envelope.", 0.0);
103 addParam(
"ARICHTOPProductionCut", m_arichtopProductionCut,
"[cm] Secondary production threshold in ARICH and TOP envelopes.", 0.0);
104 addParam(
"ECLProductionCut", m_eclProductionCut,
"[cm] Secondary production threshold in ECL envelope.", 0.0);
105 addParam(
"KLMProductionCut", m_klmProductionCut,
"[cm] Secondary production threshold in BKLM and EKLM envelopes.", 0.0);
106 addParam(
"MaxNumberSteps", m_maxNumberSteps,
107 "The maximum number of steps before the track transportation is stopped and the track is killed.", 100000);
108 addParam(
"PhotonFraction", m_photonFraction,
"The fraction of Cerenkov photons which will be kept and propagated.", 0.5);
109 addParam(
"EnableVisualization", m_EnableVisualization,
"If set to True, the Geant4 visualization support is enabled.",
false);
111 addParam(
"StoreOpticalPhotons", m_storeOpticalPhotons,
"If set to True, optical photons are stored in MCParticles.",
false);
112 addParam(
"StoreAllSecondaries", m_storeSecondaries,
113 "If set to True, all secondaries produced by Geant4 over a kinetic energy cut are stored in MCParticles. Otherwise do not store them.",
115 addParam(
"SecondariesEnergyCut", m_secondariesEnergyCut,
"[MeV] Kinetic energy cut for storing secondaries", 1.0);
116 addParam(
"StoreBremsstrahlungPhotons", m_storeBremsstrahlungPhotons,
117 "If set to True, store BremsstrahlungPhotons over a kinetic energy cut in MCParticles. Otherwise do not store them.",
false);
118 addParam(
"BremsstrahlungPhotonsEnergyCut", m_bremsstrahlungPhotonsEnergyCut,
119 "[MeV] Kinetic energy cut for storing bremsstrahlung photons", 10.0);
120 addParam(
"StorePairConversions", m_storePairConversions,
121 "If set to True, store e+ or e- from pair conversions over a kinetic energy cut in MCParticles. Otherwise do not store them.",
123 addParam(
"PairConversionsEnergyCut", m_pairConversionsEnergyCut,
124 "[MeV] Kinetic energy cut for storing e+ or e- from pair conversions", 10.0);
126 addParam(
"magneticField", m_magneticFieldName,
127 "Chooses the magnetic field stepper used by Geant4. Possible values are: default, nystrom, expliciteuler, simplerunge",
129 addParam(
"magneticCacheDistance", m_magneticCacheDistance,
130 "Minimum distance for BField lookup in cm. If the next requested point is closer than this distance than return the flast BField value. 0 means no caching",
132 addParam(
"deltaChordInMagneticField", m_deltaChordInMagneticField,
133 "[mm] The maximum miss-distance between the trajectory curve and its linear cord(s) approximation", 0.25);
134 vector<string> defaultCommandsAtPreInit;
135 addParam(
"UICommandsAtPreInit", m_uiCommandsAtPreInit,
136 "A list of Geant4 UI commands that should be applied at PreInit state, before the simulation starts.",
137 defaultCommandsAtPreInit);
138 vector<string> defaultCommandsAtIdle;
139 addParam(
"UICommandsAtIdle", m_uiCommandsAtIdle,
140 "A list of Geant4 UI commands that should be applied at Idle state, before the simulation starts.",
141 defaultCommandsAtIdle);
142 addParam(
"trajectoryStore", m_trajectoryStore,
143 "If non-zero save the full trajectory of 1=primary, 2=non-optical or 3=all particles", 0);
144 addParam(
"trajectoryDistanceTolerance", m_trajectoryDistanceTolerance,
145 "Maximum deviation from the real trajectory points when merging "
146 "segments (in cm)", 5e-4);
147 vector<float> defaultAbsorbers;
148 addParam(
"AbsorbersRadii", m_absorbers,
149 "Radii (in cm) of absorbers across which tracks will be destroyed.", defaultAbsorbers);
152 RunManager::Instance();
153 m_magneticField = NULL;
154 m_uncachedField = NULL;
155 m_magFldEquation = NULL;
157 m_chordFinder = NULL;
159 m_stepLimiter = NULL;
163 FullSimModule::~FullSimModule()
169 void FullSimModule::initialize()
176 if (m_mcParticleInputColName.empty()) {
190 RunManager& runManager = RunManager::Instance();
201 if (m_physicsList ==
"Belle2") {
218 if (m_uiCommandsAtPreInit.size() > 0) {
219 G4UImanager* uiManager = G4UImanager::GetUIpointer();
220 for (vector<string>::iterator iter = m_uiCommandsAtPreInit.begin(); iter != m_uiCommandsAtPreInit.end(); ++iter) {
221 uiManager->ApplyCommand(*iter);
225 runManager.SetUserInitialization(physicsList);
228 G4PhysListFactory physListFactory;
229 physListFactory.SetVerbose(m_runEventVerbosity);
230 G4VModularPhysicsList* physicsList = NULL;
231 if (physListFactory.IsReferencePhysList(m_physicsList)) physicsList = physListFactory.GetReferencePhysList(m_physicsList);
232 if (physicsList == NULL) B2FATAL(
"Could not load the physics list " << m_physicsList);
234 if (m_optics) physicsList->RegisterPhysics(
new G4OpticalPhysics);
241 physicsList->SetDefaultCutValue((m_productionCut / Unit::mm) * CLHEP::mm);
244 if (m_uiCommandsAtPreInit.size() > 0) {
245 G4UImanager* uiManager = G4UImanager::GetUIpointer();
246 for (vector<string>::iterator iter = m_uiCommandsAtPreInit.begin(); iter != m_uiCommandsAtPreInit.end(); ++iter) {
247 uiManager->ApplyCommand(*iter);
253 G4ParticleTable::G4PTblDicIterator* myParticleIterator = G4ParticleTable::GetParticleTable()->GetIterator();
254 myParticleIterator->reset();
255 while ((*myParticleIterator)()) {
256 G4ParticleDefinition* particle = myParticleIterator->value();
257 if (particle->GetParticleName().compare(0, 4,
"g4e_") == 0) {
258 physicsList->SetParticleCuts(1.0E+9 * CLHEP::cm, particle);
261 runManager.SetUserInitialization(physicsList);
265 if (m_magneticFieldName !=
"none") {
267 if (m_magneticCacheDistance > 0) {
268 m_uncachedField = m_magneticField;
269 m_magneticField =
new G4CachedMagneticField(m_uncachedField, m_magneticCacheDistance);
271 G4FieldManager* fieldManager = G4TransportationManager::GetTransportationManager()->GetFieldManager();
272 fieldManager->SetDetectorField(m_magneticField);
273 if (m_magneticFieldName !=
"default") {
276 m_magFldEquation =
new G4Mag_UsualEqRhs(m_magneticField);
277 if (m_magneticFieldName ==
"nystrom") {
278 m_stepper =
new G4NystromRK4(m_magFldEquation);
279 }
else if (m_magneticFieldName ==
"expliciteuler") {
280 m_stepper =
new G4HelixExplicitEuler(m_magFldEquation);
281 }
else if (m_magneticFieldName ==
"simplerunge") {
282 m_stepper =
new G4HelixSimpleRunge(m_magFldEquation);
284 B2FATAL(
"Unknown magnetic field option: " << m_magneticFieldName);
289 m_chordFinder =
new G4ChordFinder(m_magneticField, 1e-2 * CLHEP::mm, m_stepper);
290 fieldManager->SetChordFinder(m_chordFinder);
292 fieldManager->CreateChordFinder(m_magneticField);
296 G4ChordFinder* chordFinder = fieldManager->GetChordFinder();
297 B2DEBUG(1,
"Geant4 default deltaChord = " << chordFinder->GetDeltaChord());
298 chordFinder->SetDeltaChord(m_deltaChordInMagneticField * CLHEP::mm);
299 B2DEBUG(1,
"DeltaChord after reset = " << chordFinder->GetDeltaChord());
305 G4VUserPrimaryGeneratorAction* generatorAction =
new PrimaryGeneratorAction(m_mcParticleInputColName, m_mcParticleGraph);
306 runManager.SetUserAction(generatorAction);
311 runManager.SetUserAction(eventAction);
323 runManager.SetUserAction(trackingAction);
329 for (
auto& rAbsorber : m_absorbers) {
330 B2INFO(
"An absorber found at R = " << rAbsorber <<
" cm");
332 runManager.SetUserAction(steppingAction);
337 runManager.SetUserAction(stackingAction);
346 m_stepLimiter =
new G4StepLimiter();
347 G4ParticleTable::G4PTblDicIterator* partIter = G4ParticleTable::GetParticleTable()->GetIterator();
349 while ((*partIter)()) {
350 G4ParticleDefinition* currParticle = partIter->value();
351 G4ProcessVector& currProcList = *currParticle->GetProcessManager()->GetProcessList();
352 assert(currProcList.size() < INT_MAX);
353 for (
int iProcess = 0; iProcess < static_cast<int>(currProcList.size()); ++iProcess) {
354 G4Transportation* transport =
dynamic_cast<G4Transportation*
>(currProcList[iProcess]);
355 if (transport !=
nullptr) {
357 transport->SetThresholdImportantEnergy(m_thresholdImportantEnergy / Unit::MeV * CLHEP::MeV);
358 transport->SetThresholdTrials(m_thresholdTrials);
363 double zeroChargeTol = 0.01 * Unit::e;
364 if (fabs(currParticle->GetPDGCharge()) > zeroChargeTol) {
365 currParticle->GetProcessManager()->AddDiscreteProcess(m_stepLimiter);
366 B2DEBUG(100,
"Added StepLimiter process for " << currParticle->GetParticleName());
374 while ((*partIter)()) {
375 G4ParticleDefinition* currParticle = partIter->value();
376 if (currParticle->GetParticleName().compare(0, 4,
"g4e_") == 0) {
377 G4ProcessManager* processManager = currParticle->GetProcessManager();
378 if (processManager) {
379 G4ProcessVector* processList = processManager->GetProcessList();
380 assert(processList->size() < INT_MAX);
381 for (
int i = 0; i < static_cast<int>(processList->size()); ++i) {
382 if (((*processList)[i]->GetProcessName() ==
"Cerenkov") ||
383 ((*processList)[i]->GetProcessName() ==
"Scintillation") ||
384 ((*processList)[i]->GetProcessName() ==
"hFritiofCaptureAtRest")) {
385 processManager->SetProcessActivation(i,
false);
403 G4EventManager::GetEventManager()->SetVerboseLevel(m_runEventVerbosity);
404 G4RunManager::GetRunManager()->SetVerboseLevel(m_runEventVerbosity);
405 G4EventManager::GetEventManager()->GetTrackingManager()->SetVerboseLevel(
406 m_trackingVerbosity);
407 G4HadronicProcessStore::Instance()->SetVerbose(m_hadronProcessVerbosity);
408 G4EmParameters::Instance()->SetVerbose(m_emProcessVerbosity);
411 if (m_EnableVisualization) {
412 m_visManager =
new G4VisExecutive;
413 m_visManager->Initialize();
417 if (m_uiCommandsAtIdle.size() > 0) {
418 G4UImanager* uiManager = G4UImanager::GetUIpointer();
419 for (vector<string>::iterator iter = m_uiCommandsAtIdle.begin(); iter != m_uiCommandsAtIdle.end(); ++iter) {
420 uiManager->ApplyCommand(*iter);
425 if (m_trajectoryStore) {
434 B2INFO(
"Perform Geant4 final initialization: Geometry optimization, PhysicsList calculations...");
435 RunManager::Instance().beginRun(0);
436 B2INFO(
"done, Geant4 ready");
443 void FullSimModule::beginRun()
449 void FullSimModule::event()
455 RunManager::Instance().processEvent(eventMetaDataPtr->getEvent());
459 void FullSimModule::endRun()
464 void FullSimModule::terminate()
467 RunManager::Instance().endRun();
469 if (m_visManager !=
nullptr)
delete m_visManager;
470 RunManager::Instance().destroy();
472 delete m_stepLimiter;
474 if (m_chordFinder)
delete m_chordFinder;
475 if (m_stepper)
delete m_stepper;
476 if (m_magFldEquation)
delete m_magFldEquation;
477 if (m_uncachedField)
delete m_uncachedField;
478 if (m_magneticField)
delete m_magneticField;