[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 | // * * |
---|
| 21 | // * Parts of this code which have been developed by QinetiQ Ltd * |
---|
| 22 | // * under contract to the European Space Agency (ESA) are the * |
---|
| 23 | // * intellectual property of ESA. Rights to use, copy, modify and * |
---|
| 24 | // * redistribute this software for general public use are granted * |
---|
| 25 | // * in compliance with any licensing, distribution and development * |
---|
| 26 | // * policy adopted by the Geant4 Collaboration. This code has been * |
---|
| 27 | // * written by QinetiQ Ltd for the European Space Agency, under ESA * |
---|
| 28 | // * contract 17191/03/NL/LvH (Aurora Programme). * |
---|
| 29 | // * * |
---|
| 30 | // * By using, copying, modifying or distributing the software (or * |
---|
| 31 | // * any work based on the software) you agree to acknowledge its * |
---|
| 32 | // * use in resulting scientific publications, and indicate your * |
---|
| 33 | // * acceptance of all terms of the Geant4 Software license. * |
---|
| 34 | // ******************************************************************** |
---|
| 35 | // |
---|
[1340] | 36 | // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
[819] | 37 | // |
---|
| 38 | // MODULE: G4EMDissociationCrossSection.cc |
---|
| 39 | // |
---|
| 40 | // Version: B.1 |
---|
| 41 | // Date: 15/04/04 |
---|
| 42 | // Author: P R Truscott |
---|
| 43 | // Organisation: QinetiQ Ltd, UK |
---|
| 44 | // Customer: ESA/ESTEC, NOORDWIJK |
---|
| 45 | // Contract: 17191/03/NL/LvH |
---|
| 46 | // |
---|
[1340] | 47 | // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
[819] | 48 | // |
---|
| 49 | // CHANGE HISTORY |
---|
| 50 | // -------------- |
---|
| 51 | // |
---|
| 52 | // 17 October 2003, P R Truscott, QinetiQ Ltd, UK |
---|
| 53 | // Created. |
---|
| 54 | // |
---|
| 55 | // 15 March 2004, P R Truscott, QinetiQ Ltd, UK |
---|
| 56 | // Beta release |
---|
| 57 | // |
---|
[1347] | 58 | // 30 May 2005, J.P. Wellisch removed a compilation warning on gcc 3.4 for |
---|
[1340] | 59 | // geant4 7.1. |
---|
[1347] | 60 | // 09 November 2010, V.Ivanchenko make class applicable for Hydrogen but |
---|
| 61 | // set cross section for Hydrogen to zero |
---|
[819] | 62 | // |
---|
[1340] | 63 | // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
| 64 | ////////////////////////////////////////////////////////////////////////////// |
---|
[819] | 65 | // |
---|
| 66 | #include "G4EMDissociationCrossSection.hh" |
---|
| 67 | #include "G4PhysicsFreeVector.hh" |
---|
| 68 | #include "G4ParticleTable.hh" |
---|
| 69 | #include "G4IonTable.hh" |
---|
[1340] | 70 | #include "G4HadTmpUtil.hh" |
---|
[819] | 71 | #include "globals.hh" |
---|
[1340] | 72 | |
---|
| 73 | |
---|
[819] | 74 | G4EMDissociationCrossSection::G4EMDissociationCrossSection () |
---|
| 75 | { |
---|
| 76 | // |
---|
| 77 | // This function makes use of the class which can sample the virtual photon |
---|
| 78 | // spectrum, G4EMDissociationSpectrum. |
---|
| 79 | // |
---|
| 80 | thePhotonSpectrum = new G4EMDissociationSpectrum(); |
---|
| 81 | // |
---|
| 82 | // |
---|
| 83 | // Define other constants. |
---|
| 84 | // |
---|
| 85 | r0 = 1.18 * fermi; |
---|
| 86 | J = 36.8 * MeV; |
---|
| 87 | Qprime = 17.0 * MeV; |
---|
| 88 | epsilon = 0.0768; |
---|
| 89 | xd = 0.25; |
---|
| 90 | } |
---|
[1340] | 91 | ////////////////////////////////////////////////////////////////////////////// |
---|
[819] | 92 | // |
---|
| 93 | G4EMDissociationCrossSection::~G4EMDissociationCrossSection() |
---|
| 94 | { |
---|
| 95 | delete thePhotonSpectrum; |
---|
| 96 | } |
---|
[1340] | 97 | ///////////////////////////////////////////////////////////////////////////// |
---|
[819] | 98 | // |
---|
[1340] | 99 | G4bool |
---|
| 100 | G4EMDissociationCrossSection::IsIsoApplicable(const G4DynamicParticle* theDynamicParticle, |
---|
[1347] | 101 | G4int /*ZZ*/, G4int/* AA*/) |
---|
[819] | 102 | { |
---|
| 103 | // |
---|
| 104 | // The condition for the applicability of this class is that the projectile |
---|
| 105 | // must be an ion and the target must have more than one nucleon. In reality |
---|
| 106 | // the value of A for either the projectile or target could be much higher, |
---|
| 107 | // since for cases where both he projectile and target are medium to small |
---|
| 108 | // Z, the probability of the EMD process is, I think, VERY small. |
---|
| 109 | // |
---|
| 110 | if (G4ParticleTable::GetParticleTable()->GetIonTable()-> |
---|
[1347] | 111 | IsIon(theDynamicParticle->GetDefinition()) /*&& AA > 1*/) |
---|
[819] | 112 | return true; |
---|
| 113 | else |
---|
| 114 | return false; |
---|
| 115 | } |
---|
| 116 | |
---|
| 117 | G4bool G4EMDissociationCrossSection::IsApplicable |
---|
| 118 | (const G4DynamicParticle* theDynamicParticle, const G4Element* theElement) |
---|
| 119 | { |
---|
[1340] | 120 | return IsIsoApplicable(theDynamicParticle, 0, G4lrint(theElement->GetN())); |
---|
[819] | 121 | } |
---|
| 122 | |
---|
| 123 | ////////////////////////////////////////////////////////////////////////////// |
---|
| 124 | // |
---|
| 125 | G4double G4EMDissociationCrossSection::GetCrossSection |
---|
| 126 | (const G4DynamicParticle* theDynamicParticle, const G4Element* theElement, |
---|
| 127 | G4double temperature) |
---|
| 128 | { |
---|
| 129 | G4int nIso = theElement->GetNumberOfIsotopes(); |
---|
| 130 | G4double crossSection = 0; |
---|
[1347] | 131 | |
---|
| 132 | // VI protection for Hydrogen |
---|
| 133 | if(theElement->GetZ() < 1.5) { return crossSection; } |
---|
[819] | 134 | |
---|
| 135 | if (nIso) { |
---|
| 136 | G4double sig; |
---|
| 137 | G4IsotopeVector* isoVector = theElement->GetIsotopeVector(); |
---|
| 138 | G4double* abundVector = theElement->GetRelativeAbundanceVector(); |
---|
[1340] | 139 | G4int ZZ; |
---|
| 140 | G4int AA; |
---|
[819] | 141 | |
---|
| 142 | for (G4int i = 0; i < nIso; i++) { |
---|
[1340] | 143 | ZZ = (*isoVector)[i]->GetZ(); |
---|
| 144 | AA = (*isoVector)[i]->GetN(); |
---|
| 145 | sig = GetZandACrossSection(theDynamicParticle, ZZ, AA, temperature); |
---|
[819] | 146 | crossSection += sig*abundVector[i]; |
---|
| 147 | } |
---|
| 148 | |
---|
| 149 | } else { |
---|
[1340] | 150 | G4int ZZ = G4lrint(theElement->GetZ()); |
---|
| 151 | G4int AA = G4lrint(theElement->GetN()); |
---|
| 152 | crossSection = GetZandACrossSection(theDynamicParticle, ZZ, AA, temperature); |
---|
[819] | 153 | } |
---|
| 154 | |
---|
| 155 | return crossSection; |
---|
| 156 | } |
---|
| 157 | |
---|
| 158 | |
---|
[1340] | 159 | G4double |
---|
| 160 | G4EMDissociationCrossSection::GetZandACrossSection(const G4DynamicParticle *theDynamicParticle, |
---|
| 161 | G4int ZZ, G4int AA, G4double /*temperature*/) |
---|
[819] | 162 | { |
---|
[1347] | 163 | // VI protection for Hydrogen |
---|
| 164 | if(ZZ <= 1) { return 0.0; } |
---|
| 165 | |
---|
[819] | 166 | // |
---|
| 167 | // Get relevant information about the projectile and target (A, Z) and |
---|
| 168 | // velocity of the projectile. |
---|
| 169 | // |
---|
| 170 | G4ParticleDefinition *definitionP = theDynamicParticle->GetDefinition(); |
---|
| 171 | G4double AP = definitionP->GetBaryonNumber(); |
---|
| 172 | G4double ZP = definitionP->GetPDGCharge(); |
---|
| 173 | G4double b = theDynamicParticle->Get4Momentum().beta(); |
---|
| 174 | // G4double bsq = b * b; |
---|
| 175 | |
---|
| 176 | G4double AT = AA; |
---|
| 177 | G4double ZT = ZZ; |
---|
| 178 | G4double bmin = thePhotonSpectrum->GetClosestApproach(AP, ZP, AT, ZT, b); |
---|
| 179 | // |
---|
| 180 | // |
---|
| 181 | // Calculate the cross-section for the projectile and then the target. The |
---|
| 182 | // information is returned in a G4PhysicsFreeVector, which separates out the |
---|
| 183 | // cross-sections for the E1 and E2 moments of the virtual photon field, and |
---|
| 184 | // the energies (GDR and GQR). |
---|
| 185 | // |
---|
| 186 | G4PhysicsFreeVector *theProjectileCrossSections = |
---|
| 187 | GetCrossSectionForProjectile (AP, ZP, AT, ZT, b, bmin); |
---|
| 188 | G4double crossSection = |
---|
| 189 | (*theProjectileCrossSections)[0]+(*theProjectileCrossSections)[1]; |
---|
| 190 | delete theProjectileCrossSections; |
---|
| 191 | G4PhysicsFreeVector *theTargetCrossSections = |
---|
| 192 | GetCrossSectionForTarget (AP, ZP, AT, ZT, b, bmin); |
---|
| 193 | crossSection += |
---|
| 194 | (*theTargetCrossSections)[0]+(*theTargetCrossSections)[1]; |
---|
| 195 | delete theTargetCrossSections; |
---|
| 196 | |
---|
| 197 | return crossSection; |
---|
| 198 | } |
---|
| 199 | //////////////////////////////////////////////////////////////////////////////// |
---|
| 200 | // |
---|
| 201 | G4PhysicsFreeVector * |
---|
| 202 | G4EMDissociationCrossSection::GetCrossSectionForProjectile (G4double AP, |
---|
| 203 | G4double ZP, G4double /* AT */, G4double ZT, G4double b, G4double bmin) |
---|
| 204 | { |
---|
| 205 | // |
---|
| 206 | // |
---|
| 207 | // Use Wilson et al's approach to calculate the cross-sections due to the E1 |
---|
| 208 | // and E2 moments of the field at the giant dipole and quadrupole resonances |
---|
| 209 | // respectively, Note that the algorithm is traditionally applied to the |
---|
| 210 | // EMD break-up of the projectile in the field of the target, as is implemented |
---|
| 211 | // here. |
---|
| 212 | // |
---|
| 213 | // Initialise variables and calculate the energies for the GDR and GQR. |
---|
| 214 | // |
---|
| 215 | G4double AProot3 = std::pow(AP,1.0/3.0); |
---|
| 216 | G4double u = 3.0 * J / Qprime / AProot3; |
---|
| 217 | G4double R0 = r0 * AProot3; |
---|
| 218 | G4double E_GDR = hbarc / std::sqrt(0.7*amu_c2*R0*R0/8.0/J* |
---|
| 219 | (1.0 + u - (1.0 + epsilon + 3.0*u)/(1.0 + epsilon + u)*epsilon)); |
---|
| 220 | G4double E_GQR = 63.0 * MeV / AProot3; |
---|
| 221 | // |
---|
| 222 | // |
---|
| 223 | // Determine the virtual photon spectra at these energies. |
---|
| 224 | // |
---|
| 225 | G4double ZTsq = ZT * ZT; |
---|
| 226 | G4double nE1 = ZTsq * |
---|
| 227 | thePhotonSpectrum->GetGeneralE1Spectrum(E_GDR, b, bmin); |
---|
| 228 | G4double nE2 = ZTsq * |
---|
| 229 | thePhotonSpectrum->GetGeneralE2Spectrum(E_GQR, b, bmin); |
---|
| 230 | // |
---|
| 231 | // |
---|
| 232 | // Now calculate the cross-section of the projectile for interaction with the |
---|
| 233 | // E1 and E2 fields. |
---|
| 234 | // |
---|
| 235 | G4double sE1 = 60.0 * millibarn * MeV * (AP-ZP)*ZP/AP; |
---|
| 236 | G4double sE2 = 0.22 * microbarn / MeV * ZP * AProot3 * AProot3; |
---|
| 237 | if (AP > 100.0) sE2 *= 0.9; |
---|
| 238 | else if (AP > 40.0) sE2 *= 0.6; |
---|
| 239 | else sE2 *= 0.3; |
---|
| 240 | // |
---|
| 241 | // |
---|
| 242 | // ... and multiply with the intensity of the virtual photon spectra to get |
---|
| 243 | // the probability of interaction. |
---|
| 244 | // |
---|
| 245 | G4PhysicsFreeVector *theCrossSectionVector = new G4PhysicsFreeVector(2); |
---|
| 246 | theCrossSectionVector->PutValue(0, E_GDR, sE1*nE1); |
---|
| 247 | theCrossSectionVector->PutValue(1, E_GQR, sE2*nE2*E_GQR*E_GQR); |
---|
| 248 | |
---|
| 249 | return theCrossSectionVector; |
---|
| 250 | } |
---|
[1340] | 251 | |
---|
[819] | 252 | //////////////////////////////////////////////////////////////////////////////// |
---|
| 253 | // |
---|
| 254 | G4PhysicsFreeVector * |
---|
| 255 | G4EMDissociationCrossSection::GetCrossSectionForTarget (G4double AP, |
---|
| 256 | G4double ZP, G4double AT, G4double ZT, G4double b, G4double bmin) |
---|
| 257 | { |
---|
| 258 | // |
---|
| 259 | // This is a cheaky little member function to calculate the probability of |
---|
| 260 | // EMD for the target in the field of the projectile ... just by reversing the |
---|
| 261 | // A and Z's for the participants. |
---|
| 262 | // |
---|
| 263 | return GetCrossSectionForProjectile (AT, ZT, AP, ZP, b, bmin); |
---|
| 264 | } |
---|
[1340] | 265 | |
---|
| 266 | |
---|
[819] | 267 | G4double |
---|
[1340] | 268 | G4EMDissociationCrossSection::GetWilsonProbabilityForProtonDissociation(G4double A, |
---|
| 269 | G4double Z) |
---|
[819] | 270 | { |
---|
| 271 | // |
---|
| 272 | // This is a simple algorithm to choose whether a proton or neutron is ejected |
---|
| 273 | // from the nucleus in the EMD interaction. |
---|
| 274 | // |
---|
| 275 | G4double p = 0.0; |
---|
| 276 | if (Z < 6.0) |
---|
| 277 | p = 0.5; |
---|
| 278 | else if (Z < 8.0) |
---|
| 279 | p = 0.6; |
---|
| 280 | else if (Z < 14.0) |
---|
| 281 | p = 0.7; |
---|
| 282 | else |
---|
| 283 | { |
---|
| 284 | G4double p1 = (G4double) Z / (G4double) A; |
---|
| 285 | G4double p2 = 1.95*std::exp(-0.075*Z); |
---|
| 286 | if (p1 < p2) p = p1; |
---|
| 287 | else p = p2; |
---|
| 288 | } |
---|
| 289 | |
---|
| 290 | return p; |
---|
| 291 | } |
---|