Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /geant4/examples/advanced/air_shower/src/UltraFresnelLens.cc was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // --------------------------------------------------------------
0028 //                 GEANT 4 - ULTRA experiment example
0029 // --------------------------------------------------------------
0030 //
0031 // Code developed by:
0032 // B. Tome, M.C. Espirito-Santo, A. Trindade, P. Rodrigues
0033 //
0034 //    ****************************************************
0035 //    *      UltraFresnelLens.cc
0036 //    ****************************************************
0037 //
0038 //    Class for definition of the Ultra Fresnel Lens.
0039 //    An  UltraFresnelLens object is created in the UltraDetectorConstruction class.
0040 //    This class makes use of the UltraFresnelLensParameterisation class.
0041 //    The lens profile is define through the GetSagita method.
0042 //
0043 #include <cmath>
0044 #include "UltraFresnelLens.hh"
0045 #include "UltraFresnelLensParameterisation.hh"
0046 
0047 #include "G4PhysicalConstants.hh"
0048 #include "G4SystemOfUnits.hh"
0049 #include "G4Material.hh"
0050 #include "G4MaterialTable.hh"
0051 #include "G4Tubs.hh"
0052 #include "G4Cons.hh"
0053 #include "G4LogicalVolume.hh"
0054 #include "G4PVPlacement.hh"
0055 #include "G4PVParameterised.hh"
0056 #include "G4ThreeVector.hh"
0057 #include "G4VisAttributes.hh"
0058 
0059 UltraFresnelLens::UltraFresnelLens(G4double Diameter, G4int nGrooves, G4Material* Material, G4VPhysicalVolume * MotherPV, G4ThreeVector Pos)
0060 {
0061 
0062  LensMaterial    = Material ;
0063  LensDiameter    = Diameter ;
0064  NumberOfGrooves = nGrooves;
0065  GrooveWidth     = (Diameter/2.0)/nGrooves ;
0066  LensPosition    = Pos ;
0067 
0068  if( GrooveWidth <= 0 ){
0069    G4Exception("UltraFresnelLens::UltraFresnelLens()","AirSh001",FatalException,
0070            "UltraFresnelLens constructor: GrooveWidth<=0");
0071    }
0072 
0073 
0074 auto  Rmin1 = (NumberOfGrooves-1)*(GrooveWidth) ;
0075 auto  Rmax1 = (NumberOfGrooves-0)*(GrooveWidth) ;
0076 LensThickness   = GetSagita(Rmax1)-GetSagita(Rmin1) ; // Height of the highest groove
0077 
0078 BuildLens(MotherPV) ;
0079 
0080 }
0081 
0082 
0083 UltraFresnelLens::~UltraFresnelLens(  )
0084 {;}
0085 
0086 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0087 
0088 void UltraFresnelLens::BuildLens(G4VPhysicalVolume *MotherPV){
0089 
0090 auto StartPhi = 0.0 ;
0091 auto DeltaPhi = twopi ;
0092 
0093 auto  LensMotherDz = LensThickness ;
0094 
0095 
0096 auto LensMotherCylinder
0097     = new G4Tubs("LensMotherCylinder",0.0*mm,LensDiameter/2.0,LensMotherDz/2.0,StartPhi,DeltaPhi);
0098 
0099 auto LensMotherMaterial = MotherPV->GetLogicalVolume()->GetMaterial() ;
0100 auto LensMotherLV
0101     = new G4LogicalVolume(LensMotherCylinder,LensMotherMaterial,"LensMotherLV",0,0,0);
0102 auto LensMotherPV
0103     = new G4PVPlacement(0,LensPosition,"LensMotherPV",LensMotherLV,MotherPV,false,0);
0104 
0105 LensMotherLV->SetVisAttributes (G4VisAttributes::GetInvisible());
0106 
0107 
0108 auto solidGroove
0109   = new G4Cons("Groove",40.0*mm,50.0*mm,40.0*mm,40.001*mm,1.0*mm,StartPhi,DeltaPhi);
0110 
0111 auto logicalGroove
0112     = new G4LogicalVolume(solidGroove,LensMaterial,"Groove_log",0,0,0);
0113 
0114 auto FresnelLensParam = new UltraFresnelLensParameterisation(this);
0115 
0116 LensPhysicalVolume =
0117  new G4PVParameterised("LensPV",logicalGroove,LensMotherPV,kZAxis,NumberOfGrooves,FresnelLensParam) ;
0118 
0119 }
0120 
0121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0122 
0123 G4double UltraFresnelLens::GetSagita(G4double radius)
0124 {
0125   auto Conic = -1.0;
0126   auto Curvature = 0.00437636761488/mm ;
0127   G4double Aspher[8] = {4.206739256e-05/(mm),
0128                         9.6440152e-10/(mm3),
0129                        -1.4884317e-15/(mm2*mm3),
0130                               0.0/(mm*mm3*mm3),
0131                               0.0/(mm3*mm3*mm3),
0132                               0.0/(mm2*mm3*mm3*mm3),
0133                               0.0/(mm*mm3*mm3*mm3*mm3),
0134                               0.0/(mm*3*mm3*mm3*mm3*mm3)
0135                             };
0136 
0137   auto TotAspher = 0.0*mm ;
0138 
0139   for(G4int k=1;k<9;k++){
0140     TotAspher += Aspher[k-1]*std::pow(radius,2*k) ;
0141   }
0142 
0143   auto ArgSqrt = 1.0-(1.0+Conic)*std::pow(Curvature,2)*std::pow(radius,2) ;
0144 
0145   if (ArgSqrt < 0.0){
0146     G4Exception("UltraFresnelLens::GetSagita()","AirSh002",
0147         FatalException,
0148         "UltraFresnelLensParameterisation::Sagita: Square Root of <0 !");
0149   }
0150   auto Sagita_value = Curvature*std::pow(radius,2)/(1.0+std::sqrt(ArgSqrt)) + TotAspher;
0151 
0152   return Sagita_value ;
0153 
0154 
0155 }