BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
A vector of BDSBeamlineElement instances - a beamline. More...
#include <BDSBeamline.hh>
Public Member Functions | |
BDSBeamline & | operator= (const BDSBeamline &)=delete |
assignment and copy constructor not implemented nor used. | |
BDSBeamline (BDSBeamline &)=delete | |
BDSBeamline (const G4ThreeVector &initialGlobalPosition=G4ThreeVector(0, 0, 0), G4RotationMatrix *initialGlobalRotation=nullptr, G4double initialSIn=0.0) | |
BDSBeamline (G4Transform3D initialTransform, G4double initialSIn=0.0) | |
Constructor with transform instance that uses other constructor. | |
void | AddComponent (BDSAcceleratorComponent *component, BDSTiltOffset *tiltOffset=nullptr, BDSSamplerInfo *samplerInfo=nullptr) |
void | ApplyTransform3D (BDSTransform3D *component) |
void | AddBeamlineElement (BDSBeamlineElement *element) |
void | PrintAllComponents (std::ostream &out) const |
void | UpdateExtents (BDSBeamlineElement *element) |
const BDSBeamlineElement * | at (int iElement) const |
Return a reference to the element at i. | |
const BDSBeamlineElement * | GetFirstItem () const |
Return a reference to the first element. | |
const BDSBeamlineElement * | GetLastItem () const |
Return a reference to the last element. | |
const BDSBeamlineElement * | GetElement (G4String acceleratorComponentName, G4int i=0) const |
Get the ith placement of an element in the beam line. Returns null pointer if not found. | |
G4Transform3D | GetTransformForElement (const G4String &acceleratorComponentName, G4int i=0) const |
G4double | GetTotalChordLength () const |
Get the total length of the beamline - the sum of the chord length of each element. | |
G4double | GetTotalArcLength () const |
Get the total ARC length for the beamline - ie the maximum s position. | |
G4double | GetTotalAngle () const |
Get the total angle of the beamline. | |
G4bool | ElementAnglesSumToCircle () const |
Check if the sum of the angle of all elements is two pi. | |
BeamlineVector::size_type | size () const |
Get the number of elements. | |
G4ThreeVector | GetMaximumExtentPositive () const |
Get the maximum positive extent in all dimensions | |
G4ThreeVector | GetMaximumExtentNegative () const |
Get the maximum negative extent in all dimensions. | |
G4ThreeVector | GetMaximumExtentAbsolute () const |
Get the maximum extent absolute in each dimension. | |
BDSExtentGlobal | GetExtentGlobal () const |
Get the global extents for this beamline. | |
G4Transform3D | GetGlobalEuclideanTransform (G4double s, G4double x=0, G4double y=0, G4int *indexOfFoundElement=nullptr) const |
const BDSBeamlineElement * | GetElementFromGlobalS (G4double S, G4int *indexOfFoundElement=nullptr) const |
Return the element in this beam line according to a given s coordinate. | |
const_iterator | FindFromS (G4double S) const |
Returns an iterator to the beamline element at s. | |
std::vector< G4double > | GetEdgeSPositions () const |
G4double | GetSMinimum () const |
G4double | GetSMaximum () const |
const BDSBeamlineElement * | GetPrevious (const BDSBeamlineElement *element) const |
const BDSBeamlineElement * | GetNext (const BDSBeamlineElement *element) const |
const BDSBeamlineElement * | GetPrevious (G4int index) const |
const BDSBeamlineElement * | GetNext (G4int index) const |
BDSBeamlineElement * | front () const |
Return a reference to the first element. | |
BDSBeamlineElement * | back () const |
Return a reference to the last element. | |
void | PrintMemoryConsumption () const |
BDSBeamlineElement * | ProvideEndPieceElementBefore (BDSSimpleComponent *endPiece, G4int index) const |
BDSBeamlineElement * | ProvideEndPieceElementAfter (BDSSimpleComponent *endPiece, G4int index, G4bool flip=false) const |
G4bool | IndexOK (G4int index) const |
Whether the supplied index will lie within the beam line vector. | |
std::vector< G4int > | GetIndicesOfElementsOfType (const G4String &type) const |
Return vector of indices for this beam line where element of type name 'type' is found. | |
std::vector< G4int > | GetIndicesOfElementsOfType (const std::set< G4String > &types) const |
std::vector< G4int > | GetIndicesOfCollimators () const |
Return indices in order of ecol, rcol, jcol and crystalcol elements. | |
Static Public Member Functions | |
static G4double | PaddingLength () |
Access the padding length between each element added to the beamline. | |
Private Types | |
typedef std::vector< BDSBeamlineElement * > | BeamlineVector |
Typedefs up first so we can declare public iterators. | |
Private Member Functions | |
void | AddSingleComponent (BDSAcceleratorComponent *component, BDSTiltOffset *tiltOffset=nullptr, BDSSamplerInfo *samplerInfo=nullptr) |
void | RegisterElement (BDSBeamlineElement *element) |
Private Attributes | |
BeamlineVector | beamline |
Vector of beam line elements - the data. | |
G4double | sInitial |
Cache the initial S so we can tell if a requested S is too low. | |
G4double | sMaximum |
G4double | totalChordLength |
Sum of all chord lengths. | |
G4double | totalArcLength |
Sum of all arc lengths. | |
G4double | totalAngle |
Sum of all angles. | |
G4ThreeVector | maximumExtentPositive |
maximum extent in the positive coordinates in each dimension | |
G4ThreeVector | maximumExtentNegative |
maximum extent in the negative coordinates in each dimension | |
G4RotationMatrix * | previousReferenceRotationEnd |
Current reference rotation at the end of the previous element. | |
G4ThreeVector | previousReferencePositionEnd |
Current reference position at the end of the previous element. | |
G4double | previousSPositionEnd |
Current s coordinate at the end of the previous element. | |
G4bool | transformHasJustBeenApplied |
std::map< G4String, BDSBeamlineElement * > | components |
Map of objects by placement name stored in this beam line. | |
std::vector< G4double > | sEnd |
Static Private Attributes | |
static G4double | paddingLength = -1 |
The gap added for padding between each component. | |
Friends | |
std::ostream & | operator<< (std::ostream &out, BDSBeamline const &bl) |
output stream | |
typedef BeamlineVector::iterator | iterator |
Iterator mechanics. | |
typedef BeamlineVector::const_iterator | const_iterator |
Iterator mechanics. | |
typedef BeamlineVector::reverse_iterator | reverse_iterator |
Iterator mechanics. | |
typedef BeamlineVector::const_reverse_iterator | const_reverse_iterator |
Iterator mechanics. | |
iterator | begin () |
Iterator mechanics. | |
iterator | end () |
Iterator mechanics. | |
const_iterator | begin () const |
Iterator mechanics. | |
const_iterator | end () const |
Iterator mechanics. | |
reverse_iterator | rbegin () |
Iterator mechanics. | |
reverse_iterator | rend () |
Iterator mechanics. | |
const_reverse_iterator | rbegin () const |
Iterator mechanics. | |
const_reverse_iterator | rend () const |
Iterator mechanics. | |
G4bool | empty () const |
Iterator mechanics. | |
A vector of BDSBeamlineElement instances - a beamline.
This will calculate and construct a beamline as BDSAcceleratorComponents are added in sequence - ie it calculates their placement positions and rotations with respect to the start of the beamline as well as their s position in curvilinear coordinates.
Note, this is not a singleton as geometry hierarchy can be introduced by placing beamline components inside parent volumes and therefore creating a new beamline of parents. It can also be used to create multiple beam lines.
Definition at line 60 of file BDSBeamline.hh.
|
private |
Typedefs up first so we can declare public iterators.
Definition at line 64 of file BDSBeamline.hh.
typedef BeamlineVector::const_iterator BDSBeamline::const_iterator |
Iterator mechanics.
Definition at line 164 of file BDSBeamline.hh.
typedef BeamlineVector::const_reverse_iterator BDSBeamline::const_reverse_iterator |
Iterator mechanics.
Definition at line 166 of file BDSBeamline.hh.
typedef BeamlineVector::iterator BDSBeamline::iterator |
Iterator mechanics.
Definition at line 163 of file BDSBeamline.hh.
typedef BeamlineVector::reverse_iterator BDSBeamline::reverse_iterator |
Iterator mechanics.
Definition at line 165 of file BDSBeamline.hh.
BDSBeamline::BDSBeamline | ( | const G4ThreeVector & | initialGlobalPosition = G4ThreeVector(0,0,0) , |
G4RotationMatrix * | initialGlobalRotation = nullptr , |
||
G4double | initialSIn = 0.0 |
||
) |
Versatile basic constructor that allows a finite position and rotation to be applied at the beginning of the beamline in global coordinates. Remember the maximum extents of the beamline will also be displaced. The default constructor is in effect achieved via defaults.
Definition at line 50 of file BDSBeamline.cc.
References BDSGlobalConstants::Instance(), maximumExtentNegative, maximumExtentPositive, paddingLength, and previousReferenceRotationEnd.
|
explicit |
Constructor with transform instance that uses other constructor.
Definition at line 84 of file BDSBeamline.cc.
BDSBeamline::~BDSBeamline | ( | ) |
Definition at line 91 of file BDSBeamline.cc.
void BDSBeamline::AddBeamlineElement | ( | BDSBeamlineElement * | element | ) |
Add a pre-assembled beam line element. In this case, the coordinates will have been calculated external to this class and as such, it's the responsibility of the developer to make sure the coordinates are correct and do not cause overlaps. This will be useful for tunnel construction for example or for a non-contiguous beamline. Subsequent components added via the AddComponent() method will be appended in the usual way to the end coordinates of this element.
Definition at line 574 of file BDSBeamline.cc.
References beamline, BDSBeamlineElement::GetAcceleratorComponent(), RegisterElement(), and UpdateExtents().
Referenced by BDS::BuildBLMs(), BDSCurvilinearBuilder::BuildCurvilinearBeamLine1To1(), BDSCurvilinearBuilder::BuildCurvilinearBridgeBeamLine(), BDS::BuildEndPieceBeamline(), BDS::BuildPlacementGeometry(), and BDSTunnelBuilder::BuildTunnelSections().
void BDSBeamline::AddComponent | ( | BDSAcceleratorComponent * | component, |
BDSTiltOffset * | tiltOffset = nullptr , |
||
BDSSamplerInfo * | samplerInfo = nullptr |
||
) |
Add a component, but check to see if it can be dynamically upcast to a line in which case, loop over it and apply AddSingleComponent(BDSAcceleratorComponent* component) to each component Returns vector of components added
Definition at line 124 of file BDSBeamline.cc.
References AddSingleComponent(), back(), BDSBeamlineElement::GetReferencePositionEnd(), BDSBeamlineElement::GetReferencePositionStart(), BDSBeamlineElement::GetReferenceRotationEnd(), BDSBeamlineElement::GetReferenceRotationStart(), BDSBeamlineElement::GetSamplerInfo(), BDSOutput::InvalidSamplerName(), BDSOutput::PrintProtectedNames(), and BDSBeamlineElement::UpdateSamplerPlacementTransform().
Referenced by BDSDetectorConstruction::BuildBeamline(), and BDSTunnelBuilder::BuildTunnelSections().
|
private |
Add a single component and calculate its position and rotation with respect to the beginning of the beamline. Returns pointer to the component added.
Definition at line 209 of file BDSBeamline.cc.
References ApplyTransform3D(), back(), beamline, BDSSamplerPlane::ChordLength(), empty(), BDSBeamlineElement::GetAcceleratorComponent(), BDSAcceleratorComponent::GetAngle(), BDSAcceleratorComponent::GetArcLength(), BDSAcceleratorComponent::GetChordLength(), BDSBeamlineElement::GetChordLength(), BDSGeometryComponent::GetExtent(), BDSGeometryComponent::GetExtentNegative(), BDSGeometryComponent::GetExtentPositive(), BDSAcceleratorComponent::GetName(), BDSTiltOffset::GetOffset(), BDSGeometryComponent::GetPlacementOffset(), GetPrevious(), BDSBeamlineElement::GetReferencePositionEnd(), BDSBeamlineElement::GetReferenceRotationEnd(), BDSBeamlineElement::GetSPositionEnd(), BDSTiltOffset::GetTilt(), BDSAcceleratorComponent::GetType(), BDSBeamlineElement::GetType(), BDSTiltOffset::HasFiniteOffset(), BDSTiltOffset::HasFiniteTilt(), BDSAcceleratorComponent::InputFaceNormal(), BDS::IsFinite(), BDSAcceleratorComponent::OutputFaceNormal(), paddingLength, previousReferencePositionEnd, previousReferenceRotationEnd, previousSPositionEnd, RegisterElement(), sEnd, totalAngle, totalArcLength, totalChordLength, transformHasJustBeenApplied, UpdateExtents(), and BDS::WillIntersect().
Referenced by AddComponent().
void BDSBeamline::ApplyTransform3D | ( | BDSTransform3D * | component | ) |
Apply a Transform3D rotation and translation to the reference coordinates. Special method for the special case of unique component Transform3D. Modifies the end rotation and position of the last element in the beamline. Note this applies the offset dx,dy,dz first, then the rotation of coordinates
Definition at line 542 of file BDSBeamline.cc.
References back(), empty(), BDSAcceleratorComponent::GetName(), BDSBeamlineElement::GetReferencePositionEnd(), BDSBeamlineElement::GetReferenceRotationEnd(), previousReferencePositionEnd, and previousReferenceRotationEnd.
Referenced by AddSingleComponent().
|
inline |
Return a reference to the element at i.
Definition at line 120 of file BDSBeamline.hh.
References beamline.
Referenced by GetGlobalEuclideanTransform(), BDSOutputStructures::PrepareCavityInformation(), and BDSOutputStructures::PrepareCollimatorInformation().
|
inline |
Return a reference to the last element.
Definition at line 214 of file BDSBeamline.hh.
References beamline.
Referenced by AddComponent(), AddSingleComponent(), ApplyTransform3D(), BDSDetectorConstruction::BuildBeamlines(), BDS::BuildEndPieceBeamline(), BDSTunnelBuilder::BuildTunnelSections(), BDS::CalculateTeleporterDelta(), and GetLastItem().
|
inline |
Iterator mechanics.
Definition at line 167 of file BDSBeamline.hh.
References beamline.
Referenced by BDSDetectorConstruction::BuildBeamline(), BDSCurvilinearBuilder::BuildCurvilinearBeamLine1To1(), BDSCurvilinearBuilder::BuildCurvilinearBridgeBeamLine(), BDSTunnelBuilder::BuildTunnelSections(), BDSOutputROOTEventModel::Fill(), FindFromS(), and GetElement().
|
inline |
G4bool BDSBeamline::ElementAnglesSumToCircle | ( | ) | const |
Check if the sum of the angle of all elements is two pi.
Definition at line 922 of file BDSBeamline.cc.
References totalAngle.
Referenced by BDSDetectorConstruction::BuildBeamlines().
|
inline |
Iterator mechanics.
Definition at line 175 of file BDSBeamline.hh.
References beamline.
Referenced by AddSingleComponent(), ApplyTransform3D(), BDSDetectorConstruction::BuildBeamline(), BDSCurvilinearBuilder::BuildCurvilinearBeamLine1To1(), BDS::BuildEndPieceBeamline(), BDSOutput::CalculateHistogramParameters(), and BDS::CalculateTeleporterDelta().
|
inline |
Iterator mechanics.
Definition at line 168 of file BDSBeamline.hh.
References beamline.
Referenced by BDSDetectorConstruction::BuildBeamline(), BDSCurvilinearBuilder::BuildCurvilinearBeamLine1To1(), BDSCurvilinearBuilder::BuildCurvilinearBridgeBeamLine(), BDSTunnelBuilder::BuildTunnelSections(), BDSOutputROOTEventModel::Fill(), GetElement(), and BDSDetectorConstruction::SideToLocalOffset().
|
inline |
BDSBeamline::const_iterator BDSBeamline::FindFromS | ( | G4double | S | ) | const |
Returns an iterator to the beamline element at s.
Definition at line 681 of file BDSBeamline.cc.
Referenced by BDSDetectorConstruction::SideToLocalOffset().
|
inline |
Return a reference to the first element.
Definition at line 212 of file BDSBeamline.hh.
References beamline.
Referenced by BDS::CalculateTeleporterDelta(), and GetFirstItem().
std::vector< G4double > BDSBeamline::GetEdgeSPositions | ( | ) | const |
Get the global s position of each element all in one - used for histograms. For convenience, s positions are converted to metres in this function.
Definition at line 909 of file BDSBeamline.cc.
References beamline, and BDSBeamlineElement::GetSPositionEnd().
Referenced by BDSOutput::CreateHistograms().
const BDSBeamlineElement * BDSBeamline::GetElement | ( | G4String | acceleratorComponentName, |
G4int | i = 0 |
||
) | const |
Get the ith placement of an element in the beam line. Returns null pointer if not found.
Definition at line 739 of file BDSBeamline.cc.
References begin(), components, end(), BDSBeamlineElement::GetPlacementName(), and BDS::StartsWith().
Referenced by BDSDetectorConstruction::CreatePlacementTransform(), and GetTransformForElement().
const BDSBeamlineElement * BDSBeamline::GetElementFromGlobalS | ( | G4double | S, |
G4int * | indexOfFoundElement = nullptr |
||
) | const |
Return the element in this beam line according to a given s coordinate.
Definition at line 670 of file BDSBeamline.cc.
References beamline, and sEnd.
Referenced by GetGlobalEuclideanTransform().
BDSExtentGlobal BDSBeamline::GetExtentGlobal | ( | ) | const |
Get the global extents for this beamline.
Definition at line 927 of file BDSBeamline.cc.
References maximumExtentNegative, and maximumExtentPositive.
Referenced by BDSDetectorConstruction::BuildWorld(), and BDSBeamlineSet::GetExtentGlobals().
|
inline |
Return a reference to the first element.
Definition at line 123 of file BDSBeamline.hh.
References front().
Referenced by BDS::BuildEndPieceBeamline(), and BDSCurvilinearBuilder::CreateBonusSectionStart().
G4Transform3D BDSBeamline::GetGlobalEuclideanTransform | ( | G4double | s, |
G4double | x = 0 , |
||
G4double | y = 0 , |
||
G4int * | indexOfFoundElement = nullptr |
||
) | const |
Get the local to global transform for curvilinear coordinates to global coordinates. 0,0 transverse position by default. Optionally returns the index of the found element in the beam line (by reference variable).
Definition at line 601 of file BDSBeamline.cc.
References at(), BDSAcceleratorComponent::GetAngle(), BDSAcceleratorComponent::GetArcLength(), BDSAcceleratorComponent::GetChordLength(), GetElementFromGlobalS(), BDSBeamlineElement::GetSPositionStart(), BDS::IsFinite(), and totalArcLength.
Referenced by BDSBunch::ApplyCurvilinearTransform(), and BDSDetectorConstruction::CreatePlacementTransform().
std::vector< G4int > BDSBeamline::GetIndicesOfCollimators | ( | ) | const |
Return indices in order of ecol, rcol, jcol and crystalcol elements.
Definition at line 957 of file BDSBeamline.cc.
References GetIndicesOfElementsOfType().
Referenced by BDSOutputStructures::PrepareCollimatorInformation(), and BDSSDCollimator::ProcessHitsOrdered().
std::vector< G4int > BDSBeamline::GetIndicesOfElementsOfType | ( | const G4String & | type | ) | const |
Return vector of indices for this beam line where element of type name 'type' is found.
Definition at line 934 of file BDSBeamline.cc.
References beamline, BDSBeamlineElement::GetIndex(), and BDSBeamlineElement::GetType().
Referenced by GetIndicesOfCollimators(), and BDSOutputStructures::PrepareCavityInformation().
std::vector< G4int > BDSBeamline::GetIndicesOfElementsOfType | ( | const std::set< G4String > & | types | ) | const |
Apply above function to get a set of types. These are in order as they appear in the beam line.
Definition at line 945 of file BDSBeamline.cc.
References beamline, BDSBeamlineElement::GetIndex(), and BDSBeamlineElement::GetType().
|
inline |
Return a reference to the last element.
Definition at line 126 of file BDSBeamline.hh.
References back().
Referenced by BDS::BuildEndPieceBeamline(), BDSOutput::CalculateHistogramParameters(), and BDSCurvilinearBuilder::CreateBonusSectionEnd().
G4ThreeVector BDSBeamline::GetMaximumExtentAbsolute | ( | ) | const |
Get the maximum extent absolute in each dimension.
Definition at line 593 of file BDSBeamline.cc.
References maximumExtentNegative, and maximumExtentPositive.
Referenced by BDSBeamlineSet::GetMaximumExtentAbsolute().
|
inline |
Get the maximum negative extent in all dimensions.
Definition at line 154 of file BDSBeamline.hh.
References maximumExtentNegative.
|
inline |
Get the maximum positive extent in all dimensions
Definition at line 151 of file BDSBeamline.hh.
References maximumExtentPositive.
const BDSBeamlineElement * BDSBeamline::GetNext | ( | const BDSBeamlineElement * | element | ) | const |
Definition at line 709 of file BDSBeamline.cc.
const BDSBeamlineElement * BDSBeamline::GetNext | ( | G4int | index | ) | const |
Definition at line 721 of file BDSBeamline.cc.
const BDSBeamlineElement * BDSBeamline::GetPrevious | ( | const BDSBeamlineElement * | element | ) | const |
Return a pointer to the previous element. First this beamline is searched for the vector. If there is no such element or no previous element because it's the beginning, then a nullptr is returned. The caller should test on this.
Definition at line 689 of file BDSBeamline.cc.
References beamline, and GetPrevious().
Referenced by AddSingleComponent(), BDS::BuildEndPieceBeamline(), and GetPrevious().
const BDSBeamlineElement * BDSBeamline::GetPrevious | ( | G4int | index | ) | const |
Definition at line 701 of file BDSBeamline.cc.
|
inline |
Definition at line 198 of file BDSBeamline.hh.
|
inline |
Definition at line 197 of file BDSBeamline.hh.
|
inline |
Get the total angle of the beamline.
Definition at line 142 of file BDSBeamline.hh.
References totalAngle.
|
inline |
Get the total ARC length for the beamline - ie the maximum s position.
Definition at line 139 of file BDSBeamline.hh.
References totalArcLength.
Referenced by BDSDetectorConstruction::BuildBeamline(), BDSRunAction::CheckTrajectoryOptions(), and BDSSurvey::Write().
|
inline |
Get the total length of the beamline - the sum of the chord length of each element.
Definition at line 136 of file BDSBeamline.hh.
References totalChordLength.
Referenced by BDSSurvey::Write().
G4Transform3D BDSBeamline::GetTransformForElement | ( | const G4String & | acceleratorComponentName, |
G4int | i = 0 |
||
) | const |
Get the transform to the centre of the ith placement of element by name. Uses GetElement(). Exits if no such element found.
Definition at line 777 of file BDSBeamline.cc.
References GetElement(), BDSBeamlineElement::GetPositionMiddle(), and BDSBeamlineElement::GetRotationMiddle().
G4bool BDSBeamline::IndexOK | ( | G4int | index | ) | const |
Whether the supplied index will lie within the beam line vector.
Definition at line 901 of file BDSBeamline.cc.
References beamline.
Referenced by ProvideEndPieceElementAfter().
|
inlinestatic |
Access the padding length between each element added to the beamline.
Definition at line 236 of file BDSBeamline.hh.
References paddingLength.
Referenced by BDSBeamlineElement::BDSBeamlineElement(), and BDS::CalculateTeleporterDelta().
void BDSBeamline::PrintAllComponents | ( | std::ostream & | out | ) | const |
Iterate over the beamline and print out the name, position, rotation and s position of each beamline element
Definition at line 102 of file BDSBeamline.cc.
void BDSBeamline::PrintMemoryConsumption | ( | ) | const |
Feedback about memory consumption for this beamline instance - container size, size of all BDSBeamlineElement() and size of all BDSAcceleratorComponent() stored.
Definition at line 108 of file BDSBeamline.cc.
References beamline.
BDSBeamlineElement * BDSBeamline::ProvideEndPieceElementAfter | ( | BDSSimpleComponent * | endPiece, |
G4int | index, | ||
G4bool | flip = false |
||
) | const |
Calculate the placements for an end piece w.r.t. a particlar beam line element The optional flip flag is used for when the 'before' piece is used again and must be rotated.
Definition at line 856 of file BDSBeamline.cc.
References beamline, BDSAcceleratorComponent::GetChordLength(), BDSBeamlineElement::GetPositionEnd(), BDSBeamlineElement::GetRotationMiddle(), BDSBeamlineElement::GetSPositionEnd(), BDSBeamlineElement::GetTiltOffset(), IndexOK(), and paddingLength.
Referenced by BDS::BuildEndPieceBeamline().
BDSBeamlineElement * BDSBeamline::ProvideEndPieceElementBefore | ( | BDSSimpleComponent * | endPiece, |
G4int | index | ||
) | const |
Definition at line 818 of file BDSBeamline.cc.
|
inline |
Iterator mechanics.
Definition at line 171 of file BDSBeamline.hh.
References beamline.
Referenced by BDSDetectorConstruction::BuildBeamline().
|
inline |
|
private |
Register the fully created element to a map of names vs element pointers. Used to look up transforms by name.
Definition at line 729 of file BDSBeamline.cc.
References components, BDSBeamlineElement::GetName(), and BDSBeamlineElement::GetPlacementName().
Referenced by AddBeamlineElement(), and AddSingleComponent().
|
inline |
Iterator mechanics.
Definition at line 172 of file BDSBeamline.hh.
References beamline.
Referenced by BDSDetectorConstruction::BuildBeamline().
|
inline |
|
inline |
Get the number of elements.
Definition at line 148 of file BDSBeamline.hh.
References beamline.
Referenced by BDSOutputROOTEventModel::Fill().
void BDSBeamline::UpdateExtents | ( | BDSBeamlineElement * | element | ) |
Once the beamline element has been constructed and all positions and rotations use these to update the world extent of this beam line.
Definition at line 794 of file BDSBeamline.cc.
References BDSExtentGlobal::AllBoundaryPointsGlobal(), BDSBeamlineElement::GetExtentGlobal(), maximumExtentNegative, and maximumExtentPositive.
Referenced by AddBeamlineElement(), and AddSingleComponent().
|
friend |
output stream
Definition at line 115 of file BDSBeamline.cc.
|
private |
Vector of beam line elements - the data.
Definition at line 67 of file BDSBeamline.hh.
Referenced by AddBeamlineElement(), AddSingleComponent(), at(), back(), begin(), empty(), end(), front(), GetEdgeSPositions(), GetElementFromGlobalS(), GetIndicesOfElementsOfType(), GetPrevious(), IndexOK(), PrintMemoryConsumption(), ProvideEndPieceElementAfter(), rbegin(), rend(), and size().
|
private |
Map of objects by placement name stored in this beam line.
Definition at line 289 of file BDSBeamline.hh.
Referenced by GetElement(), and RegisterElement().
|
private |
maximum extent in the negative coordinates in each dimension
Definition at line 270 of file BDSBeamline.hh.
Referenced by BDSBeamline(), GetExtentGlobal(), GetMaximumExtentAbsolute(), GetMaximumExtentNegative(), and UpdateExtents().
|
private |
maximum extent in the positive coordinates in each dimension
Definition at line 269 of file BDSBeamline.hh.
Referenced by BDSBeamline(), GetExtentGlobal(), GetMaximumExtentAbsolute(), GetMaximumExtentPositive(), and UpdateExtents().
|
staticprivate |
The gap added for padding between each component.
Definition at line 286 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), BDSBeamline(), PaddingLength(), and ProvideEndPieceElementAfter().
|
private |
Current reference position at the end of the previous element.
Definition at line 276 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), and ApplyTransform3D().
|
private |
Current reference rotation at the end of the previous element.
Definition at line 273 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), ApplyTransform3D(), and BDSBeamline().
|
private |
Current s coordinate at the end of the previous element.
Definition at line 279 of file BDSBeamline.hh.
Referenced by AddSingleComponent().
|
private |
Vector of s coordinates at the end of each element. This is intended so that an iterator pointing to the s position will be the correct index for the beamline element in the main BDSBeamlineVector element. This is filled in order so it's sorted by design.
Definition at line 295 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), FindFromS(), and GetElementFromGlobalS().
|
private |
Cache the initial S so we can tell if a requested S is too low.
Definition at line 259 of file BDSBeamline.hh.
|
private |
Definition at line 260 of file BDSBeamline.hh.
|
private |
Sum of all angles.
Definition at line 267 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), ElementAnglesSumToCircle(), and GetTotalAngle().
|
private |
Sum of all arc lengths.
Definition at line 265 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), GetGlobalEuclideanTransform(), and GetTotalArcLength().
|
private |
Sum of all chord lengths.
Definition at line 263 of file BDSBeamline.hh.
Referenced by AddSingleComponent(), and GetTotalChordLength().
|
private |
Flag to remember whether a BDSTransform3D which is not a real volume has been applied to make a discrete change in the beam line or not.
Definition at line 283 of file BDSBeamline.hh.
Referenced by AddSingleComponent().