BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamPipeInfo.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 BDSBEAMPIPEINFO_H
20#define BDSBEAMPIPEINFO_H
21
22#include "BDSBeamPipeType.hh"
23
24#include "globals.hh" // geant4 types / globals
25#include "G4ThreeVector.hh"
26
27class BDSExtent;
28class G4Material;
29
41{
42public:
44 BDSBeamPipeInfo() = delete;
46 BDSBeamPipeInfo(BDSBeamPipeType beamPipeTypeIn,
47 G4double aper1In,
48 G4double aper2In,
49 G4double aper3In,
50 G4double aper4In,
51 G4Material* vacuumMaterialIn,
52 G4double beamPipeThicknessIn,
53 G4Material* beamPipeMaterialIn,
54 const G4ThreeVector& inputFaceNormalIn = G4ThreeVector(0,0,-1),
55 const G4ThreeVector& outputFaceNormalIn = G4ThreeVector(0,0,1),
56 const G4String& pointsFileNameIn = "",
57 const G4String& pointsUnitIn = "mm");
58
61 BDSBeamPipeInfo(const G4String& beamPipeTypeIn,
62 G4double aper1In,
63 G4double aper2In,
64 G4double aper3In,
65 G4double aper4In,
66 const G4String& vacuumMaterialIn,
67 G4double beamPipeThicknessIn,
68 const G4String& beamPipeMaterialIn,
69 const G4ThreeVector& inputFaceNormalIn = G4ThreeVector(0,0,-1),
70 const G4ThreeVector& outputFaceNormalIn = G4ThreeVector(0,0,1));
71
74 BDSBeamPipeInfo(const BDSBeamPipeInfo* defaultInfo,
75 const G4String& beamPipeTypeIn,
76 G4double aper1In,
77 G4double aper2In,
78 G4double aper3In,
79 G4double aper4In,
80 const G4String& vacuumMaterialIn,
81 G4double beamPipeThicknessIn,
82 const G4String& beamPipeMaterialIn,
83 const G4ThreeVector& inputFaceNormalIn = G4ThreeVector(0,0,-1),
84 const G4ThreeVector& outputFaceNormalIn = G4ThreeVector(0,0,1));
85
88 void CheckApertureInfo();
89
92 BDSExtent Extent() const;
93
95 BDSExtent ExtentInner() const;
96
98 G4double IndicativeRadius() const;
99
101 G4double IndicativeRadiusInner() const;
102
105 G4double aper1;
106 G4double aper2;
107 G4double aper3;
108 G4double aper4;
109 G4double aperOffsetX;
110 G4double aperOffsetY;
111 G4Material* vacuumMaterial;
113 G4Material* beamPipeMaterial;
114 G4ThreeVector inputFaceNormal;
115 G4ThreeVector outputFaceNormal;
117 G4String pointsUnit;
119
120private:
122 void CheckAndSetPointsInfo(const G4String& beamPipeTypeIn);
123
128 void CheckRequiredParametersSet(G4bool setAper1, G4bool setAper2,
129 G4bool setAper3, G4bool setAper4);
130
132 void InfoOKForCircular();
133
135 void InfoOKForElliptical();
136
139
141 void InfoOKForLHC();
142
145
148
150 void InfoOKForRaceTrack();
151
153 void InfoOKForOctagonal();
154
156 void InfoOKForClicPCL();
157};
158
159#endif
Holder class for all information required to describe a beam pipe model.
G4double IndicativeRadiusInner() const
Return an indicative inner extent for the beam pipe vacuum.
BDSExtent Extent() const
G4Material * beamPipeMaterial
Public member for direct access.
G4double aper3
Public member for direct access.
G4double aperOffsetX
Public member for direct access.
G4String pointsUnit
Public member for direct access.
BDSBeamPipeInfo()=delete
Deleted default constructor to ensure one of supplied constructors is used.
void InfoOKForLHC()
Aperture info check for lhc aperture.
G4double aper1
Public member for direct access.
G4ThreeVector outputFaceNormal
Public member for direct access.
void InfoOKForOctagonal()
Aperture info check for octagonal aperture.
G4Material * vacuumMaterial
Public member for direct access.
void CheckAndSetPointsInfo(const G4String &beamPipeTypeIn)
Parse the type string to extract the file name and the optional units and assign to member variables.
void InfoOKForClicPCL()
Aperture info check for CLIC PCL aperture.
G4double beamPipeThickness
Public member for direct access.
G4double aperOffsetY
Public member for direct access.
void InfoOKForRaceTrack()
Aperture info check for racetrack aperture.
void InfoOKForElliptical()
Aperture info check for elliptical aperture.
BDSExtent ExtentInner() const
Return an extent for just the raw aperture.
void InfoOKForRectEllipse()
Aperture info check for rectellipse aperture.
G4String pointsFileName
Public member for direct access.
G4double aper4
Public member for direct access.
void InfoOKForLHCDetailed()
Aperture info check for lhc detailed aperture.
void InfoOKForCircular()
Aperture info check for circular aperture.
void CheckRequiredParametersSet(G4bool setAper1, G4bool setAper2, G4bool setAper3, G4bool setAper4)
void InfoOKForRectangular()
Aperture info check for rectangular aperture.
G4double IndicativeRadius() const
Return an indicative extent of the beam pipe - typically the maximum of x or y extent.
G4double aper2
Public member for direct access.
BDSBeamPipeType beamPipeType
Public member for direct access.
G4ThreeVector inputFaceNormal
Public member for direct access.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39