Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
OpNoviceDetectorConstruction Class Reference

#include <OpNoviceDetectorConstruction.hh>

Inheritance diagram for OpNoviceDetectorConstruction:
Collaboration diagram for OpNoviceDetectorConstruction:

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 42 of file OpNoviceDetectorConstruction.hh.

Constructor & Destructor Documentation

OpNoviceDetectorConstruction::OpNoviceDetectorConstruction ( )

Definition at line 49 of file OpNoviceDetectorConstruction.cc.

51 {
52  fExpHall_x = fExpHall_y = fExpHall_z = 10.0*m;
53  fTank_x = fTank_y = fTank_z = 5.0*m;
54  fBubble_x = fBubble_y = fBubble_z = 0.5*m;
55 }
static constexpr double m
Definition: G4SIunits.hh:129
OpNoviceDetectorConstruction::~OpNoviceDetectorConstruction ( )
virtual

Definition at line 59 of file OpNoviceDetectorConstruction.cc.

59 {;}

Member Function Documentation

G4VPhysicalVolume * OpNoviceDetectorConstruction::Construct ( void  )
virtual

Implements G4VUserDetectorConstruction.

Definition at line 63 of file OpNoviceDetectorConstruction.cc.

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

Here is the call graph for this function:


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