Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FissLib.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 //
27 // This software was developed by Lawrence Livermore National Laboratory.
28 //
29 // Redistribution and use in source and binary forms, with or without
30 // modification, are permitted provided that the following conditions are met:
31 //
32 // 1. Redistributions of source code must retain the above copyright notice,
33 // this list of conditions and the following disclaimer.
34 // 2. Redistributions in binary form must reproduce the above copyright notice,
35 // this list of conditions and the following disclaimer in the documentation
36 // and/or other materials provided with the distribution.
37 // 3. The name of the author may not be used to endorse or promote products
38 // derived from this software without specific prior written permission.
39 //
40 // THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS OR IMPLIED
41 // WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF
42 // MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO
43 // EVENT SHALL THE AUTHOR BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
44 // SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
45 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
46 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
47 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
48 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
49 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
50 //
51 // Copyright (c) 2006 The Regents of the University of California.
52 // All rights reserved.
53 // UCRL-CODE-224807
54 //
55 //
56 // $Id: G4FissLib.cc 67966 2013-03-13 09:38:38Z gcosmo $
57 //
58 // neutron_hp -- source file
59 // J.M. Verbeke, Jan-2007
60 // A prototype of the low energy neutron transport model.
61 //
62 #include "G4FissLib.hh"
63 #include "G4SystemOfUnits.hh"
64 
66  :xSec(0)
67 {
68  SetMinEnergy(0.0);
69  SetMaxEnergy(20.*MeV);
70  if(!getenv("G4NEUTRONHPDATA")) {
71  G4cout << "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files." << G4endl;
72  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
73  }
74  dirName = getenv("G4NEUTRONHPDATA");
75  G4String tString = "/Fission/";
76  dirName = dirName + tString;
78  theFission = new G4ParticleHPChannel[numEle];
79 
80  for (G4int i=0; i<numEle; i++)
81  {
82 // G4cout << "G4FissLib::G4FissLib(): element "<< i << " : " << (*(G4Element::GetElementTable()))[i]->GetZ()<< G4endl;
83  if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
84  {
85  theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
86  theFission[i].Register(&theLibrary);
87  }
88  }
89 }
90 
92 {
93  delete [] theFission;
94 }
95 
98 {
100 
101  const G4Material* theMaterial = aTrack.GetMaterial();
102  G4int n = theMaterial->GetNumberOfElements();
103  G4int index = theMaterial->GetElement(0)->GetIndex();
104 
105  if (n != 1) {
106  xSec = new G4double[n];
107  G4double sum = 0;
108  G4int i;
109  G4int imat;
110  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
111  G4double rWeight;
112  G4ParticleHPThermalBoost aThermalE;
113  for (i = 0; i < n; i++) {
114  imat = theMaterial->GetElement(i)->GetIndex();
115  rWeight = NumAtomsPerVolume[i];
116  xSec[i] = theFission[imat].GetXsec(aThermalE.GetThermalEnergy(aTrack,
117  theMaterial->GetElement(i),
118  theMaterial->GetTemperature()));
119  xSec[i] *= rWeight;
120  sum+=xSec[i];
121  }
122 
123  G4double random = G4UniformRand();
124  G4double running = 0;
125  for (i = 0; i < n; i++) {
126  running += xSec[i];
127  index = theMaterial->GetElement(i)->GetIndex();
128  if(random<=running/sum) break;
129  }
130  delete [] xSec;
131  }
132 
133  //return theFission[index].ApplyYourself(aTrack);
134  //Overwrite target parameters
135  G4HadFinalState* result = theFission[index].ApplyYourself(aTrack);
136  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
137  const G4Element* target_element = (*G4Element::GetElementTable())[index];
138  const G4Isotope* target_isotope=NULL;
139  G4int iele = target_element->GetNumberOfIsotopes();
140  for ( G4int j = 0 ; j != iele ; j++ ) {
141  target_isotope=target_element->GetIsotope( j );
142  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
143  }
144  aNucleus.SetIsotope( target_isotope );
146  return result;
147 }
148 
149 const std::pair<G4double, G4double> G4FissLib::GetFatalEnergyCheckLevels() const
150 {
151  // max energy non-conservation is mass of heavy nucleus (taken from G4LFission)
152  return std::pair<G4double, G4double>(5*perCent,250*GeV);
153 }
G4double G4ParticleHPJENDLHEData::G4double result
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
Definition: G4FissLib.cc:149
static constexpr double perCent
Definition: G4SIunits.hh:332
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
Definition: G4FissLib.cc:97
G4bool Register(G4ParticleHPFinalState *theFS)
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
size_t GetIndex() const
Definition: G4Element.hh:182
void Init(G4Element *theElement, const G4String dirName)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
~G4FissLib()
Definition: G4FissLib.cc:91
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetXsec(G4double energy)
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:212
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)