BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamPipeFactory.cc
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#include "BDSBeamPipeFactory.hh"
20#include "BDSBeamPipeFactoryBase.hh"
21#include "BDSBeamPipeFactoryCircular.hh"
22#include "BDSBeamPipeFactoryCircularVacuum.hh"
23#include "BDSBeamPipeFactoryClicPCL.hh"
24#include "BDSBeamPipeFactoryElliptical.hh"
25#include "BDSBeamPipeFactoryRectangular.hh"
26#include "BDSBeamPipeFactoryLHC.hh"
27#include "BDSBeamPipeFactoryLHCDetailed.hh"
28#include "BDSBeamPipeFactoryOctagonal.hh"
29#include "BDSBeamPipeFactoryPointsFile.hh"
30#include "BDSBeamPipeFactoryRaceTrack.hh"
31#include "BDSBeamPipeFactoryRectEllipse.hh"
32#include "BDSBeamPipeInfo.hh"
33#include "BDSBeamPipeType.hh"
34#include "BDSDebug.hh"
35
36#include "globals.hh" // geant4 globals / types
37
39
41{
42 if (!instance)
44 return instance;
45}
46
48{
60}
61
62BDSBeamPipeFactory::~BDSBeamPipeFactory()
63{
64 delete circular;
65 delete elliptical;
66 delete rectangular;
67 delete lhc;
68 delete lhcdetailed;
69 delete rectellipse;
70 delete racetrack;
71 delete octagonal;
72 delete circularvacuum;
73 delete clicpcl;
74 delete pointsfile;
75 instance = nullptr;
76}
77
79{
80 switch(type.underlying())
81 {
82 case BDSBeamPipeType::circular:
83 {return circular; break;}
84 case BDSBeamPipeType::elliptical:
85 {return elliptical; break;}
86 case BDSBeamPipeType::rectangular:
87 {return rectangular; break;}
88 case BDSBeamPipeType::lhc:
89 {return lhc; break;}
90 case BDSBeamPipeType::lhcdetailed:
91 {return lhcdetailed; break;}
92 case BDSBeamPipeType::rectellipse:
93 {return rectellipse; break;}
94 case BDSBeamPipeType::racetrack:
95 {return racetrack; break;}
96 case BDSBeamPipeType::octagonal:
97 {return octagonal; break;}
98 case BDSBeamPipeType::circularvacuum:
99 {return circularvacuum; break;}
100 case BDSBeamPipeType::clicpcl:
101 {return clicpcl; break;}
102 case BDSBeamPipeType::pointsfile:
103 {return pointsfile; break;}
104 default:
105#ifdef BDSDEBUG
106 G4cout << __METHOD_NAME__ << "unknown type \"" << type << "\" - circular beampipe factory by default" << G4endl;
107#endif
108 return circular;
109 break;
110 }
111}
112
113BDSBeamPipe* BDSBeamPipeFactory::CreateBeamPipe(const G4String& name,
114 G4double length,
115 BDSBeamPipeInfo* bpi)
116{
117 if ((bpi->inputFaceNormal.z() > -1) || (bpi->outputFaceNormal.z() < 1))
118 {
119 return CreateBeamPipe(bpi->beamPipeType,
120 name,
121 length,
122 bpi->inputFaceNormal,
123 bpi->outputFaceNormal,
124 bpi->aper1,
125 bpi->aper2,
126 bpi->aper3,
127 bpi->aper4,
128 bpi->vacuumMaterial,
130 bpi->beamPipeMaterial,
131 bpi->pointsFileName,
132 bpi->pointsUnit);
133 }
134 else
135 {
136 return CreateBeamPipe(bpi->beamPipeType,
137 name,
138 length,
139 bpi->aper1,
140 bpi->aper2,
141 bpi->aper3,
142 bpi->aper4,
143 bpi->vacuumMaterial,
145 bpi->beamPipeMaterial,
146 bpi->pointsFileName,
147 bpi->pointsUnit);
148 }
149}
150
151BDSBeamPipe* BDSBeamPipeFactory::CreateBeamPipe(BDSBeamPipeType beamPipeType,
152 const G4String& name,
153 G4double length,
154 G4double aper1,
155 G4double aper2,
156 G4double aper3,
157 G4double aper4,
158 G4Material* vacuumMaterial,
159 G4double beamPipeThickness,
160 G4Material* beamPipeMaterial,
161 const G4String& pointsFileIn,
162 const G4String& pointsUnitIn)
163{
164 BDSBeamPipeFactoryBase* factory = GetAppropriateFactory(beamPipeType);
165 return factory->CreateBeamPipe(name,length,aper1,aper2,aper3,aper4,
166 vacuumMaterial,beamPipeThickness,beamPipeMaterial,
167 pointsFileIn,pointsUnitIn);
168}
169
170BDSBeamPipe* BDSBeamPipeFactory::CreateBeamPipe(BDSBeamPipeType beamPipeType,
171 const G4String& name,
172 G4double length,
173 const G4ThreeVector& inputFaceNormal,
174 const G4ThreeVector& outputFaceNormal,
175 G4double aper1,
176 G4double aper2,
177 G4double aper3,
178 G4double aper4,
179 G4Material* vacuumMaterial,
180 G4double beamPipeThickness,
181 G4Material* beamPipeMaterial,
182 const G4String& pointsFileIn,
183 const G4String& pointsUnitIn)
184{
185 BDSBeamPipeFactoryBase* factory = GetAppropriateFactory(beamPipeType);
186 return factory->CreateBeamPipe(name,length,inputFaceNormal,outputFaceNormal,aper1,
187 aper2,aper3,aper4,vacuumMaterial,beamPipeThickness,
188 beamPipeMaterial,pointsFileIn,pointsUnitIn);
189}
Abstract base class for beampipe factory classes.
virtual BDSBeamPipe * CreateBeamPipe(const G4String &nameIn, G4double lengthIn, G4double aper1=0, G4double aper2=0, G4double aper3=0, G4double aper4=0, G4Material *vacuumMaterialIn=nullptr, G4double beamPipeThicknessIn=0, G4Material *beamPipeMaterialIn=nullptr, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")=0
create a flat ended beampipe
Factory for vacuum-only circular volumes.
Factory for circular beam pipes.
Factory for CLIC post collision line beam pipes.
Factory for elliptical beam pipes.
Factory for detailed lhc aperture model beam pipes.
Factory for simple lhc aperture model beam pipes.
Factory for octagonal aperture model beampipes.
Factory for beam pipe made from an list of x,y points in a file.
Factory for racetrack aperture model beampipes.
Factory for rectellipse aperture model beampipes.
Factory for rectangular beam pipes.
The main interface for using the beam pipe factories.
BDSBeamPipeFactoryBase * racetrack
Factory instance.
BDSBeamPipeFactoryBase * elliptical
Factory instance.
static BDSBeamPipeFactory * instance
Singleton instance pointer.
BDSBeamPipeFactoryBase * circularvacuum
Factory instance.
BDSBeamPipeFactoryBase * GetAppropriateFactory(BDSBeamPipeType beamPipeTypeIn)
Return the appropriate factory singleton pointer given a type.
static BDSBeamPipeFactory * Instance()
Singleton accessor.
BDSBeamPipeFactoryBase * rectellipse
Factory instance.
BDSBeamPipeFactoryBase * lhc
Factory instance.
BDSBeamPipeFactoryBase * clicpcl
Factory instance.
BDSBeamPipeFactoryBase * pointsfile
Factory instance.
BDSBeamPipeFactoryBase * lhcdetailed
Factory instance.
BDSBeamPipeFactory()
Private constructor as singleton pattern.
BDSBeamPipeFactoryBase * rectangular
Factory instance.
BDSBeamPipeFactoryBase * octagonal
Factory instance.
BDSBeamPipeFactoryBase * circular
Factory instance.
Holder class for all information required to describe a beam pipe model.
G4Material * beamPipeMaterial
Public member for direct access.
G4double aper3
Public member for direct access.
G4String pointsUnit
Public member for direct access.
G4double aper1
Public member for direct access.
G4ThreeVector outputFaceNormal
Public member for direct access.
G4Material * vacuumMaterial
Public member for direct access.
G4double beamPipeThickness
Public member for direct access.
G4String pointsFileName
Public member for direct access.
G4double aper4
Public member for direct access.
G4double aper2
Public member for direct access.
BDSBeamPipeType beamPipeType
Public member for direct access.
G4ThreeVector inputFaceNormal
Public member for direct access.
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45
type underlying() const
return underlying value (can be used in switch statement)