Geant4  10.01.p02
G4NeutronHPFission.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 // 08-08-06 delete unnecessary and harmed declaration; Bug Report[857]
32 //
33 #include "G4NeutronHPFission.hh"
34 #include "G4SystemOfUnits.hh"
35 
36 #include "G4NeutronHPManager.hh"
37 
39  :G4HadronicInteraction("NeutronHPFission")
40  ,theFission(NULL)
41  ,numEle(0)
42  {
43  SetMinEnergy( 0.0 );
44  SetMaxEnergy( 20.*MeV );
45 /*
46  if(!getenv("G4NEUTRONHPDATA"))
47  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
48  dirName = getenv("G4NEUTRONHPDATA");
49  G4String tString = "/Fission";
50  dirName = dirName + tString;
51  numEle = G4Element::GetNumberOfElements();
52  //theFission = new G4NeutronHPChannel[numEle];
53 
54  //for (G4int i=0; i<numEle; i++)
55  //{
56  //if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
57  // if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
58  // {
59  // theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
60  // theFission[i].Register(&theFS);
61  // }
62  //}
63 
64  for ( G4int i = 0 ; i < numEle ; i++ )
65  {
66  theFission.push_back( new G4NeutronHPChannel );
67  if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
68  {
69  (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
70  (*theFission[i]).Register(&theFS);
71  }
72  }
73 */
74  }
75 
77  {
78  //delete [] theFission;
79  for ( std::vector<G4NeutronHPChannel*>::iterator
80  it = (*theFission).begin() ; it != (*theFission).end() ; it++ )
81  {
82  delete *it;
83  }
84  (*theFission).clear();
85  }
86 
89  {
90 
91  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
92 
94  const G4Material * theMaterial = aTrack.GetMaterial();
95  G4int n = theMaterial->GetNumberOfElements();
96  G4int index = theMaterial->GetElement(0)->GetIndex();
97  if(n!=1)
98  {
99  xSec = new G4double[n];
100  G4double sum=0;
101  G4int i;
102  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
103  G4double rWeight;
104  G4NeutronHPThermalBoost aThermalE;
105  for (i=0; i<n; i++)
106  {
107  index = theMaterial->GetElement(i)->GetIndex();
108  rWeight = NumAtomsPerVolume[i];
109  xSec[i] = (*(*theFission)[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
110  theMaterial->GetElement(i),
111  theMaterial->GetTemperature()));
112  xSec[i] *= rWeight;
113  sum+=xSec[i];
114  }
115  G4double random = G4UniformRand();
116  G4double running = 0;
117  for (i=0; i<n; i++)
118  {
119  running += xSec[i];
120  index = theMaterial->GetElement(i)->GetIndex();
121  //if(random<=running/sum) break;
122  if( sum == 0 || random <= running/sum ) break;
123  }
124  delete [] xSec;
125  }
126  //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
127  G4HadFinalState* result = (*(*theFission)[index]).ApplyYourself(aTrack,-2);
128 
129  //Overwrite target parameters
130  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
131  const G4Element* target_element = (*G4Element::GetElementTable())[index];
132  const G4Isotope* target_isotope=NULL;
133  G4int iele = target_element->GetNumberOfIsotopes();
134  for ( G4int j = 0 ; j != iele ; j++ ) {
135  target_isotope=target_element->GetIsotope( j );
136  if ( target_isotope->GetN() == G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
137  }
138  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
139  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
140  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
141  aNucleus.SetIsotope( target_isotope );
142 
144  return result;
145  }
146 
147 const std::pair<G4double, G4double> G4NeutronHPFission::GetFatalEnergyCheckLevels() const
148 {
149  // max energy non-conservation is mass of heavy nucleus
150  //return std::pair<G4double, G4double>(5*perCent,250*GeV);
151  return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
152 }
153 
154 
155 /*
156 void G4NeutronHPFission::addChannelForNewElement()
157 {
158  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
159  {
160  theFission.push_back( new G4NeutronHPChannel );
161  if ( (*(G4Element::GetElementTable()))[i]->GetZ() > 87 ) //TK modified for ENDF-VII
162  {
163  G4cout << "G4NeutronHPFission Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
164  (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
165  (*theFission[i]).Register(&theFS);
166  }
167  }
168  numEle = (G4int)G4Element::GetNumberOfElements();
169 }
170 */
171 
173 {
175 }
177 {
179 }
181 {
183 
184  theFission = hpmanager->GetFissionFinalStates();
185 
186  if ( !G4Threading::IsWorkerThread() ) {
187 
188  if ( theFission == NULL ) theFission = new std::vector<G4NeutronHPChannel*>;
189 
190  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
191 
192  if ( theFission->size() == G4Element::GetNumberOfElements() ) {
194  return;
195  }
196 
197  if ( !getenv("G4NEUTRONHPDATA") )
198  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
199  dirName = getenv("G4NEUTRONHPDATA");
200  G4String tString = "/Fission";
201  dirName = dirName + tString;
202 
203  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
204  (*theFission).push_back( new G4NeutronHPChannel );
205  if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
206  (*(*theFission)[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
207  (*(*theFission)[i]).Register(&theFS);
208  }
209  }
210 
212 
213  }
215 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
static const double MeV
Definition: G4SIunits.hh:193
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
void RegisterFissionFinalStates(std::vector< G4NeutronHPChannel * > *val)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
G4int GetVerboseLevel() const
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
void Register(T *inst)
Definition: G4AutoDelete.hh:65
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:93
std::vector< G4NeutronHPChannel * > * theFission
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
size_t GetIndex() const
Definition: G4Element.hh:181
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
G4bool IsWorkerThread()
Definition: G4Threading.cc:128
void SetVerboseLevel(G4int i)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
std::vector< G4NeutronHPChannel * > * GetFissionFinalStates()
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4NeutronHPFissionFS theFS
G4double GetTemperature() const
Definition: G4Material.hh:182
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
void BuildPhysicsTable(const G4ParticleDefinition &)