Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleHPInelasticURR Class Reference

#include <G4ParticleHPInelasticURR.hh>

+ Inheritance diagram for G4ParticleHPInelasticURR:

Public Member Functions

 G4ParticleHPInelasticURR ()
 
 ~G4ParticleHPInelasticURR ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 61 of file G4ParticleHPInelasticURR.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPInelasticURR()

G4ParticleHPInelasticURR::G4ParticleHPInelasticURR ( )

Definition at line 58 of file G4ParticleHPInelasticURR.cc.

58 : G4HadronicInteraction( "NeutronHPInelasticURR" ) {
59 SetMinEnergy( 0.0 * CLHEP::eV );
60 SetMaxEnergy( 20.0 * CLHEP::MeV );
61 particleHPinelastic = new G4ParticleHPInelastic( G4Neutron::Neutron(), "NeutronHPInelastic" );
62}
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101

◆ ~G4ParticleHPInelasticURR()

G4ParticleHPInelasticURR::~G4ParticleHPInelasticURR ( )

Definition at line 65 of file G4ParticleHPInelasticURR.cc.

65{}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ParticleHPInelasticURR::ApplyYourself ( const G4HadProjectile & aTrack,
G4Nucleus & aTargetNucleus )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 68 of file G4ParticleHPInelasticURR.cc.

68 {
69 if ( doNOTusePTforInelastic ) {
70 return particleHPinelastic->ApplyYourself( aTrack, aNucleus );
71 }
72 const G4Material* theMaterial = aTrack.GetMaterial();
73 G4double kineticEnergy = aTrack.GetKineticEnergy();
74 G4HadFinalState* theFinalState = nullptr;
75 if ( kineticEnergy < (*URRlimits).back().first || kineticEnergy > (*URRlimits).back().second ) {
76 return particleHPinelastic->ApplyYourself( aTrack, aNucleus );
77 }
78 G4int elementI = -1;
79 G4int isotopeJ = -1;
80 G4int A = aNucleus.GetA_asInt();
81 G4int Z = aNucleus.GetZ_asInt();
83 // finds the element and isotope of the selected target aNucleus
84 for ( G4int i = 0; i < (G4int)theMaterial->GetNumberOfElements(); ++i ) {
85 if ( Z == theMaterial->GetElement(i)->GetZasInt() ) {
86 for ( G4int j = 0; j < (G4int)theMaterial->GetElement(i)->GetNumberOfIsotopes(); ++j ) {
87 if ( A == theMaterial->GetElement(i)->GetIsotope(j)->GetN() ) {
88 isotopeJ = j;
89 break;
90 }
91 }
92 // the loop cannot be ended here because the material can have two elements with same Z but different isotopic composition
93 if ( isotopeJ != -1 ) {
94 // isotope was found and for loop is ended
95 elementI = (G4int)theMaterial->GetElement(i)->GetIndex();
96 break;
97 }
98 } // end if find element
99 } // end element loop
100 // Check whether the energy is out of the URR limits for the given element
101 if ( kineticEnergy < (*URRlimits).at(elementI).first || kineticEnergy > (*URRlimits).at(elementI).second ) {
102 // Call inelastic final state in G4ParicleHPChannel and SELECT ISOTOPE (to be improved in the future)
103 const G4Element* target_element = (*G4Element::GetElementTable())[elementI];
104 theFinalState = (*G4ParticleHPManager::GetInstance()->GetInelasticFinalStates( aTrack.GetDefinition() ))[elementI]
105 ->ApplyYourself( target_element, aTrack );
106 // Update target nucleus information according to the selected isotope
108 aNucleus.SetParameters( selectedIsotope_A, Z );
109 const G4Isotope* target_isotope = nullptr;
110 // Find the selected isotope among in the element
111 for ( G4int j = 0; j < (G4int)target_element->GetNumberOfIsotopes(); ++j ) {
112 target_isotope = target_element->GetIsotope(j);
113 if ( target_isotope->GetN() == selectedIsotope_A ) break;
114 }
115 aNucleus.SetIsotope( target_isotope );
116 } else {
117 // the energy is inside the limits of the URR and the isotope has to be found, calls the final state for the found element and isotope
118 theFinalState = (*G4ParticleHPManager::GetInstance()->GetInelasticFinalStates( aTrack.GetDefinition() ))[elementI]
119 ->ApplyYourself( isotopeJ, Z, A, aTrack );
120 }
122 return theFinalState;
123}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
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
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetN() const
Definition G4Isotope.hh:83
const G4Element * GetElement(G4int iel) const
std::size_t GetNumberOfElements() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
std::vector< G4ParticleHPChannelList * > * GetInelasticFinalStates(const G4ParticleDefinition *part) const
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()

Referenced by ApplyYourself().

◆ BuildPhysicsTable()

void G4ParticleHPInelasticURR::BuildPhysicsTable ( const G4ParticleDefinition & )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 126 of file G4ParticleHPInelasticURR.cc.

126 {
127 particleHPinelastic->BuildPhysicsTable( *(G4Neutron::Neutron()) );
128 if ( G4HadronicParameters::Instance()->GetTypeTablePT() == "njoy" ) {
129 doNOTusePTforInelastic = true;
130 } else if ( G4HadronicParameters::Instance()->GetTypeTablePT() == "calendf" ) {
131 doNOTusePTforInelastic = false;
132 // in the case of calendf probability tables, it sets the limits of the URR
134 if ( URRlimits == nullptr ) {
138 }
139 }
140}
static G4HadronicParameters * Instance()
std::vector< std::pair< G4double, G4double > > * GetURRlimits() const
void RegisterURRlimits(std::vector< std::pair< G4double, G4double > > *val)
static G4ParticleHPProbabilityTablesStore * GetInstance()
std::vector< std::pair< G4double, G4double > > * GetURRlimits()

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4ParticleHPInelasticURR::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 143 of file G4ParticleHPInelasticURR.cc.

143 {
144 // max energy non-conservation is mass of heavy nucleus
145 return std::pair< G4double, G4double >( 10.0 * perCent, 350.0 * CLHEP::GeV );
146}

◆ GetVerboseLevel()

G4int G4ParticleHPInelasticURR::GetVerboseLevel ( ) const

Definition at line 149 of file G4ParticleHPInelasticURR.cc.

149 {
151}

◆ ModelDescription()

void G4ParticleHPInelasticURR::ModelDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 159 of file G4ParticleHPInelasticURR.cc.

159 {
160 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for Inelastic reaction of neutrons in the unresolved resonance region.";
161}

◆ SetVerboseLevel()

void G4ParticleHPInelasticURR::SetVerboseLevel ( G4int newValue)

Definition at line 154 of file G4ParticleHPInelasticURR.cc.

154 {
156}

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