[968] | 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 | // |
---|
| 27 | // |
---|
| 28 | // =========================================================================== |
---|
| 29 | // GEANT4 class |
---|
| 30 | // |
---|
| 31 | // Class: G4IonParametrisedLossModel |
---|
| 32 | // |
---|
| 33 | // Base class: G4VEmModel (utils) |
---|
| 34 | // |
---|
| 35 | // Author: Anton Lechner (Anton.Lechner@cern.ch) |
---|
| 36 | // |
---|
| 37 | // First implementation: 10. 11. 2008 |
---|
| 38 | // |
---|
[1055] | 39 | // Modifications: 03. 02. 2009 - Bug fix iterators (AL) |
---|
| 40 | // 11. 03. 2009 - Introduced new table handler (G4IonDEDXHandler) |
---|
| 41 | // and modified method to add/remove tables |
---|
| 42 | // (tables are now built in initialisation phase), |
---|
| 43 | // Minor bug fix in ComputeDEDXPerVolume (AL) |
---|
[1196] | 44 | // 20. 11. 2009 - Added set-method for energy loss limit (AL) |
---|
[968] | 45 | // |
---|
| 46 | // Class description: |
---|
| 47 | // Model for computing the energy loss of ions by employing a |
---|
| 48 | // parameterisation of dE/dx tables (default ICRU 73 tables). For |
---|
| 49 | // ion-material combinations and/or projectile energies not covered |
---|
| 50 | // by this model, the G4BraggIonModel and G4BetheBloch models are |
---|
| 51 | // employed. |
---|
| 52 | // |
---|
| 53 | // Comments: |
---|
| 54 | // |
---|
| 55 | // =========================================================================== |
---|
| 56 | |
---|
| 57 | |
---|
| 58 | inline G4double G4IonParametrisedLossModel::DeltaRayMeanEnergyTransferRate( |
---|
| 59 | const G4Material* material, |
---|
| 60 | const G4ParticleDefinition* particle, |
---|
| 61 | G4double kineticEnergy, |
---|
| 62 | G4double cutEnergy) { |
---|
| 63 | |
---|
| 64 | // ############## Mean energy transferred to delta-rays ################### |
---|
| 65 | // Computes the mean energy transfered to delta-rays per unit length, |
---|
| 66 | // considering only delta-rays with energies above the energy threshold |
---|
| 67 | // (energy cut) |
---|
| 68 | // |
---|
| 69 | // The mean energy transfer rate is derived by using the differential |
---|
| 70 | // cross section given in the references below. |
---|
| 71 | // |
---|
| 72 | // See Geant4 physics reference manual (version 9.1), section 9.1.3 |
---|
| 73 | // |
---|
| 74 | // Ref.: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1. |
---|
| 75 | // B. Rossi, High energy particles, New York, NY: Prentice-Hall (1952). |
---|
| 76 | // |
---|
| 77 | // (Implementation adapted from G4BraggIonModel) |
---|
| 78 | |
---|
| 79 | |
---|
| 80 | // *** Variables: |
---|
| 81 | // kineticEnergy = kinetic energy of projectile |
---|
| 82 | // totEnergy = total energy of projectile, i.e. kinetic energy |
---|
| 83 | // plus rest energy (Mc^2) |
---|
| 84 | // betaSquared = beta of projectile squared, calculated as |
---|
| 85 | // beta^2 = 1 - 1 / (E/Mc^2)^2 |
---|
| 86 | // = T * ( E + Mc^2 ) / E^2 |
---|
| 87 | // where T = kineticEnergy, E = totEnergy |
---|
| 88 | // cutEnergy = energy threshold for secondary particle production |
---|
| 89 | // i.e. energy cut, below which energy transfered to |
---|
| 90 | // electrons is treated as continuous loss of projectile |
---|
| 91 | // maxKinEnergy = maximum energy transferable to secondary electrons |
---|
| 92 | // meanRate = mean kinetic energy of delta ray (per unit length) |
---|
| 93 | // (above cutEnergy) |
---|
| 94 | |
---|
| 95 | G4double meanRate = 0.0; |
---|
| 96 | |
---|
| 97 | G4double maxKinEnergy = MaxSecondaryEnergy(particle, kineticEnergy); |
---|
| 98 | |
---|
| 99 | if (cutEnergy < maxKinEnergy) { |
---|
| 100 | |
---|
| 101 | G4double totalEnergy = kineticEnergy + cacheMass; |
---|
| 102 | G4double betaSquared = kineticEnergy * |
---|
| 103 | (totalEnergy + cacheMass) / (totalEnergy * totalEnergy); |
---|
| 104 | |
---|
| 105 | G4double cutMaxEnergyRatio = cutEnergy / maxKinEnergy; |
---|
| 106 | |
---|
| 107 | meanRate = |
---|
| 108 | (- std::log(cutMaxEnergyRatio) - (1.0 - cutMaxEnergyRatio) * betaSquared) * |
---|
| 109 | twopi_mc2_rcl2 * |
---|
| 110 | (material->GetTotNbOfElectPerVolume()) / betaSquared; |
---|
| 111 | |
---|
| 112 | meanRate *= GetChargeSquareRatio(particle, material, kineticEnergy); |
---|
| 113 | } |
---|
| 114 | |
---|
| 115 | return meanRate; |
---|
| 116 | } |
---|
| 117 | |
---|
| 118 | |
---|
| 119 | inline |
---|
| 120 | G4double G4IonParametrisedLossModel::MaxSecondaryEnergy( |
---|
| 121 | const G4ParticleDefinition* particle, |
---|
| 122 | G4double kineticEnergy) { |
---|
| 123 | |
---|
| 124 | // ############## Maximum energy of secondaries ########################## |
---|
| 125 | // Function computes maximum energy of secondary electrons which are |
---|
| 126 | // released by an ion |
---|
| 127 | // |
---|
| 128 | // See Geant4 physics reference manual (version 9.1), section 9.1.1 |
---|
| 129 | // |
---|
| 130 | // Ref.: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1. |
---|
| 131 | // C.Caso et al. (Part. Data Group), Europ. Phys. Jour. C 3 1 (1998). |
---|
| 132 | // B. Rossi, High energy particles, New York, NY: Prentice-Hall (1952). |
---|
| 133 | // |
---|
| 134 | // (Implementation adapted from G4BraggIonModel) |
---|
| 135 | |
---|
| 136 | if(particle != cacheParticle) UpdateCache(particle); |
---|
| 137 | |
---|
| 138 | G4double tau = kineticEnergy/cacheMass; |
---|
| 139 | G4double tmax = 2.0 * electron_mass_c2 * tau * (tau + 2.) / |
---|
| 140 | (1. + 2.0 * (tau + 1.) * cacheElecMassRatio + |
---|
| 141 | cacheElecMassRatio * cacheElecMassRatio); |
---|
| 142 | |
---|
| 143 | return tmax; |
---|
| 144 | } |
---|
| 145 | |
---|
| 146 | |
---|
| 147 | inline |
---|
| 148 | void G4IonParametrisedLossModel::UpdateCache( |
---|
| 149 | const G4ParticleDefinition* particle) { |
---|
| 150 | |
---|
| 151 | cacheParticle = particle; |
---|
| 152 | cacheMass = particle -> GetPDGMass(); |
---|
| 153 | cacheElecMassRatio = electron_mass_c2 / cacheMass; |
---|
| 154 | G4double q = particle -> GetPDGCharge() / eplus; |
---|
| 155 | cacheChargeSquare = q * q; |
---|
| 156 | } |
---|
| 157 | |
---|
| 158 | |
---|
| 159 | inline |
---|
| 160 | G4double G4IonParametrisedLossModel::GetChargeSquareRatio( |
---|
| 161 | const G4ParticleDefinition* particle, |
---|
| 162 | const G4Material* material, |
---|
| 163 | G4double kineticEnergy) { // Kinetic energy |
---|
| 164 | |
---|
| 165 | G4double chargeSquareRatio = corrections -> |
---|
| 166 | EffectiveChargeSquareRatio(particle, |
---|
| 167 | material, |
---|
| 168 | kineticEnergy); |
---|
| 169 | corrFactor = chargeSquareRatio * |
---|
| 170 | corrections -> EffectiveChargeCorrection(particle, |
---|
| 171 | material, |
---|
| 172 | kineticEnergy); |
---|
| 173 | return corrFactor; |
---|
| 174 | } |
---|
| 175 | |
---|
| 176 | |
---|
| 177 | inline |
---|
| 178 | G4double G4IonParametrisedLossModel::GetParticleCharge( |
---|
| 179 | const G4ParticleDefinition* particle, |
---|
| 180 | const G4Material* material, |
---|
| 181 | G4double kineticEnergy) { // Kinetic energy |
---|
| 182 | |
---|
| 183 | return corrections -> GetParticleCharge(particle, material, kineticEnergy); |
---|
| 184 | } |
---|
| 185 | |
---|
| 186 | |
---|
| 187 | inline |
---|
| 188 | LossTableList::iterator G4IonParametrisedLossModel::IsApplicable( |
---|
| 189 | const G4ParticleDefinition* particle, // Projectile (ion) |
---|
| 190 | const G4Material* material) { // Target material |
---|
| 191 | |
---|
[1055] | 192 | LossTableList::iterator iter = lossTableList.end(); |
---|
[968] | 193 | LossTableList::iterator iterTables = lossTableList.begin(); |
---|
| 194 | LossTableList::iterator iterTables_end = lossTableList.end(); |
---|
| 195 | |
---|
| 196 | for(;iterTables != iterTables_end; iterTables++) { |
---|
| 197 | G4bool isApplicable = (*iterTables) -> |
---|
| 198 | IsApplicable(particle, material); |
---|
| 199 | if(isApplicable) { |
---|
| 200 | iter = iterTables; |
---|
| 201 | break; |
---|
| 202 | } |
---|
| 203 | } |
---|
| 204 | |
---|
| 205 | return iter; |
---|
| 206 | } |
---|
[1196] | 207 | |
---|
| 208 | |
---|
| 209 | inline |
---|
| 210 | void G4IonParametrisedLossModel::SetEnergyLossLimit( |
---|
| 211 | G4double ionEnergyLossLimit) { |
---|
| 212 | |
---|
| 213 | if(ionEnergyLossLimit > 0 && ionEnergyLossLimit <=1) { |
---|
| 214 | |
---|
| 215 | energyLossLimit = ionEnergyLossLimit; |
---|
| 216 | } |
---|
| 217 | } |
---|