Geant4 11.4.0
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 () override
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus) override
const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const override
G4int GetVerboseLevel () const
void SetVerboseLevel (G4int)
void BuildPhysicsTable (const G4ParticleDefinition &) override
void ModelDescription (std::ostream &outFile) const override
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 52 of file G4ParticleHPFission.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFission()

G4ParticleHPFission::G4ParticleHPFission ( )

Definition at line 43 of file G4ParticleHPFission.cc.

43 : G4HadronicInteraction("NeutronHPFission")
44{
45 SetMinEnergy(0.0);
46 SetMaxEnergy(20. * MeV);
47}
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4ParticleHPFission()

G4ParticleHPFission::~G4ParticleHPFission ( )
override

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()

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 63 of file G4ParticleHPFission.cc.

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

Referenced by ApplyYourself().

◆ BuildPhysicsTable()

void G4ParticleHPFission::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 133 of file G4ParticleHPFission.cc.

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

◆ GetFatalEnergyCheckLevels()

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

Reimplemented from G4HadronicInteraction.

Definition at line 117 of file G4ParticleHPFission.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPFission::GetVerboseLevel ( ) const

Definition at line 123 of file G4ParticleHPFission.cc.

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 169 of file G4ParticleHPFission.cc.

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

◆ SetVerboseLevel()

void G4ParticleHPFission::SetVerboseLevel ( G4int newValue)

Definition at line 128 of file G4ParticleHPFission.cc.


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