BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSAcceleratorModel.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 "BDSAcceleratorComponentRegistry.hh"
20#include "BDSAcceleratorModel.hh"
21#include "BDSApertureInfo.hh"
22#include "BDSBeamline.hh"
23#include "BDSBeamlineSet.hh"
24#include "BDSDebug.hh"
25#include "BDSException.hh"
26#include "BDSFieldObjects.hh"
27#include "BDSLinkComponent.hh"
28#include "BDSPhysicalVolumeInfoRegistry.hh"
29#include "BDSRegion.hh"
30#include "BDSScorerHistogramDef.hh"
31#include "BDSUtilities.hh"
32
33#include "globals.hh"
34#include "G4LogicalVolume.hh"
35#include "G4VPhysicalVolume.hh"
36#include "G4VSolid.hh"
37#include "G4VUserParallelWorld.hh"
38
39#include <cstdio>
40#include <map>
41#include <set>
42#include <vector>
43
44class G4Region;
45
46BDSAcceleratorModel* BDSAcceleratorModel::instance = nullptr;
47
48BDSAcceleratorModel* BDSAcceleratorModel::Instance()
49{
50 if (!instance)
51 {instance = new BDSAcceleratorModel();}
52 return instance;
53}
54
56 worldPV(nullptr),
57 worldLV(nullptr),
58 worldSolid(nullptr),
59 tunnelBeamline(nullptr),
60 placementBeamline(nullptr),
61 blmsBeamline(nullptr)
62{
65}
66
67BDSAcceleratorModel::~BDSAcceleratorModel()
68{
69 // User feedback as deletion can take some time
70 G4cout << "BDSAcceleratorModel> Deleting model" << G4endl;
71
72 delete worldPV;
73 delete worldLV;
74 delete worldSolid;
75
76 delete tunnelBeamline;
77 delete placementBeamline;
78 delete blmsBeamline;
79
80 for (auto world : parallelWorlds)
81 {delete world;}
82
83 mainBeamlineSet.DeleteContents();
84
85 for (auto& bl : extraBeamlines)
86 {bl.second.DeleteContents();}
87
88 for (auto lc : linkComponents)
89 {delete lc;}
90
93
94 for (auto f : fields)
95 {delete f;}
96 for (auto r : regionStorage)
97 {delete r;}
98 for (auto& a : apertures)
99 {delete a.second;}
100
101 for (auto& vr : volumeRegistries)
102 {delete vr.second;}
103
104 G4cout << "BDSAcceleratorModel> Deletion complete" << G4endl;
105
106 instance = nullptr;
107}
108
110{
111 mainBeamlineSet = setIn;
112 MapBeamlineSet(setIn);
113}
114
116 const BDSBeamlineSet& setIn)
117{
118 auto search = extraBeamlines.find(name);
119 if (search != extraBeamlines.end()) // already exists!
120 {search->second.DeleteContents();} // delete pre-existing one for replacement
121 extraBeamlines[name] = setIn;
122 MapBeamlineSet(setIn);
123}
124
125const BDSBeamlineSet& BDSAcceleratorModel::BeamlineSet(const G4String& name) const
126{
127 if (name == "main")
128 {return mainBeamlineSet;}
129
130 const auto search = extraBeamlines.find(name);
131 if (search == extraBeamlines.end())
132 {throw BDSException(__METHOD_NAME__, "No such beam line set \"" + name + "\"");}
133 else
134 {return search->second;}
135}
136
138{
139 regions[region->name] = region;
140 regionStorage.insert(region);
141}
142
143void BDSAcceleratorModel::RegisterApertures(const std::map<G4String, BDSApertureInfo*>& aperturesIn)
144{
145 apertures.insert(aperturesIn.begin(), aperturesIn.end());
146}
147
149{
150 auto result = apertures.find(name);
151 if (result != apertures.end())
152 {return result->second;}
153 else
154 {throw BDSException(__METHOD_NAME__, "Invalid aperture name \"" + name + "\"");}
155}
156
157G4Region* BDSAcceleratorModel::Region(const G4String& name) const
158{
159 auto result = regions.find(name);
160 if (result != regions.end())
161 {return result->second->g4region;}
162 else
163 {
164 G4cerr << "Invalid region name \"" << name << "\"" << G4endl;
165 G4cout << "Available regions are: " << G4endl;
166 for (const auto& r : regions)
167 {G4cout << r.first << " ";}
168 throw BDSException(__METHOD_NAME__, "invalid region name.");
169 }
170}
171
172std::set<G4LogicalVolume*>* BDSAcceleratorModel::VolumeSet(const G4String& name)
173{
174 if (volumeRegistries.find(name) == volumeRegistries.end()) // no such registry -> create one
175 {volumeRegistries[name] = new std::set<G4LogicalVolume*>();}
176 return volumeRegistries[name];
177}
178
179G4bool BDSAcceleratorModel::VolumeInSet(G4LogicalVolume* volume,
180 const G4String& registryName)
181{
182 if (volumeRegistries.find(registryName) == volumeRegistries.end())
183 {return false;} // no such registry
184 else
185 {
186 auto registry = volumeRegistries[registryName];
187 return registry->find(volume) != registry->end();
188 }
189}
190
192{
193 auto result = clToMassWorldMap.find(bl);
194 if (result != clToMassWorldMap.end())
195 {return result->second;}
196 else
197 {return nullptr;}
198}
199
201{
202 clToMassWorldMap[setIn.massWorld] = setIn.massWorld;
203 clToMassWorldMap[setIn.curvilinearWorld] = setIn.massWorld;
204 clToMassWorldMap[setIn.curvilinearBridgeWorld] = setIn.massWorld;
205
206 massWorldMapTF[setIn.massWorld] = true;
207 massWorldMapTF[setIn.curvilinearWorld] = false;
208 massWorldMapTF[setIn.curvilinearBridgeWorld] = false;
209 if (setIn.endPieces) // optional 'beam line'
210 {massWorldMapTF[setIn.endPieces] = false;}
211}
212
214{
215 return BDS::MapGetWithDefault(massWorldMapTF, bl, false);
216}
217
219 G4int& index) const
220{
221 if (!BeamlineIsMassWorld(bl))
222 {
223 // update pointer to mass world beam line
224 bl = clToMassWorldMap.at(bl);
225 // for both the curvilinear and curvilinear bridge worlds, there is one
226 // extra element, so we decrement the index by one here to get the index
227 // in the mass world.
228 if (index > 0) // shouldn't happen, but we shouldn't make it < 0
229 {index--;}
230 } // else leave unchanged..
231}
232
234{
235 scorerHistogramDefs.push_back(def);
236 // we use emplace instead of map[key] = value as that method requires a
237 // default constructor which we've deleted for this class
238 scorerHistogramDefsMap.emplace(def.uniqueName, def);
239}
240
242 const G4Transform3D& placement)
243{
244 scorerMeshPlacements[meshName] = placement;
245}
246
248{
249 const auto result = scorerHistogramDefsMap.find(name);
250 if (result != scorerHistogramDefsMap.end())
251 {return &result->second;}
252 else
253 {return nullptr;}
254}
static BDSAcceleratorComponentRegistry * Instance()
Singleton accessor.
A holder class for all representations of the accelerator model created in BDSIM.
BDSBeamline * placementBeamline
Placement beam line.
std::vector< BDSScorerHistogramDef > scorerHistogramDefs
Scorer histogram definitions cached from construction here to be used in output creation.
void RegisterBeamlineSetExtra(const G4String &name, const BDSBeamlineSet &setIn)
Register a set of beam lines to be managed and cleared up at the end of the simulation.
std::map< G4String, BDSScorerHistogramDef > scorerHistogramDefsMap
Scorer histogram definitions cached from construction here to be used in output creation.
BDSBeamline * blmsBeamline
BLMs beam line.
std::set< G4LogicalVolume * > * VolumeSet(const G4String &name)
Returns pointer to a set of logical volumes. If no set by that name exits, create it.
std::map< BDSBeamline *, G4bool > massWorldMapTF
std::map< G4String, BDSBeamlineSet > extraBeamlines
Extra beamlines.
void MassWorldBeamlineAndIndex(BDSBeamline *&bl, G4int &index) const
BDSApertureInfo * Aperture(G4String name) const
const BDSScorerHistogramDef * ScorerHistogramDef(const G4String &name)
Access all scorer histogram definitions.
std::map< G4String, std::set< G4LogicalVolume * > * > volumeRegistries
All volume registries.
BDSBeamline * CorrespondingMassWorldBeamline(BDSBeamline *bl) const
std::map< G4String, BDSApertureInfo * > apertures
All apertures.
BDSAcceleratorModel()
Default constructor is private as singleton.
std::vector< BDSFieldObjects * > fields
All field objects.
G4bool BeamlineIsMassWorld(BDSBeamline *bl) const
std::set< BDSRegion * > regionStorage
Unique storage of regions.
void MapBeamlineSet(const BDSBeamlineSet &setIn)
Utility function to apply mapping.
void RegisterApertures(const std::map< G4String, BDSApertureInfo * > &aperturesIn)
Register a map of apertures with associated names.
void RegisterScorerPlacement(const G4String &meshName, const G4Transform3D &placement)
Register a copy of the scorer placement so it can be used in the output.
G4bool VolumeInSet(G4LogicalVolume *volume, const G4String &registryName)
void RegisterBeamlineSetMain(const BDSBeamlineSet &setIn)
Register the main beam line set.
G4VPhysicalVolume * worldPV
Physical volume of the mass world.
const BDSBeamlineSet & BeamlineSet(const G4String &name) const
Accessor.
std::set< G4VUserParallelWorld * > parallelWorlds
Parallel worlds not use with beam lines.
void RegisterScorerHistogramDefinition(const BDSScorerHistogramDef &def)
std::map< BDSBeamline *, BDSBeamline * > clToMassWorldMap
Mapping from any curvilinear beam line to the corresponding mass world beam line.
G4Region * Region(const G4String &name) const
Access region information. Will exit if not found.
BDSBeamline * tunnelBeamline
Tunnel segments beam line.
void RegisterRegion(BDSRegion *region)
Register a region.
Holder class for all information required to describe an aperture.
Simple struct to return a beamline plus associated beam lines.
void DeleteContents()
Destroy objects pointed to by this instance.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
General exception with possible name of object and message.
Definition: BDSException.hh:35
static BDSPhysicalVolumeInfoRegistry * Instance()
Singleton accessor.
Range cuts for a region. Help with defaults.
Definition: BDSRegion.hh:41
G4String name
Public members for simplicity.
Definition: BDSRegion.hh:57
Definition for a scorer histogram.
G4String uniqueName
Unique name of mesh/scorer -> slash required by Geant4.
V MapGetWithDefault(const std::map< K, V > &m, const K &key, const V &defaultValue)