Geant4 Cross Reference |
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 header file 30 // 31 // Class: G4IonDEDXScalingICRU73 32 // 33 // Base class: G4VIonDEDXScalingAlgorithm 34 // 35 // Author: Anton Lechner (Anton.Lechner@cern.ch) 36 // 37 // First implementation: 10. 05. 2009 38 // 39 // Modifications: 12. 11. 2009 - Moved all decision logic concerning ICRU 73 40 // scaling for heavy ions into this class. 41 // Adapting ScalingFactorEnergy class according 42 // to changes in base class (AL). 43 // 44 // Class description: 45 // dE/dx scaling algorithm applied on top of ICRU 73 data (for ions not 46 // covered by the ICRU 73 report) 47 // 48 // Comments: 49 // 50 // =========================================================================== 51 52 #ifndef G4IONDEDXSCALINGICRU73_HH 53 #define G4IONDEDXSCALINGICRU73_HH 54 55 #include "globals.hh" 56 #include "G4VIonDEDXScalingAlgorithm.hh" 57 #include "G4Material.hh" 58 #include "G4ParticleDefinition.hh" 59 #include <vector> 60 #include "G4Exp.hh" 61 62 class G4IonDEDXScalingICRU73 : public G4VIonDEDXScalingAlgorithm { 63 64 public: 65 explicit G4IonDEDXScalingICRU73(G4int minAtomicNumberIon = 19, 66 G4int maxAtomicNumberIon = 102); 67 ~G4IonDEDXScalingICRU73(); 68 69 // Function for scaling the kinetic energy (no scaling by default). 70 // Returns scaling factor for a given ion. 71 G4double ScalingFactorEnergy( 72 const G4ParticleDefinition* particle, // Projectile (ion) 73 const G4Material* material) override; // Target material 74 75 76 // Function for scaling the dE/dx value (no scaling by default). 77 // Returns scaling factor for a given ion-material couple and 78 // a given kinetic energy. 79 G4double ScalingFactorDEDX( 80 const G4ParticleDefinition* particle, // Projectile (ion) 81 const G4Material*, // Target material 82 G4double kineticEnergy) override; // Kinetic energy 83 84 85 // Function for defining a base particle for dE/dx calculation. 86 // (no base particle by default). Returns atomic number of base 87 // particle. 88 G4int AtomicNumberBaseIon( 89 G4int atomicNumberIon, // Atomic number of ion 90 const G4Material*) override; // Target material 91 92 private: 93 void UpdateCacheParticle( 94 const G4ParticleDefinition* particle); // Projectile (ion) 95 96 void UpdateCacheMaterial( 97 const G4Material* material); // Target material 98 99 void CreateReferenceParticles(); 100 101 G4double EquilibriumCharge( 102 G4double mass, // Ion mass 103 G4double charge, // Ion charge 104 G4double atomicNumberPow, // Power of atomic nmb 105 G4double kineticEnergy); // Kinetic energy 106 107 // Scaling is only applied for ions with atomic numbers in the range 108 // defined by the following parameters: 109 G4int minAtomicNumber; 110 G4int maxAtomicNumber; 111 112 G4bool referencePrepared; 113 114 // Some properties of reference particle (Fe) are stored for faster access 115 ///////////////////////////G4ParticleDefinition* referenceFe; 116 G4int atomicNumberRefFe; 117 G4int massNumberRefFe; 118 G4double atomicNumberRefPow23Fe; 119 G4double chargeRefFe; 120 G4double massRefFe; 121 122 // Some properties of reference particle (Ar) are stored for faster access 123 ///////////////////////////G4ParticleDefinition* referenceAr; 124 G4int atomicNumberRefAr; 125 G4int massNumberRefAr; 126 G4double atomicNumberRefPow23Ar; 127 G4double chargeRefAr; 128 G4double massRefAr; 129 130 // Flag indicating the use of Fe ions as reference particles 131 G4bool useFe; 132 133 // Some properties of projectiles are stored for faster access 134 const G4ParticleDefinition* cacheParticle; 135 G4int cacheMassNumber; 136 G4int cacheAtomicNumber; 137 G4double cacheAtomicNumberPow23; 138 G4double cacheCharge; 139 G4double cacheMass; 140 141 // Material pointer 142 const G4Material* cacheMaterial; 143 }; 144 145 // ########################################################################### 146 147 inline void G4IonDEDXScalingICRU73::UpdateCacheParticle ( 148 const G4ParticleDefinition* particle) { // Projectile (ion) 149 150 if(particle != cacheParticle) { 151 152 cacheParticle = particle; 153 cacheAtomicNumber = particle -> GetAtomicNumber(); 154 cacheMassNumber = particle -> GetAtomicMass(); 155 cacheCharge = particle -> GetPDGCharge(); 156 cacheMass = particle -> GetPDGMass(); 157 cacheAtomicNumberPow23 = std::pow(G4double(cacheAtomicNumber), 2./3.); 158 } 159 } 160 161 // ########################################################################### 162 163 inline void G4IonDEDXScalingICRU73::UpdateCacheMaterial ( 164 const G4Material* material) { // Target material 165 166 if(cacheMaterial != material) { 167 168 cacheMaterial = material; 169 170 useFe = true; 171 172 size_t nmbElements = material -> GetNumberOfElements(); 173 if( nmbElements > 1 ) useFe = false; 174 175 if( material -> GetName() == "G4_WATER" ) useFe = true; 176 } 177 } 178 179 // ########################################################################### 180 181 inline G4double G4IonDEDXScalingICRU73::EquilibriumCharge( 182 G4double mass, 183 G4double charge, 184 G4double atomicNumberPow, 185 G4double kineticEnergy) { 186 187 G4double totalEnergy = kineticEnergy + mass; 188 G4double betaSquared = kineticEnergy * 189 (totalEnergy + mass) / (totalEnergy * totalEnergy); 190 191 G4double beta = std::sqrt( betaSquared ); 192 193 G4double velOverBohrVel = beta / CLHEP::fine_structure_const; 194 195 G4double q1 = 1.0 - G4Exp(-velOverBohrVel / atomicNumberPow); 196 197 return q1 * charge; 198 } 199 200 // ########################################################################### 201 202 #endif 203