Loading...
Searching...
No Matches
LXeWLSFiber.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/// \file optical/LXe/src/LXeWLSFiber.cc
28/// \brief Implementation of the LXeWLSFiber class
29//
30//
31#include "LXeWLSFiber.hh"
32
33#include "globals.hh"
34#include "G4Box.hh"
35#include "G4LogicalVolume.hh"
36#include "G4Material.hh"
37#include "G4SystemOfUnits.hh"
38#include "G4Tubs.hh"
39
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
43
44LXeWLSFiber::LXeWLSFiber(G4RotationMatrix* pRot, const G4ThreeVector& tlate,
45 G4LogicalVolume* pMotherLogical, G4bool pMany,
46 G4int pCopyNo, LXeDetectorConstruction* c)
47 : G4PVPlacement(pRot, tlate,
48 new G4LogicalVolume(new G4Box("temp", 1., 1., 1.),
49 G4Material::GetMaterial("Vacuum"), "temp"),
50 "Cladding2", pMotherLogical, pMany, pCopyNo)
51 , fConstructor(c)
52{
53 CopyValues();
54
55 // The Fiber
56 //
57 auto fiber_tube = new G4Tubs("Fiber", fFiber_rmin, fFiber_rmax, fFiber_z,
59
60 auto fiber_log = new G4LogicalVolume(
61 fiber_tube, G4Material::GetMaterial("PMMA"), "Fiber");
62
63 // Cladding (first layer)
64 //
65 auto clad1_tube = new G4Tubs("Cladding1", fClad1_rmin, fClad1_rmax,
67
68 auto clad1_log = new G4LogicalVolume(
69 clad1_tube, G4Material::GetMaterial("Pethylene1"), "Cladding1");
70
71 // Cladding (second layer)
72 //
73 auto clad2_tube = new G4Tubs("Cladding2", fClad2_rmin, fClad2_rmax,
75
77 clad2_tube, G4Material::GetMaterial("Pethylene2"), "Cladding2");
78
79 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), fiber_log, "Fiber",
80 clad1_log, false, 0);
81 new G4PVPlacement(nullptr, G4ThreeVector(0., 0., 0.), clad1_log, "Cladding1",
82 fClad2_log, false, 0);
83
84 SetLogicalVolume(fClad2_log);
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90{
91 fFiber_rmin = 0.0 * cm;
92 fFiber_rmax = 0.1 * cm;
94 fFiber_sphi = 0.0 * deg;
95 fFiber_ephi = 360. * deg;
96
97 fClad1_rmin = 0.; // fFiber_rmax;
99
103
104 fClad2_rmin = 0.; // fClad1_rmax;
106
110}
Definition of the LXeWLSFiber class.
G4double fFiber_ephi
G4double fClad1_ephi
G4double fFiber_rmin
G4double fClad2_ephi
G4double fFiber_z
G4double fClad1_rmax
G4double fClad2_sphi
void CopyValues()
G4double fClad1_rmin
G4double fFiber_sphi
G4double fClad1_z
G4double fClad2_z
G4double fClad1_sphi
LXeWLSFiber(G4RotationMatrix *pRot, const G4ThreeVector &tlate, G4LogicalVolume *pMotherLogical, G4bool pMany, G4int pCopyNo, LXeDetectorConstruction *c)
LXeDetectorConstruction * fConstructor
G4double fClad2_rmin
G4double fClad2_rmax
G4double fFiber_rmax
static G4LogicalVolume * fClad2_log

Applications | User Support | Publications | Collaboration