Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPIsoProbabilityTable Class Reference

#include <G4ParticleHPIsoProbabilityTable.hh>

+ Inheritance diagram for G4ParticleHPIsoProbabilityTable:

Public Member Functions

 G4ParticleHPIsoProbabilityTable ()=default
 
virtual ~G4ParticleHPIsoProbabilityTable ()
 
virtual void Init (G4int, G4int, G4int, G4double, const G4String &)
 
virtual G4double GetCorrelatedIsoCrossSectionPT (const G4DynamicParticle *, G4int, const G4Element *, G4double &, G4double &, std::thread::id &)
 
virtual G4double GetIsoCrossSectionPT (const G4DynamicParticle *, G4int, const G4Element *, G4double &, std::map< std::thread::id, G4double > &, std::thread::id &)
 

Protected Member Functions

G4double GetDopplerBroadenedElasticXS (const G4DynamicParticle *, G4int, G4int)
 
G4double GetDopplerBroadenedCaptureXS (const G4DynamicParticle *, G4int, G4int)
 
G4double GetDopplerBroadenedFissionXS (const G4DynamicParticle *, G4int, G4int)
 
G4double GetDopplerBroadenedInelasticXS (const G4DynamicParticle *, G4int, G4int)
 

Protected Attributes

G4int Z = 0
 
G4int A = 0
 
G4int m = -1
 
G4double T = -1.
 
G4double Emin = DBL_MAX
 
G4double Emax = 0.
 
G4int nEnergies = 0
 
std::map< std::thread::id, G4doubleenergy_cache
 
std::map< std::thread::id, G4doublexsela_cache
 
std::map< std::thread::id, G4doublexscap_cache
 
std::map< std::thread::id, G4doublexsfiss_cache
 
G4ParticleHPVectortheEnergies = nullptr
 
std::vector< std::vector< G4double > * > * theProbabilities = nullptr
 
std::vector< std::vector< G4double > * > * theElasticData = nullptr
 
std::vector< std::vector< G4double > * > * theCaptureData = nullptr
 
std::vector< std::vector< G4double > * > * theFissionData = nullptr
 
std::vector< std::vector< G4double > * > * theInelasticData = nullptr
 
G4String filename
 

Detailed Description

Definition at line 61 of file G4ParticleHPIsoProbabilityTable.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPIsoProbabilityTable()

G4ParticleHPIsoProbabilityTable::G4ParticleHPIsoProbabilityTable ( )
default

◆ ~G4ParticleHPIsoProbabilityTable()

G4ParticleHPIsoProbabilityTable::~G4ParticleHPIsoProbabilityTable ( )
virtual

Definition at line 68 of file G4ParticleHPIsoProbabilityTable.cc.

69{
70 for ( auto it = theProbabilities->cbegin(); it != theProbabilities->cend(); ++it ) {
71 delete* it;
72 }
73 for ( auto it = theElasticData->cbegin(); it != theElasticData->cend(); ++it ) {
74 delete* it;
75 }
76 for ( auto it = theCaptureData->cbegin(); it != theCaptureData->cend(); ++it ) {
77 delete* it;
78 }
79 for ( auto it = theFissionData->cbegin(); it != theFissionData->cend(); ++it ) {
80 delete* it;
81 }
82 delete theEnergies;
83 delete theProbabilities;
84 delete theElasticData;
85 delete theCaptureData;
86 delete theFissionData;
87}
std::vector< std::vector< G4double > * > * theElasticData
std::vector< std::vector< G4double > * > * theFissionData
std::vector< std::vector< G4double > * > * theCaptureData
std::vector< std::vector< G4double > * > * theProbabilities

Member Function Documentation

◆ GetCorrelatedIsoCrossSectionPT()

G4double G4ParticleHPIsoProbabilityTable::GetCorrelatedIsoCrossSectionPT ( const G4DynamicParticle * ,
G4int ,
const G4Element * ,
G4double & ,
G4double & ,
std::thread::id &  )
virtual

Reimplemented in G4ParticleHPIsoProbabilityTable_CALENDF, and G4ParticleHPIsoProbabilityTable_NJOY.

