Belle II Software  release-06-02-00
FullSimModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <simulation/modules/fullsim/FullSimModule.h>
10 #include <simulation/kernel/RunManager.h>
11 #include <simulation/kernel/DetectorConstruction.h>
12 //- #include <simulation/kernel/PhysicsList.h>
13 #include <simulation/physicslist/Belle2PhysicsList.h>
14 #include <simulation/kernel/ExtPhysicsConstructor.h>
15 #include <simulation/kernel/MagneticField.h>
16 #include <simulation/kernel/PrimaryGeneratorAction.h>
17 #include <simulation/kernel/EventAction.h>
18 #include <simulation/kernel/TrackingAction.h>
19 #include <simulation/kernel/SteppingAction.h>
20 #include <simulation/kernel/StackingAction.h>
21 
22 #include <mdst/dataobjects/MCParticle.h>
23 #include <framework/datastore/StoreObjPtr.h>
24 #include <framework/datastore/StoreArray.h>
25 #include <framework/dataobjects/EventMetaData.h>
26 #include <framework/gearbox/Unit.h>
27 
28 #include <CLHEP/Units/SystemOfUnits.h>
29 
30 #include <simulation/monopoles/G4MonopolePhysics.h>
31 #include <simulation/longlivedneutral/G4LongLivedNeutralPhysics.h>
32 
33 #include <G4TransportationManager.hh>
34 #include <G4Transportation.hh>
35 #include <G4PhysListFactory.hh>
36 #include <G4ProcessVector.hh>
37 #include <G4OpticalPhysics.hh>
38 #include <G4ParticleDefinition.hh>
39 #include <G4ParticleTable.hh>
40 #include <G4EventManager.hh>
41 #include <G4RunManager.hh>
42 #include <G4UImanager.hh>
43 #include <G4VisExecutive.hh>
44 #include <G4StepLimiter.hh>
45 #include <G4EmParameters.hh>
46 #include <G4HadronicProcessStore.hh>
47 #include <G4InuclParticleNames.hh>
48 
49 #include <G4Mag_UsualEqRhs.hh>
50 #include <G4NystromRK4.hh>
51 #include <G4HelixExplicitEuler.hh>
52 #include <G4HelixSimpleRunge.hh>
53 #include <G4CachedMagneticField.hh>
54 
55 using namespace std;
56 using namespace Belle2;
57 using namespace Belle2::Simulation;
58 using namespace Belle2::Monopoles;
59 using namespace G4InuclParticleNames;
60 
61 //-----------------------------------------------------------------
62 // Register the Module
63 //-----------------------------------------------------------------
64 REG_MODULE(FullSim)
65 
66 //-----------------------------------------------------------------
67 // Implementation
68 //-----------------------------------------------------------------
69 
70 FullSimModule::FullSimModule() : Module(), m_useNativeGeant4(true)
71 {
72  //Set module properties and the description
73  setDescription("Performs the full Geant4 detector simulation. Requires a valid geometry in memory.");
74  setPropertyFlags(c_ParallelProcessingCertified | c_TerminateInAllProcesses);
75 
76  //Parameter definition
77  addParam("InputMCParticleCollection", m_mcParticleInputColName, "The name of the input MCParticle collection.", string(""));
78  addParam("ThresholdImportantEnergy", m_thresholdImportantEnergy,
79  "[GeV] A particle which got 'stuck' and has less than this energy will be killed after 'ThresholdTrials' trials.", 0.250);
80  addParam("ThresholdTrials", m_thresholdTrials,
81  "Geant4 will try 'ThresholdTrials' times to move a particle which got 'stuck' and has an energy less than 'ThresholdImportantEnergy'.",
82  10);
83  addParam("RunEventVerbosity", m_runEventVerbosity, "Geant4 run/event verbosity: 0=silent; 1=info level; 2=debug level", 0);
84  addParam("TrackingVerbosity", m_trackingVerbosity,
85  "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.",
86  0);
87  addParam("HadronProcessVerbosity", m_hadronProcessVerbosity, "Hadron Process verbosity: 0=Silent; 1=info level; 2=debug level", 0);
88  addParam("EmProcessVerbosity", m_emProcessVerbosity, "Em Process verbosity: 0=Silent; 1=info level; 2=debug level", 0);
89  addParam("PhysicsList", m_physicsList, "The name of the physics list which is used for the simulation.", string("Belle2"));
90  addParam("StandardEM", m_standardEM, "If true, replaces fast EM physics with standard EM physics.", false);
91  addParam("RegisterOptics", m_optics, "If true, G4OpticalPhysics is registered in Geant4 PhysicsList.", true);
92  addParam("UseHighPrecisionNeutrons", m_HPneutrons, "If true, high precision neutron models used below 20 MeV.", false);
93  addParam("RegisterMonopoles", m_monopoles, "If set to true, G4MonopolePhysics is registered in Geant4 PhysicsList.", false);
94  addParam("MonopoleMagCharge", m_monopoleMagneticCharge, "The value of monopole magnetic charge in units of e+.", 1.0);
95  addParam("ProductionCut", m_productionCut,
96  "[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.",
97  0.07);
98  addParam("PXDProductionCut", m_pxdProductionCut, "[cm] Secondary production threshold in PXD envelope.", 0.0);
99  addParam("SVDProductionCut", m_svdProductionCut, "[cm] Secondary production threshold in SVD envelope.", 0.0);
100  addParam("CDCProductionCut", m_cdcProductionCut, "[cm] Secondary production threshold in CDC envelope.", 0.0);
101  addParam("ARICHTOPProductionCut", m_arichtopProductionCut, "[cm] Secondary production threshold in ARICH and TOP envelopes.", 0.0);
102  addParam("ECLProductionCut", m_eclProductionCut, "[cm] Secondary production threshold in ECL envelope.", 0.0);
103  addParam("KLMProductionCut", m_klmProductionCut, "[cm] Secondary production threshold in BKLM and EKLM envelopes.", 0.0);
104  addParam("MaxNumberSteps", m_maxNumberSteps,
105  "The maximum number of steps before the track transportation is stopped and the track is killed.", 100000);
106  addParam("PhotonFraction", m_photonFraction, "The fraction of Cerenkov photons which will be kept and propagated.", 0.5);
107  addParam("EnableVisualization", m_EnableVisualization, "If set to True, the Geant4 visualization support is enabled.", false);
108 
109  addParam("StoreOpticalPhotons", m_storeOpticalPhotons, "If set to True, optical photons are stored in MCParticles.", false);
110  addParam("StoreAllSecondaries", m_storeSecondaries,
111  "If set to True, all secondaries produced by Geant4 over a kinetic energy cut are stored in MCParticles. Otherwise do not store them.",
112  false);
113  addParam("SecondariesEnergyCut", m_secondariesEnergyCut, "[MeV] Kinetic energy cut for storing secondaries", 1.0);
114  addParam("StoreBremsstrahlungPhotons", m_storeBremsstrahlungPhotons,
115  "If set to True, store BremsstrahlungPhotons over a kinetic energy cut in MCParticles. Otherwise do not store them.", false);
116  addParam("BremsstrahlungPhotonsEnergyCut", m_bremsstrahlungPhotonsEnergyCut,
117  "[MeV] Kinetic energy cut for storing bremsstrahlung photons", 10.0);
118  addParam("StorePairConversions", m_storePairConversions,
119  "If set to True, store e+ or e- from pair conversions over a kinetic energy cut in MCParticles. Otherwise do not store them.",
120  false);
121  addParam("PairConversionsEnergyCut", m_pairConversionsEnergyCut,
122  "[MeV] Kinetic energy cut for storing e+ or e- from pair conversions", 10.0);
123 
124  addParam("magneticField", m_magneticFieldName,
125  "Chooses the magnetic field stepper used by Geant4. Possible values are: default, nystrom, expliciteuler, simplerunge",
126  string("default"));
127  addParam("magneticCacheDistance", m_magneticCacheDistance,
128  "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",
129  0.0);
130  addParam("deltaChordInMagneticField", m_deltaChordInMagneticField,
131  "[mm] The maximum miss-distance between the trajectory curve and its linear cord(s) approximation", 0.25);
132  vector<string> defaultCommandsAtPreInit;
133  addParam("UICommandsAtPreInit", m_uiCommandsAtPreInit,
134  "A list of Geant4 UI commands that should be applied at PreInit state, before the simulation starts.",
135  defaultCommandsAtPreInit);
136  vector<string> defaultCommandsAtIdle;
137  addParam("UICommandsAtIdle", m_uiCommandsAtIdle,
138  "A list of Geant4 UI commands that should be applied at Idle state, before the simulation starts.",
139  defaultCommandsAtIdle);
140  addParam("trajectoryStore", m_trajectoryStore,
141  "If non-zero save the full trajectory of 1=primary, 2=non-optical or 3=all particles", 0);
142  addParam("trajectoryDistanceTolerance", m_trajectoryDistanceTolerance,
143  "Maximum deviation from the real trajectory points when merging "
144  "segments (in cm)", 5e-4);
145  vector<float> defaultAbsorbers;
146  addParam("AbsorbersRadii", m_absorbers,
147  "Radii (in cm) of absorbers across which tracks will be destroyed.", defaultAbsorbers);
148 
149  //Make sure the instance of the run manager is created now to initialize some stuff we need for geometry
150  RunManager::Instance();
151  m_magneticField = NULL;
152  m_uncachedField = NULL;
153  m_magFldEquation = NULL;
154  m_stepper = NULL;
155  m_chordFinder = NULL;
156  m_visManager = NULL;
157  m_stepLimiter = NULL;
158 }
159 
160 
161 FullSimModule::~FullSimModule()
162 {
163 
164 }
165 
166 
167 void FullSimModule::initialize()
168 {
169  // MCParticles input and output collections can be different.
170  // Output collection is always the default one.
171  // In case we simulate only beam background events using BG mixing or BG overlay
172  // there is no input collection.
173 
174  if (m_mcParticleInputColName.empty()) {
175  // input and output collections are the same
176  // register in datastore because the input collection may not exist (case: only BG)
178  } else {
179  // input and output collections are different
180  StoreArray<MCParticle>().isRequired(m_mcParticleInputColName); // input collection
181  StoreArray<MCParticle>().registerInDataStore(); // output collection
182  }
183 
184  //Make sure the EventMetaData already exists.
186 
187  //Get the instance of the run manager.
188  RunManager& runManager = RunManager::Instance();
189 
190  //Add Geometry
191  runManager.SetUserInitialization(new DetectorConstruction());
192 
193  //Create the Physics list
194  //- PhysicsList* physicsList = new PhysicsList(m_physicsList);
195  //- physicsList->setProductionCutValue(m_productionCut);
196  //- if (m_optics) physicsList->registerOpticalPhysicsList();
197  //- runManager.SetUserInitialization(physicsList);
198 
199  if (m_physicsList == "Belle2") {
200  // Use Belle2PhysicsList
201  Belle2PhysicsList* physicsList = new Belle2PhysicsList(m_physicsList);
202  physicsList->SetVerbosity(m_runEventVerbosity);
203  physicsList->UseStandardEMPhysics(m_standardEM);
204  physicsList->UseOpticalPhysics(m_optics);
205  physicsList->UseHighPrecisionNeutrons(m_HPneutrons);
206  physicsList->SetProductionCutValue(m_productionCut);
207  physicsList->SetPXDProductionCutValue(m_pxdProductionCut);
208  physicsList->SetSVDProductionCutValue(m_svdProductionCut);
209  physicsList->SetCDCProductionCutValue(m_cdcProductionCut);
210  physicsList->SetARICHTOPProductionCutValue(m_arichtopProductionCut);
211  physicsList->SetECLProductionCutValue(m_eclProductionCut);
212  physicsList->SetKLMProductionCutValue(m_klmProductionCut);
213  physicsList->UseLongLivedNeutralParticles();
214 
215  //Apply the Geant4 UI commands in PreInit State - before initialization
216  if (m_uiCommandsAtPreInit.size() > 0) {
217  G4UImanager* uiManager = G4UImanager::GetUIpointer();
218  for (vector<string>::iterator iter = m_uiCommandsAtPreInit.begin(); iter != m_uiCommandsAtPreInit.end(); ++iter) {
219  uiManager->ApplyCommand(*iter);
220  }
221  }
222 
223  runManager.SetUserInitialization(physicsList);
224 
225  } else {
226  G4PhysListFactory physListFactory;
227  physListFactory.SetVerbose(m_runEventVerbosity);
228  G4VModularPhysicsList* physicsList = NULL;
229  if (physListFactory.IsReferencePhysList(m_physicsList)) physicsList = physListFactory.GetReferencePhysList(m_physicsList);
230  if (physicsList == NULL) B2FATAL("Could not load the physics list " << m_physicsList);
231  physicsList->RegisterPhysics(new ExtPhysicsConstructor);
232  if (m_optics) physicsList->RegisterPhysics(new G4OpticalPhysics);
233  if (m_monopoles) {
234  physicsList->RegisterPhysics(new G4MonopolePhysics(m_monopoleMagneticCharge));
235  }
236 
237  physicsList->RegisterPhysics(new G4LongLivedNeutralPhysics());
238 
239  physicsList->SetDefaultCutValue((m_productionCut / Unit::mm) * CLHEP::mm); // default is 0.7 mm
240 
241  //Apply the Geant4 UI commands in PreInit State - before initialization
242  if (m_uiCommandsAtPreInit.size() > 0) {
243  G4UImanager* uiManager = G4UImanager::GetUIpointer();
244  for (vector<string>::iterator iter = m_uiCommandsAtPreInit.begin(); iter != m_uiCommandsAtPreInit.end(); ++iter) {
245  uiManager->ApplyCommand(*iter);
246  }
247  }
248 
249  // LEP: For geant4e-specific particles, set a big step so that AlongStep computes
250  // all the energy (as is done in G4ErrorPhysicsList)
251  G4ParticleTable::G4PTblDicIterator* myParticleIterator = G4ParticleTable::GetParticleTable()->GetIterator();
252  myParticleIterator->reset();
253  while ((*myParticleIterator)()) {
254  G4ParticleDefinition* particle = myParticleIterator->value();
255  if (particle->GetParticleName().compare(0, 4, "g4e_") == 0) {
256  physicsList->SetParticleCuts(1.0E+9 * CLHEP::cm, particle);
257  }
258  }
259  runManager.SetUserInitialization(physicsList);
260  }
261 
262  //Create the magnetic field for the Geant4 simulation
263  if (m_magneticFieldName != "none") {
264  m_magneticField = new Belle2::Simulation::MagneticField();
265  if (m_magneticCacheDistance > 0) {
266  m_uncachedField = m_magneticField;
267  m_magneticField = new G4CachedMagneticField(m_uncachedField, m_magneticCacheDistance);
268  }
269  G4FieldManager* fieldManager = G4TransportationManager::GetTransportationManager()->GetFieldManager();
270  fieldManager->SetDetectorField(m_magneticField);
271  if (m_magneticFieldName != "default") {
272 
273  //We only use Magnetic field so let's try the specialized steppers
274  m_magFldEquation = new G4Mag_UsualEqRhs(m_magneticField);
275  if (m_magneticFieldName == "nystrom") {
276  m_stepper = new G4NystromRK4(m_magFldEquation);
277  } else if (m_magneticFieldName == "expliciteuler") {
278  m_stepper = new G4HelixExplicitEuler(m_magFldEquation);
279  } else if (m_magneticFieldName == "simplerunge") {
280  m_stepper = new G4HelixSimpleRunge(m_magFldEquation);
281  } else {
282  B2FATAL("Unknown magnetic field option: " << m_magneticFieldName);
283  }
284 
285  //Set a minimum stepsize (stepMinimum): The chordfinder should not attempt to limit
286  //the stepsize to something less than 10µm (which is the default value of Geant4).
287  m_chordFinder = new G4ChordFinder(m_magneticField, 1e-2 * CLHEP::mm, m_stepper);
288  fieldManager->SetChordFinder(m_chordFinder);
289  } else {
290  fieldManager->CreateChordFinder(m_magneticField);
291  }
292 
293  //Change DeltaCord (the max. miss-distance between the trajectory curve and its linear chord(s) approximation, if asked.
294  G4ChordFinder* chordFinder = fieldManager->GetChordFinder();
295  B2DEBUG(1, "Geant4 default deltaChord = " << chordFinder->GetDeltaChord());
296  chordFinder->SetDeltaChord(m_deltaChordInMagneticField * CLHEP::mm);
297  B2DEBUG(1, "DeltaChord after reset = " << chordFinder->GetDeltaChord());
298 
299  //This might be a good place to optimize the Integration parameters (DeltaOneStep, DeltaIntersection, MinEpsilon, MaxEpsilon)
300  }
301 
302  //Create the generator action which takes the MCParticle list and converts it to Geant4 primary vertices.
303  G4VUserPrimaryGeneratorAction* generatorAction = new PrimaryGeneratorAction(m_mcParticleInputColName, m_mcParticleGraph);
304  runManager.SetUserAction(generatorAction);
305 
306  //Add the event action which creates the final MCParticle list and the Relation list.
307  //The output collection name will be always "MCParticles".
308  EventAction* eventAction = new EventAction("", m_mcParticleGraph);
309  runManager.SetUserAction(eventAction);
310 
311  //Add the tracking action which handles the secondary particles created by Geant4.
312  TrackingAction* trackingAction = new TrackingAction(m_mcParticleGraph);
313  trackingAction->setIgnoreOpticalPhotons(!m_storeOpticalPhotons);
314  trackingAction->setIgnoreSecondaries(!m_storeSecondaries);
315  trackingAction->setSecondariesEnergyCut(m_secondariesEnergyCut);
316  trackingAction->setIgnoreBremsstrahlungPhotons(!m_storeBremsstrahlungPhotons);
317  trackingAction->setBremsstrahlungPhotonsEnergyCut(m_bremsstrahlungPhotonsEnergyCut);
318  trackingAction->setIgnorePairConversions(!m_storePairConversions);
319  trackingAction->setPairConversionsEnergyCut(m_pairConversionsEnergyCut);
320 
321  runManager.SetUserAction(trackingAction);
322 
323  //Add the stepping action which provides additional security checks
324  SteppingAction* steppingAction = new SteppingAction();
325  steppingAction->setMaxNumberSteps(m_maxNumberSteps);
326  steppingAction->setAbsorbersR(m_absorbers);
327  for (auto& rAbsorber : m_absorbers) {
328  B2INFO("An absorber found at R = " << rAbsorber << " cm");
329  }
330  runManager.SetUserAction(steppingAction);
331 
332  //Add the stacking action which provides performance speed ups for the handling of optical photons
333  StackingAction* stackingAction = new StackingAction();
334  stackingAction->setPropagatedPhotonFraction(m_photonFraction);
335  runManager.SetUserAction(stackingAction);
336 
337  //Initialize G4 kernel
338  runManager.Initialize();
339 
340  //Set the parameters for the G4Transportation system.
341  //To make sure we really change all G4Transportation classes, we loop over all particles
342  //even if the pointer to the G4Transportation object seems to be the same for all particles.
343  //Only one instance of G4StepLimiter is needed: see G4StepLimiterBuilder(), for example.
344  m_stepLimiter = new G4StepLimiter();
345  G4ParticleTable::G4PTblDicIterator* partIter = G4ParticleTable::GetParticleTable()->GetIterator();
346  partIter->reset();
347  while ((*partIter)()) {
348  G4ParticleDefinition* currParticle = partIter->value();
349  G4ProcessVector& currProcList = *currParticle->GetProcessManager()->GetProcessList();
350  assert(currProcList.size() < INT_MAX);
351  for (int iProcess = 0; iProcess < static_cast<int>(currProcList.size()); ++iProcess) {
352  G4Transportation* transport = dynamic_cast<G4Transportation*>(currProcList[iProcess]);
353  if (transport != nullptr) {
354  //Geant4 energy unit is MeV
355  transport->SetThresholdImportantEnergy(m_thresholdImportantEnergy / Unit::MeV * CLHEP::MeV);
356  transport->SetThresholdTrials(m_thresholdTrials);
357  break;
358  }
359  }
360  // Add StepLimiter process for charged tracks.
361  double zeroChargeTol = 0.01 * Unit::e;
362  if (fabs(currParticle->GetPDGCharge()) > zeroChargeTol) {
363  currParticle->GetProcessManager()->AddDiscreteProcess(m_stepLimiter);
364  B2DEBUG(100, "Added StepLimiter process for " << currParticle->GetParticleName());
365  }
366  }
367 
368  // Inactivate all secondary-generating processes for g4e particles. This comprises
369  // Cerenkov and Scintillation that were inserted by G4OpticalPhysics and the
370  // CaptureAtRest process for g4e anti-deuteron.
371  partIter->reset();
372  while ((*partIter)()) {
373  G4ParticleDefinition* currParticle = partIter->value();
374  if (currParticle->GetParticleName().compare(0, 4, "g4e_") == 0) {
375  G4ProcessManager* processManager = currParticle->GetProcessManager();
376  if (processManager) {
377  G4ProcessVector* processList = processManager->GetProcessList();
378  assert(processList->size() < INT_MAX);
379  for (int i = 0; i < static_cast<int>(processList->size()); ++i) {
380  if (((*processList)[i]->GetProcessName() == "Cerenkov") ||
381  ((*processList)[i]->GetProcessName() == "Scintillation") ||
382  ((*processList)[i]->GetProcessName() == "hFritiofCaptureAtRest")) {
383  processManager->SetProcessActivation(i, false);
384  }
385  }
386  }
387  }
388  }
389 
390  //Set the verbosity level of Geant4 according to the logging settings of the module
391  //int g4VerboseLevel = 0;
392  //switch (LogSystem::Instance().getCurrentLogLevel()) {
393  // case LogConfig::c_Debug : g4VerboseLevel = 2;
394  // break;
395  // case LogConfig::c_Info : g4VerboseLevel = 1;
396  // break;
397  // default: g4VerboseLevel = 0;
398  //}
399  //G4EventManager::GetEventManager()->SetVerboseLevel(g4VerboseLevel);
400  //G4RunManager::GetRunManager()->SetVerboseLevel(g4VerboseLevel);
401  G4EventManager::GetEventManager()->SetVerboseLevel(m_runEventVerbosity);
402  G4RunManager::GetRunManager()->SetVerboseLevel(m_runEventVerbosity);
403  G4EventManager::GetEventManager()->GetTrackingManager()->SetVerboseLevel(
404  m_trackingVerbosity); //turned out to be more useful as a parameter.
405  G4HadronicProcessStore::Instance()->SetVerbose(m_hadronProcessVerbosity);
406  G4EmParameters::Instance()->SetVerbose(m_emProcessVerbosity);
407 
408 
409  if (m_EnableVisualization) {
410  m_visManager = new G4VisExecutive;
411  m_visManager->Initialize();
412  }
413 
414  //Apply the Geant4 UI commands at Idle state - after initilization
415  if (m_uiCommandsAtIdle.size() > 0) {
416  G4UImanager* uiManager = G4UImanager::GetUIpointer();
417  for (vector<string>::iterator iter = m_uiCommandsAtIdle.begin(); iter != m_uiCommandsAtIdle.end(); ++iter) {
418  uiManager->ApplyCommand(*iter);
419  }
420  }
421 
422  //Store Trajectories?
423  if (m_trajectoryStore) {
424  trackingAction->setStoreTrajectories(m_trajectoryStore, m_trajectoryDistanceTolerance);
425  steppingAction->setStoreTrajectories(true);
426  }
427 
428  //Physics tables are build in run initialization. We have run independent
429  //geometry at the moment so there is no need to do this in begin run. Instead
430  //we use one Geant4 run for all Belle2 runs we might encounter. So let's do
431  //run initialization now to save memory when doing parallel processing
432  B2INFO("Perform Geant4 final initialization: Geometry optimization, PhysicsList calculations...");
433  RunManager::Instance().beginRun(0);
434  B2INFO("done, Geant4 ready");
435  //Otherwise we could use a fake run to do this and move RunManager::beginRun
436  //back to beginRun()
437  //runManager.BeamOn(0);
438 }
439 
440 
441 void FullSimModule::beginRun()
442 {
443  //Nothing to do: geometry and physics are run independent
444 }
445 
446 
447 void FullSimModule::event()
448 {
449  //Get the event meta data
450  StoreObjPtr<EventMetaData> eventMetaDataPtr;
451 
452  //Process the event
453  RunManager::Instance().processEvent(eventMetaDataPtr->getEvent());
454 }
455 
456 
457 void FullSimModule::endRun()
458 {
459  //Nothing to do: geometry and physics are run independent
460 }
461 
462 void FullSimModule::terminate()
463 {
464  //We used one Geant4 run for all Belle2 runs so end the geant4 run here
465  RunManager::Instance().endRun();
466  //And clean up the run manager
467  if (m_visManager != nullptr) delete m_visManager;
468  RunManager::Instance().destroy();
469  // Delete the step limiter process
470  delete m_stepLimiter;
471  // Delete the objects associated with transport in magnetic field
472  if (m_chordFinder) delete m_chordFinder;
473  if (m_stepper) delete m_stepper;
474  if (m_magFldEquation) delete m_magFldEquation;
475  if (m_uncachedField) delete m_uncachedField;
476  if (m_magneticField) delete m_magneticField;
477 }
Class responsible to connect to geometry to simulation.
The full Geant4 simulation module.
Definition: FullSimModule.h:42
LongLivedNeutral physics Class – to be registered in the physics list.
Base class for Modules.
Definition: Module.h:72
Monopole physics class to register on the physics list.
Custom Geant4 physics list for Belle II with options to add optical physics, standard EM physics and ...
void SetSVDProductionCutValue(G4double)
Set cut value for SVD envelope.
void SetProductionCutValue(G4double)
Use parameter to set global cut value.
void SetARICHTOPProductionCutValue(G4double)
Set cut value for ARICH and TOP envelopes.
void SetECLProductionCutValue(G4double)
Set cut value for ECL barrel, forward and backward envelopes.
void UseHighPrecisionNeutrons(G4bool)
Use high precision neutron models below 20 MeV.
void SetPXDProductionCutValue(G4double)
Set cut value for PXD envelope.
void UseLongLivedNeutralParticles()
Simulate neutral long-lived particles with given pdg and mass value.
void UseOpticalPhysics(G4bool)
Add optical photon physics.
void SetVerbosity(G4int verb)
Run/event verbosity level.
void SetCDCProductionCutValue(G4double)
Set cut value for CDC envelope.
void SetKLMProductionCutValue(G4double)
Set cut value for BKLM and EKLM envelopes.
void UseStandardEMPhysics(G4bool)
Use standard EM physics instead of EM option1.
The Event Action class.
Definition: EventAction.h:31
The Class for the Belle2 magnetic field implementation for Geant4.
Definition: MagneticField.h:28
The PrimaryGeneratorAction class inherits from G4VuserPrimaryGeneratorAction and specifies how a prim...
The run manager controls the flow of the Geant4 program and manages the event loop(s) within a run.
Definition: RunManager.h:30
void Initialize()
Initialize the Kernel.
Definition: RunManager.cc:35
The basf2 stacking action.
void setPropagatedPhotonFraction(double fraction)
Set fraction of Cerenkov photons that are actually propagated.
The Class for the stepping action.
void setAbsorbersR(const std::vector< float > &vec)
Sets the radii of absorbers for killing tracks across them.
void setStoreTrajectories(bool store)
Sets the trajectory option to enable storing of the simulated particle trajectories.
void setMaxNumberSteps(int maxSteps)
Sets the maximum number of steps before a track is stopped and killed.
The Tracking Action class.
void setIgnorePairConversions(bool ignore=true)
Set ignore flag for e+ or e- coming from gamma conversions into a pair if set to true,...
void setIgnoreBremsstrahlungPhotons(bool ignore=true)
Set ignore flag for low energy breamsstrahlung photons if set to true, breamsstrahlung photons with k...
void setSecondariesEnergyCut(double cut_MeV)
Set kinetic energy cut for secondaries.
void setStoreTrajectories(int store, double distanceTolerance)
Sets the trajectory option to enable storing of the simulated particle trajectories.
void setBremsstrahlungPhotonsEnergyCut(double cut_MeV)
Set kinetic energy cut for bremsstrahlung photons.
void setPairConversionsEnergyCut(double cut_MeV)
Set kinetic energy cut for e+ e- pair conversions.
void setIgnoreSecondaries(bool ignore=true)
Set ignore flag for low energy Geant-produced secondary particles if set to true, secondaries with ki...
void setIgnoreOpticalPhotons(bool ignore=true)
Set ignore flag for optical photons if set to true, optical photons will not be stored in MCParticles...
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.