Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPChannel.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// P. Arce, June-2014 Conversion neutron_hp to particle_hp
35//
36// June-2019 - E. Mendoza --> Modification to allow using an incomplete
37// data library if the G4NEUTRONHP_SKIP_MISSING_ISOTOPES environmental
38// flag is defined. The missing XS are set to 0.
39
40
41#include <stdlib.h>
42
44#include "globals.hh"
45#include "G4SystemOfUnits.hh"
47#include "G4HadTmpUtil.hh"
50
52{
53 return std::max(0., theChannelData->GetXsec(energy));
54}
55
57{
58 return theIsotopeWiseData[isoNumber].GetXsec(energy);
59}
60
62{
63 return theFinalStates[isoNumber]->GetXsec(energy);
64}
65
67 Init(G4Element * anElement, const G4String dirName, const G4String aFSType)
68 {
69 theFSType = aFSType;
70 Init(anElement, dirName);
71 }
72
73 void G4ParticleHPChannel::Init(G4Element * anElement, const G4String dirName)
74 {
75 theDir = dirName;
76 theElement = anElement;
77 }
78
80 {
81 registerCount++;
82 G4int Z = G4lrint(theElement->GetZ());
83
84 Z = Z-registerCount;
85 if ( registerCount > 5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material"); // for Elastic, Capture, Fission case
86 if ( Z < 1 ) return false;
87/*
88 if(registerCount<5)
89 {
90 Z = Z-registerCount;
91 }
92*/
93 //if(Z=theElement->GetZ()-5) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
94 // Bug fix by TK on behalf of AH
95 //if ( Z <=theElement->GetZ()-5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
96 G4int count = 0;
97 if(registerCount==0) count = (G4int)theElement->GetNumberOfIsotopes();
98 if(count == 0||registerCount!=0) count +=
99 theStableOnes.GetNumberOfIsotopes(Z);
100 niso = count;
101 delete [] theIsotopeWiseData;
102 theIsotopeWiseData = new G4ParticleHPIsoData [niso];
103 delete [] active;
104 active = new G4bool[niso];
105
106 delete [] theFinalStates;
107 theFinalStates = new G4ParticleHPFinalState * [niso];
108 delete theChannelData;
109 theChannelData = new G4ParticleHPVector;
110 for(G4int i=0; i<niso; ++i)
111 {
112 theFinalStates[i] = theFS->New();
113 theFinalStates[i]->SetProjectile(theProjectile);
114 }
115 count = 0;
116 G4int nIsos = niso;
117 if(theElement->GetNumberOfIsotopes()!=0&&registerCount==0)
118 {
119 for (G4int i1=0; i1<nIsos; ++i1)
120 {
121 // G4cout <<" Init: normal case"<<G4endl;
122 G4int A = theElement->GetIsotope(i1)->GetN();
123 G4int M = theElement->GetIsotope(i1)->Getm();
124 G4double frac = theElement->GetRelativeAbundanceVector()[i1]/perCent;
125 //theFinalStates[i1]->SetA_Z(A, Z);
126 //UpdateData(A, Z, count++, frac);
127 theFinalStates[i1]->SetA_Z(A, Z, M);
128 UpdateData(A, Z, M, count++, frac, theProjectile);
129 }
130 } else {
131 //G4cout <<" Init: mean case: "
132 // <<theStableOnes.GetNumberOfIsotopes(Z)<<" "
133 // <<Z<<" "<<theElement
134 // << G4endl;
135 G4int first = theStableOnes.GetFirstIsotope(Z);
136 for(G4int i1=0; i1<theStableOnes.GetNumberOfIsotopes(Z); ++i1)
137 {
138 G4int A = theStableOnes.GetIsotopeNucleonCount(first+i1);
139 G4double frac = theStableOnes.GetAbundance(first+i1);
140 theFinalStates[i1]->SetA_Z(A, Z);
141 UpdateData(A, Z, count++, frac, theProjectile);
142 }
143 }
145
146 //To avoid issuing hash by worker threads
147 if ( result ) theChannelData->Hash();
148
149 return result;
150 }
151
153 {
154 // Initialze the G4FissionFragment generator for this isomer if needed
155 if(wendtFissionGenerator)
156 {
157 wendtFissionGenerator->InitializeANucleus(A, Z, M, theDir);
158 }
159
160 theFinalStates[index]->Init(A, Z, M, theDir, theFSType, projectile);
161 if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
162
163 // the above has put the X-sec into the FS
164 theBuffer = 0;
165 if(theFinalStates[index]->HasXsec())
166 {
167 theBuffer = theFinalStates[index]->GetXsec();
168 theBuffer->Times(abundance/100.);
169 theIsotopeWiseData[index].FillChannelData(theBuffer);
170 }
171 else // get data from CrossSection directory
172 {
173 G4String tString = "/CrossSection";
174 //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
175 active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
176 if(active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
177 }
178 if(theBuffer != 0) Harmonise(theChannelData, theBuffer);
179 }
180
182 {
183 G4int s_tmp = 0, n=0, m_tmp=0;
184 G4ParticleHPVector * theMerge = new G4ParticleHPVector;
185 G4ParticleHPVector * anActive = theStore;
186 G4ParticleHPVector * aPassive = theNew;
187 G4ParticleHPVector * tmp;
188 G4int a = s_tmp, p = n, t;
189 while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
190 {
191 if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
192 {
193 G4double xa = anActive->GetEnergy(a);
194 theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
195 m_tmp++;
196 a++;
197 G4double xp = aPassive->GetEnergy(p);
198 if( std::abs(std::abs(xp-xa)/xa)<0.001 )
199 {
200 ++p;
201 }
202 } else {
203 tmp = anActive; t=a;
204 anActive = aPassive; a=p;
205 aPassive = tmp; p=t;
206 }
207 }
208 while (a!=anActive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
209 {
210 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
211 ++a;
212 }
213 while (p!=aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
214 {
215 if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
216 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
217 ++p;
218 }
219 delete theStore;
220 theStore = theMerge;
221 }
222
224 ApplyYourself(const G4HadProjectile & theTrack, G4int anIsotope)
225 {
226// G4cout << "G4ParticleHPChannel::ApplyYourself+"<<niso<<G4endl;
227 if ( anIsotope != -1 && anIsotope != -2 )
228 {
229 //Inelastic Case
230 //G4cout << "G4ParticleHPChannel Inelastic Case"
231 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
234 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
235 }
236 G4double sum=0;
237 G4int it=0;
238 G4double * xsec = new G4double[niso];
239 G4ParticleHPThermalBoost aThermalE;
240 for (G4int i=0; i<niso; i++)
241 {
242 if(theFinalStates[i]->HasAnyData())
243 {
244 xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
245 theFinalStates[i]->GetN(),
246 theFinalStates[i]->GetZ(),
247 theTrack.GetMaterial()->GetTemperature()));
248 sum += xsec[i];
249 }
250 else
251 {
252 xsec[i]=0;
253 }
254 }
255 if(sum == 0)
256 {
257// G4cout << "G4ParticleHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
258// G4cout << "G4ParticleHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
259 it = static_cast<G4int>(niso*G4UniformRand());
260 }
261 else
262 {
263// G4cout << "Are we still here? "<<sum<<G4endl;
264// G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
265 G4double random = G4UniformRand();
266 G4double running=0;
267// G4cout << "G4ParticleHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
268// G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
269 for (G4int ix=0; ix<niso; ix++)
270 {
271 running += xsec[ix];
272 //if(random<=running/sum)
273 if( sum == 0 || random <= running/sum )
274 {
275 it = ix;
276 break;
277 }
278 }
279 if(it==niso) it--;
280 }
281 delete [] xsec;
282 G4HadFinalState * theFinalState=0;
283 const G4int A = (G4int)this->GetN(it);
284 const G4int Z = (G4int)this->GetZ(it);
285 const G4int M = (G4int)this->GetM(it);
286
287 //-2:Marker for Fission
288 if(wendtFissionGenerator&&anIsotope==-2)
289 {
290 theFinalState = wendtFissionGenerator->ApplyYourself(theTrack, Z, A);
291 }
292
293 // Use the standard procedure if the G4FissionFragmentGenerator model fails
294 if (!theFinalState)
295 {
296
297 G4int icounter=0;
298 G4int icounter_max=1024;
299 while(theFinalState==0) // Loop checking, 11.05.2015, T. Koi
300 {
301 icounter++;
302 if ( icounter > icounter_max ) {
303 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
304 break;
305 }
306// G4cout << "TESTHP 24 it="<<it<<G4endl;
307 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
308 }
309 }
310
311 //G4cout <<"THE IMPORTANT RETURN"<<G4endl;
312 //G4cout << "TK G4ParticleHPChannel Elastic, Capture and Fission Cases "
313 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
317
318 return theFinalState;
319 }
320
321
323
324 G4cout<<" Element: "<<theElement->GetName()<<G4endl;
325 G4cout<<" Directory name: "<<theDir<<G4endl;
326 G4cout<<" FS name: "<<theFSType<<G4endl;
327 G4cout<<" Number of Isotopes: "<<niso<<G4endl;
328 G4cout<<" Have cross sections: "<<G4endl;
329 for(int i=0;i<niso;i++){
330 G4cout<<theFinalStates[i]->HasXsec()<<" ";
331 }
332 G4cout<<G4endl;
333 if(theChannelData){
334 G4cout<<" Cross Section (total for this channel):"<<G4endl;
335 int np=theChannelData->GetVectorLength();
336 G4cout<<np<<G4endl;
337 for(int i=0;i<np;i++){
338 G4cout<<theChannelData->GetEnergy(i)/eV<<" "<<theChannelData->GetXsec(i)<<G4endl;
339 }
340 }
341
342}
343
344
345
#define M(row, col)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4double GetZ() const
Definition: G4Element.hh:131
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4String & GetName() const
Definition: G4Element.hh:127
const G4Material * GetMaterial() const
G4int Getm() const
Definition: G4Isotope.hh:99
G4int GetN() const
Definition: G4Isotope.hh:93
G4double GetTemperature() const
Definition: G4Material.hh:177
G4bool HasAnyData(G4int isoNumber)
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
void Init(G4Element *theElement, const G4String dirName)
void Harmonise(G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4double GetM(G4int i)
G4double GetN(G4int i)
G4bool Register(G4ParticleHPFinalState *theFS)
G4double GetXsec(G4double energy)
G4double GetZ(G4int i)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
virtual G4double GetXsec(G4double)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile)
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
virtual G4ParticleHPFinalState * New()=0
void SetProjectile(G4ParticleDefinition *projectile)
G4ParticleHPVector * MakeChannelData()
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
void FillChannelData(G4ParticleHPVector *aBuffer)
G4double GetXsec(G4double energy)
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void SetData(G4int i, G4double x, G4double y)
void Times(G4double factor)
G4double GetXsec(G4int i)
G4double GetEnergy(G4int i) const
G4int GetVectorLength() const
G4double GetAbundance(G4int number)
G4int GetFirstIsotope(G4int Z)
G4int GetNumberOfIsotopes(G4int Z)
G4int GetIsotopeNucleonCount(G4int number)
void InitializeANucleus(const G4int A, const G4int Z, const G4int M, const G4String &dataDirectory)
G4HadFinalState * ApplyYourself(const G4HadProjectile &projectile, G4int Z, G4int A)
int G4lrint(double ad)
Definition: templates.hh:134