BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSSamplerRegistry.hh
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#ifndef BDSSAMPLERREGISTRY_H
20#define BDSSAMPLERREGISTRY_H
21#include "BDSDebug.hh"
22#include "BDSException.hh"
23#include "BDSSamplerPlacementRecord.hh"
24#include "BDSSamplerType.hh"
25
26#include "globals.hh" // geant4 types / globals
27#include "G4Transform3D.hh"
28
29#include <map>
30#include <set>
31#include <utility>
32#include <vector>
33
35class BDSSampler;
36
59{
60private:
62 typedef std::vector<BDSSamplerPlacementRecord> InfoVector;
63
66
67public:
70
72
86 G4int RegisterSampler(const G4String& name,
87 BDSSampler* sampler,
88 const G4Transform3D& transform = G4Transform3D(),
89 G4double S = -1000,
90 const BDSBeamlineElement* element = nullptr,
91 BDSSamplerType type = BDSSamplerType::plane,
92 G4double radius = 0);
93
95
97 typedef InfoVector::iterator iterator;
98 typedef InfoVector::const_iterator const_iterator;
99 iterator begin() {return infos.begin();}
100 iterator end() {return infos.end();}
101 const_iterator begin() const {return infos.begin();}
102 const_iterator end() const {return infos.end();}
103 G4bool empty() const {return infos.empty();}
105
107 inline G4String GetName(G4int index) const;
108 inline G4String GetNameUnique(G4int index) const;
109 inline BDSSampler* GetSampler(G4int index) const;
110 inline G4Transform3D GetTransform(G4int index) const;
111 inline G4Transform3D GetTransformInverse(G4int index) const;
112 inline G4double GetSPosition(G4int index) const;
113 inline const BDSSamplerPlacementRecord& GetInfo(G4int index) const {return infos.at(index);}
114 inline G4int GetBeamlineIndex(G4int index) const;
116
118 std::vector<G4String> GetNames() const;
119
121 std::vector<G4String> GetUniqueNames() const;
122 std::vector<G4String> GetUniqueNamesPlane() const;
123 std::vector<G4String> GetUniqueNamesCylinder() const;
124 std::vector<G4String> GetUniqueNamesSphere() const;
125
127 std::map<std::string, double> GetUniqueNameToRadiusCylinder() const;
128 std::map<std::string, double> GetUniqueNameToRadiusSphere() const;
130
131 inline std::vector<G4int> GetSamplerIDsPlane() const {return samplerIDsPerType.at(BDSSamplerType::plane);}
132 inline std::vector<G4int> GetSamplerIDsCylinder() const {return samplerIDsPerType.at(BDSSamplerType::cylinder);}
133 inline std::vector<G4int> GetSamplerIDsSphere() const {return samplerIDsPerType.at(BDSSamplerType::sphere);}
134
136 std::vector<std::pair<G4String, G4double> > GetUniquePlaneNamesAndSPosition() const;
137
139 inline G4int NumberOfExistingSamplers() const;
140 inline size_t size() const;
141
142 inline G4bool SafeIndex(G4int index) const {return index < (G4int)infos.size();}
143
144private:
147
150
154
158 std::map<G4String, G4int> existingNames;
159
161 std::set<BDSSampler*> samplerObjects;
162
164 std::map<BDSSamplerType, BDSSamplerType> samplerTypeToCategory;
165 std::map<BDSSamplerType, std::vector<G4int>> samplerIDsPerType;
166};
167
168inline G4String BDSSamplerRegistry::GetName(G4int index) const
169{return SafeIndex(index) ? infos.at(index).Name() : throw BDSException(__METHOD_NAME__, "invalid index");}
170
171inline G4String BDSSamplerRegistry::GetNameUnique(G4int index) const
172{return SafeIndex(index) ? infos.at(index).UniqueName() : throw BDSException(__METHOD_NAME__, "invalid index");}
173
175{return SafeIndex(index) ? infos.at(index).Sampler() : throw BDSException(__METHOD_NAME__, "invalid index");}
176
177inline G4Transform3D BDSSamplerRegistry::GetTransform(G4int index) const
178{return SafeIndex(index) ? infos.at(index).Transform() : throw BDSException(__METHOD_NAME__, "invalid index");}
179
180inline G4Transform3D BDSSamplerRegistry::GetTransformInverse(G4int index) const
181{return SafeIndex(index) ? infos.at(index).TransformInverse() : throw BDSException(__METHOD_NAME__, "invalid index");}
182
183inline G4double BDSSamplerRegistry::GetSPosition(G4int index) const
184{return SafeIndex(index) ? infos.at(index).SPosition() : throw BDSException(__METHOD_NAME__, "invalid index");}
185
187{return numberOfEntries;}
188
189inline size_t BDSSamplerRegistry::size() const
190{return infos.size();}
191
192inline G4int BDSSamplerRegistry::GetBeamlineIndex(G4int index) const
193{return SafeIndex(index) ? infos.at(index).BeamlineIndex() : throw BDSException(__METHOD_NAME__, "invalid index");}
194
195#endif
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
General exception with possible name of object and message.
Definition: BDSException.hh:35
Information about a registered sampler.
Associated information for the placement of a sampler.
std::vector< std::pair< G4String, G4double > > GetUniquePlaneNamesAndSPosition() const
Access all the unique names and their corresponding s position at once.
std::map< BDSSamplerType, BDSSamplerType > samplerTypeToCategory
Map to reduce 'forward' sampler types to simple sampler types for keeping a record of IDs.
G4Transform3D GetTransform(G4int index) const
Accessor.
InfoVector::iterator iterator
Iterator mechanics.
InfoVector::const_iterator const_iterator
Iterator mechanics.
std::map< std::string, double > GetUniqueNameToRadiusSphere() const
For output in standard C++ types. Also in m for units.
iterator begin()
Iterator mechanics.
iterator end()
Iterator mechanics.
BDSSampler * GetSampler(G4int index) const
Accessor.
G4int GetBeamlineIndex(G4int index) const
Accessor.
std::vector< BDSSamplerPlacementRecord > InfoVector
Typedefs up first so we can declare public iterators.
std::map< G4String, G4int > existingNames
static BDSSamplerRegistry * Instance()
Accessor for registry.
const_iterator begin() const
Iterator mechanics.
G4bool empty() const
Iterator mechanics.
G4double GetSPosition(G4int index) const
Accessor.
std::vector< G4String > GetUniqueNames() const
Access all the unique names at once.
const BDSSamplerPlacementRecord & GetInfo(G4int index) const
Accessor.
InfoVector infos
Storage of registered information.
static BDSSamplerRegistry * instance
Singleton instance.
G4int NumberOfExistingSamplers() const
Get number of registered samplers.
std::vector< G4String > GetNames() const
Access all names at once.
G4String GetNameUnique(G4int index) const
Accessor.
BDSSamplerRegistry()
Private constructor to enforce singleton pattern.
std::set< BDSSampler * > samplerObjects
Cache of unique sampler objects for memory management.
const_iterator end() const
Iterator mechanics.
G4int RegisterSampler(const G4String &name, BDSSampler *sampler, const G4Transform3D &transform=G4Transform3D(), G4double S=-1000, const BDSBeamlineElement *element=nullptr, BDSSamplerType type=BDSSamplerType::plane, G4double radius=0)
G4String GetName(G4int index) const
Accessor.
G4Transform3D GetTransformInverse(G4int index) const
Accessor.
std::map< std::string, double > GetUniqueNameToRadiusCylinder() const
For output in standard C++ types. Also in m for units.
Base class and registry of sampler instances.
Definition: BDSSampler.hh:36