Geant4 Cross Reference |
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 // neutron_hp -- source file 27 // J.P. Wellisch, Nov-1996 28 // A prototype of the low energy neutron transport model. 29 // 30 // 070618 fix memory leaking by T. Koi 31 // 071002 enable cross section dump by T. Koi 32 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties:: 33 // 081124 Protect invalid read which caused run time errors by T. Koi 34 // 100729 Add safty for 0 lenght cross sections by T. Ko 35 // P. Arce, June-2014 Conversion neutron_hp to particle_hp 36 // 37 #include "G4ParticleHPFissionData.hh" 38 39 #include "G4ElementTable.hh" 40 #include "G4HadronicParameters.hh" 41 #include "G4Neutron.hh" 42 #include "G4NucleiProperties.hh" 43 #include "G4ParticleHPData.hh" 44 #include "G4ParticleHPManager.hh" 45 #include "G4PhysicalConstants.hh" 46 #include "G4Pow.hh" 47 #include "G4SystemOfUnits.hh" 48 49 G4ParticleHPFissionData::G4ParticleHPFissionData() : G4VCrossSectionDataSet("NeutronHPFissionXS") 50 { 51 SetMinKinEnergy(0 * MeV); 52 SetMaxKinEnergy(20 * MeV); 53 54 theCrossSections = nullptr; 55 instanceOfWorker = false; 56 if (G4Threading::IsWorkerThread()) { 57 instanceOfWorker = true; 58 } 59 element_cache = nullptr; 60 material_cache = nullptr; 61 ke_cache = 0.0; 62 xs_cache = 0.0; 63 } 64 65 G4ParticleHPFissionData::~G4ParticleHPFissionData() 66 { 67 if (theCrossSections != nullptr && !instanceOfWorker) { 68 theCrossSections->clearAndDestroy(); 69 delete theCrossSections; 70 theCrossSections = nullptr; 71 } 72 } 73 74 G4bool G4ParticleHPFissionData::IsIsoApplicable(const G4DynamicParticle* dp, G4int /*Z*/, 75 G4int /*A*/, const G4Element* /*elm*/, 76 const G4Material* /*mat*/) 77 { 78 G4double eKin = dp->GetKineticEnergy(); 79 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy() 80 && dp->GetDefinition() == G4Neutron::Neutron(); 81 } 82 83 G4double G4ParticleHPFissionData::GetIsoCrossSection(const G4DynamicParticle* dp, G4int /*Z*/, 84 G4int /*A*/, const G4Isotope* /*iso*/, 85 const G4Element* element, 86 const G4Material* material) 87 { 88 if (dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache) 89 return xs_cache; 90 91 ke_cache = dp->GetKineticEnergy(); 92 element_cache = element; 93 material_cache = material; 94 G4double xs = GetCrossSection(dp, element, material->GetTemperature()); 95 xs_cache = xs; 96 return xs; 97 } 98 99 void G4ParticleHPFissionData::BuildPhysicsTable(const G4ParticleDefinition&) 100 { 101 if (G4Threading::IsWorkerThread()) { 102 theCrossSections = G4ParticleHPManager::GetInstance()->GetFissionCrossSections(); 103 return; 104 } 105 106 std::size_t numberOfElements = G4Element::GetNumberOfElements(); 107 if (theCrossSections == nullptr) 108 theCrossSections = new G4PhysicsTable(numberOfElements); 109 else 110 theCrossSections->clearAndDestroy(); 111 112 // make a PhysicsVector for each element 113 114 auto theElementTable = G4Element::GetElementTable(); 115 for (std::size_t i = 0; i < numberOfElements; ++i) { 116 G4PhysicsVector* physVec = G4ParticleHPData::Instance(G4Neutron::Neutron()) 117 ->MakePhysicsVector((*theElementTable)[i], this); 118 theCrossSections->push_back(physVec); 119 } 120 121 G4ParticleHPManager::GetInstance()->RegisterFissionCrossSections(theCrossSections); 122 } 123 124 void G4ParticleHPFissionData::DumpPhysicsTable(const G4ParticleDefinition&) 125 { 126 #ifdef G4VERBOSE 127 if (G4HadronicParameters::Instance()->GetVerboseLevel() == 0) return; 128 129 // 130 // Dump element based cross section 131 // range 10e-5 eV to 20 MeV 132 // 10 point per decade 133 // in barn 134 // 135 G4cout << G4endl; 136 G4cout << G4endl; 137 G4cout << "Fission Cross Section of Neutron HP" << G4endl; 138 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl; 139 G4cout << G4endl; 140 G4cout << "Name of Element" << G4endl; 141 G4cout << "Energy[eV] XS[barn]" << G4endl; 142 G4cout << G4endl; 143 144 std::size_t numberOfElements = G4Element::GetNumberOfElements(); 145 auto theElementTable = G4Element::GetElementTable(); 146 147 for (std::size_t i = 0; i < numberOfElements; ++i) { 148 G4cout << (*theElementTable)[i]->GetName() << G4endl; 149 150 if ((*((*theCrossSections)(i))).GetVectorLength() == 0) { 151 G4cout << "The cross-section data of the fission of this element is not available." << G4endl; 152 G4cout << G4endl; 153 continue; 154 } 155 156 for (G4int ie = 0; ie < 130; ++ie) { 157 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * eV; 158 G4bool outOfRange = false; 159 160 if (eKinetic < 20 * MeV) { 161 G4cout << eKinetic / eV << " " 162 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / barn << G4endl; 163 } 164 } 165 166 G4cout << G4endl; 167 } 168 #endif 169 } 170 171 G4double G4ParticleHPFissionData::GetCrossSection(const G4DynamicParticle* aP, const G4Element* anE, 172 G4double aT) 173 { 174 G4double result = 0; 175 if (anE->GetZ() < 88) return result; 176 G4bool outOfRange; 177 auto index = (G4int)anE->GetIndex(); 178 179 if (((*theCrossSections)(index))->GetVectorLength() == 0) return result; 180 181 // prepare neutron 182 G4double eKinetic = aP->GetKineticEnergy(); 183 G4ReactionProduct theNeutronRP(aP->GetDefinition()); 184 theNeutronRP.SetMomentum(aP->GetMomentum()); 185 theNeutronRP.SetKineticEnergy(eKinetic); 186 187 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) { 188 // NEGLECT_DOPPLER 189 G4double factor = 1.0; 190 if (eKinetic < aT * k_Boltzmann) { 191 // below 0.1 eV neutrons 192 // Have to do some, but now just igonre. 193 // Will take care after performance check. 194 // factor = factor * targetV; 195 } 196 return ((*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange)) * factor; 197 } 198 199 // prepare thermal nucleus 200 G4Nucleus aNuc; 201 G4double eps = 0.0001; 202 G4double theA = anE->GetN(); 203 G4double theZ = anE->GetZ(); 204 G4double eleMass; 205 eleMass = (G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA + eps), 206 static_cast<G4int>(theZ + eps))) 207 / G4Neutron::Neutron()->GetPDGMass(); 208 209 G4ReactionProduct boosted; 210 G4double aXsection; 211 212 // MC integration loop 213 G4int counter = 0; 214 G4double buffer = 0; 215 G4int size = G4int(std::max(10., aT / 60 * kelvin)); 216 G4ThreeVector neutronVelocity = 217 1. / G4Neutron::Neutron()->GetPDGMass() * theNeutronRP.GetMomentum(); 218 G4double neutronVMag = neutronVelocity.mag(); 219 220 while (counter == 0 221 || std::abs(buffer - result / std::max(1, counter)) 222 > 0.01 * buffer) // Loop checking, 11.05.2015, T. Koi 223 { 224 if (counter != 0) buffer = result / counter; 225 while (counter < size) // Loop checking, 11.05.2015, T. Koi 226 { 227 counter++; 228 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT); 229 boosted.Lorentz(theNeutronRP, aThermalNuc); 230 G4double theEkin = boosted.GetKineticEnergy(); 231 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange); 232 // velocity correction. 233 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum(); 234 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag; 235 result += aXsection; 236 } 237 size += size; 238 } 239 result /= counter; 240 return result; 241 } 242 243 G4int G4ParticleHPFissionData::GetVerboseLevel() const 244 { 245 return G4ParticleHPManager::GetInstance()->GetVerboseLevel(); 246 } 247 248 void G4ParticleHPFissionData::SetVerboseLevel(G4int newValue) 249 { 250 G4ParticleHPManager::GetInstance()->SetVerboseLevel(newValue); 251 } 252 253 void G4ParticleHPFissionData::CrossSectionDescription(std::ostream& outFile) const 254 { 255 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF)\n" 256 << "for induced fission reaction of neutrons below 20MeV\n"; 257 } 258