BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSScreenLayer.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 BDSSCREENLAYER_H
20#define BDSSCREENLAYER_H
21
22#include "globals.hh"
23#include "G4Colour.hh"
24#include "G4ThreeVector.hh"
25
26class G4LogicalVolume;
27class G4PVPlacement;
28class G4VSolid;
29
37{
38public:
39 BDSScreenLayer(G4ThreeVector size, // full widths in each dimension
40 G4String name,
41 G4String material,
42 G4double grooveWidth,
43 G4double grooveSpatialFrequency);
44 virtual ~BDSScreenLayer(){;}
46 inline G4LogicalVolume* GetLog() const {return log;}
47 inline G4String GetName() const {return name;}
48 inline G4ThreeVector GetSize() const {return size;}
49 inline G4PVPlacement* GetPhys() const {return phys;}
51 void SetPhys(G4PVPlacement* phys);
52 void SetColour(G4Colour col);
53 /*
54 void backInternalMirror();
55 void frontInternalMirror();
56 */
58 void AssignSampler();
59
60 //Get the AssignSampler ID.
61 G4int GetSamplerID() const {return samplerID;}
62
63protected:
66 G4ThreeVector size;
67 G4String name;
68 G4LogicalVolume* log = nullptr;
69 G4PVPlacement* phys = nullptr;
70 G4VSolid* solid = nullptr;
71
72private:
73 /*
74 class InternalMirror
75 {
76 public:
77 InternalMirror(G4int varside,
78 G4ThreeVector size,
79 G4String material,
80 G4LogicalVolume* motherLog,
81 G4PVPlacement* motherPhys);
82 ~InternalMirror();
83 void geom();
84 void compute();
85 void place();
86 void optical();
87 enum sides{BACK, FRONT};
88
89 private:
90 InternalMirror();
91 G4int side;
92 G4VSolid* solid;
93 G4LogicalVolume* log;
94 G4PVPlacement* phys;
95 G4ThreeVector motherSize;
96 G4String motherMaterial;
97 G4LogicalVolume* motherLog;
98 G4PVPlacement* motherPhys;
99 G4double thickness;
100 G4double pos;
101 };
102
103 InternalMirror* internalMirror = nullptr;
104 */
105 virtual void Build();
106 void BuildGroove();
107 virtual void BuildScreen();
108 void SetVisAttributes();
109 void CutGroove(G4double xPos);
110 void CutGrooves();
111 G4String material;
112 G4String logName;
113 G4String solidName;
114 // Geometrical objects:
115 G4LogicalVolume* grooveLog = nullptr;
116 G4VSolid* grooveSolid = nullptr;
117 G4double grooveWidth = 0.0;
118 G4double grooveSpatialFrequency = 0.0;
120 G4int nGrooves = 0;
121 G4Colour colour;
122 G4int samplerID = 0;
123};
124
125#endif
An individual screen layer for a multilayer screen.
G4LogicalVolume * GetLog() const
Accessor.
void AssignSampler()
Make this plane a sampling plane.
G4int nGrooves
Counter for the number of grooves etched into the screen.
G4String GetName() const
Accessor.
BDSScreenLayer()
Used in modules/AWAKE/include/BDSMultiFacetLayer.hh.
G4ThreeVector GetSize() const
Accessor.
G4PVPlacement * GetPhys() const
Accessor.