BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSRunManager.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSRunManager.hh"
20#include "BDSDebug.hh"
21#include "BDSDetectorConstruction.hh"
22#include "BDSExceptionHandler.hh"
23#include "BDSExtent.hh"
24#include "BDSFieldQuery.hh"
25#include "BDSPrimaryGeneratorAction.hh"
26
27#include "G4UImanager.hh"
28
29#include "CLHEP/Random/Random.h"
30
31BDSRunManager::BDSRunManager()
32{
33 // Construct an exception handler to catch Geant4 aborts.
34 // This has to be done after G4RunManager::G4RunManager() which constructs
35 // its own default exception handler which overwrites the one in G4StateManager
36 exceptionHandler = new BDSExceptionHandler();
37}
38
39BDSRunManager::~BDSRunManager()
40{
41 delete exceptionHandler;
42}
43
45{
46 G4RunManager::Initialize();
47
48 BDSExtent worldExtent;
49 if (const auto detectorConstruction = dynamic_cast<BDSDetectorConstruction*>(userDetector))
50 {
51 worldExtent = detectorConstruction->WorldExtent();
52
54 const auto& fieldQueries = detectorConstruction->FieldQueries();
55 if (!fieldQueries.empty())
56 {
57 BDSFieldQuery querier;
58 querier.QueryFields(fieldQueries);
59 }
60 }
61 if (const auto primaryGeneratorAction = dynamic_cast<BDSPrimaryGeneratorAction*>(userPrimaryGeneratorAction))
62 {primaryGeneratorAction->SetWorldExtent(worldExtent);}
63}
64
65void BDSRunManager::BeamOn(G4int n_event,const char* macroFile,G4int n_select)
66{
67 G4RunManager::BeamOn(n_event,macroFile,n_select);
68}
69
70void BDSRunManager::DoEventLoop(G4int n_event,const char* macroFile,G4int n_select)
71{
72 // save event loop state
73 if (verboseLevel > 0)
74 {
75 // Print seed to try and recreate an event in a run
76 G4cout << __METHOD_NAME__ << "Random number generator's seed="
77 << CLHEP::HepRandom::getTheSeed() << G4endl;
78 // Print generator full state to output
79 G4cout << __METHOD_NAME__ << "Random number generator's state: " << G4endl;
80 CLHEP::HepRandom::saveFullState(G4cout);
81 }
82
83 G4RunManager::DoEventLoop(n_event, macroFile, n_select);
84}
85
87{
88 // additional output
89 if (verboseLevel > 3)
90 {
91 G4cout << __METHOD_NAME__ << "Event="<<i_event<<G4endl;
92 // Print seed to try and recreate an event in a run
93 G4cout << __METHOD_NAME__ << "Random number generator's seed="
94 << CLHEP::HepRandom::getTheSeed() << G4endl;
95 // Print generator full state to output
96 G4cout << __METHOD_NAME__ << "Random number generator's state: " << G4endl;
97 CLHEP::HepRandom::saveFullState(G4cout);
98 }
99
100 //G4RunManager::ProcessOneEvent(i_event);
101
102 // This is the same as in G4RunManager, but we check the aborted event after the primary generator action
103 currentEvent = GenerateEvent(i_event);
104 if (currentEvent->IsAborted())
105 {return;}
106 eventManager->ProcessOneEvent(currentEvent);
107 AnalyzeEvent(currentEvent);
108 UpdateScoring();
109 if (i_event < n_select_msg)
110 {G4UImanager::GetUIpointer()->ApplyCommand(msgText);}
111}
112
114{
115 G4cout << "Terminate run - trying to write and close output file" << G4endl;
116 G4RunManager::AbortRun();
117}
118
Class that constructs a Geant4 model of an accelerator.
Handler that overrides Geant4's behaviour back to a normal exception.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
Class for querying the Geant4 model for field at any point.
void QueryFields(const std::vector< BDSFieldQueryInfo * > &fieldQueries)
Vector version of above function. Unique output files for each query.
Generates primary particle vertices using BDSBunch.
virtual void Initialize()
virtual void ProcessOneEvent(G4int i_event)
For additional output.
virtual void DoEventLoop(G4int n_event, const char *macroFile=nullptr, G4int n_select=-1)
For additional output.
virtual void AbortRun(G4bool)
Run G4RunManager:AbortRun(), but give some print out feedback for the user.
virtual void BeamOn(G4int n_event, const char *macroFile=nullptr, G4int n_select=-1)
Altered BeamOn function to account for Placet synchronisation.