Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
UltraFresnelLensParameterisation Class Reference

#include <UltraFresnelLensParameterisation.hh>

Inheritance diagram for UltraFresnelLensParameterisation:
G4VPVParameterisation

Public Member Functions

 UltraFresnelLensParameterisation (UltraFresnelLens *)
 
virtual ~UltraFresnelLensParameterisation ()
 
void ComputeTransformation (const G4int GrooveNo, G4VPhysicalVolume *physVol) const
 
void ComputeDimensions (G4Cons &Groove, const G4int GrooveNo, const G4VPhysicalVolume *physVol) const
 
- Public Member Functions inherited from G4VPVParameterisation
 G4VPVParameterisation ()
 
virtual ~G4VPVParameterisation ()
 
virtual G4VSolidComputeSolid (const G4int, G4VPhysicalVolume *)
 
virtual G4MaterialComputeMaterial (const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=0)
 
virtual G4bool IsNested () const
 
virtual G4VVolumeMaterialScannerGetMaterialScanner ()
 

Detailed Description

Definition at line 70 of file UltraFresnelLensParameterisation.hh.

Constructor & Destructor Documentation

UltraFresnelLensParameterisation::UltraFresnelLensParameterisation ( UltraFresnelLens Lens)

Definition at line 55 of file UltraFresnelLensParameterisation.cc.

References UltraFresnelLens::GetGrooveWidth(), UltraFresnelLens::GetNumberOfGrooves(), and UltraFresnelLens::GetSagita().

56 {
57 
58  FresnelLens = Lens ;
59  GrooveWidth = Lens->GetGrooveWidth() ;
60  NumberOfGrooves = Lens->GetNumberOfGrooves() ;
61 
62  dZOffset = Lens->GetSagita((NumberOfGrooves-0)*(GrooveWidth)) -
63  Lens->GetSagita((NumberOfGrooves-1)*(GrooveWidth)) ;
64 }
G4double GetSagita(G4double)
G4double GetGrooveWidth()
UltraFresnelLensParameterisation::~UltraFresnelLensParameterisation ( )
virtual

Definition at line 68 of file UltraFresnelLensParameterisation.cc.

69 {;}

Member Function Documentation

void UltraFresnelLensParameterisation::ComputeDimensions ( G4Cons Groove,
const G4int  GrooveNo,
const G4VPhysicalVolume physVol 
) const
virtual

Reimplemented from G4VPVParameterisation.

Definition at line 96 of file UltraFresnelLensParameterisation.cc.

References FatalException, G4cout, G4endl, G4Exception(), python.hepunit::mm, G4Cons::SetInnerRadiusMinusZ(), G4Cons::SetInnerRadiusPlusZ(), G4Cons::SetOuterRadiusMinusZ(), G4Cons::SetOuterRadiusPlusZ(), and G4Cons::SetZHalfLength().

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 }
void SetZHalfLength(G4double newDz)
void SetInnerRadiusMinusZ(G4double Rmin1)
G4double GetSagita(G4double)
void SetOuterRadiusPlusZ(G4double Rmax2)
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
void SetOuterRadiusMinusZ(G4double Rmax1)
double G4double
Definition: G4Types.hh:76
void SetInnerRadiusPlusZ(G4double Rmin2)
void UltraFresnelLensParameterisation::ComputeTransformation ( const G4int  GrooveNo,
G4VPhysicalVolume physVol 
) const
virtual

Implements G4VPVParameterisation.

Definition at line 74 of file UltraFresnelLensParameterisation.cc.

References FatalException, G4Exception(), G4VPhysicalVolume::SetRotation(), and G4VPhysicalVolume::SetTranslation().

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 }
G4double GetSagita(G4double)
void SetRotation(G4RotationMatrix *)
void SetTranslation(const G4ThreeVector &v)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double G4double
Definition: G4Types.hh:76

The documentation for this class was generated from the following files: