source: trunk/examples/extended/optical/LXe/src/LXeWLSFiber.cc @ 1009

Last change on this file since 1009 was 807, checked in by garnier, 16 years ago

update

File size: 3.8 KB
Line 
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#include "LXeWLSFiber.hh"
27#include "globals.hh"
28#include "G4LogicalSkinSurface.hh"
29#include "G4LogicalBorderSurface.hh"
30
31G4LogicalVolume* LXeWLSFiber::clad2_log=NULL;
32
33LXeWLSFiber::LXeWLSFiber(G4RotationMatrix *pRot,
34                             const G4ThreeVector &tlate,
35                             G4LogicalVolume *pMotherLogical,
36                             G4bool pMany,
37                             G4int pCopyNo,
38                             LXeDetectorConstruction* c)
39  :G4PVPlacement(pRot,tlate,
40                 new G4LogicalVolume(new G4Box("temp",1,1,1),
41                                     G4Material::GetMaterial("Vacuum"),
42                                     "temp",0,0,0),
43                 "Cladding2",pMotherLogical,pMany,pCopyNo),constructor(c)
44{
45  CopyValues();
46 
47  if(!clad2_log || updated){
48    // The Fiber
49    //
50    G4Tubs* Fiber_tube = 
51      new G4Tubs("Fiber",fiber_rmin,fiber_rmax,fiber_z,fiber_sphi,fiber_ephi);
52   
53    G4LogicalVolume* Fiber_log = 
54      new G4LogicalVolume(Fiber_tube,G4Material::GetMaterial("PMMA"),
55                          "Fiber",0,0,0);
56   
57    // Cladding (first layer)
58    //
59    G4Tubs* clad1_tube = 
60      new G4Tubs("Cladding1",clad1_rmin,clad1_rmax,clad1_z,clad1_sphi,
61                 clad1_ephi);
62   
63    G4LogicalVolume* clad1_log = 
64      new G4LogicalVolume(clad1_tube,G4Material::GetMaterial("Pethylene"),
65                          "Cladding1",0,0,0);
66     
67    // Cladding (second layer)
68    //
69    G4Tubs* clad2_tube = 
70      new G4Tubs("Cladding2",clad2_rmin,clad2_rmax,clad2_z,clad2_sphi,
71                 clad2_ephi);
72   
73    clad2_log = 
74      new G4LogicalVolume(clad2_tube,G4Material::GetMaterial("fPethylene"),
75                          "Cladding2",0,0,0);
76   
77    new G4PVPlacement(0,G4ThreeVector(0.,0.,0.),Fiber_log,
78                      "Fiber", clad1_log,false,0);
79    new G4PVPlacement(0,G4ThreeVector(0.,0.,0.),clad1_log,
80                      "Cladding1",clad2_log,false,0);
81  }
82 
83  SetLogicalVolume(clad2_log);
84}
85
86void LXeWLSFiber::CopyValues(){
87  updated=constructor->GetUpdated();
88
89  fiber_rmin = 0.00*cm;   
90  fiber_rmax = 0.10*cm;   
91  fiber_z    = constructor->GetScintX()/2;
92  fiber_sphi = 0.00*deg;
93  fiber_ephi = 360.*deg;
94 
95  clad1_rmin = 0.;// fiber_rmax;
96    clad1_rmax = fiber_rmax + 0.015*fiber_rmax;   
97 
98  clad1_z    = fiber_z;
99  clad1_sphi = fiber_sphi;
100  clad1_ephi = fiber_ephi; 
101 
102  clad2_rmin = 0.;//clad1_rmax;
103  clad2_rmax = clad1_rmax + 0.015*fiber_rmax;   
104
105  clad2_z    = fiber_z;
106  clad2_sphi = fiber_sphi;
107  clad2_ephi = fiber_ephi;
108
109}
Note: See TracBrowser for help on using the repository browser.