Geant4  10.00.p02
UltraFresnelLens.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 // --------------------------------------------------------------
28 // GEANT 4 - ULTRA experiment example
29 // --------------------------------------------------------------
30 //
31 // Code developed by:
32 // B. Tome, M.C. Espirito-Santo, A. Trindade, P. Rodrigues
33 //
34 // ****************************************************
35 // * UltraFresnelLens.cc
36 // ****************************************************
37 //
38 // Class for definition of the Ultra Fresnel Lens.
39 // An UltraFresnelLens object is created in the UltraDetectorConstruction class.
40 // This class makes use of the UltraFresnelLensParameterisation class.
41 // The lens profile is define through the GetSagita method.
42 //
43 #include <cmath>
44 #include "UltraFresnelLens.hh"
46 
47 #include "G4PhysicalConstants.hh"
48 #include "G4SystemOfUnits.hh"
49 #include "G4Material.hh"
50 #include "G4MaterialTable.hh"
51 #include "G4Tubs.hh"
52 #include "G4Cons.hh"
53 #include "G4LogicalVolume.hh"
54 #include "G4PVPlacement.hh"
55 #include "G4PVParameterised.hh"
56 #include "G4ThreeVector.hh"
57 #include "G4VisAttributes.hh"
58 
60 {
61 
62  LensMaterial = Material ;
63  LensDiameter = Diameter ;
64  NumberOfGrooves = nGrooves;
65  GrooveWidth = (Diameter/2.0)/nGrooves ;
66  LensPosition = Pos ;
67 
68  if( GrooveWidth <= 0 ){
69  G4Exception("UltraFresnelLens::UltraFresnelLens()","AirSh001",FatalException,
70  "UltraFresnelLens constructor: GrooveWidth<=0");
71  }
72 
73 
74 G4double Rmin1 = (NumberOfGrooves-1)*(GrooveWidth) ;
75 G4double Rmax1 = (NumberOfGrooves-0)*(GrooveWidth) ;
76 LensThickness = GetSagita(Rmax1)-GetSagita(Rmin1) ; // Height of the highest groove
77 
78 BuildLens(MotherPV) ;
79 
80 }
81 
82 
84 {;}
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 
90 G4double StartPhi = 0.0 ;
91 G4double DeltaPhi = twopi ;
92 
93 G4double LensMotherDz = LensThickness ;
94 
95 
96 G4Tubs *LensMotherCylinder
97  = new G4Tubs("LensMotherCylinder",0.0*mm,LensDiameter/2.0,LensMotherDz/2.0,StartPhi,DeltaPhi);
98 
99 G4Material *LensMotherMaterial = MotherPV->GetLogicalVolume()->GetMaterial() ;
100 G4LogicalVolume *LensMotherLV
101  = new G4LogicalVolume(LensMotherCylinder,LensMotherMaterial,"LensMotherLV",0,0,0);
102 G4VPhysicalVolume *LensMotherPV
103  = new G4PVPlacement(0,LensPosition,"LensMotherPV",LensMotherLV,MotherPV,false,0);
104 
106 
107 
108 G4Cons *solidGroove
109  = new G4Cons("Groove",40.0*mm,50.0*mm,40.0*mm,40.001*mm,1.0*mm,StartPhi,DeltaPhi);
110 
111 G4LogicalVolume *logicalGroove
112  = new G4LogicalVolume(solidGroove,LensMaterial,"Groove_log",0,0,0);
113 
114 G4VPVParameterisation *FresnelLensParam = new UltraFresnelLensParameterisation(this);
115 
117  new G4PVParameterised("LensPV",logicalGroove,LensMotherPV,kZAxis,NumberOfGrooves,FresnelLensParam) ;
118 
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
124 {
125 
126  G4double Conic = -1.0;
127  G4double Curvature = 0.00437636761488/mm ;
128  G4double Aspher[8] = { 4.206739256e-05/(mm),
129  9.6440152e-10/(mm3),
130  -1.4884317e-15/(mm2*mm3),
131  0.0/(mm*mm3*mm3),
132  0.0/(mm3*mm3*mm3),
133  0.0/(mm2*mm3*mm3*mm3),
134  0.0/(mm*mm3*mm3*mm3*mm3),
135  0.0/(mm*3*mm3*mm3*mm3*mm3)
136  };
137 
138  G4double TotAspher = 0.0*mm ;
139 
140  for(G4int k=1;k<9;k++){
141  TotAspher += Aspher[k-1]*std::pow(radius,2*k) ;
142  }
143 
144  G4double ArgSqrt = 1.0-(1.0+Conic)*std::pow(Curvature,2)*std::pow(radius,2) ;
145 
146  if (ArgSqrt < 0.0){
147  G4Exception("UltraFresnelLens::GetSagita()","AirSh002",
149  "UltraFresnelLensParameterisation::Sagita: Square Root of <0 !");
150  }
151  G4double Sagita_value = Curvature*std::pow(radius,2)/(1.0+std::sqrt(ArgSqrt)) + TotAspher;
152 
153  return Sagita_value ;
154 
155 
156 }
157 
158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
159 
160 
G4ThreeVector LensPosition
CLHEP::Hep3Vector G4ThreeVector
G4double GetSagita(G4double)
G4Material * GetMaterial() const
UltraFresnelLens(G4double Diameter, G4int nGrooves, G4Material *Material, G4VPhysicalVolume *MotherPV, G4ThreeVector Pos)
Definition: G4Tubs.hh:85
G4VPhysicalVolume * LensPhysicalVolume
int G4int
Definition: G4Types.hh:78
G4Material * LensMaterial
void BuildLens(G4VPhysicalVolume *)
Definition: G4Cons.hh:83
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double mm2
Definition: G4SIunits.hh:103
G4LogicalVolume * GetLogicalVolume() const
static const G4VisAttributes Invisible
double G4double
Definition: G4Types.hh:76
static const double mm
Definition: G4SIunits.hh:102
void SetVisAttributes(const G4VisAttributes *pVA)
static const double mm3
Definition: G4SIunits.hh:104