BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamline.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#ifndef BDSBEAMLINE_H
20#define BDSBEAMLINE_H
21
22#include "globals.hh" // geant4 globals / types
23#include "G4ThreeVector.hh"
24#include "G4Transform3D.hh"
25
26#include "BDSSamplerType.hh"
27#include "BDSExtentGlobal.hh"
28
29#include <iterator>
30#include <ostream>
31#include <set>
32#include <vector>
33
36class BDSSamplerInfo;
38class BDSTiltOffset;
39class BDSTransform3D;
40namespace CLHEP {
41 class HepRotation;
42}
43typedef CLHEP::HepRotation G4RotationMatrix;
44
61{
62private:
64 typedef std::vector<BDSBeamlineElement*> BeamlineVector;
65
68
69public:
72 BDSBeamline(BDSBeamline&) = delete;
73
78 BDSBeamline(const G4ThreeVector& initialGlobalPosition = G4ThreeVector(0,0,0),
79 G4RotationMatrix* initialGlobalRotation = nullptr,
80 G4double initialS = 0.0);
81
83 explicit BDSBeamline(G4Transform3D initialTransform,
84 G4double initialS = 0.0);
85
87
93 BDSTiltOffset* tiltOffset = nullptr,
94 BDSSamplerInfo* samplerInfo = nullptr);
95
101 void ApplyTransform3D(BDSTransform3D* component);
102
110
113 void PrintAllComponents(std::ostream& out) const;
114
117 void UpdateExtents(BDSBeamlineElement* element);
118
120 inline const BDSBeamlineElement* at(int iElement) const { return beamline.at(iElement);}
121
123 inline const BDSBeamlineElement* GetFirstItem() const {return front();}
124
126 inline const BDSBeamlineElement* GetLastItem() const {return back();}
127
129 const BDSBeamlineElement* GetElement(G4String acceleratorComponentName, G4int i = 0) const;
130
133 G4Transform3D GetTransformForElement(const G4String& acceleratorComponentName, G4int i = 0) const;
134
136 inline G4double GetTotalChordLength() const {return totalChordLength;}
137
139 inline G4double GetTotalArcLength() const {return totalArcLength;}
140
142 inline G4double GetTotalAngle() const {return totalAngle;}
143
145 G4bool ElementAnglesSumToCircle() const;
146
148 BeamlineVector::size_type size() const {return beamline.size();}
149
151 G4ThreeVector GetMaximumExtentPositive() const {return maximumExtentPositive;}
152
154 G4ThreeVector GetMaximumExtentNegative() const {return maximumExtentNegative;}
155
157 G4ThreeVector GetMaximumExtentAbsolute() const;
158
161
163 typedef BeamlineVector::iterator iterator;
164 typedef BeamlineVector::const_iterator const_iterator;
165 typedef BeamlineVector::reverse_iterator reverse_iterator;
166 typedef BeamlineVector::const_reverse_iterator const_reverse_iterator;
167 iterator begin() {return beamline.begin();}
168 iterator end() {return beamline.end();}
169 const_iterator begin() const {return beamline.begin();}
170 const_iterator end() const {return beamline.end();}
171 reverse_iterator rbegin() {return beamline.rbegin();}
172 reverse_iterator rend() {return beamline.rend();}
173 const_reverse_iterator rbegin() const {return beamline.rbegin();}
174 const_reverse_iterator rend() const {return beamline.rend();}
175 G4bool empty() const {return beamline.empty();}
177
181 G4Transform3D GetGlobalEuclideanTransform(G4double s,
182 G4double x = 0,
183 G4double y = 0,
184 G4int* indexOfFoundElement = nullptr) const;
185
188 G4int* indexOfFoundElement = nullptr) const;
189
191 const_iterator FindFromS(G4double S) const;
192
195 std::vector<G4double> GetEdgeSPositions() const;
196
201 const BDSBeamlineElement* GetPrevious(const BDSBeamlineElement* element) const;
202 const BDSBeamlineElement* GetNext(const BDSBeamlineElement* element) const;
203
204 const BDSBeamlineElement* GetPrevious(G4int index) const;
205 const BDSBeamlineElement* GetNext(G4int index) const;
206
207 // Accessors in a similar style to std::vector
209 inline BDSBeamlineElement* front() const {return beamline.front();}
211 inline BDSBeamlineElement* back() const {return beamline.back();}
212
214 friend std::ostream& operator<< (std::ostream &out, BDSBeamline const &bl);
215
218 void PrintMemoryConsumption() const;
219
220 BDSBeamlineElement* ProvideEndPieceElementBefore(BDSSimpleComponent* endPiece,
221 G4int index) const;
226 G4int index,
227 G4bool flip = false) const;
228
230 G4bool IndexOK(G4int index) const;
231
233 static G4double PaddingLength() {return paddingLength;}
234
236 std::vector<G4int> GetIndicesOfElementsOfType(const G4String& type) const;
237
240 std::vector<G4int> GetIndicesOfElementsOfType(const std::set<G4String>& types) const;
241
243 std::vector<G4int> GetIndicesOfCollimators() const;
244
245private:
249 BDSTiltOffset* tiltOffset = nullptr,
250 BDSSamplerInfo* samplerInfo = nullptr);
251
254 void RegisterElement(BDSBeamlineElement* element);
255
261 G4double totalAngle;
262
263 G4ThreeVector maximumExtentPositive;
264 G4ThreeVector maximumExtentNegative;
265
268
271
274
278
280 static G4double paddingLength;
281
283 std::map<G4String, BDSBeamlineElement*> components;
284
289 std::vector<G4double> sEnd;
290};
291
292#endif
Abstract class that represents a component of an accelerator.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
std::vector< BDSBeamlineElement * > BeamlineVector
Typedefs up first so we can declare public iterators.
Definition: BDSBeamline.hh:64
G4double totalChordLength
Sum of all chord lengths.
Definition: BDSBeamline.hh:257
std::vector< G4double > GetEdgeSPositions() const
Definition: BDSBeamline.cc:906
G4ThreeVector previousReferencePositionEnd
Current reference position at the end of the previous element.
Definition: BDSBeamline.hh:270
const BDSBeamlineElement * GetPrevious(const BDSBeamlineElement *element) const
Definition: BDSBeamline.cc:686
void PrintAllComponents(std::ostream &out) const
Definition: BDSBeamline.cc:100
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.
Definition: BDSBeamline.cc:736
friend std::ostream & operator<<(std::ostream &out, BDSBeamline const &bl)
output stream
Definition: BDSBeamline.cc:113
BDSBeamlineElement * front() const
Return a reference to the first element.
Definition: BDSBeamline.hh:209
BDSBeamlineElement * back() const
Return a reference to the last element.
Definition: BDSBeamline.hh:211
G4double GetTotalAngle() const
Get the total angle of the beamline.
Definition: BDSBeamline.hh:142
G4ThreeVector maximumExtentPositive
maximum extent in the positive coordinates in each dimension
Definition: BDSBeamline.hh:263
const BDSBeamlineElement * GetLastItem() const
Return a reference to the last element.
Definition: BDSBeamline.hh:126
BDSBeamline & operator=(const BDSBeamline &)=delete
assignment and copy constructor not implemented nor used.
G4double GetTotalChordLength() const
Get the total length of the beamline - the sum of the chord length of each element.
Definition: BDSBeamline.hh:136
G4double previousSPositionEnd
Current s coordinate at the end of the previous element.
Definition: BDSBeamline.hh:273
G4bool ElementAnglesSumToCircle() const
Check if the sum of the angle of all elements is two pi.
Definition: BDSBeamline.cc:919
G4double totalAngle
Sum of all angles.
Definition: BDSBeamline.hh:261
BDSBeamlineElement * ProvideEndPieceElementAfter(BDSSimpleComponent *endPiece, G4int index, G4bool flip=false) const
Definition: BDSBeamline.cc:853
G4RotationMatrix * previousReferenceRotationEnd
Current reference rotation at the end of the previous element.
Definition: BDSBeamline.hh:267
G4bool transformHasJustBeenApplied
Definition: BDSBeamline.hh:277
void RegisterElement(BDSBeamlineElement *element)
Definition: BDSBeamline.cc:726
const_iterator begin() const
Iterator mechanics.
Definition: BDSBeamline.hh:169
const_reverse_iterator rend() const
Iterator mechanics.
Definition: BDSBeamline.hh:174
void UpdateExtents(BDSBeamlineElement *element)
Definition: BDSBeamline.cc:791
G4ThreeVector GetMaximumExtentAbsolute() const
Get the maximum extent absolute in each dimension.
Definition: BDSBeamline.cc:590
void PrintMemoryConsumption() const
Definition: BDSBeamline.cc:106
void ApplyTransform3D(BDSTransform3D *component)
Definition: BDSBeamline.cc:539
static G4double PaddingLength()
Access the padding length between each element added to the beamline.
Definition: BDSBeamline.hh:233
static G4double paddingLength
The gap added for padding between each component.
Definition: BDSBeamline.hh:280
G4ThreeVector GetMaximumExtentNegative() const
Get the maximum negative extent in all dimensions.
Definition: BDSBeamline.hh:154
const_reverse_iterator rbegin() const
Iterator mechanics.
Definition: BDSBeamline.hh:173
const BDSBeamlineElement * at(int iElement) const
Return a reference to the element at i.
Definition: BDSBeamline.hh:120
reverse_iterator rend()
Iterator mechanics.
Definition: BDSBeamline.hh:172
const BDSBeamlineElement * GetFirstItem() const
Return a reference to the first element.
Definition: BDSBeamline.hh:123
G4bool IndexOK(G4int index) const
Whether the supplied index will lie within the beam line vector.
Definition: BDSBeamline.cc:898
BeamlineVector::iterator iterator
Iterator mechanics.
Definition: BDSBeamline.hh:163
BeamlineVector::size_type size() const
Get the number of elements.
Definition: BDSBeamline.hh:148
void AddBeamlineElement(BDSBeamlineElement *element)
Definition: BDSBeamline.cc:571
G4Transform3D GetTransformForElement(const G4String &acceleratorComponentName, G4int i=0) const
Definition: BDSBeamline.cc:774
const_iterator FindFromS(G4double S) const
Returns an iterator to the beamline element at s.
Definition: BDSBeamline.cc:678
G4double totalArcLength
Sum of all arc lengths.
Definition: BDSBeamline.hh:259
BeamlineVector::const_reverse_iterator const_reverse_iterator
Iterator mechanics.
Definition: BDSBeamline.hh:166
BeamlineVector beamline
Vector of beam line elements - the data.
Definition: BDSBeamline.hh:67
const BDSBeamlineElement * GetElementFromGlobalS(G4double S, G4int *indexOfFoundElement=nullptr) const
Return the element in this beam line according to a given s coordinate.
Definition: BDSBeamline.cc:667
G4bool empty() const
Iterator mechanics.
Definition: BDSBeamline.hh:175
void AddSingleComponent(BDSAcceleratorComponent *component, BDSTiltOffset *tiltOffset=nullptr, BDSSamplerInfo *samplerInfo=nullptr)
Definition: BDSBeamline.cc:207
G4ThreeVector GetMaximumExtentPositive() const
Get the maximum positive extent in all dimensions
Definition: BDSBeamline.hh:151
G4Transform3D GetGlobalEuclideanTransform(G4double s, G4double x=0, G4double y=0, G4int *indexOfFoundElement=nullptr) const
Definition: BDSBeamline.cc:598
BDSExtentGlobal GetExtentGlobal() const
Get the global extents for this beamline.
Definition: BDSBeamline.cc:924
reverse_iterator rbegin()
Iterator mechanics.
Definition: BDSBeamline.hh:171
iterator begin()
Iterator mechanics.
Definition: BDSBeamline.hh:167
const_iterator end() const
Iterator mechanics.
Definition: BDSBeamline.hh:170
std::map< G4String, BDSBeamlineElement * > components
Map of objects by placement name stored in this beam line.
Definition: BDSBeamline.hh:283
void AddComponent(BDSAcceleratorComponent *component, BDSTiltOffset *tiltOffset=nullptr, BDSSamplerInfo *samplerInfo=nullptr)
Definition: BDSBeamline.cc:122
iterator end()
Iterator mechanics.
Definition: BDSBeamline.hh:168
std::vector< G4int > GetIndicesOfElementsOfType(const G4String &type) const
Return vector of indices for this beam line where element of type name 'type' is found.
Definition: BDSBeamline.cc:931
G4double GetTotalArcLength() const
Get the total ARC length for the beamline - ie the maximum s position.
Definition: BDSBeamline.hh:139
std::vector< G4double > sEnd
Definition: BDSBeamline.hh:289
BeamlineVector::const_iterator const_iterator
Iterator mechanics.
Definition: BDSBeamline.hh:164
std::vector< G4int > GetIndicesOfCollimators() const
Return indices in order of ecol, rcol, jcol and crystalcol elements.
Definition: BDSBeamline.cc:954
BeamlineVector::reverse_iterator reverse_iterator
Iterator mechanics.
Definition: BDSBeamline.hh:165
G4ThreeVector maximumExtentNegative
maximum extent in the negative coordinates in each dimension
Definition: BDSBeamline.hh:264
Holder for +- extents in 3 dimensions with a rotation and translation.
All info required to build a sampler but not place it.
A BDSAcceleratorComponent wrapper for BDSGeometryComponent.
A holder for any placement offsets and rotations for a BDSAcceleratorComponent.
Transform in beam line that takes up no space.