BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
BDSAcceleratorComponent Class Referenceabstract

Abstract class that represents a component of an accelerator. More...

#include <BDSAcceleratorComponent.hh>

Inheritance diagram for BDSAcceleratorComponent:
Inheritance graph
Collaboration diagram for BDSAcceleratorComponent:
Collaboration graph

Public Member Functions

 BDSAcceleratorComponent (const G4String &name, G4double arcLength, G4double angle, const G4String &type, BDSBeamPipeInfo *beamPipeInfo=nullptr, const G4ThreeVector &inputFaceNormalIn=G4ThreeVector(0, 0,-1), const G4ThreeVector &outputFaceNormalIn=G4ThreeVector(0, 0, 1), BDSFieldInfo *fieldInfoIn=nullptr)
 
virtual void Initialise ()
 
virtual void SetRegion (const G4String &regionIn)
 Set the region name for this component. More...
 
void SetField (BDSFieldInfo *fieldInfoIn)
 Set the field definition for the whole component. More...
 
virtual void SetMinimumKineticEnergy (G4double)
 
virtual G4String GetName () const
 The name of the component without modification. More...
 
G4double GetAngle () const
 
G4double Sagitta () const
 Calculate the sagitta - ie the distance between the chord and the arc at the centre. More...
 
G4String GetType () const
 Get a string describing the type of the component. More...
 
G4String GetRegion () const
 Get the region name for this component. More...
 
virtual G4bool HasAField () const
 Whether this component has a field or not (ie is active). Implicit cast of pointer to bool. More...
 
virtual void SetFieldUsePlacementWorldTransform ()
 
BDSBeamPipeInfoGetBeamPipeInfo () const
 
virtual G4String Material () const
 Return the name of a material associated with the component - ie the primary material. More...
 
virtual std::set< G4LogicalVolume * > GetAcceleratorVacuumLogicalVolumes () const
 Access the 'vacuum' volume(s) in this component if any. Default is empty set. More...
 
virtual std::set< G4LogicalVolume * > GetAcceleratorMaterialLogicalVolumes () const
 Return a set of logical volumes excluding the ones in the 'vacuum' set. More...
 
void IncrementCopyNumber ()
 Increment (+1) the number of times this component has been copied. More...
 
G4int GetCopyNumber () const
 Get the number of times this component has been copied. More...
 
BDSSimpleComponentEndPieceBefore () const
 
BDSSimpleComponentEndPieceAfter () const
 
virtual void SetBiasVacuumList (const std::list< std::string > &biasVacuumListIn)
 Copy the bias list to this element. More...
 
virtual void SetBiasMaterialList (const std::list< std::string > &biasMaterialListIn)
 Copy the bias list to this element. More...
 
virtual G4double GetArcLength () const
 
virtual G4double GetChordLength () const
 
G4ThreeVector InputFaceNormal () const
 
G4ThreeVector OutputFaceNormal () const
 
G4bool AngledInputFace () const
 Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory. More...
 
G4bool AngledOutputFace () const
 Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory. More...
 
std::list< std::string > GetBiasVacuumList () const
 Access the bias list copied from parser. More...
 
std::list< std::string > GetBiasMaterialList () const
 Access the bias list copied from parser. More...
 
virtual void SetInputFaceNormal (const G4ThreeVector &input)
 Allow updating of face normals. Virtual so derived class may apply it to daughters. More...
 
virtual void SetOutputFaceNormal (const G4ThreeVector &output)
 Allow updating of face normals. Virtual so derived class may apply it to daughters. More...
 
- Public Member Functions inherited from BDSGeometryComponent
 BDSGeometryComponent (G4VSolid *containerSolidIn, G4LogicalVolume *containerLVIn, const BDSExtent &extentIn=BDSExtent(), const BDSExtent &innerExtentIn=BDSExtent(), const G4ThreeVector &placementOffsetIn=G4ThreeVector(0, 0, 0), G4RotationMatrix *placementRotationIn=nullptr)
 
 BDSGeometryComponent (G4AssemblyVolume *containerAssemblyIn, const BDSExtent &extentIn=BDSExtent(), const BDSExtent &innerExtentIn=BDSExtent(), const G4ThreeVector &placementOffsetIn=G4ThreeVector(0, 0, 0), G4RotationMatrix *placementRotationIn=nullptr)
 
 BDSGeometryComponent (const BDSGeometryComponent &component)
 Copy constructor (no copying of registered objects) More...
 
BDSGeometryComponentoperator= (const BDSGeometryComponent &)=delete
 Assignment operator not used.
 
G4bool ContainerIsAssembly () const
 Whether the container is an assembly. If not, it's a logical volume. More...
 
void SetPlacementOffset (const G4ThreeVector &offsetIn)
 Set the offset from 0,0,0 that the object should ideally be placed in its parent. More...
 
G4ThreeVector GetExtentPositive () const
 Get the extent of the object in the positive direction in all dimensions. More...
 
G4ThreeVector GetExtentNegative () const
 Get the extent of the object in the negative direction in all dimensions. More...
 
void InheritExtents (BDSGeometryComponent const *const anotherComponent)
 Update the extents of this object with those of another object. More...
 
void InheritExtents (BDSGeometryComponent const *const anotherComponent, const G4ThreeVector &offset)
 
void RegisterDaughter (BDSGeometryComponent *anotherComponent)
 
void RegisterSolid (G4VSolid *solid)
 
void RegisterSolid (const std::set< G4VSolid * > &solids)
 Similarly for a set of unique solids. More...
 
void RegisterLogicalVolume (G4LogicalVolume *logicalVolume)
 
void RegisterLogicalVolume (const std::set< G4LogicalVolume * > &localVolumes)
 Apply RegisterLogicalVolume() to a set of logical volumes. More...
 
void RegisterSensitiveVolume (G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)
 
