Geant4  10.01.p03
G4NeutronHPChannelList.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 //
33 #include "G4Element.hh"
34 #include "G4HadFinalState.hh"
35 #include "G4HadProjectile.hh"
36 #include "G4NeutronHPFinalState.hh"
37 #include "G4IonTable.hh"
38 
40 
42  {
43  nChannels = n;
45  allChannelsCreated = false;
46  theInitCount = 0;
47  }
48 
50  {
51  nChannels = 0;
52  theChannels = 0;
53  allChannelsCreated = false;
54  theInitCount = 0;
55  }
56 
58  {
59  if(theChannels!=0)
60  {
61  for(G4int i=0;i<nChannels; i++)
62  {
63  delete theChannels[i];
64  }
65  delete [] theChannels;
66  }
67  }
68 
70  #include "G4NeutronHPManager.hh"
72  {
73  G4NeutronHPThermalBoost aThermalE;
74  G4int i, ii;
75  // decide on the isotope
76  G4int numberOfIsos(0);
77  for(ii=0; ii<nChannels; ii++)
78  {
79  numberOfIsos = theChannels[ii]->GetNiso();
80  if(numberOfIsos!=0) break;
81  }
82  G4double * running= new G4double [numberOfIsos];
83  running[0] = 0;
84  for(i=0;i<numberOfIsos; i++)
85  {
86  if(i!=0) running[i] = running[i-1];
87  for(ii=0; ii<nChannels; ii++)
88  {
89  if(theChannels[ii]->HasAnyData(i))
90  {
91  running[i] +=theChannels[ii]->GetWeightedXsec(aThermalE.GetThermalEnergy(aTrack,
92  theChannels[ii]->GetN(i),
93  theChannels[ii]->GetZ(i),
94  aTrack.GetMaterial()->GetTemperature()),
95  i);
96  }
97  }
98  }
99  G4int isotope=nChannels-1;
100  G4double random=G4UniformRand();
101  for(i=0;i<numberOfIsos; i++)
102  {
103  isotope = i;
104  //if(random<running[i]/running[numberOfIsos-1]) break;
105  if(running[numberOfIsos-1] != 0) if(random<running[i]/running[numberOfIsos-1]) break;
106  }
107  delete [] running;
108 
109  // decide on the channel
110  running = new G4double[nChannels];
111  running[0]=0;
112  G4int targA=-1; // For production of unChanged
113  G4int targZ=-1;
114  for(i=0; i<nChannels; i++)
115  {
116  if(i!=0) running[i] = running[i-1];
117  if(theChannels[i]->HasAnyData(isotope))
118  {
119  running[i] += theChannels[i]->GetFSCrossSection(aThermalE.GetThermalEnergy(aTrack,
120  theChannels[i]->GetN(isotope),
121  theChannels[i]->GetZ(isotope),
122  aTrack.GetMaterial()->GetTemperature()),
123  isotope);
124  targA=(G4int)theChannels[i]->GetN(isotope); //Will be simply used the last valid value
125  targZ=(G4int)theChannels[i]->GetZ(isotope);
126  }
127  }
128 
129  //TK120607
130  if ( running[nChannels-1] == 0 )
131  {
132  //This happened usually by the miss match between the cross section data and model
133  if ( targA == -1 && targZ == -1 ) {
134  throw G4HadronicException(__FILE__, __LINE__, "NeutronHP model encounter lethal discrepancy with cross section data");
135  }
136 
137  //TK121106
138  G4cout << "Warning from NeutronHP: could not find proper reaction channel. This may cause by inconsistency between cross section and model. Unchanged final states are returned." << G4endl;
139  unChanged.Clear();
140 
141  //For Ep Check create unchanged final state including rest target
142  G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon ( targZ , targA , 0.0 );
143  G4DynamicParticle* targ_dp = new G4DynamicParticle( targ_pd , G4ThreeVector(1,0,0), 0.0 );
145  unChanged.SetMomentumChange(aTrack.Get4Momentum().vect() );
146  unChanged.AddSecondary(targ_dp);
147  //TK121106
150  delete [] running;
151  return &unChanged;
152  }
153  //TK120607
154 
155 
156  G4int lChan=0;
157  random=G4UniformRand();
158  for(i=0; i<nChannels; i++)
159  {
160  lChan = i;
161  if(running[nChannels-1] != 0) if(random<running[i]/running[nChannels-1]) break;
162  }
163  delete [] running;
164  return theChannels[lChan]->ApplyYourself(aTrack, isotope);
165  }
166 
167  void G4NeutronHPChannelList::Init(G4Element * anElement, const G4String & dirName)
168  {
169  theDir = dirName;
170 // G4cout << theDir << G4endl;
171  theElement = anElement;
172 // G4cout << theElement << G4endl;
173  ;
174  }
175 
177  const G4String & aName)
178  {
179  if(!allChannelsCreated)
180  {
181  if(nChannels!=0)
182  {
183  G4NeutronHPChannel ** theBuffer = new G4NeutronHPChannel * [nChannels+1];
184  G4int i;
185  for(i=0; i<nChannels; i++)
186  {
187  theBuffer[i] = theChannels[i];
188  }
189  delete [] theChannels;
190  theChannels = theBuffer;
191  }
192  else
193  {
195  }
196  G4String name;
197  name = aName+"/";
200  nChannels++;
201  }
202 
203  //110527TKDB Unnessary codes, Detected by gcc4.6 compiler
204  //G4bool result;
205  //result = theChannels[theInitCount]->Register(theFS);
207  theInitCount++;
208  }
G4bool Register(G4NeutronHPFinalState *theFS)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
CLHEP::Hep3Vector G4ThreeVector
G4double GetZ(G4int i)
void Init(G4Element *theElement, const G4String dirName)
static G4NeutronHPManager * GetInstance()
G4String name
Definition: TRTMaterials.hh:40
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:463
G4NeutronHPChannel ** theChannels
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
G4double GetN(G4int i)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
G4double GetKineticEnergy() const
const G4int n
const G4LorentzVector & Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
static G4ThreadLocal G4int trycounter
G4HadFinalState * ApplyYourself(const G4Element *theElement, const G4HadProjectile &aTrack)
void SetEnergyChange(G4double anEnergy)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
void Register(G4NeutronHPFinalState *theFS, const G4String &aName)
void SetMomentumChange(const G4ThreeVector &aV)