Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4NeutronHPorLFission.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// 05-11-21 NeutronHP or Low Energy Parameterization Models
28// Implemented by T. Koi (SLAC/SCCS)
29// If NeutronHP data do not available for an element, then Low Energy
30// Parameterization models handle the interactions of the element.
31// 080319 Compilation warnings - gcc-4.3.0 fix by T. Koi
32//
33
34// neutron_hp -- source file
35// J.P. Wellisch, Nov-1996
36// A prototype of the low energy neutron transport model.
37//
39#include "G4SystemOfUnits.hh"
41
43 :G4HadronicInteraction("NeutronHPorLFission")
44{
45 SetMinEnergy(0.*eV);
46 SetMaxEnergy(20.*MeV);
47
48 if( !getenv("G4NEUTRONHPDATA") )
49 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
50
51 dirName = getenv("G4NEUTRONHPDATA");
52 G4String tString = "/Fission/";
53 dirName = dirName + tString;
54// G4cout <<"G4NeutronHPorLFission::G4NeutronHPorLFission testit "<<dirName<<G4endl;
55 unavailable_elements.clear();
56
58 theFission = new G4NeutronHPChannel[numEle];
59
60 for ( G4int i = 0; i < numEle ; i++)
61 {
62 if ( (*(G4Element::GetElementTable()))[i]-> GetZ() > 87 ) //TK modified for ENDF-VII
63 {
64 theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
65 try { while(!theFission[i].Register(&theFS)) ; }
66 catch ( G4HadronicException )
67 {
68 unavailable_elements.insert ( (*(G4Element::GetElementTable()))[i]->GetName() );
69 }
70 }
71 }
72
73 if ( unavailable_elements.size() > 0 )
74 {
75 std::set< G4String>::iterator it;
76 G4cout << "HP Fission data are not available for thess elements "<< G4endl;
77 for ( it = unavailable_elements.begin() ; it != unavailable_elements.end() ; it++ )
78 {
79 G4cout << *it << G4endl;
80 }
81 G4cout << "Low Energy Parameterization Models will be used."<< G4endl;
82 }
83
84 createXSectionDataSet();
85}
86
88{
89 delete [] theFission;
90 delete theDataSet;
91}
92
94
96{
97 const G4Material * theMaterial = aTrack.GetMaterial();
98 G4int n = theMaterial->GetNumberOfElements();
99 G4int index = theMaterial->GetElement(0)->GetIndex();
100 if(n!=1)
101 {
102 G4int i;
103 xSec = new G4double[n];
104 G4double sum=0;
105 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
106 G4double rWeight;
107 G4NeutronHPThermalBoost aThermalE;
108 for (i=0; i<n; i++)
109 {
110 index = theMaterial->GetElement(i)->GetIndex();
111 rWeight = NumAtomsPerVolume[i];
112 G4double x = aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i), theMaterial->GetTemperature());
113
114 //xSec[i] = theFission[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
115 // theMaterial->GetElement(i),
116 // theMaterial->GetTemperature()));
117 xSec[i] = theFission[index].GetXsec(x);
118
119 xSec[i] *= rWeight;
120 sum+=xSec[i];
121 }
122 G4double random = G4UniformRand();
123 G4double running = 0;
124 for (i=0; i<n; i++)
125 {
126 running += xSec[i];
127 index = theMaterial->GetElement(i)->GetIndex();
128 if(random<=running/sum) break;
129 }
130 delete [] xSec;
131 // it is element-wise initialised.
132 }
133 return theFission[index].ApplyYourself(aTrack);
134}
135
136
137
139{
140 if ( unavailable_elements.find( name ) == unavailable_elements.end() )
141 return TRUE;
142 else
143 return FALSE;
144}
145
146
147
148void G4NeutronHPorLFission::createXSectionDataSet()
149{
150 theDataSet = new G4NeutronHPorLFissionData ( theFission , &unavailable_elements );
151}
152const std::pair<G4double, G4double> G4NeutronHPorLFission::GetFatalEnergyCheckLevels() const
153{
154 //return std::pair<G4double, G4double>(10*perCent,10*GeV);
155 return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
156}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
size_t GetIndex() const
Definition: G4Element.hh:182
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
const G4Material * GetMaterial() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
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 &theTrack, G4int isoNumber=-1)
G4double GetXsec(G4double energy)
void Init(G4Element *theElement, const G4String dirName)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
G4bool IsThisElementOK(G4String)
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52
#define DBL_MAX
Definition: templates.hh:83