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