void RegisterSensitiveVolume (const std::set< G4LogicalVolume * > &sensitiveVolumes, BDSSDType sensitivityType)
 
void RegisterSensitiveVolume (const std::map< G4LogicalVolume *, BDSSDType > &sensitiveVolumes)
 Copy in a mapping of sensitive volumes to sensitive detectors. More...
 
void RegisterPhysicalVolume (G4VPhysicalVolume *physicalVolume)
 
void RegisterPhysicalVolume (const std::set< G4VPhysicalVolume * > &physicalVolumes)
 Apply RegisterPhysicalVolume() to a set of physical volumes. More...
 
void RegisterRotationMatrix (G4RotationMatrix *rotationMatrix)
 
void RegisterRotationMatrix (const std::set< G4RotationMatrix * > &rotationMatrices)
 Apply Register RotationMatrix() to a set of rotation matrices. More...
 
void RegisterVisAttributes (G4VisAttributes *visAttribute)
 
void RegisterVisAttributes (const std::set< G4VisAttributes * > &visAttributes)
 Apply RegisterVisAttribute() to a set of visualisation attributes. More...
 
void RegisterUserLimits (G4UserLimits *userLimit)
 
void RegisterUserLimits (const std::set< G4UserLimits * > &userLimits)
 Apply RegisterUserLimit to a set of user limits. More...
 
void InheritObjects (BDSGeometryComponent *component)
 
virtual std::set< G4LogicalVolume * > GetAllLogicalVolumes () const
 Access all logical volumes belonging to this component. More...
 
virtual std::set< G4LogicalVolume * > GetAllBiasingVolumes () const
 Return all logical volumes that should be used for biasing minus any that are in the excluded set. More...
 
virtual std::map< G4LogicalVolume *, BDSSDTypeGetAllSensitiveVolumes () const
 Access all sensitive volumes belonging to this component. More...
 
void MakeAllVolumesSensitive (BDSSDType stype=BDSSDType::energydep)
 
virtual void AttachSensitiveDetectors ()
 Attach a sensitive detector class to all registered sensitive volumes in this component. More...
 
virtual void ExcludeLogicalVolumeFromBiasing (G4LogicalVolume *lv)
 
void StripOuterAndMakeAssemblyVolume ()
 Change from a container logical volume to an assembly volume. More...
 
G4VSolid * GetContainerSolid () const
 Accessor - see member for more info. More...
 
G4LogicalVolume * GetContainerLogicalVolume () const
 Accessor - see member for more info. More...
 
G4AssemblyVolume * GetContainerAssemblyVolume () const
 Accessor - see member for more info. More...
 
G4Transform3D GetPlacementTransform () const
 Accessor - see member for more info. More...
 
G4ThreeVector GetPlacementOffset () const
 Accessor - see member for more info. More...
 
G4RotationMatrix * GetPlacementRotation () const
 Accessor - see member for more info. More...
 
BDSExtent GetExtent () const
 Accessor - see member for more info. More...
 
BDSExtent GetInnerExtent () const
 Accessor - see member for more info. More...
 
std::pair< G4double, G4double > GetExtentX () const
 Accessor - see member for more info. More...
 
std::pair< G4double, G4double > GetExtentY () const
 Accessor - see member for more info. More...
 
std::pair< G4double, G4double > GetExtentZ () const
 Accessor - see member for more info. More...
 
std::pair< G4double, G4double > GetInnerExtentX () const
 Accessor - see member for more info. More...
 
std::pair< G4double, G4double > GetInnerExtentY () const
 Accessor - see member for more info. More...
 
std::pair< G4double, G4double > GetInnerExtentZ () const
 Accessor - see member for more info. More...
 
virtual std::set< G4VPhysicalVolume * > GetAllPhysicalVolumes () const
 Accessor - see member for more info. More...
 
virtual std::set< G4RotationMatrix * > GetAllRotationMatrices () const
 Accessor - see member for more info. More...
 
virtual std::set< G4VisAttributes * > GetAllVisAttributes () const
 Accessor - see member for more info. More...
 
virtual std::set< G4UserLimits * > GetAllUserLimits () const
 Accessor - see member for more info. More...
 
virtual std::set< BDSGeometryComponent * > GetAllDaughters () const
 Accessor - see member for more info. More...
 
virtual std::set< G4VSolid * > GetAllSolids () const
 Accessor - see member for more info. More...
 
void SetExtent (const BDSExtent &extIn)
 Set extent. More...
 
void SetInnerExtent (const BDSExtent &extIn)
 Set extent. More...
 

Static Public Attributes

static G4double lengthSafetyLarge = 0
 

Protected Member Functions

virtual void Build ()
 
virtual void BuildContainerLogicalVolume ()=0
 
void SetAcceleratorVacuumLogicalVolume (G4LogicalVolume *accVacLVIn)
 
void SetAcceleratorVacuumLogicalVolume (const std::set< G4LogicalVolume * > &accVacLVIn)
 
virtual void BuildUserLimits ()
 
virtual void AttachUserLimits () const
 Doesn't change member variables, but may change their contents. More...
 

Protected Attributes

BDSBeamPipeInfobeamPipeInfo
 Optional beam pipe recipe that is written out to the survey if it exists. More...
 
std::set< G4LogicalVolume * > acceleratorVacuumLV
 A set of logical volumes we classify as 'vacuum' for biasing purposes. More...
 
BDSSimpleComponentendPieceBefore
 
BDSSimpleComponentendPieceAfter
 
G4UserLimits * userLimits
 Cache of user limits. More...
 
BDSFieldInfofieldInfo
 Recipe for field that could overlay this whole component. More...
 
const G4String name
 Const protected member variable that may not be changed by derived classes. More...
 
const G4double arcLength
 Const protected member variable that may not be changed by derived classes. More...
 
const G4String type
 Const protected member variable that may not be changed by derived classes. More...
 
G4double chordLength
 Protected member variable that can be modified by derived classes. More...
 
