Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
OpNoviceDetectorConstruction Class Reference

#include <OpNoviceDetectorConstruction.hh>

Inheritance diagram for OpNoviceDetectorConstruction:
G4VUserDetectorConstruction

Public Member Functions

 OpNoviceDetectorConstruction ()
 
virtual ~OpNoviceDetectorConstruction ()
 
virtual G4VPhysicalVolumeConstruct ()
 
- Public Member Functions inherited from G4VUserDetectorConstruction
 G4VUserDetectorConstruction ()
 
virtual ~G4VUserDetectorConstruction ()
 
virtual void ConstructSDandField ()
 
virtual void CloneSD ()
 
virtual void CloneF ()
 
void RegisterParallelWorld (G4VUserParallelWorld *)
 
G4int ConstructParallelGeometries ()
 
void ConstructParallelSD ()
 
G4int GetNumberOfParallelWorld () const
 
G4VUserParallelWorldGetParallelWorld (G4int i) const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VUserDetectorConstruction
void SetSensitiveDetector (const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
 
void SetSensitiveDetector (G4LogicalVolume *logVol, G4VSensitiveDetector *aSD)
 

Detailed Description

Definition at line 39 of file OpNoviceDetectorConstruction.hh.

Constructor & Destructor Documentation

OpNoviceDetectorConstruction::OpNoviceDetectorConstruction ( )

Definition at line 46 of file OpNoviceDetectorConstruction.cc.

References python.hepunit::m.

48 {
49  fExpHall_x = fExpHall_y = fExpHall_z = 10.0*m;
50  fTank_x = fTank_y = fTank_z = 5.0*m;
51  fBubble_x = fBubble_y = fBubble_z = 0.5*m;
52 }
OpNoviceDetectorConstruction::~OpNoviceDetectorConstruction ( )
virtual

Definition at line 56 of file OpNoviceDetectorConstruction.cc.

56 {;}

Member Function Documentation

G4VPhysicalVolume * OpNoviceDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 60 of file OpNoviceDetectorConstruction.cc.

References test::a, G4MaterialPropertiesTable::AddConstProperty(), G4Material::AddElement(), G4MaterialPropertiesTable::AddProperty(), python.hepunit::cm3, density, dielectric_dielectric, G4OpticalSurface::DumpInfo(), python.hepunit::eV, g(), G4Material::GetIonisation(), G4LogicalSkinSurface::GetSurface(), G4LogicalSurface::GetSurfaceProperty(), glisur, ground, python.hepunit::m, python.hepunit::MeV, python.hepunit::mg, python.hepunit::mm, python.hepunit::mole, N, ns, python.hepunit::perCent, polished, G4IonisParamMat::SetBirksConstant(), G4OpticalSurface::SetFinish(), G4OpticalSurface::SetMaterialPropertiesTable(), G4Material::SetMaterialPropertiesTable(), G4OpticalSurface::SetModel(), G4PhysicsVector::SetSpline(), G4OpticalSurface::SetType(), unified, and z.

61 {
62 
63 // ------------- Materials -------------
64 
65  G4double a, z, density;
66  G4int nelements;
67 
68 // Air
69 //
70  G4Element* N = new G4Element("Nitrogen", "N", z=7 , a=14.01*g/mole);
71  G4Element* O = new G4Element("Oxygen" , "O", z=8 , a=16.00*g/mole);
72 
73  G4Material* air = new G4Material("Air", density=1.29*mg/cm3, nelements=2);
74  air->AddElement(N, 70.*perCent);
75  air->AddElement(O, 30.*perCent);
76 
77 // Water
78 //
79  G4Element* H = new G4Element("Hydrogen", "H", z=1 , a=1.01*g/mole);
80 
81  G4Material* water = new G4Material("Water", density= 1.0*g/cm3, nelements=2);
82  water->AddElement(H, 2);
83  water->AddElement(O, 1);
84 
85 //
86 // ------------ Generate & Add Material Properties Table ------------
87 //
88  const G4int nEntries = 32;
89 
90  G4double photonEnergy[nEntries] =
91  { 2.034*eV, 2.068*eV, 2.103*eV, 2.139*eV,
92  2.177*eV, 2.216*eV, 2.256*eV, 2.298*eV,
93  2.341*eV, 2.386*eV, 2.433*eV, 2.481*eV,
94  2.532*eV, 2.585*eV, 2.640*eV, 2.697*eV,
95  2.757*eV, 2.820*eV, 2.885*eV, 2.954*eV,
96  3.026*eV, 3.102*eV, 3.181*eV, 3.265*eV,
97  3.353*eV, 3.446*eV, 3.545*eV, 3.649*eV,
98  3.760*eV, 3.877*eV, 4.002*eV, 4.136*eV };
99 //
100 // Water
101 //
102  G4double refractiveIndex1[nEntries] =
103  { 1.3435, 1.344, 1.3445, 1.345, 1.3455,
104  1.346, 1.3465, 1.347, 1.3475, 1.348,
105  1.3485, 1.3492, 1.35, 1.3505, 1.351,
106  1.3518, 1.3522, 1.3530, 1.3535, 1.354,
107  1.3545, 1.355, 1.3555, 1.356, 1.3568,
108  1.3572, 1.358, 1.3585, 1.359, 1.3595,
109  1.36, 1.3608};
110 
111  G4double absorption[nEntries] =
112  {3.448*m, 4.082*m, 6.329*m, 9.174*m, 12.346*m, 13.889*m,
113  15.152*m, 17.241*m, 18.868*m, 20.000*m, 26.316*m, 35.714*m,
114  45.455*m, 47.619*m, 52.632*m, 52.632*m, 55.556*m, 52.632*m,
115  52.632*m, 47.619*m, 45.455*m, 41.667*m, 37.037*m, 33.333*m,
116  30.000*m, 28.500*m, 27.000*m, 24.500*m, 22.000*m, 19.500*m,
117  17.500*m, 14.500*m };
118 
119  G4double scintilFast[nEntries] =
120  { 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
121  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
122  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
123  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
124  1.00, 1.00, 1.00, 1.00 };
125  G4double scintilSlow[nEntries] =
126  { 0.01, 1.00, 2.00, 3.00, 4.00, 5.00, 6.00,
127  7.00, 8.00, 9.00, 8.00, 7.00, 6.00, 4.00,
128  3.00, 2.00, 1.00, 0.01, 1.00, 2.00, 3.00,
129  4.00, 5.00, 6.00, 7.00, 8.00, 9.00, 8.00,
130  7.00, 6.00, 5.00, 4.00 };
131 
133 
134  myMPT1->AddProperty("RINDEX", photonEnergy, refractiveIndex1,nEntries)
135  ->SetSpline(true);
136  myMPT1->AddProperty("ABSLENGTH", photonEnergy, absorption, nEntries)
137  ->SetSpline(true);
138  myMPT1->AddProperty("FASTCOMPONENT",photonEnergy, scintilFast, nEntries)
139  ->SetSpline(true);
140  myMPT1->AddProperty("SLOWCOMPONENT",photonEnergy, scintilSlow, nEntries)
141  ->SetSpline(true);
142 
143  myMPT1->AddConstProperty("SCINTILLATIONYIELD",50./MeV);
144  myMPT1->AddConstProperty("RESOLUTIONSCALE",1.0);
145  myMPT1->AddConstProperty("FASTTIMECONSTANT", 1.*ns);
146  myMPT1->AddConstProperty("SLOWTIMECONSTANT",10.*ns);
147  myMPT1->AddConstProperty("YIELDRATIO",0.8);
148 
149  const G4int numentries_water = 60;
150 
151  G4double energy_water[numentries_water] = {
152  1.56962*eV, 1.58974*eV, 1.61039*eV, 1.63157*eV,
153  1.65333*eV, 1.67567*eV, 1.69863*eV, 1.72222*eV,
154  1.74647*eV, 1.77142*eV, 1.7971 *eV, 1.82352*eV,
155  1.85074*eV, 1.87878*eV, 1.90769*eV, 1.93749*eV,
156  1.96825*eV, 1.99999*eV, 2.03278*eV, 2.06666*eV,
157  2.10169*eV, 2.13793*eV, 2.17543*eV, 2.21428*eV,
158  2.25454*eV, 2.29629*eV, 2.33962*eV, 2.38461*eV,
159  2.43137*eV, 2.47999*eV, 2.53061*eV, 2.58333*eV,
160  2.63829*eV, 2.69565*eV, 2.75555*eV, 2.81817*eV,
161  2.88371*eV, 2.95237*eV, 3.02438*eV, 3.09999*eV,
162  3.17948*eV, 3.26315*eV, 3.35134*eV, 3.44444*eV,
163  3.54285*eV, 3.64705*eV, 3.75757*eV, 3.87499*eV,
164  3.99999*eV, 4.13332*eV, 4.27585*eV, 4.42856*eV,
165  4.59258*eV, 4.76922*eV, 4.95999*eV, 5.16665*eV,
166  5.39129*eV, 5.63635*eV, 5.90475*eV, 6.19998*eV
167  };
168 
169  //assume 100 times larger than the rayleigh scattering for now.
170  G4double mie_water[numentries_water] = {
171  167024.4*m, 158726.7*m, 150742 *m,
172  143062.5*m, 135680.2*m, 128587.4*m,
173  121776.3*m, 115239.5*m, 108969.5*m,
174  102958.8*m, 97200.35*m, 91686.86*m,
175  86411.33*m, 81366.79*m, 76546.42*m,
176  71943.46*m, 67551.29*m, 63363.36*m,
177  59373.25*m, 55574.61*m, 51961.24*m,
178  48527.00*m, 45265.87*m, 42171.94*m,
179  39239.39*m, 36462.50*m, 33835.68*m,
180  31353.41*m, 29010.30*m, 26801.03*m,
181  24720.42*m, 22763.36*m, 20924.88*m,
182  19200.07*m, 17584.16*m, 16072.45*m,
183  14660.38*m, 13343.46*m, 12117.33*m,
184  10977.70*m, 9920.416*m, 8941.407*m,
185  8036.711*m, 7202.470*m, 6434.927*m,
186  5730.429*m, 5085.425*m, 4496.467*m,
187  3960.210*m, 3473.413*m, 3032.937*m,
188  2635.746*m, 2278.907*m, 1959.588*m,
189  1675.064*m, 1422.710*m, 1200.004*m,
190  1004.528*m, 833.9666*m, 686.1063*m
191  };
192 
193  // gforward, gbackward, forward backward ratio
194  G4double mie_water_const[3]={0.99,0.99,0.8};
195 
196  myMPT1->AddProperty("MIEHG",energy_water,mie_water,numentries_water)
197  ->SetSpline(true);
198  myMPT1->AddConstProperty("MIEHG_FORWARD",mie_water_const[0]);
199  myMPT1->AddConstProperty("MIEHG_BACKWARD",mie_water_const[1]);
200  myMPT1->AddConstProperty("MIEHG_FORWARD_RATIO",mie_water_const[2]);
201 
202  water->SetMaterialPropertiesTable(myMPT1);
203 
204  // Set the Birks Constant for the Water scintillator
205 
206  water->GetIonisation()->SetBirksConstant(0.126*mm/MeV);
207 
208 //
209 // Air
210 //
211  G4double refractiveIndex2[nEntries] =
212  { 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
213  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
214  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
215  1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00,
216  1.00, 1.00, 1.00, 1.00 };
217 
219  myMPT2->AddProperty("RINDEX", photonEnergy, refractiveIndex2, nEntries);
220 
221  air->SetMaterialPropertiesTable(myMPT2);
222 
223 //
224 // ------------- Volumes --------------
225 
226 // The experimental Hall
227 //
228  G4Box* expHall_box = new G4Box("World",fExpHall_x,fExpHall_y,fExpHall_z);
229 
230  G4LogicalVolume* expHall_log
231  = new G4LogicalVolume(expHall_box,air,"World",0,0,0);
232 
233  G4VPhysicalVolume* expHall_phys
234  = new G4PVPlacement(0,G4ThreeVector(),expHall_log,"World",0,false,0);
235 
236 // The Water Tank
237 //
238  G4Box* waterTank_box = new G4Box("Tank",fTank_x,fTank_y,fTank_z);
239 
240  G4LogicalVolume* waterTank_log
241  = new G4LogicalVolume(waterTank_box,water,"Tank",0,0,0);
242 
243  G4VPhysicalVolume* waterTank_phys
244  = new G4PVPlacement(0,G4ThreeVector(),waterTank_log,"Tank",
245  expHall_log,false,0);
246 
247 // The Air Bubble
248 //
249  G4Box* bubbleAir_box = new G4Box("Bubble",fBubble_x,fBubble_y,fBubble_z);
250 
251  G4LogicalVolume* bubbleAir_log
252  = new G4LogicalVolume(bubbleAir_box,air,"Bubble",0,0,0);
253 
254 //G4VPhysicalVolume* bubbleAir_phys =
255  new G4PVPlacement(0,G4ThreeVector(0,2.5*m,0),bubbleAir_log,"Bubble",
256  waterTank_log,false,0);
257 
258 // ------------- Surfaces --------------
259 //
260 // Water Tank
261 //
262  G4OpticalSurface* opWaterSurface = new G4OpticalSurface("WaterSurface");
263  opWaterSurface->SetType(dielectric_dielectric);
264  opWaterSurface->SetFinish(ground);
265  opWaterSurface->SetModel(unified);
266 
267  new G4LogicalBorderSurface("WaterSurface",
268  waterTank_phys,expHall_phys,opWaterSurface);
269 
270 // Air Bubble
271 //
272  G4OpticalSurface* opAirSurface = new G4OpticalSurface("AirSurface");
273  opAirSurface->SetType(dielectric_dielectric);
274  opAirSurface->SetFinish(polished);
275  opAirSurface->SetModel(glisur);
276 
277  G4LogicalSkinSurface* airSurface =
278  new G4LogicalSkinSurface("AirSurface", bubbleAir_log, opAirSurface);
279 
280  G4OpticalSurface* opticalSurface = dynamic_cast <G4OpticalSurface*>
281  (airSurface->GetSurface(bubbleAir_log)->GetSurfaceProperty());
282 
283  if (opticalSurface) opticalSurface->DumpInfo();
284 
285 //
286 // Generate & Add Material Properties Table attached to the optical surfaces
287 //
288  const G4int num = 2;
289  G4double ephoton[num] = {2.034*eV, 4.136*eV};
290 
291  //OpticalWaterSurface
292  G4double refractiveIndex[num] = {1.35, 1.40};
293  G4double specularLobe[num] = {0.3, 0.3};
294  G4double specularSpike[num] = {0.2, 0.2};
295  G4double backScatter[num] = {0.2, 0.2};
296 
298 
299  myST1->AddProperty("RINDEX", ephoton, refractiveIndex, num);
300  myST1->AddProperty("SPECULARLOBECONSTANT", ephoton, specularLobe, num);
301  myST1->AddProperty("SPECULARSPIKECONSTANT", ephoton, specularSpike, num);
302  myST1->AddProperty("BACKSCATTERCONSTANT", ephoton, backScatter, num);
303 
304  opWaterSurface->SetMaterialPropertiesTable(myST1);
305 
306  //OpticalAirSurface
307  G4double reflectivity[num] = {0.3, 0.5};
308  G4double efficiency[num] = {0.8, 1.0};
309 
311 
312  myST2->AddProperty("REFLECTIVITY", ephoton, reflectivity, num);
313  myST2->AddProperty("EFFICIENCY", ephoton, efficiency, num);
314 
315  opAirSurface->SetMaterialPropertiesTable(myST2);
316 
317 //always return the physical World
318  return expHall_phys;
319 }
void SetFinish(const G4OpticalSurfaceFinish)
G4IonisParamMat * GetIonisation() const
Definition: G4Material.hh:224
CLHEP::Hep3Vector G4ThreeVector
G4double z
Definition: TRTMaterials.hh:39
Definition: G4Box.hh:63
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
Definition: G4Material.hh:247
void SetBirksConstant(G4double value)
void DumpInfo() const
int G4int
Definition: G4Types.hh:78
G4MaterialPropertyVector * AddProperty(const char *key, G4double *PhotonEnergies, G4double *PropertyValues, G4int NumEntries)
void SetSpline(G4bool)
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
void AddConstProperty(const char *key, G4double PropertyValue)
float perCent
Definition: hepunit.py:239
**D E S C R I P T I O N
G4SurfaceProperty * GetSurfaceProperty() const
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:345
double G4double
Definition: G4Types.hh:76
void SetModel(const G4OpticalSurfaceModel model)
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
#define ns
Definition: xmlparse.cc:597
void SetType(const G4SurfaceType &type)
static G4LogicalSkinSurface * GetSurface(const G4LogicalVolume *vol)

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