Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPorLEInelastic.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//
32
33// neutron_hp -- source file
34// J.P. Wellisch, Nov-1996
35// A prototype of the low energy neutron transport model.
36//
38#include "G4SystemOfUnits.hh"
39//#include "G4NeutronHPInelasticFS.hh"
40
42 :G4HadronicInteraction("NeutronHPorLEInelastic")
43{
44 SetMinEnergy(0.*eV);
45 SetMaxEnergy(20.*MeV);
46
47// G4NeutronHPInelasticFS * theFS = new G4NeutronHPInelasticFS;
48 if(!getenv("G4NEUTRONHPDATA"))
49 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
50 dirName = getenv("G4NEUTRONHPDATA");
51 G4String tString = "/Inelastic/";
52 dirName = dirName + tString;
53// G4cout <<"G4NeutronHPorLEInelastic::G4NeutronHPorLEInelastic testit "<<dirName<<G4endl;
55 theInelastic = new G4NeutronHPChannelList[numEle];
56 unavailable_elements.clear();
57 for (G4int i=0; i<numEle; i++)
58 {
59 theInelastic[i].Init( (*(G4Element::GetElementTable()))[i] , dirName );
60 do
61 {
62
63 try
64 {
65 theInelastic[i].Register(&theNFS, "F01"); // has
66 theInelastic[i].Register(&theNXFS, "F02");
67 theInelastic[i].Register(&the2NDFS, "F03");
68 theInelastic[i].Register(&the2NFS, "F04"); // has, E Done
69 theInelastic[i].Register(&the3NFS, "F05"); // has, E Done
70 theInelastic[i].Register(&theNAFS, "F06");
71 theInelastic[i].Register(&theN3AFS, "F07");
72 theInelastic[i].Register(&the2NAFS, "F08");
73 theInelastic[i].Register(&the3NAFS, "F09");
74 theInelastic[i].Register(&theNPFS, "F10");
75 theInelastic[i].Register(&theN2AFS, "F11");
76 theInelastic[i].Register(&the2N2AFS, "F12");
77 theInelastic[i].Register(&theNDFS, "F13");
78 theInelastic[i].Register(&theNTFS, "F14");
79 theInelastic[i].Register(&theNHe3FS, "F15");
80 theInelastic[i].Register(&theND2AFS, "F16");
81 theInelastic[i].Register(&theNT2AFS, "F17");
82 theInelastic[i].Register(&the4NFS, "F18"); // has, E Done
83 theInelastic[i].Register(&the2NPFS, "F19");
84 theInelastic[i].Register(&the3NPFS, "F20");
85 theInelastic[i].Register(&theN2PFS, "F21");
86 theInelastic[i].Register(&theNPAFS, "F22");
87 theInelastic[i].Register(&thePFS, "F23");
88 theInelastic[i].Register(&theDFS, "F24");
89 theInelastic[i].Register(&theTFS, "F25");
90 theInelastic[i].Register(&theHe3FS, "F26");
91 theInelastic[i].Register(&theAFS, "F27");
92 theInelastic[i].Register(&the2AFS, "F28");
93 theInelastic[i].Register(&the3AFS, "F29");
94 theInelastic[i].Register(&the2PFS, "F30");
95 theInelastic[i].Register(&thePAFS, "F31");
96 theInelastic[i].Register(&theD2AFS, "F32");
97 theInelastic[i].Register(&theT2AFS, "F33");
98 theInelastic[i].Register(&thePDFS, "F34");
99 theInelastic[i].Register(&thePTFS, "F35");
100 theInelastic[i].Register(&theDAFS, "F36");
101 }
102
103 catch ( G4HadronicException )
104 {
105 unavailable_elements.insert ( (*(G4Element::GetElementTable()))[i]->GetName() );
106 }
107 theInelastic[i].RestartRegistration();
108 }
109 while( !theInelastic[i].HasDataInAnyFinalState());
110
111 }
112
113// delete theFS;
114 if ( unavailable_elements.size() > 0 )
115 {
116 std::set< G4String>::iterator it;
117 G4cout << "HP Inelastic data are not available for thess elements "<< G4endl;
118 for ( it = unavailable_elements.begin() ; it != unavailable_elements.end() ; it++ )
119 G4cout << *it << G4endl;
120 G4cout << "Low Energy Parameterization Models will be used."<< G4endl;
121 }
122
123 createXSectionDataSet();
124}
125
127{
128 delete [] theInelastic;
129 delete theDataSet;
130}
131
133
135{
136 G4int it=0;
137 const G4Material * theMaterial = aTrack.GetMaterial();
138 G4int n = theMaterial->GetNumberOfElements();
139 G4int index = theMaterial->GetElement(0)->GetIndex();
140 if(n!=1)
141 {
142 G4int i;
143 xSec = new G4double[n];
144 G4double sum=0;
145 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
146 G4double rWeight;
147 G4NeutronHPThermalBoost aThermalE;
148 for (i=0; i<n; i++)
149 {
150 index = theMaterial->GetElement(i)->GetIndex();
151 rWeight = NumAtomsPerVolume[i];
152 G4double x = aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i), theMaterial->GetTemperature());
153
154 //xSec[i] = theInelastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
155 // theMaterial->GetElement(i),
156 // theMaterial->GetTemperature()));
157 xSec[i] = theInelastic[index].GetXsec(x);
158
159 xSec[i] *= rWeight;
160 sum+=xSec[i];
161 }
162 G4double random = G4UniformRand();
163 G4double running = 0;
164 for (i=0; i<n; i++)
165 {
166 running += xSec[i];
167 index = theMaterial->GetElement(i)->GetIndex();
168 it = i;
169 if(random<=running/sum) break;
170 }
171 delete [] xSec;
172 // it is element-wise initialised.
173 }
174 //return theInelastic[index].ApplyYourself(aTrack);
175 return theInelastic[index].ApplyYourself( theMaterial->GetElement(it) , aTrack );
176}
177
178
179
181{
182 if ( unavailable_elements.find( name ) == unavailable_elements.end() )
183 return TRUE;
184 else
185 return FALSE;
186}
187
188
189
190void G4NeutronHPorLEInelastic::createXSectionDataSet()
191{
192 theDataSet = new G4NeutronHPorLEInelasticData ( theInelastic , &unavailable_elements );
193}
194const std::pair<G4double, G4double> G4NeutronHPorLEInelastic::GetFatalEnergyCheckLevels() const
195{
196 //return std::pair<G4double, G4double>(10*perCent,10*GeV);
197 return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
198}
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
G4double GetXsec(G4double anEnergy)
void Register(G4NeutronHPFinalState *theFS, const G4String &aName)
G4HadFinalState * ApplyYourself(const G4Element *theElement, const G4HadProjectile &aTrack)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52
#define DBL_MAX
Definition: templates.hh:83