G4ICRU73QOModel.hh

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 header 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 // Class Description:
00043 //
00044 // Quantum Harmonic Oscillator Model for energy loss using atomic shell 
00045 // structure of atoms taking into account Q^2 (main for projectile charge Q), 
00046 // Q^3 and Q^4 terms for computation of energy loss due to binary collisions. 
00047 // Can be applied on heavy negatively charged particles for the energy interval 
00048 // 10 keV - 10 MeV scaled to the proton mass.
00049 //
00050 // Used data and formula of 
00051 // 1. G4QAOLowEnergyLoss class, S.Chauvie, P.Nieminen, M.G.Pia. IEEE Trans. 
00052 //    Nucl. Sci. 54 (2007) 578.
00053 // 2. ShellStrength and ShellEnergy from ICRU'73 Report 2005,
00054 // 3. Data for Ta (Z=73) from P.Sigmund, A.Shinner. Eur. Phys. J. D15 (2001) 
00055 //    165-172
00056 //
00057 // -------------------------------------------------------------------
00058 //
00059 
00060 #ifndef G4ICRU73QOModel_h
00061 #define G4ICRU73QOModel_h 1
00062 
00063 #include <CLHEP/Units/PhysicalConstants.h>
00064 
00065 #include "G4VEmModel.hh"
00066 #include "G4AtomicShells.hh"
00067 #include "G4DensityEffectData.hh"
00068 
00069 class G4ParticleChangeForLoss;
00070 
00071 class G4ICRU73QOModel : public G4VEmModel
00072 {
00073 
00074 public:
00075 
00076   G4ICRU73QOModel(const G4ParticleDefinition* p = 0,
00077                   const G4String& nam = "ICRU73QO");
00078 
00079   virtual ~G4ICRU73QOModel();
00080 
00081   virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
00082 
00083   virtual G4double ComputeCrossSectionPerElectron(
00084                                  const G4ParticleDefinition*,
00085                                  G4double kineticEnergy,
00086                                  G4double cutEnergy,
00087                                  G4double maxEnergy);
00088                                  
00089   virtual G4double ComputeCrossSectionPerAtom(
00090                                  const G4ParticleDefinition*,
00091                                  G4double kineticEnergy,
00092                                  G4double Z, G4double A,
00093                                  G4double cutEnergy,
00094                                  G4double maxEnergy);
00095                                                                  
00096   virtual G4double CrossSectionPerVolume(const G4Material*,
00097                                  const G4ParticleDefinition*,
00098                                  G4double kineticEnergy,
00099                                  G4double cutEnergy,
00100                                  G4double maxEnergy);
00101                                  
00102   virtual G4double ComputeDEDXPerVolume(const G4Material*,
00103                                         const G4ParticleDefinition*,
00104                                         G4double kineticEnergy,
00105                                         G4double);
00106 
00107   virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
00108                                  const G4MaterialCutsCouple*,
00109                                  const G4DynamicParticle*,
00110                                  G4double tmin,
00111                                  G4double maxEnergy);
00112 
00113   // add correction to energy loss and compute non-ionizing energy loss
00114   virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
00115                                     const G4DynamicParticle*,
00116                                     G4double& eloss,
00117                                     G4double& niel,
00118                                     G4double length);
00119 
00120 protected:
00121 
00122   virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
00123                                       G4double kinEnergy);
00124 
00125 private:
00126 
00127   inline void SetParticle(const G4ParticleDefinition* p);
00128   inline void SetLowestKinEnergy(const G4double val);
00129 
00130   G4double DEDX(const G4Material* material, G4double kineticEnergy);
00131 
00132   G4double DEDXPerElement(G4int Z, G4double kineticEnergy);
00133 
00134   // hide assignment operator
00135   G4ICRU73QOModel & operator=(const  G4ICRU73QOModel &right);
00136   G4ICRU73QOModel(const  G4ICRU73QOModel&);
00137 
00138   const G4ParticleDefinition* particle;
00139   G4ParticleDefinition*       theElectron;   
00140   G4ParticleChangeForLoss*    fParticleChange;
00141   G4DensityEffectData*        denEffData;
00142 
00143   G4double mass;
00144   G4double charge;
00145   G4double chargeSquare;
00146   G4double massRate;
00147   G4double ratio;
00148   G4double lowestKinEnergy;
00149 
00150   G4bool   isInitialised;
00151 
00152   // get number of shell, energy and oscillator strenghts for material
00153   G4int GetNumberOfShells(G4int Z) const;
00154 
00155   G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const; 
00156   G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const; 
00157   G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const;
00158 
00159   // calculate stopping number for L's term
00160   G4double GetL0(G4double normEnergy) const;
00161   // terms in Z^2
00162   G4double GetL1(G4double normEnergy) const;
00163   // terms in Z^3
00164   G4double GetL2(G4double normEnergy) const;
00165   // terms in Z^4
00166   
00167 
00168   // Z of element at now avaliable for the model
00169   static const G4int NQOELEM  = 26;
00170   static const G4int NQODATA  = 130;
00171   static const G4int ZElementAvailable[NQOELEM];
00172   
00173   // number, energy and oscillator strenghts
00174   // for an harmonic oscillator model of material
00175   static const G4int startElemIndex[NQOELEM];
00176   static const G4int nbofShellsForElement[NQOELEM];
00177   static const G4double ShellEnergy[NQODATA];
00178   static const G4double SubShellOccupation[NQODATA];  // Z * ShellStrength
00179 
00180   G4int indexZ[100];
00181 
00182   //  variable for calculation of stopping number of L's term
00183   static const G4double L0[67][2];
00184   static const G4double L1[22][2];
00185   static const G4double L2[14][2];
00186   
00187   G4int sizeL0;
00188   G4int sizeL1;
00189   G4int sizeL2;
00190 
00191   static const G4double factorBethe[99];
00192   
00193 };
00194 
00195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00196 
00197 inline void G4ICRU73QOModel::SetParticle(const G4ParticleDefinition* p)
00198 {
00199   particle = p;
00200   mass = particle->GetPDGMass();
00201   charge = particle->GetPDGCharge()/CLHEP::eplus;
00202   chargeSquare = charge*charge;
00203   massRate     = mass/CLHEP::proton_mass_c2;
00204   ratio = CLHEP::electron_mass_c2/mass;
00205 }
00206 
00207 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00208 
00209 inline G4int G4ICRU73QOModel::GetNumberOfShells(G4int Z) const
00210 {
00211   G4int nShell = 0;
00212 
00213   if(indexZ[Z] >= 0) { nShell = nbofShellsForElement[indexZ[Z]]; 
00214   } else { nShell = G4AtomicShells::GetNumberOfShells(Z); }
00215 
00216   return nShell;
00217 }
00218 
00219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00220 
00221 inline G4double 
00222 G4ICRU73QOModel::GetShellEnergy(G4int Z, G4int nbOfTheShell) const
00223 {
00224   G4double shellEnergy = 0.;
00225 
00226   G4int idx = indexZ[Z];
00227 
00228   if(idx >= 0) { shellEnergy = ShellEnergy[startElemIndex[idx] + nbOfTheShell]*CLHEP::eV; 
00229   } else { shellEnergy = GetOscillatorEnergy(Z, nbOfTheShell); }
00230 
00231   return  shellEnergy;
00232 }
00233 
00234 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00235 
00236 inline G4double 
00237 G4ICRU73QOModel::GetShellStrength(G4int Z, G4int nbOfTheShell) const
00238 {
00239   G4double shellStrength = 0.;
00240 
00241   G4int idx = indexZ[Z];
00242 
00243   if(idx >= 0) { shellStrength = SubShellOccupation[startElemIndex[idx] + nbOfTheShell] / Z; 
00244   } else { shellStrength = G4double(G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell))/Z; }
00245   
00246   return shellStrength;
00247 }
00248 
00249 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00250 
00251 inline void G4ICRU73QOModel::SetLowestKinEnergy(const G4double val)
00252 {
00253   lowestKinEnergy = val;
00254 }
00255 
00256 #endif

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