| 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 | // 18-Sep-2003 First version is written by T. Koi
|
|---|
| 27 | // 10-Nov-2003 Bug fix at Cal. ke_per_n and D T. Koi
|
|---|
| 28 | // 12-Nov-2003 Add energy check at lower side T. Koi
|
|---|
| 29 | // 26-Dec-2006 Add isotope dependence D. Wright
|
|---|
| 30 |
|
|---|
| 31 | #include "G4IonsKoxCrossSection.hh"
|
|---|
| 32 | #include "G4ParticleTable.hh"
|
|---|
| 33 | #include "G4IonTable.hh"
|
|---|
| 34 | #include "G4HadTmpUtil.hh"
|
|---|
| 35 |
|
|---|
| 36 | G4double G4IonsKoxCrossSection::
|
|---|
| 37 | GetZandACrossSection(const G4DynamicParticle* aParticle, G4int ZZ,
|
|---|
| 38 | G4int AA, G4double /*temperature*/)
|
|---|
| 39 | {
|
|---|
| 40 | G4double xsection = 0.0;
|
|---|
| 41 |
|
|---|
| 42 | G4int Ap = aParticle->GetDefinition()->GetBaryonNumber();
|
|---|
| 43 | G4int Zp = G4int(aParticle->GetDefinition()->GetPDGCharge() / eplus + 0.5);
|
|---|
| 44 | G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
|
|---|
| 45 |
|
|---|
| 46 | // Apply energy check, if less than lower limit then 0 value is returned
|
|---|
| 47 | // if ( ke_per_N < lowerLimit ) return xsection;
|
|---|
| 48 |
|
|---|
| 49 | G4int At = AA;
|
|---|
| 50 | G4int Zt = ZZ;
|
|---|
| 51 |
|
|---|
| 52 | G4double one_third = 1.0 / 3.0;
|
|---|
| 53 |
|
|---|
| 54 | G4double cubicrAt = std::pow ( G4double(At) , G4double(one_third) );
|
|---|
| 55 | G4double cubicrAp = std::pow ( G4double(Ap) , G4double(one_third) );
|
|---|
| 56 |
|
|---|
| 57 | // rc divide fermi
|
|---|
| 58 | G4double Bc = Zt * Zp / ( (rc/fermi) * (cubicrAp+cubicrAt) );
|
|---|
| 59 |
|
|---|
| 60 | G4double targ_mass =
|
|---|
| 61 | G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(Zt, At);
|
|---|
| 62 | G4double proj_mass = aParticle->GetMass();
|
|---|
| 63 | G4double proj_momentum = aParticle->GetMomentum().mag();
|
|---|
| 64 |
|
|---|
| 65 | G4double Ecm = calEcm ( proj_mass , targ_mass , proj_momentum );
|
|---|
| 66 | if( Ecm <= Bc) return xsection;
|
|---|
| 67 |
|
|---|
| 68 | G4double Rvol = r0 * ( cubicrAp + cubicrAt );
|
|---|
| 69 |
|
|---|
| 70 | // G4double ke_per_N = aParticle->GetKineticEnergy() / Ap;
|
|---|
| 71 | G4double c = calCeValue ( ke_per_N / MeV );
|
|---|
| 72 |
|
|---|
| 73 | G4double a = 1.85;
|
|---|
| 74 | G4double Rsurf = r0 * (a*cubicrAp * cubicrAt/(cubicrAp + cubicrAt) - c);
|
|---|
| 75 | G4double D = 5.0 * ( At - 2 * Zt ) * Zp / ( Ap * At );
|
|---|
| 76 | Rsurf = Rsurf + D * fermi; // multiply D by fermi
|
|---|
| 77 |
|
|---|
| 78 | G4double Rint = Rvol + Rsurf;
|
|---|
| 79 | xsection = pi * Rint * Rint * ( 1 - Bc / ( Ecm / MeV ) );
|
|---|
| 80 |
|
|---|
| 81 | return xsection;
|
|---|
| 82 | }
|
|---|
| 83 |
|
|---|
| 84 |
|
|---|
| 85 | G4double G4IonsKoxCrossSection::
|
|---|
| 86 | GetCrossSection(const G4DynamicParticle* aParticle,
|
|---|
| 87 | const G4Element* anElement, G4double temperature)
|
|---|
| 88 | {
|
|---|
| 89 | G4int nIso = anElement->GetNumberOfIsotopes();
|
|---|
| 90 | G4double xsection = 0;
|
|---|
| 91 |
|
|---|
| 92 | if (nIso) {
|
|---|
| 93 | G4double sig;
|
|---|
| 94 | G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
|
|---|
| 95 | G4double* abundVector = anElement->GetRelativeAbundanceVector();
|
|---|
| 96 | G4int ZZ;
|
|---|
| 97 | G4int AA;
|
|---|
| 98 |
|
|---|
| 99 | for (G4int i = 0; i < nIso; i++) {
|
|---|
| 100 | ZZ = (*isoVector)[i]->GetZ();
|
|---|
| 101 | AA = (*isoVector)[i]->GetN();
|
|---|
| 102 | sig = GetZandACrossSection(aParticle, ZZ, AA, temperature);
|
|---|
| 103 | xsection += sig*abundVector[i];
|
|---|
| 104 | }
|
|---|
| 105 |
|
|---|
| 106 | } else {
|
|---|
| 107 | G4int ZZ = G4lrint(anElement->GetZ());
|
|---|
| 108 | G4int AA = G4lrint(anElement->GetN());
|
|---|
| 109 | xsection = GetIsoZACrossSection(aParticle, ZZ, AA, temperature);
|
|---|
| 110 | }
|
|---|
| 111 |
|
|---|
| 112 | return xsection;
|
|---|
| 113 | }
|
|---|
| 114 |
|
|---|
| 115 |
|
|---|
| 116 | G4double
|
|---|
| 117 | G4IonsKoxCrossSection::calEcm(G4double mp, G4double mt, G4double Plab)
|
|---|
| 118 | {
|
|---|
| 119 | G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
|
|---|
| 120 | G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
|
|---|
| 121 | G4double Pcm = Plab * mt / Ecm;
|
|---|
| 122 | G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
|
|---|
| 123 | return KEcm;
|
|---|
| 124 | }
|
|---|
| 125 |
|
|---|
| 126 |
|
|---|
| 127 | G4double G4IonsKoxCrossSection::calCeValue(const G4double ke)
|
|---|
| 128 | {
|
|---|
| 129 | // Calculate c value
|
|---|
| 130 | // This value is indepenent from projectile and target particle
|
|---|
| 131 | // ke is projectile kinetic energy per nucleon in the Lab system with MeV unit
|
|---|
| 132 | // fitting function is made by T. Koi
|
|---|
| 133 | // There are no data below 30 MeV/n in Kox et al.,
|
|---|
| 134 |
|
|---|
| 135 | G4double Ce;
|
|---|
| 136 | G4double log10_ke = std::log10 ( ke );
|
|---|
| 137 | if (log10_ke > 1.5)
|
|---|
| 138 | {
|
|---|
| 139 | Ce = - 10.0 / std::pow ( G4double(log10_ke) , G4double(5) ) + 2.0;
|
|---|
| 140 | }
|
|---|
| 141 | else
|
|---|
| 142 | {
|
|---|
| 143 | Ce = (-10.0/std::pow(G4double(1.5), G4double(5) ) + 2.0) /
|
|---|
| 144 | std::pow(G4double(1.5), G4double(3)) * std::pow(G4double(log10_ke), G4double(3) );
|
|---|
| 145 |
|
|---|
| 146 | }
|
|---|
| 147 | return Ce;
|
|---|
| 148 | }
|
|---|