Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPElastic.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 // 080319 Compilation warnings - gcc-4.3.0 fix by T. Koi
32 //
33 #include "G4NeutronHPElastic.hh"
34 #include "G4SystemOfUnits.hh"
35 #include "G4NeutronHPElasticFS.hh"
36 #include "G4NeutronHPManager.hh"
37 
39  :G4HadronicInteraction("NeutronHPElastic")
40  {
41  overrideSuspension = false;
43  if(!getenv("G4NEUTRONHPDATA"))
44  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
45  dirName = getenv("G4NEUTRONHPDATA");
46  G4String tString = "/Elastic";
47  dirName = dirName + tString;
48 // G4cout <<"G4NeutronHPElastic::G4NeutronHPElastic testit "<<dirName<<G4endl;
50  //theElastic = new G4NeutronHPChannel[numEle];
51  //for (G4int i=0; i<numEle; i++)
52  //{
53  // theElastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
54  // while(!theElastic[i].Register(theFS)) ;
55  //}
56  for ( G4int i = 0 ; i < numEle ; i++ )
57  {
58  theElastic.push_back( new G4NeutronHPChannel );
59  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
60  while(!(*theElastic[i]).Register(theFS)) ;
61  }
62  delete theFS;
63  SetMinEnergy(0.*eV);
64  SetMaxEnergy(20.*MeV);
65  }
66 
68  {
69  //delete [] theElastic;
70  for ( std::vector<G4NeutronHPChannel*>::iterator
71  it = theElastic.begin() ; it != theElastic.end() ; it++ )
72  {
73  delete *it;
74  }
75  theElastic.clear();
76  }
77 
79 
81  {
82 
83  if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
84 
86  const G4Material * theMaterial = aTrack.GetMaterial();
87  G4int n = theMaterial->GetNumberOfElements();
88  G4int index = theMaterial->GetElement(0)->GetIndex();
89  if(n!=1)
90  {
91  G4int i;
92  xSec = new G4double[n];
93  G4double sum=0;
94  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
95  G4double rWeight;
96  G4NeutronHPThermalBoost aThermalE;
97  for (i=0; i<n; i++)
98  {
99  index = theMaterial->GetElement(i)->GetIndex();
100  rWeight = NumAtomsPerVolume[i];
101  //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
102  xSec[i] = (*theElastic[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
103  theMaterial->GetElement(i),
104  theMaterial->GetTemperature()));
105  xSec[i] *= rWeight;
106  sum+=xSec[i];
107  }
108  G4double random = G4UniformRand();
109  G4double running = 0;
110  for (i=0; i<n; i++)
111  {
112  running += xSec[i];
113  index = theMaterial->GetElement(i)->GetIndex();
114  //if(random<=running/sum) break;
115  if( sum == 0 || random <= running/sum ) break;
116  }
117  delete [] xSec;
118  // it is element-wise initialised.
119  }
120  //G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
121  G4HadFinalState* finalState = (*theElastic[index]).ApplyYourself(aTrack);
122  if (overrideSuspension) finalState->SetStatusChange(isAlive);
123  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
125  return finalState;
126  }
127 
128 const std::pair<G4double, G4double> G4NeutronHPElastic::GetFatalEnergyCheckLevels() const
129 {
130  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
131  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
132 }
133 
134 void G4NeutronHPElastic::addChannelForNewElement()
135 {
137  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
138  {
139  G4cout << "G4NeutronHPElastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
140  theElastic.push_back( new G4NeutronHPChannel );
141  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
142  while(!(*theElastic[i]).Register(theFS)) ;
143  }
144  delete theFS;
146 }
147 
149 {
151 }
153 {
155 }
void Init()
Definition: G4IonTable.cc:89
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static G4NeutronHPManager * GetInstance()
G4int GetVerboseLevel() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:200
int G4int
Definition: G4Types.hh:78
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
void SetStatusChange(G4HadFinalStateStatus aS)
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
size_t GetIndex() const
Definition: G4Element.hh:181
const G4int n
void SetVerboseLevel(G4int i)
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