| 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 | // $Id: G4ScreenedNuclearRecoil.cc,v 1.5 2008/01/14 12:11:39 vnivanch Exp $
|
|---|
| 28 | // GEANT4 tag $Name: $
|
|---|
| 29 | //
|
|---|
| 30 | //
|
|---|
| 31 | // Class Description
|
|---|
| 32 | // Process for screened electromagnetic nuclear elastic scattering;
|
|---|
| 33 | // Physics comes from:
|
|---|
| 34 | // Marcus H. Mendenhall and Robert A. Weller,
|
|---|
| 35 | // "Algorithms for the rapid computation of classical cross
|
|---|
| 36 | // sections for screened Coulomb collisions "
|
|---|
| 37 | // Nuclear Instruments and Methods in Physics Research B58 (1991) 11-17
|
|---|
| 38 | // The only input required is a screening function phi(r/a) which is the ratio
|
|---|
| 39 | // of the actual interatomic potential for two atoms with atomic numbers Z1 and Z2,
|
|---|
| 40 | // to the unscreened potential Z1*Z2*e^2/r where e^2 is elm_coupling in Geant4 units
|
|---|
| 41 | //
|
|---|
| 42 | // First version, April 2004, Marcus H. Mendenhall, Vanderbilt University
|
|---|
| 43 | //
|
|---|
| 44 | // 5 May, 2004, Marcus Mendenhall
|
|---|
| 45 | // Added an option for enhancing hard collisions statistically, to allow
|
|---|
| 46 | // backscattering calculations to be carried out with much improved event rates,
|
|---|
| 47 | // without distorting the multiple-scattering broadening too much.
|
|---|
| 48 | // the method SetCrossSectionHardening(G4double fraction, G4double HardeningFactor)
|
|---|
| 49 | // sets what fraction of the events will be randomly hardened,
|
|---|
| 50 | // and the factor by which the impact area is reduced for such selected events.
|
|---|
| 51 | //
|
|---|
| 52 | // 21 November, 2004, Marcus Mendenhall
|
|---|
| 53 | // added static_nucleus to IsApplicable
|
|---|
| 54 | //
|
|---|
| 55 | // 7 December, 2004, Marcus Mendenhall
|
|---|
| 56 | // changed mean free path of stopping particle from 0.0 to 1.0*nanometer
|
|---|
| 57 | // to avoid new verbose warning about 0 MFP in 4.6.2p02
|
|---|
| 58 | //
|
|---|
| 59 | // 17 December, 2004, Marcus Mendenhall
|
|---|
| 60 | // added code to permit screening out overly close collisions which are
|
|---|
| 61 | // expected to be hadronic, not Coulombic
|
|---|
| 62 | //
|
|---|
| 63 | // 19 December, 2004, Marcus Mendenhall
|
|---|
| 64 | // massive rewrite to add modular physics stages and plug-in cross section table
|
|---|
| 65 | // computation. This allows one to select (e.g.) between the normal external python
|
|---|
| 66 | // process and an embedded python interpreter (which is much faster) for generating
|
|---|
| 67 | // the tables.
|
|---|
| 68 | // It also allows one to switch between sub-sampled scattering (event biasing) and
|
|---|
| 69 | // normal scattering, and between non-relativistic kinematics and relativistic
|
|---|
| 70 | // kinematic approximations, without having a class for every combination. Further, one can
|
|---|
| 71 | // add extra stages to the scattering, which can implement various book-keeping processes.
|
|---|
| 72 | //
|
|---|
| 73 | // January 2007, Marcus Mendenhall
|
|---|
| 74 | // Reorganized heavily for inclusion in Geant4 Core. All modules merged into
|
|---|
| 75 | // one source and header, all historic code removed.
|
|---|
| 76 | //
|
|---|
| 77 | // Class Description - End
|
|---|
| 78 |
|
|---|
| 79 |
|
|---|
| 80 | #include <stdio.h>
|
|---|
| 81 |
|
|---|
| 82 | #include "globals.hh"
|
|---|
| 83 |
|
|---|
| 84 | #include "G4ScreenedNuclearRecoil.hh"
|
|---|
| 85 |
|
|---|
| 86 | #include "G4ParticleTypes.hh"
|
|---|
| 87 | #include "G4ParticleTable.hh"
|
|---|
| 88 | #include "G4VParticleChange.hh"
|
|---|
| 89 | #include "G4ParticleChange.hh"
|
|---|
| 90 | #include "G4DataVector.hh"
|
|---|
| 91 | #include "G4Track.hh"
|
|---|
| 92 | #include "G4Step.hh"
|
|---|
| 93 |
|
|---|
| 94 | #include "G4Material.hh"
|
|---|
| 95 | #include "G4Element.hh"
|
|---|
| 96 | #include "G4Isotope.hh"
|
|---|
| 97 | #include "G4MaterialCutsCouple.hh"
|
|---|
| 98 | #include "G4ElementVector.hh"
|
|---|
| 99 | #include "G4IsotopeVector.hh"
|
|---|
| 100 |
|
|---|
| 101 | #include "G4RangeTest.hh"
|
|---|
| 102 | #include "G4ParticleDefinition.hh"
|
|---|
| 103 | #include "G4DynamicParticle.hh"
|
|---|
| 104 | #include "G4ProcessManager.hh"
|
|---|
| 105 | #include "G4StableIsotopes.hh"
|
|---|
| 106 |
|
|---|
| 107 | #include "Randomize.hh"
|
|---|
| 108 |
|
|---|
| 109 | #include "CLHEP/Units/PhysicalConstants.h"
|
|---|
| 110 |
|
|---|
| 111 | #include <iostream>
|
|---|
| 112 | #include <iomanip>
|
|---|
| 113 |
|
|---|
| 114 | G4ScreenedCoulombCrossSection::~G4ScreenedCoulombCrossSection()
|
|---|
| 115 | {
|
|---|
| 116 | ScreeningMap::iterator tables=screeningData.begin();
|
|---|
| 117 | for (;tables != screeningData.end(); tables++) {
|
|---|
| 118 | delete (*tables).second.EMphiData;
|
|---|
| 119 | }
|
|---|
| 120 | screeningData.clear();
|
|---|
| 121 |
|
|---|
| 122 | std::map<G4int, c2_function<G4double> *>::iterator mfpit=MFPTables.begin();
|
|---|
| 123 | for (;mfpit != MFPTables.end(); mfpit++) {
|
|---|
| 124 | delete (*mfpit).second;
|
|---|
| 125 | }
|
|---|
| 126 | MFPTables.clear();
|
|---|
| 127 |
|
|---|
| 128 | }
|
|---|
| 129 |
|
|---|
| 130 | const G4double G4ScreenedCoulombCrossSection::massmap[nMassMapElements+1]={
|
|---|
| 131 | 0, 1.007940, 4.002602, 6.941000, 9.012182, 10.811000, 12.010700,
|
|---|
| 132 | 14.006700, 15.999400, 18.998403, 20.179700, 22.989770, 24.305000, 26.981538, 28.085500,
|
|---|
| 133 | 30.973761, 32.065000, 35.453000, 39.948000, 39.098300, 40.078000, 44.955910, 47.867000,
|
|---|
| 134 | 50.941500, 51.996100, 54.938049, 55.845000, 58.933200, 58.693400, 63.546000, 65.409000,
|
|---|
| 135 | 69.723000, 72.640000, 74.921600, 78.960000, 79.904000, 83.798000, 85.467800, 87.620000,
|
|---|
| 136 | 88.905850, 91.224000, 92.906380, 95.940000, 98.000000, 101.070000, 102.905500, 106.420000,
|
|---|
| 137 | 107.868200, 112.411000, 114.818000, 118.710000, 121.760000, 127.600000, 126.904470, 131.293000,
|
|---|
| 138 | 132.905450, 137.327000, 138.905500, 140.116000, 140.907650, 144.240000, 145.000000, 150.360000,
|
|---|
| 139 | 151.964000, 157.250000, 158.925340, 162.500000, 164.930320, 167.259000, 168.934210, 173.040000,
|
|---|
| 140 | 174.967000, 178.490000, 180.947900, 183.840000, 186.207000, 190.230000, 192.217000, 195.078000,
|
|---|
| 141 | 196.966550, 200.590000, 204.383300, 207.200000, 208.980380, 209.000000, 210.000000, 222.000000,
|
|---|
| 142 | 223.000000, 226.000000, 227.000000, 232.038100, 231.035880, 238.028910, 237.000000, 244.000000,
|
|---|
| 143 | 243.000000, 247.000000, 247.000000, 251.000000, 252.000000, 257.000000, 258.000000, 259.000000,
|
|---|
| 144 | 262.000000, 261.000000, 262.000000, 266.000000, 264.000000, 277.000000, 268.000000, 281.000000,
|
|---|
| 145 | 272.000000, 285.000000, 282.500000, 289.000000, 287.500000, 292.000000};
|
|---|
| 146 |
|
|---|
| 147 | G4ParticleDefinition* G4ScreenedCoulombCrossSection::SelectRandomUnweightedTarget(const G4MaterialCutsCouple* couple)
|
|---|
| 148 | {
|
|---|
| 149 | // Select randomly an element within the material, according to number density only
|
|---|
| 150 | const G4Material* material = couple->GetMaterial();
|
|---|
| 151 | G4int nMatElements = material->GetNumberOfElements();
|
|---|
| 152 | const G4ElementVector* elementVector = material->GetElementVector();
|
|---|
| 153 | const G4Element *element=0;
|
|---|
| 154 | G4ParticleDefinition*target=0;
|
|---|
| 155 |
|
|---|
| 156 | // Special case: the material consists of one element
|
|---|
| 157 | if (nMatElements == 1)
|
|---|
| 158 | {
|
|---|
| 159 | element= (*elementVector)[0];
|
|---|
| 160 | }
|
|---|
| 161 | else
|
|---|
| 162 | {
|
|---|
| 163 | // Composite material
|
|---|
| 164 | G4double random = G4UniformRand() * material->GetTotNbOfAtomsPerVolume();
|
|---|
| 165 | G4double nsum=0.0;
|
|---|
| 166 | const G4double *atomDensities=material->GetVecNbOfAtomsPerVolume();
|
|---|
| 167 |
|
|---|
| 168 | for (G4int k=0 ; k < nMatElements ; k++ )
|
|---|
| 169 | {
|
|---|
| 170 | nsum+=atomDensities[k];
|
|---|
| 171 | element= (*elementVector)[k];
|
|---|
| 172 | if (nsum >= random) break;
|
|---|
| 173 | }
|
|---|
| 174 | }
|
|---|
| 175 |
|
|---|
| 176 | G4int N=0;
|
|---|
| 177 | G4int Z=(G4int)std::floor(element->GetZ()+0.5);
|
|---|
| 178 |
|
|---|
| 179 | G4int nIsotopes=element->GetNumberOfIsotopes();
|
|---|
| 180 | if(!nIsotopes) {
|
|---|
| 181 | if(Z<=92) {
|
|---|
| 182 | // we have no detailed material isotopic info available,
|
|---|
| 183 | // so use G4StableIsotopes table up to Z=92
|
|---|
| 184 | static G4StableIsotopes theIso; // get a stable isotope table for default results
|
|---|
| 185 | nIsotopes=theIso.GetNumberOfIsotopes(Z);
|
|---|
| 186 | G4double random = 100.0*G4UniformRand(); // values are expressed as percent, sum is 100
|
|---|
| 187 | G4int tablestart=theIso.GetFirstIsotope(Z);
|
|---|
| 188 | G4double asum=0.0;
|
|---|
| 189 | for(G4int i=0; i<nIsotopes; i++) {
|
|---|
| 190 | asum+=theIso.GetAbundance(i+tablestart);
|
|---|
| 191 | N=theIso.GetIsotopeNucleonCount(i+tablestart);
|
|---|
| 192 | if(asum >= random) break;
|
|---|
| 193 | }
|
|---|
| 194 | } else {
|
|---|
| 195 | // too heavy for stable isotope table, just use mean mass
|
|---|
| 196 | N=(G4int)std::floor(element->GetN()+0.5);
|
|---|
| 197 | }
|
|---|
| 198 | } else {
|
|---|
| 199 | G4int i;
|
|---|
| 200 | const G4IsotopeVector *isoV=element->GetIsotopeVector();
|
|---|
| 201 | G4double random = G4UniformRand();
|
|---|
| 202 | G4double *abundance=element->GetRelativeAbundanceVector();
|
|---|
| 203 | G4double asum=0.0;
|
|---|
| 204 | for(i=0; i<nIsotopes; i++) {
|
|---|
| 205 | asum+=abundance[i];
|
|---|
| 206 | N=(*isoV)[i]->GetN();
|
|---|
| 207 | if(asum >= random) break;
|
|---|
| 208 | }
|
|---|
| 209 | }
|
|---|
| 210 |
|
|---|
| 211 | // get the official definition of this nucleus, to get the correct value of A
|
|---|
| 212 | // note that GetIon is very slow, so we will cache ones we have already found ourselves.
|
|---|
| 213 | ParticleCache::iterator p=targetMap.find(Z*1000+N);
|
|---|
| 214 | if (p != targetMap.end()) {
|
|---|
| 215 | target=(*p).second;
|
|---|
| 216 | } else{
|
|---|
| 217 | target=G4ParticleTable::GetParticleTable()->GetIon(Z, N, 0.0);
|
|---|
| 218 | targetMap[Z*1000+N]=target;
|
|---|
| 219 | }
|
|---|
| 220 | return target;
|
|---|
| 221 | }
|
|---|
| 222 |
|
|---|
| 223 | void G4ScreenedCoulombCrossSection::BuildMFPTables()
|
|---|
| 224 | {
|
|---|
| 225 | const G4int nmfpvals=200;
|
|---|
| 226 |
|
|---|
| 227 | std::vector<G4double> evals(nmfpvals), mfpvals(nmfpvals);
|
|---|
| 228 |
|
|---|
| 229 | // sum up inverse MFPs per element for each material
|
|---|
| 230 | const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
|
|---|
| 231 | if (materialTable == 0)
|
|---|
| 232 | G4Exception("G4ScreenedCoulombCrossSection::BuildMFPTables - no MaterialTable found)");
|
|---|
| 233 |
|
|---|
| 234 | G4int nMaterials = G4Material::GetNumberOfMaterials();
|
|---|
| 235 |
|
|---|
| 236 | for (G4int matidx=0; matidx < nMaterials; matidx++) {
|
|---|
| 237 |
|
|---|
| 238 | const G4Material* material= (*materialTable)[matidx];
|
|---|
| 239 | const G4ElementVector &elementVector = *(material->GetElementVector());
|
|---|
| 240 | const G4int nMatElements = material->GetNumberOfElements();
|
|---|
| 241 |
|
|---|
| 242 | const G4Element *element=0;
|
|---|
| 243 | const G4double *atomDensities=material->GetVecNbOfAtomsPerVolume();
|
|---|
| 244 |
|
|---|
| 245 | G4double emin=0, emax=0; // find innermost range of cross section functions
|
|---|
| 246 | for (G4int kel=0 ; kel < nMatElements ; kel++ )
|
|---|
| 247 | {
|
|---|
| 248 | element=elementVector[kel];
|
|---|
| 249 | G4int Z=(G4int)std::floor(element->GetZ()+0.5);
|
|---|
| 250 | c2_function<G4double> &ifunc=*sigmaMap[Z];
|
|---|
| 251 | if(!kel || ifunc.xmin() > emin) emin=ifunc.xmin();
|
|---|
| 252 | if(!kel || ifunc.xmax() < emax) emax=ifunc.xmax();
|
|---|
| 253 | }
|
|---|
| 254 |
|
|---|
| 255 | G4double logint=std::log(emax/emin) / (nmfpvals-1) ; // logarithmic increment for tables
|
|---|
| 256 |
|
|---|
| 257 | // compute energy scale for interpolator. Force exact values at both ends to avoid range errors
|
|---|
| 258 | for (G4int i=1; i<nmfpvals-1; i++) evals[i]=emin*std::exp(logint*i);
|
|---|
| 259 | evals.front()=emin;
|
|---|
| 260 | evals.back()=emax;
|
|---|
| 261 |
|
|---|
| 262 | // zero out the inverse mfp sums to start
|
|---|
| 263 | for (G4int eidx=0; eidx < nmfpvals; eidx++) mfpvals[eidx] = 0.0;
|
|---|
| 264 |
|
|---|
| 265 | // sum inverse mfp for each element in this material and for each energy
|
|---|
| 266 | for (G4int kel=0 ; kel < nMatElements ; kel++ )
|
|---|
| 267 | {
|
|---|
| 268 | element=elementVector[kel];
|
|---|
| 269 | G4int Z=(G4int)std::floor(element->GetZ()+0.5);
|
|---|
| 270 | c2_function<G4double> &sigma=*sigmaMap[Z];
|
|---|
| 271 | G4double ndens = atomDensities[kel]; // compute atom fraction for this element in this material
|
|---|
| 272 |
|
|---|
| 273 | for (G4int eidx=0; eidx < nmfpvals; eidx++) {
|
|---|
| 274 | mfpvals[eidx] += ndens*sigma(evals[eidx]);
|
|---|
| 275 | }
|
|---|
| 276 | }
|
|---|
| 277 |
|
|---|
| 278 | // convert inverse mfp to regular mfp
|
|---|
| 279 | for (G4int eidx=0; eidx < nmfpvals; eidx++) {
|
|---|
| 280 | mfpvals[eidx] = 1.0/mfpvals[eidx];
|
|---|
| 281 | }
|
|---|
| 282 | // and make a new interpolating function out of the sum
|
|---|
| 283 | MFPTables[matidx] = static_cast<c2_function<G4double> *>(new log_log_interpolating_function<G4double>(
|
|---|
| 284 | evals, mfpvals));
|
|---|
| 285 | }
|
|---|
| 286 |
|
|---|
| 287 | #ifdef DEBUG
|
|---|
| 288 | for (G4int matidx=0; matidx < nMaterials; matidx++) {
|
|---|
| 289 | const G4Material* material= (*materialTable)[matidx];
|
|---|
| 290 | G4cout << "***** MFP (1MeV) ***** " << material->GetName() << " " << (*MFPTables[matidx])(1.0) << G4endl;
|
|---|
| 291 | }
|
|---|
| 292 | #endif
|
|---|
| 293 |
|
|---|
| 294 | }
|
|---|
| 295 |
|
|---|
| 296 | G4ScreenedNuclearRecoil::
|
|---|
| 297 | G4ScreenedNuclearRecoil(const G4String& processName,
|
|---|
| 298 | const G4String &ScreeningKey,
|
|---|
| 299 | G4bool GenerateRecoils,
|
|---|
| 300 | G4double RecoilCutoff, G4double PhysicsCutoff) :
|
|---|
| 301 | G4VDiscreteProcess(processName),
|
|---|
| 302 | screeningKey(ScreeningKey),
|
|---|
| 303 | generateRecoils(GenerateRecoils), avoidReactions(1),
|
|---|
| 304 | recoilCutoff(RecoilCutoff), physicsCutoff(PhysicsCutoff),
|
|---|
| 305 | hardeningFraction(0.0), hardeningFactor(1.0),
|
|---|
| 306 | externalCrossSectionConstructor(0)
|
|---|
| 307 | {
|
|---|
| 308 | highEnergyLimit=100.0*MeV;
|
|---|
| 309 | lowEnergyLimit=physicsCutoff;
|
|---|
| 310 | registerDepositedEnergy=1; // by default, don't hide NIEL
|
|---|
| 311 | MFPScale=1.0;
|
|---|
| 312 | // SetVerboseLevel(2);
|
|---|
| 313 | AddStage(new G4ScreenedCoulombClassicalKinematics);
|
|---|
| 314 | AddStage(new G4SingleScatter);
|
|---|
| 315 | }
|
|---|
| 316 |
|
|---|
| 317 | void G4ScreenedNuclearRecoil::ResetTables()
|
|---|
| 318 | {
|
|---|
| 319 | std::map<G4int, c2_function<G4double>*>::iterator xh=meanFreePathTables.begin();
|
|---|
| 320 | for(;xh != meanFreePathTables.end(); xh++) {
|
|---|
| 321 | delete (*xh).second;
|
|---|
| 322 | }
|
|---|
| 323 | meanFreePathTables.clear();
|
|---|
| 324 |
|
|---|
| 325 | std::map<G4int, G4ScreenedCoulombCrossSection*>::iterator xt=crossSectionHandlers.begin();
|
|---|
| 326 | for(;xt != crossSectionHandlers.end(); xt++) {
|
|---|
| 327 | delete (*xt).second;
|
|---|
| 328 | }
|
|---|
| 329 | crossSectionHandlers.clear();
|
|---|
| 330 | }
|
|---|
| 331 |
|
|---|
| 332 | void G4ScreenedNuclearRecoil::ClearStages()
|
|---|
| 333 | {
|
|---|
| 334 | // I don't think I like deleting the processes here... they are better abandoned
|
|---|
| 335 | // if the creator doesn't get rid of them
|
|---|
| 336 | // std::vector<G4ScreenedCollisionStage *>::iterator stage=collisionStages.begin();
|
|---|
| 337 | //for(; stage != collisionStages.end(); stage++) delete (*stage);
|
|---|
| 338 |
|
|---|
| 339 | collisionStages.clear();
|
|---|
| 340 | }
|
|---|
| 341 |
|
|---|
| 342 | G4ScreenedNuclearRecoil::~G4ScreenedNuclearRecoil()
|
|---|
| 343 | {
|
|---|
| 344 | ResetTables();
|
|---|
| 345 | }
|
|---|
| 346 |
|
|---|
| 347 | // returns true if it appears the nuclei collided, and we are interested in checking
|
|---|
| 348 | G4bool G4ScreenedNuclearRecoil::CheckNuclearCollision(
|
|---|
| 349 | G4double A, G4double a1, G4double apsis) {
|
|---|
| 350 | return avoidReactions && (apsis < (1.1*(std::pow(A,1.0/3.0)+std::pow(a1,1.0/3.0)) + 1.4)*fermi);
|
|---|
| 351 | // nuclei are within 1.4 fm (reduced pion Compton wavelength) of each other at apsis,
|
|---|
| 352 | // this is hadronic, skip it
|
|---|
| 353 | }
|
|---|
| 354 |
|
|---|
| 355 | G4ScreenedCoulombCrossSection *G4ScreenedNuclearRecoil::GetNewCrossSectionHandler(void) {
|
|---|
| 356 | G4ScreenedCoulombCrossSection *xc;
|
|---|
| 357 | if(!externalCrossSectionConstructor) xc=new G4NativeScreenedCoulombCrossSection;
|
|---|
| 358 | else xc=externalCrossSectionConstructor->create();
|
|---|
| 359 | xc->SetVerbosity(verboseLevel);
|
|---|
| 360 | return xc;
|
|---|
| 361 | }
|
|---|
| 362 |
|
|---|
| 363 | G4double G4ScreenedNuclearRecoil::GetMeanFreePath(const G4Track& track,
|
|---|
| 364 | G4double,
|
|---|
| 365 | G4ForceCondition*)
|
|---|
| 366 | {
|
|---|
| 367 | const G4DynamicParticle* incoming = track.GetDynamicParticle();
|
|---|
| 368 | G4double energy = incoming->GetKineticEnergy();
|
|---|
| 369 | G4double a1=incoming->GetDefinition()->GetPDGMass()/amu_c2;
|
|---|
| 370 |
|
|---|
| 371 | G4double meanFreePath;
|
|---|
| 372 |
|
|---|
| 373 | if (energy < lowEnergyLimit || energy < recoilCutoff) return 1.0*nanometer; /* stop slow particles! */
|
|---|
| 374 | else if (energy > highEnergyLimit*a1) energy=highEnergyLimit*a1; /* constant MFP at high energy */
|
|---|
| 375 |
|
|---|
| 376 | G4double fz1=incoming->GetDefinition()->GetPDGCharge();
|
|---|
| 377 | G4int z1=(G4int)(fz1/eplus + 0.5);
|
|---|
| 378 |
|
|---|
| 379 | std::map<G4int, G4ScreenedCoulombCrossSection*>::iterator xh=
|
|---|
| 380 | crossSectionHandlers.find(z1);
|
|---|
| 381 | G4ScreenedCoulombCrossSection *xs;
|
|---|
| 382 |
|
|---|
| 383 | if (xh==crossSectionHandlers.end()) {
|
|---|
| 384 | xs =crossSectionHandlers[z1]=GetNewCrossSectionHandler();
|
|---|
| 385 | xs->LoadData(screeningKey, z1, a1, physicsCutoff);
|
|---|
| 386 | xs->BuildMFPTables();
|
|---|
| 387 | } else xs=(*xh).second;
|
|---|
| 388 |
|
|---|
| 389 | const G4MaterialCutsCouple* materialCouple = track.GetMaterialCutsCouple();
|
|---|
| 390 | size_t materialIndex = materialCouple->GetMaterial()->GetIndex();
|
|---|
| 391 |
|
|---|
| 392 | c2_function<G4double> &mfp=*(*xs)[materialIndex];
|
|---|
| 393 |
|
|---|
| 394 | // make absolutely certain we don't get an out-of-range energy
|
|---|
| 395 | meanFreePath = mfp(std::min(std::max(energy, mfp.xmin()), mfp.xmax()));
|
|---|
| 396 |
|
|---|
| 397 | // G4cout << "MFP: " << meanFreePath << " index " << materialIndex << " energy " << energy << " MFPScale " << MFPScale << G4endl;
|
|---|
| 398 |
|
|---|
| 399 | return meanFreePath*MFPScale;
|
|---|
| 400 | }
|
|---|
| 401 |
|
|---|
| 402 | G4VParticleChange* G4ScreenedNuclearRecoil::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
|
|---|
| 403 | {
|
|---|
| 404 | validCollision=1;
|
|---|
| 405 | aParticleChange.Initialize(aTrack);
|
|---|
| 406 | NIEL=0.0; // default is no NIEL deposited
|
|---|
| 407 |
|
|---|
| 408 | // do universal setup
|
|---|
| 409 |
|
|---|
| 410 | const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
|
|---|
| 411 | G4ParticleDefinition *baseParticle=aTrack.GetDefinition();
|
|---|
| 412 |
|
|---|
| 413 | G4double fz1=baseParticle->GetPDGCharge()/eplus;
|
|---|
| 414 | G4int z1=(G4int)(fz1+0.5);
|
|---|
| 415 | G4double incidentEnergy = incidentParticle->GetKineticEnergy();
|
|---|
| 416 |
|
|---|
| 417 | // Select randomly one element and (possibly) isotope in the current material.
|
|---|
| 418 | const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
|
|---|
| 419 |
|
|---|
| 420 | if(incidentEnergy < GetRecoilCutoff()) { // check energy sanity on entry
|
|---|
| 421 | if(!baseParticle->GetProcessManager()->
|
|---|
| 422 | GetAtRestProcessVector()->size())
|
|---|
| 423 | aParticleChange.ProposeTrackStatus(fStopAndKill);
|
|---|
| 424 | else
|
|---|
| 425 | aParticleChange.ProposeTrackStatus(fStopButAlive);
|
|---|
| 426 |
|
|---|
| 427 | AddToNIEL(incidentEnergy);
|
|---|
| 428 | aParticleChange.ProposeEnergy(0.0);
|
|---|
| 429 | // stop the particle and bail out
|
|---|
| 430 | validCollision=0;
|
|---|
| 431 | }
|
|---|
| 432 |
|
|---|
| 433 | const G4Material* mat = couple->GetMaterial();
|
|---|
| 434 | G4double numberDensity=mat->GetTotNbOfAtomsPerVolume();
|
|---|
| 435 | G4double lattice=0.5/std::pow(numberDensity,1.0/3.0); // typical lattice half-spacing
|
|---|
| 436 | G4double length=GetCurrentInteractionLength();
|
|---|
| 437 | G4double sigopi=1.0/(CLHEP::pi*numberDensity*length); // this is sigma0/pi
|
|---|
| 438 |
|
|---|
| 439 | // compute the impact parameter very early, so if is rejected as too far away, little effort is wasted
|
|---|
| 440 | // this is the TRIM method for determining an impact parameter based on the flight path
|
|---|
| 441 | // this gives a cumulative distribution of N(P)= 1-exp(-pi P^2 n l)
|
|---|
| 442 | // which says the probability of NOT hitting a disk of area sigma= pi P^2 =exp(-sigma N l)
|
|---|
| 443 | // which may be reasonable
|
|---|
| 444 | G4double P;
|
|---|
| 445 | if(sigopi < lattice*lattice) {
|
|---|
| 446 | // normal long-flight approximation
|
|---|
| 447 | P = std::sqrt(-std::log(G4UniformRand()) *sigopi);
|
|---|
| 448 | } else {
|
|---|
| 449 | // short-flight limit
|
|---|
| 450 | P = std::sqrt(G4UniformRand())*lattice;
|
|---|
| 451 | }
|
|---|
| 452 |
|
|---|
| 453 | G4double fraction=GetHardeningFraction();
|
|---|
| 454 | if(fraction && G4UniformRand() < fraction) {
|
|---|
| 455 | // pick out some events, and increase the central cross section
|
|---|
| 456 | // by reducing the impact parameter
|
|---|
| 457 | P /= std::sqrt(GetHardeningFactor());
|
|---|
| 458 | }
|
|---|
| 459 |
|
|---|
| 460 |
|
|---|
| 461 | // check if we are far enough away that the energy transfer must be below cutoff,
|
|---|
| 462 | // and leave everything alone if so, saving a lot of time.
|
|---|
| 463 | if(P*P > sigopi) {
|
|---|
| 464 | if(GetVerboseLevel() > 1)
|
|---|
| 465 | printf("ScreenedNuclear impact reject: length=%.3f P=%.4f limit=%.4f\n",
|
|---|
| 466 | length/angstrom, P/angstrom,std::sqrt(sigopi)/angstrom);
|
|---|
| 467 | // no collision, don't follow up with anything
|
|---|
| 468 | validCollision=0;
|
|---|
| 469 | }
|
|---|
| 470 |
|
|---|
| 471 | // find out what we hit, and record it in our kinematics block.
|
|---|
| 472 | if(validCollision) {
|
|---|
| 473 | G4ScreenedCoulombCrossSection *xsect=GetCrossSectionHandlers()[z1];
|
|---|
| 474 | G4ParticleDefinition *recoilIon=
|
|---|
| 475 | xsect->SelectRandomUnweightedTarget(couple);
|
|---|
| 476 | kinematics.crossSection=xsect;
|
|---|
| 477 | kinematics.recoilIon=recoilIon;
|
|---|
| 478 | kinematics.impactParameter=P;
|
|---|
| 479 | kinematics.a1=baseParticle->GetPDGMass()/amu_c2;
|
|---|
| 480 | kinematics.a2=recoilIon->GetPDGMass()/amu_c2;
|
|---|
| 481 | } else {
|
|---|
| 482 | kinematics.recoilIon=0;
|
|---|
| 483 | kinematics.impactParameter=0;
|
|---|
| 484 | kinematics.a1=baseParticle->GetPDGMass()/amu_c2;
|
|---|
| 485 | kinematics.a2=0;
|
|---|
| 486 | }
|
|---|
| 487 |
|
|---|
| 488 | std::vector<G4ScreenedCollisionStage *>::iterator stage=collisionStages.begin();
|
|---|
| 489 |
|
|---|
| 490 | for(; stage != collisionStages.end(); stage++)
|
|---|
| 491 | (*stage)->DoCollisionStep(this,aTrack, aStep);
|
|---|
| 492 |
|
|---|
| 493 | if(registerDepositedEnergy) {
|
|---|
| 494 | aParticleChange.ProposeLocalEnergyDeposit(NIEL);
|
|---|
| 495 | aParticleChange.ProposeNonIonizingEnergyDeposit(NIEL);
|
|---|
| 496 | }
|
|---|
| 497 | return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
|
|---|
| 498 | }
|
|---|
| 499 |
|
|---|
| 500 | G4bool G4ScreenedCoulombClassicalKinematics::DoScreeningComputation(G4ScreenedNuclearRecoil *master,
|
|---|
| 501 | const G4ScreeningTables *screen, G4double eps, G4double beta)
|
|---|
| 502 | {
|
|---|
| 503 | G4double au=screen->au;
|
|---|
| 504 | G4CoulombKinematicsInfo &kin=master->GetKinematics();
|
|---|
| 505 | G4double A=kin.a2;
|
|---|
| 506 | G4double a1=kin.a1;
|
|---|
| 507 |
|
|---|
| 508 | G4double xx0; // first estimate of closest approach
|
|---|
| 509 | if(eps < 5.0) {
|
|---|
| 510 | G4double y=std::log(eps);
|
|---|
| 511 | G4double mlrho4=((((3.517e-4*y+1.401e-2)*y+2.393e-1)*y+2.734)*y+2.220);
|
|---|
| 512 | G4double rho4=std::exp(-mlrho4); // W&M eq. 18
|
|---|
| 513 | G4double bb2=0.5*beta*beta;
|
|---|
| 514 | xx0=std::sqrt(bb2+std::sqrt(bb2*bb2+rho4)); // W&M eq. 17
|
|---|
| 515 | } else {
|
|---|
| 516 | G4double ee=1.0/(2.0*eps);
|
|---|
| 517 | xx0=ee+std::sqrt(ee*ee+beta*beta); // W&M eq. 15 (Rutherford value)
|
|---|
| 518 | if(master->CheckNuclearCollision(A, a1, xx0*au)) return 0; // nuclei too close
|
|---|
| 519 |
|
|---|
| 520 | }
|
|---|
| 521 |
|
|---|
| 522 | c2_function<G4double> &phiData=*(screen->EMphiData);
|
|---|
| 523 | // instantiate all the needed functions statically, so no allocation is done at run time
|
|---|
| 524 | // we will be solving x^2 - x phi(x*au)/eps - beta^2 == 0.0
|
|---|
| 525 | // or, for easier scaling, x'^2 - x' au phi(x')/eps - beta^2 au^2
|
|---|
| 526 | static c2_plugin_function<G4double> phifunc;
|
|---|
| 527 | static c2_quadratic<G4double> xsq(0., 0., 0., 1.); // x^2
|
|---|
| 528 | static c2_linear<G4double> xovereps(0., 0., 0.); // will fill this in with the right slope at run time
|
|---|
| 529 | static c2_function<G4double> &xphi=xovereps*phifunc;
|
|---|
| 530 | static c2_function<G4double> &diff=xsq-xphi;
|
|---|
| 531 |
|
|---|
| 532 | xovereps.reset(0., 0.0, au/eps); // slope of x*au/eps term
|
|---|
| 533 | phifunc.set_function(phiData); // install interpolating table
|
|---|
| 534 |
|
|---|
| 535 | G4double xx1, phip, phip2;
|
|---|
| 536 | G4int root_error;
|
|---|
| 537 |
|
|---|
| 538 | xx1=diff.find_root(phiData.xmin(), std::min(10*xx0*au,phiData.xmax()),
|
|---|
| 539 | std::min(xx0*au, phiData.xmax()), beta*beta*au*au, &root_error, &phip, &phip2)/au;
|
|---|
| 540 |
|
|---|
| 541 | if(root_error) {
|
|---|
| 542 | G4cout << "Screened Coulomb Root Finder Error" << G4endl;
|
|---|
| 543 | G4cout << "au " << au << " A " << A << " a1 " << a1 << " xx1 " << xx1 << " eps " << eps << " beta " << beta << G4endl;
|
|---|
| 544 | G4cout << " xmin " << phiData.xmin() << " xmax " << std::min(10*xx0*au,phiData.xmax()) ;
|
|---|
| 545 | G4cout << " f(xmin) " << phifunc(phiData.xmin()) << " f(xmax) " << phifunc(std::min(10*xx0*au,phiData.xmax())) ;
|
|---|
| 546 | G4cout << " xstart " << std::min(xx0*au, phiData.xmax()) << " target " << beta*beta*au*au ;
|
|---|
| 547 | G4cout << G4endl;
|
|---|
| 548 | throw c2_exception("Failed root find");
|
|---|
| 549 | }
|
|---|
| 550 |
|
|---|
| 551 | phifunc.unset_function(); // throws an exception if used without setting again
|
|---|
| 552 | // phiprime is scaled by one factor of au because phi is evaluated at (xx0*au),
|
|---|
| 553 | G4double phiprime=phip*au;
|
|---|
| 554 |
|
|---|
| 555 | //lambda0 is from W&M 19
|
|---|
| 556 | G4double lambda0=1.0/std::sqrt(0.5+beta*beta/(2.0*xx1*xx1)-phiprime/(2.0*eps));
|
|---|
| 557 |
|
|---|
| 558 | //compute the 6-term Lobatto integral alpha (per W&M 21, with different coefficients)
|
|---|
| 559 | // this is probably completely un-needed but gives the highest quality results,
|
|---|
| 560 | G4double alpha=(1.0+ lambda0)/30.0;
|
|---|
| 561 | G4double xvals[]={0.98302349, 0.84652241, 0.53235309, 0.18347974};
|
|---|
| 562 | G4double weights[]={0.03472124, 0.14769029, 0.23485003, 0.18602489};
|
|---|
| 563 | for(G4int k=0; k<4; k++) {
|
|---|
| 564 | G4double x, ff;
|
|---|
| 565 | x=xx1/xvals[k];
|
|---|
| 566 | ff=1.0/std::sqrt(1.0-phiData(x*au)/(x*eps)-beta*beta/(x*x));
|
|---|
| 567 | alpha+=weights[k]*ff;
|
|---|
| 568 | }
|
|---|
| 569 |
|
|---|
| 570 | G4double thetac1=CLHEP::pi*beta*alpha/xx1; // complement of CM scattering angle
|
|---|
| 571 | G4double sintheta=std::sin(thetac1); //note sin(pi-theta)=sin(theta)
|
|---|
| 572 | G4double costheta=-std::cos(thetac1); // note cos(pi-theta)=-cos(theta)
|
|---|
| 573 | // G4double psi=std::atan2(sintheta, costheta+a1/A); // lab scattering angle (M&T 3rd eq. 8.69)
|
|---|
| 574 |
|
|---|
| 575 | // numerics note: because we checked above for reasonable values of beta which give real recoils,
|
|---|
| 576 | // we don't have to look too closely for theta -> 0 here (which would cause sin(theta)
|
|---|
| 577 | // and 1-cos(theta) to both vanish and make the atan2 ill behaved).
|
|---|
| 578 | G4double zeta=std::atan2(sintheta, 1-costheta); // lab recoil angle (M&T 3rd eq. 8.73)
|
|---|
| 579 | G4double coszeta=std::cos(zeta);
|
|---|
| 580 | G4double sinzeta=std::sin(zeta);
|
|---|
| 581 |
|
|---|
| 582 | kin.sinTheta=sintheta;
|
|---|
| 583 | kin.cosTheta=costheta;
|
|---|
| 584 | kin.sinZeta=sinzeta;
|
|---|
| 585 | kin.cosZeta=coszeta;
|
|---|
| 586 | return 1; // all OK, collision is valid
|
|---|
| 587 | }
|
|---|
| 588 |
|
|---|
| 589 | void G4ScreenedCoulombClassicalKinematics::DoCollisionStep(G4ScreenedNuclearRecoil *master,
|
|---|
| 590 | const G4Track& aTrack, const G4Step&) {
|
|---|
| 591 |
|
|---|
| 592 | if(!master->GetValidCollision()) return;
|
|---|
| 593 |
|
|---|
| 594 | G4ParticleChange &aParticleChange=master->GetParticleChange();
|
|---|
| 595 | G4CoulombKinematicsInfo &kin=master->GetKinematics();
|
|---|
| 596 |
|
|---|
| 597 | const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
|
|---|
| 598 | G4ParticleDefinition *baseParticle=aTrack.GetDefinition();
|
|---|
| 599 |
|
|---|
| 600 | G4double incidentEnergy = incidentParticle->GetKineticEnergy();
|
|---|
| 601 |
|
|---|
| 602 | // this adjustment to a1 gives the right results for soft (constant gamma)
|
|---|
| 603 | // relativistic collisions. Hard collisions are wrong anyway, since the
|
|---|
| 604 | // Coulombic and hadronic terms interfere and cannot be added.
|
|---|
| 605 | G4double gamma=(1.0+incidentEnergy/baseParticle->GetPDGMass());
|
|---|
| 606 | G4double a1=kin.a1*gamma; // relativistic gamma correction
|
|---|
| 607 |
|
|---|
| 608 | G4ParticleDefinition *recoilIon=kin.recoilIon;
|
|---|
| 609 | G4double A=recoilIon->GetPDGMass()/amu_c2;
|
|---|
| 610 | G4int Z=(G4int)((recoilIon->GetPDGCharge()/eplus)+0.5);
|
|---|
| 611 |
|
|---|
| 612 | G4double Ec = incidentEnergy*(A/(A+a1)); // energy in CM frame (non-relativistic!)
|
|---|
| 613 | const G4ScreeningTables *screen=kin.crossSection->GetScreening(Z);
|
|---|
| 614 | G4double au=screen->au; // screening length
|
|---|
| 615 |
|
|---|
| 616 | G4double beta = kin.impactParameter/au; // dimensionless impact parameter
|
|---|
| 617 | G4double eps = Ec/(screen->z1*Z*elm_coupling/au); // dimensionless energy
|
|---|
| 618 |
|
|---|
| 619 | G4bool ok=DoScreeningComputation(master, screen, eps, beta);
|
|---|
| 620 | if(!ok) {
|
|---|
| 621 | master->SetValidCollision(0); // flag bad collision
|
|---|
| 622 | return; // just bail out without setting valid flag
|
|---|
| 623 | }
|
|---|
| 624 |
|
|---|
| 625 | G4double eRecoil=4*incidentEnergy*a1*A*kin.cosZeta*kin.cosZeta/((a1+A)*(a1+A));
|
|---|
| 626 | kin.eRecoil=eRecoil;
|
|---|
| 627 |
|
|---|
| 628 | if(incidentEnergy-eRecoil < master->GetRecoilCutoff()) {
|
|---|
| 629 | if(!baseParticle->GetProcessManager()->
|
|---|
| 630 | GetAtRestProcessVector()->size())
|
|---|
| 631 | aParticleChange.ProposeTrackStatus(fStopAndKill);
|
|---|
| 632 | else
|
|---|
| 633 | aParticleChange.ProposeTrackStatus(fStopButAlive);
|
|---|
| 634 | aParticleChange.ProposeEnergy(0.0);
|
|---|
| 635 | master->AddToNIEL(incidentEnergy-eRecoil);
|
|---|
| 636 | }
|
|---|
| 637 |
|
|---|
| 638 | if(master->GetEnableRecoils() && eRecoil > master->GetRecoilCutoff()) {
|
|---|
| 639 | kin.recoilIon=recoilIon;
|
|---|
| 640 | } else {
|
|---|
| 641 | kin.recoilIon=0; // this flags no recoil to be generated
|
|---|
| 642 | master->AddToNIEL(eRecoil) ;
|
|---|
| 643 | }
|
|---|
| 644 | }
|
|---|
| 645 |
|
|---|
| 646 | void G4SingleScatter::DoCollisionStep(G4ScreenedNuclearRecoil *master,
|
|---|
| 647 | const G4Track& aTrack, const G4Step&) {
|
|---|
| 648 |
|
|---|
| 649 | if(!master->GetValidCollision()) return;
|
|---|
| 650 |
|
|---|
| 651 | G4CoulombKinematicsInfo &kin=master->GetKinematics();
|
|---|
| 652 | G4ParticleChange &aParticleChange=master->GetParticleChange();
|
|---|
| 653 |
|
|---|
| 654 | const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
|
|---|
| 655 | G4double incidentEnergy = incidentParticle->GetKineticEnergy();
|
|---|
| 656 | G4double eRecoil=kin.eRecoil;
|
|---|
| 657 |
|
|---|
| 658 | G4double azimuth=G4UniformRand()*(2.0*CLHEP::pi);
|
|---|
| 659 | G4double sa=std::sin(azimuth);
|
|---|
| 660 | G4double ca=std::cos(azimuth);
|
|---|
| 661 |
|
|---|
| 662 | G4ThreeVector recoilMomentumDirection(kin.sinZeta*ca, kin.sinZeta*sa, kin.cosZeta);
|
|---|
| 663 | G4ParticleMomentum incidentDirection = incidentParticle->GetMomentumDirection();
|
|---|
| 664 | recoilMomentumDirection=recoilMomentumDirection.rotateUz(incidentDirection);
|
|---|
| 665 | G4ThreeVector recoilMomentum=recoilMomentumDirection*std::sqrt(2.0*eRecoil*kin.a2*amu_c2);
|
|---|
| 666 |
|
|---|
| 667 | if(aParticleChange.GetEnergy() != 0.0) { // DoKinematics hasn't stopped it!
|
|---|
| 668 | G4ThreeVector beamMomentum=incidentParticle->GetMomentum()-recoilMomentum;
|
|---|
| 669 | aParticleChange.ProposeMomentumDirection(beamMomentum.unit()) ;
|
|---|
| 670 | aParticleChange.ProposeEnergy(incidentEnergy-eRecoil);
|
|---|
| 671 | }
|
|---|
| 672 |
|
|---|
| 673 | if(kin.recoilIon) {
|
|---|
| 674 | G4DynamicParticle* recoil = new G4DynamicParticle (kin.recoilIon,
|
|---|
| 675 | recoilMomentumDirection,eRecoil) ;
|
|---|
| 676 |
|
|---|
| 677 | aParticleChange.SetNumberOfSecondaries(1);
|
|---|
| 678 | aParticleChange.AddSecondary(recoil);
|
|---|
| 679 | }
|
|---|
| 680 | }
|
|---|
| 681 |
|
|---|
| 682 | G4bool G4ScreenedNuclearRecoil::
|
|---|
| 683 | IsApplicable(const G4ParticleDefinition& aParticleType)
|
|---|
| 684 | {
|
|---|
| 685 | return aParticleType == *(G4Proton::Proton()) ||
|
|---|
| 686 | aParticleType.GetParticleType() == "nucleus" ||
|
|---|
| 687 | aParticleType.GetParticleType() == "static_nucleus";
|
|---|
| 688 | }
|
|---|
| 689 |
|
|---|
| 690 |
|
|---|
| 691 | void
|
|---|
| 692 | G4ScreenedNuclearRecoil::
|
|---|
| 693 | DumpPhysicsTable(const G4ParticleDefinition&)
|
|---|
| 694 | {
|
|---|
| 695 | }
|
|---|
| 696 |
|
|---|
| 697 | // This used to be the file mhmScreenedNuclearRecoil_native.cc
|
|---|
| 698 | // it has been included here to collect this file into a smaller number of packages
|
|---|
| 699 |
|
|---|
| 700 | #include "G4DataVector.hh"
|
|---|
| 701 | #include "G4Material.hh"
|
|---|
| 702 | #include "G4Element.hh"
|
|---|
| 703 | #include "G4Isotope.hh"
|
|---|
| 704 | #include "G4MaterialCutsCouple.hh"
|
|---|
| 705 | #include "G4ElementVector.hh"
|
|---|
| 706 | #include <vector>
|
|---|
| 707 |
|
|---|
| 708 | static c2_function<G4double> &ZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
|
|---|
| 709 | {
|
|---|
| 710 | static const size_t ncoef=4;
|
|---|
| 711 | static G4double scales[ncoef]={-3.2, -0.9432, -0.4028, -0.2016};
|
|---|
| 712 | static G4double coefs[ncoef]={0.1818,0.5099,0.2802,0.0281};
|
|---|
| 713 |
|
|---|
| 714 | G4double au=0.8854*angstrom*0.529/(std::pow(z1, 0.23)+std::pow(z2,0.23));
|
|---|
| 715 | std::vector<G4double> r(npoints), phi(npoints);
|
|---|
| 716 |
|
|---|
| 717 | for(size_t i=0; i<npoints; i++) {
|
|---|
| 718 | G4double rr=(float)i/(float)(npoints-1);
|
|---|
| 719 | r[i]=rr*rr*rMax; // use quadratic r scale to make sampling fine near the center
|
|---|
| 720 | G4double sum=0.0;
|
|---|
| 721 | for(size_t j=0; j<ncoef; j++) sum+=coefs[j]*std::exp(scales[j]*r[i]/au);
|
|---|
| 722 | phi[i]=sum;
|
|---|
| 723 | }
|
|---|
| 724 |
|
|---|
| 725 | // compute the derivative at the origin for the spline
|
|---|
| 726 | G4double phiprime0=0.0;
|
|---|
| 727 | for(size_t j=0; j<ncoef; j++) phiprime0+=scales[j]*coefs[j]*std::exp(scales[j]*r[0]/au);
|
|---|
| 728 | phiprime0*=(1.0/au); // put back in natural units;
|
|---|
| 729 |
|
|---|
| 730 | *auval=au;
|
|---|
| 731 | return *static_cast<c2_function<G4double> *>(new lin_log_interpolating_function<G4double>(r, phi, false, phiprime0));
|
|---|
| 732 | }
|
|---|
| 733 |
|
|---|
| 734 | static c2_function<G4double> &MoliereScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
|
|---|
| 735 | {
|
|---|
| 736 | static const size_t ncoef=3;
|
|---|
| 737 | static G4double scales[ncoef]={-6.0, -1.2, -0.3};
|
|---|
| 738 | static G4double coefs[ncoef]={0.10, 0.55, 0.35};
|
|---|
| 739 |
|
|---|
| 740 | G4double au=0.8853*0.529*angstrom/std::sqrt(std::pow(z1, 0.6667)+std::pow(z2,0.6667));
|
|---|
| 741 | std::vector<G4double> r(npoints), phi(npoints);
|
|---|
| 742 |
|
|---|
| 743 | for(size_t i=0; i<npoints; i++) {
|
|---|
| 744 | G4double rr=(float)i/(float)(npoints-1);
|
|---|
| 745 | r[i]=rr*rr*rMax; // use quadratic r scale to make sampling fine near the center
|
|---|
| 746 | G4double sum=0.0;
|
|---|
| 747 | for(size_t j=0; j<ncoef; j++) sum+=coefs[j]*std::exp(scales[j]*r[i]/au);
|
|---|
| 748 | phi[i]=sum;
|
|---|
| 749 | }
|
|---|
| 750 |
|
|---|
| 751 | // compute the derivative at the origin for the spline
|
|---|
| 752 | G4double phiprime0=0.0;
|
|---|
| 753 | for(size_t j=0; j<ncoef; j++) phiprime0+=scales[j]*coefs[j]*std::exp(scales[j]*r[0]/au);
|
|---|
| 754 | phiprime0*=(1.0/au); // put back in natural units;
|
|---|
| 755 |
|
|---|
| 756 | *auval=au;
|
|---|
| 757 | return *static_cast<c2_function<G4double> *>(new lin_log_interpolating_function<G4double>(r, phi, false, phiprime0));
|
|---|
| 758 | }
|
|---|
| 759 |
|
|---|
| 760 | static c2_function<G4double> &LJScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
|
|---|
| 761 | {
|
|---|
| 762 | //from Loftager, Besenbacher, Jensen & Sorensen
|
|---|
| 763 | //PhysRev A20, 1443++, 1979
|
|---|
| 764 | G4double au=0.8853*0.529*angstrom/std::sqrt(std::pow(z1, 0.6667)+std::pow(z2,0.6667));
|
|---|
| 765 | std::vector<G4double> r(npoints), phi(npoints);
|
|---|
| 766 |
|
|---|
| 767 | for(size_t i=0; i<npoints; i++) {
|
|---|
| 768 | G4double rr=(float)i/(float)(npoints-1);
|
|---|
| 769 | r[i]=rr*rr*rMax; // use quadratic r scale to make sampling fine near the center
|
|---|
| 770 |
|
|---|
| 771 | G4double y=std::sqrt(9.67*r[i]/au);
|
|---|
| 772 | G4double ysq=y*y;
|
|---|
| 773 | G4double phipoly=1+y+0.3344*ysq+0.0485*y*ysq+0.002647*ysq*ysq;
|
|---|
| 774 | phi[i]=phipoly*std::exp(-y);
|
|---|
| 775 | // G4cout << r[i] << " " << phi[i] << G4endl;
|
|---|
| 776 | }
|
|---|
| 777 |
|
|---|
| 778 | // compute the derivative at the origin for the spline
|
|---|
| 779 | G4double logphiprime0=(9.67/2.0)*(2*0.3344-1.0); // #avoid 0/0 on first element
|
|---|
| 780 | logphiprime0 *= (1.0/au); // #put back in natural units
|
|---|
| 781 |
|
|---|
| 782 | *auval=au;
|
|---|
| 783 | return *static_cast<c2_function<G4double> *>(new lin_log_interpolating_function<G4double>(r, phi, false, logphiprime0*phi[0]));
|
|---|
| 784 | }
|
|---|
| 785 |
|
|---|
| 786 | G4NativeScreenedCoulombCrossSection::~G4NativeScreenedCoulombCrossSection() {
|
|---|
| 787 | }
|
|---|
| 788 |
|
|---|
| 789 | G4NativeScreenedCoulombCrossSection::G4NativeScreenedCoulombCrossSection() {
|
|---|
| 790 | AddScreeningFunction("zbl", ZBLScreening);
|
|---|
| 791 | AddScreeningFunction("lj", LJScreening);
|
|---|
| 792 | AddScreeningFunction("mol", MoliereScreening);
|
|---|
| 793 | }
|
|---|
| 794 |
|
|---|
| 795 | std::vector<G4String> G4NativeScreenedCoulombCrossSection::GetScreeningKeys() const {
|
|---|
| 796 | std::vector<G4String> keys;
|
|---|
| 797 | // find the available screening keys
|
|---|
| 798 | std::map<std::string, ScreeningFunc>::const_iterator sfunciter=phiMap.begin();
|
|---|
| 799 | for(; sfunciter != phiMap.end(); sfunciter++) keys.push_back((*sfunciter).first);
|
|---|
| 800 | return keys;
|
|---|
| 801 | }
|
|---|
| 802 |
|
|---|
| 803 | static inline G4double cm_energy(G4double a1, G4double a2, G4double t0) {
|
|---|
| 804 | // "relativistically correct energy in CM frame"
|
|---|
| 805 | G4double m1=a1*amu_c2, m2=a2*amu_c2;
|
|---|
| 806 | G4double mc2=(m1+m2);
|
|---|
| 807 | G4double f=2.0*m2*t0/(mc2*mc2);
|
|---|
| 808 | // old way: return (f < 1e-6) ? 0.5*mc2*f : mc2*(std::sqrt(1.0+f)-1.0);
|
|---|
| 809 | // formally equivalent to previous, but numerically stable for all f without conditional
|
|---|
| 810 | // uses identity (sqrt(1+x) - 1)(sqrt(1+x) + 1) = x
|
|---|
| 811 | return mc2*f/(std::sqrt(1.0+f)+1.0);
|
|---|
| 812 | }
|
|---|
| 813 |
|
|---|
| 814 | static inline G4double thetac(G4double m1, G4double m2, G4double eratio) {
|
|---|
| 815 | G4double s2th2=eratio*( (m1+m2)*(m1+m2)/(4.0*m1*m2) );
|
|---|
| 816 | G4double sth2=std::sqrt(s2th2);
|
|---|
| 817 | return 2.0*std::asin(sth2);
|
|---|
| 818 | }
|
|---|
| 819 |
|
|---|
| 820 | void G4NativeScreenedCoulombCrossSection::LoadData(G4String screeningKey, G4int z1, G4double a1, G4double recoilCutoff)
|
|---|
| 821 | {
|
|---|
| 822 | static const size_t sigLen=200; // since sigma doesn't matter much, a very coarse table will do
|
|---|
| 823 | G4DataVector energies(sigLen);
|
|---|
| 824 | G4DataVector data(sigLen);
|
|---|
| 825 |
|
|---|
| 826 | a1=standardmass(z1); // use standardized values for mass for building tables
|
|---|
| 827 |
|
|---|
| 828 | const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
|
|---|
| 829 | if (materialTable == 0)
|
|---|
| 830 | G4Exception("mhmNativeCrossSection::LoadData - no MaterialTable found)");
|
|---|
| 831 |
|
|---|
| 832 | G4int nMaterials = G4Material::GetNumberOfMaterials();
|
|---|
| 833 |
|
|---|
| 834 | for (G4int m=0; m<nMaterials; m++)
|
|---|
| 835 | {
|
|---|
| 836 | const G4Material* material= (*materialTable)[m];
|
|---|
| 837 | const G4ElementVector* elementVector = material->GetElementVector();
|
|---|
| 838 | const G4int nMatElements = material->GetNumberOfElements();
|
|---|
| 839 |
|
|---|
| 840 | for (G4int iEl=0; iEl<nMatElements; iEl++)
|
|---|
| 841 | {
|
|---|
| 842 | G4Element* element = (*elementVector)[iEl];
|
|---|
| 843 | G4int Z = (G4int) element->GetZ();
|
|---|
| 844 | G4double a2=element->GetA()*(mole/gram);
|
|---|
| 845 |
|
|---|
| 846 | if(sigmaMap.find(Z)!=sigmaMap.end()) continue; // we've already got this element
|
|---|
| 847 |
|
|---|
| 848 | // find the screening function generator we need
|
|---|
| 849 | std::map<std::string, ScreeningFunc>::iterator sfunciter=phiMap.find(screeningKey);
|
|---|
| 850 | if(sfunciter==phiMap.end()) {
|
|---|
| 851 | G4cout << "no such screening key " << screeningKey << G4endl; // FIXME later
|
|---|
| 852 | exit(1);
|
|---|
| 853 | }
|
|---|
| 854 | ScreeningFunc sfunc=(*sfunciter).second;
|
|---|
| 855 |
|
|---|
| 856 | G4double au;
|
|---|
| 857 | c2_function<G4double> &screen=sfunc(z1, Z, 200, 50.0*angstrom, &au); // generate the screening data
|
|---|
| 858 |
|
|---|
| 859 | G4ScreeningTables st;
|
|---|
| 860 | st.EMphiData=&screen; // this is our phi table
|
|---|
| 861 | st.z1=z1; st.m1=a1; st.z2=Z; st.m2=a2; st.emin=recoilCutoff;
|
|---|
| 862 | st.au=au;
|
|---|
| 863 |
|
|---|
| 864 | // now comes the hard part... build the total cross section tables from the phi table
|
|---|
| 865 | //based on (pi-thetac) = pi*beta*alpha/x0, but noting that alpha is very nearly unity, always
|
|---|
| 866 | //so just solve it wth alpha=1, which makes the solution much easier
|
|---|
| 867 | //this function returns an approximation to (beta/x0)^2=phi(x0)/(eps*x0)-1 ~ ((pi-thetac)/pi)^2
|
|---|
| 868 | //Since we don't need exact sigma values, this is good enough (within a factor of 2 almost always)
|
|---|
| 869 | //this rearranges to phi(x0)/(x0*eps) = 2*theta/pi - theta^2/pi^2
|
|---|
| 870 |
|
|---|
| 871 | c2_linear<G4double> c2au(0.0, 0.0, au);
|
|---|
| 872 | c2_composed_function<G4double> phiau(screen, c2au); // build phi(x*au) for dimensionless phi
|
|---|
| 873 | c2_linear<G4double> c2eps(0.0, 0.0, 0.0); // will store an appropriate eps inside this in loop
|
|---|
| 874 | c2_ratio<G4double> x0func(phiau, c2eps); // this will be phi(x)/(x*eps) when c2eps is correctly set
|
|---|
| 875 |
|
|---|
| 876 | G4double m1c2=a1*amu_c2;
|
|---|
| 877 | G4double escale=z1*Z*elm_coupling/au; // energy at screening distance
|
|---|
| 878 | G4double emax=m1c2; // model is doubtful in very relativistic range
|
|---|
| 879 | G4double eratkin=0.9999*(4*a1*a2)/((a1+a2)*(a1+a2)); // #maximum kinematic ratio possible at 180 degrees
|
|---|
| 880 | G4double cmfact0=st.emin/cm_energy(a1, a2, st.emin);
|
|---|
| 881 | G4double l1=std::log(emax);
|
|---|
| 882 | G4double l0=std::log(st.emin*cmfact0/eratkin);
|
|---|
| 883 |
|
|---|
| 884 | if(verbosity >=1)
|
|---|
| 885 | G4cout << "Native Screening: " << screeningKey << " " << z1 << " " << a1 << " " <<
|
|---|
| 886 | Z << " " << a2 << " " << recoilCutoff << G4endl;
|
|---|
| 887 |
|
|---|
| 888 | for(size_t idx=0; idx<sigLen; idx++) {
|
|---|
| 889 | G4double ee=std::exp(idx*((l1-l0)/sigLen)+l0);
|
|---|
| 890 | G4double gamma=1.0+ee/m1c2;
|
|---|
| 891 | G4double eratio=(cmfact0*st.emin)/ee; // factor by which ee needs to be reduced to get emin
|
|---|
| 892 | G4double theta=thetac(gamma*a1, a2, eratio);
|
|---|
| 893 |
|
|---|
| 894 | G4double eps=cm_energy(a1, a2, ee)/escale; // #make sure lab energy is converted to CM for these calculations
|
|---|
| 895 | c2eps.reset(0.0, 0.0, eps); // set correct slope in this function
|
|---|
| 896 |
|
|---|
| 897 | G4double q=theta/pi;
|
|---|
| 898 | // G4cout << ee << " " << m1c2 << " " << gamma << " " << eps << " " << theta << " " << q << G4endl;
|
|---|
| 899 | G4double x0= x0func.find_root(1e-6*angstrom/au, 0.9999*screen.xmax()/au, 1.0, 2*q-q*q);
|
|---|
| 900 | G4double betasquared=x0*x0 - x0*phiau(x0)/eps;
|
|---|
| 901 | G4double sigma=pi*betasquared*au*au;
|
|---|
| 902 | energies[idx]=ee;
|
|---|
| 903 | data[idx]=sigma;
|
|---|
| 904 | }
|
|---|
| 905 |
|
|---|
| 906 | screeningData[Z]=st;
|
|---|
| 907 | sigmaMap[Z] = static_cast<c2_function<G4double> *>(new log_log_interpolating_function<G4double>(
|
|---|
| 908 | energies, data));
|
|---|
| 909 | }
|
|---|
| 910 | }
|
|---|
| 911 | }
|
|---|