WCSim
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
WCSimConstructMaterials.cc
Go to the documentation of this file.
3 
4 #include "G4Material.hh"
5 #include "G4Element.hh"
6 #include "globals.hh"
7 #include "G4UnitsTable.hh"
8 
9 #include "G4LogicalBorderSurface.hh"
10 #include "G4LogicalSkinSurface.hh"
11 #include "G4OpBoundaryProcess.hh"
12 
13 #include "G4PhysicalConstants.hh"
14 #include "G4SystemOfUnits.hh"
15 
17 {
18  //****Materials Definitions****
19 
20  G4double density;
21  G4double a;
22 
23  //---Vacuum
24 
25  density = universe_mean_density; //from PhysicalConstants.h
26  G4double pressure = 1.e-19*pascal;
27  G4double temperature = 0.1*kelvin;
28  a = 1.01*g/mole;
29  G4Material* Vacuum =
30  new G4Material("Vacuum", 1., a, density,
31  kStateGas,temperature,pressure);
32 
33  //---Water
34 
35  a = 1.01*g/mole;
36  G4Element* elH
37  = new G4Element("Hydrogen","H", 1,a);
38 
39  a = 16.00*g/mole;
40  G4Element* elO
41  = new G4Element("Oxygen","O", 8,a);
42 
43  density = 1.00*g/cm3;
44  G4Material* Water
45  = new G4Material("Water",density,2);
46  Water->AddElement(elH, 2);
47  Water->AddElement(elO, 1);
48 
49  //---Gd doped Water
50 
51  a = 157.25*g/mole;
52  G4Element* Gd
53  = new G4Element("Gadolinium","Gd", 64,a);
54 
55  density = 1.00*g/cm3;
56  G4Material* DopedWater
57  = new G4Material("Doped Water",density,2);
58  DopedWater->AddMaterial(Water,99.9*perCent);
59  DopedWater->AddElement(Gd,0.1*perCent);
60 
61  //---Ice
62 
63  density = 0.92*g/cm3;//Ice
64  G4Material* Ice = new G4Material("Ice",density,2);
65  Ice->AddElement(elH, 2);
66  Ice->AddElement(elO, 1);
67 
68  //---Steel
69 
70  a= 12.01*g/mole;
71  G4Element* elC
72  = new G4Element("Carbon","C", 6,a);
73 
74  a = 55.85*g/mole;
75  G4Element* elFe
76  = new G4Element("Iron","Fe", 26,a);
77 
78  density = 7.8*g/cm3;
79  G4Material* Steel
80  = new G4Material("Steel",density,2);
81  Steel->AddElement(elC, 1.*perCent);
82  Steel->AddElement(elFe, 99.*perCent);
83 
84  //---Stainless steel 304L (simple example, particular alloy can be different)
85 
86  a = 51.9961*g/mole;
87  G4Element* elCr = new G4Element("Chromium", "Cr", 24., a);
88 
89  a = 58.6934*g/mole;
90  G4Element* elNi = new G4Element("Nickel", "Ni", 28., a);
91 
92  a = 54.938*g/mole;
93  G4Element* elMn = new G4Element("Manganese", "Mn", 25., a);
94 
95  a = 30.974*g/mole;
96  G4Element* elP = new G4Element("Phosphore", "P", 15., a);
97 
98  a = 28.09*g/mole;
99  G4Element* elSi = new G4Element("Silicon", "Si", 14., a);
100 
101  a = 32.065*g/mole;
102  G4Element* elS = new G4Element("Sulphur", "S", 16., a);
103 
104  density = 7.81*g/cm3;
105  G4Material* StainlessSteel = new G4Material("StainlessSteel", density, 8);
106 
107  StainlessSteel->AddElement(elFe, 70.44*perCent);
108  StainlessSteel->AddElement(elCr, 18*perCent);
109  StainlessSteel->AddElement(elC, 0.08*perCent);
110  StainlessSteel->AddElement(elNi, 8*perCent);
111  StainlessSteel->AddElement(elP, 0.45*perCent);
112  StainlessSteel->AddElement(elSi, 1*perCent);
113  StainlessSteel->AddElement(elMn, 2*perCent);
114  StainlessSteel->AddElement(elS, 0.03*perCent);
115 
116  G4MaterialPropertiesTable *mpt = new G4MaterialPropertiesTable();
117 
118  const G4int nEntries = 2;
119  G4double photonEnergy[nEntries] = {1.*eV , 7.*eV};
120 
121  //G4double rindex_Steel[nEntries] = {1.462 , 1.462}; // No I haven't gone mad
122  G4double abslength_Steel[nEntries] = {.001*mm , .001*mm};
123  //mpt->AddProperty("RINDEX", photonEnergy, rindex_Steel, nEntries);
124  mpt->AddProperty("ABSLENGTH", photonEnergy, abslength_Steel, nEntries);
125 
126  StainlessSteel->SetMaterialPropertiesTable(mpt);
127 
128  //---Solid Dry Ice
129 
130  density = 1.563*g/cm3;
131  G4Material* DryIce = new G4Material("SolidDryIce", density, 2);
132  DryIce->AddElement(elC, 1);
133  DryIce->AddElement(elO, 2);
134 
135 
136 
137  //---Air
138 
139  a = 14.01*g/mole;
140  G4Element* elN
141  = new G4Element("Nitrogen","N", 7,a);
142 
143  density = 1.290*mg/cm3;
144  G4Material* Air
145  = new G4Material("Air",density,2);
146  Air->AddElement(elN, 70.*perCent);
147  Air->AddElement(elO, 30.*perCent);
148 
149  //---Plastic
150 
151  density = 0.95*g/cm3;
152  G4Material* Plastic
153  = new G4Material("Plastic",density,2);
154  Plastic->AddElement(elC, 1);
155  Plastic->AddElement(elH, 2);
156 
157  //---Aluminum (element and material)
158 
159  a = 26.98*g/mole;
160  G4Element* elAl = new G4Element("Aluminum", "Al", 13, a);
161 
162  density = 2.7*g/cm3;
163  G4Material* Aluminum
164  = new G4Material("Aluminum",density,1);
165  Aluminum->AddElement(elAl, 1);
166 
167  //---Black sheet
168 
169  density = 0.95*g/cm3;
170  G4Material* Blacksheet
171  = new G4Material("Blacksheet",density,2);
172  Blacksheet->AddElement(elC, 1);
173  Blacksheet->AddElement(elH, 2);
174 
175  //---Tyvek - jl145
176 
177  density = 0.38*g/cm3; //cf. DuPont product handbook
178  G4Material* Tyvek
179  = new G4Material("Tyvek",density,2);
180  Tyvek->AddElement(elC, 1); //polyethylene
181  Tyvek->AddElement(elH, 2);
182 
183  //---Glass
184 
185  density = 2.20*g/cm3;
186  G4Material* SiO2 = new G4Material("SiO2",density,2);
187  SiO2->AddElement(elSi, 1);
188  SiO2->AddElement(elO , 2);
189 
190  a = 10.81*g/mole;
191  G4Element* elB = new G4Element("Boron", "B", 5, a);
192 
193  density = 2.46*g/cm3;
194  G4Material* B2O3 = new G4Material("B2O3",density,2);
195  B2O3->AddElement(elB, 2);
196  B2O3->AddElement(elO, 3);
197 
198  a = 22.99*g/mole;
199  G4Element* elNa = new G4Element("Sodium", "Na", 11, a);
200 
201  density = 2.27*g/cm3;
202  G4Material* Na2O = new G4Material("Na2O",density,2);
203  Na2O->AddElement(elNa, 2);
204  Na2O->AddElement(elO, 1);
205 
206  density = 4.00*g/cm3;
207  G4Material* Al2O3 = new G4Material("Al2O3",density,2);
208  Al2O3->AddElement(elAl, 2);
209  Al2O3->AddElement(elO, 3);
210 
211 // G4Material* blackAcryl
212 // = new G4Material("blackAcryl", density, 3);
213 // blackAcryl -> AddElement(elH, 6);
214 // blackAcryl -> AddElement(elC, 4);
215 // blackAcryl -> AddElement(elO, 2);
216 
217  density = 2.23*g/cm3;
218  G4Material* Glass
219  = new G4Material("Glass",density,4);
220  //G4Material* Glass
221  //= new G4Material("Glass",density,8); //Put in 8 materials later
222 
223  Glass->AddMaterial(SiO2, 80.6*perCent);
224  Glass->AddMaterial(B2O3, 13.0*perCent);
225  Glass->AddMaterial(Na2O, 4.0*perCent);
226  Glass->AddMaterial(Al2O3, 2.4*perCent);
227  //Glass->AddMaterial(Al2O3, 2.3*perCent);
228  //Put in 2.3 percent if the other 4 materials = 0.1 percent
229 
230  //---Rock
231 
232  // a = 16.00*g/mole; G4Element* elO = new G4Element("Oxygen","O", 8,a);
233  // a = 28.09*g/mole; G4Element* elSi = new G4Element("Silicon", "Si", 14., a);
234  // a = 26.98*g/mole; G4Element* elAl = new G4Element("Aluminum", "Al", 13, a);
235  // a = 55.85*g/mole; G4Element* elFe = new G4Element("Iron","Fe", 26,a);
236  a = 40.08*g/mole; G4Element* elCa = new G4Element("Calcium","Ca", 20,a);
237  // a = 22.99*g/mole; G4Element* elNa = new G4Element("Sodium", "Na", 11, a);
238  a = 39.10*g/mole; G4Element* elK = new G4Element("Potassium","K", 19,a);
239  a = 24.30*g/mole; G4Element* elMg = new G4Element("Magnesium","Mg",12,a);
240 
241  density = 2.7*g/cm3;
242  G4Material* Rock = new G4Material("Rock", density, 8);
243 
244  //From Daya-Bay
245  Rock->AddElement(elO, 48.50*perCent);
246  Rock->AddElement(elSi, 34.30*perCent);
247  Rock->AddElement(elAl, 8.00*perCent);
248  Rock->AddElement(elFe, 2.00*perCent);
249  Rock->AddElement(elCa, 0.20*perCent);
250  Rock->AddElement(elNa, 2.40*perCent);
251  Rock->AddElement(elK, 4.50*perCent);
252  Rock->AddElement(elMg, 0.10*perCent);
253 
254 
255 
256  // Potential alternative material definitions for HK
257  // Get nist material manager
258  //G4NistManager* nist = G4NistManager::Instance();
259  //nist->FindOrBuildMaterial("G4_AIR");
260  //nist->FindOrBuildMaterial("G4_WATER");
261  //nist->FindOrBuildMaterial("G4_Galactic");
262  //nist->FindOrBuildMaterial("G4_Pyrex_Glass");
263  //nist->FindOrBuildMaterial("G4_POLYETHYLENE");
264 
265 
266 
267 
268 
269 
270 
271 // -------------------------------------------------------------
272 // Generate & Add Material Properties Table
273 // -------------------------------------------------------------
274 
275 
276  const G4int NUMENTRIES = 32;
277 
278  G4double PPCKOV[NUMENTRIES] =
279  { 2.034E-9*GeV, 2.068E-9*GeV, 2.103E-9*GeV, 2.139E-9*GeV,
280  2.177E-9*GeV, 2.216E-9*GeV, 2.256E-9*GeV, 2.298E-9*GeV,
281  2.341E-9*GeV, 2.386E-9*GeV, 2.433E-9*GeV, 2.481E-9*GeV,
282  2.532E-9*GeV, 2.585E-9*GeV, 2.640E-9*GeV, 2.697E-9*GeV,
283  2.757E-9*GeV, 2.820E-9*GeV, 2.885E-9*GeV, 2.954E-9*GeV,
284  3.026E-9*GeV, 3.102E-9*GeV, 3.181E-9*GeV, 3.265E-9*GeV,
285  3.353E-9*GeV, 3.446E-9*GeV, 3.545E-9*GeV, 3.649E-9*GeV,
286  3.760E-9*GeV, 3.877E-9*GeV, 4.002E-9*GeV, 4.136E-9*GeV };
287 
288 
289  // default values
290  /*
291  G4double RINDEX1[NUMENTRIES] =
292  { 1.3435, 1.344, 1.3445, 1.345, 1.3455,
293  1.346, 1.3465, 1.347, 1.3475, 1.348,
294  1.3485, 1.3492, 1.35, 1.3505, 1.351,
295  1.3518, 1.3522, 1.3530, 1.3535, 1.354,
296  1.3545, 1.355, 1.3555, 1.356, 1.3568,
297  1.3572, 1.358, 1.3585, 1.359, 1.3595,
298  1.36, 1.3608};
299  */
300 
301  // M Fechner, values from refsg.F in apdetsim
302  /* G4double RINDEX1[NUMENTRIES] =
303  { 1.33332, 1.333364, 1.333396, 1.3343, 1.33465,
304  1.33502, 1.3354, 1.33579, 1.3362, 1.33663, 1.33709,
305  1.33756, 1.33806, 1.3386, 1.33915, 1.33974,
306  1.34037, 1.34105, 1.34176, 1.34253, 1.34336,
307  1.34425, 1.34521, 1.34626, 1.3474, 1.34864,
308  1.35002, 1.35153, 1.35321, 1.35507, 1.35717, 1.35955 };
309  */
310 
311  //From SFDETSIM water absorption
312  const G4int NUMENTRIES_water=60;
313 
314  G4double ENERGY_water[NUMENTRIES_water] =
315  { 1.56962e-09*GeV, 1.58974e-09*GeV, 1.61039e-09*GeV, 1.63157e-09*GeV,
316  1.65333e-09*GeV, 1.67567e-09*GeV, 1.69863e-09*GeV, 1.72222e-09*GeV,
317  1.74647e-09*GeV, 1.77142e-09*GeV,1.7971e-09*GeV, 1.82352e-09*GeV,
318  1.85074e-09*GeV, 1.87878e-09*GeV, 1.90769e-09*GeV, 1.93749e-09*GeV,
319  1.96825e-09*GeV, 1.99999e-09*GeV, 2.03278e-09*GeV, 2.06666e-09*GeV,
320  2.10169e-09*GeV, 2.13793e-09*GeV, 2.17543e-09*GeV, 2.21428e-09*GeV,
321  2.25454e-09*GeV, 2.29629e-09*GeV, 2.33962e-09*GeV, 2.38461e-09*GeV,
322  2.43137e-09*GeV, 2.47999e-09*GeV, 2.53061e-09*GeV, 2.58333e-09*GeV,
323  2.63829e-09*GeV, 2.69565e-09*GeV, 2.75555e-09*GeV, 2.81817e-09*GeV,
324  2.88371e-09*GeV, 2.95237e-09*GeV, 3.02438e-09*GeV, 3.09999e-09*GeV,
325  3.17948e-09*GeV, 3.26315e-09*GeV, 3.35134e-09*GeV, 3.44444e-09*GeV,
326  3.54285e-09*GeV, 3.64705e-09*GeV, 3.75757e-09*GeV, 3.87499e-09*GeV,
327  3.99999e-09*GeV, 4.13332e-09*GeV, 4.27585e-09*GeV, 4.42856e-09*GeV,
328  4.59258e-09*GeV, 4.76922e-09*GeV, 4.95999e-09*GeV, 5.16665e-09*GeV,
329  5.39129e-09*GeV, 5.63635e-09*GeV, 5.90475e-09*GeV, 6.19998e-09*GeV };
330 
331 
332 
333  // Air
334  G4double RINDEX_air[NUMENTRIES_water] =
335  { 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
336  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
337  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
338  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
339  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
340  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
341  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
342  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
343  1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
344  1.0, 1.0, 1.0, 1.0, 1.0, 1.0};
345 
346 
347  // M Fechner : new ; define the water refraction index using refsg.F
348  //from skdetsim using the whole range.
349  G4double RINDEX1[NUMENTRIES_water] =
350  {1.32885, 1.32906, 1.32927, 1.32948, 1.3297, 1.32992, 1.33014,
351  1.33037, 1.3306, 1.33084, 1.33109, 1.33134, 1.3316, 1.33186, 1.33213,
352  1.33241, 1.3327, 1.33299, 1.33329, 1.33361, 1.33393, 1.33427, 1.33462,
353  1.33498, 1.33536, 1.33576, 1.33617, 1.3366, 1.33705, 1.33753, 1.33803,
354  1.33855, 1.33911, 1.3397, 1.34033, 1.341, 1.34172, 1.34248, 1.34331,
355  1.34419, 1.34515, 1.3462, 1.34733, 1.34858, 1.34994, 1.35145, 1.35312,
356  1.35498, 1.35707, 1.35943, 1.36211, 1.36518, 1.36872, 1.37287, 1.37776,
357  1.38362, 1.39074, 1.39956, 1.41075, 1.42535};
358 
359 
360  G4double ABWFF = 1.0;
361 
362  // Get from the tuning parameters
363  ABWFF = WCSimTuningParams->GetAbwff();
364 
365  //T. Akiri: Values from Skdetsim
366  G4double ABSORPTION_water[NUMENTRIES_water] =
367  {
368  16.1419*cm*ABWFF, 18.278*cm*ABWFF, 21.0657*cm*ABWFF, 24.8568*cm*ABWFF, 30.3117*cm*ABWFF,
369  38.8341*cm*ABWFF, 54.0231*cm*ABWFF, 81.2306*cm*ABWFF, 120.909*cm*ABWFF, 160.238*cm*ABWFF,
370  193.771*cm*ABWFF, 215.017*cm*ABWFF, 227.747*cm*ABWFF, 243.85*cm*ABWFF, 294.036*cm*ABWFF,
371  321.647*cm*ABWFF, 342.81*cm*ABWFF, 362.827*cm*ABWFF, 378.041*cm*ABWFF, 449.378*cm*ABWFF,
372  739.434*cm*ABWFF, 1114.23*cm*ABWFF, 1435.56*cm*ABWFF, 1611.06*cm*ABWFF, 1764.18*cm*ABWFF,
373  2100.95*cm*ABWFF, 2292.9*cm*ABWFF, 2431.33*cm*ABWFF, 3053.6*cm*ABWFF, 4838.23*cm*ABWFF,
374  6539.65*cm*ABWFF, 7682.63*cm*ABWFF, 9137.28*cm*ABWFF, 12220.9*cm*ABWFF, 15270.7*cm*ABWFF,
375  19051.5*cm*ABWFF, 23671.3*cm*ABWFF, 29191.1*cm*ABWFF, 35567.9*cm*ABWFF, 42583*cm*ABWFF,
376  49779.6*cm*ABWFF, 56465.3*cm*ABWFF, 61830*cm*ABWFF, 65174.6*cm*ABWFF, 66143.7*cm*ABWFF,
377  64820*cm*ABWFF, 61635*cm*ABWFF, 57176.2*cm*ABWFF, 52012.1*cm*ABWFF, 46595.7*cm*ABWFF,
378  41242.1*cm*ABWFF, 36146.3*cm*ABWFF, 31415.4*cm*ABWFF, 27097.8*cm*ABWFF, 23205.7*cm*ABWFF,
379  19730.3*cm*ABWFF, 16651.6*cm*ABWFF, 13943.6*cm*ABWFF, 11578.1*cm*ABWFF, 9526.13*cm*ABWFF
380  };
381 
382  //Xin Qian: proposed value for absorption length
383  // G4double ABSORPTION_water[NUMENTRIES_water] =
384 // {22.8154*cm*ABWFF, 28.6144*cm*ABWFF, 35.9923*cm*ABWFF, 45.4086*cm*ABWFF, 57.4650*cm*ABWFF,
385 // 72.9526*cm*ABWFF, 75*cm*ABWFF, 81.2317*cm*ABWFF, 120.901*cm*ABWFF, 160.243*cm*ABWFF,
386 // 193.797*cm*ABWFF, 215.045*cm*ABWFF, 227.786*cm*ABWFF, 243.893*cm*ABWFF, 294.113*cm*ABWFF,
387 // 321.735*cm*ABWFF, 342.931*cm*ABWFF, 362.967*cm*ABWFF, 378.212*cm*ABWFF, 449.602*cm*ABWFF,
388 // 740.143*cm*ABWFF, 1116.06*cm*ABWFF, 1438.78*cm*ABWFF, 1615.48*cm*ABWFF, 1769.86*cm*ABWFF,
389 // 2109.67*cm*ABWFF, 2304.13*cm*ABWFF, 2444.97*cm*ABWFF, 3076.83*cm*ABWFF, 4901.5*cm*ABWFF,
390 // 6666.57*cm*ABWFF, 7873.95*cm*ABWFF, 9433.81*cm*ABWFF, 10214.5*cm*ABWFF, 10845.8*cm*ABWFF,
391 // 15746.9*cm*ABWFF, 20201.8*cm*ABWFF, 22025.8*cm*ABWFF, 21142.2*cm*ABWFF, 15083.9*cm*ABWFF,
392 // 11751*cm*ABWFF, 8795.34*cm*ABWFF, 8741.23*cm*ABWFF, 7102.37*cm*ABWFF, 6060.68*cm*ABWFF,
393 // 4498.56*cm*ABWFF, 3039.56*cm*ABWFF, 2232.2*cm*ABWFF, 1938*cm*ABWFF, 1811.58*cm*ABWFF,
394 // 1610.32*cm*ABWFF, 1338.7*cm*ABWFF, 1095.3*cm*ABWFF, 977.525*cm*ABWFF, 965.258*cm*ABWFF,
395 // 1082.86*cm*ABWFF, 876.434*cm*ABWFF, 633.723*cm*ABWFF, 389.87*cm*ABWFF, 142.011*cm*ABWFF
396 // };
397 
398  /*G4double ABSORPTION_water[NUMENTRIES_water] = //old
399  { 22.8154*cm*ABWFF, 28.6144*cm*ABWFF, 35.9923*cm*ABWFF, 45.4086*cm*ABWFF, 57.4650*cm*ABWFF,
400  72.9526*cm*ABWFF, 92.9156*cm*ABWFF, 118.737*cm*ABWFF, 152.255*cm*ABWFF, 195.925*cm*ABWFF,
401  202.429*cm*ABWFF, 224.719*cm*ABWFF, 236.407*cm*ABWFF, 245.700*cm*ABWFF, 289.017*cm*ABWFF,
402  305.810*cm*ABWFF, 316.456*cm*ABWFF, 326.797*cm*ABWFF, 347.222*cm*ABWFF, 414.938*cm*ABWFF,
403  636.943*cm*ABWFF, 934.579*cm*ABWFF, 1245.33*cm*ABWFF, 1402.52*cm*ABWFF, 1550.39*cm*ABWFF,
404  1745.20*cm*ABWFF, 1883.24*cm*ABWFF, 2016.13*cm*ABWFF, 2442.18*cm*ABWFF, 3831.28*cm*ABWFF,
405  4652.89*cm*ABWFF, 5577.04*cm*ABWFF, 6567.08*cm*ABWFF, 7559.88*cm*ABWFF, 8470.06*cm*ABWFF,
406  9205.54*cm*ABWFF, 9690.95*cm*ABWFF, 9888.36*cm*ABWFF, 9804.50*cm*ABWFF, 9482.17*cm*ABWFF,
407  8982.77*cm*ABWFF, 8369.39*cm*ABWFF, 7680.31*cm*ABWFF, 6902.11*cm*ABWFF, 6183.84*cm*ABWFF,
408  5522.27*cm*ABWFF, 4914.33*cm*ABWFF, 4357.09*cm*ABWFF, 3847.72*cm*ABWFF, 3383.51*cm*ABWFF,
409  2961.81*cm*ABWFF, 2580.08*cm*ABWFF, 2235.83*cm*ABWFF, 1926.66*cm*ABWFF, 1650.21*cm*ABWFF,
410  1404.21*cm*ABWFF, 1186.44*cm*ABWFF, 994.742*cm*ABWFF, 827.027*cm*ABWFF, 681.278*cm*ABWFF};
411  */
412 
413  /*
414  G4double ABSORPTION_water[NUMENTRIES_water] = //new
415  {25.3504*cm, 31.7938*cm, 39.9915*cm, 50.454*cm, 63.85*cm,
416  81.0584*cm, 103.24*cm, 131.93*cm, 169.172*cm, 217.694*cm,
417  224.921*cm, 249.688*cm, 262.674*cm, 273*cm, 321.13*cm, 339.789*cm,
418  351.617*cm, 363.108*cm, 385.802*cm, 461.042*cm, 707.714*cm,
419  1038.42*cm, 1383.7*cm, 1558.36*cm, 1722.65*cm, 1939.11*cm,
420  2092.49*cm, 2240.14*cm, 2962.96*cm, 4967.03*cm, 6368.58*cm,
421  8207.56*cm, 10634.2*cm, 13855.3*cm, 18157.3*cm, 23940.2*cm,
422  31766*cm, 42431.6*cm, 57074.9*cm, 77335.9*cm, 105598*cm,
423  145361*cm, 192434*cm, 183898*cm, 176087*cm, 168913*cm, 162301*cm,
424  156187*cm, 150516*cm, 145243*cm, 140327*cm, 135733*cm, 131430*cm,
425  127392*cm, 123594*cm, 120016*cm, 116640*cm, 113448*cm, 110426*cm,
426  107562*cm};
427  */
428  // M Fechner: Rayleigh scattering -- as of version 4.6.2 of GEANT,
429  // one may use one's own Rayleigh scattering lengths (the buffer is no
430  // longer overwritten for "water", see 4.6.2 release notes)
431 
432  // RAYFF = 1/ARAS, for those of you who know SKdetsim...
433  // actually that's not quite right because the scattering models
434  // are different; in G4 there is no scattering depolarization
435  // std value at SK = 0.6. But Mie scattering is implemented
436  // in SKdetsim and not in G4
437 
438 
439  // april 2005 : reduced reflections, let's increase scattering...
440  // sep 09: for the large detector like superK the old values are muc better
441  //G4double RAYFF = 1.0/1.65; //old
442  // G4double RAYFF = 1.0/1.5;
443 
444  G4double RAYFF = 0.625;
445 
446  // Get from the tuning parameters
447  RAYFF = WCSimTuningParams->GetRayff();
448  // printf("RAYFF: %f\n",RAYFF);
449 
450  //T. Akiri: Values from Skdetsim
451  G4double RAYLEIGH_water[NUMENTRIES_water] = {
452  386929*cm*RAYFF, 366249*cm*RAYFF, 346398*cm*RAYFF, 327355*cm*RAYFF, 309097*cm*RAYFF,
453  291603*cm*RAYFF, 274853*cm*RAYFF, 258825*cm*RAYFF, 243500*cm*RAYFF, 228856*cm*RAYFF,
454  214873*cm*RAYFF, 201533*cm*RAYFF, 188816*cm*RAYFF, 176702*cm*RAYFF, 165173*cm*RAYFF,
455  154210*cm*RAYFF, 143795*cm*RAYFF, 133910*cm*RAYFF, 124537*cm*RAYFF, 115659*cm*RAYFF,
456  107258*cm*RAYFF, 99318.2*cm*RAYFF, 91822.2*cm*RAYFF, 84754*cm*RAYFF, 78097.3*cm*RAYFF,
457  71836.5*cm*RAYFF, 65956*cm*RAYFF, 60440.6*cm*RAYFF, 55275.4*cm*RAYFF, 50445.6*cm*RAYFF,
458  45937*cm*RAYFF, 41735.2*cm*RAYFF, 37826.6*cm*RAYFF, 34197.6*cm*RAYFF, 30834.9*cm*RAYFF,
459  27725.4*cm*RAYFF, 24856.6*cm*RAYFF, 22215.9*cm*RAYFF, 19791.3*cm*RAYFF, 17570.9*cm*RAYFF,
460  15543*cm*RAYFF, 13696.6*cm*RAYFF, 12020.5*cm*RAYFF, 10504.1*cm*RAYFF, 9137.15*cm*RAYFF,
461  7909.45*cm*RAYFF, 6811.3*cm*RAYFF, 5833.25*cm*RAYFF, 4966.2*cm*RAYFF, 4201.36*cm*RAYFF,
462  3530.28*cm*RAYFF, 2944.84*cm*RAYFF, 2437.28*cm*RAYFF, 2000.18*cm*RAYFF, 1626.5*cm*RAYFF,
463  1309.55*cm*RAYFF, 1043.03*cm*RAYFF, 821.016*cm*RAYFF, 637.97*cm*RAYFF, 488.754*cm*RAYFF
464  };
465 
466  /*G4double RAYLEIGH_water[NUMENTRIES_water] = {
467  167024.4*cm*RAYFF, 158726.7*cm*RAYFF, 150742*cm*RAYFF,
468  143062.5*cm*RAYFF, 135680.2*cm*RAYFF, 128587.4*cm*RAYFF,
469  121776.3*cm*RAYFF, 115239.5*cm*RAYFF, 108969.5*cm*RAYFF,
470  102958.8*cm*RAYFF, 97200.35*cm*RAYFF, 91686.86*cm*RAYFF,
471  86411.33*cm*RAYFF, 81366.79*cm*RAYFF, 76546.42*cm*RAYFF,
472  71943.46*cm*RAYFF, 67551.29*cm*RAYFF, 63363.36*cm*RAYFF,
473  59373.25*cm*RAYFF, 55574.61*cm*RAYFF, 51961.24*cm*RAYFF,
474  48527.00*cm*RAYFF, 45265.87*cm*RAYFF, 42171.94*cm*RAYFF,
475  39239.39*cm*RAYFF, 36462.50*cm*RAYFF, 33835.68*cm*RAYFF,
476  31353.41*cm*RAYFF, 29010.30*cm*RAYFF, 26801.03*cm*RAYFF,
477  24720.42*cm*RAYFF, 22763.36*cm*RAYFF, 20924.88*cm*RAYFF,
478  19200.07*cm*RAYFF, 17584.16*cm*RAYFF, 16072.45*cm*RAYFF,
479  14660.38*cm*RAYFF, 13343.46*cm*RAYFF, 12117.33*cm*RAYFF,
480  10977.70*cm*RAYFF, 9920.416*cm*RAYFF, 8941.407*cm*RAYFF,
481  8036.711*cm*RAYFF, 7202.470*cm*RAYFF, 6434.927*cm*RAYFF,
482  5730.429*cm*RAYFF, 5085.425*cm*RAYFF, 4496.467*cm*RAYFF,
483  3960.210*cm*RAYFF, 3473.413*cm*RAYFF, 3032.937*cm*RAYFF,
484  2635.746*cm*RAYFF, 2278.907*cm*RAYFF, 1959.588*cm*RAYFF,
485  1675.064*cm*RAYFF, 1422.710*cm*RAYFF, 1200.004*cm*RAYFF,
486  1004.528*cm*RAYFF, 833.9666*cm*RAYFF, 686.1063*cm*RAYFF
487  };*/
488 
489 
490  // Get from the tuning parameters
491  G4double MIEFF = WCSimTuningParams->GetMieff();
492  //G4double MIEFF = 0.0;
493  // printf("MIEFF: %f\n",MIEFF);
494 
495  //Values extracted from Skdetsim
496  G4double MIE_water[NUMENTRIES_water] = {
497  7790020*cm*MIEFF, 7403010*cm*MIEFF, 7030610*cm*MIEFF, 6672440*cm*MIEFF, 6328120*cm*MIEFF,
498  5997320*cm*MIEFF, 5679650*cm*MIEFF, 5374770*cm*MIEFF, 5082340*cm*MIEFF, 4802000*cm*MIEFF,
499  4533420*cm*MIEFF, 4276280*cm*MIEFF, 4030220*cm*MIEFF, 3794950*cm*MIEFF, 3570120*cm*MIEFF,
500  3355440*cm*MIEFF, 3150590*cm*MIEFF, 2955270*cm*MIEFF, 2769170*cm*MIEFF, 2592000*cm*MIEFF,
501  2423470*cm*MIEFF, 2263300*cm*MIEFF, 2111200*cm*MIEFF, 1966900*cm*MIEFF, 1830120*cm*MIEFF,
502  1700610*cm*MIEFF, 1578100*cm*MIEFF, 1462320*cm*MIEFF, 1353040*cm*MIEFF, 1250000*cm*MIEFF,
503  1152960*cm*MIEFF, 1061680*cm*MIEFF, 975936*cm*MIEFF, 895491*cm*MIEFF, 820125*cm*MIEFF,
504  749619*cm*MIEFF, 683760*cm*MIEFF, 622339*cm*MIEFF, 565152*cm*MIEFF, 512000*cm*MIEFF,
505  462688*cm*MIEFF, 417027*cm*MIEFF, 374832*cm*MIEFF, 335923*cm*MIEFF, 300125*cm*MIEFF,
506  267267*cm*MIEFF, 237184*cm*MIEFF, 209715*cm*MIEFF, 184704*cm*MIEFF, 162000*cm*MIEFF,
507  141456*cm*MIEFF, 122931*cm*MIEFF, 106288*cm*MIEFF, 91395.2*cm*MIEFF, 78125*cm*MIEFF,
508  66355.2*cm*MIEFF, 55968.2*cm*MIEFF, 46851.2*cm*MIEFF, 38896.2*cm*MIEFF, 32000*cm*MIEFF
509  };
510 
511  //Mie scattering length values when assuming 10 times larger than Rayleigh scattering.
512  /*G4double MIE_water[NUMENTRIES_water] = {
513  3869290*cm*MIEFF, 3662490*cm*MIEFF, 3463980*cm*MIEFF, 3273550*cm*MIEFF, 3090970*cm*MIEFF,
514  2916030*cm*MIEFF, 2748530*cm*MIEFF, 2588250*cm*MIEFF, 2435000*cm*MIEFF, 2288560*cm*MIEFF,
515  2148730*cm*MIEFF, 2015330*cm*MIEFF, 1888160*cm*MIEFF, 1767020*cm*MIEFF, 1651730*cm*MIEFF,
516  1542100*cm*MIEFF, 1437950*cm*MIEFF, 1339100*cm*MIEFF, 1245370*cm*MIEFF, 1156590*cm*MIEFF,
517  1072580*cm*MIEFF, 993182*cm*MIEFF, 918222*cm*MIEFF, 847540*cm*MIEFF, 780973*cm*MIEFF,
518  718365*cm*MIEFF, 65956*cm*MIEFF, 604406*cm*MIEFF, 552754*cm*MIEFF, 504456*cm*MIEFF,
519  459370*cm*MIEFF, 417352*cm*MIEFF, 378266*cm*MIEFF, 341976*cm*MIEFF, 308349*cm*MIEFF,
520  277254*cm*MIEFF, 248566*cm*MIEFF, 222159*cm*MIEFF, 197913*cm*MIEFF, 175709*cm*MIEFF,
521  155430*cm*MIEFF, 136966*cm*MIEFF, 120205*cm*MIEFF, 105041*cm*MIEFF, 91371.5*cm*MIEFF,
522  79094.5*cm*MIEFF, 68113*cm*MIEFF, 58332.5*cm*MIEFF, 49662*cm*MIEFF, 42013.6*cm*MIEFF,
523  35302.8*cm*MIEFF, 29448.4*cm*MIEFF, 24372.8*cm*MIEFF, 20001.8*cm*MIEFF, 16265*cm*MIEFF,
524  13095.5*cm*MIEFF, 10430.3*cm*MIEFF, 8210.16*cm*MIEFF, 6379.7*cm*MIEFF, 4887.54*cm*MIEFF
525  };
526  */
527 
528  G4double MIE_water_const[3]={0.4,0.,1};// gforward, gbackward, forward backward ratio
529 
530 
531  //From SFDETSIM
532  /*
533  G4double RINDEX_glass[NUMENTRIES] =
534  { 1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
535  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
536  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
537  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
538  1.600, 1.600, 1.600, 1.600 };
539  */
540  // M Fechner : unphysical, I want to reduce reflections
541  // back to the old value 1.55
542 
543  G4double RINDEX_glass[NUMENTRIES_water] =
544  { 1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
545  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
546  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
547  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
548  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
549  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
550  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
551  1.600, 1.600, 1.600, 1.600, 1.600, 1.600, 1.600,
552  1.600, 1.600 };
553 
554  //G4double RINDEX_blacksheet[NUMENTRIES] =
555  //{ 2.500, 2.500, 2.500, 2.500, 2.500, 2.500, 2.500,
556  // 2.500, 2.500, 2.500, 2.500, 2.500, 2.500, 2.500,
557  // 2.500, 2.500, 2.500, 2.500, 2.500, 2.500, 2.500,
558  // 2.500, 2.500, 2.500, 2.500, 2.500, 2.500, 2.500,
559  // 2.500, 2.500, 2.500, 2.500 };
560 
561 
562  //G4double ABSORPTION1[NUMENTRIES] =
563  //{344.8*cm, 408.2*cm, 632.9*cm, 917.4*cm, 1234.6*cm, 1388.9*cm,
564  // 1515.2*cm, 1724.1*cm, 1886.8*cm, 2000.0*cm, 2631.6*cm, 3571.4*cm,
565  // 4545.5*cm, 4761.9*cm, 5263.2*cm, 5263.2*cm, 5555.6*cm, 5263.2*cm,
566  // 5263.2*cm, 4761.9*cm, 4545.5*cm, 4166.7*cm, 3703.7*cm, 3333.3*cm,
567  // 3000.0*cm, 2850.0*cm, 2700.0*cm, 2450.0*cm, 2200.0*cm, 1950.0*cm,
568  // 1750.0*cm, 1450.0*cm };
569 
570  /*
571  G4double ABSORPTION_glass[NUMENTRIES] =
572  { 100.0*cm, 110.0*cm, 120.0*cm, 130.0*cm, 140.0*cm, 150.0*cm, 160.0*cm,
573  165.0*cm, 170.0*cm, 175.0*cm, 180.0*cm, 185.0*cm, 190.0*cm, 195.0*cm,
574  200.0*cm, 200.0*cm, 200.0*cm, 200.0*cm, 200.0*cm, 195.0*cm, 190.0*cm,
575  185.0*cm, 180.0*cm, 175.0*cm, 170.0*cm, 160.0*cm, 150.0*cm, 140.0*cm,
576  130.0*cm, 120.0*cm, 110.0*cm, 100.0*cm };
577  */
578  // M Fechner : the quantum efficiency already takes glass abs into account
579 
580  G4double ABSORPTION_glass[NUMENTRIES_water]=
581  { 1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
582  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
583  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
584  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
585  1.0e9*cm, 1.0e9*cm,
586  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
587  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
588  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
589  1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,1.0e9*cm,
590  1.0e9*cm, 1.0e9*cm };
591 
592  G4double BLACKABS_blacksheet[NUMENTRIES_water] =
593  { 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
594  1.0e-9*cm,
595  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
596  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
597  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
598  1.0e-9*cm, 1.0e-9*cm,
599  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
600  1.0e-9*cm,
601  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
602  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
603  1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm, 1.0e-9*cm,
604  1.0e-9*cm, 1.0e-9*cm};
605 
606 
607  //The following reflectivity for blacksheet is obtained from skdetsim
608  //There is a fudge factor of 2.7 for the reflectivity of blacksheet
609  //depending on whether SK1 or SK2 simulation is used.
610  //The BlackSheetFudgeFactor is set to true if you want to use the
611  //SK2 values, false if not.
612  //G4double SK1SK2FF = 1.0;
613  //G4bool BlackSheetFudgeFactor=false;
614  //G4bool BlackSheetFudgeFactor=true;
615  // if (BlackSheetFudgeFactor) SK1SK2FF=SK1SK2FF*2.7;
616  //if (BlackSheetFudgeFactor) SK1SK2FF=SK1SK2FF*1.55;
617 
618  //July 1, 2010, F. Beroz: changed SK1SK2FF to BSRFF to avoid confusion,
619  // since this parameter is not from SK.
620 
621  G4double BSRFF = 1.0;
622 
623  // Get from the tuning parameters
624  BSRFF = WCSimTuningParams->GetBsrff();
625 
626  G4double REFLECTIVITY_blacksheet[NUMENTRIES_water] =
627  { 0.055*BSRFF, 0.055*BSRFF,
628  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
629  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
630  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
631  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
632  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
633  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
634  0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF, 0.055*BSRFF,
635  0.055*BSRFF, 0.057*BSRFF, 0.059*BSRFF, 0.060*BSRFF,
636  0.059*BSRFF, 0.058*BSRFF, 0.057*BSRFF, 0.055*BSRFF,
637  0.050*BSRFF, 0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF,
638  0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF,
639  0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF,
640  0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF,
641  0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF, 0.045*BSRFF,
642  0.045*BSRFF, 0.045*BSRFF };
643 
644  /*
645  G4double REFLECTIVITY_blacksheet[NUMENTRIES] =
646  { 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
647  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
648  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
649  0.055*SK1SK2FF, 0.057*SK1SK2FF, 0.059*SK1SK2FF, 0.060*SK1SK2FF,
650  0.059*SK1SK2FF, 0.058*SK1SK2FF, 0.057*SK1SK2FF, 0.055*SK1SK2FF,
651  0.050*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF,
652  0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF,
653  0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF };
654  */
655 
656  /*
657  G4double REFLECTIVITY_blacksheet[NUMENTRIES_water] =
658  { 0.055*SK1SK2FF, 0.055*SK1SK2FF,
659  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
660  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
661  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
662  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
663  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
664  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
665  0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF, 0.055*SK1SK2FF,
666  0.055*SK1SK2FF, 0.057*SK1SK2FF, 0.059*SK1SK2FF, 0.060*SK1SK2FF,
667  0.059*SK1SK2FF, 0.058*SK1SK2FF, 0.057*SK1SK2FF, 0.055*SK1SK2FF,
668  0.050*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF,
669  0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF,
670  0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF,
671  0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF,
672  0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF, 0.045*SK1SK2FF ,
673  0.045*SK1SK2FF, 0.045*SK1SK2FF };
674  */
675 
676  //utter fiction at this stage
677  G4double EFFICIENCY[NUMENTRIES_water] =
678  { 0.001*m };
679 
680  //utter fiction at this stage, does not matter
681  G4double RAYLEIGH_air[NUMENTRIES_water] =
682  { 0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
683  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
684  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
685  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
686  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
687  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
688  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
689  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
690  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,
691  0.001*m,0.001*m,0.001*m,0.001*m,0.001*m,0.001*m};
692 
693  //utter fiction at this stage, does not matter
694  G4double MIE_air[NUMENTRIES_water] =
695  { 0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
696  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
697  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
698  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
699  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
700  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
701  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
702  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
703  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,
704  0.1*m,0.1*m,0.1*m,0.1*m,0.1*m,0.1*m};
705 
706  G4double MIE_air_const[3]={0.99,0.99,0.8};// gforward, gbackward, forward backward ratio
707 
708 
709  //Not used yet, fictional values
710  //G4double SPECULARLOBECONSTANT1[NUMENTRIES] =
711  //{ 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
712  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
713  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
714  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
715  // 0.001, 0.001, 0.001, 0.001 };
716 
717  //Not used yet, fictional values
718  //G4double SPECULARSPIKECONSTANT1[NUMENTRIES] =
719  //{ 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
720  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
721  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
722  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
723  // 0.001, 0.001, 0.001, 0.001 };
724 
725  //Not used yet, fictional values
726  //G4double BACKSCATTERCONSTANT1[NUMENTRIES] =
727  //{ 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
728  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
729  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
730  // 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
731  // 0.001, 0.001, 0.001, 0.001 };
732 
733  G4double EFFICIENCY_blacksheet[NUMENTRIES_water] =
734  { 0.0 };
735 
736  // ------------- Surfaces --------------
737 
739  new G4OpticalSurface("WaterBSCellSurface");
740 
741  OpWaterBSSurface->SetType(dielectric_dielectric);
742  OpWaterBSSurface->SetModel(unified);
743  OpWaterBSSurface->SetFinish(groundfrontpainted);
744  OpWaterBSSurface->SetSigmaAlpha(0.1);
745 
746  const G4int NUM = 2;
747  // G4double PP[NUM] =
748  //{ 2.038E-9*GeV, 4.144E-9*GeV };
749 
750  G4double PP[NUM] = { 1.4E-9*GeV,6.2E-9*GeV};
751  G4double RINDEX_blacksheet[NUM] =
752  { 1.6, 1.6 };
753 
754  G4double SPECULARLOBECONSTANT[NUM] =
755  { 0.3, 0.3 };
756  G4double SPECULARSPIKECONSTANT[NUM] =
757  { 0.2, 0.2 };
758  G4double BACKSCATTERCONSTANT[NUM] =
759  { 0.2, 0.2 };
760 
762  new G4OpticalSurface("GlassCathodeSurface");
763  OpGlassCathodeSurface->SetType(dielectric_dielectric);
764  OpGlassCathodeSurface->SetModel(unified);
765  // OpGlassCathodeSurface->SetFinish(groundbackpainted);
766  OpGlassCathodeSurface->SetFinish(polished);
767  //OpGlassCathodeSurface->SetSigmaAlpha(0.002);
768  // was 1.0
769  // totally unphysical anyway
770  G4double RINDEX_cathode[NUM] =
771  { 1.0, 1.0 };
772  /*
773  G4double SPECULARLOBECONSTANT_glasscath[NUM] =
774  { 1.0, 1.0 };
775  // { 0.3, 0.3 };
776  G4double SPECULARSPIKECONSTANT_glasscath[NUM] =
777  { 0.0, 0.0 };
778  // { 0.2, 0.2 };
779  G4double BACKSCATTERCONSTANT_glasscath[NUM] =
780  {0.0, 0.0};
781  // { 0.2, 0.2 };
782  */
783 
784 
785  G4double RGCFF = 0.0;
786  RGCFF = WCSimTuningParams->GetRgcff();
787 
788  G4double REFLECTIVITY_glasscath[NUM] =
789  //{ 0.0+RGCFF, 0.0+RGCFF };
790  //{ RGCFF, RGCFF };
791  //{ 0.0, 0.0 };
792  { 1.0*RGCFF, 1.0*RGCFF };
793 
794 
795  G4double EFFICIENCY_glasscath[NUM] =
796  { 0.0, 0.0 };
797 
798 
799  // jl145 ----
800  //
802  new G4OpticalSurface("WaterTyCellSurface");
803 
804  OpWaterTySurface->SetType(dielectric_metal); // Only absorption and reflection
805  OpWaterTySurface->SetModel(unified);
806  OpWaterTySurface->SetFinish(ground); // ground surface with tyvek
807  OpWaterTySurface->SetSigmaAlpha(0.2);
808 
809 
810  G4double RINDEX_tyvek[NUM] =
811  { 1.5, 1.5 }; // polyethylene permittivity is ~2.25
812  G4double TySPECULARLOBECONSTANT[NUM] =
813  { 0.75, 0.75 }; // crudely estimated from A. Chavarria's thesis
814  G4double TySPECULARSPIKECONSTANT[NUM] =
815  { 0.0, 0.0 };
816  G4double TyBACKSCATTERCONSTANT[NUM] =
817  { 0.0, 0.0 };
818  // Lambertian prob is therefore 0.25
819 
820 #define NUMENTRIES_TY 33 // Number of bins of wavelength to be used for the Tyvek reflectivity
821 
822  G4double PP_TyREFLECTIVITY[NUMENTRIES_TY] = //Tyvek reflectivity wavelength bins
823  { 2.06642*eV,
824  2.10144*eV, 2.13768*eV, 2.17518*eV, 2.21402*eV, 2.25428*eV,
825  2.29602*eV, 2.33934*eV, 2.38433*eV, 2.43108*eV, 2.4797*eV,
826  2.53031*eV, 2.58302*eV, 2.63798*eV, 2.69533*eV, 2.75523*eV,
827  2.81784*eV, 2.88338*eV, 2.95203*eV, 3.02403*eV, 3.09963*eV,
828  3.17911*eV, 3.26277*eV, 3.35095*eV, 3.44403*eV, 3.54243*eV,
829  3.64662*eV, 3.75713*eV, 3.87454*eV, 3.99952*eV, 4.13284*eV,
830  4.27535*eV, 4.42804*eV};
831 
832  G4double TyREFLECTIVITY[NUMENTRIES_TY] = // Tyvek refelctivity
833  { 0.97,
834  0.97, 0.97, 0.97, 0.97, 0.97,
835  0.97, 0.97, 0.97, 0.97, 0.97,
836  0.97, 0.97, 0.97, 0.97, 0.97,
837  0.97, 0.97, 0.97, 0.97, 0.97,
838  0.96, 0.96, 0.95, 0.95, 0.95,
839  0.94, 0.93, 0.92, 0.91, 0.90,
840  0.89, 0.86};
841  //
842  // ----
843 
844 
845  G4MaterialPropertiesTable *myMPT1 = new G4MaterialPropertiesTable();
846  // M Fechner : new ; wider range for lambda
847  myMPT1->AddProperty("RINDEX", ENERGY_water, RINDEX1, NUMENTRIES_water);
848  myMPT1->AddProperty("ABSLENGTH",ENERGY_water, ABSORPTION_water, NUMENTRIES_water);
849  // M Fechner: new, don't let G4 compute it.
850  myMPT1->AddProperty("RAYLEIGH",ENERGY_water,RAYLEIGH_water,NUMENTRIES_water);
851 
852  // myMPT1->AddProperty("MIEHG",ENERGY_water,MIE_water,NUMENTRIES_water);
853 // myMPT1->AddConstProperty("MIEHG_FORWARD",MIE_water_const[0]);
854 // myMPT1->AddConstProperty("MIEHG_BACKWARD",MIE_water_const[1]);
855 // myMPT1->AddConstProperty("MIEHG_FORWARD_RATIO",MIE_water_const[2]);
856 
857 
858 
859  Water->SetMaterialPropertiesTable(myMPT1);
860  //Gd doped water has the same optical properties as pure water
861  DopedWater->SetMaterialPropertiesTable(myMPT1);
862  // myMPT1->DumpTable();
863 
864  G4MaterialPropertiesTable *myMPT2 = new G4MaterialPropertiesTable();
865  myMPT2->AddProperty("RINDEX", ENERGY_water, RINDEX_air, NUMENTRIES_water);
866  // M Fechner : what is that ?????
867  myMPT2->AddProperty("ABSLENGTH", ENERGY_water, BLACKABS_blacksheet, NUMENTRIES_water);
868  myMPT2->AddProperty("RAYLEIGH",ENERGY_water, RAYLEIGH_air, NUMENTRIES_water);
869 
870  // myMPT2->AddProperty("MIEHG",ENERGY_water,MIE_air,NUMENTRIES_water);
871 // myMPT2->AddConstProperty("MIEHG_FORWARD",MIE_air_const[0]);
872 // myMPT2->AddConstProperty("MIEHG_BACKWARD",MIE_air_const[1]);
873 // myMPT2->AddConstProperty("MIEHG_FORWARD_RATIO",MIE_air_const[2]);
874 
875  Air->SetMaterialPropertiesTable(myMPT2);
876 
877  G4MaterialPropertiesTable *myMPT3 = new G4MaterialPropertiesTable();
878  myMPT3->AddProperty("ABSLENGTH", ENERGY_water, BLACKABS_blacksheet, NUMENTRIES_water);
879  myMPT3->AddProperty("REFLECTIVITY", ENERGY_water, REFLECTIVITY_blacksheet, NUMENTRIES_water);
880  myMPT3->AddProperty("EFFICIENCY", ENERGY_water, EFFICIENCY, NUMENTRIES_water);
881  Plastic->SetMaterialPropertiesTable(myMPT3);
882 
883  G4MaterialPropertiesTable *myMPT4 = new G4MaterialPropertiesTable();
884  myMPT4->AddProperty("ABSLENGTH", ENERGY_water, BLACKABS_blacksheet, NUMENTRIES_water);
885  Blacksheet->SetMaterialPropertiesTable(myMPT4);
886 
887  G4MaterialPropertiesTable *myMPT5 = new G4MaterialPropertiesTable();
888  myMPT5->AddProperty("RINDEX", ENERGY_water, RINDEX_glass, NUMENTRIES_water);
889  myMPT5->AddProperty("ABSLENGTH",ENERGY_water, ABSORPTION_glass, NUMENTRIES_water);
890  Glass->SetMaterialPropertiesTable(myMPT5);
891 
892  // jl145 ----
893  // Abs legnth is same as blacksheet, very small.
894  G4MaterialPropertiesTable *myMPT6 = new G4MaterialPropertiesTable();
895  myMPT6->AddProperty("ABSLENGTH", ENERGY_water, BLACKABS_blacksheet, NUMENTRIES_water);
896  Tyvek->SetMaterialPropertiesTable(myMPT6);
897 
898 
899  // ------------- Surfaces --------------
900 
901  // Blacksheet
902  G4MaterialPropertiesTable *myST1 = new G4MaterialPropertiesTable();
903  myST1->AddProperty("RINDEX", ENERGY_water, RINDEX_blacksheet, NUMENTRIES_water);
904  myST1->AddProperty("SPECULARLOBECONSTANT", PP, SPECULARLOBECONSTANT, NUM);
905  myST1->AddProperty("SPECULARSPIKECONSTANT", PP, SPECULARSPIKECONSTANT, NUM);
906  myST1->AddProperty("BACKSCATTERCONSTANT", PP, BACKSCATTERCONSTANT, NUM);
907  myST1->AddProperty("REFLECTIVITY", ENERGY_water, REFLECTIVITY_blacksheet, NUMENTRIES_water);
908  myST1->AddProperty("EFFICIENCY", ENERGY_water, EFFICIENCY_blacksheet, NUMENTRIES_water);
909  OpWaterBSSurface->SetMaterialPropertiesTable(myST1);
910 
911  //Glass to Cathode surface inside PMTs
912  G4MaterialPropertiesTable *myST2 = new G4MaterialPropertiesTable();
913  myST2->AddProperty("RINDEX", PP, RINDEX_cathode, NUM);
914  // myST2->AddProperty("SPECULARLOBECONSTANT", PP, SPECULARLOBECONSTANT_glasscath, NUM);
915  // myST2->AddProperty("SPECULARSPIKECONSTANT", PP, SPECULARSPIKECONSTANT_glasscath, NUM);
916  // myST2->AddProperty("BACKSCATTERCONSTANT", PP, BACKSCATTERCONSTANT_glasscath, NUM);
917  myST2->AddProperty("REFLECTIVITY", PP, REFLECTIVITY_glasscath, NUM);
918  myST2->AddProperty("EFFICIENCY", PP, EFFICIENCY_glasscath, NUM);
919  //myST2->AddProperty("ABSLENGTH", PP, abslength_paint , NUM);
920  OpGlassCathodeSurface->SetMaterialPropertiesTable(myST2);
921 
922  //Tyvek - jl145
923  G4MaterialPropertiesTable *myST3 = new G4MaterialPropertiesTable();
924  // myST3->AddProperty("RINDEX", PP, RINDEX_tyvek, NUM);
925  // myST3->AddProperty("SPECULARLOBECONSTANT", PP, TySPECULARLOBECONSTANT, NUM);
926  // myST3->AddProperty("SPECULARSPIKECONSTANT", PP, TySPECULARSPIKECONSTANT, NUM);
927  // myST3->AddProperty("BACKSCATTERCONSTANT", PP, TyBACKSCATTERCONSTANT, NUM);
928  myST3->AddProperty("REFLECTIVITY", PP_TyREFLECTIVITY, TyREFLECTIVITY, NUM);
929  myST3->AddProperty("EFFICIENCY", PP, EFFICIENCY_blacksheet, NUM);
930  //use same efficiency as blacksheet, which is 0
931  OpWaterTySurface->SetMaterialPropertiesTable(myST3);
932 
933 }
WCSimTuningParameters * WCSimTuningParams
#define NUMENTRIES_TY