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

#include <G4ParticleHPFission.hh>

+ Inheritance diagram for G4ParticleHPFission:

Public Member Functions

 G4ParticleHPFission ()
 
 ~G4ParticleHPFission ()
 
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 G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
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 52 of file G4ParticleHPFission.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFission()

G4ParticleHPFission::G4ParticleHPFission ( )

Definition at line 42 of file G4ParticleHPFission.cc.

43 : G4HadronicInteraction("NeutronHPFission"), theFission(nullptr), numEle(0)
44 {
45 SetMinEnergy( 0.0 );
46 SetMaxEnergy( 20.*MeV );
47 }
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4ParticleHPFission()

G4ParticleHPFission::~G4ParticleHPFission ( )

Definition at line 49 of file G4ParticleHPFission.cc.

50 {
51 // Vector is shared, only master deletes it
52 // delete [] theFission;
54 if ( theFission != nullptr ) {
55 for ( auto it=theFission->cbegin(); it!=theFission->cend(); ++it ) {
56 delete *it;
57 }
58 theFission->clear();
59 }
60 }
61 }
G4bool IsMasterThread()
Definition: G4Threading.cc:124

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 63 of file G4ParticleHPFission.cc.

64 {
65
67 const G4Material * theMaterial = aTrack.GetMaterial();
68 G4int n = (G4int)theMaterial->GetNumberOfElements();
69 std::size_t index = theMaterial->GetElement(0)->GetIndex();
70 if(n!=1)
71 {
72 G4double* xSec = new G4double[n];
73 G4double sum=0;
74 G4int i;
75 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
76 G4double rWeight;
78 for (i=0; i<n; ++i)
79 {
80 index = theMaterial->GetElement(i)->GetIndex();
81 rWeight = NumAtomsPerVolume[i];
82 xSec[i] = ((*theFission)[index])
83 ->GetXsec(aThermalE.GetThermalEnergy(aTrack,
84 theMaterial->GetElement(i),
85 theMaterial->GetTemperature()));
86 xSec[i] *= rWeight;
87 sum+=xSec[i];
88 }
89 G4double random = G4UniformRand();
90 G4double running = 0;
91 for (i=0; i<n; ++i)
92 {
93 running += xSec[i];
94 index = theMaterial->GetElement(i)->GetIndex();
95 //if(random<=running/sum) break;
96 if( sum == 0 || random <= running/sum ) break;
97 }
98 delete [] xSec;
99 }
100 //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
101 G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack,-2);
102
103 //Overwrite target parameters
105 const G4Element* target_element = (*G4Element::GetElementTable())[index];
106 const G4Isotope* target_isotope=NULL;
107 G4int iele = (G4int)target_element->GetNumberOfIsotopes();
108 for ( G4int j = 0 ; j != iele ; ++j ) {
109 target_isotope=target_element->GetIsotope( j );
110 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
111 }
112 aNucleus.SetIsotope( target_isotope );
113
115 return result;
116 }
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
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
const G4Material * GetMaterial() const
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)

Referenced by ApplyYourself().

◆ BuildPhysicsTable()

void G4ParticleHPFission::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 134 of file G4ParticleHPFission.cc.

135{
136
138
139 theFission = hpmanager->GetFissionFinalStates();
140
142
143 if ( theFission == nullptr ) theFission = new std::vector<G4ParticleHPChannel*>;
144
145 if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
146
147 if ( theFission->size() == G4Element::GetNumberOfElements() ) {
149 return;
150 }
151
152 if ( !G4FindDataDir("G4NEUTRONHPDATA") )
153 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
154 dirName = G4FindDataDir("G4NEUTRONHPDATA");
155 G4String tString = "/Fission";
156 dirName = dirName + tString;
157
158 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; ++i ) {
159 theFission->push_back( new G4ParticleHPChannel );
160 if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
161 ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
162 ((*theFission)[i])->Register( new G4ParticleHPFissionFS );
163 }
164 }
165 hpmanager->RegisterFissionFinalStates( theFission );
166 }
168}
const char * G4FindDataDir(const char *)
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
void RegisterFissionFinalStates(std::vector< G4ParticleHPChannel * > *val)
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates()
void Register(T *inst)
Definition: G4AutoDelete.hh:65
void Init()
Definition: G4IonTable.cc:77

◆ GetFatalEnergyCheckLevels()

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

Reimplemented from G4HadronicInteraction.

Definition at line 118 of file G4ParticleHPFission.cc.

119{
120 // max energy non-conservation is mass of heavy nucleus
121 return std::pair<G4double, G4double>(10.0*perCent, 350.0*CLHEP::GeV);
122}

◆ GetVerboseLevel()

G4int G4ParticleHPFission::GetVerboseLevel ( ) const

Definition at line 124 of file G4ParticleHPFission.cc.

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 170 of file G4ParticleHPFission.cc.

171{
172 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF)\n"
173 << "for induced fission reaction of neutrons below 20MeV\n";
174}

◆ SetVerboseLevel()

void G4ParticleHPFission::SetVerboseLevel ( G4int  newValue)

Definition at line 129 of file G4ParticleHPFission.cc.

130{
132}

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