Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
B03DetectorConstruction Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for B03DetectorConstruction:
G4VUserDetectorConstruction

Public Member Functions

 B03DetectorConstruction ()
 
 ~B03DetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
G4VPhysicalVolumeGetWorldVolume ()
 
G4VPhysicalVolumeGetWorldVolumeAddress () const
 
void SetSensitive ()
 

Private Attributes

G4VPhysicalVolumefWorldVolume
 

Detailed Description

Definition at line 50 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ B03DetectorConstruction()

B03DetectorConstruction::B03DetectorConstruction ( )

◆ ~B03DetectorConstruction()

B03DetectorConstruction::~B03DetectorConstruction ( )

Definition at line 67 of file B03DetectorConstruction.cc.

68{;}

Member Function Documentation

◆ Construct()

G4VPhysicalVolume * B03DetectorConstruction::Construct ( )
virtual

Definition at line 72 of file B03DetectorConstruction.cc.

73{
74 G4double pos_x;
75 G4double pos_y;
76 G4double pos_z;
77
78 G4double density, pressure, temperature;
79 G4double A;
80 G4int Z;
81
82 G4String name, symbol;
83 G4double z;
84 G4double fractionmass;
85
86 A = 1.01*g/mole;
87 G4Element* elH = new G4Element(name="Hydrogen",symbol="H" , Z= 1, A);
88
89 A = 12.01*g/mole;
90 G4Element* elC = new G4Element(name="Carbon" ,symbol="C" , Z = 6, A);
91
92 A = 16.00*g/mole;
93 G4Element* elO = new G4Element(name="Oxygen" ,symbol="O" , Z= 8, A);
94
95 A = 22.99*g/mole;
96 G4Element* elNa = new G4Element(name="Natrium" ,symbol="Na" , Z=11 , A);
97
98 A = 200.59*g/mole;
99 G4Element* elHg = new G4Element(name="Hg" ,symbol="Hg" , Z=80, A);
100
101 A = 26.98*g/mole;
102 G4Element* elAl = new G4Element(name="Aluminium" ,symbol="Al" , Z=13, A);
103
104 A = 28.09*g/mole;
105 G4Element* elSi = new G4Element(name="Silicon", symbol="Si", Z=14, A);
106
107 A = 39.1*g/mole;
108 G4Element* elK = new G4Element(name="K" ,symbol="K" , Z=19 , A);
109
110 A = 69.72*g/mole;
111 G4Element* elCa = new G4Element(name="Calzium" ,symbol="Ca" , Z=31 , A);
112
113 A = 55.85*g/mole;
114 G4Element* elFe = new G4Element(name="Iron" ,symbol="Fe", Z=26, A);
115
116 density = universe_mean_density; //from PhysicalConstants.h
117 pressure = 3.e-18*pascal;
118 temperature = 2.73*kelvin;
119 G4Material *Galactic =
120 new G4Material(name="Galactic", z=1., A=1.01*g/mole, density,
121 kStateGas,temperature,pressure);
122
123 density = 2.03*g/cm3;
124 G4Material* Concrete = new G4Material("Concrete", density, 10);
125 Concrete->AddElement(elH , fractionmass= 0.01);
126 Concrete->AddElement(elO , fractionmass= 0.529);
127 Concrete->AddElement(elNa , fractionmass= 0.016);
128 Concrete->AddElement(elHg , fractionmass= 0.002);
129 Concrete->AddElement(elAl , fractionmass= 0.034);
130 Concrete->AddElement(elSi , fractionmass= 0.337);
131 Concrete->AddElement(elK , fractionmass= 0.013);
132 Concrete->AddElement(elCa , fractionmass= 0.044);
133 Concrete->AddElement(elFe , fractionmass= 0.014);
134 Concrete->AddElement(elC , fractionmass= 0.001);
135
136
137 /////////////////////////////
138 // world cylinder volume
139 ////////////////////////////
140
141 // world solid
142
143 G4double innerRadiusCylinder = 0*cm;
144 G4double outerRadiusCylinder = 101*cm; // dont't have scoring
145 // cells coinside eith world volume boundary
146 // G4double heightCylinder = 105*cm;
147 G4double heightCylinder = 100*cm;
148 G4double startAngleCylinder = 0*deg;
149 G4double spanningAngleCylinder = 360*deg;
150
151 G4Tubs *worldCylinder = new G4Tubs("worldCylinder",
152 innerRadiusCylinder,
153 outerRadiusCylinder,
154 heightCylinder,
155 startAngleCylinder,
156 spanningAngleCylinder);
157
158 // logical world
159
160 G4LogicalVolume *worldCylinder_log =
161 new G4LogicalVolume(worldCylinder, Galactic, "worldCylinder_log");
162
163 name = "shieldWorld";
164 fWorldVolume = new G4PVPlacement(0, G4ThreeVector(0,0,0), worldCylinder_log
165 ,name, 0, false, 0);
166
167 // creating 18 slobs of 10 cm thick concrete
168
169 G4double innerRadiusShield = 0*cm;
170 G4double outerRadiusShield = 100*cm;
171 G4double heightShield = 90*cm;
172 G4double startAngleShield = 0*deg;
173 G4double spanningAngleShield = 360*deg;
174
175 G4Tubs *aShield = new G4Tubs("aShield",
176 innerRadiusShield,
177 outerRadiusShield,
178 heightShield,
179 startAngleShield,
180 spanningAngleShield);
181
182 // logical shield
183
184 G4LogicalVolume *aShield_log =
185 new G4LogicalVolume(aShield, Concrete, "aShield_log");
186
187 G4VisAttributes* pShieldVis = new
188 G4VisAttributes(G4Colour(0.0,0.0,1.0));
189 pShieldVis->SetForceSolid(true);
190 aShield_log->SetVisAttributes(pShieldVis);
191
192 // physical shields
193
194 name = "concreteShield";
195
196 pos_x = 0*cm;
197 pos_y = 0*cm;
198 pos_z = 0;
199
200 new G4PVPlacement(0,
201 G4ThreeVector(pos_x, pos_y, pos_z),
202 aShield_log,
203 name,
204 worldCylinder_log,
205 false,
206 0);
207
208
209 return fWorldVolume;
210}

◆ GetWorldVolume()

G4VPhysicalVolume * B03DetectorConstruction::GetWorldVolume ( )

Definition at line 221 of file B03DetectorConstruction.cc.

221 {
222 return fWorldVolume;
223}

◆ GetWorldVolumeAddress()

G4VPhysicalVolume & B03DetectorConstruction::GetWorldVolumeAddress ( ) const

Definition at line 227 of file B03DetectorConstruction.cc.

227 {
228 return *fWorldVolume;
229}

◆ SetSensitive()

void B03DetectorConstruction::SetSensitive ( )

Member Data Documentation

◆ fWorldVolume

G4VPhysicalVolume* B03DetectorConstruction::fWorldVolume
private

Definition at line 60 of file B03DetectorConstruction.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration