BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSShield.cc
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#include "BDSBeamPipe.hh"
20#include "BDSBeamPipeFactory.hh"
21#include "BDSBeamPipeInfo.hh"
22#include "BDSDebug.hh"
23#include "BDSSDType.hh"
24#include "BDSShield.hh"
25#include "BDSUtilities.hh"
26
27#include "globals.hh"
28#include "G4Box.hh"
29#include "G4LogicalVolume.hh"
30#include "G4SubtractionSolid.hh"
31#include "G4VisAttributes.hh"
32#include "G4PVPlacement.hh"
33
34class G4Colour;
35class G4Material;
36
37BDSShield::BDSShield(const G4String& nameIn,
38 G4double lengthIn,
39 G4double horizontalWidthIn,
40 G4double xSizeIn,
41 G4double ySizeIn,
42 G4Material* materialIn,
43 G4Colour* colourIn,
44 BDSBeamPipeInfo* beamPipeInfoIn):
45 BDSAcceleratorComponent(nameIn, lengthIn, 0, "shield", beamPipeInfoIn),
46 horizontalWidth(horizontalWidthIn),
47 xSize(xSizeIn),
48 ySize(ySizeIn),
49 material(materialIn),
50 colour(colourIn)
51{;}
52
53BDSShield::~BDSShield()
54{;}
55
57{
61}
62
64{
65 containerSolid = new G4Box(name + "_container_solid",
68 chordLength*0.5);
69 containerLogicalVolume = new G4LogicalVolume(containerSolid,
71 name + "_container_lv");
72}
73
75{
76 G4VSolid* outerSolid = new G4Box(name + "_outer_solid",
80 RegisterSolid(outerSolid);
81
82 G4VSolid* shieldSolid;
83 G4LogicalVolume* shieldLV;
84
85 // only subtract inner solid if shield aperture is finite.
87 {
88 G4VSolid* innerSolid = new G4Box(name + "_inner_solid",
89 xSize,
90 ySize,
91 chordLength); // extra long for unambiguous subtraction
92 RegisterSolid(innerSolid);
93
94 shieldSolid = new G4SubtractionSolid(name + "shield_solid",
95 outerSolid, // this
96 innerSolid); // minus this
97 RegisterSolid(shieldSolid);
98
99 shieldLV = new G4LogicalVolume(shieldSolid,
100 material,
101 name+"_shield_lv");
102 }
103 else
104 {
105 shieldLV = new G4LogicalVolume(outerSolid,
106 material,
107 name+"_shield_lv");
108 }
109
110 RegisterLogicalVolume(shieldLV);
111 if (sensitiveOuter)
112 {RegisterSensitiveVolume(shieldLV, BDSSDType::energydep);}
113
114 G4VisAttributes* shieldVisAttr = new G4VisAttributes(*colour);
115 shieldVisAttr->SetVisibility(true);
116 shieldLV->SetVisAttributes(shieldVisAttr);
117 RegisterVisAttributes(shieldVisAttr);
118
119 G4PVPlacement* shieldPV = new G4PVPlacement(nullptr,
120 G4ThreeVector(),
121 shieldLV,
122 name + "_shield_pv",
123 containerLogicalVolume,
124 false,
125 0,
127
128 RegisterPhysicalVolume(shieldPV);
129}
130
132{
133 // check we should build a beam pipe
134 if (!beamPipeInfo)
135 {return;}
136
137 // check beam pipe fits
139 {
140 G4cout << __METHOD_NAME__ << "Shield will not fit around beam pipe - not building beam pipe!" << G4endl << G4endl;
141 return;
142 }
143
144 // construct and place beam pipe
145 auto bp = BDSBeamPipeFactory::Instance()->CreateBeamPipe(name,
149
150 G4PVPlacement* bpPV = new G4PVPlacement(nullptr,
151 G4ThreeVector(),
152 bp->GetContainerLogicalVolume(),
153 name+"_beampipe_pv",
154 containerLogicalVolume,
155 false,
156 0,
158
160}
Abstract class that represents a component of an accelerator.
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4bool sensitiveOuter
Useful variable often used in construction.
static G4double lengthSafety
Useful variable often used in construction.
static G4bool checkOverlaps
Useful variable often used in construction.
static G4Material * worldMaterial
Useful variable often used in construction.
G4double chordLength
Protected member variable that can be modified by derived classes.
BDSBeamPipeInfo * beamPipeInfo
Optional beam pipe recipe that is written out to the survey if it exists.
static BDSBeamPipeFactory * Instance()
Singleton accessor.
Holder class for all information required to describe a beam pipe model.
G4double aper1
Public member for direct access.
G4double beamPipeThickness
Public member for direct access.
G4double aper2
Public member for direct access.
void RegisterDaughter(BDSGeometryComponent *anotherComponent)
void RegisterLogicalVolume(G4LogicalVolume *logicalVolume)
void RegisterPhysicalVolume(G4VPhysicalVolume *physicalVolume)
void RegisterVisAttributes(G4VisAttributes *visAttribute)
void RegisterSolid(G4VSolid *solid)
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)
void BuildBeamPipe()
Build a beam pipe in the hole if required.
Definition: BDSShield.cc:131
virtual void Build()
Definition: BDSShield.cc:56
void BuildShield()
Build the outer shield geoemtry.
Definition: BDSShield.cc:74
BDSShield()=delete
Default constructor, assignment and copy constructor not used.
G4double horizontalWidth
Outer size of shield.
Definition: BDSShield.hh:73
G4Material * material
Shield material.
Definition: BDSShield.hh:76
G4double ySize
Inner vertical half width of shield.
Definition: BDSShield.hh:75
G4Colour * colour
Colour of shielding block.
Definition: BDSShield.hh:77
G4double xSize
Inner horizontal half width of shield.
Definition: BDSShield.hh:74
virtual void BuildContainerLogicalVolume()
Build a container volume for everything.
Definition: BDSShield.cc:63
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())