BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipeFactoryRaceTrack.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 "BDSBeamPipeFactoryRaceTrack.hh"
20#include "BDSBeamPipeFactoryPoints.hh"
21#include "BDSDebug.hh"
22
23#include "globals.hh" // geant4 globals / types
24#include "G4TwoVector.hh"
25
26#include <cmath>
27
28BDSBeamPipeFactoryRaceTrack::BDSBeamPipeFactoryRaceTrack()
29{;}
30
31void BDSBeamPipeFactoryRaceTrack::GenerateRaceTrack(std::vector<G4TwoVector>& vec,
32 G4double x,
33 G4double y,
34 G4double r,
35 G4int pointsPerTwoPi)
36{
37 G4int pointsPerCurve = pointsPerTwoPi / 4;
38 AppendAngle(vec, 0, 0.5*CLHEP::pi, r, pointsPerCurve, x, y);
39 AppendPoint(vec, x+r, y);
40 AppendAngle(vec, 0.5*CLHEP::pi, CLHEP::pi, r, pointsPerCurve, x, -y);
41 AppendPoint(vec, x, -y-r);
42 AppendAngle(vec, CLHEP::pi, (3./2.)*CLHEP::pi, r, pointsPerCurve, -x, -y);
43 AppendPoint(vec, -x-r, -y);
44 AppendAngle(vec, (3./2)*CLHEP::pi, CLHEP::twopi, r, pointsPerCurve, -x, y);
45 AppendPoint(vec, -x, y+r);
46}
47
48
50 G4double aper2,
51 G4double aper3,
52 G4double /*aper4*/,
53 G4double beamPipeThickness,
54 G4int pointsPerTwoPi)
55{
56 GenerateRaceTrack(vacuumEdge, aper1, aper2, aper3, pointsPerTwoPi);
57 G4double bpInner1 = aper1 + lengthSafetyLarge;
58 G4double bpInner2 = aper2 + lengthSafetyLarge;
59 G4double bpInner3 = aper3 + lengthSafetyLarge;
60 GenerateRaceTrack(beamPipeInnerEdge, bpInner1, bpInner2, bpInner3, pointsPerTwoPi);
61 G4double bpOuter1 = bpInner1 + beamPipeThickness;
62 G4double bpOuter2 = bpInner2 + beamPipeThickness;
63 G4double bpOuter3 = bpInner3 + beamPipeThickness;
64 GenerateRaceTrack(beamPipeOuterEdge, bpOuter1, bpOuter2, bpOuter3, pointsPerTwoPi);
65 G4double cont1 = bpOuter1 + lengthSafetyLarge;
66 G4double cont2 = bpOuter2 + lengthSafetyLarge;
67 G4double cont3 = bpOuter3 + lengthSafetyLarge;
68 GenerateRaceTrack(containerEdge, cont1, cont2, cont3, pointsPerTwoPi);
69 G4double contSub1 = cont1 + lengthSafetyLarge;
70 G4double contSub2 = cont2 + lengthSafetyLarge;
71 G4double contSub3 = cont3 + lengthSafetyLarge;
72 GenerateRaceTrack(containerSubtractionEdge, contSub1, contSub2, contSub3, pointsPerTwoPi);
73
74 extentX = contSub1 + contSub3;
75 extentY = contSub2 + contSub3;
76}
77
79 G4double aper2,
80 G4double aper3,
81 G4double /*aper4*/,
82 G4double beamPipeThickness)
83{
84 G4double result = std::hypot(aper1,aper2) + aper3 + beamPipeThickness;
85 result *= 1.5; // 50% margin
86#ifdef BDSDEBUG
87 G4cout << __METHOD_NAME__ << "intersection radius: " << result << G4endl;
88#endif
89 return result;
90}
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.
G4double extentY
Extents prepared by GeneratePoints function as only it knows the extents.
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.
void AppendAngle(std::vector< G4TwoVector > &vec, G4double startAngle, G4double finishAngle, G4double radius, G4int nPoints=10, G4double xOffset=0, G4double yOffset=0)
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.
virtual void GeneratePoints(G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness, G4int pointsPerTwoPi=40)
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 lengthSafetyLarge
Cache of global constants variable.