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