Belle II Software development
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 <G4FieldManager.hh>
34#include <G4TransportationManager.hh>
35#include <G4Transportation.hh>
36#include <G4PhysListFactory.hh>
37#include <G4ProcessVector.hh>
38#include <G4OpticalPhysics.hh>
39#include <G4ParticleDefinition.hh>
40#include <G4ParticleTable.hh>
41#include <G4EventManager.hh>
42#include <G4RunManager.hh>
43#include <G4UImanager.hh>
44#include <G4VisExecutive.hh>
45#include <G4StepLimiter.hh>
46#include <G4EmParameters.hh>
47#include <G4HadronicProcessStore.hh>
48#include <G4InuclParticleNames.hh>
49
50#include <G4Mag_UsualEqRhs.hh>
51#include <G4NystromRK4.hh>
52#include <G4HelixExplicitEuler.hh>
53#include <G4HelixSimpleRunge.hh>
54#include <G4CachedMagneticField.hh>
55#include <G4ChordFinder.hh>
56
57using namespace std;
58using namespace Belle2;
59using namespace Belle2::Simulation;
60using namespace Belle2::Monopoles;
61using namespace G4InuclParticleNames;
62
63//-----------------------------------------------------------------
64// Register the Module
65//-----------------------------------------------------------------
66REG_MODULE(FullSim);
67
68//-----------------------------------------------------------------
69// Implementation
70//-----------------------------------------------------------------
71
73{
74 //Set module properties and the description
75 setDescription("Performs the full Geant4 detector simulation. Requires a valid geometry in memory.");
77
78 //Parameter definition
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'.",
84 10);
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.",
88 0);
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.",
99 0.07);
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.02);
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);
110
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.",
114 false);
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.",
122 false);
123 addParam("PairConversionsEnergyCut", m_pairConversionsEnergyCut,
124 "[MeV] Kinetic energy cut for storing e+ or e- from pair conversions", 10.0);
125
126 addParam("DetailedParticleMatching", m_useDetailedParticleMatching,
127 "If true, use enhanced logic to store specific secondaries that would otherwise be ignored. "
128 "This allows to rescue particles based on their production vertex position relative to the defined region, "
129 "kinetic energy, and other properties.", false);
130 addParam("RegionZBackward", m_regionZBackward,
131 "The backward z limit [cm] of the region, used to define the region volume for particle matching.", 0.0);
132 addParam("RegionZForward", m_regionZForward,
133 "The forward z limit [cm] of the region, used to define the region volume for particle matching.", 0.0);
134 addParam("RegionRho", m_regionRho, "The rho limit [cm] of the region, used to define the region volume for particle matching.",
135 0.0);
136 addParam("KineticEnergyThreshold", m_kineticEnergyThreshold,
137 "Kinetic energy threshold [MeV]. Only particles above this energy are considered for detailed matching.", 0.0);
138 addParam("DistanceThreshold", m_distanceThreshold,
139 "Distance threshold [cm] for particles produced inside the defined region volume. "
140 "Only particles that travel further than this distance are considered.", 0.0);
141 addParam("DoNotStoreEMParticles", m_doNotStoreEMParticles,
142 "If true, electromagnetic particles (e+, e-, gamma) produced inside the defined region and passing the distance threshold "
143 "are NOT stored (they are filtered out). If false, EM particles can be stored.", false);
144 addParam("DoNotStoreNuclei", m_doNotStoreNuclei,
145 "If true, nuclei (PDG > 10^9) produced inside the defined region and passing the distance threshold "
146 "are NOT stored (they are filtered out). If false, nuclei can be stored.", false);
147 addParam("UseSeenInECL", m_useSeenInECL,
148 "If true, particles produced before the defined region (inverse of inside region) must have been 'seen' in the ECL (created hits) to be stored. "
149 "If false, all eligible particles produced before the defined region are stored.", false);
150
151 addParam("magneticField", m_magneticFieldName,
152 "Chooses the magnetic field stepper used by Geant4. Possible values are: default, nystrom, expliciteuler, simplerunge",
153 string("default"));
154 addParam("magneticCacheDistance", m_magneticCacheDistance,
155 "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",
156 0.0);
157 addParam("deltaChordInMagneticField", m_deltaChordInMagneticField,
158 "[mm] The maximum miss-distance between the trajectory curve and its linear cord(s) approximation", 0.25);
159 vector<string> defaultCommandsAtPreInit;
160 addParam("UICommandsAtPreInit", m_uiCommandsAtPreInit,
161 "A list of Geant4 UI commands that should be applied at PreInit state, before the simulation starts.",
162 defaultCommandsAtPreInit);
163 vector<string> defaultCommandsAtIdle;
164 addParam("UICommandsAtIdle", m_uiCommandsAtIdle,
165 "A list of Geant4 UI commands that should be applied at Idle state, before the simulation starts.",
166 defaultCommandsAtIdle);
167 addParam("trajectoryStore", m_trajectoryStore,
168 "If non-zero save the full trajectory of 1=primary, 2=non-optical or 3=all particles", 0);
169 addParam("trajectoryDistanceTolerance", m_trajectoryDistanceTolerance,
170 "Maximum deviation from the real trajectory points when merging "
171 "segments (in cm)", 5e-4);
172 vector<float> defaultAbsorbers;
173 addParam("AbsorbersRadii", m_absorbers,
174 "Radii (in cm) of absorbers across which tracks will be destroyed.", defaultAbsorbers);
175
176 //Make sure the instance of the run manager is created now to initialize some stuff we need for geometry
178 m_magneticField = NULL;
179 m_uncachedField = NULL;
180 m_magFldEquation = NULL;
181 m_stepper = NULL;
182 m_chordFinder = NULL;
183 m_visManager = NULL;
184 m_stepLimiter = NULL;
185}
186
187
192
193
195{
196 // MCParticles input and output collections can be different.
197 // Output collection is always the default one.
198 // In case we simulate only beam background events using BG mixing or BG overlay
199 // there is no input collection.
200
201 if (m_mcParticleInputColName.empty()) {
202 // input and output collections are the same
203 // register in datastore because the input collection may not exist (case: only BG)
205 } else {
206 // input and output collections are different
208 StoreArray<MCParticle>().registerInDataStore(); // output collection
209 }
210
211 //Make sure the EventMetaData already exists.
213
214 //Get the instance of the run manager.
215 RunManager& runManager = RunManager::Instance();
216
217 //Add Geometry
218 runManager.SetUserInitialization(new DetectorConstruction());
219
220 //Create the Physics list
221 //- PhysicsList* physicsList = new PhysicsList(m_physicsList);
222 //- physicsList->setProductionCutValue(m_productionCut);
223 //- if (m_optics) physicsList->registerOpticalPhysicsList();
224 //- runManager.SetUserInitialization(physicsList);
225
226 if (m_physicsList == "Belle2") {
227 // Use Belle2PhysicsList
229 physicsList->SetVerbosity(m_runEventVerbosity);
231 physicsList->UseOpticalPhysics(m_optics);
240 physicsList->UseLongLivedNeutralParticles();
241
242 //Apply the Geant4 UI commands in PreInit State - before initialization
243 if (m_uiCommandsAtPreInit.size() > 0) {
244 G4UImanager* uiManager = G4UImanager::GetUIpointer();
245 for (vector<string>::iterator iter = m_uiCommandsAtPreInit.begin(); iter != m_uiCommandsAtPreInit.end(); ++iter) {
246 uiManager->ApplyCommand(*iter);
247 }
248 }
249
250 runManager.SetUserInitialization(physicsList);
251
252 } else {
253 G4PhysListFactory physListFactory;
254 physListFactory.SetVerbose(m_runEventVerbosity);
255 G4VModularPhysicsList* physicsList = NULL;
256 if (physListFactory.IsReferencePhysList(m_physicsList)) physicsList = physListFactory.GetReferencePhysList(m_physicsList);
257 if (physicsList == NULL) B2FATAL("Could not load the physics list " << m_physicsList);
258 physicsList->RegisterPhysics(new ExtPhysicsConstructor);
259 if (m_optics) physicsList->RegisterPhysics(new G4OpticalPhysics);
260 if (m_monopoles) {
261 physicsList->RegisterPhysics(new G4MonopolePhysics(m_monopoleMagneticCharge));
262 }
263
264 physicsList->RegisterPhysics(new G4LongLivedNeutralPhysics());
265
266 physicsList->SetDefaultCutValue((m_productionCut / Unit::mm) * CLHEP::mm); // default is 0.7 mm
267
268 //Apply the Geant4 UI commands in PreInit State - before initialization
269 if (m_uiCommandsAtPreInit.size() > 0) {
270 G4UImanager* uiManager = G4UImanager::GetUIpointer();
271 for (vector<string>::iterator iter = m_uiCommandsAtPreInit.begin(); iter != m_uiCommandsAtPreInit.end(); ++iter) {
272 uiManager->ApplyCommand(*iter);
273 }
274 }
275
276 // // LEP: For geant4e-specific particles, set a big step so that AlongStep computes
277 // // all the energy (as is done in G4ErrorPhysicsList)
278 // G4ParticleTable::G4PTblDicIterator* myParticleIterator = G4ParticleTable::GetParticleTable()->GetIterator();
279 // myParticleIterator->reset();
280 // while ((*myParticleIterator)()) {
281 // G4ParticleDefinition* particle = myParticleIterator->value();
282 // if (particle->GetParticleName().compare(0, 4, "g4e_") == 0) {
283 // physicsList->SetParticleCuts(1.0E+9 * CLHEP::cm, particle);
284 // }
285 // }
286 runManager.SetUserInitialization(physicsList);
287 }
288
289 //Create the magnetic field for the Geant4 simulation
290 if (m_magneticFieldName != "none") {
292 if (m_magneticCacheDistance > 0) {
294 m_magneticField = new G4CachedMagneticField(m_uncachedField, m_magneticCacheDistance);
295 }
296 G4FieldManager* fieldManager = G4TransportationManager::GetTransportationManager()->GetFieldManager();
297 fieldManager->SetDetectorField(m_magneticField);
298 if (m_magneticFieldName != "default") {
299
300 //We only use Magnetic field so let's try the specialized steppers
301 m_magFldEquation = new G4Mag_UsualEqRhs(m_magneticField);
302 if (m_magneticFieldName == "nystrom") {
303 m_stepper = new G4NystromRK4(m_magFldEquation);
304 } else if (m_magneticFieldName == "expliciteuler") {
305 m_stepper = new G4HelixExplicitEuler(m_magFldEquation);
306 } else if (m_magneticFieldName == "simplerunge") {
307 m_stepper = new G4HelixSimpleRunge(m_magFldEquation);
308 } else {
309 B2FATAL("Unknown magnetic field option: " << m_magneticFieldName);
310 }
311
312 //Set a minimum stepsize (stepMinimum): The chordfinder should not attempt to limit
313 //the stepsize to something less than 10µm (which is the default value of Geant4).
314 m_chordFinder = new G4ChordFinder(m_magneticField, 1e-2 * CLHEP::mm, m_stepper);
315 fieldManager->SetChordFinder(m_chordFinder);
316 } else {
317 fieldManager->CreateChordFinder(m_magneticField);
318 }
319
320 //Change DeltaCord (the max. miss-distance between the trajectory curve and its linear chord(s) approximation, if asked.
321 G4ChordFinder* chordFinder = fieldManager->GetChordFinder();
322 B2DEBUG(1, "Geant4 default deltaChord = " << chordFinder->GetDeltaChord());
323 chordFinder->SetDeltaChord(m_deltaChordInMagneticField * CLHEP::mm);
324 B2DEBUG(1, "DeltaChord after reset = " << chordFinder->GetDeltaChord());
325
326 //This might be a good place to optimize the Integration parameters (DeltaOneStep, DeltaIntersection, MinEpsilon, MaxEpsilon)
327 }
328
329 //Create the generator action which takes the MCParticle list and converts it to Geant4 primary vertices.
330 G4VUserPrimaryGeneratorAction* generatorAction = new PrimaryGeneratorAction(m_mcParticleInputColName, m_mcParticleGraph);
331 runManager.SetUserAction(generatorAction);
332
333 //Add the event action which creates the final MCParticle list and the Relation list.
334 //The output collection name will be always "MCParticles".
335 EventAction* eventAction = new EventAction("", m_mcParticleGraph);
336 runManager.SetUserAction(eventAction);
337
338 //Add the tracking action which handles the secondary particles created by Geant4.
339 TrackingAction* trackingAction = new TrackingAction(m_mcParticleGraph);
348 trackingAction->setRegionZBackward(m_regionZBackward);
349 trackingAction->setRegionZForward(m_regionZForward);
350 trackingAction->setRegionRho(m_regionRho);
355 trackingAction->setUseSeenInECL(m_useSeenInECL);
356
357 runManager.SetUserAction(trackingAction);
358
359 //Add the stepping action which provides additional security checks
360 SteppingAction* steppingAction = new SteppingAction();
361 steppingAction->setMaxNumberSteps(m_maxNumberSteps);
362 steppingAction->setAbsorbersR(m_absorbers);
363 for (auto& rAbsorber : m_absorbers) {
364 B2INFO("An absorber found at R = " << rAbsorber << " cm");
365 }
366 runManager.SetUserAction(steppingAction);
367
368 //Add the stacking action which provides performance speed ups for the handling of optical photons
369 StackingAction* stackingAction = new StackingAction();
371 runManager.SetUserAction(stackingAction);
372
373 //Initialize G4 kernel
374 runManager.Initialize();
375
376 //Set the parameters for the G4Transportation system.
377 //To make sure we really change all G4Transportation classes, we loop over all particles
378 //even if the pointer to the G4Transportation object seems to be the same for all particles.
379 //Only one instance of G4StepLimiter is needed: see G4StepLimiterBuilder(), for example.
380 m_stepLimiter = new G4StepLimiter();
381 G4ParticleTable::G4PTblDicIterator* partIter = G4ParticleTable::GetParticleTable()->GetIterator();
382 partIter->reset();
383 while ((*partIter)()) {
384 G4ParticleDefinition* currParticle = partIter->value();
385 G4ProcessVector& currProcList = *currParticle->GetProcessManager()->GetProcessList();
386 assert(currProcList.size() < INT_MAX);
387 for (int iProcess = 0; iProcess < static_cast<int>(currProcList.size()); ++iProcess) {
388 G4Transportation* transport = dynamic_cast<G4Transportation*>(currProcList[iProcess]);
389 if (transport != nullptr) {
390 //Geant4 energy unit is MeV
391 transport->SetThresholdImportantEnergy(m_thresholdImportantEnergy / Unit::MeV * CLHEP::MeV);
392 transport->SetThresholdTrials(m_thresholdTrials);
393 break;
394 }
395 }
396 // Add StepLimiter process for charged tracks.
397 double zeroChargeTol = 0.01 * Unit::e;
398 if (fabs(currParticle->GetPDGCharge()) > zeroChargeTol) {
399 currParticle->GetProcessManager()->AddDiscreteProcess(m_stepLimiter);
400 B2DEBUG(100, "Added StepLimiter process for " << currParticle->GetParticleName());
401 }
402 }
403
404 // Deactivate all secondary-generating processes for g4e particles. This comprises
405 // Cerenkov and Scintillation that were inserted by G4OpticalPhysics and the
406 // CaptureAtRest process for g4e anti-deuteron.
407 partIter->reset();
408 while ((*partIter)()) {
409 G4ParticleDefinition* currParticle = partIter->value();
410 if (currParticle->GetParticleName().compare(0, 4, "g4e_") == 0) {
411 G4ProcessManager* processManager = currParticle->GetProcessManager();
412 if (processManager) {
413 G4ProcessVector* processList = processManager->GetProcessList();
414 assert(processList->size() < INT_MAX);
415 for (int i = 0; i < static_cast<int>(processList->size()); ++i) {
416 if (((*processList)[i]->GetProcessName() == "Cerenkov") ||
417 ((*processList)[i]->GetProcessName() == "Scintillation") ||
418 ((*processList)[i]->GetProcessName() == "hFritiofCaptureAtRest")) {
419 processManager->SetProcessActivation(i, false);
420 }
421 }
422 }
423 }
424 }
425
426 //Set the verbosity level of Geant4 according to the logging settings of the module
427 //int g4VerboseLevel = 0;
428 //switch (LogSystem::Instance().getCurrentLogLevel()) {
429 // case LogConfig::c_Debug : g4VerboseLevel = 2;
430 // break;
431 // case LogConfig::c_Info : g4VerboseLevel = 1;
432 // break;
433 // default: g4VerboseLevel = 0;
434 //}
435 //G4EventManager::GetEventManager()->SetVerboseLevel(g4VerboseLevel);
436 //G4RunManager::GetRunManager()->SetVerboseLevel(g4VerboseLevel);
437 G4EventManager::GetEventManager()->SetVerboseLevel(m_runEventVerbosity);
438 G4RunManager::GetRunManager()->SetVerboseLevel(m_runEventVerbosity);
439 G4EventManager::GetEventManager()->GetTrackingManager()->SetVerboseLevel(
440 m_trackingVerbosity); //turned out to be more useful as a parameter.
441 G4HadronicProcessStore::Instance()->SetVerbose(m_hadronProcessVerbosity);
442 G4EmParameters::Instance()->SetVerbose(m_emProcessVerbosity);
443
444
446 m_visManager = new G4VisExecutive;
447 m_visManager->Initialize();
448 }
449
450 //Apply the Geant4 UI commands at Idle state - after initilization
451 if (m_uiCommandsAtIdle.size() > 0) {
452 G4UImanager* uiManager = G4UImanager::GetUIpointer();
453 for (vector<string>::iterator iter = m_uiCommandsAtIdle.begin(); iter != m_uiCommandsAtIdle.end(); ++iter) {
454 uiManager->ApplyCommand(*iter);
455 }
456 }
457
458 //Store Trajectories?
459 if (m_trajectoryStore) {
461 steppingAction->setStoreTrajectories(true);
462 }
463
464 //Physics tables are build in run initialization. We have run independent
465 //geometry at the moment so there is no need to do this in begin run. Instead
466 //we use one Geant4 run for all Belle2 runs we might encounter. So let's do
467 //run initialization now to save memory when doing parallel processing
468 B2INFO("Perform Geant4 final initialization: Geometry optimization, PhysicsList calculations...");
470 B2INFO("done, Geant4 ready");
471 //Otherwise we could use a fake run to do this and move RunManager::beginRun
472 //back to beginRun()
473 //runManager.BeamOn(0);
474}
475
476
478{
479 //Nothing to do: geometry and physics are run independent
480}
481
482
484{
485 //Get the event meta data
486 StoreObjPtr<EventMetaData> eventMetaDataPtr;
487
488 //Process the event
489 RunManager::Instance().processEvent(eventMetaDataPtr->getEvent());
490}
491
492
494{
495 //Nothing to do: geometry and physics are run independent
496}
497
499{
500 //We used one Geant4 run for all Belle2 runs so end the geant4 run here
502 //And clean up the run manager
503 if (m_visManager != nullptr) delete m_visManager;
505 // Delete the step limiter process
506 delete m_stepLimiter;
507 // Delete the objects associated with transport in magnetic field
508 if (m_chordFinder) delete m_chordFinder;
509 if (m_stepper) delete m_stepper;
513}
Class responsible to connect to geometry to simulation.
std::vector< float > m_absorbers
The absorbers defined at given radii where tracks across them will be destroyed.
int m_emProcessVerbosity
Loss Table verbosity: 0=Silent; 1=info level; 2=debug level, default=0.
double m_trajectoryDistanceTolerance
Maximum distance to actual trajectory when merging points.
G4VisManager * m_visManager
Pointer to the visualization manager (if used)
std::vector< std::string > m_uiCommandsAtIdle
A list of Geant4 UI commands that should be applied at Idle state, after the Geant4 initialization an...
double m_monopoleMagneticCharge
The value of monopole magnetic charge in units of e+.
double m_regionZBackward
Region backward z limit for filtering secondaries.
double m_cdcProductionCut
Secondary production threshold in CDC envelope.
double m_photonFraction
The fraction of Cerenkov photons which will be kept and propagated.
bool m_doNotStoreNuclei
use is Nuclei check for filtering secondaries
bool m_HPneutrons
If true, high precision neutron models used below 20 MeV.
int m_hadronProcessVerbosity
Hadron Process verbosity: 0=Silent; 1=info level; 2=debug level, default=0.
G4Mag_UsualEqRhs * m_magFldEquation
Pointer to the equation of motion in the magnetic field (if not the default)
double m_klmProductionCut
Secondary production threshold in BKLM and EKLM envelopes.
virtual void initialize() override
Initialize the Module.
int m_thresholdTrials
Geant4 will try m_thresholdTrials times to move a particle which got 'stuck' and has an energy less t...
bool m_storePairConversions
controls storing of e+ or e- from pair conversions in MCParticles
double m_pxdProductionCut
Secondary production threshold in PXD envelope.
bool m_doNotStoreEMParticles
use is EM check for filtering secondaries
int m_trackingVerbosity
Tracking verbosity: 0=Silent; 1=Min info per step; 2=sec particles; 3=pre/post step info; 4=like 3 bu...
virtual void event() override
Performs the full Geant4 simulation.
double m_productionCut
Apply continuous energy loss to primary particle which has no longer enough energy to produce seconda...
bool m_storeOpticalPhotons
controls storing of optical photons in MCParticles
bool m_standardEM
If set to true, replaces fast EM physics with standard EM physics.
G4MagIntegratorStepper * m_stepper
Pointer to the equation-of-motion stepper (if not the default)
bool m_useNativeGeant4
If set to true, uses the Geant4 navigator and native detector construction class.
virtual void endRun() override
Called when run has ended.
double m_svdProductionCut
Secondary production threshold in SVD envelope.
virtual void terminate() override
Terminates the module.
double m_arichtopProductionCut
Secondary production threshold in ARICH and TOP envelopes.
double m_distanceThreshold
distance threshold for filtering secondaries
bool m_optics
If set to true, registers the optical physics list.
double m_pairConversionsEnergyCut
kinetic energy cut for the stored e+ or e- from pair conversions
std::vector< std::string > m_uiCommandsAtPreInit
A list of Geant4 UI commands that should be applied at PreInit state, before the Geant4 initializatio...
bool m_storeBremsstrahlungPhotons
controls storing of bremsstrahlung photons in MCParticles
bool m_useDetailedParticleMatching
If true, secondaries are ignored unless they pass additional checks (e.g.
std::string m_mcParticleInputColName
The parameter variable for the name of the input MCParticle collection.
virtual void beginRun() override
Called when a new run is started.
bool m_useSeenInECL
use seen in ECL check for filtering secondaries
double m_thresholdImportantEnergy
A particle which got 'stuck' and has less than this energy will be killed after m_thresholdTrials tri...
std::string m_physicsList
The name of the physics list which is used for the simulation.
int m_runEventVerbosity
Geant4 run/event verbosity: 0=Silent; 1=info level; 2=debug level, default=0.
double m_magneticCacheDistance
minimal distance for magnetic field lookup.
double m_eclProductionCut
Secondary production threshold in ECL envelopes.
bool m_EnableVisualization
If set to true the Geant4 visualization support is enabled.
double m_secondariesEnergyCut
kinetic energy cut for the stored Geant secondaries
double m_kineticEnergyThreshold
kinetic energy threshold for filtering secondaries
std::string m_magneticFieldName
magnetic field stepper to use
G4MagneticField * m_magneticField
Pointer to the (un)cached magnetic field.
double m_regionRho
Region rho limit for filtering secondaries.
G4MagneticField * m_uncachedField
Pointer to the uncached magnetic field (might be superseded by its cached version)
int m_maxNumberSteps
The maximum number of steps before the track transportation is stopped and the track is killed.
double m_deltaChordInMagneticField
The maximum miss-distance between the trajectory curve and its linear chord(s) approximation.
G4ChordFinder * m_chordFinder
Pointer to the equation-of-motion chord finder (if not the default)
MCParticleGraph m_mcParticleGraph
The MCParticle Graph used to manage the MCParticles before and after the simulation.
virtual ~FullSimModule()
Destructor of the module.
bool m_monopoles
If set to true, G4MonopolePhysics is registered in Geant4 PhysicsList.
double m_bremsstrahlungPhotonsEnergyCut
kinetic energy cut for the stored bremsstrahlung photons
int m_trajectoryStore
If true, store the trajectories of all primary particles.
bool m_storeSecondaries
controls storing of Geant secondaries in MCParticles
double m_regionZForward
Region forward z limit for filtering secondaries.
G4StepLimiter * m_stepLimiter
Pointer to the step limiter.
FullSimModule()
Constructor of the module.
LongLivedNeutral physics Class – to be registered in the physics list.
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition Module.cc:208
Module()
Constructor.
Definition Module.cc:30
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition Module.h:80
@ c_TerminateInAllProcesses
When using parallel processing, call this module's terminate() function in all processes().
Definition Module.h:83
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:34
The Class for the Belle2 magnetic field implementation for Geant4.
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:32
void beginRun(int runNumber)
Prepares Geant4 for a new run.
Definition RunManager.cc:50
void destroy()
Destroys the RunManager at the end of the simulation.
Definition RunManager.cc:86
void processEvent(int evtNumber)
Process a single event in Geant4.
Definition RunManager.cc:68
void endRun()
Terminates a Geant4 run.
Definition RunManager.cc:80
void Initialize()
Initialize the Kernel.
Definition RunManager.cc:36
static RunManager & Instance()
Static method to get a reference to the RunManager instance.
Definition RunManager.cc:30
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 setUseSeenInECL(bool use)
Set whether to check if particle is seen in ECL for ignoring secondaries.
void setDistanceThreshold(double threshold)
Set the distance threshold for ignoring secondaries.
void setRegionZForward(double z)
Set the forward z limit for improved matching region.
void setKineticEnergyThreshold(double threshold)
Set the kinetic energy threshold for ignoring secondaries.
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 setRegionRho(double rho)
Set the rho limit for improved matching region.
void setUseDetailedParticleMatching(bool use)
Set whether to use detailed particle matching.
void setRegionZBackward(double z)
Set the backward z limit for improved matching region.
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 setDoNotStoreNuclei(bool value)
Set whether to check if particle is Nuclei for ignoring secondaries.
void setDoNotStoreEMParticles(bool value)
Set whether to check if particle is EM for ignoring secondaries.
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.
Accessor to arrays stored in the data store.
Definition StoreArray.h:113
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
static const double mm
[millimeters]
Definition Unit.h:70
static const double e
Standard of [electric charge].
Definition Unit.h:53
static const double MeV
[megaelectronvolt]
Definition Unit.h:114
STL iterator class.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
Abstract base class for different kinds of events.
STL namespace.