FDsimG4PMT.cc
Go to the documentation of this file.
1 #include "FDsimG4Defs.hh"
2 #include "FDsimG4PMT.hh"
4 #include "FDsimG4Colours.hh"
5 #include "G4RunManager.hh"
6 #include "G4UnitsTable.hh"
7 #include "G4SDManager.hh"
8 #include "G4Material.hh"
9 #include "G4MaterialTable.hh"
10 #include "G4Polyhedra.hh"
11 #include "G4BREPSolidPolyhedra.hh"
12 #include "G4LogicalVolume.hh"
13 #include "G4PVPlacement.hh"
14 #include "G4ThreeVector.hh"
15 #include "G4VisAttributes.hh"
16 #include "G4LogicalSkinSurface.hh"
17 #include "G4OpticalSurface.hh"
18 
19 using namespace TelescopeSimulatorLX ;
20 using namespace std;
21 
22 FDsimG4PMT::FDsimG4PMT(G4String PMTname)
23 {
24  // Data from Photonis XP3062 data sheet
25 
26  fWindowSide = 21.9*mm; // Size of the hexagonal window side
27  fPhotocathodeSide = 36.0*mm/2.0*(2.0/std::sqrt(3.0)); // Size of the photocathode hexagonal side
28  fThickness = 1.*mm;
29  fWindowMaterial = G4Material::GetMaterial("LimeGlass");
30 
31  fPMT_name=PMTname;
32 
33  // Sensitive Detector Manager
34  SDmanager = G4SDManager::GetSDMpointer();
35 
36  fPMT_log = NULL;
37  fPhotocathode_log = NULL;
38 
39  CheckParameters();
40  MakeLogicalVolume();
41 }
42 
43 
45 
46 
48 {
49  if (!fWindowMaterial)
50  G4Exception("Error in FDsimG4PMT constructor: G4Material of window of PMT not specified !!!");
51 
52 
53  if ( fWindowSide < fPhotocathodeSide )
54  G4Exception("Error in FDsimG4PMT constructor: Photocathode size > PMT size !!!");
55 
56 }
57 
58 
60 {
61  G4cerr << " =================================================================== " << G4endl ;
62  G4cerr << " = PMT parameters " << G4endl ;
63  G4cerr << " = " << G4endl ;
64  G4cerr << " = PMT window side : " << fWindowSide/mm << " (mm) " << G4endl ;
65  G4cerr << " = Photocathode side : " << fPhotocathodeSide/mm << " (mm) " << G4endl ;
66  G4cerr << " = PMT thickness : " << fThickness/mm << " (mm) " << G4endl ;
67  G4cerr << " = " << G4endl ;
68  G4cerr << " = Window Material : " << fWindowMaterial->GetName() << G4endl ;
69  G4cerr << " = " << G4endl ;
70  G4cerr << " =================================================================== " << G4endl ;
71 }
72 
73 
75 {
76  const G4int numSides = 6;
77  const G4int numZplanes = 2;
78 
79  // Planar surfaces cases
80  G4double rMinVec[numZplanes] = {0.0, 0.0};
81  G4double rMaxVec[numZplanes] = {fWindowSide,fWindowSide}; // This is the radius of the hexagon corners
82  G4double zVec[numZplanes] = {-fThickness/2.0, fThickness/2.0};
83 
84  G4Polyhedra* PMT_solid
85  = new G4Polyhedra("PMT",0.,2.0*pi,numSides,numZplanes,zVec,rMinVec,rMaxVec);
86 
87 
88  fPMT_log
89  = new G4LogicalVolume(PMT_solid,(G4Material*)fWindowMaterial,"PMT",0,0,0);
90 
91 
92  fPMT_log->SetVisAttributes(new G4VisAttributes(yellow));
93 
94  // Build photocathode
95  G4double rPhotocathMinVec[numZplanes] = {0.0, 0.0};
96  G4double rPhotocathMaxVec[numZplanes] = {fPhotocathodeSide,fPhotocathodeSide}; // This is the radius of the hexagon corners
97  G4double zPhotocathVec[numZplanes] = {0.9*(fThickness/2.0), 0.95*(fThickness/2.0)};
98 
99  G4Polyhedra* Photocathode_solid
100  = new G4Polyhedra("Photocathode",0.,2.0*pi,numSides,numZplanes,zPhotocathVec,rPhotocathMinVec,rPhotocathMaxVec);
101 
102  fPhotocathode_log
103  = new G4LogicalVolume(Photocathode_solid,(G4Material*)fWindowMaterial,"Photocathode",0,0,0);
104 
105  G4String name = "Photocathode_" + fPMT_name;
106  new G4PVPlacement(0,G4ThreeVector(),fPhotocathode_log,name,fPMT_log,false,0);
107 
108  fPhotocathode_log->SetVisAttributes(new G4VisAttributes(red));
109 
110  // Photocathode surface proprieties
111 
112  G4OpticalSurface* OpticalPhotocathode = new G4OpticalSurface("PhotocathodeSurface");
113  OpticalPhotocathode->SetModel(unified);
114  OpticalPhotocathode->SetType(dielectric_metal);
115  OpticalPhotocathode->SetFinish(polished);
116 
117  G4MaterialPropertiesTable* PhotocathodeMPT = new G4MaterialPropertiesTable();
118 
119  const G4RunManager* runManager = G4RunManager::GetRunManager();
120  const FDsimG4DetectorConstruction* DetectorConstruction =
121  dynamic_cast<const FDsimG4DetectorConstruction*>(runManager->GetUserDetectorConstruction());
122 
123  const G4bool& useComplexRindex = DetectorConstruction->fPhotocathodeUseComplexRindex;
124  const G4bool& useThinFilm = DetectorConstruction->fPhotocathodeUseThinFilm;
125  const vector<G4double>& wavelength = DetectorConstruction->fPhotocathodeWavelength;
126 
127  if (useComplexRindex) {
128  const vector<G4double>& reRindex = DetectorConstruction->fPhotocathodeReRindex;
129  const vector<G4double>& imRindex = DetectorConstruction->fPhotocathodeImRindex;
130 
131  if (useThinFilm) {
132  const G4double& thickness = DetectorConstruction->fPhotocathodeThickness;
133  PhotocathodeMPT->AddConstProperty("THICKNESS",thickness) ;
134  }
135 
136  PhotocathodeMPT->AddProperty("REALRINDEX", new G4MaterialPropertyVector()) ;
137  PhotocathodeMPT->AddProperty("IMAGINARYRINDEX", new G4MaterialPropertyVector()) ;
138 
139  G4double energy;
140  unsigned nEntries = wavelength.size();
141  for (unsigned i=0; i < nEntries; i++) {
142  energy = h_Planck*c_light/wavelength[i] ;
143  PhotocathodeMPT->AddEntry("REALRINDEX",energy,reRindex[i]);
144  PhotocathodeMPT->AddEntry("IMAGINARYRINDEX",energy,imRindex[i]);
145  }
146  }
147  else {
148  const vector<G4double>& reflectivity = DetectorConstruction->fPhotocathodeReflectivity;
149  PhotocathodeMPT->AddProperty("REFLECTIVITY", new G4MaterialPropertyVector()) ;
150 
151  G4double energy;
152  unsigned nEntries = wavelength.size();
153  for (unsigned i=0; i < nEntries; i++) {
154  energy = h_Planck*c_light/wavelength[i];
155  PhotocathodeMPT->AddEntry("REFLECTIVITY",energy,reflectivity[i]);
156  }
157  }
158 
159  PhotocathodeMPT->AddProperty("EFFICIENCY", new G4MaterialPropertyVector()) ;
160  PhotocathodeMPT->AddEntry("EFFICIENCY", 0.0*eV,1.);
161  PhotocathodeMPT->AddEntry("EFFICIENCY",10.0*eV,1.);
162 
163  OpticalPhotocathode->SetMaterialPropertiesTable(PhotocathodeMPT);
164 
165  new G4LogicalSkinSurface("PhotocathodeSurface",fPhotocathode_log,OpticalPhotocathode);
166 
167  const G4int VerbosityLevel = DetectorConstruction->fVerbosityLevel;
168  if (VerbosityLevel > 0) {
169  static G4int ifirst = 1;
170  if (ifirst) {
171  G4cerr << "Material properties table for photocathode" << G4endl;
172  PhotocathodeMPT->DumpTable();
173  ifirst = 0;
174  }
175  }
176 }
const double eV
Definition: GalacticUnits.h:35
constexpr double mm
Definition: AugerUnits.h:113
static const G4Colour red(1.0, 0.0, 0.0)
FDsimG4PMT(G4String PMTname)
Definition: FDsimG4PMT.cc:22
static const G4Colour yellow(1.0, 1.0, 0.5)

, generated on Tue Sep 26 2023.