BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes
BDSIMLink Class Reference

Interface class to use BDSIM with trackers. More...

#include <BDSIMLink.hh>

Collaboration diagram for BDSIMLink:
Collaboration graph

Public Member Functions

 BDSIMLink (BDSBunch *bunchIn)
 
int Initialise (int argc, char **argv, bool usualPrintOut=true, double minimumKineticEnergy=0, bool protonsAndIonsOnly=true)
 Initialise everything given these arguments. The minimumKinetic energy should be in GeV.
 
 BDSIMLink (int argc, char **argv, bool usualPrintOut=true)
 Construct and initialise BDSIM.
 
 ~BDSIMLink ()
 The destructor opens the geometry in Geant4 and deletes everything.
 
bool Initialised () const
 
int InitialisationResult () const
 
void BeamOn (int nGenerate=-1)
 
void SelectLinkElement (const std::string &elementName, bool debug=false)
 
void SelectLinkElement (int index, bool debug=false)
 
int AddLinkCollimatorJaw (const std::string &collimatorName, const std::string &materialName, double length, double halfApertureLeft, double halfApertureRight, double rotation, double xOffset, double yOffset, double jawTiltLeft=0.0, double jawTiltRight=0.0, bool buildLeftJaw=true, bool buildRightJaw=true, bool isACrystal=false, double crystalAngle=0, bool sampleIn=false)
 Use standard C++ types as expected to be used externally.
 
BDSHitsCollectionSamplerLinkSamplerHits () const
 
void ClearSamplerHits ()
 
int GetCurrentMaximumSixTrackParticleID () const
 
void SetCurrentMaximumExternalParticleID (int currentMaximumExternalParticleID)
 
G4int NSecondariesToReturn () const
 
G4int NPrimariesToReturn () const
 
int GetLinkIndex (const std::string &elementName) const
 Get the internal index of a component by name.
 
const BDSLinkComponentGetLinkComponent (int linkID) const
 
void RegisterUserPhysicsList (G4VModularPhysicsList *userPhysicsListIn)
 Provide a physics list that will be used inplace of the BDSIM generate one.
 
G4VModularPhysicsList * UserPhysicsList () const
 Access user physics list.
 
double GetChordLengthOfLinkElement (int beamlineIndex) const
 Access the length of a component. If bad name or ID given, -1 will be returned.
 
double GetChordLengthOfLinkElement (const std::string &elementName)
 Access the length of a component. If bad name or ID given, -1 will be returned.
 
double GetArcLengthOfLinkElement (int beamlineIndex) const
 Access the length of a component. If bad name or ID given, -1 will be returned.
 
double GetArcLengthOfLinkElement (const std::string &elementName)
 Access the length of a component. If bad name or ID given, -1 will be returned.
 

Private Member Functions

int Initialise (double minimumKineticEnergy=0, bool protonsAndIonsOnly=true)
 The main function where everything is constructed.
 

Private Attributes

bool ignoreSIGINT
 For cmake testing.
 
bool usualPrintOut
 Whether to allow the usual cout output.
 
bool initialised
 Whether initialisation was completed safely.
 
int initialisationResult
 Possible to not finish initialisation but have completed ok - flag for this.
 
int argcCache
 Cache of argc.
 
char ** argvCache
 Cache of argv.
 
std::vector< BDSParticleExternal * > externalParticles
 
std::map< std::string, int > nameToElementIndex
 
std::map< int, int > linkIDToBeamlineIndex
 
int currentElementIndex
 Element to track in.
 
G4VModularPhysicsList * userPhysicsList
 Optional user registered physics list.
 
BDSParserparser
 Cache of main object in BDSIM.
 
BDSOutputbdsOutput
 Cache of main object in BDSIM.
 
BDSBunchbdsBunch
 Cache of main object in BDSIM.
 
G4RunManager * runManager
 Cache of main object in BDSIM.
 
BDSLinkDetectorConstructionconstruction
 Cache of main object in BDSIM.
 
BDSLinkRunActionrunAction
 Cache of main object in BDSIM.
 

Detailed Description

Interface class to use BDSIM with trackers.

First way to use: bds = new BDSIMLink(argc, argv); bds->TrackThin(...)

Second way to use (delayed initialisation): bds = new BDSIMLink(); // modify bds in some way bds->Initialise(argc, argv); bds->TrackThin(...)

Author
Laurie Nevay

Definition at line 57 of file BDSIMLink.hh.

Constructor & Destructor Documentation

◆ BDSIMLink() [1/2]

BDSIMLink::BDSIMLink ( BDSBunch bunchIn)
explicit

Construct an instance but don't initialise. Requires initialisation with arguments argc and arv

Definition at line 75 of file BDSIMLink.cc.

◆ BDSIMLink() [2/2]

BDSIMLink::BDSIMLink ( int  argc,
char **  argv,
bool  usualPrintOut = true 
)

Construct and initialise BDSIM.

Definition at line 92 of file BDSIMLink.cc.

References initialisationResult, and Initialise().

Here is the call graph for this function:

◆ ~BDSIMLink()

BDSIMLink::~BDSIMLink ( )

Member Function Documentation

◆ AddLinkCollimatorJaw()

int BDSIMLink::AddLinkCollimatorJaw ( const std::string &  collimatorName,
const std::string &  materialName,
double  length,
double  halfApertureLeft,
double  halfApertureRight,
double  rotation,
double  xOffset,
double  yOffset,
double  jawTiltLeft = 0.0,
double  jawTiltRight = 0.0,
bool  buildLeftJaw = true,
bool  buildRightJaw = true,
bool  isACrystal = false,
double  crystalAngle = 0,
bool  sampleIn = false 
)

Use standard C++ types as expected to be used externally.

Close the geometry in preparation for running - everything is now fixed.

Definition at line 506 of file BDSIMLink.cc.

References BDSLinkDetectorConstruction::AddLinkCollimatorJaw(), bdsOutput, construction, and BDSOutput::UpdateSamplers().

Here is the call graph for this function:

◆ BeamOn()

void BDSIMLink::BeamOn ( int  nGenerate = -1)

Generate nGenerate events. If the default argument -1 is used, the number is taken from the standard input e.g. the executable option ngenerate and then the one specified in the input gmad files as an option.

Catch aborts to close output stream/file. perhaps not all are needed.

Run in either interactive or batch mode

Definition at line 352 of file BDSIMLink.cc.

References argcCache, argvCache, BDS::HandleAborts(), ignoreSIGINT, initialisationResult, initialised, BDSGlobalConstants::Instance(), runManager, BDSVisManager::StartSession(), and BDSException::what().

Here is the call graph for this function:

◆ ClearSamplerHits()

void BDSIMLink::ClearSamplerHits ( )
inline

Definition at line 109 of file BDSIMLink.hh.

◆ GetArcLengthOfLinkElement() [1/2]

double BDSIMLink::GetArcLengthOfLinkElement ( const std::string &  elementName)

Access the length of a component. If bad name or ID given, -1 will be returned.

Definition at line 483 of file BDSIMLink.cc.

References GetArcLengthOfLinkElement(), and GetLinkIndex().

Here is the call graph for this function:

◆ GetArcLengthOfLinkElement() [2/2]

double BDSIMLink::GetArcLengthOfLinkElement ( int  beamlineIndex) const

Access the length of a component. If bad name or ID given, -1 will be returned.

Definition at line 475 of file BDSIMLink.cc.

References BDSLinkComponent::ComponentArcLength().

Referenced by GetArcLengthOfLinkElement().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetChordLengthOfLinkElement() [1/2]

double BDSIMLink::GetChordLengthOfLinkElement ( const std::string &  elementName)

Access the length of a component. If bad name or ID given, -1 will be returned.

Definition at line 468 of file BDSIMLink.cc.

References GetChordLengthOfLinkElement(), and GetLinkIndex().

Here is the call graph for this function:

◆ GetChordLengthOfLinkElement() [2/2]

double BDSIMLink::GetChordLengthOfLinkElement ( int  beamlineIndex) const

Access the length of a component. If bad name or ID given, -1 will be returned.

Definition at line 460 of file BDSIMLink.cc.

References BDSLinkComponent::ComponentChordLength().

Referenced by GetChordLengthOfLinkElement().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCurrentMaximumSixTrackParticleID()

int BDSIMLink::GetCurrentMaximumSixTrackParticleID ( ) const

Definition at line 560 of file BDSIMLink.cc.

◆ GetLinkComponent()

const BDSLinkComponent * BDSIMLink::GetLinkComponent ( int  linkID) const

Definition at line 448 of file BDSIMLink.cc.

◆ GetLinkIndex()

int BDSIMLink::GetLinkIndex ( const std::string &  elementName) const

Get the internal index of a component by name.

Definition at line 439 of file BDSIMLink.cc.

Referenced by GetArcLengthOfLinkElement(), and GetChordLengthOfLinkElement().

Here is the caller graph for this function:

