BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSDetectorConstruction.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 BDSDETECTORCONSTRUCTION_H
20#define BDSDETECTORCONSTRUCTION_H
21
22#include "BDSExtent.hh"
23#include "BDSExtentGlobal.hh"
24
25#include "globals.hh" // geant4 types / globals
26#include "G4Transform3D.hh"
27#include "G4Version.hh"
28#include "G4VUserDetectorConstruction.hh"
29
30#include <list>
31#include <map>
32#include <set>
33#include <string>
34#include <vector>
35
36class G4LogicalVolume;
37class G4Region;
38class G4VPhysicalVolume;
39
40namespace GMAD {
41 class BLMPlacement;
42 struct Element;
43 template<typename T> class FastList;
44 class Placement;
45 class Query;
46 class SamplerPlacement;
47 class ScorerMesh;
48}
49
51class BDSBeamline;
52class BDSBeamlineSet;
54class BDSFieldObjects;
57class BDSSamplerInfo;
58
59#if G4VERSION_NUMBER > 1009
61#endif
62
74class BDSDetectorConstruction: public G4VUserDetectorConstruction
75{
76public:
77 explicit BDSDetectorConstruction(BDSComponentFactoryUser* userComponentFactoryIn = nullptr);
79
83
86 virtual G4VPhysicalVolume* Construct();
87
89 virtual void ConstructSDandField();
90
94 void BuildPhysicsBias();
95
97 inline void SetDesignParticle(const BDSParticleDefinition* defIn) {designParticle = defIn;}
98
101
103 const std::vector<BDSFieldQueryInfo*>& FieldQueries() const {return fieldQueries;}
104
113 static void PlaceBeamlineInWorld(BDSBeamline* beamline,
114 G4VPhysicalVolume* containerPV,
115 G4bool checkOverlaps = false,
116 G4bool setRegions = false,
117 G4bool registerInfo = false,
118 G4bool useCLPlacementTransform = false,
119 G4bool useIncrementalCopyNumbers = false);
120
125 static G4Transform3D CreatePlacementTransform(const GMAD::Placement& placement,
126 const BDSBeamline* beamLine,
127 G4double* S = nullptr,
128 BDSExtent* placementExtent = nullptr,
129 const G4String& objectTypeForErrorMsg = "placement");
130
131 // Create a scorermesh placement transform. Turns the scorermesh into a
133 static G4Transform3D CreatePlacementTransform(const GMAD::ScorerMesh& scorerMesh,
134 const BDSBeamline* beamLine,
135 G4double* S = nullptr);
136
139 static G4Transform3D CreatePlacementTransform(const GMAD::SamplerPlacement& samplerPlacement,
140 const BDSBeamline* beamLine,
141 G4double* S = nullptr);
142
144 static G4Transform3D CreatePlacementTransform(const GMAD::BLMPlacement& blmPlacement,
145 const BDSBeamline* beamLine,
146 G4double* S = nullptr,
147 BDSExtent* blmExtent = nullptr);
148
150 static G4Transform3D CreatePlacementTransform(const GMAD::Query& queryPlacement,
151 const BDSBeamline* beamLine,
152 G4double* S = nullptr);
153
155 static G4ThreeVector SideToLocalOffset(const GMAD::Placement& placement,
156 const BDSBeamline* beamLine,
157 const BDSExtent& placementExtent);
158
161 static BDSSamplerInfo* BuildSamplerInfo(const GMAD::Element* element);
162
165 G4bool BuildSamplerWorld() const {return nSamplers > 0;}
166
167 G4bool BuildPlacementFieldsWorld() const {return buildPlacementFieldsWorld;}
168
170 static std::vector<BDSFieldQueryInfo*> PrepareFieldQueries(const BDSBeamline* mainBeamline);
171
172private:
176
180
183
185 void InitialiseRegions();
186
188 void InitialiseApertures();
189
191 void BuildBeamlines();
192
197 const G4String& name,
198 const G4Transform3D& initialTransform = G4Transform3D(),
199 G4double initialS = 0.0,
200 G4bool beamlineIsCircular = false,
201 G4bool isPlacementBeamline= false);
202
204 void BuildTunnel();
205
208 G4VPhysicalVolume* BuildWorld();
209
211 void ComponentPlacement(G4VPhysicalVolume* worldPV);
212
215 G4bool UnsuitableFirstElement(std::list<GMAD::Element>::const_iterator element);
216
219
223
226
230
231#if G4VERSION_NUMBER > 1009
233 BDSBOptrMultiParticleChangeCrossSection* BuildCrossSectionBias(const std::list<std::string>& biasList,
234 const std::list<std::string>& defaultBias,
235 const G4String& elementName);
236
239
240
242 std::vector<BDSBOptrMultiParticleChangeCrossSection*> biasObjects;
243 std::map<G4String, BDSBOptrMultiParticleChangeCrossSection*> biasSetObjects;
244#endif
245
246#ifdef BDSDEBUG
247 bool debug = true;
248#else
249 bool debug = false;
250#endif
251
253 G4bool verbose;
256
259
261 std::vector<BDSFieldObjects*> fields;
262
263 G4bool circular;
269
270 BDSComponentFactoryUser* userComponentFactory;
271
272 G4int nSamplers;
273 G4bool buildPlacementFieldsWorld;
274
276 std::set<G4LogicalVolume*> worldContentsLogicalVolumes;
277 std::set<G4LogicalVolume*> worldVacuumLogicalVolumes;
278 G4LogicalVolume* worldLogicalVolume;
279
280 std::vector<BDSFieldQueryInfo*> fieldQueries;
281
282 // for developer checks only
283#ifdef BDSCHECKUSERLIMITS
284 void PrintUserLimitsSummary(const G4VPhysicalVolume* world) const;
285 void PrintUserLimitsPV(const G4VPhysicalVolume* aPV, G4double globalMinEK) const;
286#endif
287};
288
289#endif
290
A holder class for all representations of the accelerator model created in BDSIM.
Simple struct to return a beamline plus associated beam lines.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
Factory for user specified accelerator components.
Class that constructs a Geant4 model of an accelerator.
G4bool circular
Whether or not we're building a circular machine.
BDSBeamlineSet BuildBeamline(const GMAD::FastList< GMAD::Element > &beamLine, const G4String &name, const G4Transform3D &initialTransform=G4Transform3D(), G4double initialS=0.0, G4bool beamlineIsCircular=false, G4bool isPlacementBeamline=false)
void SetDesignParticle(const BDSParticleDefinition *defIn)
Set the design particle definition.
BDSDetectorConstruction & operator=(const BDSDetectorConstruction &)=delete
assignment and copy constructor not implemented nor used
BDSBOptrMultiParticleChangeCrossSection * BuildCrossSectionBias(const std::list< std::string > &biasList, const std::list< std::string > &defaultBias, const G4String &elementName)
Function that creates physics biasing cross section.
void InitialiseRegions()
Create and set parameters for various G4Regions.
std::vector< BDSBOptrMultiParticleChangeCrossSection * > biasObjects
List of bias objects - for memory management.
BDSExtent CalculateExtentOfScorerMesh(const GMAD::ScorerMesh &sm) const
Calculate local extent of scorer mesh in 3D.
static void PlaceBeamlineInWorld(BDSBeamline *beamline, G4VPhysicalVolume *containerPV, G4bool checkOverlaps=false, G4bool setRegions=false, G4bool registerInfo=false, G4bool useCLPlacementTransform=false, G4bool useIncrementalCopyNumbers=false)
BDSAcceleratorModel * acceleratorModel
Accelerator model pointer.
static G4ThreeVector SideToLocalOffset(const GMAD::Placement &placement, const BDSBeamline *beamLine, const BDSExtent &placementExtent)
Attach component with extent2 to component with extent1 with placement.
G4bool verbose
Variable copied from global constants.
static std::vector< BDSFieldQueryInfo * > PrepareFieldQueries(const BDSBeamline *mainBeamline)
Prepare field queries from parser information.
G4VPhysicalVolume * BuildWorld()
void InitialiseApertures()
Create all aperture definitions from parser and store in BDSAcceleratorModel.
void CountPlacementFields()
Count number of fields required for placements.
G4bool canSampleAngledFaces
Whether the integrator set permits sampling elements with angled faces.
BDSExtent worldExtent
Record of the world extent.
std::set< G4LogicalVolume * > worldContentsLogicalVolumes
Cache of possibly loaded logical volumes from a world geometry file - used for biasing.
static BDSSamplerInfo * BuildSamplerInfo(const GMAD::Element *element)
BDSExtent CalculateExtentOfSamplerPlacement(const GMAD::SamplerPlacement &sp) const
Calculate local extent of custom user sampler.
void ConstructScoringMeshes()
Construct scoring meshes.
BDSExtentGlobal CalculateExtentOfScorerMeshes(const BDSBeamline *bl) const
G4bool checkOverlaps
Variable copied from global constants.
G4int nSamplers
Count of number of samplers to be built.
void ComponentPlacement(G4VPhysicalVolume *worldPV)
Place beam line, tunnel beam line, end pieces and placements in world.
const std::vector< BDSFieldQueryInfo * > & FieldQueries() const
Access vector of query objects.
void BuildBeamlines()
Build the main beam line and then any other required beam lines.
const BDSParticleDefinition * designParticle
Particle definition all components are built w.r.t. Includes rigidity etc.
void BuildTunnel()
Build the tunnel around the already constructed flat beam line.
BDSExtent WorldExtent() const
Public access to the world extent.
virtual void ConstructSDandField()
Construct sensitive detectors and fields.
std::vector< BDSFieldObjects * > fields
All fields.
virtual G4VPhysicalVolume * Construct()
static G4Transform3D CreatePlacementTransform(const GMAD::Placement &placement, const BDSBeamline *beamLine, G4double *S=nullptr, BDSExtent *placementExtent=nullptr, const G4String &objectTypeForErrorMsg="placement")
G4bool UnsuitableFirstElement(std::list< GMAD::Element >::const_iterator element)
BDSExtentGlobal CalculateExtentOfSamplerPlacements(const BDSBeamline *beamLine) const
Holder for +- extents in 3 dimensions with a rotation and translation.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
A holder for all the Geant4 field related objects.
Holder class for all information required for a field query.
Wrapper for particle definition.
All info required to build a sampler but not place it.
blm for parser.
Definition: blmplacement.h:39
Placement class for parser.
Definition: placement.h:41
Query structure class for parser.
Definition: query.h:37
Sampler placement class for parser.
ScorerMesh class for parser.
Definition: scorermesh.h:38
Parser namespace for GMAD language. Combination of Geant4 and MAD.
Element class.
Definition: element.h:43