source: trunk/examples/advanced/air_shower/src/UltraFresnelLensParameterisation.cc @ 810

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

update

File size: 4.7 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//    *      UltraFresnelLensParameterisation.cc
36//    ****************************************************
37//
38//    Class derived from G4VPVParameterisation and used to define a Fresnel lens geometry
39//    through a parameterised replication of G4Cons volumes. These volumes are frustra
40//    of cones describing the lens grooves.
41//    An  UltraFresnelLensParameterisation object is created in the UltraFresnelLens class
42//
43#include <cmath>
44#include "UltraFresnelLensParameterisation.hh"
45#include "UltraFresnelLens.hh"
46
47#include "G4VPhysicalVolume.hh"
48#include "G4ThreeVector.hh"
49#include "G4Cons.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
53UltraFresnelLensParameterisation::UltraFresnelLensParameterisation(UltraFresnelLens* Lens) 
54{
55   
56   FresnelLens     = Lens ;
57   GrooveWidth     = Lens->GetGrooveWidth() ;
58   NumberOfGrooves = Lens->GetNumberOfGrooves() ;
59
60   dZOffset = Lens->GetSagita((NumberOfGrooves-0)*(GrooveWidth)) - 
61              Lens->GetSagita((NumberOfGrooves-1)*(GrooveWidth)) ;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
66UltraFresnelLensParameterisation::~UltraFresnelLensParameterisation()
67{;}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
70
71void UltraFresnelLensParameterisation::ComputeTransformation
72(const G4int GrooveNo, G4VPhysicalVolume* physVol) const
73{
74
75  G4double  Rmin1 = (GrooveNo+0)*(GrooveWidth) ;
76  G4double  Rmax1 = (GrooveNo+1)*(GrooveWidth) ;
77
78  G4double     dZ = FresnelLens->GetSagita(Rmax1) - FresnelLens->GetSagita(Rmin1) ;
79
80  if (dZ <= 0.0){
81     G4Exception("UltraFresnelLensParameterisation::ComputeTransformation: Groove depth<0 !");
82  }
83 
84  G4ThreeVector origin(0,0,(dZ-dZOffset)/2.);
85  physVol->SetTranslation(origin);
86  physVol->SetRotation(0);
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
91void UltraFresnelLensParameterisation::ComputeDimensions
92(G4Cons& Groove, const G4int GrooveNo, const G4VPhysicalVolume*) const         
93{
94  G4double  Rmin1 = (GrooveNo+0)*(GrooveWidth) ;
95  G4double  Rmax1 = (GrooveNo+1)*(GrooveWidth) ;
96
97  G4double  Rmin2 = Rmin1 ;
98  G4double  Rmax2 = Rmin2+0.0001*mm ;
99 
100  G4double  dZ = FresnelLens->GetSagita(Rmax1) - FresnelLens->GetSagita(Rmin1) ;
101
102  if (dZ <= 0.0){
103     G4Exception("UltraFresnelLensParameterisation::ComputeDimensions: Groove depth<0 !");
104  }
105
106
107  Groove.SetInnerRadiusMinusZ(Rmin1) ;
108  Groove.SetOuterRadiusMinusZ(Rmax1) ;
109
110  Groove.SetInnerRadiusPlusZ(Rmin2) ;
111  Groove.SetOuterRadiusPlusZ(Rmax2) ;
112
113  Groove.SetZHalfLength(dZ/2.) ;
114
115#ifdef ULTRA_VERBOSE
116
117G4cout << "UltraFresnelLensParameterisation: GrooveNo " << GrooveNo+1 << 
118" Rmin1, Rmax1(mm): " << Rmin1/mm <<" "<<  Rmax1/mm << " dZ(mm) " << dZ/mm << G4endl ;
119#endif
120
121
122
123}
124
Note: See TracBrowser for help on using the repository browser.