FDsimG4DetectorConstruction.cc
Go to the documentation of this file.
1 //
2 // Description: DetectorConstruction for FDsimG4 application
3 //
4 
5 #include <string>
7 #include "FDsimG4Materials.hh"
8 #include "FDsimG4FilterSD.hh"
9 #include "FDsimG4Filter.hh"
10 #include "FDsimG4Camera.hh"
11 #include "FDsimG4Mirror.hh"
13 #include "FDsimG4MirrorSegment.hh"
14 #include "FDsimG4CorrectorRing.hh"
18 #include "FDsimG4LensSD.hh"
19 #include "FDsimG4XMLManager.hh"
20 #include "FDsimG4Colours.hh"
21 #include "FDsimG4Defs.hh"
22 
23 #include "G4SDManager.hh"
24 #include "G4Material.hh"
25 #include "G4MaterialTable.hh"
26 #include "G4Element.hh"
27 #include "G4ElementTable.hh"
28 #include "G4UserLimits.hh"
29 #include "G4Box.hh"
30 #include "G4Trap.hh"
31 #include "G4Tubs.hh"
32 #include "G4Sphere.hh"
33 #include "G4Orb.hh"
34 #include "G4UnionSolid.hh"
35 #include "G4SubtractionSolid.hh"
36 #include "G4LogicalVolume.hh"
37 #include "G4AssemblyVolume.hh"
38 #include "G4RotationMatrix.hh"
39 #include "G4ThreeVector.hh"
40 #include "G4Transform3D.hh"
41 #include "G4PVPlacement.hh"
42 #include "G4VisAttributes.hh"
43 #include "G4Colour.hh"
44 #include "G4RunManager.hh"
45 
46 #include <unconfig.h>
47 #include <utl/config.h>
48 #include <utl/ErrorLogger.h>
49 #include <det/Detector.h>
50 #include <fdet/FDetector.h>
51 #include <fdet/Eye.h>
52 #include <fdet/Telescope.h>
53 #include <fdet/Corrector.h>
54 #include <fdet/Diaphragm.h>
55 #include <fdet/Camera.h>
56 
57 #include "G4Version.hh"
58 
59 #undef __G4_V9__
60 #if (G4VERSION_NUMBER>=900)
61 #define __G4_V9__
62 #include "G4GeometryTolerance.hh"
63 #endif
64 
65 using namespace TelescopeSimulatorLX;
66 
68 
70  const FDsimG4XMLManager XMLManager) :
71  fTelescopeName(name),
72  fXMLManager(XMLManager),
73  fVerbosityLevel(0),
74  fWorld_phys(NULL),
75  fTelescopeMother_phys(NULL),
76  fLensFilterMother_phys(NULL),
77  fFilter(NULL),
78  fCorrectorRing(NULL),
79  fHasCorrectorRing(false),
80  fCorrectorRingInnerRadius(0.0*mm),
81  fCorrectorRingOuterRadius(0.0*mm),
82  fMirror(NULL),
83  fCamera(NULL),
84  fUseSensitiveDetectors(false)
85 {
86 }
87 
88 
90 {
91  delete fMaterials;
92  delete fWorld_phys;
93  delete fTelescopeMother_phys;
95  delete fFilter;
96  delete fCorrectorRing;
97  delete fMirror;
98  delete fCamera;
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
102 
103 G4VPhysicalVolume*
105 {
106  Init(); // Initialize detector geometry and material properties.
110 
112  if (fHasCorrectorRing)
116 
117  if (fVerbosityLevel > 0) {
118  G4cerr << "# #" << G4endl;
119  G4cerr << "# Done geometry #" << G4endl;
120  G4cerr << "# #" << G4endl;
121  }
122 
123  // Visualization attributes
124 
125  G4VisAttributes* UniverseVisAtt = new G4VisAttributes(G4Colour(1.0,1.0,1.0));
126  UniverseVisAtt->SetVisibility(false);
127  fWorld_phys->GetLogicalVolume()->SetVisAttributes(UniverseVisAtt);
128 
129  if (fVerbosityLevel > 0) {
130  G4cerr << *(G4Material::GetMaterialTable()) << G4endl;
131  G4cerr << "Tree of Sensitive Detectors: " << G4endl;
132  (G4SDManager::GetSDMpointer())->ListTree();
133  }
134 
135  return fWorld_phys;
136 }
137 
138 
139 G4VPhysicalVolume*
141 {
142  // The World
143  // ---------
144  G4Box * World_box
145  = new G4Box("World",fWorld_x/2.0,fWorld_y/2.0,fWorld_z/2.0);
146 
147  G4LogicalVolume * World_log
148  = new G4LogicalVolume(World_box,fWorldMaterial,"World",0,0,0);
149 
150  return new G4PVPlacement(0,G4ThreeVector(),"World",World_log,0,false,0);
151 }
152 
153 
154 G4VPhysicalVolume*
156 {
157  // Telescope mother
158  G4Box * TelMotherBox = new G4Box("TelMother",fWorld_x/4.0,fWorld_y/4.0,fWorld_z/4.0);
159 
160  G4LogicalVolume * TelMother_log
161  = new G4LogicalVolume(TelMotherBox,fWorldMaterial,"TelMother",0,0,0,true);
162 
163  TelMother_log->SetSmartless(50);
164 
165  G4ThreeVector TelMotherPos = G4ThreeVector(0.0*mm,0.0*mm,0.0*mm);
166  G4RotationMatrix *TelMotherRot = new G4RotationMatrix();
167 
168  // OZ points from the mirror to the corrector ring
169 
170  TelMotherRot->rotateZ(pi/2.);
171 
172  {
173  G4VisAttributes* tmpVisAtt = new G4VisAttributes();
174  tmpVisAtt->SetVisibility(false);
175  TelMother_log->SetVisAttributes(tmpVisAtt);
176  }
177 
178  G4VPhysicalVolume* telMother_phys =
179  new G4PVPlacement(TelMotherRot,TelMotherPos,"TelescopeMother",TelMother_log,fWorld_phys,false,0);
180 #if (G4VERSION_NUMBER>=800)
181  if (fVerbosityLevel > 2)
182  telMother_phys ->CheckOverlaps(10000);
183 #endif
184 
185  G4double xSize = fFilterFrameSize + fgTOLERANCE;
186  G4double ySize = fFilterFrameSize + fgTOLERANCE;
187  G4double zSize = 2.0*(fFilterPosition.z() +
188  2.0*fFilterFrameThickness +
190  fgTOLERANCE);
191 
192  G4Box* LensFilterSolid =
193  new G4Box("LensFilterMother",xSize/2.0,ySize/2.0,zSize/2.0);
194 
195  G4LogicalVolume * LensFilterMother_log
196  = new G4LogicalVolume(LensFilterSolid,fWorldMaterial,"LensFilterMother",0,0,0,true);
197 
198  LensFilterMother_log->SetSmartless(50);
199 
201  new G4PVPlacement(0,fCorrectorRingPosition,"LensFilterMother",LensFilterMother_log,telMother_phys,false,0);
202 
203 #if (G4VERSION_NUMBER>=800)
204  if (fVerbosityLevel > 2)
205  fLensFilterMother_phys ->CheckOverlaps(10000);
206 #endif
207 
208  {
209  G4VisAttributes* VisAtt = new G4VisAttributes(yellow);
210  VisAtt->SetVisibility(false);
211  LensFilterMother_log->SetVisAttributes(VisAtt);
212  }
213 
214  return telMother_phys;
215 }
216 
217 
219 {
221 
222  if (fVerbosityLevel > 0)
223  Filter->DumpInfo();
224 
225  G4LogicalVolume * Filter_log = Filter->GetLogicalVolume();
226 
227  if (!Filter_log)
228  G4Exception("Error !!! Filter logical volume does not exhist !!!");
229 
230 
231  // Filter sensitive detector
233  G4double FilterSDThick = 0.1*mm;
234 
235  G4double inRadius = 0.0*mm;
236  G4double ouRadius = fFilterRadius;
237  G4double dZ = FilterSDThick;
238  G4double sPhi = 0.0*deg;
239  G4double dPhi = 360.0*deg;
240 
241 
242  G4Tubs* FilterSDTubs = new G4Tubs("FilterSD",inRadius,ouRadius,dZ/2.0,sPhi,dPhi);
243  G4LogicalVolume * FilterSD_log =
244  new G4LogicalVolume(FilterSDTubs,fWorldMaterial,"FilterSD",0,0,0);
245 
246  G4ThreeVector FilterSDPos =
247  fFilterPosition - G4ThreeVector(0.,0.,fFilterThickness/2.+FilterSDThick/2.+fgTOLERANCE);
248 
249  G4PVPlacement* filterSD_phys =
250  new G4PVPlacement(0,FilterSDPos,"FilterSD",FilterSD_log,fLensFilterMother_phys,false,0);
251 
252  // Set sensitive detector
253 
254  FDsimG4FilterSD* FilterSD = new FDsimG4FilterSD("FilterSD");
255  (G4SDManager::GetSDMpointer())->AddNewDetector(FilterSD);
256  FilterSD_log->SetSensitiveDetector(FilterSD);
257 
258 #if (G4VERSION_NUMBER>=800)
259  if (fVerbosityLevel > 2)
260  filterSD_phys->CheckOverlaps(10000);
261 #endif
262 
263  {
264  G4VisAttributes* tmpVisAtt = new G4VisAttributes(yellow);
265  tmpVisAtt->SetVisibility(false);
266  FilterSD_log ->SetVisAttributes(tmpVisAtt);
267  }
268  }
269 
270  return Filter;
271 }
272 
275 {
276 #ifdef __G4_V9__
277  const double kAngTolerance = G4GeometryTolerance::GetInstance()->GetAngularTolerance();
278 #endif
279 
280  G4double inRadius = 1000.0*mm;
281  G4double ouRadius = 2900.0*mm;
282  G4double sPhi = -50.*deg;
283  G4double dPhi = 100.*deg-kAngTolerance;
284  G4double sTheta = 50.*deg;
285  G4double dTheta = 80.0*deg-kAngTolerance;
286 
287  G4Sphere *cameraSphere =
288  new G4Sphere("CameraMother",inRadius,ouRadius,sPhi,dPhi,sTheta,dTheta);
289 
290  G4LogicalVolume* mother_log
291  = new G4LogicalVolume(cameraSphere,fWorldMaterial,"CameraMother",0,0,0,true);
292 
293  mother_log->SetSmartless(50);
294 
295  G4ThreeVector motherPos(0.0*mm,0.0*mm,0.0*mm);
296  G4RotationMatrix* motherRot = new G4RotationMatrix();
297 
298  motherRot->rotateY(-pi/2.0);
299 
300  G4PVPlacement* mother_phys =
301  new G4PVPlacement(motherRot,motherPos,"CameraMother",mother_log,fTelescopeMother_phys,false,0);
302 
303  {
304  G4VisAttributes* tmpVisAtt = new G4VisAttributes();
305  tmpVisAtt->SetVisibility(false);
306  mother_log->SetVisAttributes(tmpVisAtt);
307  }
308 
309  FDsimG4Camera * Camera = new FDsimG4Camera(mother_phys);
310 
311 #if (G4VERSION_NUMBER>=800)
312  if (fVerbosityLevel > 2)
313  mother_phys->CheckOverlaps(10000);
314 #endif
315 
316  return Camera;
317 }
318 
319 
322 {
323  G4bool kUseOptimization = true;
324 
325  G4RotationMatrix *MirrorRot = new G4RotationMatrix();
326  MirrorRot->rotateY(-pi/2.0);
327 
328 #ifdef __G4_V9__
329  const double kAngTolerance = G4GeometryTolerance::GetInstance()->GetAngularTolerance();
330 #endif
331 
332  G4double inRadius = 3000.0*mm;
333  G4double ouRadius = 3500.0*mm;
334  G4double sPhi = -40.*deg;
335  G4double dPhi = 80.*deg-kAngTolerance;
336  G4double sTheta = 50.*deg;
337  G4double dTheta = 80.0*deg-kAngTolerance;
338 
339  G4Sphere* mirrorSphere = new G4Sphere("MirrorMother",
340  inRadius,ouRadius,sPhi,dPhi,sTheta,dTheta);
341 
342  G4LogicalVolume* mirrorMother_log =
343  new G4LogicalVolume(mirrorSphere,fWorldMaterial,"MirrorMother",0,0,0);
344 
345  mirrorMother_log -> SetOptimisation(kUseOptimization);
346 
347  if (kUseOptimization){
348  mirrorMother_log->SetSmartless(50);
349  }
350 
351  G4VPhysicalVolume* mirrorMother_phys =
352  new G4PVPlacement(MirrorRot, G4ThreeVector(), "MirrorMother", mirrorMother_log,
353  fTelescopeMother_phys,false, 0);
354 
355  {
356  G4VisAttributes* tmpVisAtt = new G4VisAttributes(yellow);
357  tmpVisAtt->SetVisibility(false);
358  mirrorMother_log->SetVisAttributes(tmpVisAtt);
359  }
360 
361  FDsimG4Mirror* Mirror = new FDsimG4Mirror(fTelescopeName,fMirrorDataFile,mirrorMother_phys);
362 
363 #if (G4VERSION_NUMBER>=800)
364  if (fVerbosityLevel > 2)
365  mirrorMother_phys->CheckOverlaps(10000);
366 #endif
367 
368  return Mirror;
369 }
370 
371 
374 {
375  FDsimG4VCorrectorRing *CorrectorRing = NULL;
376 
377  G4double inRadius = fCorrectorRingInnerRadius;
378  G4double ouRadius = fCorrectorRingOuterRadius;
379  G4double sPhi = -fCorrectorRingPetalSize/2.;
380  G4double dPhi = fCorrectorRingPetalSize;
381  G4double thick = fCorrectorRingThickness;
382  G4Material* material = fCorrectorRingMaterial;
383  const G4int kNumPoints = 1000;
384 
385  if (fCorrectorRingProfile == "KG"){
386  CorrectorRing =
387  new FDsimG4CorrectorRingKG(inRadius,ouRadius,sPhi,dPhi,thick,material,kNumPoints);
388  }
389  else if(fCorrectorRingProfile == "Circular"){
390  CorrectorRing =
391  new FDsimG4CorrectorRing(inRadius,ouRadius,sPhi,dPhi,thick,material,kNumPoints);
392 
393  }
394  else if(fCorrectorRingProfile == "UpperLimit"){
395  CorrectorRing =
396  new FDsimG4CorrectorRingUpperLimit(inRadius,ouRadius,sPhi,dPhi,thick,material,kNumPoints);
397 
398  }
399  else if(fCorrectorRingProfile == "LowerLimit"){
400  CorrectorRing =
401  new FDsimG4CorrectorRingLowerLimit(inRadius,ouRadius,sPhi,dPhi,thick,material,kNumPoints);
402  }
403  else{
404  G4Exception("Error in FDsimG4DetectorConstruction::ConstructCorrectorRing: unknown profile !!!");
405  }
406  G4LogicalVolume *ring_log = CorrectorRing->GetLogicalVolume();
407  G4PVPlacement * Lens_phys = NULL;
408 
409  G4String name;
410 
411  for (G4int ii=0 ; ii < fCorrectorRingNumPetals ; ii++) {
412  G4RotationMatrix * Rotation = new G4RotationMatrix();
413  G4double angle = 2.0*pi/fCorrectorRingNumPetals * ii;
414 
415  Rotation->rotateZ(angle);
416 
417  Lens_phys = new G4PVPlacement(Rotation,G4ThreeVector(),"CorrectorRing",ring_log,
418  fLensFilterMother_phys,false,ii);
419 
420 #if (G4VERSION_NUMBER>=800)
421  if (fVerbosityLevel > 2)
422  Lens_phys->CheckOverlaps(10000);
423 #endif
424  }
425 
426  // Construct lens Sensitive detector
428  ConstructLensSD();
429 
430  return CorrectorRing;
431 }
432 
433 
434 void
436 {
437  // Build the Corrector Ring Sensitive Detector
438 
439  G4double LensSDThick = 0.1*mm;
440 
441  G4double inRadius = 0.0*mm;
442  G4double ouRadius = fCorrectorRingOuterRadius;
443  G4double dZ = LensSDThick;
444  G4double sPhi = 0.0*deg;
445  G4double dPhi = 360.0*deg;
446 
447  G4Tubs* LensSDTubs = new G4Tubs("LensSD",inRadius,ouRadius,dZ/2.0,sPhi,dPhi);
448 
449  G4LogicalVolume* LensSD_log = new G4LogicalVolume(LensSDTubs,fWorldMaterial,"LensSD",0,0,0);
450 
451  G4ThreeVector LensSDPos =
452  fCorrectorRingPosition - G4ThreeVector(0.,0.,fCorrectorRingThickness+1.*mm+LensSDThick/2.+0.001*mm);
453 
454  G4PVPlacement* LensSD_phys =
455  new G4PVPlacement(0,LensSDPos,"LensSD",LensSD_log,fLensFilterMother_phys,false,0);
456 
457 #if (G4VERSION_NUMBER>=800)
458  if (fVerbosityLevel > 2)
459  LensSD_phys->CheckOverlaps(10000);
460 #endif
461 
462  // Define sensitive detector
463  FDsimG4LensSD* LensSD = new FDsimG4LensSD("LensSD");
464  (G4SDManager::GetSDMpointer())->AddNewDetector(LensSD);
465  LensSD_log->SetSensitiveDetector(LensSD);
466 
467  G4VisAttributes* tmpVisAtt = new G4VisAttributes(yellow);
468  tmpVisAtt->SetVisibility(false);
469  LensSD_log ->SetVisAttributes(tmpVisAtt);
470 }
471 
472 
473 void
475 {
476  // Definition of the Telescope Dome
477  G4double DomeThick = 20.0*cm;
478  G4ThreeVector DomeSize(5000.0*mm,5000.0*mm,5000.0*mm);
479  G4ThreeVector DomeInnerSize =
480  DomeSize-G4ThreeVector(2.0*DomeThick,2.0*DomeThick,2.0*DomeThick);
481 
482  G4Box* DomeOuterBox =
483  new G4Box("TelDome",DomeSize.x()/2.,DomeSize.y()/2.,DomeSize.z()/2.);
484 
485  G4double DomeHoleSide = fFilterFrameSize + fgTOLERANCE;
486  G4Box* DomeHole =
487  new G4Box("Hole",DomeHoleSide/2.0,DomeHoleSide/2.0,DomeThick/2.0 + 1.0*mm);
488 
489  G4ThreeVector FilterHolePosition(0.0*mm,0.0*mm,DomeInnerSize.z()/2+DomeThick/2.0);
490  G4Box* DomeInnerBox =
491  new G4Box("TelDomeIn",DomeInnerSize.x()/2.0,DomeInnerSize.y()/2.0,DomeInnerSize.z()/2.0);
492 
493  G4UnionSolid* DomeVoid =
494  new G4UnionSolid("TelescopeDomeVoid",DomeInnerBox,DomeHole,0,FilterHolePosition);
495  G4SubtractionSolid* DomeBox =
496  new G4SubtractionSolid("TelescopeDome",DomeOuterBox,DomeVoid);
497 
498  G4Material* material = G4Material::GetMaterial("Aluminum");
499  G4LogicalVolume* Dome_log = new G4LogicalVolume(DomeBox,material,"TelescopeDome",0,0,0,true);
500 
501  {
502  G4VisAttributes* DomeVisAtt = new G4VisAttributes(yellow);
503  DomeVisAtt->SetVisibility(false);
504  Dome_log->SetVisAttributes(DomeVisAtt);
505  }
506 
507  G4double domePosX = 0.0*mm;
508  G4double domePosY = 0.0*mm;
509  G4double domePosZ = -DomeSize.z()/2. + DomeThick/2.0 + fFilterPosition.z();
510  G4ThreeVector DomePosition(domePosX,domePosY,domePosZ);
511 
512  G4PVPlacement* Dome_phys =
513  new G4PVPlacement(0,DomePosition,"TelescopeDome",Dome_log,fTelescopeMother_phys,false,0);
514 
515 #if (G4VERSION_NUMBER>=800)
516  if (fVerbosityLevel > 2)
517  Dome_phys->CheckOverlaps(10000);
518 #endif
519 
520  // Construct lens supporting frame
521  G4double LensFrameX = fCorrectorRingFrameSize-fgTOLERANCE;
522  G4double LensFrameY = fCorrectorRingFrameSize-fgTOLERANCE;
523  G4double LensFrameZ = fCorrectorRingFrameThickness;
524  G4Box* LensFrameTemp = new G4Box("LensSupport",
525  LensFrameX/2.0,LensFrameY/2.0,LensFrameZ/2.0);
526 
527  // Diaphragm
528  G4double inRadius = 0.0*mm;
529  G4double ouRadius = fDiaphragmRadius + fgTOLERANCE;
531  G4double sPhi = 0.0*deg;
532  G4double dPhi = 360.0*deg;
533  G4Tubs* LensFrameHole = new G4Tubs("LensHole",inRadius,ouRadius,dZ/2.0,sPhi,dPhi);
534 
535  G4SubtractionSolid* LensFrameSolid =
536  new G4SubtractionSolid("LensFrame",LensFrameTemp,LensFrameHole);
537 
538  material = G4Material::GetMaterial("Aluminum");
539  G4LogicalVolume* LensFrame_log =
540  new G4LogicalVolume(LensFrameSolid,material,"LensFrame",0,0,0,true);
541  G4PVPlacement* LensFrame_phys =
542  new G4PVPlacement(0,G4ThreeVector(),"LensFrame",LensFrame_log,fLensFilterMother_phys,false,0);
543 #if (G4VERSION_NUMBER>=800)
544  if (fVerbosityLevel > 2)
545  LensFrame_phys->CheckOverlaps(10000);
546 #endif
547 
548  G4VisAttributes* LensVisAtt = new G4VisAttributes(gray);
549  LensVisAtt->SetVisibility(true);
550  LensFrame_log ->SetVisAttributes(LensVisAtt);
551 }
552 
553 
554 void
556 {
557 // Initialize detector geometry and material properties.
558 
559  string name = fTelescopeName.substr(0,2);
560  string eyeName;
561  if (name == "LL") eyeName = "Los Leones";
562  else
563  if (name == "LM") eyeName = "Los Morados";
564  else
565  if (name == "LA") eyeName = "Loma Amarilla";
566  else
567  if (name == "CO") eyeName = "Coihueco";
568 
569  istringstream buffer(fTelescopeName.substr(2,1));
570  int telID;
571  buffer >> telID;
572 
573  det::Detector& detector = det::Detector::GetInstance();
574 
575  const fdet::FDetector& detFD = detector.GetFDetector();
576  const fdet::Eye& eye = detFD.GetEye(eyeName);
577  const fdet::Telescope& telescope = eye.GetTelescope(telID);
578  const fdet::Corrector& corrector = telescope.GetCorrector();
579  const fdet::Diaphragm& diaphragm = telescope.GetDiaphragm();
580  const fdet::Camera& camera = telescope.GetCamera();
581 
582  fMaterials = new FDsimG4Materials(telescope);
583 
584  // Make the Geant4 world
585 
586  fWorld_x = 100.*CLHEP::m;
587  fWorld_y = 100.*CLHEP::m;
588  fWorld_z = 100.*CLHEP::m;
589  fWorldMaterial = G4Material::GetMaterial("Atmosphere");
590 
591  // We keep hardcoded the parameters of the telescope geometry for which no data
592  // seems to exhist in the Offline/FModelConfig
593 
594  fFilterRadius = 1100.0*CLHEP::mm;
597 
598  G4double posX = 0.0*CLHEP::mm;
599  G4double posY = 0.0*CLHEP::mm;
600  G4double posZ = 100.0*CLHEP::mm;
601 
602  fFilterPosition = G4ThreeVector(posX,posY,posZ);
603  fFilterMaterial = G4Material::GetMaterial("M-UG6");
604 
605  fFilterFrameSize = 2300.0*CLHEP::mm;
607 
609 
610  fHasCorrectorRing = telescope.HasCorrectorRing();
616  fCorrectorRingPosition = G4ThreeVector(0.0*CLHEP::mm,0.0*CLHEP::mm,0.0*CLHEP::mm);
617 
618  fCorrectorRingMaterial = G4Material::GetMaterial("BK7");
619 
621  fCameraNCols = telescope.GetLastColumn() - telescope.GetFirstColumn()+1; // 20
622  fCameraNRows = telescope.GetLastRow() - telescope.GetFirstRow()+1; // 22
623 
627 
628 // cout << fCameraRadiusCurvature/mm << " "
629 // << fMercedesHeight/mm << " "
630 // << fMercedesBase/mm << " "
631 // << fMercedesEfficiency << endl;
632 
633 #ifdef __TELESCOPE_SIMULATOR_LX_USE_ROOT__
634  fUseSensitiveDetectors = true;
635 #endif
636 
637  string profile = "KG";
638  fXMLManager.GetData(profile,"profile","lens");
639  fCorrectorRingProfile = profile;
640 
641  int numPetals = 1;
642  fXMLManager.GetData(numPetals,"numPetals","lens");
643  fCorrectorRingNumPetals = numPetals;
644 
645  double PetalSize = 360.0*utl::degree;
646  fXMLManager.GetData(PetalSize,"petalSize","lens");
648 
649  string MirrorDataFile;
650  if (!fXMLManager.GetData(MirrorDataFile,"mirrorDataFile","mirror"))
651  G4Exception ("Name of mirror data file not found in configuration file !");
652 
653  fMirrorDataFile = MirrorDataFile;
654 
655  double MirrorSegmentTiltSigma = 0.0*utl::degree;
656  fXMLManager.GetData(MirrorSegmentTiltSigma,"segmentTiltSigma","mirror");
657  fMirrorSegmentTiltSigma = MirrorSegmentTiltSigma/utl::degree * CLHEP::deg;
658 
659  string MirrorReflectionType = "Specular";
660  fXMLManager.GetData(MirrorReflectionType,"reflectionType","mirror");
661  fMirrorReflectionType = MirrorReflectionType;
662 
663  if (fMirrorReflectionType == "Diffusive"){
664 
665  double sigma_alpha = 0.0*utl::degree;
666  vector<double> wavelength;
667  vector<double> SpecularLobe;
668  vector<double> SpecularSpike;
669  vector<double> Backscatter;
670  vector<double> Reflectivity;
671 
673 
674  if (!fXMLManager.GetData(sigma_alpha,"sigma_alpha","mirror"))
676  if (!fXMLManager.GetData(wavelength,"wavelength","mirror"))
678  if (!fXMLManager.GetData(SpecularSpike,"specularspike","mirror"))
680  if (!fXMLManager.GetData(SpecularLobe,"specularlobe","mirror"))
682  if (!fXMLManager.GetData(Backscatter,"backscatter","mirror"))
684 
685  if (status == FDsimG4XMLManager::eNotFound)
686  G4Exception ("Mirror diffusive surface required but one or several parameters are missing.");
687 
689 
690  unsigned size = wavelength.size();
691  fMirrorRoughness.Wavelength = vector<G4double>(size);
692  fMirrorRoughness.SpecularSpike = vector<G4double>(size);
693  fMirrorRoughness.SpecularLobe = vector<G4double>(size);
694  fMirrorRoughness.Backscatter = vector<G4double>(size);
695 
696  for (unsigned i = 0 ; i < size ; i++) {
698  (wavelength[i]/utl::nanometer)*CLHEP::nm;
699  fMirrorRoughness.SpecularSpike[i] = SpecularSpike[i];
700  fMirrorRoughness.SpecularLobe[i] = SpecularLobe[i];
701  fMirrorRoughness.Backscatter[i] = Backscatter[i];
702  }
703  }
704 
705 
706  // Read filter frame surface parameters
707 
708  vector<double> wavelength;
709  vector<double> Reflectivity;
710 
712 
713  fXMLManager.GetData(Reflectivity,"reflectivity","filterframe");
714  fXMLManager.GetData(wavelength,"wavelength","filterframe");
715 
716  unsigned size = wavelength.size();
717 
718  fFilterFrameWavelength = vector<G4double>(size);
719  fFilterFrameReflectivity = vector<G4double>(size);
720  for (unsigned i = 0 ; i < size ; i++) {
722  (wavelength[i]/utl::nanometer)*CLHEP::nm;
723  fFilterFrameReflectivity[i] = Reflectivity[i];
724  }
725 
726  string FilterFrameReflectionType = "Specular";
727  fXMLManager.GetData(FilterFrameReflectionType,"reflectionType","filterframe");
728  fFilterFrameReflectionType = FilterFrameReflectionType;
729 
730  if (fFilterFrameReflectionType == "Diffusive") {
731  double sigma_alpha = 0.0*utl::degree;
732  vector<double> wavelength;
733  vector<double> SpecularLobe;
734  vector<double> SpecularSpike;
735  vector<double> Backscatter;
736 
737  status = FDsimG4XMLManager::eFound;
738 
739  if (!fXMLManager.GetData(sigma_alpha,"sigma_alpha","filterframe"))
741  if (!fXMLManager.GetData(wavelength,"wavelength","filterframe"))
743  if (!fXMLManager.GetData(SpecularSpike,"specularspike","filterframe"))
745  if (!fXMLManager.GetData(SpecularLobe,"specularlobe","filterframe"))
747  if (!fXMLManager.GetData(Backscatter,"backscatter","filterframe"))
749 
750  if (status == FDsimG4XMLManager::eNotFound)
751  G4Exception("Diffusive surface for filter frame required but one or several parameters are missing.");
752 
754 
755  size = wavelength.size();
756  fFilterFrameRoughness.Wavelength = vector<G4double>(size);
757  fFilterFrameRoughness.SpecularSpike = vector<G4double>(size);
758  fFilterFrameRoughness.SpecularLobe = vector<G4double>(size);
759  fFilterFrameRoughness.Backscatter = vector<G4double>(size);
760 
761  for (unsigned i = 0 ; i < size ; i++) {
763  (wavelength[i]/utl::nanometer)*CLHEP::nm;
764  fFilterFrameRoughness.SpecularSpike[i] = SpecularSpike[i];
765  fFilterFrameRoughness.SpecularLobe[i] = SpecularLobe[i];
766  fFilterFrameRoughness.Backscatter[i] = Backscatter[i];
767  }
768  }
769 
770  // Read photocathode parameters
771  {
772  string photocathodeMode = "reflectivity";
773  vector<double> wavelength;
774  vector<double> reflectivity;
775  vector<double> reRindex;
776  vector<double> imRindex;
777  double thickness = 0.0*utl::nanometer;
778 
781  fPhotocathodeUseThinFilm = false;
782 
783  fXMLManager.GetData(photocathodeMode,"mode","photocathode");
784 
785  if (photocathodeMode == "RefractiveIndex") {
787  }
788  else
790 
792  if (!fXMLManager.GetData(wavelength,"wavelength","photocathode") ||
793  !fXMLManager.GetData(reRindex,"realRindex","photocathode") ||
794  !fXMLManager.GetData(imRindex,"imaginaryRindex","photocathode"))
795 
796  G4Exception("RefractiveIndex mode chosen for photocathode but one or several parameters are missing.");
797  if (fXMLManager.GetData(thickness,"thickness","photocathode"))
799 
800  }
801 
803  if (!fXMLManager.GetData(reflectivity,"reflectivity","photocathode") ||
804  !fXMLManager.GetData(wavelength,"wavelength","photocathode"))
805  G4Exception("Reflectivity mode chosen for photocathode but one or several parameters are missing.");
806  }
807 
808  unsigned size = wavelength.size();
809  fPhotocathodeWavelength = vector<G4double>(size);
811  fPhotocathodeThickness = (thickness/utl::nanometer)*CLHEP::nm;
812  fPhotocathodeReRindex = vector<G4double>(size);
813  fPhotocathodeImRindex = vector<G4double>(size);
814 
815  for (unsigned i = 0 ; i < size ; i++) {
816  fPhotocathodeWavelength[i] = (wavelength[i]/utl::nanometer)*CLHEP::nm;
817  fPhotocathodeReRindex[i] = reRindex[i];
818  fPhotocathodeImRindex[i] = imRindex[i];
819  }
820  }
821  else {
822  fPhotocathodeReflectivity = vector<G4double>(size);
823  for (unsigned i = 0 ; i < size ; i++) {
824  fPhotocathodeWavelength[i] = (wavelength[i]/utl::nanometer)*CLHEP::nm;
825  fPhotocathodeReflectivity[i] = reflectivity[i];
826  }
827  }
828  }
829 }
constexpr double mm
Definition: AugerUnits.h:113
const Diaphragm & GetDiaphragm() const
Get the diaphragm that belongs to the telescope.
double GetInnerRadius() const
Inner radius of the ring.
Definition: Corrector.cc:78
const Corrector & GetCorrector() const
Get the Corrector (corrector ring) object that belongs to the telescope.
static const G4Colour gray(0.5, 0.5, 0.5)
double GetMercedesHeight() const
Height of the Mercedes.
const Eye & GetEye(const unsigned int eyeId) const
Find eye by numerical Id.
Definition: FDetector.cc:68
Detector description interface for Eye-related data.
Definition: FDetector/Eye.h:45
const Camera & GetCamera() const
Get the Camera object that belongs to the telescope.
unsigned int GetFirstColumn() const
Detector description interface for FDetector-related data.
Definition: FDetector.h:44
bool HasCorrectorRing() const
flag for corrector ring presence
double GetMercedesEfficiency() const
Average efficiency of the Mercedes.
constexpr double nanometer
Definition: AugerUnits.h:102
constexpr double deg
Definition: AugerUnits.h:140
FDsimG4DetectorConstruction(const G4String &TelescopeName, const FDsimG4XMLManager manager)
double GetOuterRadius() const
Outer radius of the ring.
Definition: Corrector.cc:86
const Telescope & GetTelescope(const unsigned int telescopeId) const
Find Telescope by numerical Id.
Top of the hierarchy of the detector description interface.
Definition: Detector.h:81
const fdet::FDetector & GetFDetector() const
Definition: Detector.cc:131
constexpr double degree
unsigned int GetLastRow() const
Description of a corrector ring.
Definition: Corrector.h:36
Status GetData(T &returnData, const std::string &componentProperty, const std::string &componentName) const
Return data as a TabulatedFunction.
G4LogicalVolume * GetLogicalVolume() const
Manager for specific FD description parameters in XML file.
double GetMercedesBase() const
Base of the Mercedes.
Detector description interface for Telescope-related data.
constexpr double cm
Definition: AugerUnits.h:117
double GetRadius() const
Radius of the diaphragm.
Definition: Diaphragm.cc:59
Description of the diaphragm.
Definition: Diaphragm.h:33
unsigned int GetLastColumn() const
double GetMeanLensThickness() const
Mean thickness of the lens.
Definition: Corrector.cc:94
double GetRadiusFocal() const
Radius of focal surface of the camera.
constexpr double m
Definition: AugerUnits.h:121
static const G4Colour yellow(1.0, 1.0, 0.5)
unsigned int GetFirstRow() const
Description of a camera.

, generated on Tue Sep 26 2023.