BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSAcceleratorComponent.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 BDSACCELERATORCOMPONENT_H
20#define BDSACCELERATORCOMPONENT_H
21
22#include "globals.hh" // geant4 globals / types
23
24#include "BDSGeometryComponent.hh"
25
26#include <list>
27#include <set>
28#include <string>
29#include <vector>
30
31class BDSBeamPipeInfo;
32class BDSFieldInfo;
34class G4LogicalVolume;
35class G4UserLimits;
36class G4VisAttributes;
37
78{
79public:
96 BDSAcceleratorComponent(const G4String& name,
97 G4double arcLength,
98 G4double angle,
99 const G4String& type,
100 BDSBeamPipeInfo* beamPipeInfo = nullptr,
101 const G4ThreeVector& inputFaceNormalIn = G4ThreeVector(0,0,-1),
102 const G4ThreeVector& outputFaceNormalIn = G4ThreeVector(0,0, 1),
103 BDSFieldInfo* fieldInfoIn = nullptr);
104
105 virtual ~BDSAcceleratorComponent();
106
111 virtual void Initialise();
112
113 // Communal constructions tasks
114
116 virtual void SetBiasVacuumList(const std::list<std::string>& biasVacuumListIn)
117 {biasVacuumList = biasVacuumListIn;}
118 virtual void SetBiasMaterialList(const std::list<std::string>& biasMaterialListIn)
119 {biasMaterialList = biasMaterialListIn;}
121
123 virtual void SetRegion(const G4String& regionIn) {region = regionIn;}
124
126 void SetField(BDSFieldInfo* fieldInfoIn);
127
130 virtual void SetMinimumKineticEnergy(G4double /*minimumKineticEnergyIn*/){;}
131
132 // Accessors
133
135 virtual inline G4String GetName() const {return name;}
136
139 virtual G4double GetArcLength() const {return arcLength;}
140 virtual G4double GetChordLength() const {return chordLength;}
142
145 inline G4double GetAngle() const {return angle;}
146
148 G4double Sagitta() const;
149
151 inline G4String GetType() const {return type;}
152
154 G4String GetRegion() const {return region;}
155
157 virtual G4bool HasAField() const {return fieldInfo;}
158
163
168
172 inline G4ThreeVector InputFaceNormal() const {return inputFaceNormal;}
173 inline G4ThreeVector OutputFaceNormal() const {return outputFaceNormal;}
175
177 G4bool AngledInputFace() const;
178 G4bool AngledOutputFace() const;
180
182 virtual G4String Material() const {return "none";}
183
185 virtual std::set<G4LogicalVolume*> GetAcceleratorVacuumLogicalVolumes() const {return acceleratorVacuumLV;}
186
188 virtual std::set<G4LogicalVolume*> GetAcceleratorMaterialLogicalVolumes() const;
189
192
194 inline G4int GetCopyNumber() const {return copyNumber;}
195
197 std::list<std::string> GetBiasVacuumList() const {return biasVacuumList;}
198 std::list<std::string> GetBiasMaterialList() const {return biasMaterialList;}
200
203 BDSSimpleComponent* EndPieceBefore() const {return endPieceBefore;}
204 BDSSimpleComponent* EndPieceAfter() const {return endPieceAfter;}
205
207 virtual void SetInputFaceNormal(const G4ThreeVector& input) {inputFaceNormal = input.unit();}
208 virtual void SetOutputFaceNormal(const G4ThreeVector& output) {outputFaceNormal = output.unit();}
210
213 static G4double lengthSafetyLarge;
214
215protected:
222 virtual void Build();
223
226 virtual void BuildContainerLogicalVolume() = 0;
227
233 inline void SetAcceleratorVacuumLogicalVolume(G4LogicalVolume* accVacLVIn)
234 {acceleratorVacuumLV.insert(accVacLVIn);}
235
236 inline void SetAcceleratorVacuumLogicalVolume(const std::set<G4LogicalVolume*>& accVacLVIn)
237 {acceleratorVacuumLV.insert(accVacLVIn.begin(), accVacLVIn.end());}
238
245 virtual void BuildUserLimits();
246
248 virtual void AttachUserLimits() const;
249
251 const G4String name;
252 const G4double arcLength;
253 const G4String type;
255
257 G4double chordLength;
258 G4double angle;
259 G4String region;
261
264
266 static G4double lengthSafety;
267 static G4Material* emptyMaterial;
268 static G4Material* worldMaterial;
269 static G4bool checkOverlaps;
270 static G4bool sensitiveOuter;
271 static G4bool sensitiveVacuum;
272 static G4VisAttributes* containerVisAttr;
274
276 std::set<G4LogicalVolume*> acceleratorVacuumLV;
277
278 BDSSimpleComponent* endPieceBefore;
279 BDSSimpleComponent* endPieceAfter;
280
281 G4UserLimits* userLimits;
282
284
285private:
290
295
296 //A vector containing the physical volumes in the accelerator component- to be used for geometric importance sampling etc.
297
302
305
307 std::list<std::string> biasVacuumList;
308 std::list<std::string> biasMaterialList;
310
311 G4ThreeVector inputFaceNormal;
312 G4ThreeVector outputFaceNormal;
313};
314
315#endif
Abstract class that represents a component of an accelerator.
BDSFieldInfo * fieldInfo
Recipe for field that could overlay this whole component.
G4UserLimits * userLimits
Cache of user limits.
G4bool AngledInputFace() const
Whether the face normal is angled at all w.r.t. the incoming / outgoing reference trajectory.
std::list< std::string > biasVacuumList
Copy of bias list from parser for this particular component.
G4double Sagitta() const
Calculate the sagitta - ie the distance between the chord and the arc at the centre.
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.
G4ThreeVector outputFaceNormal
Output face unit normal vector in outgoing reference coordinate frame.
G4String GetRegion() const
Get the region name for this component.
BDSAcceleratorComponent & operator=(const BDSAcceleratorComponent &)=delete
Assignment and copy constructor not implemented nor used.
G4ThreeVector inputFaceNormal
Input face unit normal vector in incoming reference coordinate frame.
const G4double arcLength
Const protected member variable that may not be changed by derived classes.
void SetAcceleratorVacuumLogicalVolume(G4LogicalVolume *accVacLVIn)
BDSSimpleComponent * EndPieceBefore() const
virtual G4bool HasAField() const
Whether this component has a field or not (ie is active). Implicit cast of pointer to bool.
BDSAcceleratorComponent(BDSAcceleratorComponent &)=delete
Assignment and copy constructor not implemented nor used.
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4Material * emptyMaterial
Useful variable often used in construction.
static G4bool sensitiveOuter
Useful variable often used in construction.
virtual void SetOutputFaceNormal(const G4ThreeVector &output)
Allow updating of face normals. Virtual so derived class may apply it to daughters.
static G4double lengthSafety
Useful variable often used in construction.
G4int copyNumber
Record of how many times this component has been copied.
virtual void SetMinimumKineticEnergy(G4double)
virtual G4double GetChordLength() const
static G4bool checkOverlaps
Useful variable often used in construction.
const G4String type
Const protected member variable that may not be changed by derived classes.
G4ThreeVector OutputFaceNormal() const
G4int GetCopyNumber() const
Get the number of times this component has been copied.
G4double angle
Protected member variable that can be modified by derived classes.
virtual G4double GetArcLength() const
virtual std::set< G4LogicalVolume * > GetAcceleratorVacuumLogicalVolumes() const
Access the 'vacuum' volume(s) in this component if any. Default is empty set.
G4ThreeVector InputFaceNormal() const
virtual std::set< G4LogicalVolume * > GetAcceleratorMaterialLogicalVolumes() const
Return a set of logical volumes excluding the ones in the 'vacuum' set.
static G4Material * worldMaterial
Useful variable often used in construction.
BDSAcceleratorComponent()=delete
BDSBeamPipeInfo * GetBeamPipeInfo() const
std::list< std::string > GetBiasVacuumList() 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 BuildContainerLogicalVolume()=0
G4String region
Protected member variable that can be modified by derived classes.
virtual void SetBiasVacuumList(const std::list< std::string > &biasVacuumListIn)
Copy the bias list to this element.
void SetField(BDSFieldInfo *fieldInfoIn)
Set the field definition for the whole component.
virtual void AttachUserLimits() const
Doesn't change member variables, but may change their contents.
void IncrementCopyNumber()
Increment (+1) the number of times this component has been copied.
virtual void SetBiasMaterialList(const std::list< std::string > &biasMaterialListIn)
Copy the bias list to this element.
G4double chordLength
Protected member variable that can be modified by derived classes.
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.
std::list< std::string > GetBiasMaterialList() const
Access the bias list copied from parser.
std::set< G4LogicalVolume * > acceleratorVacuumLV
A set of logical volumes we classify as 'vacuum' for biasing purposes.
std::list< std::string > biasMaterialList
Copy of bias list from parser for this particular component.
virtual void SetRegion(const G4String &regionIn)
Set the region name for this component.
static G4VisAttributes * containerVisAttr
Useful variable often used in construction.
static G4bool sensitiveVacuum
Useful variable often used in construction.
BDSBeamPipeInfo * beamPipeInfo
Optional beam pipe recipe that is written out to the survey if it exists.
virtual void SetFieldUsePlacementWorldTransform()
Holder class for all information required to describe a beam pipe model.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
A generic geometry component for a bdsim model.
A BDSAcceleratorComponent wrapper for BDSGeometryComponent.