G4TankSimulatorOG/G4TankConstruction.cc
Go to the documentation of this file.
1 #include "G4TankConstruction.h"
2 
3 #include <fwk/RunController.h>
4 #include <fwk/CentralConfig.h>
5 
6 #include <sdet/SDetector.h>
7 
8 #include <utl/CoordinateSystem.h>
9 #include <utl/Point.h>
10 #include <utl/TabulatedFunction.h>
11 #include <utl/Reader.h>
12 #include <utl/AugerUnits.h>
13 
14 #include <vector>
15 //#include <iostream>
16 #include <sstream>
17 
18 using sdet::SDetector;
20 using utl::Point;
22 using utl::Branch;
23 
24 using namespace fwk;
25 
26 using std::vector;
27 //using std::cout;
28 //using std::cerr;
29 //using std::endl;
30 
31 using std::ostringstream;
32 
33 using namespace G4TankSimulatorOG;
34 
35 // declare variables that static members above refer to.
36 double G4TankConstruction::fExpHall_x;
37 double G4TankConstruction::fExpHall_y;
38 double G4TankConstruction::fExpHall_z;
39 double G4TankConstruction::fGroundThickness;
40 
41 double G4TankConstruction::fTankRadius;
42 double G4TankConstruction::fTankHalfHeight;
43 double G4TankConstruction::fTankThickness;
44 double G4TankConstruction::fTankPos_x;
45 double G4TankConstruction::fTankPos_y;
46 double G4TankConstruction::fTankPos_z;
47 
48 utl::Point G4TankConstruction::fPmt1;
49 utl::Point G4TankConstruction::fPmt2;
50 utl::Point G4TankConstruction::fPmt3;
51 
52 double G4TankConstruction::fFaceRadius, G4TankConstruction::fFaceRadiusz;
53 double G4TankConstruction::fFaceActiveRadius, G4TankConstruction::fHeightz;
54 double G4TankConstruction::fPmtRmin, G4TankConstruction::fPmtRmax;
55 double G4TankConstruction::fPmtRzmin, G4TankConstruction::fPmtRzmax;
56 double G4TankConstruction::fInterfaceRmin, G4TankConstruction::fInterfaceRzmin;
57 double G4TankConstruction::fInterfaceRmax, G4TankConstruction::fInterfaceRzmax;
58 double G4TankConstruction::fInterfaceThickness, G4TankConstruction::fGlassThickness;
59 double G4TankConstruction::fDomeThickness;
60 double G4TankConstruction::fDomeRmin, G4TankConstruction::fDomeRmax;
61 double G4TankConstruction::fDomeRzmin, G4TankConstruction::fDomeRzmax;
62 double G4TankConstruction::fMinPhi, G4TankConstruction::fMaxPhi;
63 double G4TankConstruction::fMinTheta, G4TankConstruction::fMaxTheta;
64 double G4TankConstruction::fSolarPanelLength, G4TankConstruction::fSolarPanelWidth;
65 double G4TankConstruction::fSolarPanelThickness;
66 double G4TankConstruction::fSolarPanelX, G4TankConstruction::fSolarPanelY;
67 double G4TankConstruction::fSolarPanelZ, G4TankConstruction::fSolarPanelTiltAngle;
68 double G4TankConstruction::fElecBoxLength, G4TankConstruction::fElecBoxWidth;
69 double G4TankConstruction::fElecBoxThickness;
70 double G4TankConstruction::fElecBoxX, G4TankConstruction::fElecBoxY;
71 double G4TankConstruction::fElecBoxZ, G4TankConstruction::fElecBoxTiltAngle;
72 
73 TabulatedFunction G4TankConstruction::fPmtfaceRINDEX;
74 TabulatedFunction G4TankConstruction::fInterfaceRINDEX;
75 TabulatedFunction G4TankConstruction::fPmtdomeRINDEX;
76 TabulatedFunction G4TankConstruction::fPmtfaceABSORPTION;
77 TabulatedFunction G4TankConstruction::fInterfaceABSORPTION;
78 TabulatedFunction G4TankConstruction::fPmtdomeABSORPTION;
79 double G4TankConstruction::fSIGMA_ALPHA;
80 TabulatedFunction G4TankConstruction::fLinerREFLECTIVITY;
81 TabulatedFunction G4TankConstruction::fLinerSPECULARLOBECONSTANT;
82 TabulatedFunction G4TankConstruction::fLinerSPECULARSPIKECONSTANT;
83 TabulatedFunction G4TankConstruction::fLinerBACKSCATTERCONSTANT;
84 TabulatedFunction G4TankConstruction::fLinerTYVEK_RINDEX;
85 TabulatedFunction G4TankConstruction::fLinerABSORPTION;
86 TabulatedFunction G4TankConstruction::fWaterABSORPTION;
87 TabulatedFunction G4TankConstruction::fWaterRINDEX;
88 
89 // Older versions of Geant4 do not have the G4Ellipsoid defined
90 
91 G4TankConstruction::G4TankConstruction() : expHall_phys(0)
92 {
96 
97  // disable solar panel, electronics box, and ground for now:
98 
99  fSolarPanelEnable = false;
100  fElecBoxEnable = false;
101  fGroundEnable = false;
102 
103  Water = NULL;
104  Pyrex = NULL;
105  Pyrex1 = NULL;
106  Lucite = NULL;
107  Interface = NULL;
108  HDPE = NULL;
109  linerOpticalMPT = NULL;
110  return;
111 }
112 
114 
116 {
117  // The "experimental hall" in which the detector
118  // is contained:
119 
120  fExpHall_x = 5.0*m;
121  fExpHall_y = 5.0*m;
122  fExpHall_z = 5.0*m;
123 
124  // The optionally-simulated ground beneath the tank
125  fGroundThickness = 1.25*m;
126 
127  // The position of the center of the tank
128  // in the local Geant4 coordinate system
129 
130  fTankPos_x = 0.0*m;
131  fTankPos_y = 0.0*m;
132  fTankPos_z = 0.0*m;
133 }
134 
136 {
137  vector<double>::iterator tfIt;
138  G4double alpha;
139  // Get the pointer to the current detector
140  G4TankSimulator* const theTankSimulator =
141  dynamic_cast<G4TankSimulator*>(&RunController::GetInstance().GetModule("G4TankSimulatorOG"));
142 
143  const sdet::Station* theCurrentDetectorStation =
144  theTankSimulator->GetCurrentDetectorStation();
145  if (!theCurrentDetectorStation) return; // no point getting data yet
146 
147  // Set tank properties
148 
149  fTankRadius =
150  theCurrentDetectorStation->GetRadius() / utl::m * CLHEP::m;
151  fTankHalfHeight =
152  0.5 * theCurrentDetectorStation->GetHeight() / utl::m * CLHEP::m;
153  fTankThickness =
154  theCurrentDetectorStation->GetThickness() / utl::m * CLHEP::m;
155 
156  // Set pmt position (relative to the center of the tank)
157 
158  fPmt1 = theCurrentDetectorStation->GetPMT(1).GetPosition(); // converted to G4 units later
159  fPmt2 = theCurrentDetectorStation->GetPMT(2).GetPosition(); // converted to G4 units later
160  fPmt3 = theCurrentDetectorStation->GetPMT(3).GetPosition(); // converted to G4 units later
161 
162  // Set pmt optical properties
163  // For now, assume that all three pmts have the same properties
164 
165  fPmtfaceRINDEX =
166  theCurrentDetectorStation->GetPMT(1).GetFaceRefractionIndex();
167  for (tfIt = fPmtfaceRINDEX.XBegin(); tfIt != fPmtfaceRINDEX.XEnd(); tfIt++)
168  (*tfIt) *= eV;
170  theCurrentDetectorStation->GetPMT(1).GetRTVRefractionIndex();
171  for (tfIt = fInterfaceRINDEX.XBegin(); tfIt != fInterfaceRINDEX.XEnd();
172  tfIt++) (*tfIt) *= eV;
173  fPmtdomeRINDEX =
174  theCurrentDetectorStation->GetPMT(1).GetDomeRefractionIndex();
175  for (tfIt = fPmtdomeRINDEX.XBegin(); tfIt != fPmtdomeRINDEX.XEnd(); tfIt++)
176  (*tfIt) *= eV;
177 
178  // These parameters specify the pmt geometry ( already converted to G4 units )
179 
190  fDomeRmax = fFaceRadius + fInterfaceThickness + fDomeThickness;
191  fDomeRzmax = fFaceRadiusz + fInterfaceThickness + fDomeThickness;
192 
193  // Taking into account that the radius of the photocathode is fFaceActiveRadius
194  if ( fFaceActiveRadius/mm >= 114. )
195  {
196  fFaceActiveRadius = 113.926112814*mm;
197  fHeightz = 0.;
198  }
199  else if ( fFaceActiveRadius/mm > 97.27 )
200  {
201  alpha = acos(1./62.*(fFaceActiveRadius/mm-133.*sin(47.*deg)+62.*cos(43.*deg)));
202  fHeightz = 62.*sin(alpha)*mm;
203  }
204  else
205  {
206  alpha = asin(fFaceActiveRadius/mm/133.);
207  fHeightz = (133.*cos(alpha)-(133.-62.)*cos(47.*deg))*mm;
208  }
209 
210  // Set liner properties
211 
212  fSIGMA_ALPHA = theCurrentDetectorStation->GetLinerSigmaAlpha();
213 
214  fLinerREFLECTIVITY = theCurrentDetectorStation->GetLinerReflectivity();
215  for (tfIt = fLinerREFLECTIVITY.XBegin(); tfIt != fLinerREFLECTIVITY.XEnd();
216  tfIt++) (*tfIt) *= eV;
217  fLinerSPECULARLOBECONSTANT = theCurrentDetectorStation->GetLinerSpecularLobe();
218  for (tfIt = fLinerSPECULARLOBECONSTANT.XBegin();
219  tfIt != fLinerSPECULARLOBECONSTANT.XEnd(); tfIt++) (*tfIt) *= eV;
220  fLinerSPECULARSPIKECONSTANT = theCurrentDetectorStation->GetLinerSpecularSpike();
221  for (tfIt = fLinerSPECULARSPIKECONSTANT.XBegin();
222  tfIt != fLinerSPECULARSPIKECONSTANT.XEnd(); tfIt++) (*tfIt) *= eV;
223 
224  // Set water properties
225 
226  fWaterABSORPTION = theCurrentDetectorStation->GetWaterAbsorptionLength();
227  for (tfIt = fWaterABSORPTION.XBegin(); tfIt != fWaterABSORPTION.XEnd();
228  tfIt++) (*tfIt) *= eV;
229  for (tfIt = fWaterABSORPTION.YBegin(); tfIt != fWaterABSORPTION.YEnd();
230  tfIt++) (*tfIt) *= m;
231  fWaterRINDEX = theCurrentDetectorStation->GetWaterRefractionIndex();
232  for (tfIt = fWaterRINDEX.XBegin(); tfIt != fWaterRINDEX.XEnd(); tfIt++)
233  (*tfIt) *= eV;
234  return;
235 }
236 
238 {
239  CentralConfig* theConfig = CentralConfig::GetInstance();
240 
241  Branch TopBranch = theConfig->GetTopBranch("G4TankSimulator");
242 
243  //cout<<"Getting parameters from XML file..."<<endl;
244 
245  unsigned int size;
246 
247  vector<double> pmt_pp;
248  vector<double> pmt_face_abs;
249  vector<double> pmt_interface_abs;
250  vector<double> pmt_dome_abs;
251 
252  TopBranch.GetChild("pmt_momentum_bins").GetData(pmt_pp);
253  TopBranch.GetChild("pmtface_AbsLength").GetData(pmt_face_abs);
254  TopBranch.GetChild("interface_AbsLength").GetData(pmt_interface_abs);
255  TopBranch.GetChild("dome_AbsLength").GetData(pmt_dome_abs);
256  TopBranch.GetChild("interface_thickness").GetData(fInterfaceThickness);
257  TopBranch.GetChild("pmtface_thickness").GetData(fGlassThickness);
260 
261  size = pmt_pp.size();
262 
263  if ( pmt_face_abs.size() != size ||
264  pmt_interface_abs.size() != size ||
265  pmt_dome_abs.size() != size )
266  {
267  ostringstream errMsg;
268  errMsg<<"TabulatedFunction size mismatch for pmt properties";
269  ERROR(errMsg);
270 
271  throw utl::OutOfBoundException(errMsg.str());
272  }
273 
274  for (unsigned int i = 0; i < size; ++i)
275  {
276  fPmtfaceABSORPTION.PushBack( pmt_pp[i]*eV, pmt_face_abs[i]*m );
277  fInterfaceABSORPTION.PushBack( pmt_pp[i]*eV, pmt_interface_abs[i]*m );
278  fPmtdomeABSORPTION.PushBack( pmt_pp[i]*eV, pmt_dome_abs[i]*m );
279  }
280 
281  TopBranch.GetChild("face_radius").GetData(fFaceRadius);
282  TopBranch.GetChild("face_radius_z").GetData(fFaceRadiusz);
283  TopBranch.GetChild("face_active_radius").GetData(fFaceActiveRadius);
284  fFaceRadius*=m;
285  fFaceRadiusz*=m;
287 
288  vector<double> liner_pp;
289  vector<double> liner_rindex;
290  vector<double> liner_abs;
291  vector<double> liner_backscatter;
292 
293  TopBranch.GetChild("liner_momentum_bins").GetData(liner_pp);
294  TopBranch.GetChild("liner_Rindex").GetData(liner_rindex);
295  TopBranch.GetChild("liner_AbsLength").GetData(liner_abs);
296  TopBranch.GetChild("backScatter").GetData(liner_backscatter);
297  TopBranch.GetChild("dome_thickness").GetData(fDomeThickness);
298  fDomeThickness *= m;
299 
300  size = liner_pp.size();
301 
302  if ( liner_rindex.size() != size ||
303  liner_abs.size() != size ||
304  liner_backscatter.size() != size )
305  {
306  ostringstream errMsg;
307  errMsg<<"Tabulated function size mismatch for liner properties";
308  ERROR(errMsg);
309 
310  throw utl::OutOfBoundException(errMsg.str());
311  }
312 
313  for (unsigned int i = 0; i < size; ++i)
314  {
315  fLinerTYVEK_RINDEX.PushBack( liner_pp[i]*eV, liner_rindex[i] );
316  fLinerABSORPTION.PushBack( liner_pp[i]*eV, liner_abs[i]*m );
317  fLinerBACKSCATTERCONSTANT.PushBack( liner_pp[i]*eV, liner_backscatter[i] );
318  }
319  return;
320 }
321 
322 G4VPhysicalVolume *G4TankConstruction::Construct(void)
323 {
324  if(!expHall_phys){
325  CreateElements();
326  CreateMaterials();
327  return CreateTank();
328  }else{
329  return expHall_phys;
330  }
331 }
333 {
334  G4double a;
335  G4String name, symbol;
336  G4int z;
337 
338  // Elements
339  //---------
340 
341  a = 14.01*g/mole;
342  elN = new G4Element(name = "Nitrogen", symbol = "N", z = 7, a);
343  a = 16.00*g/mole;
344  elO = new G4Element(name = "Oxygen", symbol = "O", z = 8, a);
345  a = 29.98*g/mole;
346  elAl = new G4Element(name = "Aluminum", symbol = "Al", z = 13, a);
347  a = 55.85*g/mole;
348  elFe = new G4Element(name = "Iron", symbol = "Fe", z = 26, a);
349  a = 1.01*g/mole;
350  elH = new G4Element(name = "Hydrogen", symbol = "H", z = 1, a);
351  a = 28.09*g/mole;
352  elSi = new G4Element(name = "Silicon", symbol = "Si", z = 14, a);
353  a = 10.811*g/mole;
354  elB = new G4Element(name = "Boron", symbol = "B", z = 5, a);
355  a = 22.98977*g/mole;
356  elNa = new G4Element(name = "Sodium", symbol = "Na", z = 11, a);
357  a = 12.0107*g/mole;
358  elC = new G4Element(name = "Carbon", symbol = "C", z = 6, a);
359 }
360 
362 {
363  CreateAir();
364  CreateWater();
365  CreateVacuum();
366  CreatePyrex();
367  CreatePyrex1();
368  CreateLucite();
369  CreateInterface();
370  CreateHDPE();
371  CreateLiner();
372  CreateAluminium();
373  CreateDirt();
374  return;
375 }
376 
377 G4VPhysicalVolume *G4TankConstruction::CreateTank(void)
378 {
380  CreateHall();
381  AssembleTank();
382  return expHall_phys;
383 }
384 
386 {
387  // Air
388  //----
389  G4String name;
390  G4double density = 1.29e-03*g/cm3;
391  G4int nel;
392  Air = new G4Material(name = "Air_G4Tank", density, nel = 2);
393  Air->AddElement(elN, 0.7);
394  Air->AddElement(elO, 0.3);
395 
396  // Add material property table for Air to simulate empty tank
397  double airPP[2] = { 2.08*eV, 4.20*eV };
398  double airRINDEX[2] = { 1.000273, 1.000273 };
399  double airABSLENGTH[2] = { 10000.0*m, 10000.0*m };
400 
401  airMPT = new G4MaterialPropertiesTable();
402  airMPT->AddProperty("RINDEX", airPP, airRINDEX, 2);
403  airMPT->AddProperty("ABSLENGTH", airPP, airABSLENGTH, 2);
404 
405  Air->SetMaterialPropertiesTable(airMPT);
406  return;
407 }
408 
410 {
411 
412  // Water
413  //------
414  G4String name;
415  G4double density = 1.0*g/cm3;
416  G4int nel;
417 
418  if (Water != NULL) {
419  delete Water;
420  delete waterMPT;
421  }
422  Water = new G4Material(name = "Water", density, nel = 2);
423  Water->AddElement(elH, 2);
424  Water->AddElement(elO, 1);
425 
426  // Set water material properties
427 
428  waterMPT = new G4MaterialPropertiesTable();
429 
430  // water rindex
431 
432  int index;
433  int num = fWaterRINDEX.GetNPoints();
434 
435  double* array_ppWATER = new double[num+1];
436  double* array_waterRINDEX = new double[num+1];
437 
438  for ( index = 0; index < num; index++ )
439  {
440  array_ppWATER[index] = fWaterRINDEX[index].X();
441  array_waterRINDEX[index] = fWaterRINDEX[index].Y();
442  }
443 
444  waterMPT->AddProperty("RINDEX", array_ppWATER, array_waterRINDEX, num);
445 
446  // water abs. length
447 
449 
450  double* array_ppwaterABSORPTION = new double[num+1];
451  double* array_waterABSORPTION = new double[num+1];
452 
453  for ( index = 0; index < num; index++ )
454  {
455  array_ppwaterABSORPTION[index] = fWaterABSORPTION[index].X();
456  array_waterABSORPTION[index] = fWaterABSORPTION[index].Y();
457  }
458 
459  waterMPT->AddProperty("ABSLENGTH", array_ppwaterABSORPTION, array_waterABSORPTION, num);
460 
461  //cout<<"Water Material Properties:"<<endl;
462  //waterMPT->DumpTable();
463 
464  Water->SetMaterialPropertiesTable(waterMPT);
465 
466  delete[] array_ppWATER;
467  delete[] array_waterRINDEX;
468  delete[] array_ppwaterABSORPTION;
469  delete[] array_waterABSORPTION;
470  return;
471 }
472 
474 {
475  // Vacuum
476  //-------
477  G4String name;
478  G4double a = 1.* g / mole;
479  G4double density = universe_mean_density;
480  G4int z;
481 
482  // ATTENTION: change from STP to avoid the warning...
483  Vacuum = new G4Material(name = "Vacuum", z = 1, a, density,
484  kStateUndefined, STP_Temperature,
485  STP_Pressure);
486  return;
487 }
488 
490 {
491  // Borosilicate glass (Pyrex) for the active pmt faces
492  //----------------------------------------------------
493  /*
494  NOTE : The "face" is actually an ellipsoidal surface
495  The active area (photocathode) covers 90% of it.
496  The refractive index is the same as the dome
497  refractive index. The absorption length is set
498  so that no photons propagate through the
499  photocathode.
500  */
501 
502  // from PDG:
503  // 80% SiO2 + 13% B2O2 + 7% Na2O
504  // by fractional mass?
505 
506  G4String name;
507  G4double density = 2.23*g/cm3;
508  G4int nel, natoms, ncomponents;
509 
510  if (Pyrex != NULL) {
511  delete SiO2;
512  delete B2O2;
513  delete Na2O;
514  delete Pyrex;
515  delete pmtfaceMPT;
516  }
517 
518  SiO2 = new G4Material(name = "SiO2", density, nel = 2);
519  SiO2->AddElement(elSi, natoms = 1);
520  SiO2->AddElement(elO, natoms = 2);
521 
522  B2O2 = new G4Material(name = "B2O2", density, nel = 2);
523  B2O2->AddElement(elB, natoms = 2);
524  B2O2->AddElement(elO, natoms = 2);
525 
526  Na2O = new G4Material(name = "Na2O", density, nel = 2);
527  Na2O->AddElement(elNa, natoms = 2);
528  Na2O->AddElement(elO, natoms = 1);
529 
530  Pyrex = new G4Material(name = "Pyrex", density, ncomponents = 3);
531  Pyrex->AddMaterial(SiO2, 0.80);
532  Pyrex->AddMaterial(B2O2, 0.13);
533  Pyrex->AddMaterial(Na2O, 0.07);
534 
535  // Fill material properties table for pmt face
536 
537  pmtfaceMPT = new G4MaterialPropertiesTable();
538 
539  // pmt face rindex
540 
541  int index;
542  int num = fPmtfaceRINDEX.GetNPoints();
543 
544  double* array_pmtfaceRINDEX_PP = new double[num+1];
545  double* array_pmtfaceRINDEX = new double[num+1];
546 
547  // debug
548  //fout<<"PMT face rindex:"<<endl;
549 
550  for ( index = 0; index < num; index++ )
551  {
552  array_pmtfaceRINDEX_PP[index] = fPmtfaceRINDEX[index].X();
553  array_pmtfaceRINDEX[index] = fPmtfaceRINDEX[index].Y();
554  }
555 
556  pmtfaceMPT->AddProperty("RINDEX", array_pmtfaceRINDEX_PP, array_pmtfaceRINDEX, num);
557 
558  // Fill and add pmt face absorption length
559 
561 
562  double* array_pmtfaceABSORPTION_PP = new double[num+1];
563  double* array_pmtfaceABSORPTION = new double[num+1];
564 
565  // debug
566  //fout<<"PMT face abs. length:"<<endl;
567 
568  for ( index = 0; index < num; index++ )
569  {
570  // debug
571  //fout<< fPmtfaceABSORPTION[index].X() <<" "<< fPmtfaceABSORPTION[index].Y() <<endl;
572 
573  array_pmtfaceABSORPTION_PP[index] = fPmtfaceABSORPTION[index].X();
574  array_pmtfaceABSORPTION[index] = fPmtfaceABSORPTION[index].Y();
575  }
576 
577  pmtfaceMPT->AddProperty("ABSLENGTH", array_pmtfaceABSORPTION_PP, array_pmtfaceABSORPTION, num);
578 
579 
580  //cout<<"PMT Face Material Properties:"<<endl;
581  //pmtfaceMPT->DumpTable();
582 
583  Pyrex->SetMaterialPropertiesTable(pmtfaceMPT);
584 
585  delete[] array_pmtfaceRINDEX_PP;
586  delete[] array_pmtfaceRINDEX;
587  delete[] array_pmtfaceABSORPTION_PP;
588  delete[] array_pmtfaceABSORPTION;
589  return;
590 }
591 
593 {
594  // Borosilicate glass (Pyrex1) for the inactive part of the pmt faces
595  //-------------------------------------------------------------------
596  /*
597  NOTE : The inactive part of the "face" is about 10% of
598  the total ellipsoidal surface seen by the photons.
599  The refractive index is the same as the dome
600  refractive index. The absorption length is set
601  so that photons can propagate through this part.
602  To simplify things we use the same absorption length
603  as for the dome.
604  */
605 
606  G4String name;
607  G4double density = 2.23 * g / cm3;
608  G4int nel, natoms;
609 
610  if (Pyrex1 != NULL) {
611  delete SiO2;
612  delete B2O2;
613  delete Na2O;
614  delete Pyrex1;
615  delete pmtfaceMPT1;
616  }
617  SiO2 = new G4Material(name = "SiO2", density, nel = 2);
618  SiO2->AddElement(elSi, natoms = 1);
619  SiO2->AddElement(elO, natoms = 2);
620 
621  B2O2 = new G4Material(name = "B2O2", density, nel = 2);
622  B2O2->AddElement(elB, natoms = 2);
623  B2O2->AddElement(elO, natoms = 2);
624 
625  Na2O = new G4Material(name = "Na2O", density, nel = 2);
626  Na2O->AddElement(elNa, natoms = 2);
627  Na2O->AddElement(elO, natoms = 1);
628 
629  int ncomponents;
630  Pyrex1 = new G4Material(name = "Pyrex1", density, ncomponents = 3);
631  Pyrex1->AddMaterial(SiO2, 0.80);
632  Pyrex1->AddMaterial(B2O2, 0.13);
633  Pyrex1->AddMaterial(Na2O, 0.07);
634 
635  pmtfaceMPT1 = new G4MaterialPropertiesTable();
636 
637  // Fill material properties table for non sensitive pmt face
638 
639  // non-sensitive pmt face rindex = sensitive pmt face rindex
640 
641  int index;
642  int num = fPmtfaceRINDEX.GetNPoints();
643 
644  double* array_pmtface1RINDEX_PP = new double[num+1];
645  double* array_pmtface1RINDEX = new double[num+1];
646 
647  for ( index = 0; index < num; index++ )
648  {
649  array_pmtface1RINDEX_PP[index] = fPmtfaceRINDEX[index].X();
650  array_pmtface1RINDEX[index] = fPmtfaceRINDEX[index].Y();
651  }
652 
653  pmtfaceMPT1->AddProperty("RINDEX", array_pmtface1RINDEX_PP, array_pmtface1RINDEX, num);
654 
655 
656  // Fill and add non-sensitive pmt face absorption length = dome absorption length
657 
659 
660  double* array_pmtface1ABSORPTION_PP = new double[num+1];
661  double* array_pmtface1ABSORPTION = new double[num+1];
662 
663  for ( index = 0; index < num; index++ )
664  {
665  array_pmtface1ABSORPTION_PP[index] = fPmtdomeABSORPTION[index].X();
666  array_pmtface1ABSORPTION[index] = fPmtdomeABSORPTION[index].Y();
667  }
668 
669  pmtfaceMPT1->AddProperty("ABSLENGTH", array_pmtface1ABSORPTION_PP, array_pmtface1ABSORPTION, num);
670 
671  Pyrex1->SetMaterialPropertiesTable(pmtfaceMPT1);
672 
673  delete[] array_pmtface1RINDEX_PP;
674  delete[] array_pmtface1RINDEX;
675  delete[] array_pmtface1ABSORPTION_PP;
676  delete[] array_pmtface1ABSORPTION;
677  return;
678 }
679 
681 {
682  // Lucite (for the pmt domes)
683  G4double density = 1.20*g/cm3;
684  G4String name;
685  G4int nel, natoms;
686 
687  if (Lucite != NULL) {
688  delete CH3;
689  delete CH2;
690  delete C;
691  delete CO2;
692  delete Lucite;
693  delete pmtdomeMPT;
694  }
695 
696  CH3 = new G4Material(name = "CH3", density, nel = 2);
697  CH3->AddElement(elC, natoms = 1);
698  CH3->AddElement(elH, natoms = 3);
699 
700  CH2 = new G4Material(name = "CH2", density, nel = 2);
701  CH2->AddElement(elC, natoms = 1);
702  CH2->AddElement(elH, natoms = 2);
703 
704  C = new G4Material(name = "C", density, nel = 1);
705  C->AddElement(elC, natoms = 1);
706 
707  CO2 = new G4Material(name = "CO2", density, nel = 2);
708  CO2->AddElement(elC, natoms = 1);
709  CO2->AddElement(elO, natoms = 2);
710 
711  /*
712  NOTE : Do not know the fractional mass comp.
713  of lucite, only the molecular composition.
714  Perhaps it is a standard material defined in g4.
715  for now define by number of atoms.
716  It really doesn't matter at this point
717  anyway, since we are more concerned about
718  the optical properties.
719  */
720 
721  Lucite = new G4Material(name="Lucite", density, nel=3);
722  Lucite->AddElement(elC, natoms=4);
723  Lucite->AddElement(elH, natoms=8);
724  Lucite->AddElement(elO, natoms=2);
725 
726  pmtdomeMPT = new G4MaterialPropertiesTable();
727 
728  // pmt dome rindex
729 
730  int index;
731  int num = fPmtdomeRINDEX.GetNPoints();
732 
733  double *array_pmtdomeRINDEX_PP = new double[num+1];
734  double *array_pmtdomeRINDEX = new double[num+1];
735 
736  for ( index = 0; index < num; index++ )
737  {
738  array_pmtdomeRINDEX_PP[index] = fPmtdomeRINDEX[index].X();
739  array_pmtdomeRINDEX[index] = fPmtdomeRINDEX[index].Y();
740  }
741 
742  pmtdomeMPT->AddProperty("RINDEX", array_pmtdomeRINDEX_PP, array_pmtdomeRINDEX, num);
743 
744  // pmt dome absorption length
745 
747 
748  double *array_pmtdomeABSORPTION_PP = new double[num+1];
749  double *array_pmtdomeABSORPTION = new double[num+1];
750 
751  for ( index = 0; index < num; index++ )
752  {
753  array_pmtdomeABSORPTION_PP[index] = fPmtdomeABSORPTION[index].X();
754  array_pmtdomeABSORPTION[index] = fPmtdomeABSORPTION[index].Y();
755  }
756 
757  pmtdomeMPT->AddProperty("ABSLENGTH", array_pmtdomeABSORPTION_PP, array_pmtdomeABSORPTION, num);
758 
759  //cout<<"PMT Dome Material Properties:"<<endl;
760  //pmtdomeMPT->DumpTable();
761 
762  Lucite->SetMaterialPropertiesTable(pmtdomeMPT);
763 
764  delete[] array_pmtdomeRINDEX_PP;
765  delete[] array_pmtdomeRINDEX;
766  delete[] array_pmtdomeABSORPTION_PP;
767  delete[] array_pmtdomeABSORPTION;
768  return;
769 }
770 
772 {
773  // Interface for production tanks is Wacker SilGel 612 (Registered trademark
774  // of Wacker Group
775  G4String name;
776  G4int nel, natoms;
777  G4double density = 0.97*g/cm3; // guess - not explicitly set by Troy
778 
779  if (Interface != NULL) {
780  delete Interface;
781  delete interfaceMPT;
782  }
783  Interface = new G4Material(name = "Interface", density, nel = 3);
784  Interface->AddElement(elC, natoms = 4);
785  Interface->AddElement(elH, natoms = 8);
786  Interface->AddElement(elO, natoms = 2);
787 
788  interfaceMPT = new G4MaterialPropertiesTable();
789 
790  // interface rindex
791 
792  int index;
793  int num = fInterfaceRINDEX.GetNPoints();
794 
795  double* array_interfaceRINDEX_PP = new double[num+1];
796  double* array_interfaceRINDEX = new double[num+1];
797 
798  for ( index = 0; index < num; index++ )
799  {
800  array_interfaceRINDEX_PP[index] = fInterfaceRINDEX[index].X();
801  array_interfaceRINDEX[index] = fInterfaceRINDEX[index].Y();
802  }
803 
804  interfaceMPT->AddProperty("RINDEX", array_interfaceRINDEX_PP, array_interfaceRINDEX, num);
805 
806  // interface abs. length
807 
809 
810  double* array_interfaceABSORPTION_PP = new double[num+1];
811  double* array_interfaceABSORPTION = new double[num+1];
812 
813  for ( index = 0; index < num; index++ )
814  {
815  array_interfaceABSORPTION_PP[index] = fInterfaceABSORPTION[index].X();
816  array_interfaceABSORPTION[index] = fInterfaceABSORPTION[index].Y();
817  }
818 
819  interfaceMPT->AddProperty("ABSLENGTH", array_interfaceABSORPTION_PP, array_interfaceABSORPTION, num);
820 
821  //cout<<"Interface Material Properties:"<<endl;
822  //interfaceMPT->DumpTable();
823 
824  Interface->SetMaterialPropertiesTable(interfaceMPT);
825 
826  delete[] array_interfaceRINDEX_PP;
827  delete[] array_interfaceRINDEX;
828  delete[] array_interfaceABSORPTION_PP;
829  delete[] array_interfaceABSORPTION;
830  return;
831 }
833 {
834 
835  // HDPE (for tank walls)
836  //----------------------
837 
838  G4double density = 0.94*g/cm3;
839  G4int nel, natoms;
840  G4String name;
841  vector<double>::iterator It;
842 
843  if (HDPE != NULL) {
844  delete HDPE;
845  delete linerMPT;
846  }
847 
848  HDPE = new G4Material(name = "HDPE", density, nel = 2);
849  HDPE->AddElement(elC, natoms = 2);
850  HDPE->AddElement(elH, natoms = 4);
851 
852  linerMPT = new G4MaterialPropertiesTable();
853 
854  // liner abs. length
855 
856  int index;
857  int num = fLinerABSORPTION.GetNPoints();
858 
859  double* array_linerABSORPTION_PP = new double[num+1];
860  double* array_linerABSORPTION = new double[num+1];
861 
862  for ( index = 0; index < num; index++ )
863  {
864  array_linerABSORPTION_PP[index] = fLinerABSORPTION[index].X();
865  array_linerABSORPTION[index] = fLinerABSORPTION[index].Y();
866  }
867 
868  linerMPT->AddProperty("ABSLENGTH", array_linerABSORPTION_PP, array_linerABSORPTION, num);
869 
870  //cout<<"Liner Material Properties:"<<endl;
871  //linerMPT->DumpTable();
872 
873  HDPE->SetMaterialPropertiesTable(linerMPT);
874 
875  delete[] array_linerABSORPTION_PP;
876  delete[] array_linerABSORPTION;
877  return;
878 }
879 
881 {
882  // Optical surface properties
883  //---------------------------
884 
885 
886  // Define the optical model and set optical properties of the surfaces
887  //--------------------------------------------------------------------
888 
889  /*
890  Note that diffuse reflection is implicitly defined, since:
891  SPECULARSPIKECONSTANT+SPECULARLOBECONSTANT+BACKSCATTERCONSTANT=
892  1-DIFFUSEL0BECONSTANT
893  Set index of refraction on the small dielectric interface to 0 so that
894  all the light will undergo total internal reflection (or absorption).
895  It is done in this weird way because the only surface that works (at the
896  moment) for the unified model is the ground surface, and for such a
897  surface a transmission coefficient is always calculated if total internal
898  reflection does not occur. Thus, the only way to ensure all the
899  (non-absorbed) light is reflected back is to choose groundbackpainted
900  finish with the index of refraction for the thin dielectric material
901  set to 0. Note that groundbackpainted with some other index of
902  refraction doesn't work because light reflected off the back (painted)
903  surface undergoes only Lambertian reflection.. This will be fixed up
904  in the future..
905  */
906 
907  if (linerOpticalMPT != NULL) delete linerOpticalMPT;
908  linerOpticalMPT = new G4MaterialPropertiesTable();
909 
910  int index;
912  double* array_linerSPECULARLOBECONSTANT_PP = new double[num+1];
913  double* array_linerSPECULARLOBECONSTANT = new double[num+1];
914 
915  for ( index = 0; index < num; index++ )
916  {
917  array_linerSPECULARLOBECONSTANT_PP[index] = fLinerSPECULARLOBECONSTANT[index].X();
918  array_linerSPECULARLOBECONSTANT[index] = fLinerSPECULARLOBECONSTANT[index].Y();
919  }
920 
921  linerOpticalMPT->AddProperty("SPECULARLOBECONSTANT",
922  array_linerSPECULARLOBECONSTANT_PP,
923  array_linerSPECULARLOBECONSTANT, num);
924 
925 
927  double* array_linerSPECULARSPIKECONSTANT_PP = new double[num+1];
928  double* array_linerSPECULARSPIKECONSTANT = new double[num+1];
929 
930  for ( index = 0; index < num; index++ )
931  {
932  array_linerSPECULARSPIKECONSTANT_PP[index] = fLinerSPECULARSPIKECONSTANT[index].X();
933  array_linerSPECULARSPIKECONSTANT[index] = fLinerSPECULARSPIKECONSTANT[index].Y();
934  }
935 
936  linerOpticalMPT->AddProperty("SPECULARSPIKECONSTANT",
937  array_linerSPECULARSPIKECONSTANT_PP,
938  array_linerSPECULARSPIKECONSTANT, num);
939 
940 
942  double* array_linerREFLECTIVITY_PP = new double[num+1];
943  double* array_linerREFLECTIVITY = new double[num+1];
944 
945  for ( index = 0; index < num; index++ )
946  {
947  array_linerREFLECTIVITY_PP[index] = fLinerREFLECTIVITY[index].X();
948  array_linerREFLECTIVITY[index] = fLinerREFLECTIVITY[index].Y();
949  }
950 
951 
952  linerOpticalMPT->AddProperty("REFLECTIVITY",
953  array_linerREFLECTIVITY_PP,
954  array_linerREFLECTIVITY, num);
955 
956 
958  double* array_linerBACKSCATTERCONSTANT_PP = new double[num+1];
959  double* array_linerBACKSCATTERCONSTANT = new double[num+1];
960 
961  for ( index = 0; index < num; index++ )
962  {
963  array_linerBACKSCATTERCONSTANT_PP[index] = fLinerBACKSCATTERCONSTANT[index].X();
964  array_linerBACKSCATTERCONSTANT[index] = fLinerBACKSCATTERCONSTANT[index].Y();
965  }
966 
967  linerOpticalMPT->AddProperty("BACKSCATTERCONSTANT",
968  array_linerBACKSCATTERCONSTANT_PP,
969  array_linerBACKSCATTERCONSTANT, num);
970 
972  double* array_TYVEK_RINDEX_PP = new double[num+1];
973  double* array_TYVEK_RINDEX = new double[num+1];
974 
975  for ( index = 0; index < num; index++ )
976  {
977  array_TYVEK_RINDEX_PP[index] = fLinerTYVEK_RINDEX[index].X();
978  array_TYVEK_RINDEX[index] = fLinerTYVEK_RINDEX[index].Y();
979  }
980 
981  //cout<<"Liner Optical Properties:"<<endl;
982  //linerOpticalMPT->DumpTable();
983  linerOpticalMPT->AddProperty("RINDEX",
984  array_TYVEK_RINDEX_PP,
985  array_TYVEK_RINDEX, num);
986 
987  OpLinerSurface = new G4OpticalSurface("WallSurface");
988 
989  OpLinerSurface->SetModel(unified);
990  OpLinerSurface->SetType(dielectric_dielectric);
991  OpLinerSurface->SetFinish(groundbackpainted);
992  OpLinerSurface->SetSigmaAlpha(fSIGMA_ALPHA);
993 
994  OpLinerSurface->SetMaterialPropertiesTable(linerOpticalMPT);
995 
996  delete[] array_linerSPECULARLOBECONSTANT_PP;
997  delete[] array_linerSPECULARLOBECONSTANT;
998 
999  delete[] array_linerSPECULARSPIKECONSTANT_PP;
1000  delete[] array_linerSPECULARSPIKECONSTANT;
1001 
1002  delete[] array_linerREFLECTIVITY_PP;
1003  delete[] array_linerREFLECTIVITY;
1004 
1005  delete[] array_linerBACKSCATTERCONSTANT_PP;
1006  delete[] array_linerBACKSCATTERCONSTANT;
1007 
1008  delete[] array_TYVEK_RINDEX_PP;
1009  delete[] array_TYVEK_RINDEX;
1010 
1011  return;
1012 
1013 }
1014 
1015 
1017 {
1018  // Aluminum for electronics boxes and (possibly) solar panel
1019  //----------------------------------------------------------
1020 
1021  G4double density = 2.700*g/cm3;
1022  G4double a = 29.98*g/mole;
1023  G4String name;
1024  G4int z;
1025  Al = new G4Material(name = "Aluminum", z = 13, a, density);
1026  return;
1027 }
1028 
1030 {
1031  // Dirt for the ground
1032  //--------------------
1033  G4String name;
1034  G4int ncomponents;
1035  double fractionmass;
1036  G4double density = 1.2 * g / cm3;
1037  Dirt = new G4Material(name = "Dirt", density, ncomponents = 4);
1038  Dirt->AddElement(elO, fractionmass = 51.0*perCent);
1039  Dirt->AddElement(elSi, fractionmass = 35.2*perCent);
1040  Dirt->AddElement(elAl, fractionmass = 8.2*perCent);
1041  Dirt->AddElement(elFe, fractionmass = 5.6*perCent);
1042 }
1043 
1045 {
1046  expHall_box = new G4Box("World", fExpHall_x, fExpHall_y, fExpHall_z);
1047  ground_solid = new G4Tubs("ground_solid", 0.0*m, fExpHall_x, fGroundThickness/2.,
1048  0.0*deg, 360.0*deg);
1049  tank_solid = new G4Tubs("tank_solid", 0.0*m, fTankRadius, fTankHalfHeight,
1050  0.0*deg, 360.0*deg);
1051  top_solid = new G4Tubs("top_solid",
1053  0.0*deg, 360.0*deg);
1054  side_solid = new G4Tubs("side_solid",
1056  0.0*deg, 360.0*deg);
1057 
1058  inner_solid = new G4Ellipsoid("inner_solid", fPmtRmin, fPmtRmin, fPmtRzmin, -fPmtRzmin, 0.);
1059 
1060  // Active photocathode covered area
1061  pmt_aux = new G4Ellipsoid("pmt_aux", fPmtRmax, fPmtRmax, fPmtRzmax, -fPmtRzmax, -fHeightz);
1062 
1063  pmt_solid = new G4SubtractionSolid("pmt_solid", pmt_aux, inner_solid);
1064  pmt_aux1 = new G4Ellipsoid("pmt_aux1", fPmtRmax, fPmtRmax, fPmtRzmax, -fHeightz, 0.);
1065  pmt_solid1 = new G4SubtractionSolid("pmt_solid1", pmt_aux1, inner_solid);
1066 
1067  interface_in_aux = new G4Ellipsoid("interface_in_aux", fInterfaceRmin, fInterfaceRmin,
1069  interface_out_aux = new G4Ellipsoid("interface_out_aux", fInterfaceRmax, fInterfaceRmax,
1071  interface_solid = new G4SubtractionSolid("interface_solid", interface_out_aux, interface_in_aux);
1072  dome_in_aux = new G4Ellipsoid("dome_in_aux", fDomeRmin, fDomeRmin, fDomeRzmin, -fDomeRzmin, 0.);
1073  dome_out_aux = new G4Ellipsoid("dome_out_aux", fDomeRmax, fDomeRmax, fDomeRzmax, -fDomeRzmax, 0.0*cm);
1074  dome_solid = new G4SubtractionSolid("dome_solid", dome_out_aux, dome_in_aux);
1075 
1076  return;
1077 }
1078 
1080 {
1081  // Experimental hall
1082  //------------------
1083 
1084  expHall_log = new G4LogicalVolume(expHall_box, Air, "World", 0, 0, 0);
1085 
1086  expHall_phys = new G4PVPlacement(0, G4ThreeVector(), "World", expHall_log, 0,
1087  false, 0);
1088 
1089  expHall_log->SetVisAttributes(G4VisAttributes::Invisible);
1090 
1091  // Ground around the detector
1092  // --------------------------
1093 
1094  if (fGroundEnable)
1095  {
1096  ground_log = new G4LogicalVolume(ground_solid, Dirt, "ground_solid");
1097 
1098  ground_phys = new G4PVPlacement(0,
1099  G4ThreeVector(0, 0, -fGroundThickness/2.),
1100  "ground", ground_log, expHall_phys, false, 0);
1101  }
1102  return;
1103 }
1104 
1106 {
1107  G4SDManager *SDMan = G4SDManager::GetSDMpointer();
1108  G4String SDName;
1109 
1110  // Water tank
1111  //-----------
1112 
1113  tank_log = new G4LogicalVolume(tank_solid, Water, "tank_log", 0, 0, 0);
1114 
1115  tank_phys = new G4PVPlacement(0, G4ThreeVector(fTankPos_x, fTankPos_y, fTankPos_z +
1117  "tank", tank_log, expHall_phys, false, 0);
1118 
1119  // Top, bottom, side walls
1120  //------------------------
1121 
1122  top_log = new G4LogicalVolume(top_solid, HDPE, "top_log", 0, 0, 0);
1123  bottom_log = new G4LogicalVolume(top_solid, HDPE, "bottom_log", 0, 0, 0);
1124  side_log = new G4LogicalVolume(side_solid, HDPE, "side_log", 0, 0, 0);
1125 
1126  top_phys = new G4PVPlacement(0, G4ThreeVector(fTankPos_x, fTankPos_y,
1128  "top", top_log, expHall_phys, false, 0);
1129  bottom_phys = new G4PVPlacement(0,G4ThreeVector(fTankPos_x, fTankPos_y,
1130  fTankPos_z + 0.5*fTankThickness),
1131  "bottom", bottom_log, expHall_phys, false, 0);
1132  side_phys = new G4PVPlacement(0,G4ThreeVector(fTankPos_x, fTankPos_y,
1134  "side", side_log, expHall_phys, false, 0);
1135 
1136  //the tank liner
1137  topsurface =
1138  new G4LogicalBorderSurface("topsurface", tank_phys, top_phys,
1139  OpLinerSurface);
1140  bottomsurface =
1141  new G4LogicalBorderSurface("bottomsurface", tank_phys, bottom_phys,
1142  OpLinerSurface);
1143  sidesurface =
1144  new G4LogicalBorderSurface("sidesurface", tank_phys, side_phys,
1145  OpLinerSurface);
1146 
1147  // PMTs and such
1148  //--------------
1149  // the interior volume of the pmt
1150 
1151  // Half the tank's height is used for the z-coordinate since the PMT volumes
1152  // are placed in the logical volume of the water, not the world volume
1153  double dome1_x = fPmt1.GetX(fPmt1.GetCoordinateSystem()) / utl::m * CLHEP::m;
1154  double dome1_y = fPmt1.GetY(fPmt1.GetCoordinateSystem()) / utl::m * CLHEP::m;
1155  double dome1_z = fTankHalfHeight;
1156 
1157  double dome2_x = fPmt2.GetX(fPmt2.GetCoordinateSystem()) / utl::m * CLHEP::m;
1158  double dome2_y = fPmt2.GetY(fPmt2.GetCoordinateSystem()) / utl::m * CLHEP::m;
1159  double dome2_z = fTankHalfHeight;
1160 
1161  double dome3_x = fPmt3.GetX(fPmt3.GetCoordinateSystem()) / utl::m * CLHEP::m;
1162  double dome3_y = fPmt3.GetY(fPmt3.GetCoordinateSystem()) / utl::m * CLHEP::m;
1163  double dome3_z = fTankHalfHeight;
1164 
1165  inner_log = new G4LogicalVolume(inner_solid, Vacuum, "inner_log", 0, 0, 0);
1166 
1167  inner1_phys = new G4PVPlacement(0,
1168  G4ThreeVector(dome1_x, dome1_y, dome1_z), inner_log,
1169  "inner1", tank_log, false, 0);
1170  inner2_phys = new G4PVPlacement(0,
1171  G4ThreeVector(dome2_x, dome2_y, dome1_z), inner_log,
1172  "inner2", tank_log, false, 0);
1173  inner3_phys = new G4PVPlacement(0,
1174  G4ThreeVector(dome3_x, dome3_y, dome1_z), inner_log,
1175  "inner3", tank_log, false, 0);
1176  /*
1177  Why 3 different logical volumes?
1178  because the sensitive detectors are registered
1179  by logical volume
1180  */
1181 
1182  // Sensitive surface of pmt faces
1183 
1184  pmt1_log = new G4LogicalVolume(pmt_solid, Pyrex, "pmt1_log", 0, 0, 0);
1185  pmt2_log = new G4LogicalVolume(pmt_solid, Pyrex, "pmt2_log", 0, 0, 0);
1186  pmt3_log = new G4LogicalVolume(pmt_solid, Pyrex, "pmt3_log", 0, 0, 0);
1187 
1188  pmt1_phys = new G4PVPlacement(0, G4ThreeVector(dome1_x, dome1_y, dome1_z),
1189  pmt1_log,"pmt1", tank_log, false, 0);
1190  pmt2_phys = new G4PVPlacement(0, G4ThreeVector(dome2_x, dome2_y, dome2_z),
1191  pmt2_log,"pmt2", tank_log, false, 0);
1192  pmt3_phys = new G4PVPlacement(0, G4ThreeVector(dome3_x, dome3_y, dome3_z),
1193  pmt3_log,"pmt3", tank_log, false, 0);
1194 
1195  // Non-sensitive surface of pmt faces
1196 
1197  pmt1_log1 = new G4LogicalVolume(pmt_solid1, Pyrex1, "pmt1_log1", 0, 0, 0);
1198  pmt2_log1 = new G4LogicalVolume(pmt_solid1, Pyrex1, "pmt2_log1", 0, 0, 0);
1199  pmt3_log1 = new G4LogicalVolume(pmt_solid1, Pyrex1, "pmt3_log1", 0, 0, 0);
1200 
1201  pmt1_phys1 = new G4PVPlacement(0, G4ThreeVector(dome1_x, dome1_y, dome1_z),
1202  pmt1_log1,"pmt11", tank_log, false, 0);
1203  pmt2_phys1 = new G4PVPlacement(0, G4ThreeVector(dome2_x, dome2_y, dome2_z),
1204  pmt2_log1,"pmt21", tank_log, false, 0);
1205  pmt3_phys1 = new G4PVPlacement(0, G4ThreeVector(dome3_x, dome3_y, dome3_z),
1206  pmt3_log1,"pmt31", tank_log, false, 0);
1207 
1208  // The dome/face interface
1209 
1210  interface1_log = new G4LogicalVolume(interface_solid, Interface, "interface1_log",
1211  0, 0, 0);
1212  interface2_log = new G4LogicalVolume(interface_solid, Interface, "interface2_log",
1213  0, 0, 0);
1214  interface3_log = new G4LogicalVolume(interface_solid, Interface, "interface3_log",
1215  0, 0, 0);
1216 
1217  interface1_phys = new G4PVPlacement(0, G4ThreeVector(dome1_x, dome1_y, dome1_z),
1218  interface1_log,"interface1", tank_log, false, 0);
1219  interface2_phys = new G4PVPlacement(0, G4ThreeVector(dome2_x, dome2_y, dome2_z),
1220  interface2_log,"interface2", tank_log, false, 0);
1221  interface3_phys = new G4PVPlacement(0, G4ThreeVector(dome3_x, dome3_y, dome3_z),
1222  interface3_log,"interface3", tank_log, false, 0);
1223 
1224  // PMT domes
1225  dome1_log = new G4LogicalVolume(dome_solid, Lucite, "dome1_log", 0, 0, 0);
1226  dome2_log = new G4LogicalVolume(dome_solid, Lucite, "dome2_log", 0, 0, 0);
1227  dome3_log = new G4LogicalVolume(dome_solid, Lucite, "dome3_log", 0, 0, 0);
1228 
1229  dome1_phys = new G4PVPlacement(0, G4ThreeVector(dome1_x, dome1_y, dome1_z),
1230  dome1_log, "dome1", tank_log, false, 0);
1231  dome2_phys = new G4PVPlacement(0, G4ThreeVector(dome2_x, dome2_y, dome2_z),
1232  dome2_log, "dome2", tank_log, false, 0);
1233  dome3_phys = new G4PVPlacement(0, G4ThreeVector(dome3_x, dome3_y, dome3_z),
1234  dome3_log, "dome3", tank_log, false, 0);
1235 
1236  // Register the PMT's as sensitive detectors
1237  /*
1238  Note that each PMT must be given a unique index. If programmer
1239  screws up and gives the same index to 2 (or more) PMT's
1240  then photoelectrons will be double counted in that PMT.
1241  Perhaps there is a better way to do this.
1242  */
1243 
1244  SDName = "/Tank/pmt1";
1245  G4TankPMTAction *PMT1SD = new G4TankPMTAction(SDName, 1);
1246 
1247  SDName= "/Tank/pmt2";
1248  G4TankPMTAction *PMT2SD = new G4TankPMTAction(SDName, 2);
1249 
1250  SDName= "/Tank/pmt3";
1251  G4TankPMTAction *PMT3SD = new G4TankPMTAction(SDName, 3);
1252 
1253  SDMan->AddNewDetector(PMT1SD);
1254  SDMan->AddNewDetector(PMT2SD);
1255  SDMan->AddNewDetector(PMT3SD);
1256 
1257  pmt1_log->SetSensitiveDetector(PMT1SD);
1258  pmt2_log->SetSensitiveDetector(PMT2SD);
1259  pmt3_log->SetSensitiveDetector(PMT3SD);
1260 
1261 
1262  // Solar panel
1263  // -----------
1264 
1265  if (fSolarPanelEnable)
1266  {
1267  solarPanel_solid = new G4Box("SolarPanel", fSolarPanelThickness/2, fSolarPanelWidth/2,
1268  fSolarPanelLength/2);
1269 
1270  solarPanel_log = new G4LogicalVolume(solarPanel_solid,
1271  Al, // ** MAY BE CHANGED
1272  "solarPanel_log");
1273  G4RotationMatrix *solarPanelRot = new G4RotationMatrix();
1274  solarPanelRot->rotate(M_PI/2, G4ThreeVector(0, 0, 1));
1275  solarPanelRot->rotate(fSolarPanelTiltAngle, G4ThreeVector(0, 1, 0));
1276 
1277  solarPanel_phys = new G4PVPlacement(solarPanelRot,
1278  G4ThreeVector(fSolarPanelX, fSolarPanelY,
1279  fSolarPanelZ),
1280  "solarPanel",solarPanel_log, expHall_phys,
1281  false, 0);
1282  }
1283 
1284  // Chunk of Aluminum representing electronics box
1285  //-----------------------------------------------
1286 
1287  if (fElecBoxEnable)
1288  {
1289  elecBox_solid = new G4Box("ElecBox", fElecBoxThickness/2, fElecBoxWidth/2,
1290  fElecBoxLength/2);
1291 
1292  elecBox_log = new G4LogicalVolume(elecBox_solid,
1293  Al, // ** MAY BE CHANGED
1294  "elecBox_log");
1295  G4RotationMatrix *elecBoxRot = new G4RotationMatrix();
1296  elecBoxRot->rotate(M_PI/2, G4ThreeVector(0, 0, 1));
1297  elecBoxRot->rotate(fElecBoxTiltAngle, G4ThreeVector(0, 1, 0));
1298 
1299  elecBox_phys = new G4PVPlacement(elecBoxRot,
1300  G4ThreeVector(fElecBoxX, fElecBoxY, fElecBoxZ),
1301  "elecBox",elecBox_log, expHall_phys, false, 0);
1302  }
1303 
1304  return;
1305 }
1306 
1307 
1309 {
1310  SetDetectorParameters(); // updates G4TankConstruction data tables
1311 
1312  // clean up old geometry list
1313  G4GeometryManager::GetInstance()->OpenGeometry();
1314  G4PhysicalVolumeStore::GetInstance()->Clean();
1315  G4LogicalVolumeStore::GetInstance()->Clean();
1316  G4SolidStore::GetInstance()->Clean();
1317 
1318  // update material properties (deleting out of date materials and
1319  // replacing them with new ones).
1320  CreateWater();
1321  CreatePyrex();
1322  CreateLucite();
1323  CreateInterface();
1324  CreateHDPE();
1325  CreateLiner();
1326 
1327  // assemble objects with new sizes and materials
1328  G4RunManager::GetRunManager()->DefineWorldVolume(CreateTank(), true);
1329 }
const double eV
Definition: GalacticUnits.h:35
unsigned int GetNPoints() const
constexpr double mm
Definition: AugerUnits.h:113
constexpr double perCent
Definition: AugerUnits.h:282
Point object.
Definition: Point.h:32
const utl::TabulatedFunction & GetDomeRefractionIndex() const
Refraction index for the PMT dome.
constexpr double cm3
Definition: AugerUnits.h:119
Detector description interface for Station-related data.
ArrayIterator XEnd()
end of array of X
CoordinateSystemPtr GetCoordinateSystem() const
Get the coordinate system of the current internal representation.
Definition: BasicVector.h:234
Class to hold collection (x,y) points and provide interpolation between them.
const utl::TabulatedFunction & GetLinerSpecularLobe() const
Tyvek liner specular lobe constant as a function of photon energy.
const PMT & GetPMT(const int id) const
Get specified PMT by id.
void PushBack(const double x, const double y)
Branch GetChild(const std::string &childName) const
Get child of this Branch by child name.
Definition: Branch.cc:211
TransformerConstructor< DerivedCSPolicy > CoordinateSystem
The normal coordinate system type.
const utl::TabulatedFunction & GetWaterAbsorptionLength() const
Water absorption length as a function of photon energy.
Exception for reporting variable out of valid range.
constexpr double deg
Definition: AugerUnits.h:140
Class representing a document branch.
Definition: Branch.h:107
const utl::TabulatedFunction & GetLinerReflectivity() const
Tyvek liner reflectivity as a function of photon energy.
constexpr double mole
Definition: AugerUnits.h:262
double GetX(const CoordinateSystemPtr &coordinateSystem) const
Definition: BasicVector.h:206
const utl::Point & GetPosition() const
PMT position.
double GetRadius() const
Radius of the tank (water only)
double GetHeight() const
Height of the tank (water only)
constexpr double g
Definition: AugerUnits.h:200
void GetData(bool &b) const
Overloads of the GetData member template function.
Definition: Branch.cc:644
const utl::TabulatedFunction & GetLinerSpecularSpike() const
Tyvek liner specular spike constant as a function of photon energy.
class that handles Geant4 SD simulation
ArrayIterator YBegin()
begin of array of Y
double GetY(const CoordinateSystemPtr &coordinateSystem) const
Definition: BasicVector.h:209
static const sdet::Station * GetCurrentDetectorStation()
ArrayIterator XBegin()
begin of array of X
double GetThickness() const
Thickness of the tank walls.
constexpr double cm
Definition: AugerUnits.h:117
const utl::TabulatedFunction & GetWaterRefractionIndex() const
Water refraction index as a function of photon energy.
Detector description interface for SDetector-related data.
Definition: SDetector.h:42
Main configuration utility.
Definition: CentralConfig.h:51
ArrayIterator YEnd()
end of array of Y
const utl::TabulatedFunction & GetFaceRefractionIndex() const
Refraction index for the PMT face.
double GetLinerSigmaAlpha() const
Tyvek liner sigma_alpha roughness parameter.
double Y(const double x) const
Get or interpolate the Y value that corresponds to parameter x.
#define ERROR(message)
Macro for logging error messages.
Definition: ErrorLogger.h:165
constexpr double m
Definition: AugerUnits.h:121
const utl::TabulatedFunction & GetRTVRefractionIndex() const
Refraction index for the interface between PMT face and dome (RTV)
utl::Branch GetTopBranch(const std::string &id)
Get top branch for moduleConfigLink with given id (XML files)

, generated on Tue Sep 26 2023.