G4double angle
 Protected member variable that can be modified by derived classes. More...
 
G4String region
 Protected member variable that can be modified by derived classes. More...
 
- Protected Attributes inherited from BDSGeometryComponent
G4bool containerIsAssembly
 True if the 'container' is really an assembly; false if an LV. More...
 
G4VSolid * containerSolid
 
G4LogicalVolume * containerLogicalVolume
 
G4AssemblyVolume * containerAssembly
 
BDSExtent outerExtent
 
BDSExtent innerExtent
 
std::set< BDSGeometryComponent * > allDaughters
 registry of all daughter geometry components More...
 
std::set< G4VSolid * > allSolids
 registry of all solids belonging to this component More...
 
std::set< G4LogicalVolume * > allLogicalVolumes
 
std::map< G4LogicalVolume *, BDSSDTypesensitivity
 
G4bool overrideSensitivity
 
std::set< G4VPhysicalVolume * > allPhysicalVolumes
 registry of all physical volumes belonging to this component More...
 
std::set< G4RotationMatrix * > allRotationMatrices
 registry of all rotation matrices belonging to this component More...
 
std::set< G4VisAttributes * > allVisAttributes
 registry of all visualisation attributes belonging to this component More...
 
std::set< G4UserLimits * > allUserLimits
 registry of all user limits belonging to this component More...
 
G4ThreeVector placementOffset
 
G4RotationMatrix * placementRotation
 
std::set< G4LogicalVolume * > * lvsExcludedFromBiasing
 Volumes that should not be included when return GetAllBiasingVolumes(). More...
 

Static Protected Attributes

static G4double lengthSafety = -1
 Useful variable often used in construction. More...
 
static G4Material * emptyMaterial = nullptr
 Useful variable often used in construction. More...
 
static G4Material * worldMaterial = nullptr
 Useful variable often used in construction. More...
 
static G4bool checkOverlaps = false
 Useful variable often used in construction. More...
 
static G4bool sensitiveOuter = true
 Useful variable often used in construction. More...
 
static G4bool sensitiveVacuum = false
 Useful variable often used in construction. More...
 
static G4VisAttributes * containerVisAttr = nullptr
 Useful variable often used in construction. More...
 

Private Member Functions

 BDSAcceleratorComponent ()=delete
 
BDSAcceleratorComponentoperator= (const BDSAcceleratorComponent &)=delete
 Assignment and copy constructor not implemented nor used.
 
 BDSAcceleratorComponent (BDSAcceleratorComponent &)=delete
 Assignment and copy constructor not implemented nor used.
 

Private Attributes

G4bool initialised
 
G4int copyNumber
 Record of how many times this component has been copied. More...
 
G4ThreeVector inputFaceNormal
 Input face unit normal vector in incoming reference coordinate frame. More...
 
G4ThreeVector outputFaceNormal
 Output face unit normal vector in outgoing reference coordinate frame. More...
 
std::list< std::string > biasVacuumList
 Copy of bias list from parser for this particular component. More...
 
std::list< std::string > biasMaterialList
 Copy of bias list from parser for this particular component. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from BDSGeometryComponent
static void AttachUserLimitsToAssembly (G4AssemblyVolume *av, G4UserLimits *ul)
 Utility function to apply user limits to an assembly volume as there's not interface. More...
 

Detailed Description

Abstract class that represents a component of an accelerator.

It must be constructed with a name, length (arc), angle it induces (x,z plane in the local coordinates of the component) in the reference trajectory and a string representing its type. The class has no concept of its position in the beamline or in global coordinates. This information is contained in an instance of BDSBeamlineElement.

This is an abstract class as the derived class must provide the implementation of BuildContainerLogicalVolume() that constructs at least a single volume, which would be assumed to be a container if greater geometry hierarchy is constructed. This is the minimum required so that an instance of the derived class will operate with the rest of the placement machinery in BDSIM. The derived class should override Build() to add further geometry. If Build() is overridden, the base class BDSAcceleratorComponent::Build() should be called first. This calls BuildContainerVolume() and sets the visual attributes of it.

The class provides deferred construction through the Initialise() function to allow two stage construction if it's required. So, although a derived class is instantiated, the geometry isn't constructed until Initialise() is called.

Note, the geometry of any derived component should be nominally constructed along local z axis (beam direction) and x,y are transverse dimensions in a right-handed coordinate system.

The 'length' of the accelerator component is the arc length, which corresponds to the angle given. The chord length is internally calculated from the arc length and the angle. Therefore, the chord length must never be given.

The optional input and output face normal (unit) vectors are w.r.t. the incoming / outgoing reference trajectory and NOT the local geometry of the component.

Author
Laurie Nevay

Definition at line 77 of file BDSAcceleratorComponent.hh.

Constructor & Destructor Documentation

◆ BDSAcceleratorComponent() [1/2]

BDSAcceleratorComponent::BDSAcceleratorComponent ( const G4String &  name,
G4double  arcLength,
G4double  angle,
const G4String &  type,
BDSBeamPipeInfo beamPipeInfo = nullptr,
const G4ThreeVector &  inputFaceNormalIn = G4ThreeVector(0,0,-1),
const G4ThreeVector &  outputFaceNormalIn = G4ThreeVector(0,0, 1),
BDSFieldInfo fieldInfoIn = nullptr 
)

Constructor - this is the minimum information needed to create a BDSAcceleratorComponent instance. Methods in the class will allow the derived class to associate the appropriate volumes to the members of BDSGeometryComponent - the base class. The developer of a derived class should take care to set all members of BDSGeometryComponent in the derived class, including extents. Note, this class has arc length and chord length which are initially set to be the same, unless angle is != 0 in which case, the chord length is calculated from arc length. An associated beam pipe info instance can be attached if the component has a beam pipe. The input and output face normals should also be specified if non-zero. Additionally, a field info instance that represents a 'global' field for this component may be specified. Face normal (unit) vectors are w.r.t. the incoming / outgoing reference trajectory and NOT the local geometry of the component. The BDSBeamPipeInfo instance is associated with this class so that the survey output of BDSIM can query the aperture of any element.

