Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
37#include "G4SystemOfUnits.hh"
40#include "G4Threading.hh"
41
43 : G4HadronicInteraction("NeutronHPFission"), theFission(nullptr), numEle(0)
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 {
65
67 const G4Material * theMaterial = aTrack.GetMaterial();
68 G4int n = (G4int)theMaterial->GetNumberOfElements();
69 std::size_t index = theMaterial->GetElement(0)->GetIndex();
70 if(n!=1)
71 {
72 G4double* xSec = new G4double[n];
73 G4double sum=0;
74 G4int i;
75 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
76 G4double rWeight;
78 for (i=0; i<n; ++i)
79 {
80 index = theMaterial->GetElement(i)->GetIndex();
81 rWeight = NumAtomsPerVolume[i];
82 xSec[i] = ((*theFission)[index])
83 ->GetXsec(aThermalE.GetThermalEnergy(aTrack,
84 theMaterial->GetElement(i),
85 theMaterial->GetTemperature()));
86 xSec[i] *= rWeight;
87 sum+=xSec[i];
88 }
89 G4double random = G4UniformRand();
90 G4double running = 0;
91 for (i=0; i<n; ++i)
92 {
93 running += xSec[i];
94 index = theMaterial->GetElement(i)->GetIndex();
95 //if(random<=running/sum) break;
96 if( sum == 0 || random <= running/sum ) break;
97 }
98 delete [] xSec;
99 }
100 //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
101 G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack,-2);
102
103 //Overwrite target parameters
105 const G4Element* target_element = (*G4Element::GetElementTable())[index];
106 const G4Isotope* target_isotope=NULL;
107 G4int iele = (G4int)target_element->GetNumberOfIsotopes();
108 for ( G4int j = 0 ; j != iele ; ++j ) {
109 target_isotope=target_element->GetIsotope( j );
110 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
111 }
112 aNucleus.SetIsotope( target_isotope );
113
115 return result;
116 }
117
118const std::pair<G4double, G4double> G4ParticleHPFission::GetFatalEnergyCheckLevels() const
119{
120 // max energy non-conservation is mass of heavy nucleus
121 return std::pair<G4double, G4double>(10.0*perCent, 350.0*CLHEP::GeV);
122}
123
125{
127}
128
130{
132}
133
135{
136
138
139 theFission = hpmanager->GetFissionFinalStates();
140
142
143 if ( theFission == nullptr ) theFission = new std::vector<G4ParticleHPChannel*>;
144
145 if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
146
147 if ( theFission->size() == G4Element::GetNumberOfElements() ) {
149 return;
150 }
151
152 if ( !G4FindDataDir("G4NEUTRONHPDATA") )
153 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
154 dirName = G4FindDataDir("G4NEUTRONHPDATA");
155 G4String tString = "/Fission";
156 dirName = dirName + tString;
157
158 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; ++i ) {
159 theFission->push_back( new G4ParticleHPChannel );
160 if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
161 ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
162 ((*theFission)[i])->Register( new G4ParticleHPFissionFS );
163 }
164 }
165 hpmanager->RegisterFissionFinalStates( theFission );
166 }
168}
169
170void G4ParticleHPFission::ModelDescription(std::ostream& outFile) const
171{
172 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF)\n"
173 << "for induced fission reaction of neutrons below 20MeV\n";
174}
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:403
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetIndex() const
Definition: G4Element.hh:182
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4Material * GetMaterial() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4int GetN() const
Definition: G4Isotope.hh:93
G4double GetTemperature() const
Definition: G4Material.hh:177
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition: G4Nucleus.cc:307
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:114
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void ModelDescription(std::ostream &outFile) const
void RegisterFissionFinalStates(std::vector< G4ParticleHPChannel * > *val)
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates()
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4bool IsMasterThread()
Definition: G4Threading.cc:124