Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPFission.cc
Go to the documentation of this file.
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// 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
31// 08-08-06 delete unnecessary and harmed declaration; Bug Report[857]
32//
33// P. Arce, June-2014 Conversion neutron_hp to particle_hp
34//
36
40#include "G4SystemOfUnits.hh"
41#include "G4Threading.hh"
42
44{
45 SetMinEnergy(0.0);
46 SetMaxEnergy(20. * MeV);
47}
48
50{
51 // Vector is shared, only master deletes it
52 // delete [] theFission;
54 if (theFission != nullptr) {
55 for (auto it = theFission->cbegin(); it != theFission->cend(); ++it) {
56 delete *it;
57 }
58 theFission->clear();
59 }
60 }
61}
62
64 G4Nucleus& aNucleus)
65{
67 const G4Material* theMaterial = aTrack.GetMaterial();
68 auto n = (G4int)theMaterial->GetNumberOfElements();
69 std::size_t index = theMaterial->GetElement(0)->GetIndex();
70 if (n != 1) {
71 auto xSec = new G4double[n];
72 G4double sum = 0;
73 G4int i;
74 const G4double* NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
75 G4double rWeight;
77 for (i = 0; i < n; ++i) {
78 index = theMaterial->GetElement(i)->GetIndex();
79 rWeight = NumAtomsPerVolume[i];
80 xSec[i] = ((*theFission)[index])
81 ->GetXsec(aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i),
82 theMaterial->GetTemperature()));
83 xSec[i] *= rWeight;
84 sum += xSec[i];
85 }
86 G4double random = G4UniformRand();
87 G4double running = 0;
88 for (i = 0; i < n; ++i) {
89 running += xSec[i];
90 index = theMaterial->GetElement(i)->GetIndex();
91 // if(random<=running/sum) break;
92 if (sum == 0 || random <= running / sum) break;
93 }
94 delete[] xSec;
95 }
96 // return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
97 G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack, -2);
98
99 // Overwrite target parameters
100 aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),
101 G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
102 const G4Element* target_element = (*G4Element::GetElementTable())[index];
103 const G4Isotope* target_isotope = nullptr;
104 auto iele = (G4int)target_element->GetNumberOfIsotopes();
105 for (G4int j = 0; j != iele; ++j) {
106 target_isotope = target_element->GetIsotope(j);
107 if (target_isotope->GetN()
109 break;
110 }
111 aNucleus.SetIsotope(target_isotope);
112
114 return result;
115}
116
117const std::pair<G4double, G4double> G4ParticleHPFission::GetFatalEnergyCheckLevels() const
118{
119 // max energy non-conservation is mass of heavy nucleus
120 return std::pair<G4double, G4double>(10.0 * perCent, 350.0 * CLHEP::GeV);
121}
122
127
132
134{
136
137 theFission = hpmanager->GetFissionFinalStates();
138
140 if (theFission == nullptr) theFission = new std::vector<G4ParticleHPChannel*>;
141
142 if (numEle == (G4int)G4Element::GetNumberOfElements()) return;
143
144 if (theFission->size() == G4Element::GetNumberOfElements()) {
146 return;
147 }
148
149 if (G4FindDataDir("G4NEUTRONHPDATA") == nullptr)
151 __FILE__, __LINE__,
152 "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
153 dirName = G4FindDataDir("G4NEUTRONHPDATA");
154 G4String tString = "/Fission";
155 dirName = dirName + tString;
156
157 for (G4int i = numEle; i < (G4int)G4Element::GetNumberOfElements(); ++i) {
158 theFission->push_back(new G4ParticleHPChannel);
159 if ((*(G4Element::GetElementTable()))[i]->GetZ() > 87) { // TK modified for ENDF-VII
160 ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
161 ((*theFission)[i])->Register(new G4ParticleHPFissionFS);
162 }
163 }
164 hpmanager->RegisterFissionFinalStates(theFission);
165 }
167}
168
169void G4ParticleHPFission::ModelDescription(std::ostream& outFile) const
170{
171 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF)\n"
172 << "for induced fission reaction of neutrons below 20MeV\n";
173}
const char * G4FindDataDir(const char *)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
static size_t GetNumberOfElements()
Definition G4Element.cc:393
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
size_t GetIndex() const
Definition G4Element.hh:159
size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Material * GetMaterial() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4int GetN() const
Definition G4Isotope.hh:83
G4double GetTemperature() const
const G4Element * GetElement(G4int iel) const
const G4double * GetVecNbOfAtomsPerVolume() const
std::size_t GetNumberOfElements() const
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition G4Nucleus.cc:319
void SetIsotope(const G4Isotope *iso)
Definition G4Nucleus.hh:114
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus) override
const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const override
void ModelDescription(std::ostream &outFile) const override
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates() const
void RegisterFissionFinalStates(std::vector< G4ParticleHPChannel * > *val)
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4bool IsMasterThread()