[1316] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
| 3 | // * License and Disclaimer * |
---|
| 4 | // * * |
---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
| 7 | // * conditions of the Geant4 Software License, included in the file * |
---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
| 9 | // * include a list of copyright holders. * |
---|
| 10 | // * * |
---|
| 11 | // * Neither the authors of this software system, nor their employing * |
---|
| 12 | // * institutes,nor the agencies providing financial support for this * |
---|
| 13 | // * work make any representation or warranty, express or implied, * |
---|
| 14 | // * regarding this software system or assume any liability for its * |
---|
| 15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
| 16 | // * for the full disclaimer and the limitation of liability. * |
---|
| 17 | // * * |
---|
| 18 | // * This code implementation is the result of the scientific and * |
---|
| 19 | // * technical work of the GEANT4 collaboration. * |
---|
| 20 | // * By using, copying, modifying or distributing the software (or * |
---|
| 21 | // * any work based on the software) you agree to acknowledge its * |
---|
| 22 | // * use in resulting scientific publications, and indicate your * |
---|
| 23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 24 | // ******************************************************************** |
---|
| 25 | // |
---|
| 26 | // $Id: G4ICRU73QOModel.hh,v 1.3 2010/06/04 09:09:31 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ |
---|
| 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4ICRU73QOModel |
---|
| 35 | // |
---|
| 36 | // Author: Alexander Bagulya |
---|
| 37 | // |
---|
| 38 | // Creation date: 21.05.2010 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // |
---|
| 42 | // |
---|
| 43 | // Class Description: |
---|
| 44 | // |
---|
| 45 | // Quantum Harmonic Oscillator Model for energy loss using atomic shell |
---|
| 46 | // structure of atoms taking into account Q^2 (main for projectile charge Q), |
---|
| 47 | // Q^3 and Q^4 terms for computation of energy loss due to binary collisions. |
---|
| 48 | // Can be applied on heavy negatively charged particles for the energy interval |
---|
| 49 | // 10 keV - 10 MeV scaled to the proton mass. |
---|
| 50 | // |
---|
| 51 | // Used data and formula of |
---|
| 52 | // 1. G4QAOLowEnergyLoss class, S.Chauvie, P.Nieminen, M.G.Pia. IEEE Trans. |
---|
| 53 | // Nucl. Sci. 54 (2007) 578. |
---|
| 54 | // 2. ShellStrength and ShellEnergy from ICRU'73 Report 2005, |
---|
| 55 | // 3. Data for Ta (Z=73) from P.Sigmund, A.Shinner. Eur. Phys. J. D15 (2001) |
---|
| 56 | // 165-172 |
---|
| 57 | // |
---|
| 58 | // ------------------------------------------------------------------- |
---|
| 59 | // |
---|
| 60 | |
---|
| 61 | #ifndef G4ICRU73QOModel_h |
---|
| 62 | #define G4ICRU73QOModel_h 1 |
---|
| 63 | |
---|
| 64 | #include "G4VEmModel.hh" |
---|
| 65 | #include "G4AtomicShells.hh" |
---|
| 66 | #include "G4DensityEffectData.hh" |
---|
| 67 | |
---|
| 68 | class G4ParticleChangeForLoss; |
---|
| 69 | |
---|
| 70 | class G4ICRU73QOModel : public G4VEmModel |
---|
| 71 | { |
---|
| 72 | |
---|
| 73 | public: |
---|
| 74 | |
---|
| 75 | G4ICRU73QOModel(const G4ParticleDefinition* p = 0, |
---|
| 76 | const G4String& nam = "ICRU73QO"); |
---|
| 77 | |
---|
| 78 | virtual ~G4ICRU73QOModel(); |
---|
| 79 | |
---|
| 80 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
| 81 | |
---|
| 82 | virtual G4double MinEnergyCut(const G4ParticleDefinition*, |
---|
| 83 | const G4MaterialCutsCouple*); |
---|
| 84 | |
---|
| 85 | virtual G4double ComputeCrossSectionPerElectron( |
---|
| 86 | const G4ParticleDefinition*, |
---|
| 87 | G4double kineticEnergy, |
---|
| 88 | G4double cutEnergy, |
---|
| 89 | G4double maxEnergy); |
---|
| 90 | |
---|
| 91 | virtual G4double ComputeCrossSectionPerAtom( |
---|
| 92 | const G4ParticleDefinition*, |
---|
| 93 | G4double kineticEnergy, |
---|
| 94 | G4double Z, G4double A, |
---|
| 95 | G4double cutEnergy, |
---|
| 96 | G4double maxEnergy); |
---|
| 97 | |
---|
| 98 | virtual G4double CrossSectionPerVolume(const G4Material*, |
---|
| 99 | const G4ParticleDefinition*, |
---|
| 100 | G4double kineticEnergy, |
---|
| 101 | G4double cutEnergy, |
---|
| 102 | G4double maxEnergy); |
---|
| 103 | |
---|
| 104 | virtual G4double ComputeDEDXPerVolume(const G4Material*, |
---|
| 105 | const G4ParticleDefinition*, |
---|
| 106 | G4double kineticEnergy, |
---|
| 107 | G4double); |
---|
| 108 | |
---|
| 109 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 110 | const G4MaterialCutsCouple*, |
---|
| 111 | const G4DynamicParticle*, |
---|
| 112 | G4double tmin, |
---|
| 113 | G4double maxEnergy); |
---|
| 114 | |
---|
| 115 | // add correction to energy loss and compute non-ionizing energy loss |
---|
| 116 | virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*, |
---|
| 117 | const G4DynamicParticle*, |
---|
| 118 | G4double& eloss, |
---|
| 119 | G4double& niel, |
---|
| 120 | G4double length); |
---|
| 121 | |
---|
| 122 | protected: |
---|
| 123 | |
---|
| 124 | virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*, |
---|
| 125 | G4double kinEnergy); |
---|
| 126 | |
---|
| 127 | private: |
---|
| 128 | |
---|
| 129 | inline void SetParticle(const G4ParticleDefinition* p); |
---|
| 130 | inline void SetLowestKinEnergy(const G4double val); |
---|
| 131 | |
---|
| 132 | G4double DEDX(const G4Material* material, G4double kineticEnergy); |
---|
| 133 | |
---|
| 134 | G4double DEDXPerElement(G4int Z, G4double kineticEnergy); |
---|
| 135 | |
---|
| 136 | // hide assignment operator |
---|
| 137 | G4ICRU73QOModel & operator=(const G4ICRU73QOModel &right); |
---|
| 138 | G4ICRU73QOModel(const G4ICRU73QOModel&); |
---|
| 139 | |
---|
| 140 | const G4ParticleDefinition* particle; |
---|
| 141 | G4ParticleDefinition* theElectron; |
---|
| 142 | G4ParticleChangeForLoss* fParticleChange; |
---|
| 143 | G4DensityEffectData* denEffData; |
---|
| 144 | |
---|
| 145 | G4double mass; |
---|
| 146 | G4double q; |
---|
| 147 | G4double chargeSquare; |
---|
| 148 | G4double massRate; |
---|
| 149 | G4double ratio; |
---|
| 150 | G4double lowestKinEnergy; |
---|
| 151 | |
---|
| 152 | G4bool isInitialised; |
---|
| 153 | |
---|
| 154 | // get number of shell, energy and oscillator strenghts for material |
---|
| 155 | G4int GetNumberOfShells(G4int Z) const; |
---|
| 156 | |
---|
| 157 | G4double GetShellEnergy(G4int Z, G4int nbOfTheShell) const; |
---|
| 158 | G4double GetOscillatorEnergy(G4int Z, G4int nbOfTheShell) const; |
---|
| 159 | G4double GetShellStrength(G4int Z, G4int nbOfTheShell) const; |
---|
| 160 | |
---|
| 161 | // calculate stopping number for L's term |
---|
| 162 | G4double GetL0(G4double normEnergy) const; |
---|
| 163 | // terms in Z^2 |
---|
| 164 | G4double GetL1(G4double normEnergy) const; |
---|
| 165 | // terms in Z^3 |
---|
| 166 | G4double GetL2(G4double normEnergy) const; |
---|
| 167 | // terms in Z^4 |
---|
| 168 | |
---|
| 169 | |
---|
| 170 | // Z of element at now avaliable for the model |
---|
| 171 | static const G4int NQOELEM = 26; |
---|
| 172 | static const G4int NQODATA = 130; |
---|
| 173 | static const G4int ZElementAvailable[NQOELEM]; |
---|
| 174 | |
---|
| 175 | // number, energy and oscillator strenghts |
---|
| 176 | // for an harmonic oscillator model of material |
---|
| 177 | static const G4int startElemIndex[NQOELEM]; |
---|
| 178 | static const G4int nbofShellsForElement[NQOELEM]; |
---|
| 179 | static const G4double ShellEnergy[NQODATA]; |
---|
| 180 | static const G4double SubShellOccupation[NQODATA]; // Z * ShellStrength |
---|
| 181 | |
---|
| 182 | G4int indexZ[100]; |
---|
| 183 | |
---|
| 184 | // variable for calculation of stopping number of L's term |
---|
| 185 | static const G4double L0[67][2]; |
---|
| 186 | static const G4double L1[22][2]; |
---|
| 187 | static const G4double L2[14][2]; |
---|
| 188 | |
---|
| 189 | G4int sizeL0; |
---|
| 190 | G4int sizeL1; |
---|
| 191 | G4int sizeL2; |
---|
| 192 | |
---|
| 193 | static const G4double factorBethe[99]; |
---|
| 194 | |
---|
| 195 | }; |
---|
| 196 | |
---|
| 197 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 198 | |
---|
| 199 | inline void G4ICRU73QOModel::SetParticle(const G4ParticleDefinition* p) |
---|
| 200 | { |
---|
| 201 | particle = p; |
---|
| 202 | mass = particle->GetPDGMass(); |
---|
| 203 | q = particle->GetPDGCharge()/eplus; |
---|
| 204 | chargeSquare = q*q; |
---|
| 205 | massRate = mass/proton_mass_c2; |
---|
| 206 | ratio = electron_mass_c2/mass; |
---|
| 207 | } |
---|
| 208 | |
---|
| 209 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 210 | |
---|
| 211 | inline G4int G4ICRU73QOModel::GetNumberOfShells(G4int Z) const |
---|
| 212 | { |
---|
| 213 | G4int nShell = 0; |
---|
| 214 | |
---|
| 215 | if(indexZ[Z] >= 0) { nShell = nbofShellsForElement[indexZ[Z]]; |
---|
| 216 | } else { nShell = G4AtomicShells::GetNumberOfShells(Z); } |
---|
| 217 | |
---|
| 218 | return nShell; |
---|
| 219 | } |
---|
| 220 | |
---|
| 221 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 222 | |
---|
| 223 | inline G4double G4ICRU73QOModel::GetShellEnergy(G4int Z, |
---|
| 224 | G4int nbOfTheShell) const |
---|
| 225 | { |
---|
| 226 | G4double shellEnergy = 0.; |
---|
| 227 | |
---|
| 228 | G4int idx = indexZ[Z]; |
---|
| 229 | |
---|
| 230 | if(indexZ[Z] >= 0) { shellEnergy = ShellEnergy[startElemIndex[idx] + nbOfTheShell]*eV; |
---|
| 231 | } else { shellEnergy = GetOscillatorEnergy(Z, nbOfTheShell); } |
---|
| 232 | |
---|
| 233 | return shellEnergy; |
---|
| 234 | } |
---|
| 235 | |
---|
| 236 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 237 | |
---|
| 238 | inline G4double G4ICRU73QOModel::GetShellStrength(G4int Z, |
---|
| 239 | G4int nbOfTheShell) const |
---|
| 240 | { |
---|
| 241 | G4double shellStrength = 0.; |
---|
| 242 | |
---|
| 243 | G4int idx = indexZ[Z]; |
---|
| 244 | |
---|
| 245 | if(indexZ[Z] >= 0) { shellStrength = SubShellOccupation[startElemIndex[idx] + nbOfTheShell] / Z; |
---|
| 246 | } else { shellStrength = 1.0 / Z * G4AtomicShells::GetNumberOfElectrons(Z,nbOfTheShell); } |
---|
| 247 | |
---|
| 248 | return shellStrength; |
---|
| 249 | } |
---|
| 250 | |
---|
| 251 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 252 | |
---|
| 253 | inline void G4ICRU73QOModel::SetLowestKinEnergy(const G4double val) |
---|
| 254 | { |
---|
| 255 | lowestKinEnergy = val; |
---|
| 256 | } |
---|
| 257 | |
---|
| 258 | #endif |
---|