Definition at line 49 of file BDSAcceleratorComponent.cc.

References arcLength, checkOverlaps, chordLength, containerVisAttr, emptyMaterial, BDSMaterials::GetMaterial(), initialised, BDSGlobalConstants::Instance(), BDSMaterials::Instance(), BDS::IsFinite(), lengthSafety, lengthSafetyLarge, name, sensitiveOuter, sensitiveVacuum, and worldMaterial.

Here is the call graph for this function:

◆ ~BDSAcceleratorComponent()

BDSAcceleratorComponent::~BDSAcceleratorComponent ( )
virtual

Definition at line 109 of file BDSAcceleratorComponent.cc.

◆ BDSAcceleratorComponent() [2/2]

BDSAcceleratorComponent::BDSAcceleratorComponent ( )
privatedelete

Private default constructor to force use of provided constructors, which ensure an object meets the requirements for the rest of the construction and placement machinery in BDSIM

Member Function Documentation

◆ AngledInputFace()

G4bool BDSAcceleratorComponent::AngledInputFace ( ) const

Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.

Definition at line 159 of file BDSAcceleratorComponent.cc.

References inputFaceNormal, and BDS::IsFinite().

Referenced by BDSBeamlineElement::AngledInputFace().

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

◆ AngledOutputFace()

G4bool BDSAcceleratorComponent::AngledOutputFace ( ) const

Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.

Definition at line 167 of file BDSAcceleratorComponent.cc.

References BDS::IsFinite(), and outputFaceNormal.

Referenced by BDSBeamlineElement::AngledOutputFace().

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

◆ AttachUserLimits()

void BDSAcceleratorComponent::AttachUserLimits ( ) const
protectedvirtual

Doesn't change member variables, but may change their contents.

Definition at line 186 of file BDSAcceleratorComponent.cc.

References BDSGeometryComponent::AttachUserLimitsToAssembly(), BDSGeometryComponent::containerIsAssembly, and userLimits.

Referenced by Build(), BDSElement::Build(), and BDSCavityElement::Build().

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

◆ Build()

void BDSAcceleratorComponent::Build ( )
protectedvirtual

This calls BuildContainerLogicalVolume() and then sets the visual attributes of the container logical volume. This should be overridden by derived class to add more geometry apart from the container volume. The overridden Build() function can however, call make use of this function to call BuildContainerLogicalVolume() by calling BDSAcceleratorComponent::Build() at the beginning.

Reimplemented in BDSCollimator, BDSCollimatorCrystal, BDSDegrader, BDSDrift, BDSElement, BDSGap, BDSMagnet, BDSScreen, BDSShield, BDSSimpleComponent, BDSTeleporter, BDSTerminator, BDSUndulator, BDSWireScanner, BDSCavityElement, BDSCollimatorJaw, and BDSDump.

Definition at line 136 of file BDSAcceleratorComponent.cc.

References AttachUserLimits(), BuildContainerLogicalVolume(), BuildUserLimits(), and containerVisAttr.

Referenced by BDSCollimator::Build(), BDSDegrader::Build(), BDSGap::Build(), BDSShield::Build(), BDSUndulator::Build(), BDSWireScanner::Build(), BDSCollimatorJaw::Build(), BDSDump::Build(), and Initialise().

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

◆ BuildContainerLogicalVolume()

virtual void BDSAcceleratorComponent::BuildContainerLogicalVolume ( )
protectedpure virtual

Build the container solid and logical volume that all parts of the component will contained within - must be provided by derived class.

Implemented in BDSCollimator, BDSCollimatorCrystal, BDSDegrader, BDSDrift, BDSElement, BDSGap, BDSLaserWire, BDSLine, BDSLinkComponent, BDSMagnet, BDSShield, BDSSimpleComponent, BDSTeleporter, BDSTerminator, BDSTransform3D, BDSTunnelSection, BDSUndulator, BDSWireScanner, BDSCavityElement, BDSCollimatorJaw, and BDSDump.

Referenced by Build().

Here is the caller graph for this function:

◆ BuildUserLimits()

void BDSAcceleratorComponent::BuildUserLimits ( )
protectedvirtual

This tests to see if the length of the BDSAcceleratorComponent is shorter than the global step length in the global users limits and if so build a unique one for this component and register it (memory management). It's provided by the member userLimits. This will be nullptr until this function is called, which is called in this class's Build(). Putting it here makes the same G4UserLimits object available to all derived classes potentially saving creation of a duplicate object.

Reimplemented in BDSDump.

Definition at line 175 of file BDSAcceleratorComponent.cc.

References arcLength, chordLength, BDS::CreateUserLimits(), BDSGlobalConstants::Instance(), BDSGeometryComponent::RegisterUserLimits(), and userLimits.

Referenced by Build(), BDSElement::Build(), BDSMagnet::Build(), BDSSimpleComponent::Build(), and BDSCavityElement::Build().

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

◆ EndPieceAfter()

BDSSimpleComponent * BDSAcceleratorComponent::EndPieceAfter ( ) const
inline

Definition at line 204 of file BDSAcceleratorComponent.hh.

◆ EndPieceBefore()

BDSSimpleComponent * BDSAcceleratorComponent::EndPieceBefore ( ) const
inline

Whether this component has an optional end piece that should be placed independently or not depending on other items in the beamline.

Definition at line 203 of file BDSAcceleratorComponent.hh.

Referenced by BDS::BuildEndPieceBeamline().

Here is the caller graph for this function:

◆ GetAcceleratorMaterialLogicalVolumes()

std::set< G4LogicalVolume * > BDSAcceleratorComponent::GetAcceleratorMaterialLogicalVolumes ( ) const
virtual