Definition at line 93 of file G4ParticleHPIsoProbabilityTable.cc.

95{
96 G4Exception( "G4ParticleHPIsoProbabilityTable::GetCorrelatedIsoCrossSectionPT", "hadhp01",
97 FatalException, "The base method G4ParticleHPIsoProbabilityTable::GetCorrelatedIsoCrossSectionPT"
98 "is trying to be accessed, which is not allowed,"
99 "the derived class for NJOY or CALENDF was not properly declared." );
100 return 0.0;
101}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)

◆ GetDopplerBroadenedCaptureXS()

G4double G4ParticleHPIsoProbabilityTable::GetDopplerBroadenedCaptureXS ( const G4DynamicParticle * dP,
G4int indexEl,
G4int isotopeJ )
protected

Definition at line 150 of file G4ParticleHPIsoProbabilityTable.cc.

150 {
151 G4double result = 0.0;
152 G4ReactionProduct theNeutron( dP->GetDefinition() );
153 theNeutron.SetMomentum( dP->GetMomentum() );
154 theNeutron.SetKineticEnergy( dP->GetKineticEnergy() );
155 // prepare thermal nucleus
156 G4Nucleus aNuc;
158 G4ReactionProduct boosted;
159 G4double aXsection;
160 G4int counter = 0;
161 G4double buffer = 0.0;
162 G4int size = G4int( std::max( 10.0, T / 60.0 * kelvin));
163 G4ThreeVector neutronVelocity = 1.0 / G4Neutron::Neutron()->GetPDGMass() * theNeutron.GetMomentum();
164 G4double neutronVMag = neutronVelocity.mag();
165 while ( counter == 0 || std::abs( buffer - result / std::max( 1, counter ) ) > 0.03 * buffer ) {
166 if ( counter ) buffer = result / counter;
167 while ( counter < size ) {
168 ++counter;
169 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass, T );
170 boosted.Lorentz( theNeutron, aThermalNuc );
171 G4double theEkin = boosted.GetKineticEnergy();
172 aXsection = (*G4ParticleHPManager::GetInstance()->GetCaptureFinalStates())[indexEl]->GetWeightedXsec( theEkin, isotopeJ );
173 // velocity correction
174 G4ThreeVector targetVelocity = 1.0 / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
175 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
176 result += aXsection;
177 }
178 size += size;
179 }
180 result /= counter;
181 return result;
182}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition G4Nucleus.cc:248
std::vector< G4ParticleHPChannel * > * GetCaptureFinalStates() const
static G4ParticleHPManager * GetInstance()
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const

Referenced by G4ParticleHPIsoProbabilityTable_CALENDF::GetCorrelatedIsoCrossSectionPT(), G4ParticleHPIsoProbabilityTable_NJOY::GetCorrelatedIsoCrossSectionPT(), G4ParticleHPIsoProbabilityTable_CALENDF::GetIsoCrossSectionPT(), and G4ParticleHPIsoProbabilityTable_NJOY::GetIsoCrossSectionPT().

◆ GetDopplerBroadenedElasticXS()

G4double G4ParticleHPIsoProbabilityTable::GetDopplerBroadenedElasticXS ( const G4DynamicParticle * dP,
G4int indexEl,
G4int isotopeJ )
protected

Definition at line 115 of file G4ParticleHPIsoProbabilityTable.cc.

115 {
116 G4double result = 0.0;
117 G4ReactionProduct theNeutron( dP->GetDefinition() );
118 theNeutron.SetMomentum( dP->GetMomentum() );
119 theNeutron.SetKineticEnergy( dP->GetKineticEnergy() );
120 // prepare thermal nucleus
121 G4Nucleus aNuc;
123 G4ReactionProduct boosted;
124 G4double aXsection;
125 G4int counter = 0;
126 G4double buffer = 0.0;
127 G4int size = G4int( std::max( 10.0, T / 60.0 * kelvin ) );
128 G4ThreeVector neutronVelocity = 1.0 / G4Neutron::Neutron()->GetPDGMass() * theNeutron.GetMomentum();
129 G4double neutronVMag = neutronVelocity.mag();
130 while ( counter == 0 || std::abs( buffer - result / std::max(1, counter) ) > 0.03 * buffer ) {
131 if ( counter ) buffer = result / counter;
132 while ( counter < size ) {
133 ++counter;
134 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass, T );
135 boosted.Lorentz( theNeutron, aThermalNuc );
136 G4double theEkin = boosted.GetKineticEnergy();
137 aXsection = (*G4ParticleHPManager::GetInstance()->GetElasticFinalStates())[indexEl]->GetWeightedXsec( theEkin, isotopeJ );
138 // velocity correction
139 G4ThreeVector targetVelocity = 1.0 / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
140 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
141 result += aXsection;
142 }
143 size += size;
144 }
145 result /= counter;
146 return result;
147}
std::vector< G4ParticleHPChannel * > * GetElasticFinalStates() const

