BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSRunAction.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 "BDSAcceleratorModel.hh"
20#include "BDSAuxiliaryNavigator.hh"
21#include "BDSBeamline.hh"
22#include "BDSBunch.hh"
23#include "BDSDebug.hh"
24#include "BDSEventAction.hh"
25#include "BDSEventInfo.hh"
26#include "BDSException.hh"
27#include "BDSGlobalConstants.hh"
28#include "BDSOutput.hh"
29#include "BDSParser.hh"
30#include "BDSRunAction.hh"
31#include "BDSSamplerPlacementRecord.hh"
32#include "BDSSamplerRegistry.hh"
33
34#include "parser/beamBase.h"
35#include "parser/optionsBase.h"
36
37#include "globals.hh" // geant4 globals / types
38#include "G4ParticleTable.hh"
39#include "G4ProcessManager.hh"
40#include "G4ProcessVector.hh"
41#include "G4Run.hh"
42#include "G4Version.hh"
43
44#include "CLHEP/Random/Random.h"
45
46#include <chrono>
47#include <ctime>
48#include <sstream>
49#include <string>
50#include <vector>
51
52#if G4VERSION_NUMBER > 1049
53#include "BDSPhysicsUtilities.hh"
54#include "G4Positron.hh"
55#include "G4Electron.hh"
56#endif
57
58BDSRunAction::BDSRunAction(BDSOutput* outputIn,
59 BDSBunch* bunchGeneratorIn,
60 G4bool usingIonsIn,
61 BDSEventAction* eventActionIn,
62 G4String trajectorySamplerIDIn):
63 output(outputIn),
64 starttime(time(nullptr)),
65 seedStateAtStart(""),
66 info(nullptr),
67 bunchGenerator(bunchGeneratorIn),
68 usingIons(usingIonsIn),
69 cpuStartTime(std::clock_t()),
70 eventAction(eventActionIn),
71 trajectorySamplerID(trajectorySamplerIDIn)
72{;}
73
74BDSRunAction::~BDSRunAction()
75{
76 delete info;
77}
78
79void BDSRunAction::BeginOfRunAction(const G4Run* aRun)
80{
81 if (BDSGlobalConstants::Instance()->PrintPhysicsProcesses())
83
84 BDSAuxiliaryNavigator::ResetNavigatorStates();
85
86 // Bunch generator beginning of run action (optional mean subtraction).
87 bunchGenerator->BeginOfRunAction(aRun->GetNumberOfEventToBeProcessed());
88
91
92 info = new BDSEventInfo();
93
94 // save the random engine state
95 std::stringstream ss;
96 CLHEP::HepRandom::saveFullState(ss);
97 seedStateAtStart = ss.str();
99
100 // get the current time
101 starttime = time(nullptr);
102 info->SetStartTime(starttime);
103
104 // Output feedback
105 G4cout << __METHOD_NAME__ << "Run " << aRun->GetRunID()
106 << " start. Time is " << asctime(localtime(&starttime)) << G4endl;
107
109 output->NewFile();
110
111 // Write options now file open.
113 output->FillOptions(ob);
114
115 // Write beam
117 output->FillBeam(bb);
118
119 // Write model now file open.
120 output->FillModel();
121
122 // Write out geant4 data including particle tables.
124
125#if G4VERSION_NUMBER > 1049
126 // this apparently has to be done in the run action and doesn't work if done earlier
128 BDS::FixGeant105ThreshholdsForParticle(G4Positron::Definition());
129 BDS::FixGeant105ThreshholdsForParticle(G4Electron::Definition());
130#endif
131
132 cpuStartTime = std::clock();
133}
134
135void BDSRunAction::EndOfRunAction(const G4Run* aRun)
136{
137 // Get the current time
138 time_t stoptime = time(nullptr);
139 info->SetStopTime(stoptime);
140 // Run duration
141 G4float duration = difftime(stoptime, starttime);
142 info->SetDurationWall(G4double(duration));
143
144 // Calculate the elapsed CPU time for the event.
145 auto cpuEndTime = std::clock();
146 G4double durationCPU = static_cast<G4double>(cpuEndTime - cpuStartTime) / CLOCKS_PER_SEC;
147 info->SetDurationCPU(durationCPU);
148
149 // Output feedback
150 G4cout << G4endl << __METHOD_NAME__ << "Run " << aRun->GetRunID()
151 << " end. Time is " << asctime(localtime(&stoptime));
152
153 // Write output
154 output->FillRun(info);
155 output->CloseFile();
156 info->Flush();
157
158 // note difftime only calculates to the integer second
159 G4cout << __METHOD_NAME__ << "Run Duration >> " << (int)duration << " s" << G4endl;
160}
161
163{
164 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
165 G4ParticleTable::G4PTblDicIterator* particleIterator = particleTable->GetIterator();
166 particleIterator->reset();
167 while ((*particleIterator)())
168 {
169 G4ParticleDefinition* particle = particleIterator->value();
170 G4cout << "Particle: \"" << particle->GetParticleName() << "\", defined processes are: " << G4endl;
171 G4ProcessManager* pManager = particle->GetProcessManager();
172 if (!pManager)
173 {continue;}
174 G4ProcessVector* processList = pManager->GetProcessList();
175 if (!processList)
176 {continue;}
177 for (G4int i = 0; i < (G4int)processList->size(); i++)
178 {G4cout << "\"" << (*processList)[i]->GetProcessName() << "\"" << G4endl;}
179 G4cout << G4endl;
180 }
181}
182
184{
185 if (trajectorySamplerID.empty())
186 {return;}
187
188 std::vector<G4int> samplerIDs;
189 std::istringstream is(trajectorySamplerID);
190 G4String tok;
191 while (is >> tok)
192 {
194 G4int i = 0;
195 G4bool found = false;
196 for (const auto& samplerInfo : *samplerRegistry)
197 {
198 if (samplerInfo.UniqueName() == tok)
199 {
200 samplerIDs.push_back(i);
201 found = true;
202 break;
203 }
204 i++;
205 }
206 if (!found)
207 {throw BDSException(__METHOD_NAME__, "Error: sampler \"" + tok + "\" named in the option storeTrajectorySamplerID was not found.");}
208 }
209
211}
212
214{
215 // TODO - with multiple beam lines this will have to check for the maximum S coordinate
216 G4double maxS = 1*CLHEP::m; // 1m for margin
217 const BDSBeamline* blm = BDSAcceleratorModel::Instance()->BeamlineMain();
218 if (blm)
219 {maxS += blm->GetTotalArcLength();}
220 std::vector<std::pair<double,double>> sRangeToStore = BDSGlobalConstants::Instance()->StoreTrajectoryELossSRange();
221 for (const auto& range : sRangeToStore)
222 {
223 if (range.first > maxS)
224 {throw BDSException(__METHOD_NAME__, "S coordinate " + std::to_string(range.first / CLHEP::m) + "m in option storeTrajectoryElossSRange is beyond the length of the beam line (2m margin).");}
225 }
226}
const BDSBeamline * BeamlineMain() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
G4double GetTotalArcLength() const
Get the total ARC length for the beamline - ie the maximum s position.
Definition: BDSBeamline.hh:139
The base class for bunch distribution generators.
Definition: BDSBunch.hh:47
virtual const BDSParticleDefinition * ParticleDefinition() const
Access the beam particle definition.
Definition: BDSBunch.hh:95
virtual void BeginOfRunAction(G4int numberOfEvents)
Definition: BDSBunch.cc:252
Process information at the event level.
void SetSamplerIDsForTrajectories(const std::vector< G4int > &samplerIDsIn)
Update the vector of sampler IDs to match for trajectories.
Interface to store event information use G4 hooks.
Definition: BDSEventInfo.hh:42
void SetStopTime(const time_t &stopTimeIn)
Setters.
Definition: BDSEventInfo.hh:51
void SetStartTime(const time_t &startTimeIn)
Setters.
Definition: BDSEventInfo.hh:50
void SetDurationCPU(G4float durationCPUIn)
Setters.
Definition: BDSEventInfo.hh:53
void SetSeedStateAtStart(const G4String &seedStateAtStartIn)
Setters.
Definition: BDSEventInfo.hh:54
void SetDurationWall(G4float durationWallIn)
Setters.
Definition: BDSEventInfo.hh:52
General exception with possible name of object and message.
Definition: BDSException.hh:35
static BDSGlobalConstants * Instance()
Access method.
Output base class that defines interface for all output types.
Definition: BDSOutput.hh:73
virtual void NewFile()=0
Open a new file. This should call WriteHeader() in it.
virtual void InitialiseGeometryDependent()
Definition: BDSOutput.cc:166
void FillParticleData(G4bool writeIons)
Fill the local structure particle data with information. Also calls WriteParticleData().
Definition: BDSOutput.cc:186
void FillRun(const BDSEventInfo *info)
Copy run information to output structure.
Definition: BDSOutput.cc:371
void FillBeam(const GMAD::BeamBase *beam)
Definition: BDSOutput.cc:204
void FillModel()
Definition: BDSOutput.cc:218
void FillOptions(const GMAD::OptionsBase *options)
Definition: BDSOutput.cc:211
virtual void CloseFile()=0
const GMAD::OptionsBase * GetOptionsBase() const
Return bare options base class.
Definition: BDSParser.hh:53
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
const GMAD::BeamBase * GetBeamBase() const
Return bare beam base class.
Definition: BDSParser.hh:62
void SetTrajectorySamplerIDs() const
BDSBunch * bunchGenerator
Cache of bunch generator.
Definition: BDSRunAction.hh:77
std::clock_t cpuStartTime
Start time of run.
Definition: BDSRunAction.hh:79
BDSEventAction * eventAction
Event action for updating information at start of run.
Definition: BDSRunAction.hh:80
G4String trajectorySamplerID
Copy of option.
Definition: BDSRunAction.hh:81
std::string seedStateAtStart
Seed state at start of the run.
Definition: BDSRunAction.hh:75
void CheckTrajectoryOptions() const
void PrintAllProcessesForAllParticles() const
const G4bool usingIons
Cache of whether ions are being used (for particle table writing out).
Definition: BDSRunAction.hh:78
BDSOutput * output
Cache of output instance. Not owned by this class.
Definition: BDSRunAction.hh:73
Associated information for the placement of a sampler.
static BDSSamplerRegistry * Instance()
Accessor for registry.
Options for a beam distribution.
Definition: beamBase.h:35
Basic options class independent of Geant4.
Definition: optionsBase.h:36
void FixGeant105ThreshholdsForBeamParticle(const BDSParticleDefinition *particleDefinition)
Apply FixGeant105ThreshholdsForParticle to the beam particle definition.
void FixGeant105ThreshholdsForParticle(const G4ParticleDefinition *particleDefinition)