Return a set of logical volumes excluding the ones in the 'vacuum' set.

Reimplemented in BDSLine.

Definition at line 199 of file BDSAcceleratorComponent.cc.

References acceleratorVacuumLV, and BDSGeometryComponent::GetAllBiasingVolumes().

Referenced by BDSElement::BuildContainerLogicalVolume(), and BDSDetectorConstruction::BuildPhysicsBias().

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

◆ GetAcceleratorVacuumLogicalVolumes()

virtual std::set< G4LogicalVolume * > BDSAcceleratorComponent::GetAcceleratorVacuumLogicalVolumes ( ) const
inlinevirtual

Access the 'vacuum' volume(s) in this component if any. Default is empty set.

Reimplemented in BDSLine.

Definition at line 185 of file BDSAcceleratorComponent.hh.

References acceleratorVacuumLV.

Referenced by BDSCollimatorCrystal::Build(), BDSWireScanner::Build(), BDSCavityElement::BuildField(), and BDSDetectorConstruction::BuildPhysicsBias().

Here is the caller graph for this function:

◆ GetAngle()

G4double BDSAcceleratorComponent::GetAngle ( ) const
inline

Get the angle the component induces in the reference trajectory (rad). Note, this is 0 for h and v kickers.

Definition at line 145 of file BDSAcceleratorComponent.hh.

References angle.

Referenced by BDSBeamline::AddSingleComponent(), BDSLinkOpaqueBox::Angled(), BDSBeamlineElement::GetAngle(), BDSBeamline::GetGlobalEuclideanTransform(), BDSLinkOpaqueBox::PlaceOutputSampler(), and BDSSurvey::Write().

Here is the caller graph for this function:

◆ GetArcLength()

virtual G4double BDSAcceleratorComponent::GetArcLength ( ) const
inlinevirtual

Access the length of the component. Note there is no z length - this is chord length. Only chord OR arc makes it explicit.

Reimplemented in BDSLine.

Definition at line 139 of file BDSAcceleratorComponent.hh.

References arcLength.

Referenced by BDSBeamline::AddSingleComponent(), BDSLinkOpaqueBox::ArcLength(), BDSBeamlineElement::GetArcLength(), BDSBeamline::GetGlobalEuclideanTransform(), and BDSSurvey::Write().

Here is the caller graph for this function:

◆ GetBeamPipeInfo()

BDSBeamPipeInfo * BDSAcceleratorComponent::GetBeamPipeInfo ( ) const
inline

Access beam pipe information, which is stored in this class to provide aperture information when making a survey of the beamline consisting of accelerator components.

Definition at line 167 of file BDSAcceleratorComponent.hh.

References beamPipeInfo.

Referenced by BDSMagnet::BDSMagnet(), BDSBeamlineElement::GetBeamPipeInfo(), BDSCollimatorCrystal::Material(), BDSDrift::Material(), and BDSSurvey::Write().

Here is the caller graph for this function:

◆ GetBiasMaterialList()

std::list< std::string > BDSAcceleratorComponent::GetBiasMaterialList ( ) const
inline

Access the bias list copied from parser.

Definition at line 198 of file BDSAcceleratorComponent.hh.

References biasMaterialList.

Referenced by BDSDetectorConstruction::BuildPhysicsBias().

Here is the caller graph for this function:

◆ GetBiasVacuumList()

std::list< std::string > BDSAcceleratorComponent::GetBiasVacuumList ( ) const
inline

Access the bias list copied from parser.

Definition at line 197 of file BDSAcceleratorComponent.hh.

References biasVacuumList.

Referenced by BDSDetectorConstruction::BuildPhysicsBias().

Here is the caller graph for this function:

◆ GetChordLength()

virtual G4double BDSAcceleratorComponent::GetChordLength ( ) const
inlinevirtual

Access the length of the component. Note there is no z length - this is chord length. Only chord OR arc makes it explicit.

Reimplemented in BDSLine.

Definition at line 140 of file BDSAcceleratorComponent.hh.

References chordLength.

Referenced by BDSBeamline::AddSingleComponent(), BDS::BuildEndPieceBeamline(), BDS::BuildPlacementGeometry(), BDSLinkOpaqueBox::ChordLength(), BDSCurvilinearBuilder::CreateBridgeElementFromComponent(), BDSBeamlineElement::GetChordLength(), BDSBeamline::GetGlobalEuclideanTransform(), BDSLinkOpaqueBox::PlaceOutputSampler(), BDSBeamline::ProvideEndPieceElementAfter(), and BDSSurvey::Write().

Here is the caller graph for this function:

◆ GetCopyNumber()

G4int BDSAcceleratorComponent::GetCopyNumber ( ) const
inline

Get the number of times this component has been copied.

Definition at line 194 of file BDSAcceleratorComponent.hh.

References copyNumber.

Referenced by BDSBeamlineElement::BDSBeamlineElement().

Here is the caller graph for this function:

◆ GetName()

virtual G4String BDSAcceleratorComponent::GetName ( ) const
inlinevirtual

◆ GetRegion()

G4String BDSAcceleratorComponent::GetRegion ( ) const
inline

Get the region name for this component.

Definition at line 154 of file BDSAcceleratorComponent.hh.

References region.

◆ GetType()

G4String BDSAcceleratorComponent::GetType ( ) const
inline

Get a string describing the type of the component.

Definition at line 151 of file BDSAcceleratorComponent.hh.

References type.

Referenced by BDSBeamline::AddSingleComponent(), BDSDetectorConstruction::BuildBeamline(), BDSBeamlineElement::GetType(), BDSAcceleratorComponentRegistry::RegisterComponent(), and BDSSurvey::Write().

Here is the caller graph for this function:

◆ HasAField()

virtual G4bool BDSAcceleratorComponent::HasAField ( ) const
inlinevirtual

Whether this component has a field or not (ie is active). Implicit cast of pointer to bool.

