FDsimG4Filter.cc
Go to the documentation of this file.
1 //#warning Preliminary implementation of filter frame. Exact geometry and material optical properties still needed.
2 #include "FDsimG4Filter.hh"
4 #include "FDsimG4Colours.hh"
5 
6 #include "G4UnitsTable.hh"
7 #include "G4PVPlacement.hh"
8 #include "G4Box.hh"
9 #include "G4Tubs.hh"
10 #include "G4UnionSolid.hh"
11 #include "G4SubtractionSolid.hh"
12 #include "G4OpticalSurface.hh"
13 #include "G4LogicalBorderSurface.hh"
14 #include "G4LogicalSkinSurface.hh"
15 #include "G4Material.hh"
16 #include "G4MaterialTable.hh"
17 #include "G4LogicalVolume.hh"
18 #include "G4ThreeVector.hh"
19 #include "G4VisAttributes.hh"
20 #include "G4RunManager.hh"
21 
22 #include <utl/unconfig.h>
23 #include <utl/config.h>
24 
25 using namespace std;
26 using namespace TelescopeSimulatorLX;
27 
28 
29 FDsimG4Filter::FDsimG4Filter(G4VPhysicalVolume* const mother)
30 {
31  fDetector =
32  dynamic_cast<const FDsimG4DetectorConstruction*>(G4RunManager::GetRunManager()->GetUserDetectorConstruction());
33 
34  fMother = mother;
35  fPosition = fDetector->fFilterPosition;
36  fThickness = fDetector->fFilterThickness;
37  fRadius = 0.5 * fDetector->fFilterDiameter;
38  fMaterial = fDetector->fFilterMaterial;
39 
40  fFrameSize = fDetector->fFilterFrameSize;
41  fFrameArmWidth = 5.0*cm;
42  fFrameThickness = fDetector->fFilterFrameThickness;
43 
44  fFilter_log = nullptr;
45  MakeFilter();
46 }
47 
48 
49 void
50 FDsimG4Filter::DumpInfo()
51 {
52  G4cerr << " ================================================================\n"
53  " = Filter parameters\n"
54  " =\n"
55  " = Radius : " << fRadius/mm << " mm\n"
56  " =\n"
57  " = Thickness : " << fThickness/mm << " mm\n"
58  " =\n"
59  " = Material : " << fMaterial->GetName()
60  << " =\n"
61  " ================================================================" << G4endl;
62 }
63 
64 
65 void
66 FDsimG4Filter::MakeFilter()
67 {
68  const G4double tolerance = 0.1*mm;
69 
70  MakeFrame();
71 
72  G4Tubs* const filterSolid =
73  new G4Tubs("Filter", 0, fRadius, 0.5*fThickness, 0, 2*pi);
74 
75  fFilter_log = new G4LogicalVolume(filterSolid, (G4Material*)fMaterial, "Filter", 0,0,0);
76  fFilter_log->SetVisAttributes(new G4VisAttributes(magenta));
77 
78  const G4ThreeVector filterPos =
79  fPosition -
80  G4ThreeVector(0, 0, 0.5*(fFrameThickness + fThickness) + tolerance);
81 
82  /*G4VPhysicalVolume* const filterPhys =*/
83  new G4PVPlacement(0, filterPos, "UVFilter", fFilter_log, fMother, false, 0);
84  G4VPhysicalVolume* const framePhys =
85  new G4PVPlacement(0, fPosition, "Filter frame", fFrame_log, fMother, false, 0);
86 
87  new G4LogicalBorderSurface("Frame Surface", fMother, framePhys, fFrameOpticalSurface);
88 
89 #if (G4VERSION_NUMBER >= 900)
90  const G4int verbosityLevel = fDetector->fVerbosityLevel;
91  if (verbosityLevel > 2) {
92  filterPhys->CheckOverlaps(10000);
93  framePhys->CheckOverlaps(10000);
94  }
95 #endif
96 }
97 
98 
99 void
100 FDsimG4Filter::MakeFrame()
101 {
102  const G4double tolerance = 0.1*mm;
103 
104  const G4double frameX = fFrameSize;
105  const G4double frameY = fFrameSize;
106  const G4double frameZ = fFrameThickness;
107 
108  G4Box* const frameBox = new G4Box("Frame_tmp", 0.5*frameX, 0.5*frameY, 0.5*frameZ);
109 
110  G4Tubs* const frameHole = new G4Tubs("Hole", 0, fRadius+tolerance, 10*frameZ, 0, 2*pi);
111 
112  G4Box* const frameArmBox = new G4Box("FrameArm", 0.5*fFrameArmWidth, 0.5*frameY - 0.5*(0.5*frameY - fRadius), 0.5*frameZ);
113 
114  const G4double yOffset = 0.5*frameY/3;
115  // This offset is required for visualization purposes
116  // see http://hypernews.slac.stanford.edu:5090/HyperNews/geant4/get/visualization/558/1.html
117  const G4double zOffset = fFrameThickness;
118 
119  const G4ThreeVector frameArmPos_1 = G4ThreeVector(0, 0, zOffset);
120  const G4ThreeVector frameArmPos_2 = G4ThreeVector(0, yOffset, zOffset);
121  const G4ThreeVector frameArmPos_3 = G4ThreeVector(0, -yOffset, zOffset);
122 
123  G4SubtractionSolid* const frameTmp = new G4SubtractionSolid("Frame", frameBox, frameHole);
124 
125  G4RotationMatrix* const zRot90 = new G4RotationMatrix();
126  zRot90->rotateZ(90*deg);
127 
128  G4UnionSolid* const frame_1 = new G4UnionSolid("Frame_1", frameTmp, frameArmBox, nullptr, frameArmPos_1);
129  G4UnionSolid* const Frame_2 = new G4UnionSolid("Frame_2", frame_1, frameArmBox, zRot90, frameArmPos_2);
130  G4UnionSolid* const fFrameSolid = new G4UnionSolid("Frame_3", frame_2, frameArmBox, zRot90, frameArmPos_3);
131 
132  fFrame_log = new G4LogicalVolume(fFrameSolid, G4Material::GetMaterial("Aluminum"), "Filter Frame", 0,0,0);
133  fFrame_log->SetVisAttributes(new G4VisAttributes(gray));
134 
135  // Set reflection properties
136 
137  const G4String& frameReflectionType = fDetector->fFilterFrameReflectionType;
138  const vector<G4double>& wavelength = fDetector->fFilterFrameWavelength;
139  const vector<G4double>& reflectivity = fDetector->fFilterFrameReflectivity;
140 
141  fFrameOpticalSurface = new G4OpticalSurface("Filter frame");
142  fFrameOpticalSurface->SetModel(unified);
143  fFrameOpticalSurface->SetType(dielectric_metal);
144  fFrameOpticalSurface->SetFinish(polished);
145 
146  G4MaterialPropertiesTable* const frameMPT = new G4MaterialPropertiesTable();
147  frameMPT->AddProperty("REFLECTIVITY", new G4MaterialPropertyVector());
148 
149  for (unsigned iwl = 0, n = wavelength.size(); iwl < n; ++iwl) {
150  const G4double energy = h_Planck * c_light / wavelength[iwl];
151  frameMPT->AddEntry("REFLECTIVITY", energy, reflectivity[iwl]);
152  }
153 
154  const G4int VerbosityLevel = fDetector->fVerbosityLevel;
155  if (VerbosityLevel > 0) {
156  G4cerr << "==== Filter Frame Surface properties ====\n"
157  "==== Reflection type : " << FrameReflectionType << "====" << G4endl;
158  }
159 
160  if (FrameReflectionType == "Diffusive") {
161  const SurfaceRoughness& frameRoughness = fDetector->fFilterFrameRoughness;
162  const G4double& sigma_alpha = frameRoughness.Sigma_alpha;
163  const vector<G4double>& wavelengthRoughness = frameRoughness.Wavelength;
164  const vector<G4double>& specularSpike = frameRoughness.SpecularSpike;
165  const vector<G4double>& specularLobe = frameRoughness.SpecularLobe;
166  const vector<G4double>& backscatter = frameRoughness.Backscatter;
167 
168  fFrameOpticalSurface->SetFinish(ground);
169  fFrameOpticalSurface->SetSigmaAlpha(sigma_alpha);
170 
171  frameMPT->AddProperty("SPECULARSPIKECONSTANT", new G4MaterialPropertyVector());
172  frameMPT->AddProperty("SPECULARLOBECONSTANT", new G4MaterialPropertyVector());
173  frameMPT->AddProperty("BACKSCATTERCONSTANT", new G4MaterialPropertyVector());
174 
175  for (unsigned iwl = 0, n = wavelengthRoughness.size(); iwl < n; ++iwl) {
176  const G4double energy = h_Planck * c_light / wavelengthRoughness[iwl];
177  frameMPT->AddEntry("SPECULARSPIKECONSTANT", energy, specularSpike[iwl]);
178  frameMPT->AddEntry("SPECULARLOBECONSTANT", energy, specularLobe[iwl]);
179  frameMPT->AddEntry("BACKSCATTERCONSTANT", energy, backscatter[iwl]);
180  }
181 
182  if (VerbosityLevel > 0) {
183  G4cerr << "==== sigma_alpha = " << sigma_alpha/deg << " deg. ====\n"
184  "==== Parameters vs wavelength (Wl (nm) \t Reflectivity \t SpecularSpike \t SpecularLobe \t Backscatter) ====" << G4endl;
185 
186  for (unsigned i = 0, n = wavelengthRoughness.size(); i < n; ++i) {
187  G4cerr << wavelengthRoughness[i]/nm << '\t'
188  << reflectivity[i] << '\t'
189  << specularSpike[i] << '\t'
190  << specularLobe[i] << '\t'
191  << backscatter[i] << G4endl;
192  }
193  }
194 
195  }
196 
197  fFrameOpticalSurface->SetMaterialPropertiesTable(FrameMPT);
198 }
constexpr double mm
Definition: AugerUnits.h:113
constexpr double deg
Definition: AugerUnits.h:140
static const G4Colour magenta(1.0, 0.5, 1.0)
constexpr double gray
Definition: AugerUnits.h:269
constexpr double cm
Definition: AugerUnits.h:117

, generated on Tue Sep 26 2023.