Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
UltraDetectorConstruction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // --------------------------------------------------------------
28 // GEANT 4 - ULTRA experiment example
29 // --------------------------------------------------------------
30 //
31 // Code developed by:
32 // B. Tome, M.C. Espirito-Santo, A. Trindade, P. Rodrigues
33 //
34 // ****************************************************
35 // * UltraDetectorConstruction.cc
36 // ****************************************************
37 //
38 // Class used in the definition of the Ultra setup consisting of:
39 // - the UVscope detector
40 // - an optional reflecting surface
41 // Optical photons can reach the UVscope either directly or after reflection in the
42 // surface, which can be polished or diffusing.
43 // The main part of the UVscope definition is the Fresnel lens construction based
44 // on the UltraFresnelLens class.
45 //
46 #include <cmath>
47 
49 #include "UltraPMTSD.hh"
50 #include "UltraFresnelLens.hh"
51 
52 #include "G4PhysicalConstants.hh"
53 #include "G4SystemOfUnits.hh"
54 #include "G4Material.hh"
55 #include "G4MaterialTable.hh"
56 #include "G4Element.hh"
57 #include "G4ElementTable.hh"
59 #include "G4Box.hh"
60 #include "G4Sphere.hh"
61 #include "G4Tubs.hh"
62 #include "G4LogicalVolume.hh"
63 #include "G4RotationMatrix.hh"
64 #include "G4ThreeVector.hh"
65 #include "G4Transform3D.hh"
66 #include "G4PVPlacement.hh"
67 #include "G4OpBoundaryProcess.hh"
68 #include "G4VisAttributes.hh"
69 #include "G4Colour.hh"
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72 
74  logicalPMT(0)
75 {
76  // Define wavelength limits for materials definition
77  lambda_min = 200*nm ;
78  lambda_max = 700*nm ;
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82 
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88 {
89  ConstructTableMaterials();
90 
91 
92 
93 // The experimental Hall
94 // ---------------------
95 
96  G4double World_x = 1.*m;
97  G4double World_y = 1.*m;
98  G4double World_z = 2*m;
99 
100  G4Box * World_box = new G4Box("World",World_x,World_y,World_z);
101 
102  // Get Air pointer from static funcion - (G4Material::GetMaterial)
103 
104 G4String name;
105 G4Material *Air = G4Material::GetMaterial(name = "Air");
106 G4LogicalVolume *World_log ;
107 World_log = new G4LogicalVolume(World_box,Air,"World",0,0,0);
108 
109 G4VPhysicalVolume *World_phys ;
110 World_phys = new G4PVPlacement(0,G4ThreeVector(),"World",World_log,0,false,0);
111 
112  G4VisAttributes* UniverseVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
113  UniverseVisAtt->SetVisibility(true);
114  UniverseVisAtt->SetForceWireframe(true);
115  World_log->SetVisAttributes(UniverseVisAtt);
117 
118 
119 
120  G4cout << "\n \n \n \n \n \n \n \n \n \n \n \n \n " << G4endl ;
121 
122  G4cout << "######################################################" << G4endl ;
123  G4cout << "# #" << G4endl ;
124  G4cout << "# #" << G4endl ;
125  G4cout << "# UltraDetectorConstruction: #" << G4endl ;
126  G4cout << "# #" << G4endl ;
127  G4cout << "# #" << G4endl ;
128 
129  ConstructUVscope(World_phys);
130 
131 
132  G4cout << "# #" << G4endl ;
133  G4cout << "# #" << G4endl ;
134  G4cout << "######################################################" << G4endl ;
135 
136 
137 #ifdef ULTRA_MIRROR_USE
138 
139  G4cout << "Using mirror reflecting surface " << G4endl ;
140 
141  G4VPhysicalVolume* Mirror ;
142  Mirror = ConstructMirror(World_phys);
143 
144 #elif ULTRA_GROUND_USE
145 
146  G4cout << "Using ground reflecting surface " << G4endl ;
147 
148  G4VPhysicalVolume* Ground ;
149  Ground = ConstructGround(World_phys);
150 
151 #else
152 
153  G4cout << "No reflecting surface used" << G4endl ;
154 
155 #endif
156 
157  return World_phys;
158 }
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161 
163 {
164  UltraPMTSD* PMTSD = new UltraPMTSD("PMTSD");
165  SetSensitiveDetector(logicalPMT,PMTSD);
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169 
170 void UltraDetectorConstruction::ConstructTableMaterials()
171 {
172  G4double a, z, density;
174  G4int nel;
175 
176 
177 // ------------- Elements -------------
178  a = 1.01*g/mole;
179  G4Element* elH = new G4Element(name="Hydrogen", symbol="H", z=1., a);
180 
181  a = 12.01*g/mole;
182  G4Element* elC = new G4Element(name="Carbon", symbol="C", z=6., a);
183 
184  a = 14.01*g/mole;
185  G4Element* elN = new G4Element(name="Nitrogen", symbol="N", z=7., a);
186 
187  a = 16.00*g/mole;
188  G4Element* elO = new G4Element(name="Oxygen", symbol="O", z=8., a);
189 
190  a = 28.09*g/mole;
191  G4Element* elSi = new G4Element(name="Silicon", symbol="Si", z=14., a);
192 
193 
194 // ------------- Materials -------------
195 
196 
197 // Air
198 // ---
199  density = 1.29e-03*g/cm3;
200  G4Material* Air = new G4Material(name="Air", density, nel=2);
201  Air->AddElement(elN, .7);
202  Air->AddElement(elO, .3);
203 
204 
205 // Aluminum
206 // ---------
207  a = 26.98*g/mole;
208  density = 2.7*g/cm3;
209  new G4Material(name="Aluminum", z=13., a, density);
210 
211 
212 // Quartz
213 // -------
214 // density = 2.200*g/cm3; // fused quartz
215  density = 2.64*g/cm3; // crystalline quartz (c.f. PDG)
216  G4Material *Quartz = new G4Material(name="Quartz",density, nel=2);
217  Quartz->AddElement(elSi, 1) ;
218  Quartz->AddElement(elO , 2) ;
219 
220 
221 // PMMA C5H8O2 ( Acrylic )
222 // -------------
223  density = 1.19*g/cm3;
224  G4Material* Acrylic = new G4Material(name="Acrylic", density, nel=3);
225  Acrylic->AddElement(elC, 5);
226  Acrylic->AddElement(elH, 8);
227  Acrylic->AddElement(elO, 2);
228 
229 
230 /////////////////////////////////////////////
231 // Construct Material Properties Tables
232 /////////////////////////////////////////////
233 
234  const G4int NUMENTRIES = 2;
235 
236  // Energy bins
237  G4double X_RINDEX[NUMENTRIES] = {h_Planck*c_light/lambda_max, h_Planck*c_light/lambda_min} ;
238 
239 
240  // Air
241  G4double RINDEX_AIR[NUMENTRIES] = {1.00, 1.00} ;
242 
243 // Air refractive index at 20 oC and 1 atm (from PDG)
244  for(G4int j=0 ; j<NUMENTRIES ; j++){
245  RINDEX_AIR[j] = RINDEX_AIR[j] + 2.73*std::pow(10.0,-4) ;
246  }
247 
249  MPT_Air->AddProperty("RINDEX", X_RINDEX, RINDEX_AIR, NUMENTRIES);
250  Air->SetMaterialPropertiesTable(MPT_Air);
251 
252 //////////////////////////////////////////////////////////////////////////////////////
253 // Photomultiplier (PMT) window
254 // The refractive index is for lime glass;
255 // wavelength dependence is not included and value at 400nm is used.
256 //////////////////////////////////////////////////////////////////////////////////////
257 
258  // Refractive index
259 
260  const G4int N_RINDEX_QUARTZ = 2 ;
261  G4double X_RINDEX_QUARTZ[N_RINDEX_QUARTZ] = {h_Planck*c_light/lambda_max, h_Planck*c_light/lambda_min} ;
262  G4double RINDEX_QUARTZ[N_RINDEX_QUARTZ] = {1.54, 1.54};
263 
265  MPT_PMT->AddProperty("RINDEX", X_RINDEX_QUARTZ, RINDEX_QUARTZ, N_RINDEX_QUARTZ);
266 
267  Quartz->SetMaterialPropertiesTable(MPT_PMT);
268 
269 
270 //////////////////////////////////////////////////////////////////
271 // ACRYLIC Optical properties
272 //////////////////////////////////////////////////////////////////
273 
274 // Refractive index
275 
276  const G4int NENTRIES = 11 ;
277  G4double LAMBDA_ACRYLIC[NENTRIES] ;
278 
279 
280  G4double RINDEX_ACRYLIC[NENTRIES] ;
281  G4double ENERGY_ACRYLIC[NENTRIES] ;
282 
283 // Parameterization for refractive index of High Grade PMMA
284 
285  G4double bParam[4] = {1760.7010,-1.3687,2.4388e-3,-1.5178e-6} ;
286 
287  for(G4int i=0;i<NENTRIES; i++){
288 
289  LAMBDA_ACRYLIC[i] = lambda_min + i*(lambda_max-lambda_min)/float(NENTRIES-1) ;
290  RINDEX_ACRYLIC[i] = 0.0 ;
291 
292  for (G4int jj=0 ; jj<4 ; jj++)
293  {
294  RINDEX_ACRYLIC[i] += (bParam[jj]/1000.0)*std::pow(LAMBDA_ACRYLIC[i]/nm,jj) ;
295  }
296 
297  ENERGY_ACRYLIC[i] = h_Planck*c_light/LAMBDA_ACRYLIC[i] ; // Convert from wavelength to energy ;
298 // G4cout << ENERGY_ACRYLIC[i]/eV << " " << LAMBDA_ACRYLIC[i]/nm << " " << RINDEX_ACRYLIC[i] << G4endl ;
299 
300  }
301 
303  MPT_Acrylic->AddProperty("RINDEX", ENERGY_ACRYLIC, RINDEX_ACRYLIC, NENTRIES);
304 
305 
306 // Absorption
307  const G4int NENT = 25 ;
308  G4double LAMBDAABS[NENT] =
309  {
310  100.0,
311  246.528671, 260.605103, 263.853516, 266.019104, 268.726105,
312  271.433136, 273.598724, 276.305725, 279.554138, 300.127380,
313  320.159241, 340.191101, 360.764343, 381.337585, 399.745239,
314  421.401276, 440.891724, 460.382172, 480.414001, 500.987274,
315  520.477722, 540.509583, 559.458618,
316  700.0
317  } ;
318 
319  G4double ABS[NENT] = // Transmission (in %) of 3mm thick PMMA
320  {
321  0.0000000,
322  0.0000000, 5.295952, 9.657321, 19.937695, 29.283491,
323  39.252335, 48.598133, 58.255451, 65.109039, 79.439247,
324  85.669785, 89.719627, 91.277260, 91.588783, 91.900307,
325  91.588783, 91.277260, 91.277260, 91.588783, 91.588783,
326  91.900307, 91.900307, 91.588783,
327  91.5
328  } ;
329 
330 
331  MPT_Acrylic->AddProperty("ABSLENGTH", new G4MaterialPropertyVector()) ;
332  for(G4int i=0;i<NENT; i++){
333  G4double energy = h_Planck*c_light/(LAMBDAABS[i]*nm) ;
334  G4double abslength ;
335 
336  if (ABS[i] <= 0.0) {
337  abslength = 1.0/kInfinity ;
338  }
339  else {
340  abslength = -3.0*mm/(std::log(ABS[i]/100.0)) ;
341  }
342 
343  MPT_Acrylic->AddEntry("ABSLENGTH", energy, abslength);
344 
345  }
346 
347  Acrylic->SetMaterialPropertiesTable(MPT_Acrylic);
348 
349 
350 //////////////////////////////////////////////////////////////////
351 
353 
354 }
355 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
356 
357 G4VPhysicalVolume* UltraDetectorConstruction::ConstructMirror(G4VPhysicalVolume *World_phys){
358 
359  G4double Mirror_x = 40.0*cm;
360  G4double Mirror_y = 40.0*cm;
361  G4double Mirror_z = 1*cm;
362 
363  G4Box * boxMirror = new G4Box("Mirror",Mirror_x,Mirror_y,Mirror_z);
364 
365  // Get Air pointer from static funcion - (G4Material::GetMaterial)
366 
367 G4String name;
368 G4Material *Al = G4Material::GetMaterial(name = "Aluminum");
369 G4LogicalVolume *logMirror ;
370 logMirror = new G4LogicalVolume(boxMirror,Al,"Mirror",0,0,0);
371 
372 
373 G4ThreeVector SurfacePosition = G4ThreeVector(0*m,0*m,1.5*m) ;
374 
375 // Rotate reflecting surface by 45. degrees around the OX axis.
376 
377 G4RotationMatrix *Surfrot = new G4RotationMatrix(G4ThreeVector(1.0,0.0,0.0),-pi/4.);
378 
379 G4VPhysicalVolume *physMirror ;
380 physMirror = new G4PVPlacement(Surfrot,SurfacePosition,"MirrorPV",logMirror,World_phys,false,0);
381 
382 G4VisAttributes* SurfaceVisAtt = new G4VisAttributes(G4Colour(0.0,0.0,1.0));
383 SurfaceVisAtt->SetVisibility(true);
384 SurfaceVisAtt->SetForceWireframe(true);
385 logMirror->SetVisAttributes(SurfaceVisAtt);
386 
387 
388 //////////////////////////////////////////////////////////////////////////////////////////
389 // Optical properties of the interface between the Air and Reflective Surface
390 // For Mirror, reflectivity is set at 95% and specular reflection is assumed.
391 
392 
393 G4OpticalSurface *OpticalAirMirror = new G4OpticalSurface("AirMirrorSurface");
394 OpticalAirMirror->SetModel(unified);
395 OpticalAirMirror->SetType(dielectric_dielectric);
396 OpticalAirMirror->SetFinish(polishedfrontpainted);
397 
398 const G4int NUM = 2;
399 G4double XX[NUM] = {h_Planck*c_light/lambda_max, h_Planck*c_light/lambda_min} ;
400 G4double ICEREFLECTIVITY[NUM] = { 0.95, 0.95 };
401 
403 AirMirrorMPT->AddProperty("REFLECTIVITY", XX, ICEREFLECTIVITY,NUM);
404 OpticalAirMirror->SetMaterialPropertiesTable(AirMirrorMPT);
405 
406 
407 
408 new G4LogicalBorderSurface("Air/Mirror Surface",World_phys,physMirror,OpticalAirMirror);
409 
410  return physMirror ;
411 
412 }
413 
414 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
415 
416 G4VPhysicalVolume* UltraDetectorConstruction::ConstructGround(G4VPhysicalVolume *World_phys){
417 
418  G4double Ground_x = 40.0*cm;
419  G4double Ground_y = 40.0*cm;
420  G4double Ground_z = 1*cm;
421 
422  G4Box * boxGround = new G4Box("Ground",Ground_x,Ground_y,Ground_z);
423 
424  // Get Air pointer from static funcion - (G4Material::GetMaterial)
425 
426 G4String name;
427 G4Material *Al = G4Material::GetMaterial(name = "Aluminum");
428 G4LogicalVolume *logGround ;
429 logGround = new G4LogicalVolume(boxGround,Al,"Ground",0,0,0);
430 
431 
432 G4ThreeVector SurfacePosition = G4ThreeVector(0*m,0*m,1.5*m) ;
433 
434 // Rotate reflecting surface by 45. degrees around the OX axis.
435 
436 G4RotationMatrix *Surfrot = new G4RotationMatrix(G4ThreeVector(1.0,0.0,0.0),-pi/4.);
437 
438 G4VPhysicalVolume *physGround ;
439 physGround = new G4PVPlacement(Surfrot,SurfacePosition,"GroundPV",logGround,World_phys,false,0);
440 
441 G4VisAttributes* SurfaceVisAtt = new G4VisAttributes(G4Colour(0.0,0.0,1.0));
442 SurfaceVisAtt->SetVisibility(true);
443 SurfaceVisAtt->SetForceWireframe(true);
444 logGround->SetVisAttributes(SurfaceVisAtt);
445 
446 
447 //////////////////////////////////////////////////////////////////////////////////////////
448 // Optical properties of the interface between the Air and Reflective Surface
449 // For Ground, reflectivity is set to 95% and diffusive reflection is assumed.
450 
451 
452 G4OpticalSurface *OpticalAirGround = new G4OpticalSurface("AirGroundSurface");
453 OpticalAirGround->SetModel(unified);
454 OpticalAirGround->SetType(dielectric_dielectric);
455 OpticalAirGround->SetFinish(groundfrontpainted);
456 
457  const G4int NUM = 2;
458 G4double XX[NUM] = {h_Planck*c_light/lambda_max, h_Planck*c_light/lambda_min} ;
459 G4double ICEREFLECTIVITY[NUM] = { 0.95, 0.95 };
460 
462 AirGroundMPT->AddProperty("REFLECTIVITY", XX, ICEREFLECTIVITY,NUM);
463 OpticalAirGround->SetMaterialPropertiesTable(AirGroundMPT);
464 
465 
466 new G4LogicalBorderSurface("Air/Ground Surface",World_phys,physGround,OpticalAirGround);
467 
468  return physGround ;
469 
470 }
471 
472 
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
475 
476 G4VPhysicalVolume* UltraDetectorConstruction::ConstructUVscope(G4VPhysicalVolume *World_phys)
477 {
478 
479  // ------------- Volumes --------------
480 
481  ////////////////////////////////////////////////////////////////////////////////////////////////////////
482 
483  G4cout << "# #" << G4endl ;
484  G4cout << "# Building the Telescope ... #" << G4endl ;
485  G4cout << "# #" << G4endl ;
486 
487  /////////////////////////////////////////////////////////////
488  // UVscope housing is a cylinder made of 1 mm thick aluminum
489  /////////////////////////////////////////////////////////////
490 
491  G4double UVscopeHeight = 1030.0*mm ;
492  G4double UVscopeDiameter = 518.0*mm ;
493  G4double UVscopeThickness = 1.0*mm ;
494  G4double UVscopeBaffle = 514.0*mm ;
495 
496  G4double UVscopeInnerRadius = UVscopeDiameter/2.0-UVscopeThickness ;
497  G4double UVscopeOuterRadius = UVscopeDiameter/2.0 ;
498 
499  G4ThreeVector UVscopePosition = G4ThreeVector(0.0*m,0.0*m,-1.0*m) ;
500  G4String name;
501  G4Material* Al = G4Material::GetMaterial(name = "Aluminum");
502 
503 
504  G4Tubs *solidUVscope =
505  new G4Tubs("UVscopeSolid",UVscopeInnerRadius,UVscopeOuterRadius,UVscopeHeight/2.0,0.0,twopi) ;
506  G4LogicalVolume *logicUVscope =
507  new G4LogicalVolume(solidUVscope,Al,"UVscopeLV",0,0,0);
508  G4VPhysicalVolume *physicalUVscope =
509  new G4PVPlacement(0,UVscopePosition,"UVSCopePV",logicUVscope,World_phys,false,0);
510 
511 
512  //////////////////////////////////////
513  // Back cover of the UVscope cylinder
514  //////////////////////////////////////
515 
516  G4Tubs *solidUVscopeBack =
517  new G4Tubs("UVscopeBackSolid",0.0,UVscopeOuterRadius,UVscopeThickness/2.0,0.0,twopi) ;
518 
519  G4LogicalVolume *logicUVscopeBack =
520  new G4LogicalVolume(solidUVscopeBack,Al,"UVscopeBackLV",0,0,0);
521 
522  G4ThreeVector UVscopeBackPosition ;
523  UVscopeBackPosition = UVscopePosition+G4ThreeVector(0.0*mm,0.0*mm,-(UVscopeHeight/2.0+UVscopeThickness/2.0)) ;
524  G4VPhysicalVolume *physicalUVscopeBack =
525  new G4PVPlacement(0,UVscopeBackPosition,"UVscopeBack",logicUVscopeBack,World_phys,false,0);
526 
527 
528 
529  ////////////////////////////////////////////////////////////////////////////////////////////////////////
530 
531  G4cout << "# #" << G4endl ;
532  G4cout << "# Building the Fresnel lens ... #" << G4endl ;
533  G4cout << "# #" << G4endl ;
534 
535  G4double LensDiameter = 457*mm ; // Size of the optical active area of the lens.
536  G4int LensNumOfGrooves = 13 ;
537  //G4int LensNumOfGrooves = 129 ;
538  //G4int LensNumOfGrooves = 1287 ;
539 
540  G4double LensBorderThickness = 2.8*mm ; // Thickness of the border area.
541  G4double LensFocalLength = 441.973*mm ; // This parameter depends on the lens geometry, etc !!
542  G4Material *LensMaterial = G4Material::GetMaterial(name = "Acrylic") ;
543  G4ThreeVector LensPosition = UVscopePosition+G4ThreeVector(0.0*mm,0.0*mm,UVscopeHeight/2.0-UVscopeBaffle) ;
544 
545 
546  FresnelLens = new UltraFresnelLens(LensDiameter,LensNumOfGrooves,LensMaterial,World_phys,LensPosition) ;
547 
548 
549  ///////////////////////////////////
550  // Lens supporting ring (aluminum)
551  ///////////////////////////////////
552 
553  G4Tubs *solidLensFrame = new G4Tubs("LensFrame",LensDiameter/2.0,UVscopeInnerRadius,LensBorderThickness/2.0,0.0,twopi) ;
554  G4LogicalVolume *logicLensFrame = new G4LogicalVolume(solidLensFrame,Al,"LensFrameLV",0,0,0);
555 
556  G4ThreeVector LensFramePosition ;
557  LensFramePosition = LensPosition+G4ThreeVector(0.0*mm,0.0*mm,-((FresnelLens->GetThickness())/2.0+solidLensFrame->GetDz())) ;
558 
559  G4VPhysicalVolume *physicalLensFrame =
560  new G4PVPlacement(0,LensFramePosition,"LensFramePV",logicLensFrame,World_phys,false,0);
561 
562  ////////////////////////////////////////////////////////////////////////////////////////////////////////
563 
564 
565  G4cout << "# #" << G4endl ;
566  G4cout << "# Building the photomultiplier ... #" << G4endl ;
567  G4cout << "# #" << G4endl ;
568 
569 
570  // Photomultiplier window is a spherical section made of quartz
571 
572  G4double PMT_thick = 1.0*mm ; // Thickness of PMT window
573  G4double PMT_curv = 65.5*mm ; // Radius of curvature of PMT window
574  G4double StartTheta = (180.0-31.2)*pi/180. ;
575  G4double EndTheta = 31.2*pi/180. ;
576 
577  G4Sphere *solidPMT ;
578  solidPMT = new G4Sphere("PMT_solid",PMT_curv-PMT_thick,PMT_curv,0.0,twopi,StartTheta,EndTheta);
579 
580  G4Material* Quartz = G4Material::GetMaterial(name = "Quartz");
581  logicalPMT = new G4LogicalVolume(solidPMT,Quartz,"PMT_log",0,0,0);
582 
583 
584  // Place PMT is at Lens Focus
585 
586  G4ThreeVector PMTpos = LensPosition + G4ThreeVector(0.0*cm,0.0*cm,-(LensFocalLength+PMT_curv)) ;
587 
588  // Rotate PMT window through the axis OX by an angle = 180. degrees
589 
590  G4RotationMatrix *PMTrot = new G4RotationMatrix(G4ThreeVector(1.0,0.0,0.0),pi);
591  new G4PVPlacement(PMTrot,PMTpos,"PMT1",logicalPMT,World_phys,false,0);
592 
593 
594  G4VisAttributes* PMTVisAtt = new G4VisAttributes(true,G4Colour(0.0,0.0,1.0)) ;
595  logicalPMT->SetVisAttributes(PMTVisAtt);
596 
597  //////////////////////////////////////////////////////////////////////////////////////////
598  // Optical properties of the interface between the Air and the walls of the
599  // UVscope cylinder (5% reflectivity)
600 
601 
602  G4cout << "# Defining interface's optical properties ... #" << G4endl ;
603  G4cout << "# #" << G4endl ;
604 
605 
606  G4OpticalSurface *OpticalAirPaint = new G4OpticalSurface("AirPaintSurface");
607  OpticalAirPaint->SetModel(unified);
608  OpticalAirPaint->SetType(dielectric_dielectric);
609  OpticalAirPaint->SetFinish(groundfrontpainted);
610 
611  const G4int NUM = 2;
612  G4double XX[NUM] = {h_Planck*c_light/lambda_max, h_Planck*c_light/lambda_min} ;
613  G4double BLACKPAINTREFLECTIVITY[NUM] = { 0.05, 0.05 };
614  //G4double WHITEPAINTREFLECTIVITY[NUM] = { 0.99, 0.99 };
615 
617  AirPaintMPT->AddProperty("REFLECTIVITY", XX, BLACKPAINTREFLECTIVITY,NUM);
618  OpticalAirPaint->SetMaterialPropertiesTable(AirPaintMPT);
619 
620  //OpticalAirPaint->DumpInfo();
621 
622  new G4LogicalBorderSurface("Air/UVscope Cylinder Surface",World_phys,physicalUVscope,OpticalAirPaint);
623 
624  new G4LogicalBorderSurface("Air/LensFrame Surface",World_phys,physicalLensFrame,OpticalAirPaint);
625 
626  new G4LogicalBorderSurface("Air/UVscope Back Cover Surface",World_phys,physicalUVscopeBack,OpticalAirPaint);
627 
628 
629  /////////////////////////////////////////////////////////////////////////////////////
630 
631 
632  G4VisAttributes* LensVisAtt = new G4VisAttributes(G4Colour(1.0,0.0,0.0)) ; // Red
633  LensVisAtt ->SetVisibility(true);
634 
635 
636  if (FresnelLens){
637  FresnelLens->GetPhysicalVolume()->GetLogicalVolume()->SetVisAttributes(LensVisAtt);
638  }
639 
640  G4VisAttributes* UVscopeVisAtt = new G4VisAttributes(G4Colour(0.5,0.5,0.5)) ; // Gray
641  UVscopeVisAtt ->SetVisibility(true);
642 
643  physicalUVscope ->GetLogicalVolume()->SetVisAttributes(UVscopeVisAtt);
644  physicalUVscopeBack ->GetLogicalVolume()->SetVisAttributes(UVscopeVisAtt);
645  physicalLensFrame ->GetLogicalVolume()->SetVisAttributes(UVscopeVisAtt);
646 
647  /////////////////////////////////////////////////////////////////////////////////////
648 
649  G4cout << "# #" << G4endl ;
650  G4cout << "# UVscope is built ! ... #" << G4endl ;
651  G4cout << "# #" << G4endl ;
652 
653  return physicalUVscope;
654 }
655 
656 
657 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
658 
659 
660 
void SetFinish(const G4OpticalSurfaceFinish)
G4String symbol
Definition: TRTMaterials.hh:40
void SetForceWireframe(G4bool)
G4Material * Air
Definition: TRTMaterials.hh:57
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
float h_Planck
Definition: hepunit.py:263
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
G4double z
Definition: TRTMaterials.hh:39
Definition: G4Box.hh:63
G4VPhysicalVolume * GetPhysicalVolume()
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
Definition: G4Material.hh:247
void AddEntry(const char *key, G4double aPhotonEnergy, G4double aPropertyValue)
void SetVisibility(G4bool)
Definition: G4Tubs.hh:84
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:564
G4Element * elC
Definition: TRTMaterials.hh:48
const XML_Char * name
G4double GetDz() const
int G4int
Definition: G4Types.hh:78
G4MaterialPropertyVector * AddProperty(const char *key, G4double *PhotonEnergies, G4double *PropertyValues, G4int NumEntries)
G4Element * elN
Definition: TRTMaterials.hh:44
G4Element * elH
Definition: TRTMaterials.hh:50
G4double density
Definition: TRTMaterials.hh:39
function g(Y1, Y2, PT2)
Definition: hijing1.383.f:5205
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
G4PhysicsOrderedFreeVector G4MaterialPropertyVector
G4Element * elO
Definition: TRTMaterials.hh:46
void SetSensitiveDetector(const G4String &logVolName, G4VSensitiveDetector *aSD, G4bool multi=false)
G4LogicalVolume * GetLogicalVolume() const
static const G4VisAttributes Invisible
#define G4endl
Definition: G4ios.hh:61
void AddElement(G4Element *element, G4int nAtoms)
Definition: G4Material.cc:345
double G4double
Definition: G4Types.hh:76
void SetModel(const G4OpticalSurfaceModel model)
G4double GetThickness()
void SetMaterialPropertiesTable(G4MaterialPropertiesTable *anMPT)
void SetType(const G4SurfaceType &type)
void SetVisAttributes(const G4VisAttributes *pVA)
float c_light
Definition: hepunit.py:257
G4int nel
Definition: TRTMaterials.hh:41
G4Material * Al
Definition: TRTMaterials.hh:74