BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBLMFactory.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 "BDSBLM.hh"
20#include "BDSBLMFactory.hh"
21#include "BDSBLMType.hh"
22#include "BDSDebug.hh"
23#include "BDSException.hh"
24#include "BDSExtent.hh"
25#include "BDSGeometryExternal.hh"
26#include "BDSGeometryFactory.hh"
27#include "BDSMaterials.hh"
28#include "BDSUtilities.hh"
29
30#include "parser/blmplacement.h"
31
32#include "globals.hh"
33#include "G4Box.hh"
34#include "G4LogicalVolume.hh"
35#include "G4Orb.hh"
36#include "G4Tubs.hh"
37
38#include <map>
39
40BDSBLMFactory::BDSBLMFactory()
41{;}
42
43BDSBLMFactory::~BDSBLMFactory()
44{;}
45
46void BDSBLMFactory::PositiveFinite(G4double value,
47 const G4String& parameterName,
48 const G4String& blmName) const
49{
50 if (!BDS::IsFinite(value) || value < 0)
51 {throw BDSException(__METHOD_NAME__, "\"" + parameterName + "\" must be +ve finite for blm named \"" + blmName + "\"");}
52}
53
54BDSBLM* BDSBLMFactory::CreateBLM(const GMAD::BLMPlacement& bp,
55 G4VSensitiveDetector* sd)
56{
57 return CreateBLM(G4String(bp.name),
58 G4String(bp.geometryFile),
59 G4String(bp.geometryType),
60 G4String(bp.blmMaterial),
61 bp.blm1 * CLHEP::m,
62 bp.blm2 * CLHEP::m,
63 bp.blm3 * CLHEP::m,
64 bp.blm4 * CLHEP::m,
65 sd,
66 bp.bias);
67}
68
69BDSBLM* BDSBLMFactory::CreateBLM(const G4String& name,
70 const G4String& geometryFile,
71 const G4String& geometryType,
72 const G4String& material,
73 G4double blm1,
74 G4double blm2,
75 G4double blm3,
76 G4double /*blm4*/,
77 G4VSensitiveDetector* sd,
78 const G4String& bias)
79{
80 // if geometry file is specified then we load the external file.
81 BDSBLM* result = nullptr;
82 if (!geometryFile.empty())
83 {
84 BDSGeometryExternal* blmGeom = BDSGeometryFactory::Instance()->BuildGeometry(name, geometryFile);
85 result = new BDSBLM(blmGeom);
86 }
87 else
88 {
89 BDSBLMType shape = BDS::DetermineBLMType(geometryType);
90 switch (shape.underlying())
91 {
92 case BDSBLMType::cylinder:
93 {result = BuildBLMCylinder(name, material, blm1, blm2); break;}
94 case BDSBLMType::cube:
95 {result = BuildBLMCube(name, material, blm1, blm2, blm3); break;}
96 case BDSBLMType::sphere:
97 {result = BuildBLMSphere(name, material, blm1); break;}
98 default:
99 {break;}
100 }
101 }
102 if (!result)
103 {return result;}
104
105 // attach sensitivity
106 for (auto lv : result->GetAllLogicalVolumes())
107 {lv->SetSensitiveDetector(sd);}
108 result->GetContainerLogicalVolume()->SetSensitiveDetector(sd); // attach separately to the container
109
110 // biasing - name of objects attached here and built later
111 result->SetBias(bias);
112 return result;
113}
114
116 const G4String& material,
117 G4double halfLength,
118 G4double radius)
119{
120 PositiveFinite(halfLength, "blm1 (half length)", name);
121 PositiveFinite(radius, "blm2 (radius)", name);
122 G4Tubs* shape = new G4Tubs(name + "_solid",
123 0,
124 radius,
125 halfLength,
126 0,
127 CLHEP::twopi);
128 BDSExtent ext = BDSExtent(radius, radius, halfLength);
129 return CommonConstruction(name, material, shape, ext);
130}
131
133 const G4String& material,
134 G4double halfLengthX,
135 G4double halfLengthY,
136 G4double halfLengthZ)
137{
138 PositiveFinite(halfLengthX, "blm1 (half length in x)", name);
139 PositiveFinite(halfLengthY, "blm2 (half length in y)", name);
140 PositiveFinite(halfLengthZ, "blm3 (half length in z)", name);
141 G4Box* shape = new G4Box(name + "_solid", halfLengthX, halfLengthY, halfLengthZ);
142 BDSExtent ext = BDSExtent(halfLengthX, halfLengthY, halfLengthZ);
143 return CommonConstruction(name, material, shape, ext);
144}
145
147 const G4String& material,
148 G4double radius)
149{
150 PositiveFinite(radius, "blm1 (radius)", name);
151 G4Orb* shape = new G4Orb(name + "_solid", radius);
152 BDSExtent ext = BDSExtent(radius, radius, radius);
153 return CommonConstruction(name, material, shape, ext);
154}
155
157 const G4String& material,
158 G4VSolid* shape,
159 const BDSExtent& extent)
160{
161 G4Material* mat = BDSMaterials::Instance()->GetMaterial(material);
162 G4LogicalVolume* lv = new G4LogicalVolume(shape, mat, name + "_lv");
163
164 BDSBLM* result = new BDSBLM(shape, lv, extent);
165 return result;
166}
BDSBLM * BuildBLMCylinder(const G4String &name, const G4String &material, G4double halfLength, G4double radius)
Build the geometry for a cylinder.
BDSBLM * BuildBLMCube(const G4String &name, const G4String &material, G4double halfLengthX, G4double halfLengthY, G4double halfLengthZ)
Build the geometry for a cube.
BDSBLM * CommonConstruction(const G4String &name, const G4String &material, G4VSolid *shape, const BDSExtent &extent)
BDSBLM * BuildBLMSphere(const G4String &name, const G4String &material, G4double radius)
Build the geometry for a sphere.
A beam loss monitor.
Definition: BDSBLM.hh:42
void SetBias(const G4String &biasIn)
Definition: BDSBLM.hh:62
General exception with possible name of object and message.
Definition: BDSException.hh:35
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4LogicalVolume * GetContainerLogicalVolume() const
Accessor - see member for more info.
virtual std::set< G4LogicalVolume * > GetAllLogicalVolumes() const
Access all logical volumes belonging to this component.
A loaded piece of externally provided geometry.
static BDSGeometryFactory * Instance()
Singleton accessor.
BDSGeometryExternal * BuildGeometry(G4String componentName, const G4String &formatAndFilePath, std::map< G4String, G4Colour * > *colourMapping=nullptr, G4bool autoColour=true, G4double suggestedLength=0, G4double suggestedHorizontalWidth=0, std::vector< G4String > *namedVacuumVolumes=nullptr, G4bool makeSensitive=true, BDSSDType sensitivityType=BDSSDType::energydep, G4bool stripOuterVolumeAndMakeAssembly=false, G4UserLimits *userLimitsToAttachToAllLVs=nullptr, G4bool dontReloadGeometry=false)
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
G4Material * GetMaterial(G4String material) const
Get material by name.
Improve type-safety of native enum data type in C++.
type underlying() const
return underlying value (can be used in switch statement)
blm for parser.
Definition: blmplacement.h:39
std::string name
Name of this samplerplacement.
Definition: blmplacement.h:41
BDSBLMType DetermineBLMType(G4String blmType)
function that gives corresponding enum value for string (case-insensitive)
Definition: BDSBLMType.cc:39
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())