G4BetheHeitlerModel.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:     G4BetheHeitlerModel
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 // 24-06-05 Increase number of bins to 200 (V.Ivantchenko)
00042 // 16-11-05 replace shootBit() by G4UniformRand()  mma
00043 // 04-12-05 SetProposedKineticEnergy(0.) for the killed photon (mma)
00044 // 20-02-07 SelectRandomElement is called for any initial gamma energy 
00045 //          in order to have selected element for polarized model (VI)
00046 // 25-10-10 Removed unused table, added element selector (VI) 
00047 //
00048 // Class Description:
00049 //
00050 // -------------------------------------------------------------------
00051 //
00052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00054 
00055 #include "G4BetheHeitlerModel.hh"
00056 #include "G4PhysicalConstants.hh"
00057 #include "G4SystemOfUnits.hh"
00058 #include "G4Electron.hh"
00059 #include "G4Positron.hh"
00060 #include "G4Gamma.hh"
00061 #include "Randomize.hh"
00062 #include "G4ParticleChangeForGamma.hh"
00063 
00064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00065 
00066 using namespace std;
00067 
00068 G4BetheHeitlerModel::G4BetheHeitlerModel(const G4ParticleDefinition*,
00069                                          const G4String& nam)
00070   : G4VEmModel(nam)
00071 {
00072   fParticleChange = 0;
00073   theGamma    = G4Gamma::Gamma();
00074   thePositron = G4Positron::Positron();
00075   theElectron = G4Electron::Electron();
00076 }
00077 
00078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00079 
00080 G4BetheHeitlerModel::~G4BetheHeitlerModel()
00081 {}
00082 
00083 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00084 
00085 void G4BetheHeitlerModel::Initialise(const G4ParticleDefinition* p,
00086                                      const G4DataVector& cuts)
00087 {
00088   if(!fParticleChange) { fParticleChange = GetParticleChangeForGamma(); }
00089   InitialiseElementSelectors(p, cuts);
00090 }
00091 
00092 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00093 
00094 G4double 
00095 G4BetheHeitlerModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
00096                                                 G4double GammaEnergy, G4double Z,
00097                                                 G4double, G4double, G4double)
00098 // Calculates the microscopic cross section in GEANT4 internal units.
00099 // A parametrized formula from L. Urban is used to estimate
00100 // the total cross section.
00101 // It gives a good description of the data from 1.5 MeV to 100 GeV.
00102 // below 1.5 MeV: sigma=sigma(1.5MeV)*(GammaEnergy-2electronmass)
00103 //                                   *(GammaEnergy-2electronmass) 
00104 {
00105   static const G4double GammaEnergyLimit = 1.5*MeV;
00106   G4double xSection = 0.0 ;
00107   if ( Z < 0.9 || GammaEnergy <= 2.0*electron_mass_c2 ) { return xSection; }
00108 
00109   static const G4double
00110     a0= 8.7842e+2*microbarn, a1=-1.9625e+3*microbarn, a2= 1.2949e+3*microbarn,
00111     a3=-2.0028e+2*microbarn, a4= 1.2575e+1*microbarn, a5=-2.8333e-1*microbarn;
00112 
00113   static const G4double
00114     b0=-1.0342e+1*microbarn, b1= 1.7692e+1*microbarn, b2=-8.2381   *microbarn,
00115     b3= 1.3063   *microbarn, b4=-9.0815e-2*microbarn, b5= 2.3586e-3*microbarn;
00116 
00117   static const G4double
00118     c0=-4.5263e+2*microbarn, c1= 1.1161e+3*microbarn, c2=-8.6749e+2*microbarn,
00119     c3= 2.1773e+2*microbarn, c4=-2.0467e+1*microbarn, c5= 6.5372e-1*microbarn;
00120 
00121   G4double GammaEnergySave = GammaEnergy;
00122   if (GammaEnergy < GammaEnergyLimit) { GammaEnergy = GammaEnergyLimit; }
00123 
00124   G4double X=log(GammaEnergy/electron_mass_c2), X2=X*X, X3=X2*X, X4=X3*X, X5=X4*X;
00125 
00126   G4double F1 = a0 + a1*X + a2*X2 + a3*X3 + a4*X4 + a5*X5,
00127            F2 = b0 + b1*X + b2*X2 + b3*X3 + b4*X4 + b5*X5,
00128            F3 = c0 + c1*X + c2*X2 + c3*X3 + c4*X4 + c5*X5;     
00129 
00130   xSection = (Z + 1.)*(F1*Z + F2*Z*Z + F3);
00131 
00132   if (GammaEnergySave < GammaEnergyLimit) {
00133 
00134     X = (GammaEnergySave  - 2.*electron_mass_c2)
00135       / (GammaEnergyLimit - 2.*electron_mass_c2);
00136     xSection *= X*X;
00137   }
00138 
00139   if (xSection < 0.) { xSection = 0.; }
00140   return xSection;
00141 }
00142 
00143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00144 
00145 void G4BetheHeitlerModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
00146                                             const G4MaterialCutsCouple* couple,
00147                                             const G4DynamicParticle* aDynamicGamma,
00148                                             G4double,
00149                                             G4double)
00150 // The secondaries e+e- energies are sampled using the Bethe - Heitler
00151 // cross sections with Coulomb correction.
00152 // A modified version of the random number techniques of Butcher & Messel
00153 // is used (Nuc Phys 20(1960),15).
00154 //
00155 // GEANT4 internal units.
00156 //
00157 // Note 1 : Effects due to the breakdown of the Born approximation at
00158 //          low energy are ignored.
00159 // Note 2 : The differential cross section implicitly takes account of 
00160 //          pair creation in both nuclear and atomic electron fields.
00161 //          However triplet prodution is not generated.
00162 {
00163   const G4Material* aMaterial = couple->GetMaterial();
00164 
00165   G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
00166   G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
00167 
00168   G4double epsil ;
00169   G4double epsil0 = electron_mass_c2/GammaEnergy ;
00170   if(epsil0 > 1.0) { return; }
00171 
00172   // do it fast if GammaEnergy < 2. MeV
00173   static const G4double Egsmall=2.*MeV;
00174 
00175   // select randomly one element constituing the material
00176   const G4Element* anElement = SelectRandomAtom(aMaterial, theGamma, GammaEnergy);
00177 
00178   if (GammaEnergy < Egsmall) {
00179 
00180     epsil = epsil0 + (0.5-epsil0)*G4UniformRand();
00181 
00182   } else {
00183     // now comes the case with GammaEnergy >= 2. MeV
00184 
00185     // Extract Coulomb factor for this Element
00186     G4double FZ = 8.*(anElement->GetIonisation()->GetlogZ3());
00187     if (GammaEnergy > 50.*MeV) { FZ += 8.*(anElement->GetfCoulomb()); }
00188 
00189     // limits of the screening variable
00190     G4double screenfac = 136.*epsil0/(anElement->GetIonisation()->GetZ3());
00191     G4double screenmax = exp ((42.24 - FZ)/8.368) - 0.952 ;
00192     G4double screenmin = min(4.*screenfac,screenmax);
00193 
00194     // limits of the energy sampling
00195     G4double epsil1 = 0.5 - 0.5*sqrt(1. - screenmin/screenmax) ;
00196     G4double epsilmin = max(epsil0,epsil1) , epsilrange = 0.5 - epsilmin;
00197 
00198     //
00199     // sample the energy rate of the created electron (or positron)
00200     //
00201     //G4double epsil, screenvar, greject ;
00202     G4double  screenvar, greject ;
00203 
00204     G4double F10 = ScreenFunction1(screenmin) - FZ;
00205     G4double F20 = ScreenFunction2(screenmin) - FZ;
00206     G4double NormF1 = max(F10*epsilrange*epsilrange,0.); 
00207     G4double NormF2 = max(1.5*F20,0.);
00208 
00209     do {
00210       if ( NormF1/(NormF1+NormF2) > G4UniformRand() ) {
00211         epsil = 0.5 - epsilrange*pow(G4UniformRand(), 0.333333);
00212         screenvar = screenfac/(epsil*(1-epsil));
00213         greject = (ScreenFunction1(screenvar) - FZ)/F10;
00214               
00215       } else { 
00216         epsil = epsilmin + epsilrange*G4UniformRand();
00217         screenvar = screenfac/(epsil*(1-epsil));
00218         greject = (ScreenFunction2(screenvar) - FZ)/F20;
00219       }
00220 
00221     } while( greject < G4UniformRand() );
00222 
00223   }   //  end of epsil sampling
00224    
00225   //
00226   // fixe charges randomly
00227   //
00228 
00229   G4double ElectTotEnergy, PositTotEnergy;
00230   if (G4UniformRand() > 0.5) {
00231 
00232     ElectTotEnergy = (1.-epsil)*GammaEnergy;
00233     PositTotEnergy = epsil*GammaEnergy;
00234      
00235   } else {
00236     
00237     PositTotEnergy = (1.-epsil)*GammaEnergy;
00238     ElectTotEnergy = epsil*GammaEnergy;
00239   }
00240 
00241   //
00242   // scattered electron (positron) angles. ( Z - axis along the parent photon)
00243   //
00244   //  universal distribution suggested by L. Urban 
00245   // (Geant3 manual (1993) Phys211),
00246   //  derived from Tsai distribution (Rev Mod Phys 49,421(1977))
00247 
00248   G4double u;
00249   const G4double a1 = 0.625 , a2 = 3.*a1 , d = 27. ;
00250 
00251   if (9./(9.+d) >G4UniformRand()) u= - log(G4UniformRand()*G4UniformRand())/a1;
00252   else                            u= - log(G4UniformRand()*G4UniformRand())/a2;
00253 
00254   G4double TetEl = u*electron_mass_c2/ElectTotEnergy;
00255   G4double TetPo = u*electron_mass_c2/PositTotEnergy;
00256   G4double Phi  = twopi * G4UniformRand();
00257   G4double dxEl= sin(TetEl)*cos(Phi),dyEl= sin(TetEl)*sin(Phi),dzEl=cos(TetEl);
00258   G4double dxPo=-sin(TetPo)*cos(Phi),dyPo=-sin(TetPo)*sin(Phi),dzPo=cos(TetPo);
00259    
00260   //
00261   // kinematic of the created pair
00262   //
00263   // the electron and positron are assumed to have a symetric
00264   // angular distribution with respect to the Z axis along the parent photon.
00265 
00266   G4double ElectKineEnergy = max(0.,ElectTotEnergy - electron_mass_c2);
00267 
00268   G4ThreeVector ElectDirection (dxEl, dyEl, dzEl);
00269   ElectDirection.rotateUz(GammaDirection);   
00270 
00271   // create G4DynamicParticle object for the particle1  
00272   G4DynamicParticle* aParticle1= new G4DynamicParticle(
00273                      theElectron,ElectDirection,ElectKineEnergy);
00274   
00275   // the e+ is always created (even with Ekine=0) for further annihilation.
00276 
00277   G4double PositKineEnergy = max(0.,PositTotEnergy - electron_mass_c2);
00278 
00279   G4ThreeVector PositDirection (dxPo, dyPo, dzPo);
00280   PositDirection.rotateUz(GammaDirection);   
00281 
00282   // create G4DynamicParticle object for the particle2 
00283   G4DynamicParticle* aParticle2= new G4DynamicParticle(
00284                       thePositron,PositDirection,PositKineEnergy);
00285 
00286   // Fill output vector
00287   fvect->push_back(aParticle1);
00288   fvect->push_back(aParticle2);
00289 
00290   // kill incident photon
00291   fParticleChange->SetProposedKineticEnergy(0.);
00292   fParticleChange->ProposeTrackStatus(fStopAndKill);   
00293 }
00294 
00295 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

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