Geant4-11
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 ( std::vector<G4ParticleHPChannel*>::iterator
58 it = theElastic->begin() ; it != theElastic->end() ; it++ ) {
59 delete *it;
60 }
61 theElastic->clear();
62 }
63 }
64}
65
66
68{
69 return this->ApplyYourself(aTrack, aNucleus, 0);
70}
71
72
73//--------------------------------------------------------
74// New method added by L. Thulliez (CEA-Saclay) 2021/05/04
75//--------------------------------------------------------
77{
79 const G4Material * theMaterial = aTrack.GetMaterial();
80 G4int n = theMaterial->GetNumberOfElements();
81 G4int index = theMaterial->GetElement(0)->GetIndex();
82
83 if ( ! isFromTSL ) {
84 if ( n != 1 ) {
85 G4int i;
86 G4double* xSec = new G4double[n];
87 G4double sum=0;
88 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
89 G4double rWeight;
91 for ( i = 0; i < n; i++ ) {
92 index = theMaterial->GetElement(i)->GetIndex();
93 rWeight = NumAtomsPerVolume[i];
94 xSec[i] = ((*theElastic)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
95 theMaterial->GetElement(i),
96 theMaterial->GetTemperature()));
97 xSec[i] *= rWeight;
98 sum+=xSec[i];
99 }
100 G4double random = G4UniformRand();
101 G4double running = 0;
102 for ( i = 0; i < n; i++ ) {
103 running += xSec[i];
104 index = theMaterial->GetElement(i)->GetIndex();
105 if ( sum == 0 || random <= running/sum ) break;
106 }
107 delete [] xSec;
108 }
109 } else {
110 G4int i;
111 if ( n != 1 ) {
112 for ( i = 0; i < n; i++ ) {
113 if ( aNucleus.GetZ_asInt() == (G4int)(theMaterial->GetElement(i)->GetZ()) ) {
114 index = theMaterial->GetElement(i)->GetIndex();
115 }
116 }
117 }
118 }
119
120 G4HadFinalState* finalState = ((*theElastic)[index])->ApplyYourself(aTrack);
122
123 // Overwrite target parameters
125 const G4Element* target_element = (*G4Element::GetElementTable())[index];
126 const G4Isotope* target_isotope=nullptr;
127 G4int iele = target_element->GetNumberOfIsotopes();
128 for ( G4int j = 0 ; j != iele ; j++ ) {
129 target_isotope=target_element->GetIsotope( j );
130 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
131 }
132 aNucleus.SetIsotope( target_isotope );
133
135 return finalState;
136}
137
138
139const std::pair<G4double, G4double> G4ParticleHPElastic::GetFatalEnergyCheckLevels() const
140{
141 //return std::pair<G4double, G4double>(10*perCent,10*GeV);
142 return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
143}
144
145
147{
149}
150
151
153{
155}
156
157
159{
160
162
163 theElastic = hpmanager->GetElasticFinalStates();
164
166
167 if ( theElastic == nullptr ) theElastic = new std::vector<G4ParticleHPChannel*>;
168
169 if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
170
171 if ( theElastic->size() == G4Element::GetNumberOfElements() ) {
173 return;
174 }
175
177 if(!std::getenv("G4NEUTRONHPDATA"))
178 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
179 dirName = std::getenv("G4NEUTRONHPDATA");
180 G4String tString = "/Elastic";
181 dirName = dirName + tString;
182 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
183 theElastic->push_back( new G4ParticleHPChannel );
184 ((*theElastic)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
185 //while(!((*theElastic)[i])->Register(theFS)) ;
186 ((*theElastic)[i])->Register(theFS) ;
187 }
188 delete theFS;
190
191 }
193}
194
195
196void G4ParticleHPElastic::ModelDescription(std::ostream& outFile) const
197{
198 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for inelastic reaction of neutrons below 20MeV\n";
199}
@ isAlive
static constexpr double perCent
Definition: G4SIunits.hh:325
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double MeV
Definition: G4SIunits.hh:200
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:397
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
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:178
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
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 &)
std::vector< G4ParticleHPChannel * > * theElastic
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)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
G4bool IsMasterThread()
Definition: G4Threading.cc:124
void Init()
Definition: G4IonTable.cc:77
#define DBL_MAX
Definition: templates.hh:62