BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamPipeFactoryBase.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 "BDSBeamPipeFactoryBase.hh"
20#include "BDSColours.hh"
21#include "BDSColourFromMaterial.hh"
22#include "BDSDebug.hh"
23#include "BDSException.hh"
24#include "BDSGlobalConstants.hh"
25#include "BDSMaterials.hh"
26#include "BDSSDType.hh"
27#include "BDSUtilities.hh"
28
29#include "globals.hh" // geant4 globals / types
30#include "G4LogicalVolume.hh"
31#include "G4Material.hh"
32#include "G4PVPlacement.hh"
33#include "G4ThreeVector.hh"
34#include "G4UserLimits.hh"
35#include "G4VisAttributes.hh"
36
37#include <limits>
38#include <set>
39
41{
43 sensitiveBeamPipe = g->SensitiveBeamPipe();
44 sensitiveVacuum = g->StoreELossVacuum();
45 storeApertureImpacts = g->StoreApertureImpacts();
46 CleanUpBase(); // non-virtual call in constructor
47}
48
49void BDSBeamPipeFactoryBase::CleanUp()
50{
52}
53
55{
57 vacuumSolid = nullptr;
58 beamPipeSolid = nullptr;
59 containerSolid = nullptr;
61 vacuumLV = nullptr;
62 beamPipeLV = nullptr;
63 containerLV = nullptr;
64 vacuumPV = nullptr;
65 beamPipePV = nullptr;
66
67 inputFaceNormal = G4ThreeVector(0,0,-1);
68 outputFaceNormal = G4ThreeVector(0,0, 1);
69}
70
72 G4Material* vacuumMaterialIn,
73 G4Material* beamPipeMaterialIn,
74 G4double length)
75{
76 allSolids.insert(vacuumSolid);
77 allSolids.insert(beamPipeSolid);
79 BuildLogicalVolumes(nameIn, vacuumMaterialIn, beamPipeMaterialIn);
81 SetVisAttributes(beamPipeMaterialIn);
83 SetUserLimits(length);
85 PlaceComponents(nameIn);
86}
87
89 G4Material* vacuumMaterialIn,
90 G4Material* beamPipeMaterialIn)
91{
92 // build the logical volumes
93 vacuumLV = new G4LogicalVolume(vacuumSolid,
94 vacuumMaterialIn,
95 nameIn + "_vacuum_lv");
96
97 beamPipeLV = new G4LogicalVolume(beamPipeSolid,
98 beamPipeMaterialIn,
99 nameIn + "_beampipe_lv");
100
101 G4Material* emptyMaterial = BDSMaterials::Instance()->GetMaterial(BDSGlobalConstants::Instance()->EmptyMaterial());
102 containerLV = new G4LogicalVolume(containerSolid,
103 emptyMaterial,
104 nameIn + "_container_lv");
105 allLogicalVolumes.insert(vacuumLV);
106 allLogicalVolumes.insert(beamPipeLV);
107}
108
109void BDSBeamPipeFactoryBase::SetVisAttributes(G4Material* beamPipeMaterialIn)
110{
111 G4Colour* c = BDSColourFromMaterial::Instance()->GetColourWithDefault(beamPipeMaterialIn,
112 BDSColours::Instance()->GetColour("beampipe"));
113 G4VisAttributes* pipeVisAttr = new G4VisAttributes(*c);
114 pipeVisAttr->SetVisibility(true);
115 pipeVisAttr->SetForceLineSegmentsPerCircle(nSegmentsPerCircle);
116 allVisAttributes.insert(pipeVisAttr);
117 beamPipeLV->SetVisAttributes(pipeVisAttr);
118 vacuumLV->SetVisAttributes(BDSGlobalConstants::Instance()->ContainerVisAttr());
119 containerLV->SetVisAttributes(BDSGlobalConstants::Instance()->ContainerVisAttr());
120}
121
123{
124 auto defaultUL = BDSGlobalConstants::Instance()->DefaultUserLimits();
125 //copy the default and update with the length of the object rather than the default 1m
126 G4UserLimits* ul = BDS::CreateUserLimits(defaultUL, length);
127 if (ul != defaultUL) // if it's not the default register it
128 {allUserLimits.insert(ul);}
129 vacuumLV->SetUserLimits(ul);
130 containerLV->SetUserLimits(ul);
131
132 G4UserLimits* beamPipeUL = ul;
133 if (BDSGlobalConstants::Instance()->BeamPipeIsInfiniteAbsorber())
134 {// new beam pipe user limits, copy from updated default user limits above
135 beamPipeUL = new G4UserLimits(*ul);
136 beamPipeUL->SetUserMinEkine(std::numeric_limits<double>::max());
137 allUserLimits.insert(beamPipeUL);
138 }
139 beamPipeLV->SetUserLimits(beamPipeUL);
140}
141
142void BDSBeamPipeFactoryBase::PlaceComponents(const G4String& nameIn)
143{
144 // PLACEMENT
145 // place the components inside the container
146 // note we don't need the pointer for anything - it's registered upon construction with g4
147 vacuumPV = new G4PVPlacement(nullptr, // no rotation
148 G4ThreeVector(), // position
149 vacuumLV, // lv to be placed
150 nameIn + "_vacuum_pv", // name
151 containerLV, // mother lv to be placed in
152 false, // no boolean operation
153 0, // copy number
154 checkOverlaps); // whether to check overlaps
155
156 beamPipePV = new G4PVPlacement(nullptr, // no rotation
157 G4ThreeVector(), // position
158 beamPipeLV, // lv to be placed
159 nameIn + "_beampipe_pipe_pv", // name
160 containerLV, // mother lv to be placed in
161 false, // no boolean operation
162 0, // copy number
163 checkOverlaps); // whether to check overlaps
164 allPhysicalVolumes.insert(vacuumPV);
165 allPhysicalVolumes.insert(beamPipePV);
166}
167
169 G4double containerRadius,
170 G4bool containerIsCircular)
171{
172 // build the BDSBeamPipe instance and return it
173 BDSBeamPipe* aPipe = new BDSBeamPipe(containerSolid,containerLV,extent,
175 vacuumLV,containerIsCircular,containerRadius,
177
178 // register objects
179 aPipe->RegisterSolid(allSolids);
180 aPipe->RegisterLogicalVolume(allLogicalVolumes); //using geometry component base class method
181 aPipe->RegisterPhysicalVolume(allPhysicalVolumes);
182 if (sensitiveVacuum)
183 {aPipe->RegisterSensitiveVolume(vacuumLV, BDSSDType::energydepvacuum);}
184 if (beamPipeLV)
185 {
186 if (sensitiveBeamPipe && storeApertureImpacts)// in the case of the circular vacuum, there isn't a beampipeLV
187 {aPipe->RegisterSensitiveVolume(beamPipeLV, BDSSDType::aperturecomplete);}
188 else if (storeApertureImpacts)
189 {aPipe->RegisterSensitiveVolume(beamPipeLV, BDSSDType::apertureimpacts);}
190 else
191 {aPipe->RegisterSensitiveVolume(beamPipeLV, BDSSDType::energydep);}
192 }
193 aPipe->RegisterUserLimits(allUserLimits);
194 aPipe->RegisterVisAttributes(allVisAttributes);
195
196 return aPipe;
197}
198
200 const G4ThreeVector &inputface,
201 const G4ThreeVector &outputface,
202 G4double beampipeRadius,
203 const G4String& name)
204{
205 G4bool intersects = BDS::WillIntersect(inputface.angle(), outputface.angle(), beampipeRadius*2, length);
206
207 if (intersects)
208 {throw BDSException(__METHOD_NAME__, "a volume in the \"" + name + "\" beam pipe geometry cannot be constructed as it's angled faces intersect within it's extent");}
209}
virtual void SetVisAttributes(G4Material *beamPipeMaterialIn)
Set visual attributes.
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
G4ThreeVector inputFaceNormal
For recording the face normals in the finished pipe component.
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 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
static BDSColourFromMaterial * Instance()
Singleton pattern.
G4Colour * GetColourWithDefault(const G4Material *material, G4Colour *defaultIn) const
Get colour from name - if not found return the supplied default.
static BDSColours * Instance()
singleton pattern
Definition: BDSColours.cc:33
General exception with possible name of object and message.
Definition: BDSException.hh:35
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4bool checkOverlaps
Cache of global constants variable.
G4double nSegmentsPerCircle
Cache of global constants variable.
virtual void FactoryBaseCleanUp()
Empty containers for next use - factories are never deleted so can't rely on scope.
void RegisterLogicalVolume(G4LogicalVolume *logicalVolume)
void RegisterPhysicalVolume(G4VPhysicalVolume *physicalVolume)
void RegisterUserLimits(G4UserLimits *userLimit)
void RegisterVisAttributes(G4VisAttributes *visAttribute)
void RegisterSolid(G4VSolid *solid)
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
G4Material * GetMaterial(G4String material) const
Get material by name.
G4bool WillIntersect(const G4ThreeVector &incomingNormal, const G4ThreeVector &outgoingNormal, const G4double &zSeparation, const BDSExtent &incomingExtent, const BDSExtent &outgoingExtent)
G4UserLimits * CreateUserLimits(G4UserLimits *defaultUL, G4double length, G4double fraction=1.6)