G4ICRU73QOModel.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:   G4ICRU73QOModel
00034 //
00035 // Author:        Alexander Bagulya
00036 //
00037 // Creation date: 21.05.2010
00038 //
00039 // Modifications: 
00040 //
00041 //
00042 // -------------------------------------------------------------------
00043 //
00044 
00045 
00046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00047 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00048 
00049 #include "G4ICRU73QOModel.hh"
00050 #include "G4PhysicalConstants.hh"
00051 #include "G4SystemOfUnits.hh"
00052 #include "Randomize.hh"
00053 #include "G4Electron.hh"
00054 #include "G4ParticleChangeForLoss.hh"
00055 #include "G4LossTableManager.hh"
00056 #include "G4AntiProton.hh"
00057 
00058 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00059 
00060 using namespace std;
00061 
00062 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00063 
00064 G4ICRU73QOModel::G4ICRU73QOModel(const G4ParticleDefinition* p, const G4String& nam)
00065   : G4VEmModel(nam),
00066     particle(0),
00067     isInitialised(false)
00068 {
00069   mass = charge = chargeSquare = massRate = ratio = 0.0;
00070   if(p) { SetParticle(p); }
00071   SetHighEnergyLimit(10.0*MeV);
00072 
00073   lowestKinEnergy  = 5.0*keV;
00074 
00075   sizeL0 = 67;
00076   sizeL1 = 22;
00077   sizeL2 = 14;
00078 
00079   theElectron = G4Electron::Electron();
00080 
00081   for (G4int i = 0; i < 100; ++i)
00082     {
00083       indexZ[i] = -1;
00084     }
00085   for(G4int i = 0; i < NQOELEM; ++i) 
00086     {
00087       if(ZElementAvailable[i] > 0) {
00088         indexZ[ZElementAvailable[i]] = i;
00089       }
00090     }
00091   fParticleChange = 0;
00092   denEffData = 0;
00093 }
00094 
00095 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00096 
00097 G4ICRU73QOModel::~G4ICRU73QOModel()
00098 {}
00099 
00100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00101 
00102 void G4ICRU73QOModel::Initialise(const G4ParticleDefinition* p,
00103                                  const G4DataVector&)
00104 {
00105   if(p != particle) SetParticle(p);
00106 
00107   // always false before the run
00108   SetDeexcitationFlag(false);
00109 
00110   if(!isInitialised) {
00111     isInitialised = true;
00112 
00113     G4String pname = particle->GetParticleName();
00114     fParticleChange = GetParticleChangeForLoss();
00115     const G4MaterialTable* mtab = G4Material::GetMaterialTable(); 
00116     denEffData = (*mtab)[0]->GetIonisation()->GetDensityEffectData();
00117   }
00118 }
00119 
00120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00121 
00122 G4double G4ICRU73QOModel::ComputeCrossSectionPerElectron(
00123                                            const G4ParticleDefinition* p,
00124                                                  G4double kineticEnergy,
00125                                                  G4double cutEnergy,
00126                                                  G4double maxKinEnergy)
00127 {
00128   G4double cross     = 0.0;
00129   G4double tmax      = MaxSecondaryEnergy(p, kineticEnergy);
00130   G4double maxEnergy = std::min(tmax,maxKinEnergy);
00131   if(cutEnergy < maxEnergy) {
00132 
00133     G4double energy  = kineticEnergy + mass;
00134     G4double energy2 = energy*energy;
00135     G4double beta2   = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
00136     cross = 1.0/cutEnergy - 1.0/maxEnergy - beta2*log(maxEnergy/cutEnergy)/tmax;
00137 
00138     cross *= CLHEP::twopi_mc2_rcl2*chargeSquare/beta2;
00139   }
00140  
00141   return cross;
00142 }
00143 
00144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00145 
00146 G4double G4ICRU73QOModel::ComputeCrossSectionPerAtom(
00147                                            const G4ParticleDefinition* p,
00148                                                  G4double kineticEnergy,
00149                                                  G4double Z, G4double,
00150                                                  G4double cutEnergy,
00151                                                  G4double maxEnergy)
00152 {
00153   G4double cross = Z*ComputeCrossSectionPerElectron
00154                                          (p,kineticEnergy,cutEnergy,maxEnergy);
00155   return cross;
00156 }
00157 
00158 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00159 
00160 G4double G4ICRU73QOModel::CrossSectionPerVolume(
00161                                            const G4Material* material,
00162                                            const G4ParticleDefinition* p,
00163                                                  G4double kineticEnergy,
00164                                                  G4double cutEnergy,
00165                                                  G4double maxEnergy)
00166 {
00167   G4double eDensity = material->GetElectronDensity();
00168   G4double cross = eDensity*ComputeCrossSectionPerElectron
00169                                          (p,kineticEnergy,cutEnergy,maxEnergy);
00170   return cross;
00171 }
00172 
00173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00174 
00175 G4double G4ICRU73QOModel::ComputeDEDXPerVolume(const G4Material* material,
00176                                                const G4ParticleDefinition* p,
00177                                                G4double kineticEnergy,
00178                                                G4double cutEnergy)
00179 {
00180   SetParticle(p);
00181   G4double tmax  = MaxSecondaryEnergy(p, kineticEnergy);
00182   G4double tkin  = kineticEnergy/massRate;
00183   G4double dedx  = 0.0;
00184   if(tkin > lowestKinEnergy) { dedx = DEDX(material, tkin); }
00185   else { dedx = DEDX(material, lowestKinEnergy)*sqrt(tkin/lowestKinEnergy); }
00186 
00187   if (cutEnergy < tmax) {
00188 
00189     G4double tau   = kineticEnergy/mass;
00190     G4double gam   = tau + 1.0;
00191     G4double bg2   = tau * (tau+2.0);
00192     G4double beta2 = bg2/(gam*gam);
00193     G4double x     = cutEnergy/tmax;
00194 
00195     dedx += chargeSquare*( log(x) + (1.0 - x)*beta2 ) * twopi_mc2_rcl2
00196           * material->GetElectronDensity()/beta2;
00197   }
00198   if(dedx < 0.0) { dedx = 0.0; }
00199   return dedx;
00200 }
00201 
00202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00203 
00204 G4double G4ICRU73QOModel::DEDX(const G4Material* material,
00205                                G4double kineticEnergy) 
00206 {
00207   G4double eloss = 0.0;
00208   const G4int numberOfElements = material->GetNumberOfElements();
00209   const G4double* theAtomicNumDensityVector =
00210                                  material->GetAtomicNumDensityVector();
00211   
00212   // Bragg's rule calculation
00213   const G4ElementVector* theElementVector =
00214                            material->GetElementVector() ;
00215   
00216   //  loop for the elements in the material
00217   for (G4int i=0; i<numberOfElements; ++i)
00218     {
00219       const G4Element* element = (*theElementVector)[i] ;
00220       eloss   += DEDXPerElement(G4int(element->GetZ()), kineticEnergy)
00221                                  * theAtomicNumDensityVector[i] * G4int(element->GetZ());
00222     }      
00223   return eloss;
00224 }
00225 
00226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00227 
00228 G4double G4ICRU73QOModel::DEDXPerElement(G4int AtomicNumber,
00229                                          G4double kineticEnergy)
00230 {
00231   G4int Z = AtomicNumber;
00232   if(Z > 97) { Z = 97; }
00233   G4int nbOfShells = GetNumberOfShells(Z);
00234   if(nbOfShells < 1) { nbOfShells = 1; }
00235 
00236   G4double v = CLHEP::c_light * std::sqrt( 2.0*kineticEnergy/proton_mass_c2 );
00237 
00238   G4double fBetheVelocity = CLHEP::fine_structure_const*CLHEP::c_light/v;
00239 
00240   G4double tau   = kineticEnergy/proton_mass_c2;
00241   G4double gam   = tau + 1.0;
00242   G4double bg2   = tau * (tau+2.0);
00243   G4double beta2 = bg2/(gam*gam);
00244 
00245   G4double l0Term = 0, l1Term = 0, l2Term = 0;
00246   
00247   for (G4int nos = 0; nos < nbOfShells; ++nos){
00248     
00249     G4double NormalizedEnergy = (2.0*CLHEP::electron_mass_c2*beta2) / 
00250       GetShellEnergy(Z,nos);
00251 
00252     G4double shStrength = GetShellStrength(Z,nos);
00253 
00254     G4double l0 = GetL0(NormalizedEnergy);
00255     l0Term += shStrength  * l0; 
00256 
00257     G4double l1 = GetL1(NormalizedEnergy);
00258     l1Term += shStrength * l1; 
00259 
00260     G4double l2 = GetL2(NormalizedEnergy);
00261     l2Term += shStrength * l2;
00262 
00263   }
00264   G4double dedx  = 2*CLHEP::twopi_mc2_rcl2*chargeSquare*factorBethe[Z]*
00265     (l0Term + charge*fBetheVelocity*l1Term 
00266      + chargeSquare*fBetheVelocity*fBetheVelocity*l2Term)/beta2;
00267   return dedx;
00268 }
00269 
00270 
00271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00272 
00273 G4double G4ICRU73QOModel::GetOscillatorEnergy(G4int Z,
00274                                               G4int nbOfTheShell) const
00275 { 
00276   G4int idx = denEffData->GetElementIndex(Z, kStateUndefined);
00277   if(idx == -1) { idx = denEffData->GetElementIndex(Z-1, kStateUndefined); }
00278   G4double PlasmaEnergy = denEffData->GetPlasmaEnergy(idx);
00279  
00280   G4double PlasmaEnergy2 = PlasmaEnergy * PlasmaEnergy;
00281 
00282   G4double plasmonTerm = 0.66667 * G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell)  
00283     * PlasmaEnergy2 / (Z*Z) ; 
00284 
00285   G4double ionTerm = std::exp(0.5) * (G4AtomicShells::GetBindingEnergy(Z,nbOfTheShell)) ;
00286   G4double ionTerm2 = ionTerm*ionTerm ;
00287    
00288   G4double oscShellEnergy = std::sqrt( ionTerm2 + plasmonTerm );
00289 
00290   return  oscShellEnergy;
00291 }
00292 
00293 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00294 
00295 G4double G4ICRU73QOModel::GetL0(G4double normEnergy) const 
00296 {
00297   G4int n;
00298   
00299   for(n = 0; n < sizeL0; n++) {
00300     if( normEnergy < L0[n][0] ) break;
00301   }
00302   if(0 == n) { n = 1; }
00303   if(n >= sizeL0) { n = sizeL0 - 1; }
00304 
00305   G4double l0    = L0[n][1];
00306   G4double l0p   = L0[n-1][1];
00307   G4double bethe = l0p + (l0 - l0p) * ( normEnergy - L0[n-1][0]) / 
00308                   (L0[n][0] - L0[n-1][0]);
00309 
00310   return bethe ;
00311 }
00312 
00313 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00314 
00315 G4double G4ICRU73QOModel::GetL1(G4double normEnergy) const
00316 {
00317   G4int n;
00318 
00319   for(n = 0; n < sizeL1; n++) {
00320     if( normEnergy < L1[n][0] ) break;
00321   }
00322   if(0 == n) n = 1 ;
00323   if(n >= sizeL1) n = sizeL1 - 1 ;
00324 
00325   G4double l1    = L1[n][1];
00326   G4double l1p   = L1[n-1][1];
00327   G4double barkas= l1p + (l1 - l1p) * ( normEnergy - L1[n-1][0]) / 
00328                   (L1[n][0] - L1[n-1][0]);
00329 
00330   return barkas;
00331 }
00332 
00333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00334 
00335 G4double G4ICRU73QOModel::GetL2(G4double normEnergy) const
00336 {
00337   G4int n;
00338   for(n = 0; n < sizeL2; n++) {
00339     if( normEnergy < L2[n][0] ) break;
00340   }
00341   if(0 == n) n = 1 ;
00342   if(n >= sizeL2) n = sizeL2 - 1 ;
00343 
00344   G4double l2    = L2[n][1];
00345   G4double l2p   = L2[n-1][1];
00346   G4double bloch = l2p + (l2 - l2p) * ( normEnergy - L2[n-1][0]) / 
00347                   (L2[n][0] - L2[n-1][0]);
00348 
00349   return bloch;
00350 }
00351 
00352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00353 
00354 void G4ICRU73QOModel::CorrectionsAlongStep(const G4MaterialCutsCouple*,
00355                                            const G4DynamicParticle*,
00356                                            G4double&,
00357                                            G4double&,
00358                                            G4double)
00359 {}
00360 
00361 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00362 
00363 void G4ICRU73QOModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
00364                                         const G4MaterialCutsCouple*,
00365                                         const G4DynamicParticle* dp,
00366                                         G4double xmin,
00367                                         G4double maxEnergy)
00368 {
00369   G4double tmax = MaxSecondaryKinEnergy(dp);
00370   G4double xmax = std::min(tmax, maxEnergy);
00371   if(xmin >= xmax) { return; }
00372 
00373   G4double kineticEnergy = dp->GetKineticEnergy();
00374   G4double energy  = kineticEnergy + mass;
00375   G4double energy2 = energy*energy;
00376   G4double beta2   = kineticEnergy*(kineticEnergy + 2.0*mass)/energy2;
00377   G4double grej    = 1.0;
00378   G4double deltaKinEnergy, f;
00379 
00380   G4ThreeVector direction = dp->GetMomentumDirection();
00381 
00382   // sampling follows ...
00383   do {
00384     G4double x = G4UniformRand();
00385     deltaKinEnergy = xmin*xmax/(xmin*(1.0 - x) + xmax*x);
00386 
00387     f = 1.0 - beta2*deltaKinEnergy/tmax;
00388 
00389     if(f > grej) {
00390         G4cout << "G4ICRU73QOModel::SampleSecondary Warning! "
00391                << "Majorant " << grej << " < "
00392                << f << " for e= " << deltaKinEnergy
00393                << G4endl;
00394     }
00395 
00396   } while( grej*G4UniformRand() >= f );
00397 
00398   G4double deltaMomentum =
00399            sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
00400   G4double totMomentum = energy*sqrt(beta2);
00401   G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
00402                                    (deltaMomentum * totMomentum);
00403   if(cost > 1.0) { cost = 1.0; }
00404   G4double sint = sqrt((1.0 - cost)*(1.0 + cost));
00405 
00406   G4double phi = twopi * G4UniformRand() ;
00407 
00408   G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
00409   deltaDirection.rotateUz(direction);
00410 
00411   // Change kinematics of primary particle
00412   kineticEnergy       -= deltaKinEnergy;
00413   G4ThreeVector finalP = direction*totMomentum - deltaDirection*deltaMomentum;
00414   finalP               = finalP.unit();
00415   
00416   fParticleChange->SetProposedKineticEnergy(kineticEnergy);
00417   fParticleChange->SetProposedMomentumDirection(finalP);
00418 
00419   // create G4DynamicParticle object for delta ray
00420   G4DynamicParticle* delta = new G4DynamicParticle(theElectron,deltaDirection,
00421                                                    deltaKinEnergy);
00422 
00423   vdp->push_back(delta);
00424 }
00425 
00426 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00427 
00428 G4double G4ICRU73QOModel::MaxSecondaryEnergy(const G4ParticleDefinition* pd,
00429                                              G4double kinEnergy)
00430 {
00431   if(pd != particle) { SetParticle(pd); }
00432   G4double tau  = kinEnergy/mass;
00433   G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
00434                   (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
00435   return tmax;
00436 }
00437 
00438 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00439 
00440 const G4int G4ICRU73QOModel::ZElementAvailable[NQOELEM] = {1,2,4,6,7,8,10,13,14,-18,
00441                                                   22,26,28,29,32,36,42,47,
00442                                                   50,54,73,74,78,79,82,92};
00443 
00444 const G4int G4ICRU73QOModel::nbofShellsForElement[NQOELEM] = {1,1,2,3,3,3,3,4,5,4,
00445                                                          5,5,5,5,6,4,6,6,
00446                                                          7,6,6,8,7,7,9,9};
00447 
00448 const G4int G4ICRU73QOModel::startElemIndex[NQOELEM] = {0,1,2,4,7,10,13,16,20,25,
00449                                                   29,34,39,44,49,55,59,65,
00450                                                   71,78,84,90,98,105,112,121};
00451 
00452 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00453 
00454 // SubShellOccupation = Z * ShellStrength
00455 const G4double G4ICRU73QOModel::SubShellOccupation[NQODATA] = 
00456   {
00457     1.000, // H 0
00458     2.000, // He 1
00459     1.930, 2.070, // Be 2-3
00460     1.992, 1.841, 2.167, // C 4-6    
00461     1.741, 1.680, 3.579, // N 7-9
00462     1.802, 1.849, 4.349, // O 10-12
00463     1.788, 2.028, 6.184, // Ne 13-15
00464     1.623, 2.147, 6.259, 2.971, // Al 16-19
00465     1.631, 2.094, 6.588, 2.041, 1.646, // Si 20-24
00466     1.535, 8.655, 1.706, 6.104, // Ar 25-28
00467     1.581, 8.358, 8.183, 2.000, 1.878, // Ti 29-33
00468     1.516, 8.325, 8.461, 6.579, 1.119, // Fe 34-38
00469     1.422, 7.81, 8.385, 8.216, 2.167, // Ni 39-43
00470     1.458, 8.049, 8.79, 9.695, 1.008, // Cu 44-48
00471     1.442, 7.791, 7.837, 10.122, 2.463, 2.345, // Ge 49-54
00472     1.645, 7.765, 19.192, 7.398, // Kr 55-58
00473     1.313, 6.409, 19.229, 8.633, 5.036, 1.380, // Mo 59-64
00474     1.295, 6.219, 18.751, 8.748, 10.184, 1.803, // Ag 65-70
00475     1.277, 6.099, 20.386, 8.011, 10.007, 2.272, 1.948, // Sn 71-77
00476     1.563, 6.312, 21.868, 5.762, 11.245, 7.250, // Xe 78-83
00477     0.9198, 6.5408, 18.9727, 24.9149, 15.0161, 6.6284, // Ta 84-89
00478     1.202, 5.582, 19.527, 18.741, 8.411, 14.387, 4.042, 2.108, // W 90-97
00479     1.159, 5.467, 18.802, 33.905, 8.300, 9.342, 1.025, // Pt 98-104
00480     1.124, 5.331, 18.078, 34.604, 8.127, 10.414, 1.322, // Au 105-111
00481     2.000, 8.000, 18.000, 18.000, 14.000, 8.000, 10.000, 2.000, 2.000, // Pb 112-120
00482     2.000, 8.000, 18.000, 32.000, 18.000, 8.000, 2.000, 1.000, 3.000  // U 121-129
00483 };
00484 
00485 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00486 
00487 // ShellEnergy in eV
00488 const G4double G4ICRU73QOModel::ShellEnergy[NQODATA] = 
00489   {
00490     19.2, // H
00491     41.8, // He
00492     209.11, 21.68, // Be
00493     486.2, 60.95, 23.43, // C
00494     732.61, 100.646, 23.550, // N    
00495     965.1, 129.85, 31.60, // O
00496     1525.9, 234.9, 56.18, // Ne
00497     2701, 476.5, 150.42, 16.89, // Al
00498     3206.1, 586.4, 186.8, 23.52, 14.91, // Si
00499     5551.6, 472.43, 124.85, 22.332, // Ar
00500     8554.6, 850.58, 93.47, 39.19, 19.46, // Ti
00501     12254.7, 1279.29, 200.35, 49.19, 17.66, // Fe
00502     14346.9, 1532.28, 262.71, 74.37, 23.03, // Ni
00503     15438.5, 1667.96, 294.1, 70.69, 16.447, // Cu
00504     19022.1, 2150.79, 455.79, 179.87, 57.89, 20.95, // Ge
00505     24643, 2906.4, 366.85, 22.24, // Kr
00506     34394, 4365.3, 589.36, 129.42, 35.59, 18.42, // Mo
00507     43664.3, 5824.91, 909.79, 175.47, 54.89, 19.63, // Ag
00508     49948, 6818.2, 1036.1, 172.65, 70.89, 33.87, 14.54, // Sn
00509     58987, 8159, 1296.6, 356.75, 101.03, 16.52, // Xe
00510     88926, 18012, 3210, 575, 108.7, 30.8, // Ta
00511     115025.9, 17827.44, 3214.36, 750.41, 305.21, 105.50, 38.09, 21.25, // W
00512     128342, 20254, 3601.8, 608.1, 115.0, 42.75, 17.04, // Pt
00513     131872, 20903, 3757.4, 682.1, 105.2, 44.89, 17.575, // Au
00514     154449, 25067, 5105.0, 987.44, 247.59, 188.1, 40.61, 19.2, 15.17, // Pb
00515     167282, 27868, 6022.7, 1020.4, 244.81, 51.33, 13, 11.06, 14.43  // U
00516 };
00517 
00518 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00519 
00520 // Data for L0 from: Sigmund P., Haagerup U. Phys. Rev. A34 (1986) 892-910
00521 const G4double G4ICRU73QOModel::L0[67][2] =
00522 {
00523   {0.00,        0.000001},
00524   {0.10,        0.000001},
00525   {0.12,        0.00001},
00526   {0.14,        0.00005},
00527   {0.16,        0.00014},
00528   {0.18,        0.00030},
00529   {0.20,        0.00057},
00530   {0.25,        0.00189},
00531   {0.30,        0.00429},
00532   {0.35,        0.00784},
00533   {0.40,        0.01248},
00534   {0.45,        0.01811},
00535   {0.50,        0.02462},
00536   {0.60,        0.03980},
00537   {0.70,        0.05731},
00538   {0.80,        0.07662},
00539   {0.90,        0.09733},
00540   {1.00,        0.11916},
00541   {1.20,        0.16532},
00542   {1.40,        0.21376},
00543   {1.60,        0.26362},
00544   {1.80,        0.31428},
00545   {2.00,        0.36532},
00546   {2.50,        0.49272},
00547   {3.00,        0.61765},
00548   {3.50,        0.73863},
00549   {4.00,        0.85496},
00550   {4.50,        0.96634},
00551   {5.00,        1.07272},
00552   {6.00,        1.27086},
00553   {7.00,        1.45075},
00554   {8.00,        1.61412},
00555   {9.00,        1.76277},
00556   {10.00,       1.89836},
00557   {12.00,       2.13625},
00558   {14.00,       2.33787},
00559   {16.00,       2.51093},
00560   {18.00,       2.66134},
00561   {20.00,       2.79358},
00562   {25.00,       3.06539},
00563   {30.00,       3.27902},
00564   {35.00,       3.45430},
00565   {40.00,       3.60281},
00566   {45.00,       3.73167},
00567   {50.00,       3.84555},
00568   {60.00,       4.04011},
00569   {70.00,       4.20264},
00570   {80.00,       4.34229},
00571   {90.00,       4.46474},
00572   {100.00,      4.57378},
00573   {120.00,      4.76155},
00574   {140.00,      4.91953},
00575   {160.00,      5.05590},
00576   {180.00,      5.17588},
00577   {200.00,      5.28299},
00578   {250.00,      5.50925},
00579   {300.00,      5.69364},
00580   {350.00,      5.84926},
00581   {400.00,      5.98388},
00582   {450.00,      6.10252},
00583   {500.00,      6.20856},
00584   {600.00,      6.39189},
00585   {700.00,      6.54677},
00586   {800.00,      6.68084},
00587   {900.00,      6.79905},
00588   {1000.00,     6.90474}
00589 };
00590 
00591 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00592 
00593 // Data for L1 from: Mikkelsen H.H., Sigmund P. Phys. Rev. A40 (1989) 101-116
00594 const G4double G4ICRU73QOModel::L1[22][2] =
00595 {
00596   {0.00,       -0.000001},
00597   {0.10,       -0.00001},
00598   {0.20,       -0.00049},
00599   {0.30,       -0.00084},
00600   {0.40,        0.00085},
00601   {0.50,        0.00519},
00602   {0.60,        0.01198},
00603   {0.70,        0.02074},
00604   {0.80,        0.03133},
00605   {0.90,        0.04369},
00606   {1.00,        0.06035},
00607   {2.00,        0.24023},
00608   {3.00,        0.44284},
00609   {4.00,        0.62012},
00610   {5.00,        0.77031},
00611   {6.00,        0.90390},
00612   {7.00,        1.02705},
00613   {8.00,        1.10867},
00614   {9.00,        1.17546},
00615   {10.00,       1.21599},
00616   {15.00,       1.24349},
00617   {20.00,       1.16752}
00618 };
00619 
00620 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00621 
00622 // Data for L2 from: Mikkelsen H.H. Nucl. Instr. Meth. B58 (1991) 136-148
00623 const G4double G4ICRU73QOModel::L2[14][2] =
00624 {
00625   {0.00,        0.000001},
00626   {0.10,        0.00001},
00627   {0.20,        0.00000},
00628   {0.40,       -0.00120},
00629   {0.60,       -0.00036},
00630   {0.80,        0.00372},
00631   {1.00,        0.01298},
00632   {2.00,        0.08296},
00633   {4.00,        0.21953},
00634   {6.00,        0.23903},
00635   {8.00,        0.20893},
00636   {10.00,       0.10879},
00637   {20.00,      -0.88409},       
00638   {40.00,      -1.13902}
00639 };
00640 
00641 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00642 
00643 // Correction obtained by V.Ivanchenko using G4BetheBlochModel
00644 const G4double G4ICRU73QOModel::factorBethe[99] = { 1.0, 
00645 0.9637, 0.9872, 0.9469, 0.9875, 0.91, 0.989, 0.9507, 0.9773, 0.8621, 0.979,   // 1 - 10
00646 0.8357, 0.868, 0.9417, 0.9466, 0.8911, 0.905, 0.944, 0.9607, 0.928, 0.96,   // 11 - 20
00647 0.9098, 0.976, 0.8425, 0.8099, 0.7858, 0.947, 0.7248, 0.9106, 0.9246, 0.6821,   // 21 - 30
00648 0.7223, 0.9784, 0.774, 0.7953, 0.829, 0.9405, 0.8318, 0.8583, 0.8563, 0.8481,   // 31 - 40
00649 0.8207, 0.9033, 0.8063, 0.7837, 0.7818, 0.744, 0.875, 0.7693, 0.7871, 0.8459,   // 41 - 50
00650 0.8231, 0.8462, 0.853, 0.8736, 0.856, 0.8762, 0.8629, 0.8323, 0.8064, 0.7828,   // 51 - 60
00651 0.7533, 0.7273, 0.7093, 0.7157, 0.6823, 0.6612, 0.6418, 0.6395, 0.6323, 0.6221,   // 61 - 70
00652 0.6497, 0.6746, 0.8568, 0.8541, 0.6958, 0.6962, 0.7051, 0.863, 0.8588, 0.7226,   // 71 - 80
00653 0.7454, 0.78, 0.7783, 0.7996, 0.8216, 0.8632, 0.8558, 0.8792, 0.8745, 0.8676,   // 81 - 90
00654 0.8321, 0.8272, 0.7999, 0.7934, 0.7787, 0.7851, 0.7692, 0.7598}; 

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