Referenced by G4ParticleHPIsoProbabilityTable_CALENDF::GetCorrelatedIsoCrossSectionPT(), G4ParticleHPIsoProbabilityTable_NJOY::GetCorrelatedIsoCrossSectionPT(), G4ParticleHPIsoProbabilityTable_CALENDF::GetIsoCrossSectionPT(), and G4ParticleHPIsoProbabilityTable_NJOY::GetIsoCrossSectionPT().

◆ GetDopplerBroadenedFissionXS()

G4double G4ParticleHPIsoProbabilityTable::GetDopplerBroadenedFissionXS ( const G4DynamicParticle * dP,
G4int indexEl,
G4int isotopeJ )
protected

Definition at line 185 of file G4ParticleHPIsoProbabilityTable.cc.

185 {
186 G4double result = 0.0;
187 G4ReactionProduct theNeutron( dP->GetDefinition() );
188 theNeutron.SetMomentum( dP->GetMomentum() );
189 theNeutron.SetKineticEnergy( dP->GetKineticEnergy() );
190 // prepare thermal nucleus
191 G4Nucleus aNuc;
192 G4double eleMass;
194 G4ReactionProduct boosted;
195 G4double aXsection;
196 G4int counter = 0;
197 G4double buffer = 0.0;
198 G4int size = G4int( std::max( 10.0, T / 60.0 * kelvin ) );
199 G4ThreeVector neutronVelocity = 1.0 / G4Neutron::Neutron()->GetPDGMass() * theNeutron.GetMomentum();
200 G4double neutronVMag = neutronVelocity.mag();
201 while ( counter == 0 || std::abs( buffer - result / std::max( 1, counter ) ) > 0.01 * buffer ) {
202 if ( counter ) buffer = result / counter;
203 while ( counter < size ) {
204 ++counter;
205 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass, T );
206 boosted.Lorentz( theNeutron, aThermalNuc );
207 G4double theEkin = boosted.GetKineticEnergy();
208 aXsection = (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[indexEl]->GetWeightedXsec( theEkin, isotopeJ );
209 // velocity correction
210 G4ThreeVector targetVelocity = 1.0 / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
211 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
212 result += aXsection;
213 }
214 size += size;
215 }
216 result /= counter;
217 return result;
218}
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates() const

Referenced by G4ParticleHPIsoProbabilityTable_CALENDF::GetCorrelatedIsoCrossSectionPT(), G4ParticleHPIsoProbabilityTable_NJOY::GetCorrelatedIsoCrossSectionPT(), G4ParticleHPIsoProbabilityTable_CALENDF::GetIsoCrossSectionPT(), and G4ParticleHPIsoProbabilityTable_NJOY::GetIsoCrossSectionPT().

◆ GetDopplerBroadenedInelasticXS()

G4double G4ParticleHPIsoProbabilityTable::GetDopplerBroadenedInelasticXS ( const G4DynamicParticle * dP,
G4int indexEl,
G4int isotopeJ )
protected

Definition at line 221 of file G4ParticleHPIsoProbabilityTable.cc.