◆ InitialisationResult()

int BDSIMLink::InitialisationResult ( ) const
inline

Definition at line 81 of file BDSIMLink.hh.

◆ Initialise() [1/2]

int BDSIMLink::Initialise ( double  minimumKineticEnergy = 0,
bool  protonsAndIonsOnly = true 
)
private

The main function where everything is constructed.

Initialize executable command line options reader object

Parse lattice file

Update options generated by parser with those from executable options.

Check options for consistency

Explicitly initialise materials to construct required materials before global constants.

No longer needed. Everything can safely use BDSGlobalConstants from now on.

Force construction of global constants after parser has been initialised (requires materials first). This uses the options and beam from BDSParser. Non-const as we'll update the particle definition.

Initialize random number generator

Construct output

Check geant4 exists in the current environment

Construct mandatory run manager (the G4 kernel) and register mandatory initialization classes.

Register the geometry and parallel world construction methods with run manager.

Here the geometry isn't actually constructed - this is called by the runManager->Initialize()

Instantiate the specific type of bunch distribution.

We no longer need beamParticle so delete it to avoid confusion. The definition is held inside bdsBunch (can be updated dynamically).

Print the geometry tolerance

Set user action classes

Initialize G4 kernel

Implement bias operations on all volumes only after G4RunManager::Initialize()

Set verbosity levels at run and G4 event level. Per event and stepping are controlled in event, tracking and stepping action. These have to be done here due to the order of construction in Geant4.

Close the geometry in preparation for running - everything is now fixed.

Definition at line 124 of file BDSIMLink.cc.

References BDSParser::AmalgamateBeam(), BDSParser::AmalgamateOptions(), argcCache, argvCache, bdsBunch, bdsOutput, BDSExecOptions::Beam(), BDS::BuildPhysics(), BDSParser::CheckOptions(), BDS::ConstructAndRegisterParallelWorlds(), BDS::ConstructDesignAndBeamParticle(), construction, BDS::ConstructParallelWorldPhysics(), BDSBunchFactory::CreateBunch(), BDSOutputFactory::CreateOutput(), currentElementIndex, GMAD::BeamBase::distrType, BDSGeometryWriter::ExportGeometry(), BDS::Geant4EnvironmentIsSet(), BDSParser::GetBeam(), BDSParser::GetOptions(), BDSExecOptions::IgnoreSIGINT(), ignoreSIGINT, initialised, BDSOutput::InitialiseGeometryDependent(), BDSExecOptions::InputFileName(), BDSGlobalConstants::Instance(), BDSMaterials::Instance(), BDSParser::Instance(), BDSParticleDefinition::Name(), BDSOutput::NewFile(), BDSExecOptions::Options(), parser, BDSBunch::ParticleDefinition(), GMAD::OptionsBase::physicsList, BDSMaterials::PrepareRequiredMaterials(), BDSExecOptions::Print(), BDSExecOptions::PrintCopyright(), BDS::PrintDefinedParticles(), BDS::PrintPrimaryParticleProcesses(), GMAD::OptionsBase::recreate, BDS::RegisterSamplerPhysics(), runAction, runManager, BDSLinkDetectorConstruction::SetDesignParticle(), BDSBunch::SetOptions(), userPhysicsList, usualPrintOut, GMAD::OptionsBase::verbose, and BDS::VerboseEventStop().

Here is the call graph for this function:

◆ Initialise() [2/2]

int BDSIMLink::Initialise ( int  argc,
char **  argv,
bool  usualPrintOut = true,
double  minimumKineticEnergy = 0,
bool  protonsAndIonsOnly = true 
)

Initialise everything given these arguments. The minimumKinetic energy should be in GeV.

Definition at line 111 of file BDSIMLink.cc.

References argcCache, argvCache, initialisationResult, Initialise(), and usualPrintOut.

Referenced by BDSIMLink(), and Initialise().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialised()

bool BDSIMLink::Initialised ( ) const
inline

Accessor as to whether BDSIM kernel is initialised - ie all geometry and physics constructed.

Definition at line 79 of file BDSIMLink.hh.

References initialised.

◆ NPrimariesToReturn()

G4int BDSIMLink::NPrimariesToReturn ( ) const
inline

Definition at line 115 of file BDSIMLink.hh.

◆ NSecondariesToReturn()

G4int BDSIMLink::NSecondariesToReturn ( ) const
inline

Definition at line 114 of file BDSIMLink.hh.

◆ RegisterUserPhysicsList()

