[819] | 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 | // |
---|
[1055] | 26 | // $Id: G4eCoulombScatteringModel.hh,v 1.43 2009/05/10 16:09:29 vnivanch Exp $ |
---|
| 27 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[819] | 28 | // |
---|
| 29 | // ------------------------------------------------------------------- |
---|
| 30 | // |
---|
| 31 | // GEANT4 Class header file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4eCoulombScatteringModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko |
---|
| 37 | // |
---|
| 38 | // Creation date: 19.02.2006 |
---|
| 39 | // |
---|
| 40 | // Modifications: |
---|
| 41 | // 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the |
---|
| 42 | // logic of building - only elements from G4ElementTable |
---|
| 43 | // 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim |
---|
| 44 | // 19.08.06 V.Ivanchenko add inline function ScreeningParameter and |
---|
| 45 | // make some members protected |
---|
| 46 | // 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e- |
---|
[961] | 47 | // 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion |
---|
[819] | 48 | // |
---|
| 49 | // Class Description: |
---|
| 50 | // |
---|
| 51 | // Implementation of eCoulombScattering of pointlike charge particle |
---|
| 52 | // on Atomic Nucleus for interval of scattering anles in Lab system |
---|
| 53 | // thetaMin - ThetaMax, nucleus recoil is neglected. |
---|
| 54 | // The model based on analysis of J.M.Fernandez-Varea et al. |
---|
| 55 | // NIM B73(1993)447 originated from G.Wentzel Z.Phys. 40(1927)590 with |
---|
| 56 | // screening parameter from H.A.Bethe Phys. Rev. 89 (1953) 1256. |
---|
| 57 | // |
---|
| 58 | |
---|
| 59 | // ------------------------------------------------------------------- |
---|
| 60 | // |
---|
| 61 | |
---|
| 62 | #ifndef G4eCoulombScatteringModel_h |
---|
| 63 | #define G4eCoulombScatteringModel_h 1 |
---|
| 64 | |
---|
| 65 | #include "G4VEmModel.hh" |
---|
| 66 | #include "G4PhysicsTable.hh" |
---|
[1055] | 67 | #include "globals.hh" |
---|
[819] | 68 | #include "G4NistManager.hh" |
---|
| 69 | |
---|
| 70 | class G4ParticleChangeForGamma; |
---|
| 71 | class G4ParticleDefinition; |
---|
| 72 | |
---|
| 73 | class G4eCoulombScatteringModel : public G4VEmModel |
---|
| 74 | { |
---|
| 75 | |
---|
| 76 | public: |
---|
| 77 | |
---|
[961] | 78 | G4eCoulombScatteringModel(const G4String& nam = "eCoulombScattering"); |
---|
[819] | 79 | |
---|
| 80 | virtual ~G4eCoulombScatteringModel(); |
---|
| 81 | |
---|
| 82 | virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&); |
---|
| 83 | |
---|
| 84 | virtual G4double ComputeCrossSectionPerAtom( |
---|
| 85 | const G4ParticleDefinition*, |
---|
| 86 | G4double kinEnergy, |
---|
| 87 | G4double Z, |
---|
| 88 | G4double A, |
---|
| 89 | G4double cut, |
---|
| 90 | G4double emax); |
---|
| 91 | |
---|
| 92 | virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*, |
---|
| 93 | const G4MaterialCutsCouple*, |
---|
| 94 | const G4DynamicParticle*, |
---|
| 95 | G4double tmin, |
---|
| 96 | G4double maxEnergy); |
---|
| 97 | |
---|
[961] | 98 | inline void SetRecoilThreshold(G4double eth); |
---|
| 99 | |
---|
[819] | 100 | protected: |
---|
| 101 | |
---|
[961] | 102 | G4double CrossSectionPerAtom(); |
---|
[819] | 103 | |
---|
[961] | 104 | G4double SampleCosineTheta(); |
---|
[819] | 105 | |
---|
[961] | 106 | inline void DefineMaterial(const G4MaterialCutsCouple*); |
---|
| 107 | |
---|
[819] | 108 | inline void SetupParticle(const G4ParticleDefinition*); |
---|
| 109 | |
---|
[961] | 110 | inline void SetupKinematic(G4double kinEnergy, G4double cut); |
---|
[819] | 111 | |
---|
[961] | 112 | inline void SetupTarget(G4double Z, G4double kinEnergy); |
---|
[819] | 113 | |
---|
| 114 | private: |
---|
| 115 | |
---|
[961] | 116 | void ComputeMaxElectronScattering(G4double cut); |
---|
| 117 | |
---|
[819] | 118 | // hide assignment operator |
---|
| 119 | G4eCoulombScatteringModel & operator=(const G4eCoulombScatteringModel &right); |
---|
| 120 | G4eCoulombScatteringModel(const G4eCoulombScatteringModel&); |
---|
| 121 | |
---|
| 122 | protected: |
---|
[961] | 123 | |
---|
[819] | 124 | const G4ParticleDefinition* theProton; |
---|
| 125 | const G4ParticleDefinition* theElectron; |
---|
| 126 | const G4ParticleDefinition* thePositron; |
---|
| 127 | |
---|
[961] | 128 | G4ParticleTable* theParticleTable; |
---|
[819] | 129 | G4ParticleChangeForGamma* fParticleChange; |
---|
| 130 | G4NistManager* fNistManager; |
---|
[961] | 131 | const G4DataVector* currentCuts; |
---|
[819] | 132 | |
---|
[961] | 133 | const G4MaterialCutsCouple* currentCouple; |
---|
| 134 | const G4Material* currentMaterial; |
---|
| 135 | const G4Element* currentElement; |
---|
| 136 | G4int currentMaterialIndex; |
---|
| 137 | |
---|
[819] | 138 | G4double coeff; |
---|
| 139 | G4double cosThetaMin; |
---|
| 140 | G4double cosThetaMax; |
---|
[961] | 141 | G4double cosTetMinNuc; |
---|
[819] | 142 | G4double cosTetMaxNuc; |
---|
[961] | 143 | G4double cosTetMaxNuc2; |
---|
[819] | 144 | G4double cosTetMaxElec; |
---|
[961] | 145 | G4double cosTetMaxElec2; |
---|
[819] | 146 | G4double q2Limit; |
---|
[961] | 147 | G4double recoilThreshold; |
---|
[819] | 148 | G4double elecXSection; |
---|
| 149 | G4double nucXSection; |
---|
| 150 | G4double ecut; |
---|
| 151 | |
---|
| 152 | // projectile |
---|
| 153 | const G4ParticleDefinition* particle; |
---|
| 154 | |
---|
| 155 | G4double chargeSquare; |
---|
| 156 | G4double spin; |
---|
| 157 | G4double mass; |
---|
| 158 | G4double tkin; |
---|
| 159 | G4double mom2; |
---|
| 160 | G4double invbeta2; |
---|
[1055] | 161 | G4double kinFactor; |
---|
[961] | 162 | G4double etag; |
---|
| 163 | G4double lowEnergyLimit; |
---|
[819] | 164 | |
---|
| 165 | // target |
---|
| 166 | G4double targetZ; |
---|
[1055] | 167 | G4double targetMass; |
---|
[819] | 168 | G4double screenZ; |
---|
| 169 | G4double formfactA; |
---|
[961] | 170 | G4int idxelm; |
---|
[1055] | 171 | G4int iz; |
---|
[819] | 172 | |
---|
| 173 | private: |
---|
| 174 | |
---|
| 175 | G4double alpha2; |
---|
| 176 | G4double faclim; |
---|
| 177 | |
---|
[1055] | 178 | static G4double ScreenRSquare[100]; |
---|
| 179 | static G4double FormFactor[100]; |
---|
| 180 | |
---|
[819] | 181 | G4bool isInitialised; |
---|
| 182 | }; |
---|
| 183 | |
---|
| 184 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 185 | |
---|
[961] | 186 | inline |
---|
| 187 | void G4eCoulombScatteringModel::DefineMaterial(const G4MaterialCutsCouple* cup) |
---|
| 188 | { |
---|
| 189 | if(cup != currentCouple) { |
---|
| 190 | currentCouple = cup; |
---|
| 191 | currentMaterial = cup->GetMaterial(); |
---|
| 192 | currentMaterialIndex = currentCouple->GetIndex(); |
---|
| 193 | } |
---|
| 194 | } |
---|
| 195 | |
---|
| 196 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 197 | |
---|
[819] | 198 | inline |
---|
| 199 | void G4eCoulombScatteringModel::SetupParticle(const G4ParticleDefinition* p) |
---|
| 200 | { |
---|
| 201 | // Initialise mass and charge |
---|
| 202 | if(p != particle) { |
---|
| 203 | particle = p; |
---|
| 204 | mass = particle->GetPDGMass(); |
---|
| 205 | spin = particle->GetPDGSpin(); |
---|
| 206 | G4double q = particle->GetPDGCharge()/eplus; |
---|
| 207 | chargeSquare = q*q; |
---|
| 208 | tkin = 0.0; |
---|
| 209 | } |
---|
| 210 | } |
---|
| 211 | |
---|
| 212 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 213 | |
---|
[961] | 214 | inline void G4eCoulombScatteringModel::SetupKinematic(G4double ekin, |
---|
| 215 | G4double cut) |
---|
[819] | 216 | { |
---|
[961] | 217 | if(ekin != tkin || ecut != cut) { |
---|
| 218 | tkin = ekin; |
---|
| 219 | mom2 = tkin*(tkin + 2.0*mass); |
---|
[819] | 220 | invbeta2 = 1.0 + mass*mass/mom2; |
---|
[961] | 221 | cosTetMinNuc = cosThetaMin; |
---|
| 222 | cosTetMaxNuc = cosThetaMax; |
---|
[1055] | 223 | if(mass < MeV && cosThetaMin < 1.0 && ekin <= 10.*cut) { |
---|
[961] | 224 | cosTetMinNuc = ekin*(cosThetaMin + 1.0)/(10.*cut) - 1.0; |
---|
| 225 | } |
---|
| 226 | ComputeMaxElectronScattering(cut); |
---|
| 227 | } |
---|
[819] | 228 | } |
---|
| 229 | |
---|
| 230 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 231 | |
---|
[961] | 232 | inline void G4eCoulombScatteringModel::SetupTarget(G4double Z, G4double e) |
---|
[819] | 233 | { |
---|
[961] | 234 | if(Z != targetZ || e != etag) { |
---|
| 235 | etag = e; |
---|
[819] | 236 | targetZ = Z; |
---|
[1055] | 237 | iz= G4int(Z); |
---|
[961] | 238 | if(iz > 99) iz = 99; |
---|
[1055] | 239 | targetMass = fNistManager->GetAtomicMassAmu(iz)*amu_c2; |
---|
| 240 | G4double m12 = mass*mass; |
---|
| 241 | G4double x = 1.0 + mass/targetMass; |
---|
| 242 | kinFactor = (1.0 + m12/mom2)/mom2; |
---|
| 243 | |
---|
| 244 | screenZ = ScreenRSquare[iz]/mom2; |
---|
| 245 | if(iz > 1) { |
---|
| 246 | screenZ *=(1.13 + 3.76*Z*Z*alpha2); |
---|
| 247 | kinFactor /= (x*x); |
---|
[961] | 248 | } |
---|
[1055] | 249 | // if(iz > 1) screenZ *=(1.13 + std::min(0.5,3.76*Z*Z*invbeta2*alpha2)); |
---|
| 250 | formfactA = FormFactor[iz]*mom2; |
---|
[961] | 251 | cosTetMaxNuc2 = cosTetMaxNuc; |
---|
[1055] | 252 | if(1 == iz && particle == theProton && cosTetMaxNuc2 < 0.0) { |
---|
[961] | 253 | cosTetMaxNuc2 = 0.0; |
---|
| 254 | } |
---|
| 255 | cosTetMaxElec2 = cosTetMaxElec; |
---|
[819] | 256 | } |
---|
| 257 | } |
---|
| 258 | |
---|
| 259 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 260 | |
---|
[961] | 261 | inline void G4eCoulombScatteringModel::SetRecoilThreshold(G4double eth) |
---|
| 262 | { |
---|
| 263 | recoilThreshold = eth; |
---|
| 264 | } |
---|
| 265 | |
---|
| 266 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 267 | |
---|
[819] | 268 | #endif |
---|