BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSMultilayerScreen.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 BDSMULTILAYERSCREEN_H
20#define BDSMULTILAYERSCREEN_H
21
22#include <map>
23#include <vector>
24
25#include "globals.hh"
26#include "G4RotationMatrix.hh"
27#include "G4TwoVector.hh"
28#include "G4ThreeVector.hh"
29
30class G4LogicalVolume;
31class G4PVPlacement;
32class G4VSolid;
33
34class BDSColourWheel;
35class BDSScreenLayer;
36
44{
45public:
46 BDSMultilayerScreen(const G4TwoVector& xysizeIn,
47 const G4String& nameIn);
48 virtual ~BDSMultilayerScreen();
49
51 void AddScreenLayer(G4double thickness,
52 G4String material,
53 G4String name,
54 G4int isSampler = 0,
55 G4double grooveWidth = 0,
56 G4double grooveSpatialFrequency = 0);
58 void AddScreenLayer(BDSScreenLayer* layer, G4int isSampler=0);
59
61 inline const G4ThreeVector& GetSize() const {return size;}
63 inline BDSScreenLayer* ScreenLayer(G4int layer) {return screenLayers[layer];}
65 BDSScreenLayer* ScreenLayer(const G4String& layerName);
67 BDSScreenLayer* LastLayer() const {return screenLayers.back();}
69 inline G4int NLayers() const {return (G4int)screenLayers.size();}
70
72 inline void SetPhys(G4PVPlacement* physIn) {phys = physIn;}
73
75 void Build();
76
78 virtual void Place(G4RotationMatrix* rot,
79 const G4ThreeVector& pos,
80 G4LogicalVolume* motherVol);
81
83 void ReflectiveSurface(G4int layer1, G4int layer2);
85 void RoughSurface(G4int layer1, G4int layer2);
86
87private:
94 void ComputeDimensions();
95
97 void BuildMotherVolume();
98
100 void PlaceLayers();
101
102 G4TwoVector xysize;
103 G4String name;
104 G4ThreeVector size;
105
107 G4LogicalVolume* log;
108 G4PVPlacement* phys;
109 G4VSolid* solid;
111 std::vector<BDSScreenLayer*> screenLayers;
112 std::map<G4String, BDSScreenLayer*> screenLayerNames;
113 std::vector<G4double> screenLayerZPos;
114
118};
119
120#endif
Three colours that are supplied sequentially.
An accelerator component for diagnostics screens e.g. OTR. Screen inside beam pipe.
void AddScreenLayer(G4double thickness, G4String material, G4String name, G4int isSampler=0, G4double grooveWidth=0, G4double grooveSpatialFrequency=0)
Construct and add a screen layer.
const G4ThreeVector & GetSize() const
Full size of screen in x,y,z.
std::vector< BDSScreenLayer * > screenLayers
Main storage of all layers.
void RoughSurface(G4int layer1, G4int layer2)
Construct a rough surface between two layers.
BDSMultilayerScreen & operator=(const BDSMultilayerScreen &)=delete
Private default constructor to force use of provided one.
BDSMultilayerScreen(const BDSMultilayerScreen &)=delete
Private default constructor to force use of provided one.
void SetPhys(G4PVPlacement *physIn)
Set a physical placement to member variable.
BDSMultilayerScreen()=delete
Private default constructor to force use of provided one.
G4PVPlacement * phys
Geometrical objects:
G4ThreeVector size
Extent.
G4VSolid * solid
Geometrical objects:
void BuildMotherVolume()
Build a container volume.
std::map< G4String, BDSScreenLayer * > screenLayerNames
Map of names for looking up.
std::vector< G4double > screenLayerZPos
Cache of calculated z locations for layers.
virtual void Place(G4RotationMatrix *rot, const G4ThreeVector &pos, G4LogicalVolume *motherVol)
Place the member logical volume 'log' with a given transform in a given mother volume.
BDSScreenLayer * ScreenLayer(G4int layer)
Get a particular screen layer by index.
G4int NLayers() const
Size of the screen.
BDSColourWheel * colourWheel
void Build()
Construct container, compute dimensions then place layers.
G4LogicalVolume * log
Geometrical objects:
G4TwoVector xysize
X,Y size of multilayer screen.
void ReflectiveSurface(G4int layer1, G4int layer2)
Construct a reflective surface between two layers.
G4String name
Name of multilayer screen.
BDSScreenLayer * LastLayer() const
Get the last layer of the screen.
void PlaceLayers()
Place each layer in order in the container volume.
An individual screen layer for a multilayer screen.