Reimplemented in BDSMagnet.

Definition at line 157 of file BDSAcceleratorComponent.hh.

References fieldInfo.

Referenced by BDS::BuildPlacementGeometry().

Here is the caller graph for this function:

◆ IncrementCopyNumber()

void BDSAcceleratorComponent::IncrementCopyNumber ( )
inline

Increment (+1) the number of times this component has been copied.

Definition at line 191 of file BDSAcceleratorComponent.hh.

References copyNumber.

Referenced by BDSBeamlineElement::BDSBeamlineElement().

Here is the caller graph for this function:

◆ Initialise()

void BDSAcceleratorComponent::Initialise ( )
virtual

Two stage construction - first instantiate class, and then second, call this method to run Build() which constructs geometry. This allows common construction tasks to be done in one place in BDSComponentFactory rather than pass as arguments through the constructors of all derived classes. Also builds read out geometry.

Reimplemented in BDSLine, and BDSTunnelSection.

Definition at line 116 of file BDSAcceleratorComponent.cc.

References Build(), fieldInfo, initialised, BDSFieldBuilder::Instance(), and BDSFieldBuilder::RegisterFieldForConstruction().

Referenced by BDSDetectorConstruction::BuildBeamline(), BDS::BuildBLMs(), BDS::BuildPlacementGeometry(), and BDSComponentFactory::CreateComponent().

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

◆ InputFaceNormal()

G4ThreeVector BDSAcceleratorComponent::InputFaceNormal ( ) const
inline

Access face normal unit vector. This is w.r.t. the incoming / outgoing reference trajectory and NOT the local geometry of the component. Ie for an SBend with no pole face rotation this is incoming (0,0,-1). Does not account for tilt.

Definition at line 172 of file BDSAcceleratorComponent.hh.

References inputFaceNormal.

Referenced by BDSBeamline::AddSingleComponent(), BDSBeamlineElement::InputFaceNormal(), and BDSParallelWorldSampler::Place().

Here is the caller graph for this function:

◆ Material()

virtual G4String BDSAcceleratorComponent::Material ( ) const
inlinevirtual

Return the name of a material associated with the component - ie the primary material.

Reimplemented in BDSCollimator, BDSCollimatorCrystal, BDSDrift, BDSMagnet, BDSWireScanner, and BDSDump.

Definition at line 182 of file BDSAcceleratorComponent.hh.

Referenced by BDSBeamlineElement::GetMaterial(), BDSCollimatorCrystal::Material(), BDSDrift::Material(), and BDSMagnet::Material().

Here is the caller graph for this function:

◆ OutputFaceNormal()

G4ThreeVector BDSAcceleratorComponent::OutputFaceNormal ( ) const
inline

Access face normal unit vector. This is w.r.t. the incoming / outgoing reference trajectory and NOT the local geometry of the component. Ie for an SBend with no pole face rotation this is incoming (0,0,-1). Does not account for tilt.

Definition at line 173 of file BDSAcceleratorComponent.hh.

References outputFaceNormal.

Referenced by BDSBeamline::AddSingleComponent(), BDS::BuildEndPieceBeamline(), BDSCurvilinearBuilder::CreateAngledBridgeComponent(), BDSBeamlineElement::OutputFaceNormal(), and BDSParallelWorldSampler::Place().

Here is the caller graph for this function:

◆ Sagitta()

G4double BDSAcceleratorComponent::Sagitta ( ) const

Calculate the sagitta - ie the distance between the chord and the arc at the centre.

Definition at line 150 of file BDSAcceleratorComponent.cc.

References angle, arcLength, chordLength, and BDS::IsFinite().

Referenced by BDSLinkOpaqueBox::PlaceOutputSampler().

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

◆ SetAcceleratorVacuumLogicalVolume() [1/2]

void BDSAcceleratorComponent::SetAcceleratorVacuumLogicalVolume ( const std::set< G4LogicalVolume * > &  accVacLVIn)
inlineprotected

Definition at line 236 of file BDSAcceleratorComponent.hh.

◆ SetAcceleratorVacuumLogicalVolume() [2/2]

void BDSAcceleratorComponent::SetAcceleratorVacuumLogicalVolume ( G4LogicalVolume *  accVacLVIn)
inlineprotected

Assign the accelerator tracking volume - only callable by derived classes - ie not public. This is just setting a reference to the accelerator volume and it is not deleted by this class (BDSAcceleratorComponent) - therefore, the derived class should also deal with memory management of this volume - whether this is by using the inherited (from BDSGeometryComponent) RegisterLogicalVolume() or by manually deleting itself.

Definition at line 233 of file BDSAcceleratorComponent.hh.

References acceleratorVacuumLV.

Referenced by BDSCollimator::Build(), BDSCollimatorCrystal::Build(), BDSDrift::Build(), BDSUndulator::Build(), BDSCavityElement::Build(), BDSCollimatorJaw::Build(), BDSMagnet::BuildBeampipe(), BDSElement::BuildContainerLogicalVolume(), BDSLaserWire::BuildContainerLogicalVolume(), and BDSWireScanner::BuildContainerLogicalVolume().

Here is the caller graph for this function:

◆ SetBiasMaterialList()

virtual void BDSAcceleratorComponent::SetBiasMaterialList ( const std::list< std::string > &  biasMaterialListIn)
inlinevirtual

Copy the bias list to this element.

Reimplemented in BDSLine.

Definition at line 118 of file BDSAcceleratorComponent.hh.

References biasMaterialList.

Referenced by BDSComponentFactory::CreateComponent(), and BDSLine::SetBiasMaterialList().

Here is the caller graph for this function:

◆ SetBiasVacuumList()

virtual void BDSAcceleratorComponent::SetBiasVacuumList ( const std::list< std::string > &  biasVacuumListIn)
inlinevirtual

Copy the bias list to this element.

Reimplemented in BDSLine.

