BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSLinkDetectorConstruction.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 BDSLINKDETECTORCONSTRUCTION_H
20#define BDSLINKDETECTORCONSTRUCTION_H
21#include "BDSBeamline.hh"
22#include "BDSCollimatorJaw.hh"
23#include "BDSExtent.hh"
24
25#include "G4ThreeVector.hh"
26#include "G4Types.hh"
27#include "G4Version.hh"
28#include "G4VUserDetectorConstruction.hh"
29
30#include <string>
31
32class BDSBeamline;
35class BDSLinkRegistry;
37class G4Box;
38class G4ChannelingOptrMultiParticleChangeCrossSection;
39class G4VPhysicalVolume;
40
47class BDSLinkDetectorConstruction: public G4VUserDetectorConstruction
48{
49public:
52
54 virtual G4VPhysicalVolume* Construct();
55
57 G4int AddLinkCollimatorJaw(const std::string& collimatorName,
58 const std::string& materialName,
59 G4double length,
60 G4double halfApertureLeft,
61 G4double halfApertureRight,
62 G4double rotation,
63 G4double xOffset,
64 G4double yOffset,
65 G4bool buildLeftJaw = true,
66 G4bool buildRightJaw = true,
67 G4bool isACrystal = false,
68 G4double crystalAngle = 0,
69 G4bool sampleIn = false);
70
72 inline void SetDesignParticle(const BDSParticleDefinition* defIn) {designParticle = defIn;}
73 inline void SetPrimaryGeneratorAction(BDSLinkPrimaryGeneratorAction* pgIn) {primaryGeneratorAction = pgIn;}
74
75
77 inline BDSExtent WorldExtent() const {return worldExtent;}
78 inline BDSLinkRegistry* LinkRegistry() const {return linkRegistry;}
80
81 void BuildPhysicsBias();
82
83 inline const std::map<std::string, G4int>& NameToElementIndex() const {return nameToElementIndex;}
84 inline const std::map<int, int>& LinkIDToBeamlineIndex() const {return linkIDToBeamlineIndex;}
85 inline G4int NumberOfElements() const {return linkBeamline ? (G4int)linkBeamline->size() : 0;}
86 inline void SetSamplerWorldID(G4int samplerWorldIDIn) {samplerWorldID = samplerWorldIDIn;}
87 inline const BDSBeamline* LinkBeamline() const {return linkBeamline;}
88
89 private:
92 void UpdateWorldSolid();
93
95 G4int PlaceOneComponent(const BDSBeamlineElement* element, const G4String& originalName);
96
97 G4Box* worldSolid;
98 G4VPhysicalVolume* worldPV;
99 BDSExtent worldExtent;
100 BDSBeamline* linkBeamline;
101 BDSLinkRegistry* linkRegistry;
102 BDSLinkPrimaryGeneratorAction* primaryGeneratorAction;
103
106
107#if G4VERSION_NUMBER > 1039
108 G4ChannelingOptrMultiParticleChangeCrossSection* crystalBiasing;
109#endif
110
113
114 std::map<std::string, G4int> nameToElementIndex;
115 std::map<G4int, G4int> linkIDToBeamlineIndex;
116};
117
118#endif
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
BeamlineVector::size_type size() const
Get the number of elements.
Definition: BDSBeamline.hh:148
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
Construction of the geometry in the case of a link model.
BDSExtent WorldExtent() const
Accessor.
G4int samplerWorldID
Cache of the index to which parallel world the sampler one is.
G4int PlaceOneComponent(const BDSBeamlineElement *element, const G4String &originalName)
Place a beam line element in the world.
const BDSParticleDefinition * designParticle
Particle definition all components are built w.r.t. Includes rigidity etc.
BDSLinkRegistry * LinkRegistry() const
Accessor.
BDSLinkDetectorConstruction()
Default constructor.
std::map< std::string, G4int > nameToElementIndex
Build up a copy here too.
std::map< G4int, G4int > linkIDToBeamlineIndex
Special linkID to linkBeamline index.
G4int AddLinkCollimatorJaw(const std::string &collimatorName, const std::string &materialName, G4double length, G4double halfApertureLeft, G4double halfApertureRight, G4double rotation, G4double xOffset, G4double yOffset, G4bool buildLeftJaw=true, G4bool buildRightJaw=true, G4bool isACrystal=false, G4double crystalAngle=0, G4bool sampleIn=false)
Interface to append a collimator of jaw style to the linking.
void SetDesignParticle(const BDSParticleDefinition *defIn)
Set the design particle definition.
Generates primary particle vertices using BDSBunch.
Registry / map of components for tracker linkage.
Wrapper for particle definition.