BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
A holder class for all representations of the accelerator model created in BDSIM. More...
#include <BDSAcceleratorModel.hh>
Public Member Functions | |
G4VPhysicalVolume * | WorldPV () const |
Access the physical volume of the world. | |
G4LogicalVolume * | WorldLV () const |
Access the logical volume of the world. | |
void | RegisterBeamlineSetMain (const BDSBeamlineSet &setIn) |
Register the main beam line set. | |
void | RegisterBeamlineSetExtra (const G4String &name, const BDSBeamlineSet &setIn) |
Register a set of beam lines to be managed and cleared up at the end of the simulation. | |
void | RegisterParallelWorld (G4VUserParallelWorld *world) |
void | RegisterPlacementBeamline (BDSBeamline *placementBLIn) |
Register the beam line of arbitrary placements. | |
BDSBeamline * | PlacementBeamline () const |
Access the beam line of arbitrary placements. | |
void | RegisterPlacementFieldPlacements (const std::vector< BDSPlacementToMake > &pIn) |
Register complete placements to make for field volumes in parallel world. | |
const std::vector< BDSPlacementToMake > & | PlacementFieldPWPlacements () const |
Access field volumes for parallel world. | |
void | RegisterBLMs (BDSBeamline *blmsIn) |
Register a beam line of blm placements. | |
BDSBeamline * | BLMsBeamline () const |
Access the beam line of blm placements. | |
void | RegisterTunnelBeamline (BDSBeamline *beamlineIn) |
Register the beam line containing all the tunnel segments. | |
BDSBeamline * | TunnelBeamline () const |
Access the beam line containing all the tunnel segments. | |
void | RegisterFields (std::vector< BDSFieldObjects * > &fieldsIn) |
Register all field objects. | |
void | RegisterRegion (BDSRegion *region) |
Register a region. | |
void | RegisterAperture (const G4String &name, BDSApertureInfo *apertureIn) |
Register a single aperture. | |
void | RegisterApertures (const std::map< G4String, BDSApertureInfo * > &aperturesIn) |
Register a map of apertures with associated names. | |
BDSApertureInfo * | Aperture (G4String name) const |
G4Region * | Region (const G4String &name) const |
Access region information. Will exit if not found. | |
std::set< G4LogicalVolume * > * | VolumeSet (const G4String &name) |
Returns pointer to a set of logical volumes. If no set by that name exits, create it. | |
G4bool | VolumeInSet (G4LogicalVolume *volume, const G4String ®istryName) |
BDSBeamline * | CorrespondingMassWorldBeamline (BDSBeamline *bl) const |
G4bool | BeamlineIsMassWorld (BDSBeamline *bl) const |
void | MassWorldBeamlineAndIndex (BDSBeamline *&bl, G4int &index) const |
void | RegisterScorerHistogramDefinition (const BDSScorerHistogramDef &def) |
void | RegisterScorerPlacement (const G4String &meshName, const G4Transform3D &placement) |
Register a copy of the scorer placement so it can be used in the output. | |
void | RegisterLinkComponent (BDSLinkComponent *linkComponentIn) |
const std::set< BDSLinkComponent * > & | LinkComponents () const |
void | RegisterWorldPV (G4VPhysicalVolume *worldIn) |
Register constituent of world. | |
void | RegisterWorldLV (G4LogicalVolume *worldIn) |
Register constituent of world. | |
void | RegisterWorldSolid (G4VSolid *worldIn) |
Register constituent of world. | |
const BDSBeamlineSet & | BeamlineSetMain () const |
Accessor. | |
const BDSBeamlineSet & | BeamlineSet (const G4String &name) const |
Accessor. | |
const std::map< G4String, BDSBeamlineSet > & | ExtraBeamlines () const |
Accessor. | |
const BDSBeamline * | BeamlineMain () const |
Accessor. | |
const std::vector< BDSScorerHistogramDef > & | ScorerHistogramDefinitions () const |
Access all scorer histogram definitions. | |
const std::map< G4String, BDSScorerHistogramDef > & | ScorerHistogramDefinitionsMap () const |
Access all scorer histogram definitions. | |
const BDSScorerHistogramDef * | ScorerHistogramDef (const G4String &name) |
Access all scorer histogram definitions. | |
const std::map< G4String, G4Transform3D > & | ScorerMeshPlacementsMap () const |
Access all scorer histogram definitions. | |
Static Public Member Functions | |
static BDSAcceleratorModel * | Instance () |
Private Member Functions | |
BDSAcceleratorModel () | |
Default constructor is private as singleton. | |
void | MapBeamlineSet (const BDSBeamlineSet &setIn) |
Utility function to apply mapping. | |
Private Attributes | |
G4VPhysicalVolume * | worldPV |
Physical volume of the mass world. | |
G4LogicalVolume * | worldLV |
G4VSolid * | worldSolid |
BDSBeamlineSet | mainBeamlineSet |
std::map< G4String, BDSBeamlineSet > | extraBeamlines |
Extra beamlines. | |
std::set< G4VUserParallelWorld * > | parallelWorlds |
Parallel worlds not use with beam lines. | |
std::map< BDSBeamline *, BDSBeamline * > | clToMassWorldMap |
Mapping from any curvilinear beam line to the corresponding mass world beam line. | |
std::map< BDSBeamline *, G4bool > | massWorldMapTF |
BDSBeamline * | tunnelBeamline |
Tunnel segments beam line. | |
BDSBeamline * | placementBeamline |
Placement beam line. | |
BDSBeamline * | blmsBeamline |
BLMs beam line. | |
std::vector< BDSFieldObjects * > | fields |
All field objects. | |
std::map< G4String, BDSRegion * > | regions |
std::set< BDSRegion * > | regionStorage |
Unique storage of regions. | |
std::map< G4String, BDSApertureInfo * > | apertures |
All apertures. | |
std::vector< BDSPlacementToMake > | placementFieldPlacements |
std::map< G4String, std::set< G4LogicalVolume * > * > | volumeRegistries |
All volume registries. | |
std::map< G4String, G4Transform3D > | scorerMeshPlacements |
std::set< BDSLinkComponent * > | linkComponents |
std::vector< BDSScorerHistogramDef > | scorerHistogramDefs |
Scorer histogram definitions cached from construction here to be used in output creation. | |
std::map< G4String, BDSScorerHistogramDef > | scorerHistogramDefsMap |
Scorer histogram definitions cached from construction here to be used in output creation. | |
Static Private Attributes | |
static BDSAcceleratorModel * | instance = nullptr |
A holder class for all representations of the accelerator model created in BDSIM.
This can be extended to allow inspection of the model. Holds the readout geometry physical world in a location independent of detector construction.
In future, there may be several more beamlines than just a flat one, and perhaps grouped into a more hierarchical version. These can be contained here and this class can be extended as required.
Definition at line 60 of file BDSAcceleratorModel.hh.
BDSAcceleratorModel::~BDSAcceleratorModel | ( | ) |
Definition at line 67 of file BDSAcceleratorModel.cc.
|
private |
Default constructor is private as singleton.
Definition at line 55 of file BDSAcceleratorModel.cc.
References BDSAcceleratorComponentRegistry::Instance(), and BDSPhysicalVolumeInfoRegistry::Instance().
BDSApertureInfo * BDSAcceleratorModel::Aperture | ( | G4String | name | ) | const |
Access an aperture definition. Will exit if not found. Note, we use pointers as we purposively don't provide a default constructor for BDSApertureInfo.
Definition at line 148 of file BDSAcceleratorModel.cc.
References apertures.
Referenced by BDSParallelWorldSampler::BuildSampler(), and BDSDetectorConstruction::CalculateExtentOfSamplerPlacement().
G4bool BDSAcceleratorModel::BeamlineIsMassWorld | ( | BDSBeamline * | bl | ) | const |
Return whether a beam line is a mass world beam line. If in the unlikely event the beam line isn't registered, false is returned by default.
Definition at line 213 of file BDSAcceleratorModel.cc.
References BDS::MapGetWithDefault(), and massWorldMapTF.
Referenced by MassWorldBeamlineAndIndex().
|
inline |
Accessor.
Definition at line 86 of file BDSAcceleratorModel.hh.
Referenced by BDSBunch::ApplyCurvilinearTransform(), BDSDetectorConstruction::BuildTunnel(), BDSOutput::CalculateHistogramParameters(), BDSRunAction::CheckTrajectoryOptions(), BDSParallelWorldSampler::Construct(), BDSDetectorConstruction::ConstructScoringMeshes(), BDSOutput::CreateHistograms(), BDSOutputROOTEventModel::Fill(), BDSOutputStructures::PrepareCavityInformation(), and BDSOutputStructures::PrepareCollimatorInformation().
const BDSBeamlineSet & BDSAcceleratorModel::BeamlineSet | ( | const G4String & | name | ) | const |
Accessor.
Definition at line 125 of file BDSAcceleratorModel.cc.
References extraBeamlines.
Referenced by BDSParallelWorldCurvilinear::Construct(), and BDSParallelWorldCurvilinearBridge::Construct().
|
inline |
Accessor.
Definition at line 83 of file BDSAcceleratorModel.hh.
Referenced by BDSDetectorConstruction::BuildWorld(), BDSDetectorConstruction::ComponentPlacement(), and BDSDetectorConstruction::Construct().
|
inline |
Access the beam line of blm placements.
Definition at line 107 of file BDSAcceleratorModel.hh.
References blmsBeamline.
BDSBeamline * BDSAcceleratorModel::CorrespondingMassWorldBeamline | ( | BDSBeamline * | bl | ) | const |
Find a corresponding mass world beam line for a curvilinear (or bridge) beam line from the registered beam line sets.
Definition at line 191 of file BDSAcceleratorModel.cc.
References clToMassWorldMap.
Referenced by BDSSDCollimator::ProcessHitsOrdered().
|
inline |
Accessor.
Definition at line 85 of file BDSAcceleratorModel.hh.
References extraBeamlines.
Referenced by BDSDetectorConstruction::BuildWorld(), and BDSDetectorConstruction::ComponentPlacement().
|
static |
Definition at line 48 of file BDSAcceleratorModel.cc.
|
inline |
Definition at line 172 of file BDSAcceleratorModel.hh.
|
private |
Utility function to apply mapping.
Definition at line 200 of file BDSAcceleratorModel.cc.
References clToMassWorldMap, and massWorldMapTF.
Referenced by RegisterBeamlineSetExtra(), and RegisterBeamlineSetMain().
void BDSAcceleratorModel::MassWorldBeamlineAndIndex | ( | BDSBeamline *& | bl, |
G4int & | index | ||
) | const |
Update a beam line pointer and index if required for the equivalent ones in the mass world beam line. If the beam line supplied is a mass world one, nothing is done.
Definition at line 218 of file BDSAcceleratorModel.cc.
References BeamlineIsMassWorld(), and clToMassWorldMap.
|
inline |
Access the beam line of arbitrary placements.
Definition at line 95 of file BDSAcceleratorModel.hh.
References placementBeamline.
Referenced by BDSDetectorConstruction::BuildWorld().
|
inline |
Access field volumes for parallel world.
Definition at line 101 of file BDSAcceleratorModel.hh.
References placementFieldPlacements.
Referenced by BDSParallelWorldPlacementFields::Construct().
G4Region * BDSAcceleratorModel::Region | ( | const G4String & | name | ) | const |
Access region information. Will exit if not found.
Definition at line 157 of file BDSAcceleratorModel.cc.
Referenced by BDSDetectorConstruction::PlaceBeamlineInWorld().
|
inline |
Register a single aperture.
Definition at line 122 of file BDSAcceleratorModel.hh.
References apertures.
void BDSAcceleratorModel::RegisterApertures | ( | const std::map< G4String, BDSApertureInfo * > & | aperturesIn | ) |
Register a map of apertures with associated names.
Definition at line 143 of file BDSAcceleratorModel.cc.
References apertures.
Referenced by BDSDetectorConstruction::InitialiseApertures().
void BDSAcceleratorModel::RegisterBeamlineSetExtra | ( | const G4String & | name, |
const BDSBeamlineSet & | setIn | ||
) |
Register a set of beam lines to be managed and cleared up at the end of the simulation.
Definition at line 115 of file BDSAcceleratorModel.cc.
References extraBeamlines, and MapBeamlineSet().
Referenced by BDSDetectorConstruction::BuildBeamlines().
void BDSAcceleratorModel::RegisterBeamlineSetMain | ( | const BDSBeamlineSet & | setIn | ) |
Register the main beam line set.
Definition at line 109 of file BDSAcceleratorModel.cc.
References MapBeamlineSet().
Referenced by BDSDetectorConstruction::BuildBeamlines().
|
inline |
Register a beam line of blm placements.
Definition at line 104 of file BDSAcceleratorModel.hh.
References blmsBeamline.
Referenced by BDSDetectorConstruction::Construct().
|
inline |
Register all field objects.
Definition at line 116 of file BDSAcceleratorModel.hh.
References fields.
Referenced by BDSDetectorConstruction::ConstructSDandField().
|
inline |
Definition at line 171 of file BDSAcceleratorModel.hh.
|
inline |
Definition at line 89 of file BDSAcceleratorModel.hh.
|
inline |
Register the beam line of arbitrary placements.
Definition at line 92 of file BDSAcceleratorModel.hh.
References placementBeamline.
Referenced by BDSDetectorConstruction::Construct().
|
inline |
Register complete placements to make for field volumes in parallel world.
Definition at line 98 of file BDSAcceleratorModel.hh.
References placementFieldPlacements.
Referenced by BDS::BuildPlacementGeometry().
void BDSAcceleratorModel::RegisterRegion | ( | BDSRegion * | region | ) |
Register a region.
Definition at line 137 of file BDSAcceleratorModel.cc.
References BDSRegion::name, and regionStorage.
Referenced by BDSDetectorConstruction::InitialiseRegions().
void BDSAcceleratorModel::RegisterScorerHistogramDefinition | ( | const BDSScorerHistogramDef & | def | ) |
Register a scorer histogram definition so it can be used in the output. The definition is stored both in a vector and a map. Note, repeated entries will exist in the vector but be overwritten in the map.
Definition at line 233 of file BDSAcceleratorModel.cc.
References scorerHistogramDefs, scorerHistogramDefsMap, and BDSScorerHistogramDef::uniqueName.
Referenced by BDSDetectorConstruction::ConstructScoringMeshes().
void BDSAcceleratorModel::RegisterScorerPlacement | ( | const G4String & | meshName, |
const G4Transform3D & | placement | ||
) |
Register a copy of the scorer placement so it can be used in the output.
Definition at line 241 of file BDSAcceleratorModel.cc.
Referenced by BDSDetectorConstruction::ConstructScoringMeshes().
|
inline |
Register the beam line containing all the tunnel segments.
Definition at line 110 of file BDSAcceleratorModel.hh.
References tunnelBeamline.
Referenced by BDSDetectorConstruction::BuildTunnel().
|
inline |
Register constituent of world.
Definition at line 68 of file BDSAcceleratorModel.hh.
Referenced by BDSDetectorConstruction::BuildWorld().
|
inline |
Register constituent of world.
Definition at line 67 of file BDSAcceleratorModel.hh.
References worldPV.
Referenced by BDSDetectorConstruction::BuildWorld().
|
inline |
Register constituent of world.
Definition at line 69 of file BDSAcceleratorModel.hh.
Referenced by BDSDetectorConstruction::BuildWorld().
const BDSScorerHistogramDef * BDSAcceleratorModel::ScorerHistogramDef | ( | const G4String & | name | ) |
Access all scorer histogram definitions.
Definition at line 247 of file BDSAcceleratorModel.cc.
References scorerHistogramDefsMap.
|
inline |
Access all scorer histogram definitions.
Definition at line 165 of file BDSAcceleratorModel.hh.
References scorerHistogramDefs.
|
inline |
Access all scorer histogram definitions.
Definition at line 166 of file BDSAcceleratorModel.hh.
References scorerHistogramDefsMap.
Referenced by BDSOutput::CreateHistograms().
|
inline |
Access all scorer histogram definitions.
Definition at line 168 of file BDSAcceleratorModel.hh.
Referenced by BDSOutput::FillModel().
|
inline |
Access the beam line containing all the tunnel segments.
Definition at line 113 of file BDSAcceleratorModel.hh.
References tunnelBeamline.
Referenced by BDSDetectorConstruction::BuildWorld(), BDSDetectorConstruction::ComponentPlacement(), and BDSOutput::CreateHistograms().
G4bool BDSAcceleratorModel::VolumeInSet | ( | G4LogicalVolume * | volume, |
const G4String & | registryName | ||
) |
Check whether a volume is in a registry of volumes (a set). If no such registry exists then return false.
Definition at line 179 of file BDSAcceleratorModel.cc.
References volumeRegistries.
std::set< G4LogicalVolume * > * BDSAcceleratorModel::VolumeSet | ( | const G4String & | name | ) |
Returns pointer to a set of logical volumes. If no set by that name exits, create it.
Definition at line 172 of file BDSAcceleratorModel.cc.
References volumeRegistries.
Referenced by BDSCollimator::Build(), BDSCollimatorJaw::Build(), BDSElement::BuildContainerLogicalVolume(), BDSDetectorConstruction::BuildPhysicsBias(), BDSCrystalFactory::CommonConstruction(), BDSLinkTrackingAction::PostUserTrackingAction(), BDSTrackingAction::PostUserTrackingAction(), and BDSCollimatorCrystal::RegisterCrystalLVs().
|
inline |
Access the logical volume of the world.
Definition at line 73 of file BDSAcceleratorModel.hh.
Referenced by BDSDetectorConstruction::ConstructScoringMeshes().
|
inline |
Access the physical volume of the world.
Definition at line 72 of file BDSAcceleratorModel.hh.
References worldPV.
|
private |
All apertures.
Definition at line 205 of file BDSAcceleratorModel.hh.
Referenced by Aperture(), RegisterAperture(), and RegisterApertures().
|
private |
BLMs beam line.
Definition at line 200 of file BDSAcceleratorModel.hh.
Referenced by BLMsBeamline(), and RegisterBLMs().
|
private |
Mapping from any curvilinear beam line to the corresponding mass world beam line.
Definition at line 189 of file BDSAcceleratorModel.hh.
Referenced by CorrespondingMassWorldBeamline(), MapBeamlineSet(), and MassWorldBeamlineAndIndex().
|
private |
Extra beamlines.
Definition at line 184 of file BDSAcceleratorModel.hh.
Referenced by BeamlineSet(), ExtraBeamlines(), and RegisterBeamlineSetExtra().
|
private |
All field objects.
Definition at line 202 of file BDSAcceleratorModel.hh.
Referenced by RegisterFields().
|
staticprivate |
Definition at line 177 of file BDSAcceleratorModel.hh.
|
private |
Definition at line 219 of file BDSAcceleratorModel.hh.
|
private |
Definition at line 183 of file BDSAcceleratorModel.hh.
|
private |
Map from beam line pointer to whether that beam line object is a mass world one, i.e. not a curvilinear or bridge one. 'TF' for true false.
Definition at line 193 of file BDSAcceleratorModel.hh.
Referenced by BeamlineIsMassWorld(), and MapBeamlineSet().
|
private |
Parallel worlds not use with beam lines.
Definition at line 186 of file BDSAcceleratorModel.hh.
|
private |
Placement beam line.
Definition at line 199 of file BDSAcceleratorModel.hh.
Referenced by PlacementBeamline(), and RegisterPlacementBeamline().
|
private |
Placements for parallel world field volumes for geometry placements (the placement beam line). These contain lvs and transforms to place in a parallel world when it's built for the coordinate look ups.
Definition at line 209 of file BDSAcceleratorModel.hh.
Referenced by PlacementFieldPWPlacements(), and RegisterPlacementFieldPlacements().
|
private |
Definition at line 203 of file BDSAcceleratorModel.hh.
|
private |
Unique storage of regions.
Definition at line 204 of file BDSAcceleratorModel.hh.
Referenced by RegisterRegion().
|
private |
Scorer histogram definitions cached from construction here to be used in output creation.
Definition at line 214 of file BDSAcceleratorModel.hh.
Referenced by RegisterScorerHistogramDefinition(), and ScorerHistogramDefinitions().
|
private |
Scorer histogram definitions cached from construction here to be used in output creation.
Definition at line 215 of file BDSAcceleratorModel.hh.
Referenced by RegisterScorerHistogramDefinition(), ScorerHistogramDef(), and ScorerHistogramDefinitionsMap().
|
private |
Definition at line 217 of file BDSAcceleratorModel.hh.
|
private |
Tunnel segments beam line.
Definition at line 198 of file BDSAcceleratorModel.hh.
Referenced by RegisterTunnelBeamline(), and TunnelBeamline().
|
private |
All volume registries.
Definition at line 211 of file BDSAcceleratorModel.hh.
Referenced by VolumeInSet(), and VolumeSet().
|
private |
Definition at line 180 of file BDSAcceleratorModel.hh.
|
private |
Physical volume of the mass world.
Definition at line 179 of file BDSAcceleratorModel.hh.
Referenced by RegisterWorldPV(), and WorldPV().
|
private |
Definition at line 181 of file BDSAcceleratorModel.hh.