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

#include <G4NeutronHPFission.hh>

+ Inheritance diagram for G4NeutronHPFission:

Public Member Functions

 G4NeutronHPFission ()
 
 ~G4NeutronHPFission ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
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)
 
const G4HadronicInteractionGetMyPointer () const
 
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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) 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
 

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 53 of file G4NeutronHPFission.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPFission()

G4NeutronHPFission::G4NeutronHPFission ( )

Definition at line 38 of file G4NeutronHPFission.cc.

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

◆ ~G4NeutronHPFission()

G4NeutronHPFission::~G4NeutronHPFission ( )

Definition at line 72 of file G4NeutronHPFission.cc.

73 {
74 //delete [] theFission;
75 for ( std::vector<G4NeutronHPChannel*>::iterator
76 it = theFission.begin() ; it != theFission.end() ; it++ )
77 {
78 delete *it;
79 }
80 theFission.clear();
81 }

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 84 of file G4NeutronHPFission.cc.

85 {
86
87 if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
88
90 const G4Material * theMaterial = aTrack.GetMaterial();
91 G4int n = theMaterial->GetNumberOfElements();
92 G4int index = theMaterial->GetElement(0)->GetIndex();
93 if(n!=1)
94 {
95 xSec = new G4double[n];
96 G4double sum=0;
97 G4int i;
98 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
99 G4double rWeight;
100 G4NeutronHPThermalBoost aThermalE;
101 for (i=0; i<n; i++)
102 {
103 index = theMaterial->GetElement(i)->GetIndex();
104 rWeight = NumAtomsPerVolume[i];
105 xSec[i] = (*theFission[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
106 theMaterial->GetElement(i),
107 theMaterial->GetTemperature()));
108 xSec[i] *= rWeight;
109 sum+=xSec[i];
110 }
111 G4double random = G4UniformRand();
112 G4double running = 0;
113 for (i=0; i<n; i++)
114 {
115 running += xSec[i];
116 index = theMaterial->GetElement(i)->GetIndex();
117 //if(random<=running/sum) break;
118 if( sum == 0 || random <= running/sum ) break;
119 }
120 delete [] xSec;
121 }
122 //return theFission[index].ApplyYourself(aTrack);
123 G4HadFinalState* result = (*theFission[index]).ApplyYourself(aTrack);
126 return result;
127 }
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
size_t GetIndex() const
Definition: G4Element.hh:182
const G4Material * GetMaterial() const
G4double GetTemperature() const
Definition: G4Material.hh:181
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:201
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)

Referenced by ApplyYourself().

◆ GetFatalEnergyCheckLevels()

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

Reimplemented from G4HadronicInteraction.

Definition at line 129 of file G4NeutronHPFission.cc.

130{
131 // max energy non-conservation is mass of heavy nucleus
132 //return std::pair<G4double, G4double>(5*perCent,250*GeV);
133 return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
134}
#define DBL_MAX
Definition: templates.hh:83

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