BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipeFactoryClicPCL.cc
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#include "BDSBeamPipeFactoryBase.hh"
20#include "BDSBeamPipeFactoryClicPCL.hh"
21#include "BDSBeamPipe.hh"
22#include "BDSDebug.hh"
23#include "BDSExtent.hh"
24
25#include "globals.hh" // geant4 globals / types
26#include "G4EllipticalTube.hh"
27#include "G4CutTubs.hh"
28#include "G4IntersectionSolid.hh"
29#include "G4SubtractionSolid.hh"
30#include "G4ThreeVector.hh"
31#include "G4VSolid.hh"
32
33#include <cmath> // sin, cos, fabs
34
35BDSBeamPipeFactoryClicPCL::BDSBeamPipeFactoryClicPCL():
36 extentYLow(0),
37 extentYHigh(0)
38{;}
39
40void BDSBeamPipeFactoryClicPCL::GenerateClicPCL(std::vector<G4TwoVector>& vec,
41 G4double aper1,
42 G4double aper2,
43 G4double aper3,
44 G4double aper4,
45 G4int pointsPerTwoPi,
46 G4double margin)
47{
48 G4double a1 = aper1 + margin;
49 G4double a2 = aper2 + margin;
50 G4double a3 = aper3 + margin;
51 G4double a4 = aper4;
52
53 G4int nPoints = 0.5*pointsPerTwoPi;
54 AppendAngleEllipse(vec, -CLHEP::halfpi, CLHEP::halfpi, a1, a2, nPoints, 0, a4);
55 AppendPoint(vec, a1, a4);
56 AppendAngleEllipse(vec, CLHEP::halfpi, CLHEP::halfpi + CLHEP::pi, a1, a3, nPoints);
57 AppendPoint(vec, -a1, 0);
58}
59
61 G4double aper2,
62 G4double aper3,
63 G4double aper4,
64 G4double beamPipeThickness,
65 G4int pointsPerTwoPi)
66{
67 G4double bpInMargin = lengthSafetyLarge;
68 G4double bpOutMargin = bpInMargin + beamPipeThickness;
69 G4double cont1Margin = bpOutMargin + lengthSafetyLarge;
70 G4double cont2Margin = cont1Margin + lengthSafetyLarge;
71 GenerateClicPCL(vacuumEdge, aper1, aper2, aper3, aper4, pointsPerTwoPi);
72 GenerateClicPCL(beamPipeInnerEdge, aper1, aper2, aper3, aper4, pointsPerTwoPi, bpInMargin);
73 GenerateClicPCL(beamPipeOuterEdge, aper1, aper2, aper3, aper4, pointsPerTwoPi, bpOutMargin);
74 GenerateClicPCL(containerEdge, aper1, aper2, aper3, aper4, pointsPerTwoPi, cont1Margin);
75 GenerateClicPCL(containerSubtractionEdge, aper1, aper2, aper3, aper4, pointsPerTwoPi, cont2Margin);
76
77 extentX = aper1 + cont1Margin;
78 extentYLow = -(std::abs(aper3) + cont1Margin);
79 extentYHigh = aper2 + aper4 + cont1Margin;
80}
81
83 G4double aper2,
84 G4double aper3,
85 G4double aper4,
86 G4double beamPipeThickness)
87{
88 G4double bottom = std::max(std::abs(aper3), aper1);
89 G4double top = std::max(aper1, aper2);
90 G4double hyp = std::hypot(aper1, aper4);
91
92 G4double result = std::max(std::max(bottom, top), hyp) + beamPipeThickness;
93 result *= 1.2; // 20% margin
94 return result;
95}
96
98{
99 extentYLow = 0;
100 extentYHigh = 0;
101}
102
104{
107}
108
110 G4Material* vacuumMaterialIn,
111 G4Material* beamPipeMaterialIn,
112 G4double lengthIn)
113{
114 BDSBeamPipeFactoryBase::CommonConstruction(nameIn, vacuumMaterialIn,
115 beamPipeMaterialIn, lengthIn);
116
117 // record extents
120 -lengthIn*0.5,lengthIn*0.5);
121
122 // calculate radius if a tube were to be place around it
123 G4double containerRadius = ext.MaximumAbsTransverse();
124
125 BDSBeamPipe* aPipe = BuildBeamPipeAndRegisterVolumes(ext,containerRadius);
126
127 return aPipe;
128}
void CommonConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double length)
finalise beampipe construction
BDSBeamPipe * BuildBeamPipeAndRegisterVolumes(BDSExtent extent, G4double containerRadius, G4bool containerIsCircular=false)
build beampipe and register logical volumes
virtual BDSBeamPipe * CommonFinalConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double lengthIn)
Overloads BDSBeamPipeFactoryPoints to make asymmetric extents, otherwise the same.
virtual void GeneratePoints(G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness, G4int pointsPerTwoPi=40)
G4double extentYLow
Cache of extent to pass around.
virtual void CleanUp()
Clear member vectors and run base class clean up to clear pointers between runs.
virtual G4double CalculateIntersectionRadius(G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness)
Calculate the radius of the solid used for intersection for angled faces.
G4double extentYHigh
Cache of extent to pass around.
void CleanUpClicPCL()
Clear member vectors - used for both initialisation and virtual CleanUp.
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.
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.
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.
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
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double MaximumAbsTransverse() const
Return the maximum absolute value considering only x,y.
Definition: BDSExtent.cc:171
G4double lengthSafetyLarge
Cache of global constants variable.