Geant4 10.7.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 ()
 
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 41 of file G4ParticleHPFission.cc.

42 :G4HadronicInteraction("NeutronHPFission")
43 ,theFission(NULL)
44 ,numEle(0)
45 {
46 SetMinEnergy( 0.0 );
47 SetMaxEnergy( 20.*MeV );
48/*
49 if(!std::getenv("G4NEUTRONHPDATA"))
50 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
51 dirName = std::getenv("G4NEUTRONHPDATA");
52 G4String tString = "/Fission";
53 dirName = dirName + tString;
54 numEle = G4Element::GetNumberOfElements();
55 //theFission = new G4ParticleHPChannel[numEle];
56
57 //for (G4int i=0; i<numEle; i++)
58 //{
59 //if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
60 // if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
61 // {
62 // theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
63 // theFission[i].Register(&theFS);
64 // }
65 //}
66
67 for ( G4int i = 0 ; i < numEle ; i++ )
68 {
69 theFission.push_back( new G4ParticleHPChannel );
70 if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
71 {
72 (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
73 (*theFission[i]).Register(&theFS);
74 }
75 }
76*/
77 }
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)

◆ ~G4ParticleHPFission()

G4ParticleHPFission::~G4ParticleHPFission ( )

Definition at line 79 of file G4ParticleHPFission.cc.

80 {
81 //Vector is shared, only master deletes it
82 //delete [] theFission;
84 if ( theFission != NULL ) {
85 for ( std::vector<G4ParticleHPChannel*>::iterator
86 it = theFission->begin() ; it != theFission->end() ; it++ ) {
87 delete *it;
88 }
89 theFission->clear();
90 }
91 }
92 }
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 95 of file G4ParticleHPFission.cc.

96 {
97
99 const G4Material * theMaterial = aTrack.GetMaterial();
100 G4int n = theMaterial->GetNumberOfElements();
101 G4int index = theMaterial->GetElement(0)->GetIndex();
102 if(n!=1)
103 {
104 G4double* xSec = new G4double[n];
105 G4double sum=0;
106 G4int i;
107 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
108 G4double rWeight;
109 G4ParticleHPThermalBoost aThermalE;
110 for (i=0; i<n; i++)
111 {
112 index = theMaterial->GetElement(i)->GetIndex();
113 rWeight = NumAtomsPerVolume[i];
114 xSec[i] = ((*theFission)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
115 theMaterial->GetElement(i),
116 theMaterial->GetTemperature()));
117 xSec[i] *= rWeight;
118 sum+=xSec[i];
119 }
120 G4double random = G4UniformRand();
121 G4double running = 0;
122 for (i=0; i<n; i++)
123 {
124 running += xSec[i];
125 index = theMaterial->GetElement(i)->GetIndex();
126 //if(random<=running/sum) break;
127 if( sum == 0 || random <= running/sum ) break;
128 }
129 delete [] xSec;
130 }
131 //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
132 G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack,-2);
133
134 //Overwrite target parameters
136 const G4Element* target_element = (*G4Element::GetElementTable())[index];
137 const G4Isotope* target_isotope=NULL;
138 G4int iele = target_element->GetNumberOfIsotopes();
139 for ( G4int j = 0 ; j != iele ; j++ ) {
140 target_isotope=target_element->GetIsotope( j );
141 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
142 }
143 //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
144 //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
145 //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
146 aNucleus.SetIsotope( target_isotope );
147
149 return result;
150 }
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetIndex() const
Definition: G4Element.hh:181
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
const G4Material * GetMaterial() const
G4int GetN() const
Definition: G4Isotope.hh:93
G4double GetTemperature() const
Definition: G4Material.hh:180
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
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 183 of file G4ParticleHPFission.cc.

184{
185
187
188 theFission = hpmanager->GetFissionFinalStates();
189
191
192 if ( theFission == NULL ) theFission = new std::vector<G4ParticleHPChannel*>;
193
194 if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
195
196 if ( theFission->size() == G4Element::GetNumberOfElements() ) {
198 return;
199 }
200
201 if ( !std::getenv("G4NEUTRONHPDATA") )
202 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
203 dirName = std::getenv("G4NEUTRONHPDATA");
204 G4String tString = "/Fission";
205 dirName = dirName + tString;
206
207 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
208 theFission->push_back( new G4ParticleHPChannel );
209 if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
210 ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
211 ((*theFission)[i])->Register( new G4ParticleHPFissionFS );
212 }
213 }
214
215 hpmanager->RegisterFissionFinalStates( theFission );
216
217 }
218
220}
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
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 152 of file G4ParticleHPFission.cc.

153{
154 // max energy non-conservation is mass of heavy nucleus
155 //return std::pair<G4double, G4double>(5*perCent,250*GeV);
156 return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
157}
#define DBL_MAX
Definition: templates.hh:62

◆ GetVerboseLevel()

G4int G4ParticleHPFission::GetVerboseLevel ( ) const

Definition at line 175 of file G4ParticleHPFission.cc.

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 222 of file G4ParticleHPFission.cc.

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

◆ SetVerboseLevel()

void G4ParticleHPFission::SetVerboseLevel ( G4int  newValue)

Definition at line 179 of file G4ParticleHPFission.cc.

180{
182}

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