BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamlineElement.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 BDSBEAMLINEELEMENT_H
20#define BDSBEAMLINEELEMENT_H
21
22#include "BDSAcceleratorComponent.hh"
23#include "BDSExtent.hh"
24#include "BDSExtentGlobal.hh"
25#include "BDSSamplerInfo.hh"
26#include "BDSSamplerType.hh"
27#include "BDSTiltOffset.hh"
28
29#include "globals.hh" // geant4 globals / types
30#include "G4ThreeVector.hh"
31
32#include <ostream>
33#include <set>
34
35namespace HepGeom {
36 class Transform3D;
37}
38typedef HepGeom::Transform3D G4Transform3D;
39namespace CLHEP {
40 class HepRotation;
41}
42typedef CLHEP::HepRotation G4RotationMatrix;
43class G4VPhysicalVolume;
44
64{
65public:
66 BDSBeamlineElement() = delete;
68 BDSBeamlineElement& operator=(const BDSBeamlineElement&) = delete;
70 const G4ThreeVector& positionStartIn,
71 const G4ThreeVector& positionMiddleIn,
72 const G4ThreeVector& positionEndIn,
73 G4RotationMatrix* rotationStartIn,
74 G4RotationMatrix* rotationMiddleIn,
75 G4RotationMatrix* rotationEndIn,
76 const G4ThreeVector& referencePositionStartIn,
77 const G4ThreeVector& referencePositionMiddleIn,
78 const G4ThreeVector& referencePositionEndIn,
79 G4RotationMatrix* referenceRotationStartIn,
80 G4RotationMatrix* referenceRotationMiddleIn,
81 G4RotationMatrix* referenceRotationEndIn,
82 G4double sPositionStartIn,
83 G4double sPositionMiddleIn,
84 G4double sPositionEndIn,
85 BDSTiltOffset* tiltOffsetIn = nullptr,
86 BDSSamplerInfo* samplerInfoIn = nullptr,
87 G4int indexIn = -1);
88
90
93 std::set<G4VPhysicalVolume*> PlaceElement(const G4String& pvName,
94 G4VPhysicalVolume* containerPV,
95 G4bool useCLPlacementTransform,
96 G4int copyNumber,
97 G4bool checkOverlaps) const;
98
100 static std::set<G4VPhysicalVolume*> GetPVsFromAssembly(G4AssemblyVolume* av);
101
104 inline G4String GetName() const {return component->GetName();}
105 inline G4String GetType() const {return component->GetType();}
106 inline G4double GetArcLength() const {return component->GetArcLength();}
107 inline G4double GetChordLength() const {return component->GetChordLength();}
108 inline G4double GetAngle() const {return component->GetAngle();}
110 inline BDSExtent GetExtent() const {return component->GetExtent();}
111 inline G4String GetPlacementName() const {return placementName;}
112 inline G4int GetCopyNo() const {return copyNumber;}
113 inline G4ThreeVector GetPositionStart() const {return positionStart;}
114 inline G4ThreeVector GetPositionMiddle() const {return positionMiddle;}
115 inline G4ThreeVector GetPositionEnd() const {return positionEnd;}
116 inline G4RotationMatrix* GetRotationStart() const {return rotationStart;}
117 inline G4RotationMatrix* GetRotationMiddle() const {return rotationMiddle;}
118 inline G4RotationMatrix* GetRotationEnd() const {return rotationEnd;}
119 inline G4ThreeVector GetReferencePositionStart() const {return referencePositionStart;}
120 inline G4ThreeVector GetReferencePositionMiddle() const {return referencePositionMiddle;}
121 inline G4ThreeVector GetReferencePositionEnd() const {return referencePositionEnd;}
122 inline G4RotationMatrix* GetReferenceRotationStart() const {return referenceRotationStart;}
123 inline G4RotationMatrix* GetReferenceRotationMiddle() const {return referenceRotationMiddle;}
124 inline G4RotationMatrix* GetReferenceRotationEnd() const {return referenceRotationEnd;}
125 inline G4double GetSPositionStart() const {return sPositionStart;}
126 inline G4double GetSPositionMiddle() const {return sPositionMiddle;}
127 inline G4double GetSPositionEnd() const {return sPositionEnd;}
128 inline BDSTiltOffset* GetTiltOffset() const {return tiltOffset;}
129 inline G4Transform3D* GetPlacementTransform() const {return placementTransform;}
130 inline G4Transform3D* GetPlacementTransformCL() const {return placementTransformCL;}
131 inline BDSSamplerInfo* GetSamplerInfo() const {return samplerInfo;}
132 inline BDSSamplerType GetSamplerType() const {return samplerInfo ? samplerInfo->samplerType : BDSSamplerType::none;}
133 inline G4Transform3D* GetSamplerPlacementTransform() const {return samplerPlacementTransform;}
134 inline G4int GetIndex() const {return index;}
135 inline G4String GetMaterial() const {return component->Material();}
137
140
142 G4ThreeVector InputFaceNormal() const;
143 G4ThreeVector OutputFaceNormal() const;
145
147 G4bool AngledInputFace() const {return component->AngledInputFace();}
148 G4bool AngledOutputFace() const {return component->AngledOutputFace();}
150
152 inline G4double GetTilt() const {return tiltOffset ? tiltOffset->GetTilt() : 0.0;}
153
155 inline void SetReferencePositionEnd(G4ThreeVector newReferencePositionEnd);
156 inline void SetReferenceRotationEnd(G4RotationMatrix* newReferenceRotatonEnd);
158
161 void UpdateSamplerPlacementTransform(const G4Transform3D& tranfsormIn);
162
164 friend std::ostream& operator<< (std::ostream& out, BDSBeamlineElement const &element);
165
167 G4bool Overlaps(const BDSBeamlineElement* otherElement) const;
168
169private:
172
178
181
183 G4ThreeVector positionStart;
184 G4ThreeVector positionMiddle;
185 G4ThreeVector positionEnd;
187
189 G4RotationMatrix* rotationStart;
190 G4RotationMatrix* rotationMiddle;
191 G4RotationMatrix* rotationEnd;
193
198 G4ThreeVector referencePositionEnd;
200
203 G4RotationMatrix* referenceRotationStart;
204 G4RotationMatrix* referenceRotationMiddle;
205 G4RotationMatrix* referenceRotationEnd;
207
211 G4double sPositionEnd;
213
216
219 G4Transform3D* placementTransform;
220
225 G4Transform3D* placementTransformCL;
226
227 BDSSamplerInfo* samplerInfo;
228
236
238 G4int index;
239};
240
241#endif
Abstract class that represents a component of an accelerator.
G4bool AngledInputFace() const
Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.
virtual G4String Material() const
Return the name of a material associated with the component - ie the primary material.
virtual G4String GetName() const
The name of the component without modification.
virtual G4double GetChordLength() const
virtual G4double GetArcLength() const
BDSBeamPipeInfo * GetBeamPipeInfo() const
G4bool AngledOutputFace() const
Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.
G4String GetType() const
Get a string describing the type of the component.
Holder class for all information required to describe a beam pipe model.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4double GetChordLength() const
Accessor.
G4ThreeVector positionEnd
Global coordinates for the start, middle and end of this beamline element.
G4ThreeVector OutputFaceNormal() const
Return face normal accounting for placement tilt of this component.
G4double GetAngle() const
Accessor.
G4ThreeVector GetPositionMiddle() const
Accessor.
G4ThreeVector referencePositionStart
G4RotationMatrix * referenceRotationMiddle
BDSAcceleratorComponent * component
The accelerator component.
void SetReferencePositionEnd(G4ThreeVector newReferencePositionEnd)
Reassign the end variable as required when applying a transform.
G4ThreeVector GetReferencePositionEnd() const
Accessor.
G4int index
Index of this item in the beamline - saves keeping track of iterators and conversion.
G4ThreeVector GetPositionEnd() const
Accessor.
G4ThreeVector positionStart
Global coordinates for the start, middle and end of this beamline element.
G4double sPositionEnd
S Positions for the start, middle and end of this beamline element.
friend std::ostream & operator<<(std::ostream &out, BDSBeamlineElement const &element)
Output stream.
BDSTiltOffset * tiltOffset
The tilt and offset this element was constructed with. Default is nullptr.
BDSSamplerType GetSamplerType() const
Accessor.
G4RotationMatrix * GetReferenceRotationStart() const
Accessor.
G4Transform3D * placementTransform
G4RotationMatrix * GetRotationEnd() const
Accessor.
std::set< G4VPhysicalVolume * > PlaceElement(const G4String &pvName, G4VPhysicalVolume *containerPV, G4bool useCLPlacementTransform, G4int copyNumber, G4bool checkOverlaps) const
G4ThreeVector GetPositionStart() const
Accessor.
BDSExtent GetExtent() const
Accessor.
G4String GetPlacementName() const
Accessor.
BDSTiltOffset * GetTiltOffset() const
Accessor.
void UpdateSamplerPlacementTransform(const G4Transform3D &tranfsormIn)
G4ThreeVector referencePositionEnd
G4ThreeVector InputFaceNormal() const
Return face normal accounting for placement tilt of this component.
G4bool AngledOutputFace() const
Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.
G4Transform3D * samplerPlacementTransform
G4Transform3D * GetPlacementTransform() const
Accessor.
G4RotationMatrix * GetReferenceRotationEnd() const
Accessor.
G4String GetMaterial() const
Accessor.
G4String GetName() const
Accessor.
G4Transform3D * placementTransformCL
G4double sPositionStart
S Positions for the start, middle and end of this beamline element.
BDSSamplerInfo * GetSamplerInfo() const
Accessor.
G4int GetIndex() const
Accessor.
G4ThreeVector GetReferencePositionMiddle() const
Accessor.
G4double GetSPositionEnd() const
Accessor.
G4RotationMatrix * GetRotationStart() const
Accessor.
G4ThreeVector GetReferencePositionStart() const
Accessor.
BDSAcceleratorComponent * GetAcceleratorComponent() const
Accessor.
G4Transform3D * GetPlacementTransformCL() const
Accessor.
G4RotationMatrix * GetReferenceRotationMiddle() const
Accessor.
G4bool Overlaps(const BDSBeamlineElement *otherElement) const
Whether this beam line element will overlaps in 3D Cartesian coordinates with another.
G4ThreeVector positionMiddle
Global coordinates for the start, middle and end of this beamline element.
G4int copyNumber
identification number of AcceleratorComponent (0 for first volume of given type)
G4RotationMatrix * referenceRotationStart
G4Transform3D * GetSamplerPlacementTransform() const
Accessor.
G4RotationMatrix * rotationStart
Global rotation matrices for the start, middle and end of this beamline element.
G4RotationMatrix * rotationMiddle
Global rotation matrices for the start, middle and end of this beamline element.
G4double GetSPositionMiddle() const
Accessor.
G4double sPositionMiddle
S Positions for the start, middle and end of this beamline element.
void SetReferenceRotationEnd(G4RotationMatrix *newReferenceRotatonEnd)
Reassign the end variable as required when applying a transform.
G4int GetCopyNo() const
Accessor.
G4RotationMatrix * rotationEnd
Global rotation matrices for the start, middle and end of this beamline element.
G4RotationMatrix * referenceRotationEnd
G4double GetArcLength() const
Accessor.
G4RotationMatrix * GetRotationMiddle() const
Accessor.
G4double GetTilt() const
Convenience accessor.
G4String GetType() const
Accessor.
BDSBeamPipeInfo * GetBeamPipeInfo() const
Accessor.
static std::set< G4VPhysicalVolume * > GetPVsFromAssembly(G4AssemblyVolume *av)
Utility method to account for the interface in G4AssemblyVolume.
G4bool AngledInputFace() const
Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.
G4ThreeVector referencePositionMiddle
BDSExtentGlobal GetExtentGlobal() const
Create a global extent object from the extent of the component.
G4double GetSPositionStart() const
Accessor.
Holder for +- extents in 3 dimensions with a rotation and translation.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
BDSExtent GetExtent() const
Accessor - see member for more info.
All info required to build a sampler but not place it.
A holder for any placement offsets and rotations for a BDSAcceleratorComponent.
G4double GetTilt() const
Accessor.