BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSIMLink.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 "BDSIMLink.hh"
20
21#include <csignal>
22#include <cstdlib> // standard headers
23#include <cstdio>
24
25#include "G4EventManager.hh" // Geant4 includes
26#include "G4GeometryManager.hh"
27#include "G4GeometryTolerance.hh"
28#include "G4Version.hh"
29#include "G4VModularPhysicsList.hh"
30
31#include "BDSAcceleratorModel.hh"
32#include "BDSAperturePointsLoader.hh"
33#include "BDSBeamPipeFactory.hh"
34#include "BDSBunch.hh"
35#include "BDSBunchFactory.hh"
36#include "BDSBunchSixTrackLink.hh"
37#include "BDSCavityFactory.hh"
38#include "BDSColours.hh"
39#include "BDSComponentFactoryUser.hh"
40#include "BDSDebug.hh"
41#include "BDSException.hh"
42#include "BDSExecOptions.hh"
43#include "BDSFieldFactory.hh"
44#include "BDSFieldLoader.hh"
45#include "BDSGeometryFactory.hh"
46#include "BDSGeometryFactorySQL.hh"
47#include "BDSGeometryWriter.hh"
48#include "BDSGlobalConstants.hh"
49#include "BDSLinkComponent.hh"
50#include "BDSLinkDetectorConstruction.hh"
51#include "BDSLinkEventAction.hh"
52#include "BDSLinkPrimaryGeneratorAction.hh"
53#include "BDSLinkRunAction.hh"
54#include "BDSLinkRunManager.hh"
55#include "BDSLinkStackingAction.hh"
56#include "BDSLinkTrackingAction.hh"
57#include "BDSMaterials.hh"
58#include "BDSOutput.hh"
59#include "BDSOutputFactory.hh"
60#include "BDSParallelWorldUtilities.hh"
61#include "BDSParser.hh"
62#include "BDSParticleExternal.hh"
63#include "BDSParticleDefinition.hh"
64#include "BDSPhysicsUtilities.hh"
65#include "BDSRandom.hh"
66#include "BDSSamplerRegistry.hh"
67#include "BDSSDManager.hh"
68#include "BDSTemporaryFiles.hh"
69#include "BDSUtilities.hh"
70#include "BDSVisManager.hh"
71
72#include <map>
73#include <set>
74
76 ignoreSIGINT(false),
77 usualPrintOut(true),
78 initialised(false),
79 initialisationResult(1),
80 argcCache(0),
81 argvCache(nullptr),
82 parser(nullptr),
83 bdsOutput(nullptr),
84 bdsBunch(bunchIn),
85 runManager(nullptr),
86 construction(nullptr),
87 runAction(nullptr),
88 currentElementIndex(0),
89 userPhysicsList(nullptr)
90{;}
91
92BDSIMLink::BDSIMLink(int argc, char** argv, bool usualPrintOutIn):
93 ignoreSIGINT(false),
94 usualPrintOut(usualPrintOutIn),
95 initialised(false),
96 initialisationResult(1),
97 argcCache(argc),
98 argvCache(argv),
99 parser(nullptr),
100 bdsOutput(nullptr),
101 bdsBunch(nullptr),
102 runManager(nullptr),
103 construction(nullptr),
104 runAction(nullptr),
105 currentElementIndex(0),
106 userPhysicsList(nullptr)
107{
109}
110
112 char** argv,
113 bool usualPrintOutIn,
114 double minimumKineticEnergy,
115 bool protonsAndIonsOnly)
116{
117 argcCache = argc;
118 argvCache = argv;
119 usualPrintOut = usualPrintOutIn;
120 initialisationResult = Initialise(minimumKineticEnergy, protonsAndIonsOnly);
122}
123
124int BDSIMLink::Initialise(double minimumKineticEnergy,
125 bool protonsAndIonsOnly)
126{
127 minimumKineticEnergy *= CLHEP::GeV; // units
128 G4bool trackerDebug = false;
129
131 const BDSExecOptions* execOptions = new BDSExecOptions(argcCache,argvCache);
132 if (usualPrintOut)
133 {execOptions->Print();}
134 ignoreSIGINT = execOptions->IgnoreSIGINT(); // different sig catching for cmake
135 execOptions->PrintCopyright();
136#ifdef BDSDEBUG
137 G4cout << __METHOD_NAME__ << "DEBUG mode is on." << G4endl;
138#endif
139
141 parser = BDSParser::Instance(execOptions->InputFileName());
143 parser->AmalgamateOptions(execOptions->Options());
144 parser->AmalgamateBeam(execOptions->Beam(), execOptions->Options().recreate);
147
150
152 delete execOptions;
153
158
160 BDSRandom::CreateRandomNumberGenerator(globalConstants->RandomEngine());
161 BDSRandom::SetSeed(); // set the seed from options
162
164 bdsOutput = BDSOutputFactory::CreateOutput(globalConstants->OutputFormat(),
165 globalConstants->OutputFileName());
166
169 {
170 G4cerr << "No Geant4 environmental variables found - please source geant4.sh environment" << G4endl;
171 G4cout << "A common fault is the wrong Geant4 environment as compared to the one BDSIM was compiled with." << G4endl;
172 return 1;
173 }
174
178
181
183 auto parallelWorldsRequiringPhysics = BDS::ConstructAndRegisterParallelWorlds(construction, true, false);
184
185 runManager->SetUserInitialization(construction);
186
187 // Set filters used in sensitive detectors that transfer particles back
188 BDSSDManager::Instance()->SetLinkMinimumEK(minimumKineticEnergy);
189 BDSSDManager::Instance()->SetLinkProtonsAndIonsOnly(protonsAndIonsOnly);
190
191#ifdef BDSDEBUG
192 G4cout << __METHOD_NAME__ << "> Constructing physics processes" << G4endl;
193#endif
194 G4String physicsListName = parser->GetOptions().physicsList;
195
196#if G4VERSION_NUMBER > 1049
197 // from 10.5 onwards they have a looping particle killer that warnings and kills particles
198 // deemed to be looping that are <100 MeV. This is unrelated to the primary energy so troublesome.
199 // set to the 'low' limits here ~10keV. This must be done before any physics is created as the
200 // parameters are copied into the transportation physics process for each particle and it's very
201 // hard to sift through and fix afterwards
202 G4PhysicsListHelper::GetPhysicsListHelper()->UseLowLooperThresholds();
203#endif
204 auto parallelWorldPhysics = BDS::ConstructParallelWorldPhysics(parallelWorldsRequiringPhysics);
205 G4int physicsVerbosity = BDSGlobalConstants::Instance()->PhysicsVerbosity();
206 G4VModularPhysicsList* physList;
207 if (userPhysicsList)
208 {
209 G4cout << "Using externally registered user defined physics list" << G4endl;
210 physList = userPhysicsList;
211 }
212 else
213 {physList = BDS::BuildPhysics(physicsListName, physicsVerbosity);}
214
215 BDS::RegisterSamplerPhysics(parallelWorldPhysics, physList);
216 // Construction of the physics lists defines the necessary particles and therefore
217 // we can calculate the beam rigidity for the particle the beam is designed w.r.t. This
218 // must happen before the geometry is constructed (which is called by
219 // runManager->Initialize()).
220 BDSParticleDefinition* designParticle = nullptr;
221 BDSParticleDefinition* beamParticle = nullptr;
222 G4bool beamDifferentFromDesignParticle = false;
224 globalConstants->FFact(),
225 designParticle,
226 beamParticle,
227 beamDifferentFromDesignParticle);
228 if (usualPrintOut)
229 {
230 G4cout << "Design particle properties: " << G4endl << *designParticle;
231 if (beamDifferentFromDesignParticle)
232 {G4cout << "Beam particle properties: " << G4endl << *beamParticle;}
233 }
234 // update rigidity where needed
235 construction->SetDesignParticle(designParticle);
236 //BDSFieldFactory::SetDesignParticle(designParticle);
237
238 //auto biasPhysics = BDS::BuildAndAttachBiasWrapper(parser->GetBiasing());
239 //if (biasPhysics)//could be nullptr and can't be passed to geant4 like this
240 // {physList->RegisterPhysics(biasPhysics);}
241 runManager->SetUserInitialization(physList);
242
244 GMAD::Beam b;
245 b.distrType = "sixtracklink";
246 if (!bdsBunch)
247 {
249 parser->GetBeam(),
250 globalConstants->BeamlineTransform(),
251 globalConstants->BeamlineS(),
252 globalConstants->GeneratePrimariesOnly());
253 }
254 else
255 {
256 bdsBunch->SetOptions(beamParticle, b, BDSBunchType::reference); // update particle definition only
257 }
260 delete beamParticle;
261
263 G4GeometryTolerance* theGeometryTolerance = G4GeometryTolerance::GetInstance();
264 if (usualPrintOut)
265 {
266 G4cout << __METHOD_NAME__ << "Geometry Tolerances: " << G4endl;
267 G4cout << __METHOD_NAME__ << std::setw(12) << "Surface: " << std::setw(7) << theGeometryTolerance->GetSurfaceTolerance() << " mm" << G4endl;
268 G4cout << __METHOD_NAME__ << std::setw(12) << "Angular: " << std::setw(7) << theGeometryTolerance->GetAngularTolerance() << " rad" << G4endl;
269 G4cout << __METHOD_NAME__ << std::setw(12) << "Radial: " << std::setw(7) << theGeometryTolerance->GetRadialTolerance() << " mm" << G4endl;
270 }
273 BDSLinkEventAction* eventAction = new BDSLinkEventAction(bdsOutput, runAction, trackerDebug);
274 runManager->SetUserAction(eventAction);
275 runManager->SetUserAction(runAction);
276 G4int verboseSteppingEventStart = globalConstants->VerboseSteppingEventStart();
277 G4int verboseSteppingEventStop = BDS::VerboseEventStop(verboseSteppingEventStart,
278 globalConstants->VerboseSteppingEventContinueFor());
279 runManager->SetUserAction(new BDSLinkTrackingAction(globalConstants->Batch(),
280 eventAction,
281 verboseSteppingEventStart,
282 verboseSteppingEventStop,
283 globalConstants->VerboseSteppingPrimaryOnly(),
284 globalConstants->VerboseSteppingLevel()));
285 runManager->SetUserAction(new BDSLinkStackingAction(globalConstants, std::set<G4int>(), protonsAndIonsOnly, minimumKineticEnergy));
286
287 /*
288 // Only add steppingaction if it is actually used, so do check here (for performance reasons)
289 G4int verboseSteppingEventStart = globalConstants->VerboseSteppingEventStart();
290 G4int verboseSteppingEventStop = BDS::VerboseEventStop(verboseSteppingEventStart,
291 globalConstants->VerboseSteppingEventContinueFor());
292 if (globalConstants->VerboseSteppingBDSIM())
293 {
294 runManager->SetUserAction(new BDSSteppingAction(true,
295 verboseSteppingEventStart,
296 verboseSteppingEventStop));
297 }
298 */
299
300 auto primaryGeneratorAction = new BDSLinkPrimaryGeneratorAction(bdsBunch, &currentElementIndex, construction, trackerDebug);
301 construction->SetPrimaryGeneratorAction(primaryGeneratorAction);
302 runManager->SetUserAction(primaryGeneratorAction);
303 //BDSFieldFactory::SetPrimaryGeneratorAction(primaryGeneratorAction);
304
306 runManager->Initialize();
307
309 construction->BuildPhysicsBias();
310
311 if (usualPrintOut && BDSGlobalConstants::Instance()->PhysicsVerbose())
312 {
315 }
316
320 runManager->SetVerboseLevel(std::min(globalConstants->VerboseRunLevel(), globalConstants->PhysicsVerbosity()));
321 G4EventManager::GetEventManager()->SetVerboseLevel(globalConstants->VerboseEventLevel());
322 G4EventManager::GetEventManager()->GetTrackingManager()->SetVerboseLevel(globalConstants->VerboseTrackingLevel());
323
325 G4bool bCloseGeometry = G4GeometryManager::GetInstance()->CloseGeometry();
326 if (!bCloseGeometry)
327 {
328 G4cerr << __METHOD_NAME__ << "error - geometry not closed." << G4endl;
329 return 1;
330 }
331
332 if (globalConstants->ExportGeometry())
333 {
334 BDSGeometryWriter geometrywriter;
335 geometrywriter.ExportGeometry(globalConstants->ExportType(),
336 globalConstants->ExportFileName());
337 }
338
339 const auto& nameInds = construction->NameToElementIndex();
340 nameToElementIndex.insert(nameInds.begin(), nameInds.end());
341
342 if (bdsOutput)
343 {
346 }
347
348 initialised = true;
349 return 0;
350}
351
352void BDSIMLink::BeamOn(int nGenerate)
353{
355 {return;} // a mode where we don't do anything
356
357 G4cout.precision(10);
359 struct sigaction act;
360 act.sa_handler = &BDS::HandleAborts;
361 sigemptyset(&act.sa_mask);
362 act.sa_flags = 0;
363 if (!ignoreSIGINT)
364 {sigaction(SIGINT, &act, nullptr);}
365 sigaction(SIGABRT, &act, nullptr);
366 sigaction(SIGTERM, &act, nullptr);
367 sigaction(SIGSEGV, &act, nullptr);
368
370 try
371 {
372 if (!BDSGlobalConstants::Instance()->Batch()) // Interactive mode
373 {
374 BDSVisManager visManager = BDSVisManager(BDSGlobalConstants::Instance()->VisMacroFileName(),
375 BDSGlobalConstants::Instance()->Geant4MacroFileName());
376 visManager.StartSession(argcCache, argvCache);
377 }
378 else
379 {// batch mode
380 if (nGenerate < 0)
381 {runManager->BeamOn(BDSGlobalConstants::Instance()->NGenerate());}
382 else
383 {runManager->BeamOn(nGenerate);}
384 }
385 }
386 catch (const BDSException& exception)
387 {
388 // don't do this for now in case it's dangerous and we try tracking with open geometry
389 //G4GeometryManager::GetInstance()->OpenGeometry();
390 G4cout << exception.what() << G4endl;
391 exit(1);
392 }
393}
394
396{
398 G4GeometryManager::GetInstance()->OpenGeometry();
399
400#ifdef BDSDEBUG
401 G4cout << __METHOD_NAME__ << "deleting..." << G4endl;
402#endif
403 if (bdsOutput)
404 {bdsOutput->CloseFile();}
405 delete bdsOutput;
406
407 try
408 {
409 // order important here because of singletons relying on each other
413 delete BDSAcceleratorModel::Instance();
415 delete BDSFieldFactory::Instance(); // this uses BDSGlobalConstants which uses BDSMaterials
417 delete BDSMaterials::Instance();
418
419 // instances not used in this file, but no other good location for deletion
420 if (initialisationResult < 2)
421 {
422 delete BDSColours::Instance();
424 //delete BDSSDManager::Instance();
427 }
428 }
429 catch (...)
430 {;} // ignore any exception as this is a destructor
431
432 delete runManager;
433 delete parser;
434
435 if (usualPrintOut)
436 {G4cout << __METHOD_NAME__ << "End of Run. Thank you for using BDSIM!" << G4endl;}
437}
438
439int BDSIMLink::GetLinkIndex(const std::string& elementName) const
440{
441 int result = -1;
442 auto search = nameToElementIndex.find(elementName);
443 if (search != nameToElementIndex.end())
444 {result = search->second;}
445 return result;
446}
447
448const BDSLinkComponent* BDSIMLink::GetLinkComponent(int linkID) const
449{
450 const BDSBeamline* bl = construction->LinkBeamline();
451 if (!bl)
452 {return nullptr;}
453 if (linkID > (int)bl->size())
454 {return nullptr;}
455 const auto rawAccComponent = bl->at(linkID)->GetAcceleratorComponent();
456 const auto linkComponent = dynamic_cast<const BDSLinkComponent*>(rawAccComponent);
457 return linkComponent;
458}
459
460double BDSIMLink::GetChordLengthOfLinkElement(int beamlineIndex) const
461{
462 const BDSLinkComponent* component = GetLinkComponent(beamlineIndex);
463 if (!component)
464 {return -1.0;} // play it safe
465 return component->ComponentChordLength();
466}
467
468double BDSIMLink::GetChordLengthOfLinkElement(const std::string& elementName)
469{
470 int linkID = GetLinkIndex(elementName);
471 int beamlineIndex = linkIDToBeamlineIndex[linkID];
472 return GetChordLengthOfLinkElement(beamlineIndex);
473}
474
475double BDSIMLink::GetArcLengthOfLinkElement(int beamlineIndex) const
476{
477 const BDSLinkComponent* component = GetLinkComponent(beamlineIndex);
478 if (!component)
479 {return -1.0;} // play it safe
480 return component->ComponentArcLength();
481}
482
483double BDSIMLink::GetArcLengthOfLinkElement(const std::string& elementName)
484{
485 int linkID = GetLinkIndex(elementName);
486 int beamlineIndex = linkIDToBeamlineIndex[linkID];
487 return GetArcLengthOfLinkElement(beamlineIndex);
488}
489
490void BDSIMLink::SelectLinkElement(const std::string& elementName, G4bool debug)
491{
492 if (debug)
493 {G4cout << "Searching for " << elementName;}
494 currentElementIndex = GetLinkIndex(elementName);
495 if (debug)
496 {G4cout << ", Index " << currentElementIndex << G4endl;}
497}
498
499void BDSIMLink::SelectLinkElement(int index, G4bool debug)
500{
501 if (debug)
502 {G4cout << "Searching for " << index << G4endl;}
503 currentElementIndex = index;
504}
505
506int BDSIMLink::AddLinkCollimatorJaw(const std::string& collimatorName,
507 const std::string& materialName,
508 double length,
509 double halfApertureLeft,
510 double halfApertureRight,
511 double rotation,
512 double xOffset,
513 double yOffset,
514 bool buildLeftJaw,
515 bool buildRightJaw,
516 bool isACrystal,
517 double crystalAngle,
518 bool sampleIn)
519{
520 G4GeometryManager* gm = G4GeometryManager::GetInstance();
521 if (gm->IsGeometryClosed())
522 {gm->OpenGeometry();}
523
524 G4int linkID = construction->AddLinkCollimatorJaw(collimatorName,
525 materialName,
526 length,
527 halfApertureLeft,
528 halfApertureRight,
529 rotation,
530 xOffset,
531 yOffset,
532 buildLeftJaw,
533 buildRightJaw,
534 isACrystal,
535 crystalAngle,
536 sampleIn);
537 // update this class's nameToElementIndex map
538 nameToElementIndex = construction->NameToElementIndex();
539 linkIDToBeamlineIndex = construction->LinkIDToBeamlineIndex();
540
541 if (bdsOutput)
543
545 G4bool bCloseGeometry = gm->CloseGeometry();
546 if (!bCloseGeometry)
547 {throw BDSException(__METHOD_NAME__, "error - geometry not closed.");}
548 return (int)linkID;
549}
550
551BDSHitsCollectionSamplerLink* BDSIMLink::SamplerHits() const
552{
553 return runAction ? runAction->SamplerHits() : nullptr;
554}
555
556int BDSIMLink::GetCurrentMaximumSixTrackParticleID() const
557{
558 return runAction ? runAction->MaximumExternalParticleID() : 0;
559}
560
561void BDSIMLink::SetCurrentMaximumExternalParticleID(int currentMaximumExternalParticleID)
562{
563 if (runAction)
564 {runAction->SetMaximumExternalParticleID(currentMaximumExternalParticleID);}
565}
static BDSAperturePointsCache * Instance()
Access the singleton instance.
void ClearCachedFiles()
Delete all cached points from memory and clear the map of files loaded.
static BDSBeamPipeFactory * Instance()
Singleton accessor.
BDSAcceleratorComponent * GetAcceleratorComponent() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
const BDSBeamlineElement * at(int iElement) const
Return a reference to the element at i.
Definition: BDSBeamline.hh:120
BeamlineVector::size_type size() const
Get the number of elements.
Definition: BDSBeamline.hh:148
static BDSBunch * CreateBunch(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const G4Transform3D &beamlineTransform=G4Transform3D::Identity, G4double beamlineS=0, G4bool generatePrimariesOnlyIn=false)
factory method
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 SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
Definition: BDSBunch.cc:78
static BDSCavityFactory * Instance()
Singleton accessor.
static BDSColours * Instance()
singleton pattern
Definition: BDSColours.cc:33
General exception with possible name of object and message.
Definition: BDSException.hh:35
const char * what() const noexcept override
Override message in std::exception.
Definition: BDSException.hh:55
Executable option processing for BDSIM.
const GMAD::Options & Options() const
Accessor for options generated by command line parsing.
void Print() const
Print out the commands and their set values.
G4bool IgnoreSIGINT() const
Accessor.
const GMAD::Beam & Beam() const
Accessor for beam options generate by command line parsing.
void PrintCopyright() const
Print out the copyright information (no exit).
G4String InputFileName() const
Acessor for convenience for the one thing that's needed before the parser options.
static BDSFieldFactory * Instance()
Public accessor method for singleton pattern.
static BDSFieldLoader * Instance()
Singleton accessor.
static BDSGeometryFactory * Instance()
Singleton accessor.
A class for writing fully constructed geometry from BDSIM out to other formats.
void ExportGeometry(G4String geometryType, G4String geometryFileName)
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
A BDSAcceleratorComponent wrapper for BDSLinkOpaqueBox.
G4double ComponentChordLength() const
Accessor.
G4double ComponentArcLength() const
Accessor.
Construction of the geometry in the case of a link model.
G4int AddLinkCollimatorJaw(const std::string &collimatorName, const std::string &materialName, G4double length, G4double halfApertureLeft, G4double halfApertureRight, G4double rotation, G4double xOffset, G4double yOffset, G4bool buildLeftJaw=true, G4bool buildRightJaw=true, G4bool isACrystal=false, G4double crystalAngle=0, G4bool sampleIn=false)
Interface to append a collimator of jaw style to the linking.
void SetDesignParticle(const BDSParticleDefinition *defIn)
Set the design particle definition.
Process information at the event level for Link to trackers.
Generates primary particle vertices using BDSBunch.
Simplified run action to hold link hits.
Wrapper from G4RunManager that provides our exception handler.
BDSIM's Geant4 stacking action.
Action to decide whether or not to store trajectory information.
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
void PrepareRequiredMaterials(G4bool verbose=false)
converts parser material list
static BDSOutput * CreateOutput(BDSOutputType format, const G4String &fileName, G4int fileNumberOffset=-1, G4int compressionLevel=-1)
factory method
virtual void NewFile()=0
Open a new file. This should call WriteHeader() in it.
virtual void InitialiseGeometryDependent()
Definition: BDSOutput.cc:166
virtual void UpdateSamplers()
Interface to allow updating samplers with dynamic construction. Only for link - not for regular use.
Definition: BDSOutput.hh:94
virtual void CloseFile()=0
void AmalgamateBeam(const GMAD::Beam &beamIn, bool recreate)
Amalgamate the input beam definition with the ones stored in teh parser.
Definition: BDSParser.cc:65
void CheckOptions()
Check options for consistency. This also checks the beam options.
Definition: BDSParser.cc:84
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
const GMAD::Beam & GetBeam() const
Return beam.
Definition: BDSParser.hh:56
const GMAD::Options & GetOptions() const
Return options.
Definition: BDSParser.hh:50
void AmalgamateOptions(const GMAD::Options &optionsIn)
Amalgamate the input options with the ones stored in the parser.
Definition: BDSParser.cc:76
Wrapper for particle definition.
G4String Name() const
Accessor.
static BDSSamplerRegistry * Instance()
Accessor for registry.
static BDSTemporaryFiles * Instance()
Singleton accessor.
The BDS Visualisation Manager.
void StartSession(int argc, char **argv)
Start interactive mode.
std::string distrType
beam parameters
Definition: beamBase.h:45
Beam class.
Definition: beam.h:44
std::string physicsList
list of physics processes
Definition: optionsBase.h:113
bool verbose
General verbosity.
Definition: optionsBase.h:64
bool recreate
Whether to recreate from a file or not.
Definition: optionsBase.h:92
std::vector< G4ParallelWorldPhysics * > ConstructParallelWorldPhysics(const std::vector< G4VUserParallelWorld * > &worlds)
Construct the parallel physics process for each sampler world.
void PrintDefinedParticles()
void RegisterSamplerPhysics(const std::vector< G4ParallelWorldPhysics * > &processes, G4VModularPhysicsList *physicsList)
Register each parallel physics process to the main physics list.
std::vector< G4VUserParallelWorld * > ConstructAndRegisterParallelWorlds(G4VUserDetectorConstruction *massWorld, G4bool buildSamplerWorld, G4bool buildPlacementFieldsWorld)
void PrintPrimaryParticleProcesses(const G4String &primaryParticleName)
G4bool Geant4EnvironmentIsSet()
Check if the geant4 environmental variables necessary for a run are set.
void HandleAborts(int signal_number)
void ConstructDesignAndBeamParticle(const GMAD::Beam &beamDefinition, G4double ffact, BDSParticleDefinition *&designParticle, BDSParticleDefinition *&beamParticle, G4bool &beamDifferentFromDesignParticle)
G4int VerboseEventStop(G4int verboseEventStart, G4int verboseEventContinueFor)
G4VModularPhysicsList * BuildPhysics(const G4String &physicsList, G4int verbosity=1)