BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSGlobalConstants.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 "BDSBeamPipeInfo.hh"
20#include "BDSDebug.hh"
21#include "BDSException.hh"
22#include "BDSGlobalConstants.hh"
23#include "BDSFieldInfo.hh"
24#include "BDSIntegratorSetType.hh"
25#include "BDSOutputType.hh"
26#include "BDSPTCOneTurnMap.hh"
27#include "BDSParser.hh"
28#include "BDSSamplerCustom.hh"
29#include "BDSSamplerPlane.hh"
30#include "BDSSamplerRegistry.hh"
31#include "BDSTrajectoryFilter.hh"
32#include "BDSTrajectoryOptions.hh"
33#include "BDSTunnelInfo.hh"
34#include "BDSUtilities.hh"
35#include "BDSWarning.hh"
36
37#include "globals.hh"
38#include "G4Colour.hh"
39#include "G4RotationMatrix.hh"
40#include "G4ThreeVector.hh"
41#include "G4Transform3D.hh"
42#include "G4UserLimits.hh"
43#include "G4Version.hh"
44#include "G4VisAttributes.hh"
45
46#include "CLHEP/Units/SystemOfUnits.h"
47#include "CLHEP/Vector/EulerAngles.h"
48
49#include <cstdlib>
50#include <limits>
51#include <sstream>
52#include <string>
53#include <stdexcept>
54#include <utility>
55
57
59{
60 if (!instance)
61 {instance = new BDSGlobalConstants(BDSParser::Instance()->GetOptions());}
62 return instance;
63}
64
66 options(opt),
67 turnsTaken(0)
68{
71
73
74 samplerDiameter = G4double(options.samplerDiameter)*CLHEP::m;
75 curvilinearDiameter = 5*CLHEP::m;
76 curvilinearDiameterShrunkForBends = false;
77
78 // beam pipe
80 options.aper1 * CLHEP::m,
81 options.aper2 * CLHEP::m,
82 options.aper3 * CLHEP::m,
83 options.aper4 * CLHEP::m,
85 options.beampipeThickness * CLHEP::m,
87
88 // magnet geometry
89 G4double horizontalWidth = options.horizontalWidth * CLHEP::m;
91 {
92 G4cerr << __METHOD_NAME__ << "Error: option \"horizontalWidth\" " << horizontalWidth
93 << " must be greater than 2x (\"aper1\" + \"beamPipeThickness\") ("
94 << defaultBeamPipeModel->aper1 << " + " << defaultBeamPipeModel->beamPipeThickness << ")" << G4endl;
95 throw BDSException(__METHOD_NAME__,"error in beam pipe defaults");
96 }
98
99 // tunnel
101 options.tunnelThickness * CLHEP::m,
102 options.tunnelSoilThickness * CLHEP::m,
106 options.tunnelFloorOffset * CLHEP::m,
107 options.tunnelAper1 * CLHEP::m,
108 options.tunnelAper2 * CLHEP::m,
109 options.storeElossTunnel,
111
112 // defaults - parameters of the laserwire process
113 itsLaserwireWavelength = 0.532 * CLHEP::micrometer;
114 itsLaserwireDir = G4ThreeVector(1,0,0);
115 itsLaserwireTrackPhotons = true;
116 itsLaserwireTrackElectrons = true;
117
118 // initialise the default vis attributes and user limits that
119 // can be copied by various bits of geometry
120 InitVisAttributes();
121 InitDefaultUserLimits();
122
124
126
127 BDSSamplerPlane::chordLength = 10*LengthSafety();
128 BDSSamplerCustom::chordLength = 10*BDSSamplerPlane::chordLength;
129
131
132 // mark which trajectory filters have been set
133 trajectoryFiltersSet[BDSTrajectoryFilter::depth] = options.HasBeenSet("storeTrajectoryDepth");
134 trajectoryFiltersSet[BDSTrajectoryFilter::particle] = options.HasBeenSet("storeTrajectoryParticle") ||
135 options.HasBeenSet("storeTrajectoryParticleID");
136 trajectoryFiltersSet[BDSTrajectoryFilter::energyThreshold] = options.HasBeenSet("storeTrajectoryEnergyThreshold");
137 trajectoryFiltersSet[BDSTrajectoryFilter::sampler] = options.HasBeenSet("storeTrajectorySamplerID");
138 trajectoryFiltersSet[BDSTrajectoryFilter::elossSRange] = options.HasBeenSet("storeTrajectoryElossSRange");
139 trajectoryFiltersSet[BDSTrajectoryFilter::minimumZ] = options.HasBeenSet("trajCutGTZ");
140 trajectoryFiltersSet[BDSTrajectoryFilter::maximumR] = options.HasBeenSet("trajCutLTR");
141
142 if (StoreMinimalData())
143 {
144 G4cout << "\nGlobal option> storing minimal data\n" << G4endl;
145 // these options are made with respect to the defaults in parser/optionsBase.cc - i.e. no need ot set false
146 // for a default that is false -> saves code
147 // the one on the right MUST be the one returned by this class in the getter function in the header
148 auto& o = options;
149 std::map<std::string, bool*> otc = {
150 {"storeApertureImpacts", &o.storeApertureImpacts},
151 {"storeApertureImpactsHistograms", &o.storeApertureImpactsHistograms},
152 {"storeCollimatorInfo", &o.storeCollimatorInfo},
153 {"storeCollimatorHits", &o.storeCollimatorHits},
154 {"storeCollimatorHitsLinks", &o.storeCollimatorHitsLinks},
155 {"storeCollimatorHitsIons", &o.storeCollimatorHitsIons},
156 {"storeCollimatorHitsAll", &o.storeCollimatorHitsAll},
157 {"storePrimaryHistograms", &o.storePrimaryHistograms},
158 {"storeTrajectoryTransportationSteps", &o.storeTrajectoryTransportationSteps},
159 {"storeModel", &o.storeModel}
160 };
161 for (auto& no : otc)
162 {
163 if (!options.HasBeenSet(no.first))
164 {*no.second = false;}
165 }
166 // try again for options that have multiple versions and check if any are set
167 // even though there's only one member bool we turn on / off in the options
168 std::map<std::set<std::string>, bool*> otcMultiple = {
169 {{"storeELoss", "storeEloss"}, &o.storeEloss},
170 {{"storeELoss","storeEloss"}, &o.storeEloss},
171 {{"storeELossHistograms", "storeElossHistograms"}, &o.storeElossHistograms},
172 {{"storeTrajectory", "storeTrajectories"}, &o.storeTrajectory},
173 {{"storeParticleData","storeGeant4Data"}, &o.storeParticleData},
174 {{"storePrimaries", "writePrimaries"}, &o.storePrimaries},
175 };
176 for (auto& no : otcMultiple)
177 {
178 bool hasBeenSet = false;
179 for (auto& op : no.first)
180 {hasBeenSet = hasBeenSet || options.HasBeenSet(op);}
181 if (!hasBeenSet)
182 {*no.second = false;}
183 }
184 }
185#if G4VERSION_NUMBER > 1079
186 if (options.HasBeenSet("scintYieldFactor"))
187 {BDS::Warning("The option \"scintYieldFactor\" has no effect with Geant4 11.0 onwards");}
188#endif
189}
190
192{
193 G4ThreeVector offset = G4ThreeVector(options.beamlineX * CLHEP::m,
194 options.beamlineY * CLHEP::m,
195 options.beamlineZ * CLHEP::m);
196
197 G4RotationMatrix rm;
199 {
200 G4ThreeVector axis = G4ThreeVector(options.beamlineAxisX,
203
204 rm = G4RotationMatrix(axis, options.beamlineAngle * CLHEP::rad);
205 }
206 else
207 {
208 G4double phi = options.beamlinePhi * CLHEP::rad;
209 G4double theta = options.beamlineTheta * CLHEP::rad;
210 G4double psi = options.beamlinePsi * CLHEP::rad;
211 CLHEP::HepEulerAngles ang = CLHEP::HepEulerAngles(phi, theta, psi);
212 rm = G4RotationMatrix(ang);
213 }
214
215 beamlineTransform = G4Transform3D(rm, offset);
216}
217
218void BDSGlobalConstants::InitVisAttributes()
219{
220 //for vacuum volumes
221 invisibleVisAttr = new G4VisAttributes(G4Colour::Black());
222 invisibleVisAttr->SetVisibility(false);
223 invisibleVisAttr->SetForceLineSegmentsPerCircle(options.nSegmentsPerCircle);
224
225 //for normally invisible volumes like marker / container volumes in debug mode
226 visibleDebugVisAttr = new G4VisAttributes(); //green
227 visibleDebugVisAttr->SetColour(0,0.6,0,0.1);
228 visibleDebugVisAttr->SetVisibility(true);
229 visibleDebugVisAttr->SetForceLineSegmentsPerCircle(options.nSegmentsPerCircle);
230}
231
232void BDSGlobalConstants::InitDefaultUserLimits()
233{
234 auto pteAsVector = BDS::SplitOnWhiteSpace(ParticlesToExcludeFromCuts());
235 // construct the set of PDG IDs
236 for (G4int i = 0; i < (G4int)pteAsVector.size(); i++)
237 {
238 try
239 {
240 G4int pdgID = std::stoi(pteAsVector[i]);
241 particlesToExcludeFromCutsAsSet.insert(pdgID);
242 }
243 catch (std::logic_error& e)
244 {
245 G4String msg = "Particle ID " + pteAsVector[i] + " at index " + std::to_string(i);
246 msg += " in the option particlesToExcludeFromCutsAsSet cannot be converted to an integer";
247 throw BDSException(__METHOD_NAME__, msg);
248 }
249 }
250
251 defaultUserLimits = new G4UserLimits("default_cuts");
252 const G4double maxTime = MaxTime();
253 if (maxTime > 0)
254 {
255 G4cout << __METHOD_NAME__ << "Setting maximum tracking time to " << maxTime << " ns" << G4endl;
256 defaultUserLimits->SetUserMaxTime(maxTime);
257 }
258 defaultUserLimits->SetMaxAllowedStep(MaxStepLength());
259 defaultUserLimits->SetUserMaxTrackLength(MaxTrackLength());
260 G4double minEK = MinimumKineticEnergy();
261 defaultUserLimits->SetUserMinEkine(minEK);
262 if (minEK > 0)
263 {G4cout << __METHOD_NAME__ << "Default minimum kinetic energy for model: " << minEK/CLHEP::GeV << " GeV" << G4endl;}
264 G4double minR = MinimumRange();
265 defaultUserLimits->SetUserMinRange(minR);
266 if (minR > 0)
267 {G4cout << __METHOD_NAME__ << "Default minimum range for user limits: " << minR/CLHEP::mm << " mm" << G4endl;}
268
269
270 BDSFieldInfo::defaultUL = defaultUserLimits; // update static member for field definitions
271
272 defaultUserLimitsTunnel = new G4UserLimits(*defaultUserLimits);
273 defaultUserLimitsTunnel->SetType("default_cuts_tunnel");
274 if (BDS::IsFinite(MinimumKineticEnergyTunnel()))
275 {defaultUserLimitsTunnel->SetUserMinEkine(MinimumKineticEnergyTunnel());}
276 if (TunnelIsInfiniteAbsorber())
277 {defaultUserLimitsTunnel->SetUserMinEkine(std::numeric_limits<double>::max());}
278}
279
280G4int BDSGlobalConstants::PrintModuloEvents() const
281{
282 G4int nGenerate = NGenerate();
283 G4double fraction = PrintFractionEvents();
284 G4int printModulo = (G4int)ceil(nGenerate * fraction);
285 if (printModulo < 0)
286 {printModulo = 1;}
287
288 if (!Batch())
289 {printModulo = 1;} // interactive -> print every event
290 return printModulo;
291}
292
293G4int BDSGlobalConstants::PrintModuloTurns() const
294{
295 G4int nTurns = TurnsToTake();
296 G4double fraction = PrintFractionTurns();
297 G4int printModulo = (G4int)ceil(nTurns * fraction);
298 if (printModulo < 0)
299 {printModulo = 1;}
300
301 if (!Batch())
302 {printModulo = 1;} // interactive -> print every turn
303 return printModulo;
304}
305
306BDSGlobalConstants::~BDSGlobalConstants()
307{
309 delete tunnelInfo;
310 delete defaultUserLimits;
311 delete defaultUserLimitsTunnel;
312 delete invisibleVisAttr;
313 delete visibleDebugVisAttr;
314
315 instance = nullptr;
316}
317
319{
320 if (options.storeTrajectoryELossSRange.empty())
321 {return;}
322 std::istringstream is(options.storeTrajectoryELossSRange);
323 std::string tok;
324 while (is >> tok)
325 {
326 std::size_t loc = tok.find(':',0);
327 if (loc == std::string::npos)
328 {throw BDSException(__METHOD_NAME__, "Error: no ':' character found in option storeTrajectoryELossSRange \"" + options.storeTrajectoryELossSRange + "\" - invalid range.");}
329 G4double rstart = 0;
330 G4double rend = 0;
331 try
332 {
333 rstart = std::stod(tok.substr(0, loc));
334 rend = std::stod(tok.substr(loc+1,tok.size()));
335 }
336 catch (const std::invalid_argument& e)
337 {
338 G4cerr << "Invalid value \"" << tok << "\" in option storeTrajectoryELossSRange." << G4endl;
339 throw BDSException(__METHOD_NAME__, "Error: can't convert string to number for option storeTrajectoryELossSRange.");
340 }
341 if (rend < rstart)
342 {
343 G4String message = "Error in option storeTrajectoryElossSRange - end point "
344 + std::to_string(rend) + " is less than beginning " + std::to_string(rstart) + ".";
345 throw BDSException(__METHOD_NAME__, message);
346 }
347 elossSRange.emplace_back(rstart*CLHEP::m, rend*CLHEP::m);
348 }
349}
350
352{
353 // only if true we let this option take precedence
355 {return false;}
356 else
357 {return options.storeTrajectoryTransportationSteps;}
358}
359
361{
363 StoreTrajectoryKineticEnergy(),
364 StoreTrajectoryMomentumVector(),
365 StoreTrajectoryProcesses(),
366 StoreTrajectoryTime(),
367 StoreTrajectoryLocal(),
368 StoreTrajectoryLinks(),
369 StoreTrajectoryIon(),
370 StoreTrajectoryMaterial()};
371
372 if (StoreTrajectoryAllVariables())
373 {
374 result.storeKineticEnergy = true;
375 result.storeMomentumVector = true;
376 result.storeProcesses = true;
377 result.storeTime = true;
378 result.storeLocal = true;
379 result.storeLinks = true;
380 result.storeIon = true;
381 result.storeMaterial = true;
382 }
383 return result;
384}
Holder class for all information required to describe a beam pipe model.
G4double aper1
Public member for direct access.
G4double beamPipeThickness
Public member for direct access.
General exception with possible name of object and message.
Definition: BDSException.hh:35
static G4UserLimits * defaultUL
Cache of default user limits.
A class that holds global options and constants.
GMAD::Options options
Options instance that this is largely based on and extends.
void ProcessTrajectoryELossSRange()
Process the option string and fill the below vector.
G4double curvilinearDiameter
Curvilinear diameter for CL volumes - defaults to samplerDiameter.
BDSTunnelInfo * tunnelInfo
Tunnel model.
std::bitset< BDS::NTrajectoryFilters > trajectoryFiltersSet
Which filters were used in the options.
static BDSGlobalConstants * Instance()
Access method.
BDSOutputType outputType
Output type enum for output format to be used.
BDSIntegratorSetType integratorSet
Integrator type enum for integrator set to be used.
std::vector< std::pair< G4double, G4double > > elossSRange
Pairs of S ranges to link trajectories to.
BDS::TrajectoryOptions StoreTrajectoryOptions() const
options that require some implementation.
BDSBeamPipeInfo * defaultBeamPipeModel
Default beam pipe model information.
void InitialiseBeamlineTransform()
Prepare the G4Transform3D instance from the options for the initial beam line transform.
BDSGlobalConstants()=delete
Unused default constructors.
BDSMagnetGeometryType magnetGeometryType
Magnet geometry.
G4double samplerDiameter
Cache of sampler diameter in this class so it can be updated.
G4Transform3D beamlineTransform
Transform for start of beam line.
G4int numberToGenerate
Number of particles to generate can be set from outside (by e.g. BDSBunchPtc)
G4bool StoreTrajectoryTransportationSteps() const
options that require some implementation.
void ResetTurnNumber()
Setter.
static BDSGlobalConstants * instance
Singleton instance.
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
static G4double chordLength
The chord length for all is fixed and can be static.
Holder struct of all information required to create a section of tunnel.
std::string magnetGeometryType
default magnet geometry parameters
Definition: optionsBase.h:152
double beamlineZ
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:132
int nGenerate
The number of primary events to simulate.
Definition: optionsBase.h:91
double beamlineAngle
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:139
double tunnelSoilThickness
tunnel geometry parameters
Definition: optionsBase.h:207
std::string outputFormat
Parameter for output format.
Definition: optionsBase.h:49
bool beamlineAxisAngle
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:140
double tunnelFloorOffset
tunnel geometry parameters
Definition: optionsBase.h:211
std::string tunnelMaterial
tunnel geometry parameters
Definition: optionsBase.h:208
int nSegmentsPerCircle
Number of facets per 2pi in visualisation.
Definition: optionsBase.h:390
std::string soilMaterial
tunnel geometry parameters
Definition: optionsBase.h:209
double aper1
default beampipe parameters
Definition: optionsBase.h:181
double aper2
default beampipe parameters
Definition: optionsBase.h:182
double beamlinePhi
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:133
std::string vacMaterial
vacuum material
Definition: optionsBase.h:189
bool buildTunnelFloor
tunnel geometry parameters
Definition: optionsBase.h:210
double beamlineAxisY
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:137
double tunnelAper1
tunnel geometry parameters
Definition: optionsBase.h:212
std::string beampipeMaterial
default beampipe parameters
Definition: optionsBase.h:185
double beampipeThickness
default beampipe parameters
Definition: optionsBase.h:179
double beamlineX
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:130
std::string apertureType
default beampipe parameters
Definition: optionsBase.h:180
double tunnelAper2
tunnel geometry parameters
Definition: optionsBase.h:213
double beamlineY
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:131
double beamlineAxisX
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:136
double tunnelThickness
tunnel geometry parameters
Definition: optionsBase.h:206
double aper4
default beampipe parameters
Definition: optionsBase.h:184
double aper3
default beampipe parameters
Definition: optionsBase.h:183
std::string tunnelType
tunnel geometry parameters
Definition: optionsBase.h:205
double beamlineAxisZ
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:138
double beamlinePsi
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:135
bool trajNoTransportation
kept only for backwards compatibility.
Definition: optionsBase.h:351
bool tunnelVisible
tunnel geometry parameters
Definition: optionsBase.h:214
double beamlineTheta
Initial beam line transform w.r.t. the world coordinate frame.
Definition: optionsBase.h:134
Options class.
Definition: options.h:44
bool HasBeenSet(const std::string &name) const
Whether a parameter has been set using the set_value method or not.
Definition: options.cc:139
BDSMagnetGeometryType DetermineMagnetGeometryType(G4String geometryType)
function to determine the enum type of the magnet geometry (case-insensitive)
BDSIntegratorSetType DetermineIntegratorSetType(G4String integratorSet)
Function that gives corresponding enum value for string (case-insensitive)
BDSOutputType DetermineOutputType(G4String outputType)
Determine the output format to be used from the input string.
std::vector< G4String > SplitOnWhiteSpace(const G4String &input)
Split a string on whitespace and return a vector of these 'words'.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())