File indexing completed on 2025-02-23 09:22:22
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031 #include "LXeWLSFiber.hh"
0032
0033 #include "G4Box.hh"
0034 #include "G4LogicalVolume.hh"
0035 #include "G4Material.hh"
0036 #include "G4SystemOfUnits.hh"
0037 #include "G4Tubs.hh"
0038 #include "globals.hh"
0039
0040 G4LogicalVolume* LXeWLSFiber::fClad2_log = nullptr;
0041
0042
0043
0044 LXeWLSFiber::LXeWLSFiber(G4RotationMatrix* pRot, const G4ThreeVector& tlate,
0045 G4LogicalVolume* pMotherLogical, G4bool pMany, G4int pCopyNo,
0046 LXeDetectorConstruction* c)
0047 : G4PVPlacement(
0048 pRot, tlate,
0049 new G4LogicalVolume(new G4Box("temp", 1., 1., 1.), G4Material::GetMaterial("Vacuum"), "temp"),
0050 "Cladding2", pMotherLogical, pMany, pCopyNo),
0051 fConstructor(c)
0052 {
0053 CopyValues();
0054
0055
0056
0057 auto fiber_tube =
0058 new G4Tubs("Fiber", fFiber_rmin, fFiber_rmax, fFiber_z, fFiber_sphi, fFiber_ephi);
0059
0060 auto fiber_log = new G4LogicalVolume(fiber_tube, G4Material::GetMaterial("PMMA"), "Fiber");
0061
0062
0063
0064 auto clad1_tube =
0065 new G4Tubs("Cladding1", fClad1_rmin, fClad1_rmax, fClad1_z, fClad1_sphi, fClad1_ephi);
0066
0067 auto clad1_log =
0068 new G4LogicalVolume(clad1_tube, G4Material::GetMaterial("Pethylene1"), "Cladding1");
0069
0070
0071
0072 auto clad2_tube =
0073 new G4Tubs("Cladding2", fClad2_rmin, fClad2_rmax, fClad2_z, fClad2_sphi, fClad2_ephi);
0074
0075 fClad2_log = new G4LogicalVolume(clad2_tube, G4Material::GetMaterial("Pethylene2"), "Cladding2");
0076
0077 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), fiber_log, "Fiber", clad1_log, false, 0);
0078 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), clad1_log, "Cladding1", fClad2_log, false,
0079 0);
0080
0081 SetLogicalVolume(fClad2_log);
0082 }
0083
0084
0085
0086 void LXeWLSFiber::CopyValues()
0087 {
0088 fFiber_rmin = 0.0 * cm;
0089 fFiber_rmax = 0.1 * cm;
0090 fFiber_z = fConstructor->GetScintX() / 2.;
0091 fFiber_sphi = 0.0 * deg;
0092 fFiber_ephi = 360. * deg;
0093
0094 fClad1_rmin = 0.;
0095 fClad1_rmax = fFiber_rmax + 0.015 * fFiber_rmax;
0096
0097 fClad1_z = fFiber_z;
0098 fClad1_sphi = fFiber_sphi;
0099 fClad1_ephi = fFiber_ephi;
0100
0101 fClad2_rmin = 0.;
0102 fClad2_rmax = fClad1_rmax + 0.015 * fFiber_rmax;
0103
0104 fClad2_z = fFiber_z;
0105 fClad2_sphi = fFiber_sphi;
0106 fClad2_ephi = fFiber_ephi;
0107 }