Definition at line 116 of file BDSAcceleratorComponent.hh.

References biasVacuumList.

Referenced by BDSComponentFactory::CreateComponent(), and BDSLine::SetBiasVacuumList().

Here is the caller graph for this function:

◆ SetField()

void BDSAcceleratorComponent::SetField ( BDSFieldInfo fieldInfoIn)

Set the field definition for the whole component.

Definition at line 145 of file BDSAcceleratorComponent.cc.

References fieldInfo.

Referenced by BDS::BuildPlacementGeometry(), and BDSComponentFactory::SetFieldDefinitions().

Here is the caller graph for this function:

◆ SetFieldUsePlacementWorldTransform()

void BDSAcceleratorComponent::SetFieldUsePlacementWorldTransform ( )
virtual

For when an accelerator component is used in a placement, we need to flag the field transform should come from a different parallel world and not the regularly curvilinear ones. This should be used before Initialise() when the field is registered for construction.

Reimplemented in BDSMagnet.

Definition at line 209 of file BDSAcceleratorComponent.cc.

References fieldInfo.

Referenced by BDS::BuildPlacementGeometry(), and BDSMagnet::SetFieldUsePlacementWorldTransform().

Here is the caller graph for this function:

◆ SetInputFaceNormal()

virtual void BDSAcceleratorComponent::SetInputFaceNormal ( const G4ThreeVector &  input)
inlinevirtual

Allow updating of face normals. Virtual so derived class may apply it to daughters.

Reimplemented in BDSMagnet.

Definition at line 207 of file BDSAcceleratorComponent.hh.

References inputFaceNormal.

Referenced by BDSCollimatorCrystal::Build(), BDSDrift::Build(), BDSUndulator::Build(), BDSWireScanner::BuildContainerLogicalVolume(), BDSMagnetOuter::SetInputFaceNormal(), and BDSMagnetOuter::SetOutputFaceNormal().

Here is the caller graph for this function:

◆ SetMinimumKineticEnergy()

virtual void BDSAcceleratorComponent::SetMinimumKineticEnergy ( G4double  )
inlinevirtual

Set a minimum kinetic energy for a component by component level cut. By default does nothing - only certain components implement a response or user limits to this.

Reimplemented in BDSCollimator.

Definition at line 130 of file BDSAcceleratorComponent.hh.

Referenced by BDSComponentFactory::CreateComponent().

Here is the caller graph for this function:

◆ SetOutputFaceNormal()

virtual void BDSAcceleratorComponent::SetOutputFaceNormal ( const G4ThreeVector &  output)
inlinevirtual

Allow updating of face normals. Virtual so derived class may apply it to daughters.

Reimplemented in BDSMagnet.

Definition at line 208 of file BDSAcceleratorComponent.hh.

References outputFaceNormal.

Referenced by BDSCollimatorCrystal::Build(), BDSDrift::Build(), BDSUndulator::Build(), and BDSWireScanner::BuildContainerLogicalVolume().

Here is the caller graph for this function:

◆ SetRegion()

virtual void BDSAcceleratorComponent::SetRegion ( const G4String &  regionIn)
inlinevirtual

Set the region name for this component.

Reimplemented in BDSLine.

Definition at line 123 of file BDSAcceleratorComponent.hh.

References region.

Referenced by BDSComponentFactory::CreateComponent(), and BDSLine::SetRegion().

Here is the caller graph for this function:

Field Documentation

◆ acceleratorVacuumLV

std::set<G4LogicalVolume*> BDSAcceleratorComponent::acceleratorVacuumLV
protected

A set of logical volumes we classify as 'vacuum' for biasing purposes.

Definition at line 276 of file BDSAcceleratorComponent.hh.

Referenced by GetAcceleratorMaterialLogicalVolumes(), GetAcceleratorVacuumLogicalVolumes(), and SetAcceleratorVacuumLogicalVolume().

◆ angle

G4double BDSAcceleratorComponent::angle
protected

Protected member variable that can be modified by derived classes.

Definition at line 258 of file BDSAcceleratorComponent.hh.

Referenced by BDSMagnet::BuildBeampipe(), BDSElement::BuildContainerLogicalVolume(), BDSMagnet::BuildOuter(), GetAngle(), and Sagitta().

◆ arcLength

const G4double BDSAcceleratorComponent::arcLength
protected

Const protected member variable that may not be changed by derived classes.

Definition at line 252 of file BDSAcceleratorComponent.hh.

Referenced by BDSAcceleratorComponent(), BDSDegrader::Build(), BDSElement::BuildContainerLogicalVolume(), BuildUserLimits(), GetArcLength(), BDSDegrader::PlaceWedge(), and Sagitta().

◆ beamPipeInfo

BDSBeamPipeInfo* BDSAcceleratorComponent::beamPipeInfo
protected

◆ biasMaterialList

std::list<std::string> BDSAcceleratorComponent::biasMaterialList
private

Copy of bias list from parser for this particular component.

Definition at line 308 of file BDSAcceleratorComponent.hh.

Referenced by GetBiasMaterialList(), and SetBiasMaterialList().

◆ biasVacuumList

std::list<std::string> BDSAcceleratorComponent::biasVacuumList
private

Copy of bias list from parser for this particular component.

Definition at line 307 of file BDSAcceleratorComponent.hh.

Referenced by GetBiasVacuumList(), and SetBiasVacuumList().

◆ checkOverlaps

G4bool BDSAcceleratorComponent::checkOverlaps = false
staticprotected

◆ chordLength

G4double BDSAcceleratorComponent::chordLength
protected

◆ containerVisAttr

G4VisAttributes * BDSAcceleratorComponent::containerVisAttr = nullptr
staticprotected

◆ copyNumber

G4int BDSAcceleratorComponent::copyNumber
private

Record of how many times this component has been copied.

Definition at line 304 of file BDSAcceleratorComponent.hh.

