Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleHPFission Class Reference

#include <G4ParticleHPFission.hh>

Inheritance diagram for G4ParticleHPFission:
Collaboration diagram for G4ParticleHPFission:

Public Member Functions

 G4ParticleHPFission ()
 
 ~G4ParticleHPFission ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () 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,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void InitialiseModel ()
 

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

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(!getenv("G4NEUTRONHPDATA"))
50  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
51  dirName = 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)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)
static constexpr double MeV
Definition: G4SIunits.hh:214

Here is the call graph for this function:

G4ParticleHPFission::~G4ParticleHPFission ( )

Definition at line 79 of file G4ParticleHPFission.cc.

80  {
81  //delete [] theFission;
82  if ( theFission != NULL ) {
83  for ( std::vector<G4ParticleHPChannel*>::iterator
84  it = theFission->begin() ; it != theFission->end() ; it++ ) {
85  delete *it;
86  }
87  theFission->clear();
88  }
89  }

Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 92 of file G4ParticleHPFission.cc.

93  {
94 
96  const G4Material * theMaterial = aTrack.GetMaterial();
97  G4int n = theMaterial->GetNumberOfElements();
98  G4int index = theMaterial->GetElement(0)->GetIndex();
99  if(n!=1)
100  {
101  G4double* xSec = new G4double[n];
102  G4double sum=0;
103  G4int i;
104  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
105  G4double rWeight;
106  G4ParticleHPThermalBoost aThermalE;
107  for (i=0; i<n; i++)
108  {
109  index = theMaterial->GetElement(i)->GetIndex();
110  rWeight = NumAtomsPerVolume[i];
111  xSec[i] = ((*theFission)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
112  theMaterial->GetElement(i),
113  theMaterial->GetTemperature()));
114  xSec[i] *= rWeight;
115  sum+=xSec[i];
116  }
117  G4double random = G4UniformRand();
118  G4double running = 0;
119  for (i=0; i<n; i++)
120  {
121  running += xSec[i];
122  index = theMaterial->GetElement(i)->GetIndex();
123  //if(random<=running/sum) break;
124  if( sum == 0 || random <= running/sum ) break;
125  }
126  delete [] xSec;
127  }
128  //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
129  G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack,-2);
130 
131  //Overwrite target parameters
132  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
133  const G4Element* target_element = (*G4Element::GetElementTable())[index];
134  const G4Isotope* target_isotope=NULL;
135  G4int iele = target_element->GetNumberOfIsotopes();
136  for ( G4int j = 0 ; j != iele ; j++ ) {
137  target_isotope=target_element->GetIsotope( j );
138  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
139  }
140  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
141  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
142  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
143  aNucleus.SetIsotope( target_isotope );
144 
146  return result;
147  }
G4double G4ParticleHPJENDLHEData::G4double result
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
size_t GetIndex() const
Definition: G4Element.hh:182
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
G4double GetTemperature() const
Definition: G4Material.hh:182
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()

Here is the call graph for this function:

void G4ParticleHPFission::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 180 of file G4ParticleHPFission.cc.

181 {
182 
184 
185  theFission = hpmanager->GetFissionFinalStates();
186 
187  if ( G4Threading::IsMasterThread() ) {
188 
189  if ( theFission == NULL ) theFission = new std::vector<G4ParticleHPChannel*>;
190 
191  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
192 
193  if ( theFission->size() == G4Element::GetNumberOfElements() ) {
195  return;
196  }
197 
198  if ( !getenv("G4NEUTRONHPDATA") )
199  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
200  dirName = getenv("G4NEUTRONHPDATA");
201  G4String tString = "/Fission";
202  dirName = dirName + tString;
203 
204  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
205  theFission->push_back( new G4ParticleHPChannel );
206  if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
207  ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
208  ((*theFission)[i])->Register( new G4ParticleHPFissionFS );
209  }
210  }
211 
212  hpmanager->RegisterFissionFinalStates( theFission );
213 
214  }
215 
217 }
static G4ParticleHPManager * GetInstance()
void Init()
Definition: G4IonTable.cc:90
void RegisterFissionFinalStates(std::vector< G4ParticleHPChannel * > *val)
int G4int
Definition: G4Types.hh:78
void Register(T *inst)
Definition: G4AutoDelete.hh:65
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates()
G4bool IsMasterThread()
Definition: G4Threading.cc:146
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398

Here is the call graph for this function:

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

Reimplemented from G4HadronicInteraction.

Definition at line 149 of file G4ParticleHPFission.cc.

150 {
151  // max energy non-conservation is mass of heavy nucleus
152  //return std::pair<G4double, G4double>(5*perCent,250*GeV);
153  return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
154 }
static constexpr double perCent
Definition: G4SIunits.hh:332
#define DBL_MAX
Definition: templates.hh:83
G4int G4ParticleHPFission::GetVerboseLevel ( ) const

Definition at line 172 of file G4ParticleHPFission.cc.

173 {
175 }
static G4ParticleHPManager * GetInstance()

Here is the call graph for this function:

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

Reimplemented from G4HadronicInteraction.

Definition at line 219 of file G4ParticleHPFission.cc.

220 {
221  outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n";
222 }
void G4ParticleHPFission::SetVerboseLevel ( G4int  newValue)

Definition at line 176 of file G4ParticleHPFission.cc.

177 {
179 }
static G4ParticleHPManager * GetInstance()

Here is the call graph for this function:


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