[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.cc,v 1.69 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 file |
---|
| 32 | // |
---|
| 33 | // |
---|
| 34 | // File name: G4eCoulombScatteringModel |
---|
| 35 | // |
---|
| 36 | // Author: Vladimir Ivanchenko |
---|
| 37 | // |
---|
| 38 | // Creation date: 22.08.2005 |
---|
| 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 |
---|
| 45 | // 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e- |
---|
[961] | 46 | // 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion |
---|
[819] | 47 | // |
---|
| 48 | // Class Description: |
---|
| 49 | // |
---|
| 50 | // ------------------------------------------------------------------- |
---|
| 51 | // |
---|
| 52 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 53 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 54 | |
---|
| 55 | #include "G4eCoulombScatteringModel.hh" |
---|
| 56 | #include "Randomize.hh" |
---|
| 57 | #include "G4DataVector.hh" |
---|
| 58 | #include "G4ElementTable.hh" |
---|
| 59 | #include "G4PhysicsLogVector.hh" |
---|
| 60 | #include "G4ParticleChangeForGamma.hh" |
---|
| 61 | #include "G4Electron.hh" |
---|
| 62 | #include "G4Positron.hh" |
---|
| 63 | #include "G4Proton.hh" |
---|
[961] | 64 | #include "G4ParticleTable.hh" |
---|
[819] | 65 | |
---|
| 66 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 67 | |
---|
[1055] | 68 | G4double G4eCoulombScatteringModel::ScreenRSquare[] = {0.0}; |
---|
| 69 | G4double G4eCoulombScatteringModel::FormFactor[] = {0.0}; |
---|
| 70 | |
---|
[819] | 71 | using namespace std; |
---|
| 72 | |
---|
[961] | 73 | G4eCoulombScatteringModel::G4eCoulombScatteringModel(const G4String& nam) |
---|
[819] | 74 | : G4VEmModel(nam), |
---|
[961] | 75 | cosThetaMin(1.0), |
---|
| 76 | cosThetaMax(-1.0), |
---|
| 77 | q2Limit(TeV*TeV), |
---|
[819] | 78 | alpha2(fine_structure_const*fine_structure_const), |
---|
| 79 | faclim(100.0), |
---|
| 80 | isInitialised(false) |
---|
| 81 | { |
---|
| 82 | fNistManager = G4NistManager::Instance(); |
---|
[961] | 83 | theParticleTable = G4ParticleTable::GetParticleTable(); |
---|
[819] | 84 | theElectron = G4Electron::Electron(); |
---|
| 85 | thePositron = G4Positron::Positron(); |
---|
| 86 | theProton = G4Proton::Proton(); |
---|
[961] | 87 | currentMaterial = 0; |
---|
| 88 | currentElement = 0; |
---|
[1055] | 89 | lowEnergyLimit = keV; |
---|
[819] | 90 | G4double p0 = electron_mass_c2*classic_electr_radius; |
---|
| 91 | coeff = twopi*p0*p0; |
---|
[961] | 92 | tkin = targetZ = mom2 = DBL_MIN; |
---|
[819] | 93 | elecXSection = nucXSection = 0.0; |
---|
[1055] | 94 | recoilThreshold = 100.*keV; |
---|
[819] | 95 | ecut = DBL_MAX; |
---|
| 96 | particle = 0; |
---|
[961] | 97 | currentCouple = 0; |
---|
[1055] | 98 | |
---|
| 99 | // Thomas-Fermi screening radii |
---|
| 100 | // Formfactors from A.V. Butkevich et al., NIM A 488 (2002) 282 |
---|
| 101 | |
---|
| 102 | if(0.0 == ScreenRSquare[0]) { |
---|
| 103 | G4double a0 = electron_mass_c2/0.88534; |
---|
| 104 | G4double constn = 6.937e-6/(MeV*MeV); |
---|
| 105 | |
---|
| 106 | ScreenRSquare[0] = alpha2*a0*a0; |
---|
| 107 | for(G4int j=1; j<100; j++) { |
---|
| 108 | G4double x = a0*fNistManager->GetZ13(j); |
---|
| 109 | ScreenRSquare[j] = alpha2*x*x; |
---|
| 110 | x = fNistManager->GetA27(j); |
---|
| 111 | FormFactor[j] = constn*x*x; |
---|
| 112 | } |
---|
| 113 | } |
---|
[819] | 114 | } |
---|
| 115 | |
---|
| 116 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 117 | |
---|
| 118 | G4eCoulombScatteringModel::~G4eCoulombScatteringModel() |
---|
[961] | 119 | {} |
---|
[819] | 120 | |
---|
| 121 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 122 | |
---|
| 123 | void G4eCoulombScatteringModel::Initialise(const G4ParticleDefinition* p, |
---|
[961] | 124 | const G4DataVector& cuts) |
---|
[819] | 125 | { |
---|
[961] | 126 | SetupParticle(p); |
---|
| 127 | currentCouple = 0; |
---|
| 128 | elecXSection = nucXSection = 0.0; |
---|
| 129 | tkin = targetZ = mom2 = DBL_MIN; |
---|
| 130 | ecut = etag = DBL_MAX; |
---|
| 131 | cosThetaMin = cos(PolarAngleLimit()); |
---|
| 132 | currentCuts = &cuts; |
---|
| 133 | //G4cout << "!!! G4eCoulombScatteringModel::Initialise for " |
---|
| 134 | // << p->GetParticleName() << " cos(TetMin)= " << cosThetaMin |
---|
| 135 | // << " cos(TetMax)= " << cosThetaMax <<G4endl; |
---|
[819] | 136 | if(!isInitialised) { |
---|
| 137 | isInitialised = true; |
---|
[1055] | 138 | fParticleChange = GetParticleChangeForGamma(); |
---|
[819] | 139 | } |
---|
[961] | 140 | if(mass < GeV && particle->GetParticleType() != "nucleus") { |
---|
| 141 | InitialiseElementSelectors(p,cuts); |
---|
| 142 | } |
---|
| 143 | } |
---|
[819] | 144 | |
---|
[961] | 145 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
[819] | 146 | |
---|
[961] | 147 | void G4eCoulombScatteringModel::ComputeMaxElectronScattering(G4double cutEnergy) |
---|
| 148 | { |
---|
| 149 | ecut = cutEnergy; |
---|
| 150 | G4double tmax = tkin; |
---|
| 151 | cosTetMaxElec = 1.0; |
---|
| 152 | if(mass > MeV) { |
---|
| 153 | G4double ratio = electron_mass_c2/mass; |
---|
| 154 | G4double tau = tkin/mass; |
---|
| 155 | tmax = 2.0*electron_mass_c2*tau*(tau + 2.)/ |
---|
| 156 | (1.0 + 2.0*ratio*(tau + 1.0) + ratio*ratio); |
---|
| 157 | cosTetMaxElec = 1.0 - std::min(cutEnergy, tmax)*electron_mass_c2/mom2; |
---|
| 158 | } else { |
---|
| 159 | |
---|
| 160 | if(particle == theElectron) tmax *= 0.5; |
---|
| 161 | G4double t = std::min(cutEnergy, tmax); |
---|
| 162 | G4double mom21 = t*(t + 2.0*electron_mass_c2); |
---|
| 163 | G4double t1 = tkin - t; |
---|
| 164 | //G4cout << "tkin= " << tkin << " t= " << t << " t1= " << t1 << G4endl; |
---|
| 165 | if(t1 > 0.0) { |
---|
| 166 | G4double mom22 = t1*(t1 + 2.0*mass); |
---|
| 167 | G4double ctm = (mom2 + mom22 - mom21)*0.5/sqrt(mom2*mom22); |
---|
| 168 | //G4cout << "ctm= " << ctm << G4endl; |
---|
[1055] | 169 | if(ctm < 1.0) cosTetMaxElec = ctm; |
---|
| 170 | if(ctm < -1.0) cosTetMaxElec = -1.0; |
---|
[961] | 171 | } |
---|
| 172 | } |
---|
[819] | 173 | } |
---|
| 174 | |
---|
| 175 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 176 | |
---|
| 177 | G4double G4eCoulombScatteringModel::ComputeCrossSectionPerAtom( |
---|
| 178 | const G4ParticleDefinition* p, |
---|
| 179 | G4double kinEnergy, |
---|
[961] | 180 | G4double Z, G4double, |
---|
[819] | 181 | G4double cutEnergy, G4double) |
---|
| 182 | { |
---|
| 183 | //G4cout << "### G4eCoulombScatteringModel::ComputeCrossSectionPerAtom for " |
---|
[961] | 184 | // << p->GetParticleName()<<" Z= "<<Z<<" e(MeV)= "<< kinEnergy/MeV << G4endl; |
---|
| 185 | G4double xsec = 0.0; |
---|
| 186 | SetupParticle(p); |
---|
| 187 | G4double ekin = std::max(lowEnergyLimit, kinEnergy); |
---|
| 188 | SetupKinematic(ekin, cutEnergy); |
---|
| 189 | if(cosTetMaxNuc < cosTetMinNuc) { |
---|
| 190 | SetupTarget(Z, ekin); |
---|
| 191 | xsec = CrossSectionPerAtom(); |
---|
| 192 | } |
---|
| 193 | /* |
---|
| 194 | G4cout << "e(MeV)= " << ekin/MeV << "cosTetMinNuc= " << cosTetMinNuc |
---|
| 195 | << " cosTetMaxNuc= " << cosTetMaxNuc |
---|
| 196 | << " cosTetMaxElec= " << cosTetMaxElec |
---|
| 197 | << " screenZ= " << screenZ |
---|
| 198 | << " formfactA= " << formfactA |
---|
| 199 | << " cosTetMaxHad= " << cosTetMaxHad << G4endl; |
---|
| 200 | */ |
---|
| 201 | return xsec; |
---|
[819] | 202 | } |
---|
| 203 | |
---|
| 204 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 205 | |
---|
[961] | 206 | G4double G4eCoulombScatteringModel::CrossSectionPerAtom() |
---|
[819] | 207 | { |
---|
[961] | 208 | // This method needs initialisation before be called |
---|
[1055] | 209 | G4double fac = coeff*targetZ*chargeSquare*kinFactor; |
---|
[819] | 210 | elecXSection = 0.0; |
---|
[961] | 211 | nucXSection = 0.0; |
---|
[819] | 212 | |
---|
[961] | 213 | G4double x = 1.0 - cosTetMinNuc; |
---|
| 214 | G4double x1 = x + screenZ; |
---|
[819] | 215 | |
---|
[961] | 216 | if(cosTetMaxElec2 < cosTetMinNuc) { |
---|
| 217 | elecXSection = fac*(cosTetMinNuc - cosTetMaxElec2)/ |
---|
| 218 | (x1*(1.0 - cosTetMaxElec2 + screenZ)); |
---|
| 219 | nucXSection = elecXSection; |
---|
[819] | 220 | } |
---|
| 221 | |
---|
[961] | 222 | //G4cout << "XS tkin(MeV)= " << tkin<<" xs= " <<nucXSection |
---|
| 223 | // << " costmax= " << cosTetMaxNuc2 |
---|
| 224 | // << " costmin= " << cosTetMinNuc << " Z= " << targetZ <<G4endl; |
---|
| 225 | if(cosTetMaxNuc2 < cosTetMinNuc) { |
---|
| 226 | G4double s = screenZ*formfactA; |
---|
| 227 | G4double z1 = 1.0 - cosTetMaxNuc2 + screenZ; |
---|
[1055] | 228 | G4double s1 = 1.0 - s; |
---|
| 229 | G4double d = s1/formfactA; |
---|
[961] | 230 | //G4cout <<"x1= "<<x1<<" z1= " <<z1<<" s= "<<s << " d= " <<d <<G4endl; |
---|
| 231 | if(d < 0.2*x1) { |
---|
| 232 | G4double x2 = x1*x1; |
---|
| 233 | G4double z2 = z1*z1; |
---|
| 234 | x = (1.0/(x1*x2) - 1.0/(z1*z2) - d*1.5*(1.0/(x2*x2) - 1.0/(z2*z2)))/ |
---|
| 235 | (3.0*formfactA*formfactA); |
---|
| 236 | } else { |
---|
| 237 | G4double x2 = x1 + d; |
---|
| 238 | G4double z2 = z1 + d; |
---|
[1055] | 239 | x = (1.0/x1 - 1.0/z1 + 1.0/x2 - 1.0/z2 - 2.0*log(z1*x2/(z2*x1))/d)/(s1*s1); |
---|
[961] | 240 | } |
---|
| 241 | nucXSection += fac*targetZ*x; |
---|
[819] | 242 | } |
---|
[961] | 243 | //G4cout<<" cross(bn)= "<<nucXSection/barn<<" xsElec(bn)= "<<elecXSection/barn |
---|
| 244 | // << " Asc= " << screenZ << G4endl; |
---|
[819] | 245 | |
---|
[961] | 246 | return nucXSection; |
---|
[819] | 247 | } |
---|
| 248 | |
---|
| 249 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
| 250 | |
---|
| 251 | void G4eCoulombScatteringModel::SampleSecondaries( |
---|
[961] | 252 | std::vector<G4DynamicParticle*>* fvect, |
---|
[819] | 253 | const G4MaterialCutsCouple* couple, |
---|
| 254 | const G4DynamicParticle* dp, |
---|
| 255 | G4double cutEnergy, |
---|
[961] | 256 | G4double) |
---|
[819] | 257 | { |
---|
| 258 | G4double kinEnergy = dp->GetKineticEnergy(); |
---|
[961] | 259 | if(kinEnergy <= DBL_MIN) return; |
---|
| 260 | DefineMaterial(couple); |
---|
| 261 | SetupParticle(dp->GetDefinition()); |
---|
| 262 | G4double ekin = std::max(lowEnergyLimit, kinEnergy); |
---|
| 263 | SetupKinematic(ekin, cutEnergy); |
---|
| 264 | //G4cout << "G4eCoulombScatteringModel::SampleSecondaries e(MeV)= " |
---|
| 265 | // << kinEnergy << " " << particle->GetParticleName() << G4endl; |
---|
| 266 | |
---|
| 267 | // Choose nucleus |
---|
| 268 | currentElement = SelectRandomAtom(couple,particle,ekin,cutEnergy,ekin); |
---|
[819] | 269 | |
---|
[961] | 270 | SetupTarget(currentElement->GetZ(),ekin); |
---|
| 271 | |
---|
| 272 | G4double cost = SampleCosineTheta(); |
---|
| 273 | G4double z1 = 1.0 - cost; |
---|
| 274 | if(z1 < 0.0) return; |
---|
[819] | 275 | |
---|
[961] | 276 | G4double sint = sqrt(z1*(1.0 + cost)); |
---|
| 277 | |
---|
| 278 | //G4cout<<"## Sampled sint= " << sint << " Z= " << targetZ |
---|
| 279 | // << " screenZ= " << screenZ << " cn= " << formfactA << G4endl; |
---|
| 280 | |
---|
| 281 | G4double phi = twopi * G4UniformRand(); |
---|
[819] | 282 | |
---|
[961] | 283 | G4ThreeVector direction = dp->GetMomentumDirection(); |
---|
| 284 | G4ThreeVector newDirection(cos(phi)*sint,sin(phi)*sint,cost); |
---|
| 285 | newDirection.rotateUz(direction); |
---|
| 286 | |
---|
| 287 | fParticleChange->ProposeMomentumDirection(newDirection); |
---|
| 288 | |
---|
| 289 | // recoil sampling assuming a small recoil |
---|
| 290 | // and first order correction to primary 4-momentum |
---|
| 291 | if(lowEnergyLimit < kinEnergy) { |
---|
| 292 | G4int ia = SelectIsotopeNumber(currentElement); |
---|
| 293 | G4double Trec = z1*mom2/(amu_c2*G4double(ia)); |
---|
| 294 | |
---|
[1055] | 295 | if(Trec > recoilThreshold) { |
---|
[961] | 296 | G4ParticleDefinition* ion = theParticleTable->FindIon(iz, ia, 0, iz); |
---|
| 297 | Trec = z1*mom2/ion->GetPDGMass(); |
---|
| 298 | if(Trec < kinEnergy) { |
---|
| 299 | G4ThreeVector dir = (direction - newDirection).unit(); |
---|
| 300 | G4DynamicParticle* newdp = new G4DynamicParticle(ion, dir, Trec); |
---|
| 301 | fvect->push_back(newdp); |
---|
| 302 | fParticleChange->SetProposedKineticEnergy(kinEnergy - Trec); |
---|
| 303 | } |
---|
| 304 | } |
---|
| 305 | } |
---|
| 306 | |
---|
| 307 | return; |
---|
| 308 | } |
---|
| 309 | |
---|
| 310 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 311 | |
---|
| 312 | G4double G4eCoulombScatteringModel::SampleCosineTheta() |
---|
| 313 | { |
---|
| 314 | G4double costm = cosTetMaxNuc2; |
---|
[819] | 315 | G4double formf = formfactA; |
---|
[961] | 316 | G4double prob = 0.0; |
---|
| 317 | G4double xs = CrossSectionPerAtom(); |
---|
| 318 | if(xs > 0.0) prob = elecXSection/xs; |
---|
| 319 | |
---|
| 320 | // scattering off e or A? |
---|
| 321 | if(G4UniformRand() < prob) { |
---|
| 322 | costm = cosTetMaxElec2; |
---|
[819] | 323 | formf = 0.0; |
---|
| 324 | } |
---|
[961] | 325 | |
---|
[819] | 326 | /* |
---|
[961] | 327 | G4cout << "SampleCost: e(MeV)= " << tkin |
---|
[819] | 328 | << " ctmin= " << cosThetaMin |
---|
| 329 | << " ctmaxN= " << cosTetMaxNuc |
---|
| 330 | << " ctmax= " << costm |
---|
[961] | 331 | << " Z= " << targetZ << " A= " << targetA |
---|
[819] | 332 | << G4endl; |
---|
| 333 | */ |
---|
[961] | 334 | if(costm >= cosTetMinNuc) return 2.0; |
---|
[819] | 335 | |
---|
[961] | 336 | G4double x1 = 1. - cosTetMinNuc + screenZ; |
---|
| 337 | G4double x2 = 1. - costm + screenZ; |
---|
| 338 | G4double x3 = cosTetMinNuc - costm; |
---|
| 339 | G4double grej, z1; |
---|
[819] | 340 | do { |
---|
[961] | 341 | z1 = x1*x2/(x1 + G4UniformRand()*x3) - screenZ; |
---|
[819] | 342 | grej = 1.0/(1.0 + formf*z1); |
---|
| 343 | } while ( G4UniformRand() > grej*grej ); |
---|
| 344 | |
---|
[961] | 345 | //G4cout << "z= " << z1 << " cross= " << nucXSection/barn |
---|
| 346 | // << " crossE= " << elecXSection/barn << G4endl; |
---|
[819] | 347 | |
---|
[961] | 348 | return 1.0 - z1; |
---|
[819] | 349 | } |
---|
| 350 | |
---|
| 351 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
| 352 | |
---|
| 353 | |
---|