source: trunk/examples/advanced/air_shower/src/UltraFresnelLens.cc @ 1303

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

update

File size: 5.6 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//
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"
45#include "UltraFresnelLensParameterisation.hh"
46
47#include "G4Material.hh"
48#include "G4MaterialTable.hh"
49#include "G4Tubs.hh"
50#include "G4Cons.hh"
51#include "G4LogicalVolume.hh"
52#include "G4PVPlacement.hh"
53#include "G4PVParameterised.hh"
54#include "G4ThreeVector.hh"
55#include "G4VisAttributes.hh"
56
57UltraFresnelLens::UltraFresnelLens( G4double Diameter, G4int nGrooves, G4Material* Material, G4VPhysicalVolume * MotherPV, G4ThreeVector Pos)
58{
59
60 LensMaterial    = Material ;
61 LensDiameter    = Diameter ;
62 NumberOfGrooves = nGrooves;
63 GrooveWidth     = (Diameter/2.0)/nGrooves ;
64 LensPosition    = Pos ;
65
66 if( GrooveWidth <= 0 ){
67       G4Exception("UltraFresnelLens constructor: GrooveWidth<=0");
68   }
69
70
71G4double  Rmin1 = (NumberOfGrooves-1)*(GrooveWidth) ;
72G4double  Rmax1 = (NumberOfGrooves-0)*(GrooveWidth) ;
73LensThickness   = GetSagita(Rmax1)-GetSagita(Rmin1) ; // Height of the highest groove
74
75BuildLens(MotherPV) ;
76
77}
78
79
80UltraFresnelLens::~UltraFresnelLens(  )
81{;}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
85void UltraFresnelLens::BuildLens(G4VPhysicalVolume *MotherPV){
86
87G4double StartPhi = 0.0 ;
88G4double DeltaPhi = twopi ;
89
90G4double  LensMotherDz = LensThickness ;
91
92
93G4Tubs *LensMotherCylinder
94    = new G4Tubs("LensMotherCylinder",0.0*mm,LensDiameter/2.0,LensMotherDz/2.0,StartPhi,DeltaPhi);
95
96G4Material *LensMotherMaterial = MotherPV->GetLogicalVolume()->GetMaterial() ;
97G4LogicalVolume *LensMotherLV
98    = new G4LogicalVolume(LensMotherCylinder,LensMotherMaterial,"LensMotherLV",0,0,0);
99G4VPhysicalVolume *LensMotherPV
100    = new G4PVPlacement(0,LensPosition,"LensMotherPV",LensMotherLV,MotherPV,false,0);
101
102LensMotherLV->SetVisAttributes (G4VisAttributes::Invisible);
103
104
105G4Cons *solidGroove
106  = new G4Cons("Groove",40.0*mm,50.0*mm,40.0*mm,40.001*mm,1.0*mm,StartPhi,DeltaPhi);
107
108G4LogicalVolume *logicalGroove
109    = new G4LogicalVolume(solidGroove,LensMaterial,"Groove_log",0,0,0);
110
111G4VPVParameterisation *FresnelLensParam = new UltraFresnelLensParameterisation(this);
112
113LensPhysicalVolume = 
114 new G4PVParameterised("LensPV",logicalGroove,LensMotherPV,kZAxis,NumberOfGrooves,FresnelLensParam) ;
115
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119
120G4double UltraFresnelLens::GetSagita(G4double radius) 
121{
122
123  G4double Conic = -1.0;
124  G4double Curvature = 0.00437636761488/mm ;
125  G4double Aspher[8] = {      4.206739256e-05/(mm),
126                                9.6440152e-10/(mm3),
127                               -1.4884317e-15/(mm2*mm3),
128                                          0.0/(mm*mm3*mm3),
129                                          0.0/(mm3*mm3*mm3),
130                                          0.0/(mm2*mm3*mm3*mm3),
131                                          0.0/(mm*mm3*mm3*mm3*mm3),
132                                          0.0/(mm*3*mm3*mm3*mm3*mm3)
133                            };
134
135  G4double TotAspher = 0.0*mm ;
136
137  for(G4int k=1;k<9;k++){
138    TotAspher += Aspher[k-1]*std::pow(radius,2*k) ;
139  }
140
141  G4double ArgSqrt = 1.0-(1.0+Conic)*std::pow(Curvature,2)*std::pow(radius,2) ; 
142
143  if (ArgSqrt < 0.0){
144     G4Exception("UltraFresnelLensParameterisation::Sagita: Square Root of <0 !");
145  }
146  G4double Sagita_value = Curvature*std::pow(radius,2)/(1.0+std::sqrt(ArgSqrt)) + TotAspher;
147
148  return Sagita_value ;
149
150                             
151}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
154
155
Note: See TracBrowser for help on using the repository browser.