BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSSDManager.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 BDSSDMANAGER_H
20#define BDSSDMANAGER_H
21
22#include "BDSSDSamplerLink.hh"
23#include "BDSSDType.hh"
24
25#include "G4String.hh"
26#include "G4Types.hh"
27#include "G4Version.hh"
28
29#include <map>
30#include <set>
31#include <vector>
32
34class BDSSDCollimator;
37class BDSLinkRegistry;
40class BDSSDSampler;
43class BDSSDTerminator;
44class BDSSDThinThing;
45class BDSSDVolumeExit;
46
47class G4VSDFilter;
48
49#if G4VERSION_NUMBER < 1030
50// In this case we use only the energy counter SD and return it
51// as a base class pointer. Include header so casting works.
52#include "BDSSDEnergyDeposition.hh"
53#include "BDSSDEnergyDepositionGlobal.hh"
54#endif
55
68{
69public:
70 static BDSSDManager* Instance();
71
72 BDSSDManager(const BDSSDManager&) = delete;
73 BDSSDManager& operator=(const BDSSDManager&) = delete;
75
82 G4VSensitiveDetector* SensitiveDetector(const BDSSDType sdType,
83 G4bool applyOptions = false) const;
84
86 inline BDSSDSampler* SamplerPlane() const {return samplerPlane;}
87
90
93
95 inline BDSSDSamplerLink* SamplerLink() const {return samplerLink;}
96
99 inline BDSSDTerminator* Terminator() const {return terminator;}
100
103
106
109
112
115
118
120 inline BDSSDVolumeExit* WorldExit() const {return worldExit;}
121
124
125#if G4VERSION_NUMBER > 1029
127 inline G4VSensitiveDetector* WorldComplete() const {return worldCompleteSD;}
128 inline G4VSensitiveDetector* ApertureComplete() const {return apertureCompleteSD;}
129#else
132 inline G4VSensitiveDetector* WorldComplete() const {return energyDepositionWorld;}
133 inline G4VSensitiveDetector* ApertureComplete() const {return energyDeposition;}
134#endif
135
137 inline BDSSDCollimator* Collimator() const {return collimatorSD;}
138
140 inline BDSMultiSensitiveDetectorOrdered* CollimatorComplete() const {return collimatorCompleteSD;}
141
144 inline BDSSDThinThing* ThinThing() const {return thinThingSD;}
145
147 inline BDSMultiSensitiveDetectorOrdered* WireComplete() const {return wireCompleteSD;}
148
151 void RegisterPrimitiveScorerName(const G4String& nameIn, G4double unit = 1.0);
152
154 void RegisterPrimitiveScorerNames(const std::vector<G4String>& namesIn,
155 const std::vector<G4double>* units = nullptr);
156
158 inline const std::vector<G4String>& PrimitiveScorerNamesComplete() const {return primitiveScorerNamesComplete;}
159
161 inline const std::vector<G4String>& PrimitiveScorerNames() const {return primitiveScorerNames;}
162
164 inline const std::map<G4String, G4double>& PrimitiveScorerUnits() const {return primitiveScorerNameToUnit;}
165
167 void SetLinkRegistry(BDSLinkRegistry* registry);
168 inline void SetLinkMinimumEK(G4double minimumEKIn) {samplerLink->SetMinimumEK(minimumEKIn);}
169 inline void SetLinkProtonsAndIonsOnly(G4bool protonsAndIonsOnlyIn) {samplerLink->SetProtonsAndIonsOnly(protonsAndIonsOnlyIn);}
170
174 void ConstructSamplerSDsForParticleSets(const std::map<int, std::set<int>>& samplerFilterIDtoPDGSetIn);
175
177 inline const std::vector<G4String>& ExtraSamplerWithFilterNamesComplete() const {return extraSamplerWithFilterNamesComplete;}
178
180 inline const std::vector<G4String>& ExtraSamplerCylinderWithFilterNamesComplete() const {return extraSamplerCylinderWithFilterNamesComplete;}
181
183 inline const std::vector<G4String>& ExtraSamplerSphereWithFilterNamesComplete() const {return extraSamplerSphereWithFilterNamesComplete;}
184
186 BDSSDSampler* SamplerPlaneWithFilter(G4int ID) const;
187
190
193
194private:
196 BDSSDManager();
197
198 static BDSSDManager* instance;
199
214#if G4VERSION_NUMBER > 1029
215 G4VSensitiveDetector* worldCompleteSD;
216 G4VSensitiveDetector* apertureCompleteSD;
217#endif
219 BDSSDCollimator* collimatorSD;
220 BDSMultiSensitiveDetectorOrdered* collimatorCompleteSD;
221 BDSSDThinThing* thinThingSD;
222 BDSMultiSensitiveDetectorOrdered* wireCompleteSD;
223
225 std::map<G4String, G4VSDFilter*> filters;
226
228 std::vector<G4String> primitiveScorerNamesComplete;
229
231 std::vector<G4String> primitiveScorerNames;
232
234 std::map<G4String, G4double> primitiveScorerNameToUnit;
235
236 std::map<G4int, BDSSDSampler*> extraSamplersWithFilters;
237 std::map<G4int, BDSSDSamplerCylinder*> extraSamplerCylindersWithFilters;
238 std::map<G4int, BDSSDSamplerSphere*> extraSamplerSpheresWithFilters;
239 std::map<G4int, BDSSDFilterPDGIDSet*> extraSamplerFilters;
240 std::vector<G4String> extraSamplerWithFilterNamesComplete;
241 std::vector<G4String> extraSamplerCylinderWithFilterNamesComplete;
242 std::vector<G4String> extraSamplerSphereWithFilterNamesComplete;
243
260};
261
262#endif
Registry / map of components for tracker linkage.
Modified G4MultiSensitiveDetector that retains order and passes hits in sequence.
Generates BDSHitsEnergyDepositions from step information - uses curvilinear coords.
The sensitive detector class that provides sensitivity to collimators instances.
Generates BDSHitsEnergyDepositionGlobal from step information.
Generates BDSHitsEnergyDepositions from step information - uses curvilinear coords.
Filter for a set of PDG IDs (ints).
A singleton class that holds all required sensitive detector class instances.
Definition: BDSSDManager.hh:68
BDSSDEnergyDeposition * EnergyDeposition() const
SD for general energy deposition.
G4double collimatorHitsMinimumKE
Cache of global constant option.
const std::map< G4String, G4double > & PrimitiveScorerUnits() const
Access the map of units for primitive scorers.
BDSSDManager()
Private default constructor for singleton.
Definition: BDSSDManager.cc:76
BDSSDSamplerSphere * samplerSphere
SD instance.
G4bool storeCollimatorHitsAll
Cache of global constant option.
std::map< G4String, G4VSDFilter * > filters
Map of all filters used. This class owns a single instance of each.
BDSSDVolumeExit * WorldExit() const
SD for world exit hits.
BDSSDTerminator * terminator
SD instance.
BDSSDSamplerCylinder * SamplerCylinderWithFilter(G4int ID) const
Access the relevant SD for a given particle filter set ID. It will return nullptr if the ID is invali...
G4bool generateCollimatorHits
Cache of global constant option.
BDSSDVolumeExit * worldExit
SD instance.
G4bool generateELossVacuumHits
Cache of global constant option.
BDSSDEnergyDepositionGlobal * energyDepositionWorldContents
SD instance.
void RegisterPrimitiveScorerName(const G4String &nameIn, G4double unit=1.0)
BDSSDApertureImpacts * apertureImpacts
SD instance.
void ConstructSamplerSDsForParticleSets(const std::map< int, std::set< int > > &samplerFilterIDtoPDGSetIn)
BDSMultiSensitiveDetectorOrdered * CollimatorComplete() const
SD for collimator impacts + energy deposition at the same time in order.
BDSSDEnergyDeposition * EnergyDepositionTunnel() const
SD for tunnel energy counter.
BDSSDTerminator * Terminator() const
Definition: BDSSDManager.hh:99
const std::vector< G4String > & ExtraSamplerCylinderWithFilterNamesComplete() const
Access a vector of names of extra samplers so we can identify the hits collections.
const std::vector< G4String > & ExtraSamplerSphereWithFilterNamesComplete() const
Access a vector of names of extra samplers so we can identify the hits collections.
G4bool generateELossWorldContents
Cache of global constant option.
G4bool generateELossTunnelHits
Cache of global constant option.
const std::vector< G4String > & PrimitiveScorerNames() const
Access a vector of the just primitive scorer part of the names.
BDSSDEnergyDeposition * EnergyDepositionFull() const
SD for general energy deposition but always include extra half of information.
BDSSDSamplerLink * samplerLink
SD instance.
G4bool storeELossWorld
Cache of global constant option.
BDSSDSamplerSphere * SamplerSphereWithFilter(G4int ID) const
Access the relevant SD for a given particle filter set ID. It will return nullptr if the ID is invali...
BDSSDEnergyDeposition * energyDeposition
SD instance.
void RegisterPrimitiveScorerNames(const std::vector< G4String > &namesIn, const std::vector< G4double > *units=nullptr)
Loop over a vector and apply above single name function.
G4VSensitiveDetector * worldCompleteSD
SD instance.
G4bool storeApertureImpactsIons
Cache of global constant option.
G4VSensitiveDetector * apertureCompleteSD
SD instance.
BDSSDSamplerLink * SamplerLink() const
SD for link samplers.
Definition: BDSSDManager.hh:95
BDSSDCollimator * Collimator() const
SD for collimator impact locations.
BDSSDThinThing * ThinThing() const
G4bool storeApertureImpactsAll
Cache of global constant option.
BDSSDApertureImpacts * ApertureImpacts() const
SD for aperture impact hits.
BDSSDEnergyDeposition * energyDepositionTunnel
SD instance.
BDSSDSampler * samplerPlane
SD instance.
G4bool generateELossHits
Cache of global constant option.
G4VSensitiveDetector * SensitiveDetector(const BDSSDType sdType, G4bool applyOptions=false) const
BDSSDSamplerCylinder * samplerCylinder
SD instance.
BDSSDSamplerCylinder * SamplerCylinder() const
SD for samplers (cylinder type).
Definition: BDSSDManager.hh:89
const std::vector< G4String > & PrimitiveScorerNamesComplete() const
Access a vector the full primitive scorer names as registered.
BDSSDManager(const BDSSDManager &)=delete
Singleton accessor.
std::map< G4String, G4double > primitiveScorerNameToUnit
Map of primitive scorer names to units.
BDSSDSampler * SamplerPlane() const
SD for samplers (plane type). See also SamplerPlaneWithFilter below.
Definition: BDSSDManager.hh:86
BDSSDSamplerSphere * SamplerSphere() const
SD for samplers (sphere type).
Definition: BDSSDManager.hh:92
BDSSDSampler * SamplerPlaneWithFilter(G4int ID) const
Access the relevant SD for a given particle filter set ID. It will return nullptr if the ID is invali...
G4bool generateApertureImpacts
Cache of global constant option.
BDSMultiSensitiveDetectorOrdered * WireComplete() const
SD for wire scanner wires that is a composite of thin things + energy deposition full.
G4bool storeCollimatorHitsIons
Cache of global constant option.
G4VSensitiveDetector * WorldComplete() const
SD for multiple SDs for world - energy loss and exit.
BDSSDEnergyDepositionGlobal * EnergyDepositionWorldContents() const
SD for energy deposition in things that were already placed in the externally provided world.
BDSSDEnergyDeposition * energyDepositionFull
SD instance.
G4double apertureImpactsMinimumKE
Cache of global constant option.
BDSSDEnergyDepositionGlobal * EnergyDepositionWorld() const
SD for energy deposition in the world volume.
std::vector< G4String > primitiveScorerNamesComplete
Vector of complete (including mesh name) primitive scorer names.
G4bool storeELossExtras
Cache of global constant option.
std::vector< G4String > primitiveScorerNames
Just the primitive scorer part of the name.
void SetLinkRegistry(BDSLinkRegistry *registry)
If samplerLink member exists, set the registry to look up links for that SD.
BDSSDEnergyDeposition * EnergyDepositionVacuum() const
SD for energy deposition in vacuum volumes.
BDSSDEnergyDepositionGlobal * energyDepositionWorld
SD instance.
BDSSDEnergyDeposition * energyDepositionVacuum
SD instance.
const std::vector< G4String > & ExtraSamplerWithFilterNamesComplete() const
Access a vector of names of extra samplers so we can identify the hits collections.
The sensitive detector class that provides sensitivity to BDSSamplerCylinder instances.
The sensitive detector class that provides sensitivity to BDSSamplerSphere instances.
The sensitive detector class that provides sensitivity to BDSSampler instances.
Definition: BDSSDSampler.hh:44
Sensitivity that measures primary particle turns for terminator.
The sensitive detector class that provides sensitivity to record thin thing hits.
Generates BDSHitEnergyDepositionGlobals if a particle is leaving a volume.
Improve type-safety of native enum data type in C++.