Geant4  10.02.p02
UltraFresnelLensParameterisation.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 // * 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 
46 #include "UltraFresnelLens.hh"
47 
48 #include "G4SystemOfUnits.hh"
49 #include "G4VPhysicalVolume.hh"
50 #include "G4ThreeVector.hh"
51 #include "G4Cons.hh"
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54 
56 {
57 
58  FresnelLens = Lens ;
59  GrooveWidth = Lens->GetGrooveWidth() ;
61 
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
69 {;}
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72 
74 (const G4int GrooveNo, G4VPhysicalVolume* physVol) const
75 {
76 
77  G4double Rmin1 = (GrooveNo+0)*(GrooveWidth) ;
78  G4double Rmax1 = (GrooveNo+1)*(GrooveWidth) ;
79 
80  G4double dZ = FresnelLens->GetSagita(Rmax1) - FresnelLens->GetSagita(Rmin1) ;
81 
82  if (dZ <= 0.0){
83  G4Exception("UltraFresnelLensParameterisation::ComputeTransformation()",
84  "AirSh003",FatalException,
85  "UltraFresnelLensParameterisation::ComputeTransformation: Groove depth<0 !");
86  }
87 
88  G4ThreeVector origin(0,0,(dZ-dZOffset)/2.);
89  physVol->SetTranslation(origin);
90  physVol->SetRotation(0);
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
96 (G4Cons& Groove, const G4int GrooveNo, const G4VPhysicalVolume*) const
97 {
98  G4double Rmin1 = (GrooveNo+0)*(GrooveWidth) ;
99  G4double Rmax1 = (GrooveNo+1)*(GrooveWidth) ;
100 
101  G4double Rmin2 = Rmin1 ;
102  G4double Rmax2 = Rmin2+0.0001*mm ;
103 
104  G4double dZ = FresnelLens->GetSagita(Rmax1) - FresnelLens->GetSagita(Rmin1) ;
105 
106  if (dZ <= 0.0){
107  G4Exception("UltraFresnelLensParameterisation::ComputeDimensions()",
108  "AirSh004",FatalException,
109  "UltraFresnelLensParameterisation::ComputeDimensions: Groove depth<0 !");
110  }
111 
112 
113  Groove.SetInnerRadiusMinusZ(Rmin1) ;
114  Groove.SetOuterRadiusMinusZ(Rmax1) ;
115 
116  Groove.SetInnerRadiusPlusZ(Rmin2) ;
117  Groove.SetOuterRadiusPlusZ(Rmax2) ;
118 
119  Groove.SetZHalfLength(dZ/2.) ;
120 
121 #ifdef ULTRA_VERBOSE
122 
123 G4cout << "UltraFresnelLensParameterisation: GrooveNo " << GrooveNo+1 <<
124 " Rmin1, Rmax1(mm): " << Rmin1/mm <<" "<< Rmax1/mm << " dZ(mm) " << dZ/mm << G4endl ;
125 #endif
126 
127 
128 
129 }
130 
void SetZHalfLength(G4double newDz)
void SetInnerRadiusMinusZ(G4double Rmin1)
void ComputeDimensions(G4Cons &Groove, const G4int GrooveNo, const G4VPhysicalVolume *physVol) const
CLHEP::Hep3Vector G4ThreeVector
G4double GetSagita(G4double)
void SetOuterRadiusPlusZ(G4double Rmax2)
int G4int
Definition: G4Types.hh:78
void SetRotation(G4RotationMatrix *)
G4GLOB_DLL std::ostream G4cout
G4double GetGrooveWidth()
Definition: G4Cons.hh:83
void SetTranslation(const G4ThreeVector &v)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void ComputeTransformation(const G4int GrooveNo, G4VPhysicalVolume *physVol) const
#define G4endl
Definition: G4ios.hh:61
void SetOuterRadiusMinusZ(G4double Rmax1)
double G4double
Definition: G4Types.hh:76
void SetInnerRadiusPlusZ(G4double Rmin2)
static const double mm
Definition: G4SIunits.hh:114