Geant4 10.7.0
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"
48
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 = 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;
137 i1<theStableOnes.GetNumberOfIsotopes(Z);
138 i1++)
139 {
140 G4int A = theStableOnes.GetIsotopeNucleonCount(first+i1);
141 G4double frac = theStableOnes.GetAbundance(first+i1);
142 theFinalStates[i1]->SetA_Z(A, Z);
143 UpdateData(A, Z, count++, frac, theProjectile);
144 }
145 }
147
148 //To avoid issuing hash by worker threads
149 if ( result ) theChannelData->Hash();
150
151 return result;
152 }
153
154 //void G4ParticleHPChannel::UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
156 {
157 // Initialze the G4FissionFragment generator for this isomer if needed
158 if(wendtFissionGenerator)
159 {
160 wendtFissionGenerator->InitializeANucleus(A, Z, M, theDir);
161 }
162
163 theFinalStates[index]->Init(A, Z, M, theDir, theFSType, projectile);
164 if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
165
166 // the above has put the X-sec into the FS
167 theBuffer = 0;
168 if(theFinalStates[index]->HasXsec())
169 {
170 theBuffer = theFinalStates[index]->GetXsec();
171 theBuffer->Times(abundance/100.);
172 theIsotopeWiseData[index].FillChannelData(theBuffer);
173 }
174 else // get data from CrossSection directory
175 {
176 G4String tString = "/CrossSection";
177 //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
178 active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
179 if(active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
180 }
181 if(theBuffer != 0) Harmonise(theChannelData, theBuffer);
182 }
183
185 {
186 G4int s_tmp = 0, n=0, m_tmp=0;
187 G4ParticleHPVector * theMerge = new G4ParticleHPVector;
188 G4ParticleHPVector * anActive = theStore;
189 G4ParticleHPVector * aPassive = theNew;
190 G4ParticleHPVector * tmp;
191 G4int a = s_tmp, p = n, t;
192 while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
193 {
194 if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
195 {
196 G4double xa = anActive->GetEnergy(a);
197 theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
198 m_tmp++;
199 a++;
200 G4double xp = aPassive->GetEnergy(p);
201 if( std::abs(std::abs(xp-xa)/xa)<0.001 )
202 {
203 p++;
204 }
205 } else {
206 tmp = anActive; t=a;
207 anActive = aPassive; a=p;
208 aPassive = tmp; p=t;
209 }
210 }
211 while (a!=anActive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
212 {
213 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
214 a++;
215 }
216 while (p!=aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
217 {
218 if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
219 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
220 p++;
221 }
222 delete theStore;
223 theStore = theMerge;
224 }
225
227
229 ApplyYourself(const G4HadProjectile & theTrack, G4int anIsotope)
230 {
231// G4cout << "G4ParticleHPChannel::ApplyYourself+"<<niso<<G4endl;
232 if ( anIsotope != -1 && anIsotope != -2 )
233 {
234 //Inelastic Case
235 //G4cout << "G4ParticleHPChannel Inelastic Case"
236 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
239 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
240 }
241 G4double sum=0;
242 G4int it=0;
243 G4double * xsec = new G4double[niso];
244 G4ParticleHPThermalBoost aThermalE;
245 for (G4int i=0; i<niso; i++)
246 {
247 if(theFinalStates[i]->HasAnyData())
248 {
249 xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
250 theFinalStates[i]->GetN(),
251 theFinalStates[i]->GetZ(),
252 theTrack.GetMaterial()->GetTemperature()));
253 sum += xsec[i];
254 }
255 else
256 {
257 xsec[i]=0;
258 }
259 }
260 if(sum == 0)
261 {
262// G4cout << "G4ParticleHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
263// G4cout << "G4ParticleHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
264 it = static_cast<G4int>(niso*G4UniformRand());
265 }
266 else
267 {
268// G4cout << "Are we still here? "<<sum<<G4endl;
269// G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
270 G4double random = G4UniformRand();
271 G4double running=0;
272// G4cout << "G4ParticleHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
273// G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
274 for (G4int ix=0; ix<niso; ix++)
275 {
276 running += xsec[ix];
277 //if(random<=running/sum)
278 if( sum == 0 || random <= running/sum )
279 {
280 it = ix;
281 break;
282 }
283 }
284 if(it==niso) it--;
285 }
286 delete [] xsec;
287 G4HadFinalState * theFinalState=0;
288 const G4int A = (G4int)this->GetN(it);
289 const G4int Z = (G4int)this->GetZ(it);
290 const G4int M = (G4int)this->GetM(it);
291
292 //-2:Marker for Fission
293 if(wendtFissionGenerator&&anIsotope==-2)
294 {
295 theFinalState = wendtFissionGenerator->ApplyYourself(theTrack, Z, A);
296 }
297
298 // Use the standard procedure if the G4FissionFragmentGenerator model fails
299 if (!theFinalState)
300 {
301
302 G4int icounter=0;
303 G4int icounter_max=1024;
304 while(theFinalState==0) // Loop checking, 11.05.2015, T. Koi
305 {
306 icounter++;
307 if ( icounter > icounter_max ) {
308 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
309 break;
310 }
311// G4cout << "TESTHP 24 it="<<it<<G4endl;
312 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
313 }
314 }
315
316 //G4cout <<"THE IMPORTANT RETURN"<<G4endl;
317 //G4cout << "TK G4ParticleHPChannel Elastic, Capture and Fission Cases "
318 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
322
323 return theFinalState;
324 }
325
326
328
329 G4cout<<" Element: "<<theElement->GetName()<<G4endl;
330 G4cout<<" Directory name: "<<theDir<<G4endl;
331 G4cout<<" FS name: "<<theFSType<<G4endl;
332 G4cout<<" Number of Isotopes: "<<niso<<G4endl;
333 G4cout<<" Have cross sections: "<<G4endl;
334 for(int i=0;i<niso;i++){
335 G4cout<<theFinalStates[i]->HasXsec()<<" ";
336 }
337 G4cout<<G4endl;
338 if(theChannelData){
339 G4cout<<" Cross Section (total for this channel):"<<G4endl;
340 int np=theChannelData->GetVectorLength();
341 G4cout<<np<<G4endl;
342 for(int i=0;i<np;i++){
343 G4cout<<theChannelData->GetEnergy(i)/eV<<" "<<theChannelData->GetXsec(i)<<G4endl;
344 }
345 }
346
347}
348
349
350
double A(double temperature)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:130
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
const G4String & GetName() const
Definition: G4Element.hh:126
const G4Material * GetMaterial() const
G4int Getm() const
Definition: G4Isotope.hh:99
G4int GetN() const
Definition: G4Isotope.hh:93
G4double GetTemperature() const
Definition: G4Material.hh:180
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