BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSScreenFrameRectangular.cc
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#include "BDSGlobalConstants.hh"
20#include "BDSMaterials.hh"
21#include "BDSScreenFrameRectangular.hh"
22
23#include "globals.hh"
24#include "G4Box.hh"
25#include "G4LogicalVolume.hh"
26#include "G4PVPlacement.hh"
27#include "G4SubtractionSolid.hh"
28#include "G4RotationMatrix.hh"
29#include "G4ThreeVector.hh"
30#include "G4TwoVector.hh"
31
32
34 G4ThreeVector sizeIn,
35 G4TwoVector windowSizeIn,
36 G4TwoVector windowOffsetIn,
37 G4Material* materialIn):
38 BDSScreenFrame(nameIn, sizeIn, windowSizeIn, windowOffsetIn, materialIn),
39 zeroRot(new G4RotationMatrix(0,0,0))
40{
41 Build();
42}
43
44BDSScreenFrameRectangular::~BDSScreenFrameRectangular()
45{
46 delete zeroRot;
47}
48
50{
51 G4double lenSaf = BDSGlobalConstants::Instance()->LengthSafety();
52 //G4double tinyLenSaf = 1e-3*CLHEP::nm;
53
54 G4Box* frameBox = new G4Box("frameBox",
55 size.x()/2.0,
56 size.y()/2.0,
57 size.z()/2.0);
58
59
60 G4Box* windowBox = new G4Box("windowBox",
61 windowSize.x()/2.0 + lenSaf,
62 windowSize.y()/2.0 + lenSaf,
63 size.z()/2.0 + lenSaf); // must be thicker than frame.
64
65 cavityPos = G4ThreeVector(windowOffset.x(), windowOffset.y(), 0);
66 G4SubtractionSolid* frameSS = new G4SubtractionSolid((G4String)"frameSS", frameBox, windowBox, zeroRot, cavityPos);
67 logVol = new G4LogicalVolume(frameSS, material, name+"_lv");
68 SetVisAtts();
69}
70
71void BDSScreenFrameRectangular::Place(G4RotationMatrix* rot,
72 G4ThreeVector pos,
73 G4LogicalVolume* motherVol)
74{
75 new G4PVPlacement(rot,
76 pos,
77 logVol,
78 "screenFrame_pv",
79 motherVol,
80 false,
81 0,
83
84}
85
87{
89}
static BDSGlobalConstants * Instance()
Access method.
virtual void SetVisAtts()
Set the visual attributes to member logVol.
void Place(G4RotationMatrix *rot, G4ThreeVector pos, G4LogicalVolume *motherVol)
Place member logVol (from base class) into argument motherVol.
virtual void Build()
Construct the geometry.
BDSScreenFrameRectangular()=delete
Remove default constructor to force use of supplied one.
A frame for the vacuum window in e.g. BDSMultilayerScreen.
G4bool checkOverlaps
Cache of checking overlaps from global constants.
virtual void SetVisAtts()
Set the visual attributes to member logVol.