void BDSIMLink::RegisterUserPhysicsList ( G4VModularPhysicsList *  userPhysicsListIn)
inline

Provide a physics list that will be used inplace of the BDSIM generate one.

Definition at line 129 of file BDSIMLink.hh.

References userPhysicsList.

◆ SamplerHits()

BDSHitsCollectionSamplerLink * BDSIMLink::SamplerHits ( ) const

Definition at line 555 of file BDSIMLink.cc.

◆ SetCurrentMaximumExternalParticleID()

void BDSIMLink::SetCurrentMaximumExternalParticleID ( int  currentMaximumExternalParticleID)

Definition at line 565 of file BDSIMLink.cc.

◆ UserPhysicsList()

G4VModularPhysicsList * BDSIMLink::UserPhysicsList ( ) const
inline

Access user physics list.

Definition at line 130 of file BDSIMLink.hh.

References userPhysicsList.

Field Documentation

◆ argcCache

int BDSIMLink::argcCache
private

Cache of argc.

Definition at line 141 of file BDSIMLink.hh.

Referenced by BeamOn(), and Initialise().

◆ argvCache

char** BDSIMLink::argvCache
private

Cache of argv.

Definition at line 142 of file BDSIMLink.hh.

Referenced by BeamOn(), and Initialise().

◆ bdsBunch

BDSBunch* BDSIMLink::bdsBunch
private

Cache of main object in BDSIM.

Definition at line 147 of file BDSIMLink.hh.

Referenced by Initialise().

◆ bdsOutput

BDSOutput* BDSIMLink::bdsOutput
private

Cache of main object in BDSIM.

Definition at line 146 of file BDSIMLink.hh.

Referenced by AddLinkCollimatorJaw(), Initialise(), and ~BDSIMLink().

◆ construction

BDSLinkDetectorConstruction* BDSIMLink::construction
private

Cache of main object in BDSIM.

Definition at line 149 of file BDSIMLink.hh.

Referenced by AddLinkCollimatorJaw(), and Initialise().

◆ currentElementIndex

int BDSIMLink::currentElementIndex
private

Element to track in.

Definition at line 156 of file BDSIMLink.hh.

Referenced by Initialise().

◆ externalParticles

std::vector<BDSParticleExternal*> BDSIMLink::externalParticles
private

Definition at line 153 of file BDSIMLink.hh.

◆ ignoreSIGINT

bool BDSIMLink::ignoreSIGINT
private

For cmake testing.

Definition at line 137 of file BDSIMLink.hh.

Referenced by BeamOn(), and Initialise().

◆ initialisationResult

int BDSIMLink::initialisationResult
private

Possible to not finish initialisation but have completed ok - flag for this.

Definition at line 140 of file BDSIMLink.hh.

Referenced by BDSIMLink(), BeamOn(), Initialise(), and ~BDSIMLink().

◆ initialised

bool BDSIMLink::initialised
private

Whether initialisation was completed safely.

Definition at line 139 of file BDSIMLink.hh.

Referenced by BeamOn(), Initialise(), and Initialised().

◆ linkIDToBeamlineIndex

std::map<int, int> BDSIMLink::linkIDToBeamlineIndex
private

Definition at line 155 of file BDSIMLink.hh.

◆ nameToElementIndex

std::map<std::string, int> BDSIMLink::nameToElementIndex
private

Definition at line 154 of file BDSIMLink.hh.

◆ parser

BDSParser* BDSIMLink::parser
private

Cache of main object in BDSIM.

Definition at line 145 of file BDSIMLink.hh.

Referenced by Initialise(), and ~BDSIMLink().

◆ runAction

BDSLinkRunAction* BDSIMLink::runAction
private

Cache of main object in BDSIM.

Definition at line 150 of file BDSIMLink.hh.

Referenced by Initialise().

◆ runManager

G4RunManager* BDSIMLink::runManager
private

Cache of main object in BDSIM.

Definition at line 148 of file BDSIMLink.hh.

Referenced by BeamOn(), Initialise(), and ~BDSIMLink().

◆ userPhysicsList

G4VModularPhysicsList* BDSIMLink::userPhysicsList
private

Optional user registered physics list.

Definition at line 157 of file BDSIMLink.hh.

Referenced by Initialise(), RegisterUserPhysicsList(), and UserPhysicsList().

◆ usualPrintOut

bool BDSIMLink::usualPrintOut
private

Whether to allow the usual cout output.

Definition at line 138 of file BDSIMLink.hh.

Referenced by Initialise(), and ~BDSIMLink().


The documentation for this class was generated from the following files: