Geant4  10.02.p01
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 "G4ParticleHPFissionFS.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4ParticleHPManager.hh"
39 #include "G4Threading.hh"
40 
42  :G4HadronicInteraction("NeutronHPFission")
43  ,theFission(NULL)
44  ,numEle(0)
45  {
46  SetMinEnergy( 0.0 );
47  SetMaxEnergy( 20.*MeV );
48 /*
49  if(!getenv("G4NEUTRONHPDATA"))
50  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
51  dirName = getenv("G4NEUTRONHPDATA");
52  G4String tString = "/Fission";
53  dirName = dirName + tString;
54  numEle = G4Element::GetNumberOfElements();
55  //theFission = new G4ParticleHPChannel[numEle];
56 
57  //for (G4int i=0; i<numEle; i++)
58  //{
59  //if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
60  // if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
61  // {
62  // theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
63  // theFission[i].Register(&theFS);
64  // }
65  //}
66 
67  for ( G4int i = 0 ; i < numEle ; i++ )
68  {
69  theFission.push_back( new G4ParticleHPChannel );
70  if((*(G4Element::GetElementTable()))[i]->GetZ()>87) //TK modified for ENDF-VII
71  {
72  (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
73  (*theFission[i]).Register(&theFS);
74  }
75  }
76 */
77  }
78 
80  {
81  //delete [] theFission;
82  for ( std::vector<G4ParticleHPChannel*>::iterator
83  it = theFission->begin() ; it != theFission->end() ; it++ )
84  {
85  delete *it;
86  }
87  theFission->clear();
88  }
89 
92  {
93 
95  const G4Material * theMaterial = aTrack.GetMaterial();
96  G4int n = theMaterial->GetNumberOfElements();
97  G4int index = theMaterial->GetElement(0)->GetIndex();
98  if(n!=1)
99  {
100  G4double* xSec = new G4double[n];
101  G4double sum=0;
102  G4int i;
103  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
104  G4double rWeight;
105  G4ParticleHPThermalBoost aThermalE;
106  for (i=0; i<n; i++)
107  {
108  index = theMaterial->GetElement(i)->GetIndex();
109  rWeight = NumAtomsPerVolume[i];
110  xSec[i] = ((*theFission)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
111  theMaterial->GetElement(i),
112  theMaterial->GetTemperature()));
113  xSec[i] *= rWeight;
114  sum+=xSec[i];
115  }
116  G4double random = G4UniformRand();
117  G4double running = 0;
118  for (i=0; i<n; i++)
119  {
120  running += xSec[i];
121  index = theMaterial->GetElement(i)->GetIndex();
122  //if(random<=running/sum) break;
123  if( sum == 0 || random <= running/sum ) break;
124  }
125  delete [] xSec;
126  }
127  //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
128  G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack,-2);
129 
130  //Overwrite target parameters
131  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
132  const G4Element* target_element = (*G4Element::GetElementTable())[index];
133  const G4Isotope* target_isotope=NULL;
134  G4int iele = target_element->GetNumberOfIsotopes();
135  for ( G4int j = 0 ; j != iele ; j++ ) {
136  target_isotope=target_element->GetIsotope( j );
137  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
138  }
139  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
140  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
141  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
142  aNucleus.SetIsotope( target_isotope );
143 
145  return result;
146  }
147 
148 const std::pair<G4double, G4double> G4ParticleHPFission::GetFatalEnergyCheckLevels() const
149 {
150  // max energy non-conservation is mass of heavy nucleus
151  //return std::pair<G4double, G4double>(5*perCent,250*GeV);
152  return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
153 }
154 /*
155 void G4ParticleHPFission::addChannelForNewElement()
156 {
157  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
158  {
159  theFission.push_back( new G4ParticleHPChannel );
160  if ( (*(G4Element::GetElementTable()))[i]->GetZ() > 87 ) //TK modified for ENDF-VII
161  {
162  G4cout << "G4ParticleHPFission Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
163  (*theFission[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
164  (*theFission[i]).Register(&theFS);
165  }
166  }
167  numEle = (G4int)G4Element::GetNumberOfElements();
168 }
169 */
170 
172 {
174 }
176 {
178 }
180 {
181 
183 
184  theFission = hpmanager->GetFissionFinalStates();
185 
186  if ( G4Threading::IsMasterThread() ) {
187 
188  if ( theFission == NULL ) theFission = new std::vector<G4ParticleHPChannel*>;
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 G4ParticleHPChannel );
205  if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
206  ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
207  ((*theFission)[i])->Register( new G4ParticleHPFissionFS );
208  }
209  }
210 
212 
213  }
214 
216 }
217 
218 void G4ParticleHPFission::ModelDescription(std::ostream& outFile) const
219 {
220  outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n";
221 }
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:86
static const double MeV
Definition: G4SIunits.hh:211
void RegisterFissionFinalStates(std::vector< G4ParticleHPChannel * > *val)
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:97
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:329
const G4int n
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates()
std::vector< G4ParticleHPChannel * > * theFission
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
const G4Material * GetMaterial() const
G4bool IsMasterThread()
Definition: G4Threading.cc:130
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
virtual void ModelDescription(std::ostream &outFile) const
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:211
G4int GetVerboseLevel() const