BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Public Member Functions | Protected Attributes
BDSBeamPipe Class Reference

A holder class for a piece of beam pipe geometry. More...

#include <BDSBeamPipe.hh>

Inheritance diagram for BDSBeamPipe:
Inheritance graph
Collaboration diagram for BDSBeamPipe:
Collaboration graph

Public Member Functions

 BDSBeamPipe (G4VSolid *containerSolidIn, G4LogicalVolume *containerLVIn, BDSExtent extentIn, G4VSolid *containerSubtractionSolidIn, G4LogicalVolume *vacuumLVIn, G4bool containerIsCircularIn=false, G4double containerRadiusIn=0.0, G4ThreeVector inputFaceNormalIn=G4ThreeVector(0, 0,-1), G4ThreeVector outputFaceNormalIn=G4ThreeVector(0, 0, 1))
 
G4VSolid * GetContainerSubtractionSolid () const
 default destructor sufficient as G4 manages solids and LVs More...
 
G4LogicalVolume * GetVacuumLogicalVolume () const
 Access the vacuum volume to set fields and limits. More...
 
G4bool ContainerIsCircular () const
 
G4double GetContainerRadius () const
 If it is circular, we need the radius. More...
 
std::set< G4LogicalVolume * > GetVolumesForField () const
 
G4ThreeVector InputFaceNormal () const
 Accessor. More...
 
G4ThreeVector OutputFaceNormal () const
 Accessor. 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...
 
virtual G4String GetName () const
 Accessor - see member for more info. 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...
 

Protected Attributes

G4VSolid * containerSubtractionSolid
 
G4LogicalVolume * vacuumLogicalVolume
 
G4bool containerIsCircular
 
G4double containerRadius
 
G4ThreeVector inputFaceNormal
 
G4ThreeVector outputFaceNormal
 
- 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...
 

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

A holder class for a piece of beam pipe geometry.

This does not implement the construction of the beampipe but merely holds all relevant objects and information.

Face normals are ensured to be unit.

Author
Laurie Nevay

Definition at line 44 of file BDSBeamPipe.hh.

Constructor & Destructor Documentation

◆ BDSBeamPipe()

BDSBeamPipe::BDSBeamPipe ( G4VSolid *  containerSolidIn,
G4LogicalVolume *  containerLVIn,
BDSExtent  extentIn,
G4VSolid *  containerSubtractionSolidIn,
G4LogicalVolume *  vacuumLVIn,
G4bool  containerIsCircularIn = false,
G4double  containerRadiusIn = 0.0,
G4ThreeVector  inputFaceNormalIn = G4ThreeVector(0,0,-1),
G4ThreeVector  outputFaceNormalIn = G4ThreeVector(0,0, 1) 
)

constructor has BDSGeometryComponent members first, then everything extra for this derived class

Parameters
containerSolidInContainer solid.
containerLVInContainer logical volume.
containerIsCircularInWhether the container is circular.

Definition at line 25 of file BDSBeamPipe.cc.

◆ ~BDSBeamPipe()

BDSBeamPipe::~BDSBeamPipe ( )
virtual

Definition at line 43 of file BDSBeamPipe.cc.

Member Function Documentation

◆ ContainerIsCircular()

G4bool BDSBeamPipe::ContainerIsCircular ( ) const
inline

Flag to tell whether the parent volume needn't use a subtraction solid and can simply use a G4Tubs for example

Definition at line 68 of file BDSBeamPipe.hh.

Referenced by BDSMagnetOuterFactoryCylindrical::CreateCylindricalSolids(), BDSMagnetOuterFactoryLHC::CreateLHCDipole(), BDSMagnetOuterFactoryLHC::CreateQuadrupole(), and BDSMagnetOuterFactoryCylindrical::TestInputParameters().

Here is the caller graph for this function:

◆ GetContainerRadius()

G4double BDSBeamPipe::GetContainerRadius ( ) const
inline

◆ GetContainerSubtractionSolid()

G4VSolid * BDSBeamPipe::GetContainerSubtractionSolid ( ) const
inline

default destructor sufficient as G4 manages solids and LVs

Access a solid for beampipe subtraction - note this is typically longer than the actual beampipe for unambiguous subtraction

Definition at line 63 of file BDSBeamPipe.hh.

Referenced by BDSMagnetOuterFactoryCylindrical::CreateCylindricalSolids(), BDSMagnetOuterFactoryLHC::CreateLHCDipole(), and BDSMagnetOuterFactoryLHC::CreateQuadrupole().

Here is the caller graph for this function:

◆ GetVacuumLogicalVolume()

G4LogicalVolume * BDSBeamPipe::GetVacuumLogicalVolume ( ) const
inline

Access the vacuum volume to set fields and limits.

Definition at line 65 of file BDSBeamPipe.hh.

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

Here is the caller graph for this function:

◆ GetVolumesForField()

std::set< G4LogicalVolume * > BDSBeamPipe::GetVolumesForField ( ) const

Get all volumes except the container logical volume as this is the optimal set of volumes for putting fields on.

Definition at line 50 of file BDSBeamPipe.cc.

References BDSGeometryComponent::GetAllLogicalVolumes(), BDSGeometryComponent::GetContainerLogicalVolume(), and BDS::IsFinite().

Referenced by BDSMagnet::BuildVacuumField().

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

◆ InputFaceNormal()

G4ThreeVector BDSBeamPipe::InputFaceNormal ( ) const
inline

◆ OutputFaceNormal()

G4ThreeVector BDSBeamPipe::OutputFaceNormal ( ) const
inline

Field Documentation

◆ containerIsCircular

G4bool BDSBeamPipe::containerIsCircular
protected

Definition at line 84 of file BDSBeamPipe.hh.

◆ containerRadius

G4double BDSBeamPipe::containerRadius
protected

Definition at line 85 of file BDSBeamPipe.hh.

◆ containerSubtractionSolid

G4VSolid* BDSBeamPipe::containerSubtractionSolid
protected

Definition at line 82 of file BDSBeamPipe.hh.

◆ inputFaceNormal

G4ThreeVector BDSBeamPipe::inputFaceNormal
protected

Definition at line 86 of file BDSBeamPipe.hh.

◆ outputFaceNormal

G4ThreeVector BDSBeamPipe::outputFaceNormal
protected

Definition at line 87 of file BDSBeamPipe.hh.

◆ vacuumLogicalVolume

G4LogicalVolume* BDSBeamPipe::vacuumLogicalVolume
protected

Definition at line 83 of file BDSBeamPipe.hh.


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