BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSBeamPipe.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 BDSBEAMPIPE_H
20#define BDSBEAMPIPE_H
21
22#include "BDSGeometryComponent.hh"
23
24#include "globals.hh"
25#include "G4ThreeVector.hh"
26
27#include <set>
28
29class BDSExtent;
30class G4LogicalVolume;
31class G4VSolid;
32
45{
46public:
49 BDSBeamPipe(G4VSolid* containerSolidIn,
50 G4LogicalVolume* containerLVIn,
51 BDSExtent extentIn,
52 G4VSolid* containerSubtractionSolidIn,
53 G4LogicalVolume* vacuumLVIn,
54 G4bool containerIsCircularIn = false,
55 G4double containerRadiusIn = 0.0,
56 G4ThreeVector inputFaceNormalIn = G4ThreeVector(0,0,-1),
57 G4ThreeVector outputFaceNormalIn = G4ThreeVector(0,0, 1));
58
59 virtual ~BDSBeamPipe();
60
63 inline G4VSolid* GetContainerSubtractionSolid() const {return containerSubtractionSolid;}
65 inline G4LogicalVolume* GetVacuumLogicalVolume() const {return vacuumLogicalVolume;};
68 inline G4bool ContainerIsCircular() const {return containerIsCircular;}
70 inline G4double GetContainerRadius() const {return containerRadius;}
71
73 inline G4ThreeVector InputFaceNormal() const {return inputFaceNormal;}
74 inline G4ThreeVector OutputFaceNormal() const {return outputFaceNormal;}
76
79 std::set<G4LogicalVolume*> GetVolumesForField() const;
80
81protected:
82 G4VSolid* containerSubtractionSolid;
83 G4LogicalVolume* vacuumLogicalVolume;
84 G4bool containerIsCircular;
85 G4double containerRadius;
86 G4ThreeVector inputFaceNormal;
87 G4ThreeVector outputFaceNormal;
88};
89
90#endif
A holder class for a piece of beam pipe geometry.
Definition: BDSBeamPipe.hh:45
G4LogicalVolume * GetVacuumLogicalVolume() const
Access the vacuum volume to set fields and limits.
Definition: BDSBeamPipe.hh:65
G4ThreeVector InputFaceNormal() const
Accessor.
Definition: BDSBeamPipe.hh:73
G4VSolid * GetContainerSubtractionSolid() const
default destructor sufficient as G4 manages solids and LVs
Definition: BDSBeamPipe.hh:63
BDSBeamPipe(G4VSolid *containerSolidIn, G4LogicalVolume *containerLVIn, BDSExtent extentIn, G4VSolid *containerSubtractionSolidIn, G4LogicalVolume *vacuumLVIn, G4bool containerIsCircularIn=false, G4double containerRadiusIn=0.0, G4ThreeVector inputFaceNormalIn=G4ThreeVector(0, 0,-1), G4ThreeVector outputFaceNormalIn=G4ThreeVector(0, 0, 1))
Definition: BDSBeamPipe.cc:25
G4double GetContainerRadius() const
If it is circular, we need the radius.
Definition: BDSBeamPipe.hh:70
G4bool ContainerIsCircular() const
Definition: BDSBeamPipe.hh:68
G4ThreeVector OutputFaceNormal() const
Accessor.
Definition: BDSBeamPipe.hh:74
std::set< G4LogicalVolume * > GetVolumesForField() const
Definition: BDSBeamPipe.cc:50
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
A generic geometry component for a bdsim model.