Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPFissionURR.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//
27// -------------------------------------------------------------------
28//
29// Geant4 source file
30//
31// File name: G4ParticleHPFissionURR.cc
32//
33// Authors: Marek Zmeskal (CTU, Czech Technical University in Prague, Czech Republic)
34// Loic Thulliez (CEA France)
35//
36// Creation date: 4 June 2024
37//
38// Description: Class to handle URR range, can be omitted once the
39// proper isotope cross-section is stored in ParticleHP.
40//
41// Modifications:
42//
43// -------------------------------------------------------------------
44//
45//
46
53#include "G4SystemOfUnits.hh"
54#include "G4Threading.hh"
55
56
58 SetMinEnergy( 0.0 * CLHEP::eV );
59 SetMaxEnergy( 20.0 * CLHEP::MeV );
60 particleHPfission = new G4ParticleHPFission;
61}
62
63
65
66
68 const G4Material* theMaterial = aTrack.GetMaterial();
69 G4double kineticEnergy = aTrack.GetKineticEnergy();
70 G4HadFinalState* theFinalState = nullptr;
71 if ( kineticEnergy < (*URRlimits).back().first || kineticEnergy > (*URRlimits).back().second ) {
72 return particleHPfission->ApplyYourself( aTrack, aNucleus );
73 }
74 G4int elementI = -1;
75 G4int isotopeJ = -1;
76 G4int A = aNucleus.GetA_asInt();
77 G4int Z = aNucleus.GetZ_asInt();
78 // finds the element and isotope of the selected target aNucleus
79 for ( G4int i = 0; i < (G4int)theMaterial->GetNumberOfElements(); ++i ) {
80 if ( Z == theMaterial->GetElement(i)->GetZasInt() ) {
81 for ( G4int j = 0; j < (G4int)theMaterial->GetElement(i)->GetNumberOfIsotopes(); ++j ) {
82 if ( A == theMaterial->GetElement(i)->GetIsotope(j)->GetN() ) {
83 isotopeJ = j;
84 break;
85 }
86 }
87 // the loop cannot be ended here because the material can have two elements with same Z but different isotopic composition
88 if ( isotopeJ != -1 ) {
89 // isotope was found and for loop is ended
90 elementI = (G4int)theMaterial->GetElement(i)->GetIndex();
91 break;
92 }
93 } // end if find element
94 } // end element loop
95 // Check whether the energy is out of the URR limits for the given element
96 if ( kineticEnergy < (*URRlimits).at(elementI).first || kineticEnergy > (*URRlimits).at(elementI).second ) {
97 // Call fission final state in G4ParicleHPChannel and SELECT ISOTOPE (to be improved in the future)
99 theFinalState = (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->ApplyYourself( aTrack, -2 );
100 // Update target nucleus information according to the selected isotope
102 aNucleus.SetParameters( selectedIsotope_A, Z );
103 const G4Element* target_element = (*G4Element::GetElementTable())[elementI];
104 const G4Isotope* target_isotope = nullptr;
105 // Find the selected isotope among in the element
106 for ( G4int j = 0; j < (G4int)target_element->GetNumberOfIsotopes(); ++j ) {
107 target_isotope = target_element->GetIsotope(j);
108 if ( target_isotope->GetN() == selectedIsotope_A ) break;
109 }
110 aNucleus.SetIsotope( target_isotope );
112 } else {
113 // the energy is inside the limits of the URR, part copied from G4ParticleHPChannel::ApplyYourself
114 if ( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ) {
115 if ( (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->GetWendtFissionGenerator() ) {
116 theFinalState = (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->GetWendtFissionGenerator()->ApplyYourself( aTrack, Z, A );
117 }
118 }
119 if ( ! theFinalState ) {
120 G4int icounter = 0;
121 G4int icounter_max = 1024;
122 while ( theFinalState == nullptr ) {
123 icounter++;
124 if ( icounter > icounter_max ) {
125 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
126 break;
127 }
128 // calls the final state for the found element and isotope
129 theFinalState = ((*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->GetFinalStates())[isotopeJ]->ApplyYourself( aTrack );
130 }
131 }
132 }
133 return theFinalState;
134}
135
136
138 particleHPfission->BuildPhysicsTable( *(G4Neutron::Neutron()) );
140 if ( URRlimits == nullptr ) {
144 }
145}
146
147
148const std::pair< G4double, G4double > G4ParticleHPFissionURR::GetFatalEnergyCheckLevels() const {
149 // max energy non-conservation is mass of heavy nucleus
150 return std::pair< G4double, G4double >( 10.0 * perCent, 350.0 * CLHEP::GeV );
151}
152
153
157
158
162
163
164void G4ParticleHPFissionURR::ModelDescription( std::ostream& outFile ) const {
165 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for fission reaction of neutrons in the unresolved resonance region.";
166}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
std::size_t GetIndex() const
Definition G4Element.hh:159
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
G4int GetZasInt() const
Definition G4Element.hh:120
const G4Material * GetMaterial() const
G4double GetKineticEnergy() const
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)
G4int GetN() const
Definition G4Isotope.hh:83
const G4Element * GetElement(G4int iel) const
std::size_t GetNumberOfElements() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
virtual void ModelDescription(std::ostream &outFile) const
void BuildPhysicsTable(const G4ParticleDefinition &)
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates() const
std::vector< std::pair< G4double, G4double > > * GetURRlimits() const
void RegisterURRlimits(std::vector< std::pair< G4double, G4double > > *val)
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
static G4ParticleHPProbabilityTablesStore * GetInstance()
std::vector< std::pair< G4double, G4double > > * GetURRlimits()