Referenced by GetCopyNumber(), and IncrementCopyNumber().

◆ emptyMaterial

G4Material * BDSAcceleratorComponent::emptyMaterial = nullptr
staticprotected

◆ endPieceAfter

BDSSimpleComponent* BDSAcceleratorComponent::endPieceAfter
protected

Definition at line 279 of file BDSAcceleratorComponent.hh.

◆ endPieceBefore

BDSSimpleComponent* BDSAcceleratorComponent::endPieceBefore
protected

Definition at line 278 of file BDSAcceleratorComponent.hh.

◆ fieldInfo

BDSFieldInfo* BDSAcceleratorComponent::fieldInfo
protected

Recipe for field that could overlay this whole component.

Definition at line 283 of file BDSAcceleratorComponent.hh.

Referenced by HasAField(), Initialise(), SetField(), and SetFieldUsePlacementWorldTransform().

◆ initialised

G4bool BDSAcceleratorComponent::initialised
private

Boolean record of whether this component has been already initialised. This check protects against duplicate initialisation and therefore the potential memory leaks that would ensue.

Definition at line 301 of file BDSAcceleratorComponent.hh.

Referenced by BDSAcceleratorComponent(), and Initialise().

◆ inputFaceNormal

G4ThreeVector BDSAcceleratorComponent::inputFaceNormal
private

Input face unit normal vector in incoming reference coordinate frame.

Definition at line 311 of file BDSAcceleratorComponent.hh.

Referenced by AngledInputFace(), InputFaceNormal(), and SetInputFaceNormal().

◆ lengthSafety

G4double BDSAcceleratorComponent::lengthSafety = -1
staticprotected

◆ lengthSafetyLarge

G4double BDSAcceleratorComponent::lengthSafetyLarge = 0
static

A larger length safety that can be used where tracking accuracy isn't required or more tolerant geometry is required (1um).

Definition at line 213 of file BDSAcceleratorComponent.hh.

Referenced by BDSAcceleratorComponent(), BDSUndulator::Build(), BDSUndulator::BuildContainerLogicalVolume(), BDSDump::BuildContainerLogicalVolume(), BDSCollimatorRectangular::BuildInnerCollimator(), and BDSShield::BuildShield().

◆ name

const G4String BDSAcceleratorComponent::name
protected

Const protected member variable that may not be changed by derived classes.

Definition at line 251 of file BDSAcceleratorComponent.hh.

Referenced by BDSScreen::AddScreenLayer(), BDSAcceleratorComponent(), BDSCollimator::Build(), BDSCollimatorCrystal::Build(), BDSDegrader::Build(), BDSDrift::Build(), BDSUndulator::Build(), BDSWireScanner::Build(), BDSCavityElement::Build(), BDSCollimatorJaw::Build(), BDSMagnet::BuildBeampipe(), BDSShield::BuildBeamPipe(), BDSCollimator::BuildContainerLogicalVolume(), BDSDegrader::BuildContainerLogicalVolume(), BDSElement::BuildContainerLogicalVolume(), BDSLaserWire::BuildContainerLogicalVolume(), BDSMagnet::BuildContainerLogicalVolume(), BDSShield::BuildContainerLogicalVolume(), BDSTeleporter::BuildContainerLogicalVolume(), BDSTerminator::BuildContainerLogicalVolume(), BDSUndulator::BuildContainerLogicalVolume(), BDSWireScanner::BuildContainerLogicalVolume(), BDSCollimatorJaw::BuildContainerLogicalVolume(), BDSDump::BuildContainerLogicalVolume(), BDSCollimatorElliptical::BuildInnerCollimator(), BDSCollimatorRectangular::BuildInnerCollimator(), BDSShield::BuildShield(), BDSCollimator::CheckParameters(), BDSCollimatorElliptical::CheckParameters(), BDSCollimatorJaw::CheckParameters(), BDSTunnelFactoryElliptical::CreateTunnelSection(), BDSTunnelFactoryRectAboveGround::CreateTunnelSection(), BDSTunnelFactoryElliptical::CreateTunnelSectionAngled(), BDSTunnelFactoryRectAboveGround::CreateTunnelSectionAngled(), GetName(), and BDSMagnet::PlaceComponents().

◆ outputFaceNormal

G4ThreeVector BDSAcceleratorComponent::outputFaceNormal
private

Output face unit normal vector in outgoing reference coordinate frame.

Definition at line 312 of file BDSAcceleratorComponent.hh.

Referenced by AngledOutputFace(), OutputFaceNormal(), and SetOutputFaceNormal().

◆ region

G4String BDSAcceleratorComponent::region
protected

Protected member variable that can be modified by derived classes.

Definition at line 259 of file BDSAcceleratorComponent.hh.

Referenced by GetRegion(), and SetRegion().

◆ sensitiveOuter

G4bool BDSAcceleratorComponent::sensitiveOuter = true
staticprotected

◆ sensitiveVacuum

G4bool BDSAcceleratorComponent::sensitiveVacuum = false
staticprotected

Useful variable often used in construction.

Definition at line 271 of file BDSAcceleratorComponent.hh.

Referenced by BDSAcceleratorComponent(), BDSCollimator::Build(), BDSCollimatorJaw::Build(), and BDSLaserWire::BuildContainerLogicalVolume().

◆ type

const G4String BDSAcceleratorComponent::type
protected

Const protected member variable that may not be changed by derived classes.

Definition at line 253 of file BDSAcceleratorComponent.hh.

Referenced by GetType().

◆ userLimits

G4UserLimits* BDSAcceleratorComponent::userLimits
protected

◆ worldMaterial

G4Material * BDSAcceleratorComponent::worldMaterial = nullptr
staticprotected

Useful variable often used in construction.

Definition at line 268 of file BDSAcceleratorComponent.hh.

Referenced by BDSAcceleratorComponent(), and BDSShield::BuildContainerLogicalVolume().


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