BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipeFactory.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 "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 "BDSBeamPipeFactoryRhombus.hh"
33#include "BDSBeamPipeInfo.hh"
34#include "BDSBeamPipeType.hh"
35#include "BDSDebug.hh"
36
37#include "globals.hh" // geant4 globals / types
38
40
42{
43 if (!instance)
45 return instance;
46}
47
49{
62}
63
64BDSBeamPipeFactory::~BDSBeamPipeFactory()
65{
66 delete circular;
67 delete elliptical;
68 delete rectangular;
69 delete lhc;
70 delete lhcdetailed;
71 delete rectellipse;
72 delete racetrack;
73 delete rhombus;
74 delete octagonal;
75 delete circularvacuum;
76 delete clicpcl;
77 delete pointsfile;
78 instance = nullptr;
79}
80
82{
83 switch(type.underlying())
84 {
85 case BDSBeamPipeType::circular:
86 {return circular; break;}
87 case BDSBeamPipeType::elliptical:
88 {return elliptical; break;}
89 case BDSBeamPipeType::rectangular:
90 {return rectangular; break;}
91 case BDSBeamPipeType::lhc:
92 {return lhc; break;}
93 case BDSBeamPipeType::lhcdetailed:
94 {return lhcdetailed; break;}
95 case BDSBeamPipeType::rectellipse:
96 {return rectellipse; break;}
97 case BDSBeamPipeType::racetrack:
98 {return racetrack; break;}
99 case BDSBeamPipeType::octagonal:
100 {return octagonal; break;}
101 case BDSBeamPipeType::circularvacuum:
102 {return circularvacuum; break;}
103 case BDSBeamPipeType::clicpcl:
104 {return clicpcl; break;}
105 case BDSBeamPipeType::pointsfile:
106 {return pointsfile; break;}
107 case BDSBeamPipeType::rhombus:
108 {return rhombus; break;}
109 default:
110#ifdef BDSDEBUG
111 G4cout << __METHOD_NAME__ << "unknown type \"" << type << "\" - circular beampipe factory by default" << G4endl;
112#endif
113 return circular;
114 break;
115 }
116}
117
118BDSBeamPipe* BDSBeamPipeFactory::CreateBeamPipe(const G4String& name,
119 G4double length,
120 BDSBeamPipeInfo* bpi)
121{
122 if ((bpi->inputFaceNormal.z() > -1) || (bpi->outputFaceNormal.z() < 1))
123 {
124 return CreateBeamPipe(bpi->beamPipeType,
125 name,
126 length,
127 bpi->inputFaceNormal,
128 bpi->outputFaceNormal,
129 bpi->aper1,
130 bpi->aper2,
131 bpi->aper3,
132 bpi->aper4,
133 bpi->vacuumMaterial,
135 bpi->beamPipeMaterial,
136 bpi->pointsFileName,
137 bpi->pointsUnit);
138 }
139 else
140 {
141 return CreateBeamPipe(bpi->beamPipeType,
142 name,
143 length,
144 bpi->aper1,
145 bpi->aper2,
146 bpi->aper3,
147 bpi->aper4,
148 bpi->vacuumMaterial,
150 bpi->beamPipeMaterial,
151 bpi->pointsFileName,
152 bpi->pointsUnit);
153 }
154}
155
156BDSBeamPipe* BDSBeamPipeFactory::CreateBeamPipe(BDSBeamPipeType beamPipeType,
157 const G4String& name,
158 G4double length,
159 G4double aper1,
160 G4double aper2,
161 G4double aper3,
162 G4double aper4,
163 G4Material* vacuumMaterial,
164 G4double beamPipeThickness,
165 G4Material* beamPipeMaterial,
166 const G4String& pointsFileIn,
167 const G4String& pointsUnitIn)
168{
169 BDSBeamPipeFactoryBase* factory = GetAppropriateFactory(beamPipeType);
170 return factory->CreateBeamPipe(name,length,aper1,aper2,aper3,aper4,
171 vacuumMaterial,beamPipeThickness,beamPipeMaterial,
172 pointsFileIn,pointsUnitIn);
173}
174
175BDSBeamPipe* BDSBeamPipeFactory::CreateBeamPipe(BDSBeamPipeType beamPipeType,
176 const G4String& name,
177 G4double length,
178 const G4ThreeVector& inputFaceNormal,
179 const G4ThreeVector& outputFaceNormal,
180 G4double aper1,
181 G4double aper2,
182 G4double aper3,
183 G4double aper4,
184 G4Material* vacuumMaterial,
185 G4double beamPipeThickness,
186 G4Material* beamPipeMaterial,
187 const G4String& pointsFileIn,
188 const G4String& pointsUnitIn)
189{
190 BDSBeamPipeFactoryBase* factory = GetAppropriateFactory(beamPipeType);
191 return factory->CreateBeamPipe(name,length,inputFaceNormal,outputFaceNormal,aper1,
192 aper2,aper3,aper4,vacuumMaterial,beamPipeThickness,
193 beamPipeMaterial,pointsFileIn,pointsUnitIn);
194}
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.
Factory for rhombus aperture model beampipes.
The main interface for using the beam pipe factories.
BDSBeamPipeFactoryBase * racetrack
Factory instance.
BDSBeamPipeFactoryBase * rhombus
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)