Geant4  10.01.p03
G4ParticleHPFission.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 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
34 //
35 #include "G4ParticleHPFission.hh"
36 #include "G4SystemOfUnits.hh"
37 
38 #include "G4ParticleHPManager.hh"
39 
41  :G4HadronicInteraction("NeutronHPFission")
42  {
43  SetMinEnergy( 0.0 );
44  SetMaxEnergy( 20.*MeV );
45  if(!getenv("G4NEUTRONHPDATA"))
46  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
47  dirName = getenv("G4NEUTRONHPDATA");
48  G4String tString = "/Fission";
49  dirName = dirName + tString;
51  //theFission = new G4ParticleHPChannel[numEle];
52 
53  //for (G4int i=0; i<numEle; i++)
54  //{
55  //if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
56  // if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
57  // {
58  // theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
59  // theFission[i].Register(&theFS);
60  // }
61  //}
62 
63  for ( G4int i = 0 ; i < numEle ; i++ )
64  {
65  theFission.push_back( new G4ParticleHPChannel );
66  if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
67  {
69  (*theFission[i]).Register(&theFS);
70  }
71  }
72  }
73 
75  {
76  //delete [] theFission;
77  for ( std::vector<G4ParticleHPChannel*>::iterator
78  it = theFission.begin() ; it != theFission.end() ; it++ )
79  {
80  delete *it;
81  }
82  theFission.clear();
83  }
84 
87  {
88 
90 
92  const G4Material * theMaterial = aTrack.GetMaterial();
93  G4int n = theMaterial->GetNumberOfElements();
94  G4int index = theMaterial->GetElement(0)->GetIndex();
95  if(n!=1)
96  {
97  xSec = new G4double[n];
98  G4double sum=0;
99  G4int i;
100  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
101  G4double rWeight;
102  G4ParticleHPThermalBoost aThermalE;
103  for (i=0; i<n; i++)
104  {
105  index = theMaterial->GetElement(i)->GetIndex();
106  rWeight = NumAtomsPerVolume[i];
107  xSec[i] = (*theFission[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
108  theMaterial->GetElement(i),
109  theMaterial->GetTemperature()));
110  xSec[i] *= rWeight;
111  sum+=xSec[i];
112  }
113  G4double random = G4UniformRand();
114  G4double running = 0;
115  for (i=0; i<n; i++)
116  {
117  running += xSec[i];
118  index = theMaterial->GetElement(i)->GetIndex();
119  //if(random<=running/sum) break;
120  if( sum == 0 || random <= running/sum ) break;
121  }
122  delete [] xSec;
123  }
124  //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
125  G4HadFinalState* result = (*theFission[index]).ApplyYourself(aTrack,-2);
126 
127  //Overwrite target parameters
128  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
129  const G4Element* target_element = (*G4Element::GetElementTable())[index];
130  const G4Isotope* target_isotope=NULL;
131  G4int iele = target_element->GetNumberOfIsotopes();
132  for ( G4int j = 0 ; j != iele ; j++ ) {
133  target_isotope=target_element->GetIsotope( j );
134  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
135  }
136  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
137  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
138  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
139  aNucleus.SetIsotope( target_isotope );
140 
142  return result;
143  }
144 
145 const std::pair<G4double, G4double> G4ParticleHPFission::GetFatalEnergyCheckLevels() const
146 {
147  // max energy non-conservation is mass of heavy nucleus
148  //return std::pair<G4double, G4double>(5*perCent,250*GeV);
149  return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
150 }
151 
152 
153 
155 {
156  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
157  {
158  theFission.push_back( new G4ParticleHPChannel );
159  if ( (*(G4Element::GetElementTable()))[i]->GetZ() > 87 ) //TK modified for ENDF-VII
160  {
161  G4cout << "G4ParticleHPFission Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
163  (*theFission[i]).Register(&theFS);
164  }
165  }
167 }
168 
170 {
172 }
174 {
176 }
static G4ParticleHPManager * GetInstance()
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
G4ParticleHPFissionFS theFS
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
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
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
G4int GetVerboseLevel() const
std::vector< G4ParticleHPChannel * > theFission