Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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  if ( theFission != NULL ) {
83  for ( std::vector<G4ParticleHPChannel*>::iterator
84  it = theFission->begin() ; it != theFission->end() ; it++ ) {
85  delete *it;
86  }
87  theFission->clear();
88  }
89  }
90 
93  {
94 
96  const G4Material * theMaterial = aTrack.GetMaterial();
97  G4int n = theMaterial->GetNumberOfElements();
98  G4int index = theMaterial->GetElement(0)->GetIndex();
99  if(n!=1)
100  {
101  G4double* xSec = new G4double[n];
102  G4double sum=0;
103  G4int i;
104  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
105  G4double rWeight;
106  G4ParticleHPThermalBoost aThermalE;
107  for (i=0; i<n; i++)
108  {
109  index = theMaterial->GetElement(i)->GetIndex();
110  rWeight = NumAtomsPerVolume[i];
111  xSec[i] = ((*theFission)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
112  theMaterial->GetElement(i),
113  theMaterial->GetTemperature()));
114  xSec[i] *= rWeight;
115  sum+=xSec[i];
116  }
117  G4double random = G4UniformRand();
118  G4double running = 0;
119  for (i=0; i<n; i++)
120  {
121  running += xSec[i];
122  index = theMaterial->GetElement(i)->GetIndex();
123  //if(random<=running/sum) break;
124  if( sum == 0 || random <= running/sum ) break;
125  }
126  delete [] xSec;
127  }
128  //return theFission[index].ApplyYourself(aTrack); //-2:Marker for Fission
129  G4HadFinalState* result = ((*theFission)[index])->ApplyYourself(aTrack,-2);
130 
131  //Overwrite target parameters
132  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
133  const G4Element* target_element = (*G4Element::GetElementTable())[index];
134  const G4Isotope* target_isotope=NULL;
135  G4int iele = target_element->GetNumberOfIsotopes();
136  for ( G4int j = 0 ; j != iele ; j++ ) {
137  target_isotope=target_element->GetIsotope( j );
138  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
139  }
140  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
141  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
142  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
143  aNucleus.SetIsotope( target_isotope );
144 
146  return result;
147  }
148 
149 const std::pair<G4double, G4double> G4ParticleHPFission::GetFatalEnergyCheckLevels() const
150 {
151  // max energy non-conservation is mass of heavy nucleus
152  //return std::pair<G4double, G4double>(5*perCent,250*GeV);
153  return std::pair<G4double, G4double>(5*perCent,DBL_MAX);
154 }
155 /*
156 void G4ParticleHPFission::addChannelForNewElement()
157 {
158  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
159  {
160  theFission.push_back( new G4ParticleHPChannel );
161  if ( (*(G4Element::GetElementTable()))[i]->GetZ() > 87 ) //TK modified for ENDF-VII
162  {
163  G4cout << "G4ParticleHPFission 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 {
182 
184 
185  theFission = hpmanager->GetFissionFinalStates();
186 
187  if ( G4Threading::IsMasterThread() ) {
188 
189  if ( theFission == NULL ) theFission = new std::vector<G4ParticleHPChannel*>;
190 
191  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
192 
193  if ( theFission->size() == G4Element::GetNumberOfElements() ) {
195  return;
196  }
197 
198  if ( !getenv("G4NEUTRONHPDATA") )
199  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
200  dirName = getenv("G4NEUTRONHPDATA");
201  G4String tString = "/Fission";
202  dirName = dirName + tString;
203 
204  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
205  theFission->push_back( new G4ParticleHPChannel );
206  if ((*(G4Element::GetElementTable()))[i]->GetZ()>87) { //TK modified for ENDF-VII
207  ((*theFission)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
208  ((*theFission)[i])->Register( new G4ParticleHPFissionFS );
209  }
210  }
211 
212  hpmanager->RegisterFissionFinalStates( theFission );
213 
214  }
215 
217 }
218 
219 void G4ParticleHPFission::ModelDescription(std::ostream& outFile) const
220 {
221  outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for induced fission reaction of neutrons below 20MeV\n";
222 }
G4double G4ParticleHPJENDLHEData::G4double result
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:90
static constexpr double perCent
Definition: G4SIunits.hh:332
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:405
size_t GetIndex() const
Definition: G4Element.hh:182
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
std::vector< G4ParticleHPChannel * > * GetFissionFinalStates()
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4double GetTemperature() const
Definition: G4Material.hh:182
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
G4bool IsMasterThread()
Definition: G4Threading.cc:146
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
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:212
G4int GetVerboseLevel() const