Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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  {
41  SetMinEnergy( 0.0 );
42  SetMaxEnergy( 20.*MeV );
43  if(!getenv("G4NEUTRONHPDATA"))
44  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
45  dirName = getenv("G4NEUTRONHPDATA");
46  G4String tString = "/Fission";
47  dirName = dirName + tString;
49  //theFission = new G4NeutronHPChannel[numEle];
50 
51  //for (G4int i=0; i<numEle; i++)
52  //{
53  //if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
54  // if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
55  // {
56  // theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
57  // theFission[i].Register(&theFS);
58  // }
59  //}
60 
61  for ( G4int i = 0 ; i < numEle ; i++ )
62  {
63  theFission.push_back( new G4NeutronHPChannel );
64  if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
65  {
66  (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
67  (*theFission[i]).Register(&theFS);
68  }
69  }
70  }
71 
73  {
74  //delete [] theFission;
75  for ( std::vector<G4NeutronHPChannel*>::iterator
76  it = theFission.begin() ; it != theFission.end() ; it++ )
77  {
78  delete *it;
79  }
80  theFission.clear();
81  }
82 
85  {
86 
87  if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
88 
90  const G4Material * theMaterial = aTrack.GetMaterial();
91  G4int n = theMaterial->GetNumberOfElements();
92  G4int index = theMaterial->GetElement(0)->GetIndex();
93  if(n!=1)
94  {
95  xSec = new G4double[n];
96  G4double sum=0;
97  G4int i;
98  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
99  G4double rWeight;
100  G4NeutronHPThermalBoost aThermalE;
101  for (i=0; i<n; i++)
102  {
103  index = theMaterial->GetElement(i)->GetIndex();
104  rWeight = NumAtomsPerVolume[i];
105  xSec[i] = (*theFission[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
106  theMaterial->GetElement(i),
107  theMaterial->GetTemperature()));
108  xSec[i] *= rWeight;
109  sum+=xSec[i];
110  }
111  G4double random = G4UniformRand();
112  G4double running = 0;
113  for (i=0; i<n; i++)
114  {
115  running += xSec[i];
116  index = theMaterial->GetElement(i)->GetIndex();
117  //if(random<=running/sum) break;
118  if( sum == 0 || random <= running/sum ) break;
119  }
120  delete [] xSec;
121  }
122  //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
123  G4HadFinalState* result = (*theFission[index]).ApplyYourself(aTrack,-2);
124  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
126  return result;
127  }
128 
129 const std::pair<G4double, G4double> G4NeutronHPFission::GetFatalEnergyCheckLevels() const
130 {
131  // max energy non-conservation is mass of heavy nucleus
132  //return std::pair<G4double, G4double>(5*perCent,250*GeV);
133  return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
134 }
135 
136 
137 
138 void G4NeutronHPFission::addChannelForNewElement()
139 {
140  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
141  {
142  theFission.push_back( new G4NeutronHPChannel );
143  if ( (*(G4Element::GetElementTable()))[i]->GetZ() > 87 ) //TK modified for ENDF-VII
144  {
145  G4cout << "G4NeutronHPFission Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
146  (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
147  (*theFission[i]).Register(&theFS);
148  }
149  }
151 }
152 
154 {
156 }
158 {
160 }
void Init()
Definition: G4IonTable.cc:89
static G4NeutronHPManager * GetInstance()
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
G4int GetVerboseLevel() const
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
void Register(T *inst)
Definition: G4AutoDelete.hh:65
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
size_t GetIndex() const
Definition: G4Element.hh:181
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
const G4int n
void SetVerboseLevel(G4int i)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void SetMaxEnergy(const G4double anEnergy)
float perCent
Definition: hepunit.py:239
G4double GetTemperature() const
Definition: G4Material.hh:180
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198