Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/air_shower/src/UltraFresnelLens.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 "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 
 59 UltraFresnelLens::UltraFresnelLens(G4double Diameter, G4int nGrooves, G4Material* Material, G4VPhysicalVolume * MotherPV, G4ThreeVector Pos)
 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 auto  Rmin1 = (NumberOfGrooves-1)*(GrooveWidth) ;
 75 auto  Rmax1 = (NumberOfGrooves-0)*(GrooveWidth) ;
 76 LensThickness   = GetSagita(Rmax1)-GetSagita(Rmin1) ; // Height of the highest groove
 77 
 78 BuildLens(MotherPV) ;
 79 
 80 }
 81 
 82 
 83 UltraFresnelLens::~UltraFresnelLens(  )
 84 {;}
 85 
 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 87 
 88 void UltraFresnelLens::BuildLens(G4VPhysicalVolume *MotherPV){
 89 
 90 auto StartPhi = 0.0 ;
 91 auto DeltaPhi = twopi ;
 92 
 93 auto  LensMotherDz = LensThickness ;
 94 
 95 
 96 auto LensMotherCylinder
 97     = new G4Tubs("LensMotherCylinder",0.0*mm,LensDiameter/2.0,LensMotherDz/2.0,StartPhi,DeltaPhi);
 98 
 99 auto LensMotherMaterial = MotherPV->GetLogicalVolume()->GetMaterial() ;
100 auto LensMotherLV
101     = new G4LogicalVolume(LensMotherCylinder,LensMotherMaterial,"LensMotherLV",0,0,0);
102 auto LensMotherPV
103     = new G4PVPlacement(0,LensPosition,"LensMotherPV",LensMotherLV,MotherPV,false,0);
104 
105 LensMotherLV->SetVisAttributes (G4VisAttributes::GetInvisible());
106 
107 
108 auto solidGroove
109   = new G4Cons("Groove",40.0*mm,50.0*mm,40.0*mm,40.001*mm,1.0*mm,StartPhi,DeltaPhi);
110 
111 auto logicalGroove
112     = new G4LogicalVolume(solidGroove,LensMaterial,"Groove_log",0,0,0);
113 
114 auto FresnelLensParam = new UltraFresnelLensParameterisation(this);
115 
116 LensPhysicalVolume =
117  new G4PVParameterised("LensPV",logicalGroove,LensMotherPV,kZAxis,NumberOfGrooves,FresnelLensParam) ;
118 
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122 
123 G4double UltraFresnelLens::GetSagita(G4double radius)
124 {
125   auto Conic = -1.0;
126   auto Curvature = 0.00437636761488/mm ;
127   G4double Aspher[8] = {4.206739256e-05/(mm),
128                         9.6440152e-10/(mm3),
129                        -1.4884317e-15/(mm2*mm3),
130                         0.0/(mm*mm3*mm3),
131                         0.0/(mm3*mm3*mm3),
132                         0.0/(mm2*mm3*mm3*mm3),
133                         0.0/(mm*mm3*mm3*mm3*mm3),
134                         0.0/(mm*3*mm3*mm3*mm3*mm3)
135                             };
136 
137   auto TotAspher = 0.0*mm ;
138 
139   for(G4int k=1;k<9;k++){
140     TotAspher += Aspher[k-1]*std::pow(radius,2*k) ;
141   }
142 
143   auto ArgSqrt = 1.0-(1.0+Conic)*std::pow(Curvature,2)*std::pow(radius,2) ;
144 
145   if (ArgSqrt < 0.0){
146     G4Exception("UltraFresnelLens::GetSagita()","AirSh002",
147     FatalException,
148     "UltraFresnelLensParameterisation::Sagita: Square Root of <0 !");
149   }
150   auto Sagita_value = Curvature*std::pow(radius,2)/(1.0+std::sqrt(ArgSqrt)) + TotAspher;
151 
152   return Sagita_value ;
153 
154 
155 }
156