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