BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipeFactoryPoints.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 BDSBEAMPIPEFACTORYPOINTS_H
20#define BDSBEAMPIPEFACTORYPOINTS_H
21
22#include "BDSBeamPipeFactoryBase.hh"
23#include "BDSBeamPipe.hh"
24
25#include "globals.hh"
26#include "G4TwoVector.hh"
27
28#include <vector>
29
30class G4Material;
31
48{
49public:
52
54 virtual BDSBeamPipe* CreateBeamPipe(const G4String& nameIn,
55 G4double lengthIn,
56 G4double aper1In,
57 G4double aper2In,
58 G4double aper3In,
59 G4double aper4In,
60 G4Material* vacuumMaterialIn,
61 G4double beamPipeThicknessIn,
62 G4Material* beamPipeMaterialIn,
63 const G4String& pointsFileIn = "",
64 const G4String& pointsUnitIn = "");
65
67 virtual BDSBeamPipe* CreateBeamPipe(const G4String& nameIn,
68 G4double lengthIn,
69 const G4ThreeVector& inputFaceNormalIn,
70 const G4ThreeVector& outputFaceNormalIn,
71 G4double aper1In,
72 G4double aper2In,
73 G4double aper3In,
74 G4double aper4In,
75 G4Material* vacuumMaterialIn,
76 G4double beamPipeThicknessIn,
77 G4Material* beamPipeMaterialIn,
78 const G4String& pointsFileIn = "",
79 const G4String& pointsUnitIn = "");
80
81protected:
85 virtual void GeneratePoints(G4double aper1,
86 G4double aper2,
87 G4double aper3,
88 G4double aper4,
89 G4double beamPipeThickness,
90 G4int pointsPerTwoPi = 40) = 0;
91
97 virtual G4double CalculateIntersectionRadius(G4double aper1,
98 G4double aper2,
99 G4double aper3,
100 G4double aper4,
101 G4double beamPipeThickness) = 0;
102
104 virtual void CleanUp();
105
106 void AppendPoint(std::vector<G4TwoVector>& vec,
107 G4double x,
108 G4double y);
109
112 void AppendAngle(std::vector<G4TwoVector>& vec,
113 G4double startAngle,
114 G4double finishAngle,
115 G4double radius,
116 G4int nPoints = 10,
117 G4double xOffset = 0,
118 G4double yOffset = 0);
119
121 void AppendAngleEllipse(std::vector<G4TwoVector>& vec,
122 G4double startAngle,
123 G4double finishAngle,
124 G4double radiusA, // radius in horizontal
125 G4double radiusB, // radius in vertical
126 G4int nPoints = 10,
127 G4double xOffset = 0,
128 G4double yOffset = 0);
129
131 std::vector<G4TwoVector> vacuumEdge;
132
134 std::vector<G4TwoVector> beamPipeInnerEdge;
135
137 std::vector<G4TwoVector> beamPipeOuterEdge;
138
140 std::vector<G4TwoVector> containerEdge;
141
143 std::vector<G4TwoVector> containerSubtractionEdge;
144
146 G4double extentX;
147 G4double extentY;
149
151 G4String pointsFile;
152 G4String pointsUnit;
154
155private:
162
165
168
170 void CleanUpPoints();
171
175 void CreateSolids(const G4String& name,
176 G4double length,
177 G4bool buildLongForIntersection = false);
178
180 void CreateSolidsAngled(const G4String& name,
181 G4double length,
182 const G4ThreeVector& inputFace,
183 const G4ThreeVector& outputFace);
184
187 virtual BDSBeamPipe* CommonFinalConstruction(const G4String& nameIn,
188 G4Material* vacuumMaterialIn,
189 G4Material* beamPipeMaterialIn,
190 G4double lengthIn);
191};
192
193#endif
Abstract base class for beampipe factory classes.
Factory for beam pipes defined by a series of x,y points that are extruded.
G4VSolid * beamPipeInnerSolid
Solid for inner surface of beam pipe.
std::vector< G4TwoVector > containerSubtractionEdge
Vector of x,y coordinates for the container subtraction volume.
std::vector< G4TwoVector > beamPipeInnerEdge
Vector of x,y coordinates for beam pipe inner edge.
virtual void CleanUp()
Clear member vectors and run base class clean up to clear pointers between runs.
G4VSolid * beamPipeOuterSolid
Solid for outer surface of beam pipe.
virtual void GeneratePoints(G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness, G4int pointsPerTwoPi=40)=0
G4double extentY
Extents prepared by GeneratePoints function as only it knows the extents.
virtual G4double CalculateIntersectionRadius(G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness)=0
void CleanUpPoints()
Actual clean up here for this class only.
virtual BDSBeamPipe * CommonFinalConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double lengthIn)
G4double extentX
Extents prepared by GeneratePoints function as only it knows the extents.
std::vector< G4TwoVector > vacuumEdge
Vector of x,y coordinates for vacuum extruded solid edge.
virtual BDSBeamPipe * CreateBeamPipe(const G4String &nameIn, G4double lengthIn, G4double aper1In, G4double aper2In, G4double aper3In, G4double aper4In, G4Material *vacuumMaterialIn, G4double beamPipeThicknessIn, G4Material *beamPipeMaterialIn, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")
Required overloaded method from BDSBeamPipeFactoryBase to build straight piece of beam pipe.
void CreateSolidsAngled(const G4String &name, G4double length, const G4ThreeVector &inputFace, const G4ThreeVector &outputFace)
Create angled solids. Uses CreateSolids() method.
void AppendAngle(std::vector< G4TwoVector > &vec, G4double startAngle, G4double finishAngle, G4double radius, G4int nPoints=10, G4double xOffset=0, G4double yOffset=0)
void CreateSolids(const G4String &name, G4double length, G4bool buildLongForIntersection=false)
G4String pointsFile
Temporary copy of input variables.
std::vector< G4TwoVector > containerEdge
Vector of x,y coordinates for the container volume.
std::vector< G4TwoVector > beamPipeOuterEdge
Vector of x,y coordinates for beam pipe outer edge.
G4String pointsUnit
Temporary copy of input variables.
void AppendAngleEllipse(std::vector< G4TwoVector > &vec, G4double startAngle, G4double finishAngle, G4double radiusA, G4double radiusB, G4int nPoints=10, G4double xOffset=0, G4double yOffset=0)
Generate 2-vector points (and append them) about an ellipse.
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45