Geant4  10.00.p02
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 <stdlib.h>
35 
36 #include "G4NeutronHPChannel.hh"
37 #include "globals.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4NeutronHPFinalState.hh"
40 #include "G4HadTmpUtil.hh"
41 
42 #include "G4NeutronHPManager.hh"
44 
46  {
47  return std::max(0., theChannelData->GetXsec(energy));
48  }
49 
51  {
52  return theIsotopeWiseData[isoNumber].GetXsec(energy);
53  }
54 
56  {
57  return theFinalStates[isoNumber]->GetXsec(energy);
58  }
59 
61  Init(G4Element * anElement, const G4String dirName, const G4String aFSType)
62  {
63  theFSType = aFSType;
64  Init(anElement, dirName);
65  }
66 
67  void G4NeutronHPChannel::Init(G4Element * anElement, const G4String dirName)
68  {
69  theDir = dirName;
70  theElement = anElement;
71  }
72 
74  {
75  registerCount++;
76  G4int Z = G4lrint(theElement->GetZ());
77 
78  Z = Z-registerCount;
79  if ( registerCount > 5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material"); // for Elastic, Capture, Fission case
80  if ( Z < 1 ) return false;
81 /*
82  if(registerCount<5)
83  {
84  Z = Z-registerCount;
85  }
86 */
87  //if(Z=theElement->GetZ()-5) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
88  // Bug fix by TK on behalf of AH
89  if ( Z <=theElement->GetZ()-5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
90  G4int count = 0;
91  if(registerCount==0) count = theElement->GetNumberOfIsotopes();
92  if(count == 0||registerCount!=0) count +=
94  niso = count;
95  delete [] theIsotopeWiseData;
97  delete [] active;
98  active = new G4bool[niso];
99 
100  delete [] theFinalStates;
102  delete theChannelData;
104  for(G4int i=0; i<niso; i++)
105  {
106  theFinalStates[i] = theFS->New();
107  }
108  count = 0;
109  G4int nIsos = niso;
110  if(theElement->GetNumberOfIsotopes()!=0&&registerCount==0)
111  {
112  for (G4int i1=0; i1<nIsos; i1++)
113  {
114  // G4cout <<" Init: normal case"<<G4endl;
115  G4int A = theElement->GetIsotope(i1)->GetN();
116  G4int M = theElement->GetIsotope(i1)->Getm();
118  //theFinalStates[i1]->SetA_Z(A, Z);
119  //UpdateData(A, Z, count++, frac);
120  theFinalStates[i1]->SetA_Z(A, Z, M);
121  UpdateData(A, Z, M, count++, frac);
122  }
123  } else {
124  //G4cout <<" Init: mean case: "
125  // <<theStableOnes.GetNumberOfIsotopes(Z)<<" "
126  // <<Z<<" "<<theElement
127  // << G4endl;
129  for(G4int i1=0;
131  i1++)
132  {
134  G4double frac = theStableOnes.GetAbundance(first+i1);
135  theFinalStates[i1]->SetA_Z(A, Z);
136  UpdateData(A, Z, count++, frac);
137  }
138  }
139  G4bool result = HasDataInAnyFinalState();
140  return result;
141  }
142 
143  //void G4NeutronHPChannel::UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
145  {
146  // Initialze the G4FissionFragment generator for this isomer if needed
148  {
150  }
151 
152  theFinalStates[index]->Init(A, Z, M, theDir, theFSType);
153  if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
154 
155  // the above has put the X-sec into the FS
156  theBuffer = 0;
157  if(theFinalStates[index]->HasXsec())
158  {
159  theBuffer = theFinalStates[index]->GetXsec();
160  theBuffer->Times(abundance/100.);
162  }
163  else // get data from CrossSection directory
164  {
165  G4String tString = "/CrossSection";
166  //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
167  active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
168  if(active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
169  }
171  }
172 
174  {
175  G4int s_tmp = 0, n=0, m_tmp=0;
176  G4NeutronHPVector * theMerge = new G4NeutronHPVector;
177  G4NeutronHPVector * anActive = theStore;
178  G4NeutronHPVector * aPassive = theNew;
179  G4NeutronHPVector * tmp;
180  G4int a = s_tmp, p = n, t;
181  while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength())
182  {
183  if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
184  {
185  G4double xa = anActive->GetEnergy(a);
186  theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
187  m_tmp++;
188  a++;
189  G4double xp = aPassive->GetEnergy(p);
190  if( std::abs(std::abs(xp-xa)/xa)<0.001 )
191  {
192  p++;
193  }
194  } else {
195  tmp = anActive; t=a;
196  anActive = aPassive; a=p;
197  aPassive = tmp; p=t;
198  }
199  }
200  while (a!=anActive->GetVectorLength())
201  {
202  theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
203  a++;
204  }
205  while (p!=aPassive->GetVectorLength())
206  {
207  if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
208  theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
209  p++;
210  }
211  delete theStore;
212  theStore = theMerge;
213  }
214 
216 
218  ApplyYourself(const G4HadProjectile & theTrack, G4int anIsotope)
219  {
220 // G4cout << "G4NeutronHPChannel::ApplyYourself+"<<niso<<G4endl;
221  if ( anIsotope != -1 && anIsotope != -2 )
222  {
223  //Inelastic Case
224  //G4cout << "G4NeutronHPChannel Inelastic Case"
225  //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
228  return theFinalStates[anIsotope]->ApplyYourself(theTrack);
229  }
230  G4double sum=0;
231  G4int it=0;
232  G4double * xsec = new G4double[niso];
233  G4NeutronHPThermalBoost aThermalE;
234  for (G4int i=0; i<niso; i++)
235  {
236  if(theFinalStates[i]->HasAnyData())
237  {
238  xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
239  theFinalStates[i]->GetN(),
240  theFinalStates[i]->GetZ(),
241  theTrack.GetMaterial()->GetTemperature()));
242  sum += xsec[i];
243  }
244  else
245  {
246  xsec[i]=0;
247  }
248  }
249  if(sum == 0)
250  {
251 // G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
252 // G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
253  it = static_cast<G4int>(niso*G4UniformRand());
254  }
255  else
256  {
257 // G4cout << "Are we still here? "<<sum<<G4endl;
258 // G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
259  G4double random = G4UniformRand();
260  G4double running=0;
261 // G4cout << "G4NeutronHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
262 // G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
263  for (G4int ix=0; ix<niso; ix++)
264  {
265  running += xsec[ix];
266  //if(random<=running/sum)
267  if( sum == 0 || random <= running/sum )
268  {
269  it = ix;
270  break;
271  }
272  }
273  if(it==niso) it--;
274  }
275  delete [] xsec;
276  G4HadFinalState * theFinalState=0;
277  const G4int A = (G4int)this->GetN(it);
278  const G4int Z = (G4int)this->GetZ(it);
279  const G4int M = (G4int)this->GetM(it);
280 
281  //-2:Marker for Fission
282  if(wendtFissionGenerator&&anIsotope==-2)
283  {
284  theFinalState = wendtFissionGenerator->ApplyYourself(theTrack, Z, A);
285  }
286 
287  // Use the standard procedure if the G4FissionFragmentGenerator model fails
288  if (!theFinalState)
289  {
290  while(theFinalState==0)
291  {
292 // G4cout << "TESTHP 24 it="<<it<<G4endl;
293  theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
294  }
295  }
296 
297  //G4cout <<"THE IMPORTANT RETURN"<<G4endl;
298  //G4cout << "TK G4NeutronHPChannel Elastic, Capture and Fission Cases "
299  //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
303 
304  return theFinalState;
305  }
306 
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double GetEnergy(G4int i) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &projectile, G4int Z, G4int A)
void Harmonise(G4NeutronHPVector *&theStore, G4NeutronHPVector *theNew)
G4int GetVectorLength() const
virtual G4NeutronHPFinalState * New()=0
G4bool Register(G4NeutronHPFinalState *theFS)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4NeutronHPFinalState ** theFinalStates
G4NeutronHPVector * theChannelData
G4double GetZ(G4int i)
void Init(G4Element *theElement, const G4String dirName)
static G4NeutronHPManager * GetInstance()
G4int GetFirstIsotope(G4int Z)
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetZ() const
Definition: G4Element.hh:131
void Times(G4double factor)
G4double a
Definition: TRTMaterials.hh:39
void SetData(G4int i, G4double x, G4double y)
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
int G4int
Definition: G4Types.hh:78
G4NeutronHPIsoData * theIsotopeWiseData
G4NeutronHPVector * theBuffer
G4double GetN(G4int i)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
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:87
G4int Getm() const
Definition: G4Isotope.hh:100
bool G4bool
Definition: G4Types.hh:79
G4StableIsotopes theStableOnes
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetXsec(G4double energy)
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
static const G4double A[nN]
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
G4double GetXsec(G4double energy)
G4int GetNumberOfIsotopes(G4int Z)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4int GetIsotopeNucleonCount(G4int number)
void FillChannelData(G4NeutronHPVector *aBuffer)
G4double GetXsec(G4int i)
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4double energy(const ThreeVector &p, const G4double m)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType)
G4bool HasAnyData(G4int isoNumber)
G4double GetM(G4int i)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4double GetTemperature() const
Definition: G4Material.hh:180
const G4Material * GetMaterial() const
G4double GetAbundance(G4int number)
double G4double
Definition: G4Types.hh:76
virtual G4double GetXsec(G4double)
G4NeutronHPVector * MakeChannelData()
G4WendtFissionFragmentGenerator *const wendtFissionGenerator