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
G4ParticleHPElastic.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// 080319 Compilation warnings - gcc-4.3.0 fix by T. Koi
32//
33// P. Arce, June-2014 Conversion neutron_hp to particle_hp
34//
36#include "G4SystemOfUnits.hh"
39#include "G4Threading.hh"
41
42
44 : G4HadronicInteraction("NeutronHPElastic"), theElastic(nullptr), numEle(0)
45{
46 overrideSuspension = false;
47 SetMinEnergy(0.*eV);
48 SetMaxEnergy(20.*MeV);
49}
50
51
53{
54 //the vectror is shared among threads, only master deletes
56 if ( theElastic != nullptr ) {
57 for ( auto it=theElastic->cbegin(); it!=theElastic->cend(); ++it ) {
58 delete *it;
59 }
60 theElastic->clear();
61 }
62 }
63}
64
65
67{
68 return this->ApplyYourself(aTrack, aNucleus, 0);
69}
70
71
72//--------------------------------------------------------
73// New method added by L. Thulliez (CEA-Saclay) 2021/05/04
74//--------------------------------------------------------
76{
78 const G4Material * theMaterial = aTrack.GetMaterial();
79 G4int n = (G4int)theMaterial->GetNumberOfElements();
80 std::size_t index = theMaterial->GetElement(0)->GetIndex();
81
82 if ( ! isFromTSL ) {
83 if ( n != 1 ) {
84 G4int i;
85 G4double* xSec = new G4double[n];
86 G4double sum=0;
87 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
88 G4double rWeight;
90 for ( i = 0; i < n; ++i ) {
91 index = theMaterial->GetElement(i)->GetIndex();
92 rWeight = NumAtomsPerVolume[i];
93 xSec[i] = ((*theElastic)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
94 theMaterial->GetElement(i),
95 theMaterial->GetTemperature()));
96 xSec[i] *= rWeight;
97 sum+=xSec[i];
98 }
99 G4double random = G4UniformRand();
100 G4double running = 0;
101 for ( i = 0; i < n; ++i ) {
102 running += xSec[i];
103 index = theMaterial->GetElement(i)->GetIndex();
104 if ( sum == 0 || random <= running/sum ) break;
105 }
106 delete [] xSec;
107 }
108 } else {
109 G4int i;
110 if ( n != 1 ) {
111 for ( i = 0; i < n; ++i ) {
112 if ( aNucleus.GetZ_asInt() == (G4int)(theMaterial->GetElement(i)->GetZ()) ) {
113 index = theMaterial->GetElement(i)->GetIndex();
114 }
115 }
116 }
117 }
118
119 G4HadFinalState* finalState = ((*theElastic)[index])->ApplyYourself(aTrack);
120 if (overrideSuspension) finalState->SetStatusChange(isAlive);
121
122 // Overwrite target parameters
124 const G4Element* target_element = (*G4Element::GetElementTable())[index];
125 const G4Isotope* target_isotope=nullptr;
126 G4int iele = (G4int)target_element->GetNumberOfIsotopes();
127 for ( G4int j = 0 ; j != iele ; ++j ) {
128 target_isotope=target_element->GetIsotope( j );
129 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
130 }
131 aNucleus.SetIsotope( target_isotope );
132
134 return finalState;
135}
136
137
138const std::pair<G4double, G4double> G4ParticleHPElastic::GetFatalEnergyCheckLevels() const
139{
140 // max energy non-conservation is mass of heavy nucleus
141 return std::pair<G4double, G4double>(10.0*perCent, 350.0*CLHEP::GeV);
142}
143
144
146{
148}
149
150
152{
154}
155
156
158{
159
161
162 theElastic = hpmanager->GetElasticFinalStates();
163
165
166 if ( theElastic == nullptr ) theElastic = new std::vector<G4ParticleHPChannel*>;
167
168 if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
169
170 if ( theElastic->size() == G4Element::GetNumberOfElements() ) {
172 return;
173 }
174
176 if(!G4FindDataDir("G4NEUTRONHPDATA"))
177 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
178 dirName = G4FindDataDir("G4NEUTRONHPDATA");
179 G4String tString = "/Elastic";
180 dirName = dirName + tString;
181 for ( G4int i = numEle; i < (G4int)G4Element::GetNumberOfElements(); ++i ) {
182 theElastic->push_back( new G4ParticleHPChannel );
183 ((*theElastic)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
184 //while(!((*theElastic)[i])->Register(theFS)) ;
185 ((*theElastic)[i])->Register(theFS) ;
186 }
187 delete theFS;
188 hpmanager->RegisterElasticFinalStates( theElastic );
189
190 }
192}
193
194
195void G4ParticleHPElastic::ModelDescription(std::ostream& outFile) const
196{
197 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for inelastic reaction of neutrons below 20MeV\n";
198}
const char * G4FindDataDir(const char *)
@ isAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4double GetZ() const
Definition: G4Element.hh:131
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
void SetStatusChange(G4HadFinalStateStatus aS)
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
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
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
void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void ModelDescription(std::ostream &outFile) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
void RegisterElasticFinalStates(std::vector< G4ParticleHPChannel * > *val)
static G4ParticleHPManager * GetInstance()
std::vector< G4ParticleHPChannel * > * GetElasticFinalStates()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
G4bool IsMasterThread()
Definition: G4Threading.cc:124