BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSGeometryFactoryBase.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 "BDSColours.hh"
20#include "BDSColourFromMaterial.hh"
21#include "BDSGeometryExternal.hh"
22#include "BDSGeometryFactoryBase.hh"
23#include "BDSGlobalConstants.hh"
24#include "BDSSDType.hh"
25#include "BDSUtilities.hh"
26
27#include "globals.hh"
28#include "G4AssemblyVolume.hh"
29#include "G4Colour.hh"
30#include "G4LogicalVolume.hh"
31#include "G4RotationMatrix.hh"
32#include "G4String.hh"
33#include "G4VisAttributes.hh"
34#include "G4VPhysicalVolume.hh"
35
36#include <algorithm>
37#include <map>
38#include <set>
39#include <vector>
40
41BDSGeometryFactoryBase::BDSGeometryFactoryBase()
42{
44}
45
46BDSGeometryFactoryBase::~BDSGeometryFactoryBase()
47{;}
48
49std::set<G4VisAttributes*> BDSGeometryFactoryBase::ApplyColourMapping(std::set<G4LogicalVolume*>& lvsIn,
50 std::map<G4String, G4Colour*>* mapping,
51 G4bool autoColour,
52 const G4String& prefixToStripFromName)
53{
54 std::set<G4VisAttributes*> visAttributes; // empty set
55
56 // no mapping, just return.
57 if (!mapping && !autoColour)
58 {return visAttributes;}
59
60 if (mapping)
61 {
62 if (mapping->size() == 1)
63 {// only one colour for all - simpler
64 G4VisAttributes* vis = new G4VisAttributes(*BDSColours::Instance()->GetColour("gdml"));
65 vis->SetVisibility(true);
66 visAttributes.insert(vis);
67 for (auto lv : lvsIn)
68 {lv->SetVisAttributes(*vis);}
69 return visAttributes;
70 }
71
72 // else iterate over all lvs and required vis attributes
73 // prepare required vis attributes
74 std::map<G4String, G4VisAttributes*> attMap;
75 for (const auto& it : *mapping)
76 {
77 G4VisAttributes* vis = new G4VisAttributes(*(it.second));
78 vis->SetVisibility(true);
79 visAttributes.insert(vis);
80 attMap[it.first] = vis;
81 }
82
83 for (auto lv : lvsIn)
84 {// iterate over all volumes
85 const G4String& name = lv->GetName();
86 for (const auto& it : attMap)
87 {// iterate over all mappings to find first one that matches substring
88 if (BDS::StrContains(name, it.first))
89 {
90 lv->SetVisAttributes(it.second);
91 break;
92 }
93 }
94 }
95 }
96
97 if (autoColour)
98 {
99 for (auto lv : lvsIn)
100 {
101 const G4VisAttributes* existingVis = lv->GetVisAttributes();
102 if (!existingVis)
103 {
104 G4Colour* c = BDSColourFromMaterial::Instance()->GetColour(lv->GetMaterial(), prefixToStripFromName);
105 G4VisAttributes* vis = new G4VisAttributes(*c);
106 vis->SetVisibility(true);
107 visAttributes.insert(vis);
108 lv->SetVisAttributes(vis);
109 }
110 }
111 }
112
113 return visAttributes;
114}
115
116void BDSGeometryFactoryBase::ApplyUserLimits(const std::set<G4LogicalVolume*>& lvsIn,
117 G4UserLimits* userLimits)
118{
119 for (auto& lv : lvsIn)
120 {lv->SetUserLimits(userLimits);}
121}
122
124 const std::set<G4LogicalVolume*>& allLogicalVolumesIn,
125 BDSSDType generalSensitivity,
126 const std::set<G4LogicalVolume*>& vacuumLogicalVolumes,
127 BDSSDType vacuumSensitivity)
128{
129 std::map<G4LogicalVolume*, BDSSDType> sensitivityMapping;
130 std::set<G4LogicalVolume*> notVacuumVolumes;
131 std::set_difference(allLogicalVolumesIn.begin(),
132 allLogicalVolumesIn.end(),
133 vacuumLogicalVolumes.begin(),
134 vacuumLogicalVolumes.end(),
135 std::inserter(notVacuumVolumes, notVacuumVolumes.begin()));
136 for (auto lv : notVacuumVolumes)
137 {sensitivityMapping[lv] = generalSensitivity;}
138 for (auto lv : vacuumLogicalVolumes)
139 {sensitivityMapping[lv] = vacuumSensitivity;}
140 result->RegisterSensitiveVolume(sensitivityMapping);
141}
142
144{
145 CleanUpBase();
146}
147
148G4String BDSGeometryFactoryBase:: PreprocessedName(const G4String& objectName,
149 const G4String& /*acceleratorComponentName*/) const
150{return objectName;}
151
152std::set<G4LogicalVolume*> BDSGeometryFactoryBase::GetVolumes(const std::set<G4LogicalVolume*>& allLVs,
153 std::vector<G4String>* volumeNames,
154 G4bool preprocessFile,
155 const G4String& componentName) const
156{
157 if (!volumeNames)
158 {return std::set<G4LogicalVolume*>();}
159
160 std::vector<G4String> expectedVolumeNames;
161 if (preprocessFile)
162 {
163 // transform the names to those the preprocessor would produce
164 expectedVolumeNames.resize(volumeNames->size());
165 std::transform(volumeNames->begin(),
166 volumeNames->end(),
167 expectedVolumeNames.begin(),
168 [&](const G4String& n){return PreprocessedName(n, componentName);});
169 }
170 else
171 {expectedVolumeNames = *volumeNames;}
172
173 std::set<G4LogicalVolume*> volsMatched;
174 for (const auto& en : expectedVolumeNames)
175 {
176 for (auto lv : allLVs)
177 {
178 if (lv->GetName() == en)
179 {volsMatched.insert(lv);}
180 }
181 }
182 return volsMatched;
183}
184
186{
188 xmin = 0;
189 xmax = 0;
190 ymin = 0;
191 ymax = 0;
192 zmin = 0;
193 zmax = 0;
194}
195
197{
198 xmin = std::min(ext.XNeg(), xmin);
199 xmax = std::max(ext.XPos(), xmax);
200 ymin = std::min(ext.YNeg(), ymin);
201 ymax = std::max(ext.YPos(), ymax);
202 zmin = std::min(ext.ZNeg(), zmin);
203 zmax = std::max(ext.ZPos(), zmax);
204}
205
206void BDSGeometryFactoryBase::ExpandExtent(G4double x0, G4double rx,
207 G4double y0, G4double ry,
208 G4double z0, G4double rz)
209{
210 xmin = std::min(x0-rx, xmin);
211 xmax = std::max(x0+rx, xmax);
212 ymin = std::min(y0-ry, ymin);
213 ymax = std::max(y0+ry, ymax);
214 zmin = std::min(z0-rz, zmin);
215 zmax = std::max(z0+rz, zmax);
216}
217
218void BDSGeometryFactoryBase::ExpandExtent(G4double x0, G4double xLow, G4double xHigh,
219 G4double y0, G4double yLow, G4double yHigh,
220 G4double z0, G4double zLow, G4double zHigh)
221{
222 xmin = std::min(x0+xLow, xmin);
223 xmax = std::max(x0+xHigh, xmax);
224 ymin = std::min(y0+yLow, ymin);
225 ymax = std::max(y0+yHigh, ymax);
226 zmin = std::min(z0+zLow, zmin);
227 zmax = std::max(z0+zHigh, zmax);
228}
static BDSColourFromMaterial * Instance()
Singleton pattern.
G4Colour * GetColour(const G4Material *material, const G4String &prefixToStripFromName="")
Get colour from name.
static BDSColours * Instance()
singleton pattern
Definition: BDSColours.cc:33
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double XPos() const
Accessor.
Definition: BDSExtent.hh:66
G4double ZPos() const
Accessor.
Definition: BDSExtent.hh:70
G4double XNeg() const
Accessor.
Definition: BDSExtent.hh:65
G4double YNeg() const
Accessor.
Definition: BDSExtent.hh:67
G4double YPos() const
Accessor.
Definition: BDSExtent.hh:68
G4double ZNeg() const
Accessor.
Definition: BDSExtent.hh:69
virtual void FactoryBaseCleanUp()
Empty containers for next use - factories are never deleted so can't rely on scope.
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)
A loaded piece of externally provided geometry.
G4double ymax
Extent in one dimension.
virtual void ApplySensitivity(BDSGeometryExternal *result, const std::set< G4LogicalVolume * > &allLogicalVolumesIn, BDSSDType generalSensitivity, const std::set< G4LogicalVolume * > &vacuumLogicalVolumes, BDSSDType vacuumSensitivity)
Attach the relevant general and vacuum sensitivity to each volume.
void ExpandExtent(const BDSExtent &extent)
Expand the acuumulated extents using a (possibly asymmetric) extent instance.
virtual G4String PreprocessedName(const G4String &objectName, const G4String &acceleratorComponentName) const
std::set< G4LogicalVolume * > GetVolumes(const std::set< G4LogicalVolume * > &allLVs, std::vector< G4String > *volumeNames, G4bool preprocessGDML, const G4String &componentName) const
Get the volumes that match the name. Volume names are matched exactly and are case sensitive.
virtual void CleanUp()
Virtual clean up that derived classes can override that calls CleanUpBase().
G4double xmin
Extent in one dimension.
G4double ymin
Extent in one dimension.
G4double zmin
Extent in one dimension.
G4double zmax
Extent in one dimension.
G4double xmax
Extent in one dimension.
virtual void ApplyUserLimits(const std::set< G4LogicalVolume * > &lvsIn, G4UserLimits *userLimits)
Attach a set of user limits to every logical volume supplied.
virtual std::set< G4VisAttributes * > ApplyColourMapping(std::set< G4LogicalVolume * > &lvs, std::map< G4String, G4Colour * > *mapping, G4bool autoColour, const G4String &preprocessPrefix="")
Improve type-safety of native enum data type in C++.
G4bool StrContains(const G4String &str, const G4String &test)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
Definition: BDSUtilities.cc:66