BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSMagnetOuter.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 BDSMAGNETOUTER_H
20#define BDSMAGNETOUTER_H
21
22#include "BDSGeometryComponent.hh"
23
24#include "globals.hh"
25#include "G4ThreeVector.hh"
26
29class G4LogicalVolume;
30class G4VSolid;
31
46{
47public:
48 BDSMagnetOuter(G4VSolid* containerSolidIn,
49 G4LogicalVolume* containerLVIn,
50 const BDSExtent& extent,
51 BDSGeometryComponent* magnetContainerIn,
52 const G4ThreeVector& placementOffset = G4ThreeVector(0,0,0),
53 BDSSimpleComponent* endPieceBeforeIn = nullptr,
54 BDSSimpleComponent* endPieceAfterIn = nullptr,
55 const G4ThreeVector& inputFaceNormalIn = G4ThreeVector(0,0,-1),
56 const G4ThreeVector& outputFaceNormalIn = G4ThreeVector(0,0, 1));
58 BDSGeometryComponent* magnetContainerIn,
59 BDSSimpleComponent* endPieceBeforeIn = nullptr,
60 BDSSimpleComponent* endPieceAfterIn = nullptr,
61 const G4ThreeVector& inputFaceNormalIn = G4ThreeVector(0,0,-1),
62 const G4ThreeVector& outputFaceNormalIn = G4ThreeVector(0,0, 1));
64 BDSGeometryComponent* magnetContainerIn);
65 virtual ~BDSMagnetOuter();
66
71 BDSGeometryComponent* GetMagnetContainer() const {return magnetContainer;}
72
74 BDSSimpleComponent* EndPieceBefore() const {return endPieceBefore;}
75 BDSSimpleComponent* EndPieceAfter() const {return endPieceAfter;}
77
78 void SetEndPieceBefore(BDSSimpleComponent* endPieceIn) {endPieceBefore = endPieceIn;}
79 void SetEndPieceAfter(BDSSimpleComponent* endPieceIn) {endPieceAfter = endPieceIn;}
80
83
85 void ClearEndPieces();
86
88 inline G4ThreeVector InputFaceNormal() const {return inputFaceNormal;}
89 inline G4ThreeVector OutputFaceNormal() const {return outputFaceNormal;}
91
93 void SetInputFaceNormal(const G4ThreeVector& input);
94 void SetOutputFaceNormal(const G4ThreeVector& output);
96
99 BDSGeometryExternal* ExternalGeometry() const {return externalGeometry;}
100
101protected:
102 BDSGeometryComponent* magnetContainer;
103 BDSSimpleComponent* endPieceBefore;
104 BDSSimpleComponent* endPieceAfter;
105 G4ThreeVector inputFaceNormal;
106 G4ThreeVector outputFaceNormal;
107 BDSGeometryExternal* externalGeometry;
108};
109
110#endif
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
A generic geometry component for a bdsim model.
A loaded piece of externally provided geometry.
An object for both the returned magnet outer body but also a tight fitting container for the whole ma...
BDSGeometryComponent * GetMagnetContainer() const
void ClearMagnetContainer()
Clear the memory of the now unneeded magnet container object.
void ClearEndPieces()
Clear the memory of the possibly unneeded end piece objects.
void SetInputFaceNormal(const G4ThreeVector &input)
Setter for face normals.
G4ThreeVector InputFaceNormal() const
Accessor.
BDSSimpleComponent * EndPieceAfter() const
Access the end piece.
BDSGeometryExternal * ExternalGeometry() const
G4ThreeVector OutputFaceNormal() const
Accessor.
void SetOutputFaceNormal(const G4ThreeVector &output)
Setter for face normals.
BDSSimpleComponent * EndPieceBefore() const
Access the end piece.
A BDSAcceleratorComponent wrapper for BDSGeometryComponent.