Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4NeutronHPCapture.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//
32#include "G4NeutronHPCapture.hh"
33#include "G4SystemOfUnits.hh"
36#include "G4ParticleTable.hh"
37#include "G4IonTable.hh"
38
39#include "G4NeutronHPManager.hh"
40
42 :G4HadronicInteraction("NeutronHPCapture")
43 {
44 SetMinEnergy( 0.0 );
45 SetMaxEnergy( 20.*MeV );
46// G4cout << "Capture : start of construction!!!!!!!!"<<G4endl;
47 if(!getenv("G4NEUTRONHPDATA"))
48 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
49 dirName = getenv("G4NEUTRONHPDATA");
50 G4String tString = "/Capture";
51 dirName = dirName + tString;
53// G4cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<<G4endl;
54// G4cout <<"Disname="<<dirName<<" numEle="<<numEle<<G4endl;
55 //theCapture = new G4NeutronHPChannel[numEle];
56// G4cout <<"G4NeutronHPChannel constructed"<<G4endl;
58 //for (G4int i=0; i<numEle; i++)
59 //{
60// // G4cout << "initializing theCapture "<<i<<" "<< numEle<<G4endl;
61 // theCapture[i].Init((*(G4Element::GetElementTable()))[i], dirName);
62 // theCapture[i].Register(theFS);
63 //}
64 for ( G4int i = 0 ; i < numEle ; i++ )
65 {
66 theCapture.push_back( new G4NeutronHPChannel );
67 (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
68 (*theCapture[i]).Register(theFS);
69 }
70 delete theFS;
71// G4cout << "-------------------------------------------------"<<G4endl;
72// G4cout << "Leaving G4NeutronHPCapture::G4NeutronHPCapture"<<G4endl;
73 }
74
76 {
77 //delete [] theCapture;
78// G4cout << "Leaving G4NeutronHPCapture::~G4NeutronHPCapture"<<G4endl;
79 for ( std::vector<G4NeutronHPChannel*>::iterator
80 ite = theCapture.begin() ; ite != theCapture.end() ; ite++ )
81 {
82 delete *ite;
83 }
84 theCapture.clear();
85 }
86
89 {
90
91 if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
92
94 if(getenv("NeutronHPCapture")) G4cout <<" ####### G4NeutronHPCapture called"<<G4endl;
95 const G4Material * theMaterial = aTrack.GetMaterial();
96 G4int n = theMaterial->GetNumberOfElements();
97 G4int index = theMaterial->GetElement(0)->GetIndex();
98 if(n!=1)
99 {
100 xSec = new G4double[n];
101 G4double sum=0;
102 G4int i;
103 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
104 G4double rWeight;
105 G4NeutronHPThermalBoost aThermalE;
106 for (i=0; i<n; i++)
107 {
108 index = theMaterial->GetElement(i)->GetIndex();
109 rWeight = NumAtomsPerVolume[i];
110 //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
111 xSec[i] = (*theCapture[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 if(i==n) i=std::max(0, n-1);
127 delete [] xSec;
128 }
129
130 //return theCapture[index].ApplyYourself(aTrack);
131 //G4HadFinalState* result = theCapture[index].ApplyYourself(aTrack);
132 G4HadFinalState* result = (*theCapture[index]).ApplyYourself(aTrack);
135 return result;
136 }
137
138const std::pair<G4double, G4double> G4NeutronHPCapture::GetFatalEnergyCheckLevels() const
139{
140 //return std::pair<G4double, G4double>(10*perCent,10*GeV);
141 return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
142}
143
144void G4NeutronHPCapture::addChannelForNewElement()
145{
147 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
148 {
149 G4cout << "G4NeutronHPCapture Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
150 theCapture.push_back( new G4NeutronHPChannel );
151 (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
152 (*theCapture[i]).Register(theFS);
153 }
154 delete theFS;
156}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
size_t GetIndex() const
Definition: G4Element.hh:182
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
const G4Material * GetMaterial() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4double GetTemperature() const
Definition: G4Material.hh:181
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:201
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
#define DBL_MAX
Definition: templates.hh:83