G4KleinNishinaCompton.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 // $Id$
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4KleinNishinaCompton
00034 //
00035 // Author:        Vladimir Ivanchenko on base of Michel Maire code
00036 //
00037 // Creation date: 15.03.2005
00038 //
00039 // Modifications:
00040 // 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
00041 // 27-03-06 Remove upper limit of cross section (V.Ivantchenko)
00042 //
00043 // Class Description:
00044 //
00045 // -------------------------------------------------------------------
00046 //
00047 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00048 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00049 
00050 #include "G4KleinNishinaCompton.hh"
00051 #include "G4PhysicalConstants.hh"
00052 #include "G4SystemOfUnits.hh"
00053 #include "G4Electron.hh"
00054 #include "G4Gamma.hh"
00055 #include "Randomize.hh"
00056 #include "G4DataVector.hh"
00057 #include "G4ParticleChangeForGamma.hh"
00058 
00059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00060 
00061 using namespace std;
00062 
00063 G4KleinNishinaCompton::G4KleinNishinaCompton(const G4ParticleDefinition*,
00064                                              const G4String& nam)
00065   : G4VEmModel(nam)
00066 {
00067   theGamma = G4Gamma::Gamma();
00068   theElectron = G4Electron::Electron();
00069   lowestGammaEnergy = 1.0*eV;
00070   fParticleChange = 0;
00071 }
00072 
00073 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00074 
00075 G4KleinNishinaCompton::~G4KleinNishinaCompton()
00076 {}
00077 
00078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00079 
00080 void G4KleinNishinaCompton::Initialise(const G4ParticleDefinition*,
00081                                        const G4DataVector&)
00082 {
00083   if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
00084 }
00085 
00086 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00087 
00088 G4double G4KleinNishinaCompton::ComputeCrossSectionPerAtom(
00089                                        const G4ParticleDefinition*,
00090                                              G4double GammaEnergy,
00091                                              G4double Z, G4double,
00092                                              G4double, G4double)
00093 {
00094   G4double xSection = 0.0 ;
00095   if ( Z < 0.9999 )                 return xSection;
00096   if ( GammaEnergy < 0.1*keV      ) return xSection;
00097   //  if ( GammaEnergy > (100.*GeV/Z) ) return xSection;
00098 
00099   static const G4double a = 20.0 , b = 230.0 , c = 440.0;
00100   
00101   static const G4double
00102     d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527   *barn, d4=-1.9798e+1*barn,
00103     e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
00104     f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
00105      
00106   G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
00107            p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
00108 
00109   G4double T0  = 15.0*keV; 
00110   if (Z < 1.5) T0 = 40.0*keV; 
00111 
00112   G4double X   = max(GammaEnergy, T0) / electron_mass_c2;
00113   xSection = p1Z*std::log(1.+2.*X)/X
00114                + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
00115                 
00116   //  modification for low energy. (special case for Hydrogen)
00117   if (GammaEnergy < T0) {
00118     G4double dT0 = 1.*keV;
00119     X = (T0+dT0) / electron_mass_c2 ;
00120     G4double sigma = p1Z*log(1.+2*X)/X
00121                     + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
00122     G4double   c1 = -T0*(sigma-xSection)/(xSection*dT0);             
00123     G4double   c2 = 0.150; 
00124     if (Z > 1.5) c2 = 0.375-0.0556*log(Z);
00125     G4double    y = log(GammaEnergy/T0);
00126     xSection *= exp(-y*(c1+c2*y));          
00127   }
00128   //  G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << xSection << G4endl;
00129   return xSection;
00130 }
00131 
00132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00133 
00134 void G4KleinNishinaCompton::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
00135                                               const G4MaterialCutsCouple*,
00136                                               const G4DynamicParticle* aDynamicGamma,
00137                                               G4double,
00138                                               G4double)
00139 {
00140   // The scattered gamma energy is sampled according to Klein - Nishina formula.
00141   // The random number techniques of Butcher & Messel are used 
00142   // (Nuc Phys 20(1960),15).
00143   // Note : Effects due to binding of atomic electrons are negliged.
00144  
00145   G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
00146 
00147   // extra protection
00148   if(gamEnergy0 < lowestGammaEnergy) {
00149     fParticleChange->ProposeTrackStatus(fStopAndKill);
00150     fParticleChange->ProposeLocalEnergyDeposit(gamEnergy0);
00151     fParticleChange->SetProposedKineticEnergy(0.0);
00152     return;
00153   }
00154 
00155   G4double E0_m = gamEnergy0 / electron_mass_c2 ;
00156 
00157   G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
00158 
00159   //
00160   // sample the energy rate of the scattered gamma 
00161   //
00162 
00163   G4double epsilon, epsilonsq, onecost, sint2, greject ;
00164 
00165   G4double eps0       = 1./(1. + 2.*E0_m);
00166   G4double epsilon0sq = eps0*eps0;
00167   G4double alpha1     = - log(eps0);
00168   G4double alpha2     = 0.5*(1.- epsilon0sq);
00169 
00170   do {
00171     if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
00172       epsilon   = exp(-alpha1*G4UniformRand());   // eps0**r
00173       epsilonsq = epsilon*epsilon; 
00174 
00175     } else {
00176       epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
00177       epsilon   = sqrt(epsilonsq);
00178     };
00179 
00180     onecost = (1.- epsilon)/(epsilon*E0_m);
00181     sint2   = onecost*(2.-onecost);
00182     greject = 1. - epsilon*sint2/(1.+ epsilonsq);
00183 
00184   } while (greject < G4UniformRand());
00185  
00186   //
00187   // scattered gamma angles. ( Z - axis along the parent gamma)
00188   //
00189 
00190   if(sint2 < 0.0) { sint2 = 0.0; }
00191   G4double cosTeta = 1. - onecost; 
00192   G4double sinTeta = sqrt (sint2);
00193   G4double Phi     = twopi * G4UniformRand();
00194 
00195   //
00196   // update G4VParticleChange for the scattered gamma
00197   //
00198    
00199   G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
00200   gamDirection1.rotateUz(gamDirection0);
00201   G4double gamEnergy1 = epsilon*gamEnergy0;
00202   if(gamEnergy1 > lowestGammaEnergy) {
00203     fParticleChange->ProposeMomentumDirection(gamDirection1);
00204     fParticleChange->SetProposedKineticEnergy(gamEnergy1);
00205   } else { 
00206     fParticleChange->ProposeTrackStatus(fStopAndKill);
00207     fParticleChange->ProposeLocalEnergyDeposit(gamEnergy1);
00208     fParticleChange->SetProposedKineticEnergy(0.0);
00209   }
00210 
00211   //
00212   // kinematic of the scattered electron
00213   //
00214 
00215   G4double eKinEnergy = gamEnergy0 - gamEnergy1;
00216 
00217   if(eKinEnergy > DBL_MIN) {
00218     G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
00219     eDirection = eDirection.unit();
00220 
00221     // create G4DynamicParticle object for the electron.
00222     G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
00223     fvect->push_back(dp);
00224   }
00225 }
00226 
00227 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00228 
00229 

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