BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSCollimatorJaw.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 BDSCOLLIMATORJAW_H
20#define BDSCOLLIMATORJAW_H
21
22#include "BDSCollimator.hh"
23
24#include "globals.hh" // geant4 types / globals
25#include "G4Material.hh"
26
27class G4Colour;
28class G4VSolid;
29
37{
38public:
39 BDSCollimatorJaw(const G4String& nameIn,
40 G4double lengthIn,
41 G4double horizontalWidthIn,
42 G4double xHalfGapIn,
43 G4double yHalfHeightIn,
44 G4double xsizeLeftIn,
45 G4double xsizeRightIn,
46 G4double leftJawTiltIn,
47 G4double rightJawTiltIn,
48 G4bool buildLeftJawIn,
49 G4bool buildRightJawIn,
50 G4Material* collimatorMaterialIn,
51 G4Material* vacuumMaterialIn,
52 G4Colour* colourIn = nullptr);
53 virtual ~BDSCollimatorJaw();
54
55 inline G4double GetJawTiltLeft() const {return jawTiltLeft;}
56 inline G4double GetJawTiltRight() const {return jawTiltRight;}
57
58protected:
61 virtual void CheckParameters() override;
62
64 virtual void Build() override;
65
67 virtual void BuildContainerLogicalVolume() override;
68
70 virtual void BuildInnerCollimator() final {;}
71
72 G4VSolid* jawSolid;
73 G4double xSizeLeft;
74 G4double xSizeRight;
75 G4double xHalfGap;
76 G4double jawTiltLeft;
77 G4double jawTiltRight;
78 G4double jawHalfWidth;
79 G4double yHalfHeight;
80 G4bool buildLeftJaw;
83
84private:
86 BDSCollimatorJaw() = delete;
87
92};
93
94#endif
Collimator with only two jaw and no top bit.
G4double jawHalfWidth
Half width of each jaw.
G4bool buildAperture
Build aperture or not.
virtual void Build() override
Override function in BDSCollimator for totally different construction.
G4double xSizeRight
Offset of jaw 2.
virtual void BuildContainerLogicalVolume() override
Override function in BDSCollimator for different size based container.
G4double xSizeLeft
Offset of jaw 1.
BDSCollimatorJaw(BDSCollimatorJaw &)=delete
Assignment and copy constructor not implemented nor used.
BDSCollimatorJaw & operator=(const BDSCollimatorJaw &)=delete
Assignment and copy constructor not implemented nor used.
G4double jawTiltLeft
Tilt of jaw 1 (angle in x-z plane)
G4double jawTiltRight
Tilt of jaw 2 (angle in x-z plane)
G4double yHalfHeight
Half height of each jaw.
G4bool buildRightJaw
Build right jaw or not.
virtual void BuildInnerCollimator() final
To fulfill inheritance but unused.
virtual void CheckParameters() override
G4VSolid * jawSolid
Jaw solid.
G4double xHalfGap
Half gap separation between jaws.
G4bool buildLeftJaw
Build left jaw or not.
BDSCollimatorJaw()=delete
Private default constructor to force the use of the supplied one.
Base class for collimators with common construction.