Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPChannel.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// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
31// 071031 bug fix T. Koi on behalf of A. Howard
32// 081203 bug fix in Register method by T. Koi
33//
34#include "G4NeutronHPChannel.hh"
35#include "globals.hh"
36#include "G4SystemOfUnits.hh"
38#include "G4HadTmpUtil.hh"
39
40#include "G4NeutronHPManager.hh"
42
44 {
45 return std::max(0., theChannelData->GetXsec(energy));
46 }
47
49 {
50 return theIsotopeWiseData[isoNumber].GetXsec(energy);
51 }
52
54 {
55 return theFinalStates[isoNumber]->GetXsec(energy);
56 }
57
59 Init(G4Element * anElement, const G4String dirName, const G4String aFSType)
60 {
61 theFSType = aFSType;
62 Init(anElement, dirName);
63 }
64
65 void G4NeutronHPChannel::Init(G4Element * anElement, const G4String dirName)
66 {
67 theDir = dirName;
68 theElement = anElement;
69 }
70
72 {
73 registerCount++;
74 G4int Z = G4lrint(theElement->GetZ());
75
76 Z = Z-registerCount;
77 if ( registerCount > 5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material"); // for Elastic, Capture, Fission case
78 if ( Z < 1 ) return false;
79/*
80 if(registerCount<5)
81 {
82 Z = Z-registerCount;
83 }
84*/
85 //if(Z=theElement->GetZ()-5) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
86 // Bug fix by TK on behalf of AH
87 if ( Z <=theElement->GetZ()-5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
88 G4int count = 0;
89 if(registerCount==0) count = theElement->GetNumberOfIsotopes();
90 if(count == 0||registerCount!=0) count +=
91 theStableOnes.GetNumberOfIsotopes(Z);
92 niso = count;
93 delete [] theIsotopeWiseData;
94 theIsotopeWiseData = new G4NeutronHPIsoData [niso];
95 delete [] active;
96 active = new G4bool[niso];
97
98 delete [] theFinalStates;
99 theFinalStates = new G4NeutronHPFinalState * [niso];
100 delete theChannelData;
101 theChannelData = new G4NeutronHPVector;
102 for(G4int i=0; i<niso; i++)
103 {
104 theFinalStates[i] = theFS->New();
105 }
106 count = 0;
107 G4int nIsos = niso;
108 if(theElement->GetNumberOfIsotopes()!=0&&registerCount==0)
109 {
110 for (G4int i1=0; i1<nIsos; i1++)
111 {
112 // G4cout <<" Init: normal case"<<G4endl;
113 G4int A = theElement->GetIsotope(i1)->GetN();
114 G4int M = theElement->GetIsotope(i1)->Getm();
115 G4double frac = theElement->GetRelativeAbundanceVector()[i1]/perCent;
116 //theFinalStates[i1]->SetA_Z(A, Z);
117 //UpdateData(A, Z, count++, frac);
118 theFinalStates[i1]->SetA_Z(A, Z, M);
119 UpdateData(A, Z, M, count++, frac);
120 }
121 } else {
122 //G4cout <<" Init: mean case: "
123 // <<theStableOnes.GetNumberOfIsotopes(Z)<<" "
124 // <<Z<<" "<<theElement
125 // << G4endl;
126 G4int first = theStableOnes.GetFirstIsotope(Z);
127 for(G4int i1=0;
128 i1<theStableOnes.GetNumberOfIsotopes(Z);
129 i1++)
130 {
131 G4int A = theStableOnes.GetIsotopeNucleonCount(first+i1);
132 G4double frac = theStableOnes.GetAbundance(first+i1);
133 theFinalStates[i1]->SetA_Z(A, Z);
134 UpdateData(A, Z, count++, frac);
135 }
136 }
138 return result;
139 }
140
141 //void G4NeutronHPChannel::UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
143 {
144 //theFinalStates[index]->Init(A, Z, theDir, theFSType);
145 theFinalStates[index]->Init(A, Z, M, theDir, theFSType);
146 if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
147
148 // the above has put the X-sec into the FS
149 theBuffer = 0;
150 if(theFinalStates[index]->HasXsec())
151 {
152 theBuffer = theFinalStates[index]->GetXsec();
153 theBuffer->Times(abundance/100.);
154 theIsotopeWiseData[index].FillChannelData(theBuffer);
155 }
156 else // get data from CrossSection directory
157 {
158 G4String tString = "/CrossSection";
159 //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
160 active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
161 if(active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
162 }
163 if(theBuffer != 0) Harmonise(theChannelData, theBuffer);
164 }
165
167 {
168 G4int s_tmp = 0, n=0, m_tmp=0;
169 G4NeutronHPVector * theMerge = new G4NeutronHPVector;
170 G4NeutronHPVector * anActive = theStore;
171 G4NeutronHPVector * aPassive = theNew;
172 G4NeutronHPVector * tmp;
173 G4int a = s_tmp, p = n, t;
174 while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength())
175 {
176 if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
177 {
178 G4double xa = anActive->GetEnergy(a);
179 theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
180 m_tmp++;
181 a++;
182 G4double xp = aPassive->GetEnergy(p);
183 if( std::abs(std::abs(xp-xa)/xa)<0.001 )
184 {
185 p++;
186 }
187 } else {
188 tmp = anActive; t=a;
189 anActive = aPassive; a=p;
190 aPassive = tmp; p=t;
191 }
192 }
193 while (a!=anActive->GetVectorLength())
194 {
195 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
196 a++;
197 }
198 while (p!=aPassive->GetVectorLength())
199 {
200 if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
201 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
202 p++;
203 }
204 delete theStore;
205 theStore = theMerge;
206 }
207
209
211 ApplyYourself(const G4HadProjectile & theTrack, G4int anIsotope)
212 {
213// G4cout << "G4NeutronHPChannel::ApplyYourself+"<<niso<<G4endl;
214 if ( anIsotope != -1 )
215 {
216 //Inelastic Case
217 //G4cout << "G4NeutronHPChannel Inelastic Case"
218 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
221 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
222 }
223 G4double sum=0;
224 G4int it=0;
225 G4double * xsec = new G4double[niso];
226 G4NeutronHPThermalBoost aThermalE;
227 for (G4int i=0; i<niso; i++)
228 {
229 if(theFinalStates[i]->HasAnyData())
230 {
231 xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
232 theFinalStates[i]->GetN(),
233 theFinalStates[i]->GetZ(),
234 theTrack.GetMaterial()->GetTemperature()));
235 sum += xsec[i];
236 }
237 else
238 {
239 xsec[i]=0;
240 }
241 }
242 if(sum == 0)
243 {
244// G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
245// G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
246 it = static_cast<G4int>(niso*G4UniformRand());
247 }
248 else
249 {
250// G4cout << "Are we still here? "<<sum<<G4endl;
251// G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
252 G4double random = G4UniformRand();
253 G4double running=0;
254// G4cout << "G4NeutronHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
255// G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
256 for (G4int ix=0; ix<niso; ix++)
257 {
258 running += xsec[ix];
259 //if(random<=running/sum)
260 if( sum == 0 || random <= running/sum )
261 {
262 it = ix;
263 break;
264 }
265 }
266 if(it==niso) it--;
267 }
268 delete [] xsec;
269 G4HadFinalState * theFinalState=0;
270 while(theFinalState==0)
271 {
272// G4cout << "TESTHP 24 it="<<it<<G4endl;
273 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
274 }
275// G4cout <<"THE IMPORTANT RETURN"<<G4endl;
276 //G4cout << "TK G4NeutronHPChannel Elastic, Capture and Fission Cases "
277 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
280 return theFinalState;
281 }
282
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:131
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
const G4Material * GetMaterial() const
G4int Getm() const
Definition: G4Isotope.hh:100
G4int GetN() const
Definition: G4Isotope.hh:94
G4double GetTemperature() const
Definition: G4Material.hh:181
G4double GetN(G4int i)
G4bool Register(G4NeutronHPFinalState *theFS)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
void Harmonise(G4NeutronHPVector *&theStore, G4NeutronHPVector *theNew)
G4bool HasAnyData(G4int isoNumber)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
G4double GetXsec(G4double energy)
G4double GetZ(G4int i)
void Init(G4Element *theElement, const G4String dirName)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
virtual G4double GetXsec(G4double)
virtual G4NeutronHPFinalState * New()=0
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType)
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
G4NeutronHPVector * MakeChannelData()
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
void FillChannelData(G4NeutronHPVector *aBuffer)
G4double GetXsec(G4double energy)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4int GetVectorLength() const
G4double GetEnergy(G4int i) const
G4double GetXsec(G4int i)
void SetData(G4int i, G4double x, G4double y)
void Times(G4double factor)
G4double GetAbundance(G4int number)
G4int GetFirstIsotope(G4int Z)
G4int GetNumberOfIsotopes(G4int Z)
G4int GetIsotopeNucleonCount(G4int number)
int G4lrint(double ad)
Definition: templates.hh:163