BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipeFactoryBase.hh
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#ifndef BDSBEAMPIPEFACTORYBASE_H
20#define BDSBEAMPIPEFACTORYBASE_H
21
22#include "BDSBeamPipe.hh"
23#include "BDSFactoryBase.hh"
24
25#include "globals.hh" // geant4 globals / types
26#include "G4RotationMatrix.hh"
27
28#include <set>
29
30class G4LogicalVolume;
31class G4Material;
32class G4PVPlacement;
33class G4UserLimits;
34class G4VSolid;
35
58{
59public:
61 virtual BDSBeamPipe* CreateBeamPipe(const G4String& nameIn, // name
62 G4double lengthIn, // length [mm]
63 G4double aper1 = 0, // aperture parameter 1
64 G4double aper2 = 0, // aperture parameter 2
65 G4double aper3 = 0, // aperture parameter 3
66 G4double aper4 = 0, // aperture parameter 4
67 G4Material* vacuumMaterialIn = nullptr, // vacuum material
68 G4double beamPipeThicknessIn = 0, // beampipe thickness [mm]
69 G4Material* beamPipeMaterialIn = nullptr,// beampipe material
70 const G4String& pointsFileIn = "",
71 const G4String& pointsUnitIn = ""
72 ) = 0;
73
77 virtual BDSBeamPipe* CreateBeamPipe(const G4String& nameIn,
78 G4double lengthIn,
79 const G4ThreeVector& inputFaceNormalIn,
80 const G4ThreeVector& outputFaceNormalIn,
81 G4double aper1 = 0,
82 G4double aper2 = 0,
83 G4double aper3 = 0,
84 G4double aper4 = 0,
85 G4Material* vacuumMaterialIn = nullptr,
86 G4double beamPipeThicknessIn = 0,
87 G4Material* beamPipeMaterialIn = nullptr,
88 const G4String& pointsFileIn = "",
89 const G4String& pointsUnitIn = "") = 0;
90
93
94protected:
97
100 void CleanUpBase();
101
102 virtual void CleanUp();
103
105 void CommonConstruction(const G4String& nameIn,
106 G4Material* vacuumMaterialIn,
107 G4Material* beamPipeMaterialIn,
108 G4double length);
109
112 G4double containerRadius,
113 G4bool containerIsCircular = false);
114
116 static void CheckAngledVolumeCanBeBuilt(G4double length,
117 const G4ThreeVector &inputfaceAngle,
118 const G4ThreeVector &outputfaceAngle,
119 G4double horizontalWidth,
120 const G4String& name);
121
122 // methods called by CommonConstruction, can be implemented by derived classes
123
125 virtual void BuildLogicalVolumes(const G4String& nameIn,
126 G4Material* vacuumMaterialIn,
127 G4Material* beamPipeMaterialIn);
129 virtual void SetVisAttributes(G4Material* beamPipeMaterialIn,
130 G4Material* vacuumMaterialIn);
131
133 virtual void SetUserLimits(G4double length);
134
136 virtual void PlaceComponents(const G4String& nameIn);
137
141 G4VSolid* vacuumSolid;
142 G4VSolid* beamPipeSolid;
143 G4VSolid* containerSolid;
146 G4LogicalVolume* vacuumLV;
147 G4LogicalVolume* beamPipeLV;
148 G4LogicalVolume* containerLV;
149 G4PVPlacement* vacuumPV;
150 G4PVPlacement* beamPipePV;
151
153 G4ThreeVector inputFaceNormal;
154 G4ThreeVector outputFaceNormal;
156};
157
158#endif
Abstract base class for beampipe factory classes.
G4ThreeVector outputFaceNormal
For recording the face normals in the finished pipe component.
G4bool storeApertureImpacts
Whether to store aperture impacts.
BDSBeamPipeFactoryBase()
base constructor
virtual void BuildLogicalVolumes(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn)
build logical volumes
virtual void SetVisAttributes(G4Material *beamPipeMaterialIn, G4Material *vacuumMaterialIn)
Set visual attributes.
virtual BDSBeamPipe * CreateBeamPipe(const G4String &nameIn, G4double lengthIn, const G4ThreeVector &inputFaceNormalIn, const G4ThreeVector &outputFaceNormalIn, 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
G4ThreeVector inputFaceNormal
For recording the face normals in the finished pipe component.
virtual ~BDSBeamPipeFactoryBase()
Virtual base destructor.
void CommonConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double length)
finalise beampipe construction
G4bool sensitiveVacuum
Whether the vacuum will record any energy deposition.
BDSBeamPipe * BuildBeamPipeAndRegisterVolumes(BDSExtent extent, G4double containerRadius, G4bool containerIsCircular=false)
build beampipe and register logical volumes
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
virtual void PlaceComponents(const G4String &nameIn)
Place volumes.
G4bool sensitiveBeamPipe
Whether the beam pipe will record energy deposition.
static void CheckAngledVolumeCanBeBuilt(G4double length, const G4ThreeVector &inputfaceAngle, const G4ThreeVector &outputfaceAngle, G4double horizontalWidth, const G4String &name)
check if a beam pipe volume with angled faces can be constructed
G4VSolid * containerSubtractionSolid
Longer (in length) version of container solid for unambiguous subtraction.
virtual void SetUserLimits(G4double length)
Set user limits.
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
Common temporary storage for all factories no matter what geometry.