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 // 12-Apr-06 fix in delayed neutron and photon emission without FS data by T. Koi 31 // 07-Sep-11 M. Kelsey -- Follow change to G4HadFinalState interface 32 // P. Arce, June-2014 Conversion neutron_hp to particle_hp 33 // 34 35 #include "G4ParticleHPFissionFS.hh" 36 37 #include "G4DynamicParticleVector.hh" 38 #include "G4Exp.hh" 39 #include "G4IonTable.hh" 40 #include "G4Nucleus.hh" 41 #include "G4ParticleHPFissionERelease.hh" 42 #include "G4PhysicalConstants.hh" 43 #include "G4PhysicsModelCatalog.hh" 44 45 G4ParticleHPFissionFS::G4ParticleHPFissionFS() 46 { 47 secID = G4PhysicsModelCatalog::GetModelID("model_NeutronHPFission"); 48 hasXsec = false; 49 produceFissionFragments = false; 50 } 51 52 void G4ParticleHPFissionFS::Init(G4double A, G4double Z, G4int M, const G4String& dirName, 53 const G4String& aFSType, G4ParticleDefinition* projectile) 54 { 55 theFS.Init(A, Z, M, dirName, aFSType, projectile); 56 theFC.Init(A, Z, M, dirName, aFSType, projectile); 57 theSC.Init(A, Z, M, dirName, aFSType, projectile); 58 theTC.Init(A, Z, M, dirName, aFSType, projectile); 59 theLC.Init(A, Z, M, dirName, aFSType, projectile); 60 61 theFF.Init(A, Z, M, dirName, aFSType, projectile); 62 if (G4ParticleHPManager::GetInstance()->GetProduceFissionFragments() && theFF.HasFSData()) { 63 G4cout << "Fission fragment production is now activated in HP package for " 64 << "Z = " << (G4int)Z << ", A = " << (G4int)A << G4endl; 65 G4cout << "As currently modeled this option precludes production of delayed neutrons from " 66 "fission fragments." 67 << G4endl; 68 produceFissionFragments = true; 69 } 70 } 71 72 G4HadFinalState* G4ParticleHPFissionFS::ApplyYourself(const G4HadProjectile& theTrack) 73 { 74 // Because it may change by UI command 75 produceFissionFragments = G4ParticleHPManager::GetInstance()->GetProduceFissionFragments(); 76 77 // prepare neutron 78 if (theResult.Get() == nullptr) theResult.Put(new G4HadFinalState); 79 theResult.Get()->Clear(); 80 G4double eKinetic = theTrack.GetKineticEnergy(); 81 const G4HadProjectile* incidentParticle = &theTrack; 82 G4ReactionProduct theNeutron( 83 const_cast<G4ParticleDefinition*>(incidentParticle->GetDefinition())); 84 theNeutron.SetMomentum(incidentParticle->Get4Momentum().vect()); 85 theNeutron.SetKineticEnergy(eKinetic); 86 87 // prepare target 88 G4Nucleus aNucleus; 89 G4ReactionProduct theTarget; 90 G4double targetMass = theFS.GetMass(); 91 G4ThreeVector neuVelo = 92 (1. / incidentParticle->GetDefinition()->GetPDGMass()) * theNeutron.GetMomentum(); 93 theTarget = 94 aNucleus.GetBiasedThermalNucleus(targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature()); 95 theTarget.SetDefinition( 96 G4IonTable::GetIonTable()->GetIon(G4int(theBaseZ), G4int(theBaseA), 0.0)); // TESTPHP 97 98 // set neutron and target in the FS classes 99 theFS.SetNeutronRP(theNeutron); 100 theFS.SetTarget(theTarget); 101 theFC.SetNeutronRP(theNeutron); 102 theFC.SetTarget(theTarget); 103 theSC.SetNeutronRP(theNeutron); 104 theSC.SetTarget(theTarget); 105 theTC.SetNeutronRP(theNeutron); 106 theTC.SetTarget(theTarget); 107 theLC.SetNeutronRP(theNeutron); 108 theLC.SetTarget(theTarget); 109 theFF.SetNeutronRP(theNeutron); 110 theFF.SetTarget(theTarget); 111 112 // boost to target rest system and decide on channel. 113 theNeutron.Lorentz(theNeutron, -1 * theTarget); 114 115 // dice the photons 116 117 G4DynamicParticleVector* thePhotons; 118 thePhotons = theFS.GetPhotons(); 119 120 // select the FS in charge 121 122 eKinetic = theNeutron.GetKineticEnergy(); 123 G4double xSec[4]; 124 xSec[0] = theFC.GetXsec(eKinetic); 125 xSec[1] = xSec[0] + theSC.GetXsec(eKinetic); 126 xSec[2] = xSec[1] + theTC.GetXsec(eKinetic); 127 xSec[3] = xSec[2] + theLC.GetXsec(eKinetic); 128 G4int it; 129 unsigned int i = 0; 130 G4double random = G4UniformRand(); 131 if (xSec[3] == 0) { 132 it = -1; 133 } 134 else { 135 for (i = 0; i < 4; i++) { 136 it = i; 137 if (random < xSec[i] / xSec[3]) break; 138 } 139 } 140 141 // dice neutron multiplicities, energies and momenta in Lab. @@ 142 // no energy conservation on an event-to-event basis. we rely on the data to be ok. @@ 143 // also for mean, we rely on the consistancy of the data. @@ 144 145 G4int Prompt = 0, delayed = 0, all = 0; 146 G4DynamicParticleVector* theNeutrons = nullptr; 147 switch (it) // check logic, and ask, if partials can be assumed to correspond to individual 148 // particles @@@ 149 { 150 case 0: 151 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0); 152 if (Prompt == 0 && delayed == 0) Prompt = all; 153 theNeutrons = theFC.ApplyYourself(Prompt); // delayed always in FS 154 // take 'U' into account explicitly (see 5.4) in the sampling of energy @@@@ 155 break; 156 case 1: 157 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 1); 158 if (Prompt == 0 && delayed == 0) Prompt = all; 159 theNeutrons = theSC.ApplyYourself(Prompt); // delayed always in FS, off done in FSFissionFS 160 break; 161 case 2: 162 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 2); 163 if (Prompt == 0 && delayed == 0) Prompt = all; 164 theNeutrons = theTC.ApplyYourself(Prompt); // delayed always in FS 165 break; 166 case 3: 167 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 3); 168 if (Prompt == 0 && delayed == 0) Prompt = all; 169 theNeutrons = theLC.ApplyYourself(Prompt); // delayed always in FS 170 break; 171 default: 172 break; 173 } 174 175 // dice delayed neutrons and photons, and fallback 176 // for Prompt in case channel had no FS data; add all paricles to FS. 177 178 if (produceFissionFragments) delayed = 0; 179 180 G4double* theDecayConstants; 181 182 if (theNeutrons != nullptr) { 183 theDecayConstants = new G4double[delayed]; 184 for (i = 0; i < theNeutrons->size(); ++i) { 185 theResult.Get()->AddSecondary(theNeutrons->operator[](i), secID); 186 } 187 delete theNeutrons; 188 189 G4DynamicParticleVector* theDelayed = nullptr; 190 theDelayed = theFS.ApplyYourself(0, delayed, theDecayConstants); 191 for (i = 0; i < theDelayed->size(); i++) { 192 G4double time = -G4Log(G4UniformRand()) / theDecayConstants[i]; 193 time += theTrack.GetGlobalTime(); 194 theResult.Get()->AddSecondary(theDelayed->operator[](i), secID); 195 theResult.Get()->GetSecondary(theResult.Get()->GetNumberOfSecondaries() - 1)->SetTime(time); 196 } 197 delete theDelayed; 198 } 199 else { 200 theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0); 201 theDecayConstants = new G4double[delayed]; 202 if (Prompt == 0 && delayed == 0) Prompt = all; 203 theNeutrons = theFS.ApplyYourself(Prompt, delayed, theDecayConstants); 204 G4int i0; 205 for (i0 = 0; i0 < Prompt; ++i0) { 206 theResult.Get()->AddSecondary(theNeutrons->operator[](i0), secID); 207 } 208 209 for (i0 = Prompt; i0 < Prompt + delayed; ++i0) { 210 // Protect against the very rare case of division by zero 211 G4double time = 0.0; 212 if (theDecayConstants[i0 - Prompt] > 1.0e-30) { 213 time = -G4Log(G4UniformRand()) / theDecayConstants[i0 - Prompt]; 214 } 215 else { 216 G4ExceptionDescription ed; 217 ed << " theDecayConstants[i0-Prompt]=" << theDecayConstants[i0 - Prompt] 218 << " -> cannot sample the time : set it to 0.0 !" << G4endl; 219 G4Exception("G4ParticleHPFissionFS::ApplyYourself ", "HAD_FISSIONHP_001", JustWarning, ed); 220 } 221 222 time += theTrack.GetGlobalTime(); 223 theResult.Get()->AddSecondary(theNeutrons->operator[](i0), secID); 224 theResult.Get()->GetSecondary(theResult.Get()->GetNumberOfSecondaries() - 1)->SetTime(time); 225 } 226 delete theNeutrons; 227 } 228 delete[] theDecayConstants; 229 230 std::size_t nPhotons = 0; 231 if (thePhotons != nullptr) { 232 nPhotons = thePhotons->size(); 233 for (i = 0; i < nPhotons; ++i) { 234 theResult.Get()->AddSecondary(thePhotons->operator[](i), secID); 235 } 236 delete thePhotons; 237 } 238 239 // finally deal with local energy depositions. 240 241 G4ParticleHPFissionERelease* theERelease = theFS.GetEnergyRelease(); 242 G4double eDepByFragments = theERelease->GetFragmentKinetic(); 243 // theResult.SetLocalEnergyDeposit(eDepByFragments); 244 if (!produceFissionFragments) theResult.Get()->SetLocalEnergyDeposit(eDepByFragments); 245 // clean up the primary neutron 246 theResult.Get()->SetStatusChange(stopAndKill); 247 248 if (produceFissionFragments) { 249 G4int fragA_Z = 0; 250 G4int fragA_A = 0; 251 G4int fragA_M = 0; 252 // System is traget rest! 253 theFF.GetAFissionFragment(eKinetic, fragA_Z, fragA_A, fragA_M); 254 if (0 == fragA_A) { return theResult.Get(); } 255 G4int fragB_Z = (G4int)theBaseZ - fragA_Z; 256 G4int fragB_A = (G4int)theBaseA - fragA_A - Prompt; 257 258 G4IonTable* pt = G4IonTable::GetIonTable(); 259 // Excitation energy is not taken into account 260 G4ParticleDefinition* pdA = pt->GetIon(fragA_Z, fragA_A, 0.0); 261 G4ParticleDefinition* pdB = pt->GetIon(fragB_Z, fragB_A, 0.0); 262 263 // Isotropic Distribution 264 G4double phi = twopi * G4UniformRand(); 265 // Bug #1745 DHW G4double theta = pi*G4UniformRand(); 266 G4double costheta = 2. * G4UniformRand() - 1.; 267 G4double theta = std::acos(costheta); 268 G4double sinth = std::sin(theta); 269 G4ThreeVector direction(sinth * std::cos(phi), sinth * std::sin(phi), costheta); 270 271 // Just use ENDF value for this 272 G4double ER = eDepByFragments; 273 G4double ma = pdA->GetPDGMass(); 274 G4double mb = pdB->GetPDGMass(); 275 G4double EA = ER / (1 + ma / mb); 276 G4double EB = ER - EA; 277 auto dpA = new G4DynamicParticle(pdA, direction, EA); 278 auto dpB = new G4DynamicParticle(pdB, -direction, EB); 279 theResult.Get()->AddSecondary(dpA, secID); 280 theResult.Get()->AddSecondary(dpB, secID); 281 } 282 283 return theResult.Get(); 284 } 285