BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamPipeFactoryCircularVacuum.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 "BDSBeamPipeFactoryCircularVacuum.hh"
20#include "BDSBeamPipe.hh"
21#include "BDSExtent.hh"
22#include "BDSGlobalConstants.hh"
23#include "BDSUtilities.hh"
24
25#include "globals.hh" // geant4 globals / types
26#include "G4CutTubs.hh"
27#include "G4LogicalVolume.hh"
28#include "G4ThreeVector.hh"
29#include "G4Tubs.hh"
30#include "G4UserLimits.hh"
31#include "G4VSolid.hh"
32
33#include <cmath>
34#include <set>
35#include <utility>
36
37class G4Material;
38
39BDSBeamPipeFactoryCircularVacuum::BDSBeamPipeFactoryCircularVacuum()
40{;}
41
43 G4double lengthIn,
44 G4double aper1In,
45 G4double /*aper2In*/,
46 G4double /*aper3In*/,
47 G4double /*aper4In*/,
48 G4Material* vacuumMaterialIn,
49 G4double /*beamPipeThicknessIn*/,
50 G4Material* /*beamPipeMaterialIn*/,
51 const G4String& /*pointsFileIn*/,
52 const G4String& /*pointsUnitIn*/)
53{
54 // clean up after last usage
55 CleanUp();
56
57 containerSolid = new G4Tubs(nameIn + "_container_solid", // name
58 0, // inner radius
59 aper1In, // outer radius
60 lengthIn*0.5, // half length
61 0, // rotation start angle
62 CLHEP::twopi); // rotation finish angle
63
64 G4double containerRadius = aper1In + lengthSafetyLarge;
65
66 return CommonFinalConstruction(nameIn, vacuumMaterialIn, lengthIn, containerRadius);
67}
68
70 G4double lengthIn,
71 const G4ThreeVector& inputFaceNormalIn,
72 const G4ThreeVector& outputFaceNormalIn,
73 G4double aper1In,
74 G4double /*aper2In*/,
75 G4double /*aper3In*/,
76 G4double /*aper4In */,
77 G4Material* vacuumMaterialIn,
78 G4double /*beamPipeThicknessIn*/,
79 G4Material* /*beamPipeMaterialIn*/,
80 const G4String& /*pointsFileIn*/,
81 const G4String& /*pointsUnitIn*/)
82{
83 // clean up after last usage
84 CleanUp();
85 inputFaceNormal = inputFaceNormalIn;
86 outputFaceNormal = outputFaceNormalIn;
87
88 // check faces of angled container volume don't intersect - if it can be built, remaining angled volumes can be built
90
91 containerSolid = new G4CutTubs(nameIn + "_container_solid", // name
92 0, // inner radius
93 aper1In, // outer radius
94 lengthIn*0.5, // half length
95 0, // rotation start angle
96 CLHEP::twopi, // rotation finish angle
97 inputFaceNormal, // input face normal
98 outputFaceNormal); // rotation finish angle
99
100 G4double containerRadius = aper1In + lengthSafetyLarge;
101
102 return CommonFinalConstruction(nameIn, vacuumMaterialIn, lengthIn, containerRadius);
103}
104
108 G4Material* vacuumMaterialIn,
109 G4double lengthIn,
110 G4double containerRadiusIn)
111{
112 // prepare a longer container subtraction solid
113 // doesn't have to be angled as it's only used for transverse subtraction
114 containerSubtractionSolid = new G4Tubs(nameIn + "_container_sub_solid",// name
115 0, // inner radius
116 containerRadiusIn, // outer radius
117 lengthIn*4, // full length for unambiguous subtraction
118 0, // rotation start angle
119 CLHEP::twopi); // rotation finish angle
120
121 // We don't use BDSBeamPipeFactoryBase::CommonConstruction as this factory doesn't
122 // produce a beam pipe solid etc - only vacuum so it wouldn't work.
123
124 containerLV = new G4LogicalVolume(containerSolid,
125 vacuumMaterialIn,
126 nameIn + "_container_lv");
127
128 containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->ContainerVisAttr());
129 vacuumLV = containerLV; // copy pointer for referencing in BuildBeamPipeAndRegisterVolumes.
130
131 // user limits
132 auto defaultUL = BDSGlobalConstants::Instance()->DefaultUserLimits();
133 //copy the default and update with the length of the object rather than the default 1m
134 G4UserLimits* ul = BDS::CreateUserLimits(defaultUL, lengthIn);
135 if (ul != defaultUL) // if it's not the default register it
136 {allUserLimits.insert(ul);}
137 containerLV->SetUserLimits(ul);
138 vacuumLV->SetUserLimits(ul);
139
140 // record extents
141 BDSExtent ext = BDSExtent(containerRadiusIn, containerRadiusIn, lengthIn*0.5);
142
143 return BDSBeamPipeFactoryBase::BuildBeamPipeAndRegisterVolumes(ext, containerRadiusIn);
144}
G4ThreeVector outputFaceNormal
For recording the face normals in the finished pipe component.
G4ThreeVector inputFaceNormal
For recording the face normals in the finished pipe component.
BDSBeamPipe * BuildBeamPipeAndRegisterVolumes(BDSExtent extent, G4double containerRadius, G4bool containerIsCircular=false)
build beampipe and register logical volumes
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 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="")
create a flat ended beampipe
BDSBeamPipe * CommonFinalConstruction(const G4String &nameIn, G4Material *vacuumMaterialIn, G4double lengthIn, G4double containerRadiusIn)
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double lengthSafetyLarge
Cache of global constants variable.
static BDSGlobalConstants * Instance()
Access method.
G4UserLimits * CreateUserLimits(G4UserLimits *defaultUL, G4double length, G4double fraction=1.6)