221 {
222 G4double result = 0.0;
223 G4ReactionProduct theNeutron( dP->GetDefinition() );
224 theNeutron.SetMomentum( dP->GetMomentum() );
225 theNeutron.SetKineticEnergy( dP->GetKineticEnergy() );
226 // prepare thermal nucleus
227 G4Nucleus aNuc;
229 G4ReactionProduct boosted;
230 G4double aXsection;
231 G4int counter = 0;
232 G4int failCount = 0;
233 G4double buffer = 0.0;
234 G4int size = G4int( std::max( 10.0, T / 60.0 * kelvin ) );
235 G4ThreeVector neutronVelocity = 1.0 / G4Neutron::Neutron()->GetPDGMass() * theNeutron.GetMomentum();
236 G4double neutronVMag = neutronVelocity.mag();
237 #ifndef G4PHP_DOPPLER_LOOP_ONCE
238 while ( counter == 0 || std::abs( buffer - result / std::max( 1, counter ) ) > 0.01 * buffer ) {
239 if ( counter ) buffer = result / counter;
240 while ( counter < size ) {
241 ++counter;
242 #endif
243 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass / G4Neutron::Neutron()->GetPDGMass(), T );
244 boosted.Lorentz( theNeutron, aThermalNuc );
245 G4double theEkin = boosted.GetKineticEnergy();
247 ->GetWeightedXsec( theEkin, isotopeJ ) ) * barn;
248 if ( aXsection < 0.0 ) {
249 if ( failCount < 1000 ) {
250 ++failCount;
251 #ifndef G4PHP_DOPPLER_LOOP_ONCE
252 --counter;
253 continue;
254 #endif
255 } else {
256 aXsection = 0.0;
257 }
258 }
259 // velocity correction
260 G4ThreeVector targetVelocity = 1.0 / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
261 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
262 result += aXsection;
263 }
264#ifndef G4PHP_DOPPLER_LOOP_ONCE
265 size += size;
266 }
267 result /= counter;
268#endif
269 return result;
270}
std::vector< G4ParticleHPChannelList * > * GetInelasticFinalStates(const G4ParticleDefinition *part) const

Referenced by G4ParticleHPIsoProbabilityTable_CALENDF::GetCorrelatedIsoCrossSectionPT(), and G4ParticleHPIsoProbabilityTable_CALENDF::GetIsoCrossSectionPT().

◆ GetIsoCrossSectionPT()

G4double G4ParticleHPIsoProbabilityTable::GetIsoCrossSectionPT ( const G4DynamicParticle * ,
G4int ,
const G4Element * ,
G4double & ,
std::map< std::thread::id, G4double > & ,
std::thread::id &  )
virtual

Reimplemented in G4ParticleHPIsoProbabilityTable_CALENDF, and G4ParticleHPIsoProbabilityTable_NJOY.

Definition at line 104 of file G4ParticleHPIsoProbabilityTable.cc.

106{
107 G4Exception( "G4ParticleHPIsoProbabilityTable::GetIsoCrossSectionPT", "hadhp01",
108 FatalException, "The base method G4ParticleHPIsoProbabilityTable::GetIsoCrossSectionPT"
109 "is trying to be accessed, which is not allowed,"
110 "the derived class for NJOY or CALENDF was not properly declared." );
111 return 0.0;
112}

◆ Init()

void G4ParticleHPIsoProbabilityTable::Init ( G4int ,
G4int ,
G4int ,
G4double ,
const G4String &  )
virtual

Member Data Documentation

◆ A

◆ Emax

◆ Emin

G4double G4ParticleHPIsoProbabilityTable::Emin = DBL_MAX
protected

◆ energy_cache

◆ filename

G4String G4ParticleHPIsoProbabilityTable::filename
protected

◆ m

G4int G4ParticleHPIsoProbabilityTable::m = -1
protected

◆ nEnergies

G4int G4ParticleHPIsoProbabilityTable::nEnergies = 0
protected

◆ T

◆ theCaptureData

◆ theElasticData

◆ theEnergies

◆ theFissionData

◆ theInelasticData

std::vector< std::vector< G4double >* >* G4ParticleHPIsoProbabilityTable::theInelasticData = nullptr
protected

Definition at line 99 of file G4ParticleHPIsoProbabilityTable.hh.

◆ theProbabilities

◆ xscap_cache

◆ xsela_cache

◆ xsfiss_cache

◆ Z


The documentation for this class was generated from the following files: