BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Undulator. More...
#include <BDSUndulator.hh>
Public Member Functions | |
BDSUndulator (const G4String &nameIn, G4double lengthIn, G4double periodIn, G4double magnetHeightIn, G4double magnetWidthIn, G4double undulatorGapIn, BDSBeamPipeInfo *beamPipeInfoIn, BDSFieldInfo *vacuumFieldInfoIn, BDSFieldInfo *outerFieldInfoIn, const G4String &materialIn="iron") | |
![]() | |
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 ®ionIn) |
Set the region name for this component. | |
void | SetField (BDSFieldInfo *fieldInfoIn) |
Set the field definition for the whole component. | |
virtual void | SetMinimumKineticEnergy (G4double) |
virtual G4String | GetName () const |
The name of the component without modification. | |
G4double | GetAngle () const |
G4double | Sagitta () const |
Calculate the sagitta - ie the distance between the chord and the arc at the centre. | |
G4String | GetType () const |
Get a string describing the type of the component. | |
G4String | GetRegion () const |
Get the region name for this component. | |
virtual G4bool | HasAField () const |
Whether this component has a field or not (ie is active). Implicit cast of pointer to bool. | |
virtual void | SetFieldUsePlacementWorldTransform () |
BDSBeamPipeInfo * | GetBeamPipeInfo () const |
virtual G4String | Material () const |
Return the name of a material associated with the component - ie the primary material. | |
virtual std::set< G4LogicalVolume * > | GetAcceleratorVacuumLogicalVolumes () const |
Access the 'vacuum' volume(s) in this component if any. Default is empty set. | |
virtual std::set< G4LogicalVolume * > | GetAcceleratorMaterialLogicalVolumes () const |
Return a set of logical volumes excluding the ones in the 'vacuum' set. | |
void | IncrementCopyNumber () |
Increment (+1) the number of times this component has been copied. | |
G4int | GetCopyNumber () const |
Get the number of times this component has been copied. | |
BDSSimpleComponent * | EndPieceBefore () const |
BDSSimpleComponent * | EndPieceAfter () const |
virtual void | SetBiasVacuumList (const std::list< std::string > &biasVacuumListIn) |
Copy the bias list to this element. | |
virtual void | SetBiasMaterialList (const std::list< std::string > &biasMaterialListIn) |
Copy the bias list to this element. | |
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. | |
G4bool | AngledOutputFace () const |
Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory. | |
std::list< std::string > | GetBiasVacuumList () const |
Access the bias list copied from parser. | |
std::list< std::string > | GetBiasMaterialList () const |
Access the bias list copied from parser. | |
virtual void | SetInputFaceNormal (const G4ThreeVector &input) |
Allow updating of face normals. Virtual so derived class may apply it to daughters. | |
virtual void | SetOutputFaceNormal (const G4ThreeVector &output) |
Allow updating of face normals. Virtual so derived class may apply it to daughters. | |
![]() | |
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) | |
BDSGeometryComponent & | operator= (const BDSGeometryComponent &)=delete |
Assignment operator not used. | |
G4bool | ContainerIsAssembly () const |
Whether the container is an assembly. If not, it's a logical volume. | |
void | SetPlacementOffset (const G4ThreeVector &offsetIn) |
Set the offset from 0,0,0 that the object should ideally be placed in its parent. | |
G4ThreeVector | GetExtentPositive () const |
Get the extent of the object in the positive direction in all dimensions. | |
G4ThreeVector | GetExtentNegative () const |
Get the extent of the object in the negative direction in all dimensions. | |
void | InheritExtents (BDSGeometryComponent const *const anotherComponent) |
Update the extents of this object with those of another object. | |
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. | |
void | RegisterLogicalVolume (G4LogicalVolume *logicalVolume) |
void | RegisterLogicalVolume (const std::set< G4LogicalVolume * > &localVolumes) |
Apply RegisterLogicalVolume() to a set of logical volumes. | |
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. | |
void | RegisterPhysicalVolume (G4VPhysicalVolume *physicalVolume) |
void | RegisterPhysicalVolume (const std::set< G4VPhysicalVolume * > &physicalVolumes) |
Apply RegisterPhysicalVolume() to a set of physical volumes. | |
void | RegisterRotationMatrix (G4RotationMatrix *rotationMatrix) |
void | RegisterRotationMatrix (const std::set< G4RotationMatrix * > &rotationMatrices) |
Apply Register RotationMatrix() to a set of rotation matrices. | |
void | RegisterVisAttributes (G4VisAttributes *visAttribute) |
void | RegisterVisAttributes (const std::set< G4VisAttributes * > &visAttributes) |
Apply RegisterVisAttribute() to a set of visualisation attributes. | |
void | RegisterUserLimits (G4UserLimits *userLimit) |
void | RegisterUserLimits (const std::set< G4UserLimits * > &userLimits) |
Apply RegisterUserLimit to a set of user limits. | |
void | InheritObjects (BDSGeometryComponent *component) |
virtual std::set< G4LogicalVolume * > | GetAllLogicalVolumes () const |
Access all logical volumes belonging to this component. | |
virtual std::set< G4LogicalVolume * > | GetAllBiasingVolumes () const |
Return all logical volumes that should be used for biasing minus any that are in the excluded set. | |
virtual std::map< G4LogicalVolume *, BDSSDType > | GetAllSensitiveVolumes () const |
Access all sensitive volumes belonging to this component. | |
void | MakeAllVolumesSensitive (BDSSDType stype=BDSSDType::energydep) |
virtual void | AttachSensitiveDetectors () |
Attach a sensitive detector class to all registered sensitive volumes in this component. | |
virtual void | ExcludeLogicalVolumeFromBiasing (G4LogicalVolume *lv) |
void | StripOuterAndMakeAssemblyVolume () |
Change from a container logical volume to an assembly volume. | |
G4VSolid * | GetContainerSolid () const |
Accessor - see member for more info. | |
G4LogicalVolume * | GetContainerLogicalVolume () const |
Accessor - see member for more info. | |
G4AssemblyVolume * | GetContainerAssemblyVolume () const |
Accessor - see member for more info. | |
G4Transform3D | GetPlacementTransform () const |
Accessor - see member for more info. | |
G4ThreeVector | GetPlacementOffset () const |
Accessor - see member for more info. | |
G4RotationMatrix * | GetPlacementRotation () const |
Accessor - see member for more info. | |
BDSExtent | GetExtent () const |
Accessor - see member for more info. | |
BDSExtent | GetInnerExtent () const |
Accessor - see member for more info. | |
std::pair< G4double, G4double > | GetExtentX () const |
Accessor - see member for more info. | |
std::pair< G4double, G4double > | GetExtentY () const |
Accessor - see member for more info. | |
std::pair< G4double, G4double > | GetExtentZ () const |
Accessor - see member for more info. | |
std::pair< G4double, G4double > | GetInnerExtentX () const |
Accessor - see member for more info. | |
std::pair< G4double, G4double > | GetInnerExtentY () const |
Accessor - see member for more info. | |
std::pair< G4double, G4double > | GetInnerExtentZ () const |
Accessor - see member for more info. | |
virtual std::set< G4VPhysicalVolume * > | GetAllPhysicalVolumes () const |
Accessor - see member for more info. | |
virtual std::set< G4RotationMatrix * > | GetAllRotationMatrices () const |
Accessor - see member for more info. | |
virtual std::set< G4VisAttributes * > | GetAllVisAttributes () const |
Accessor - see member for more info. | |
virtual std::set< G4UserLimits * > | GetAllUserLimits () const |
Accessor - see member for more info. | |
virtual std::set< BDSGeometryComponent * > | GetAllDaughters () const |
Accessor - see member for more info. | |
virtual std::set< G4VSolid * > | GetAllSolids () const |
Accessor - see member for more info. | |
void | SetExtent (const BDSExtent &extIn) |
Set extent. | |
void | SetInnerExtent (const BDSExtent &extIn) |
Set extent. | |
Protected Member Functions | |
virtual void | Build () |
virtual void | BuildContainerLogicalVolume () |
![]() | |
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. | |
Protected Attributes | |
BDSFieldInfo * | vacuumFieldInfo |
BDSFieldInfo * | outerFieldInfo |
const G4double | undulatorPeriod |
const G4double | horizontalWidth |
Element width (and height) | |
G4double | undulatorMagnetHeight |
Full magnet box height. | |
G4double | undulatorGap |
Full undulator gap. | |
G4int | numMagnets |
Total number of magnets (1 undulator period is 2 magnets) | |
G4String | material |
Undulator magnet material. | |
![]() | |
BDSBeamPipeInfo * | beamPipeInfo |
Optional beam pipe recipe that is written out to the survey if it exists. | |
std::set< G4LogicalVolume * > | acceleratorVacuumLV |
A set of logical volumes we classify as 'vacuum' for biasing purposes. | |
BDSSimpleComponent * | endPieceBefore |
BDSSimpleComponent * | endPieceAfter |
G4UserLimits * | userLimits |
Cache of user limits. | |
BDSFieldInfo * | fieldInfo |
Recipe for field that could overlay this whole component. | |
const G4String | name |
Const protected member variable that may not be changed by derived classes. | |
const G4double | arcLength |
Const protected member variable that may not be changed by derived classes. | |
const G4String | type |
Const protected member variable that may not be changed by derived classes. | |
G4double | chordLength |
Protected member variable that can be modified by derived classes. | |
G4double | angle |
Protected member variable that can be modified by derived classes. | |
G4String | region |
Protected member variable that can be modified by derived classes. | |
![]() | |
G4bool | containerIsAssembly |
True if the 'container' is really an assembly; false if an LV. | |
G4VSolid * | containerSolid |
G4LogicalVolume * | containerLogicalVolume |
G4AssemblyVolume * | containerAssembly |
BDSExtent | outerExtent |
BDSExtent | innerExtent |
std::set< BDSGeometryComponent * > | allDaughters |
registry of all daughter geometry components | |
std::set< G4VSolid * > | allSolids |
registry of all solids belonging to this component | |
std::set< G4LogicalVolume * > | allLogicalVolumes |
std::map< G4LogicalVolume *, BDSSDType > | sensitivity |
G4bool | overrideSensitivity |
std::set< G4VPhysicalVolume * > | allPhysicalVolumes |
registry of all physical volumes belonging to this component | |
std::set< G4RotationMatrix * > | allRotationMatrices |
registry of all rotation matrices belonging to this component | |
std::set< G4VisAttributes * > | allVisAttributes |
registry of all visualisation attributes belonging to this component | |
std::set< G4UserLimits * > | allUserLimits |
registry of all user limits belonging to this component | |
G4ThreeVector | placementOffset |
G4RotationMatrix * | placementRotation |
std::set< G4LogicalVolume * > * | lvsExcludedFromBiasing |
Volumes that should not be included when return GetAllBiasingVolumes(). | |
Additional Inherited Members | |
![]() | |
static void | AttachUserLimitsToAssembly (G4AssemblyVolume *av, G4UserLimits *ul) |
Utility function to apply user limits to an assembly volume as there's not interface. | |
![]() | |
static G4double | lengthSafetyLarge = 0 |
![]() | |
static G4double | lengthSafety = -1 |
Useful variable often used in construction. | |
static G4Material * | emptyMaterial = nullptr |
Useful variable often used in construction. | |
static G4Material * | worldMaterial = nullptr |
Useful variable often used in construction. | |
static G4bool | checkOverlaps = false |
Useful variable often used in construction. | |
static G4bool | sensitiveOuter = true |
Useful variable often used in construction. | |
static G4bool | sensitiveVacuum = false |
Useful variable often used in construction. | |
static G4VisAttributes * | containerVisAttr = nullptr |
Useful variable often used in construction. | |
Undulator.
Requires a valid BDSBeamPipeInfo* instance.
Definition at line 38 of file BDSUndulator.hh.
BDSUndulator::BDSUndulator | ( | const G4String & | nameIn, |
G4double | lengthIn, | ||
G4double | periodIn, | ||
G4double | magnetHeightIn, | ||
G4double | magnetWidthIn, | ||
G4double | undulatorGapIn, | ||
BDSBeamPipeInfo * | beamPipeInfoIn, | ||
BDSFieldInfo * | vacuumFieldInfoIn, | ||
BDSFieldInfo * | outerFieldInfoIn, | ||
const G4String & | materialIn = "iron" |
||
) |
Definition at line 42 of file BDSUndulator.cc.
|
virtual |
Definition at line 73 of file BDSUndulator.cc.
|
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 from BDSAcceleratorComponent.
Definition at line 128 of file BDSUndulator.cc.
References BDSAcceleratorComponent::beamPipeInfo, BDSAcceleratorComponent::Build(), BDSAcceleratorComponent::checkOverlaps, BDSAcceleratorComponent::chordLength, BDSGeometryComponent::GetContainerLogicalVolume(), BDSMaterials::GetMaterial(), BDSBeamPipe::GetVacuumLogicalVolume(), horizontalWidth, BDSBeamPipe::InputFaceNormal(), BDSBeamPipeFactory::Instance(), BDSColours::Instance(), BDSFieldBuilder::Instance(), BDSMaterials::Instance(), BDS::IsFinite(), BDSAcceleratorComponent::lengthSafetyLarge, BDSFieldInfo::MagnetStrength(), material, BDSAcceleratorComponent::name, numMagnets, BDSBeamPipe::OutputFaceNormal(), BDSGeometryComponent::RegisterDaughter(), BDSFieldBuilder::RegisterFieldForConstruction(), BDSGeometryComponent::RegisterLogicalVolume(), BDSGeometryComponent::RegisterPhysicalVolume(), BDSGeometryComponent::RegisterSensitiveVolume(), BDSGeometryComponent::RegisterSolid(), BDSGeometryComponent::RegisterVisAttributes(), BDSAcceleratorComponent::sensitiveOuter, BDSAcceleratorComponent::SetAcceleratorVacuumLogicalVolume(), BDSAcceleratorComponent::SetInputFaceNormal(), BDSAcceleratorComponent::SetOutputFaceNormal(), undulatorGap, and undulatorMagnetHeight.
|
protectedvirtual |
Build the container solid and logical volume that all parts of the component will contained within - must be provided by derived class.
Implements BDSAcceleratorComponent.
Definition at line 76 of file BDSUndulator.cc.
References BDSAcceleratorComponent::beamPipeInfo, BDSBeamPipeInfo::beamPipeThickness, BDSAcceleratorComponent::chordLength, BDSExtent::DY(), BDSAcceleratorComponent::emptyMaterial, BDSBeamPipeInfo::Extent(), horizontalWidth, BDSGlobalConstants::Instance(), BDS::IsFinite(), BDSAcceleratorComponent::lengthSafetyLarge, BDSAcceleratorComponent::name, numMagnets, BDSGeometryComponent::SetExtent(), undulatorGap, and undulatorMagnetHeight.
|
protected |
Element width (and height)
Definition at line 62 of file BDSUndulator.hh.
Referenced by Build(), and BuildContainerLogicalVolume().
|
protected |
|
protected |
Total number of magnets (1 undulator period is 2 magnets)
Definition at line 65 of file BDSUndulator.hh.
Referenced by Build(), and BuildContainerLogicalVolume().
|
protected |
Definition at line 60 of file BDSUndulator.hh.
|
protected |
Full undulator gap.
Definition at line 64 of file BDSUndulator.hh.
Referenced by Build(), and BuildContainerLogicalVolume().
|
protected |
Full magnet box height.
Definition at line 63 of file BDSUndulator.hh.
Referenced by Build(), and BuildContainerLogicalVolume().
|
protected |
Definition at line 61 of file BDSUndulator.hh.
|
protected |
Definition at line 59 of file BDSUndulator.hh.