G4hBetheBlochModel.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // -------------------------------------------------------------------
00028 //
00029 // GEANT4 Class file
00030 //
00031 //
00032 // File name:     G4hBetheBlochModel
00033 //
00034 // Author:        V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
00035 // 
00036 // Creation date: 20 July 2000
00037 //
00038 // Modifications: 
00039 // 20/07/2000  V.Ivanchenko First implementation
00040 // 03/10/2000  V.Ivanchenko clean up accoding to CodeWizard
00041 //
00042 // Class Description: 
00043 //
00044 // Bethe-Bloch ionisation model
00045 //
00046 // Class Description: End 
00047 //
00048 // -------------------------------------------------------------------
00049 //
00050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00051 
00052 #include "G4hBetheBlochModel.hh"
00053 
00054 #include "globals.hh"
00055 #include "G4PhysicalConstants.hh"
00056 #include "G4SystemOfUnits.hh"
00057 #include "G4DynamicParticle.hh"
00058 #include "G4ParticleDefinition.hh"
00059 #include "G4Material.hh"
00060 
00061 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00062 
00063 G4hBetheBlochModel::G4hBetheBlochModel(const G4String& name)
00064   : G4VLowEnergyModel(name), 
00065     lowEnergyLimit(1.*MeV),
00066     highEnergyLimit(100.*GeV),
00067     twoln10(2.*std::log(10.)),
00068     bg2lim(0.0169), 
00069     taulim(8.4146e-3)
00070 {;}
00071 
00072 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00073 
00074 G4hBetheBlochModel::~G4hBetheBlochModel() 
00075 {;}
00076 
00077 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00078 
00079 G4double G4hBetheBlochModel::TheValue(const G4DynamicParticle* particle,
00080                                       const G4Material* material) 
00081 {
00082   G4double energy = particle->GetKineticEnergy() ;
00083   G4double particleMass = particle->GetMass() ;
00084 
00085   G4double eloss  = BetheBlochFormula(material,energy,particleMass) ;
00086 
00087   return eloss ;
00088 }
00089 
00090 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00091 
00092 G4double G4hBetheBlochModel::TheValue(const G4ParticleDefinition* aParticle,
00093                                       const G4Material* material,
00094                                       G4double kineticEnergy) 
00095 {
00096   G4double particleMass = aParticle->GetPDGMass() ;
00097   G4double eloss  = BetheBlochFormula(material,kineticEnergy,particleMass) ;
00098 
00099   return eloss ;
00100 }
00101 
00102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00103 
00104 G4double G4hBetheBlochModel::HighEnergyLimit(
00105                              const G4ParticleDefinition* ,
00106                              const G4Material* ) const
00107 {
00108   return highEnergyLimit ;
00109 }
00110 
00111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00112 
00113 G4double G4hBetheBlochModel::LowEnergyLimit(
00114                              const G4ParticleDefinition* aParticle,
00115                              const G4Material* material) const
00116 {
00117   G4double taul = (material->GetIonisation()->GetTaul())*
00118                   (aParticle->GetPDGMass()) ;
00119   return taul ;
00120 }
00121 
00122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00123 
00124 G4double G4hBetheBlochModel::HighEnergyLimit(
00125                              const G4ParticleDefinition* ) const
00126 {
00127   return highEnergyLimit ;
00128 }
00129 
00130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00131 
00132 G4double G4hBetheBlochModel::LowEnergyLimit(
00133                              const G4ParticleDefinition* ) const
00134 {
00135   return lowEnergyLimit ;
00136 }
00137 
00138 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00139  
00140 G4bool G4hBetheBlochModel::IsInCharge(const G4DynamicParticle* ,
00141                                       const G4Material* ) const
00142 {
00143   return true ;
00144 }
00145 
00146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00147  
00148 G4bool G4hBetheBlochModel::IsInCharge(const G4ParticleDefinition* ,
00149                                       const G4Material* ) const
00150 {
00151   return true ;
00152 }
00153 
00154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00155 
00156 G4double G4hBetheBlochModel::BetheBlochFormula(
00157                              const G4Material* material,
00158                                    G4double kineticEnergy,
00159                                    G4double particleMass) const
00160 {
00161   // This member function is applied normally to proton/antiproton
00162   G4double ionloss ;
00163 
00164   G4double rateMass = electron_mass_c2/particleMass ;
00165 
00166   G4double taul = material->GetIonisation()->GetTaul() ;
00167   G4double tau  = kineticEnergy/particleMass ;    // tau is relative energy
00168   
00169   // It is not normal case for this function
00170   // for low energy parametrisation have to be applied
00171   if ( tau < taul ) tau = taul ; 
00172     
00173   // some local variables 
00174     
00175   G4double gamma,bg2,beta2,tmax,x,delta,sh ;
00176   G4double electronDensity = material->GetElectronDensity();
00177   G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
00178   G4double eexc2 = eexc*eexc ;
00179   G4double cden  = material->GetIonisation()->GetCdensity();
00180   G4double mden  = material->GetIonisation()->GetMdensity();
00181   G4double aden  = material->GetIonisation()->GetAdensity();
00182   G4double x0den = material->GetIonisation()->GetX0density();
00183   G4double x1den = material->GetIonisation()->GetX1density();
00184   G4double* shellCorrectionVector =
00185             material->GetIonisation()->GetShellCorrectionVector();
00186     
00187   gamma = tau + 1.0 ;
00188   bg2 = tau*(tau+2.0) ;
00189   beta2 = bg2/(gamma*gamma) ;
00190   tmax = 2.*electron_mass_c2*bg2/(1.+2.*gamma*rateMass+rateMass*rateMass) ;
00191         
00192   ionloss = std::log(2.0*electron_mass_c2*bg2*tmax/eexc2)-2.0*beta2 ;
00193     
00194   // density correction     
00195   x = std::log(bg2)/twoln10 ;
00196   if ( x < x0den ) {
00197     delta = 0.0 ;
00198 
00199   } else { 
00200     delta = twoln10*x - cden ;
00201     if ( x < x1den ) delta += aden*std::pow((x1den-x),mden) ;
00202   }
00203     
00204   // shell correction 
00205   sh = 0.0 ;      
00206   x  = 1.0 ;
00207 
00208   if ( bg2 > bg2lim ) {
00209     for (G4int k=0; k<=2; k++) {
00210         x *= bg2 ;
00211         sh += shellCorrectionVector[k]/x;
00212     }
00213 
00214   } else {
00215     for (G4int k=0; k<=2; k++) {
00216         x *= bg2lim ;
00217         sh += shellCorrectionVector[k]/x;
00218     }
00219     sh *= std::log(tau/taul)/std::log(taulim/taul) ;     
00220   }
00221     
00222   // now compute the total ionization loss
00223     
00224   ionloss -= delta + sh ;
00225   ionloss *= twopi_mc2_rcl2*electronDensity/beta2 ;
00226 
00227   if ( ionloss < 0.0) ionloss = 0.0 ;
00228   
00229   return ionloss;
00230 }
00231 
00232 

Generated on Mon May 27 17:48:27 2013 for Geant4 by  doxygen 1.4.7