Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 "G4NeutronHPManager.hh"
34 #include "G4SystemOfUnits.hh"
35 #include "G4NeutronHPCaptureFS.hh"
36 #include "G4NeutronHPDeExGammas.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4IonTable.hh"
39 
41  :G4HadronicInteraction("NeutronHPCapture")
42  {
43  SetMinEnergy( 0.0 );
44  SetMaxEnergy( 20.*MeV );
45 // G4cout << "Capture : start of construction!!!!!!!!"<<G4endl;
46  if(!getenv("G4NEUTRONHPDATA"))
47  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
48  dirName = getenv("G4NEUTRONHPDATA");
49  G4String tString = "/Capture";
50  dirName = dirName + tString;
52 // G4cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<<G4endl;
53 // G4cout <<"Disname="<<dirName<<" numEle="<<numEle<<G4endl;
54  //theCapture = new G4NeutronHPChannel[numEle];
55 // G4cout <<"G4NeutronHPChannel constructed"<<G4endl;
57  //for (G4int i=0; i<numEle; i++)
58  //{
59 // // G4cout << "initializing theCapture "<<i<<" "<< numEle<<G4endl;
60  // theCapture[i].Init((*(G4Element::GetElementTable()))[i], dirName);
61  // theCapture[i].Register(theFS);
62  //}
63  for ( G4int i = 0 ; i < numEle ; i++ )
64  {
65  theCapture.push_back( new G4NeutronHPChannel );
66  (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
67  (*theCapture[i]).Register(theFS);
68  }
69  delete theFS;
70 // G4cout << "-------------------------------------------------"<<G4endl;
71 // G4cout << "Leaving G4NeutronHPCapture::G4NeutronHPCapture"<<G4endl;
72  }
73 
75  {
76  //delete [] theCapture;
77 // G4cout << "Leaving G4NeutronHPCapture::~G4NeutronHPCapture"<<G4endl;
78  for ( std::vector<G4NeutronHPChannel*>::iterator
79  ite = theCapture.begin() ; ite != theCapture.end() ; ite++ )
80  {
81  delete *ite;
82  }
83  theCapture.clear();
84  }
85 
88  {
89 
90  if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
91 
93  if(getenv("NeutronHPCapture")) G4cout <<" ####### G4NeutronHPCapture called"<<G4endl;
94  const G4Material * theMaterial = aTrack.GetMaterial();
95  G4int n = theMaterial->GetNumberOfElements();
96  G4int index = theMaterial->GetElement(0)->GetIndex();
97  if(n!=1)
98  {
99  xSec = new G4double[n];
100  G4double sum=0;
101  G4int i;
102  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
103  G4double rWeight;
104  G4NeutronHPThermalBoost aThermalE;
105  for (i=0; i<n; i++)
106  {
107  index = theMaterial->GetElement(i)->GetIndex();
108  rWeight = NumAtomsPerVolume[i];
109  //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
110  xSec[i] = (*theCapture[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  if(i==n) i=std::max(0, n-1);
126  delete [] xSec;
127  }
128 
129  //return theCapture[index].ApplyYourself(aTrack);
130  //G4HadFinalState* result = theCapture[index].ApplyYourself(aTrack);
131  G4HadFinalState* result = (*theCapture[index]).ApplyYourself(aTrack);
132  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
134  return result;
135  }
136 
137 const std::pair<G4double, G4double> G4NeutronHPCapture::GetFatalEnergyCheckLevels() const
138 {
139  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
140  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
141 }
142 
143 void G4NeutronHPCapture::addChannelForNewElement()
144 {
146  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
147  {
148  G4cout << "G4NeutronHPCapture Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
149  theCapture.push_back( new G4NeutronHPChannel );
150  (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
151  (*theCapture[i]).Register(theFS);
152  }
153  delete theFS;
155 }
156 
158 {
160 }
162 {
164 }
void Init()
Definition: G4IonTable.cc:89
static G4NeutronHPManager * GetInstance()
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
G4int GetVerboseLevel() const
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
void Register(T *inst)
Definition: G4AutoDelete.hh:65
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
size_t GetIndex() const
Definition: G4Element.hh:181
const G4int n
void SetVerboseLevel(G4int i)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void SetMaxEnergy(const G4double anEnergy)
float perCent
Definition: hepunit.py:239
G4double GetTemperature() const
Definition: G4Material.hh:180
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198