Geant4  10.01.p02
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 
38 #include "G4Threading.hh"
39 
41  :G4HadronicInteraction("NeutronHPElastic")
42  ,theElastic(NULL)
43  ,numEle(0)
44  {
45  overrideSuspension = false;
46 /*
47  G4NeutronHPElasticFS * theFS = new G4NeutronHPElasticFS;
48  if(!getenv("G4NEUTRONHPDATA"))
49  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
50  dirName = getenv("G4NEUTRONHPDATA");
51  G4String tString = "/Elastic";
52  dirName = dirName + tString;
53 // G4cout <<"G4NeutronHPElastic::G4NeutronHPElastic testit "<<dirName<<G4endl;
54  numEle = G4Element::GetNumberOfElements();
55  //theElastic = new G4NeutronHPChannel[numEle];
56  //for (G4int i=0; i<numEle; i++)
57  //{
58  // theElastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
59  // while(!theElastic[i].Register(theFS)) ;
60  //}
61  for ( G4int i = 0 ; i < numEle ; i++ )
62  {
63  theElastic.push_back( new G4NeutronHPChannel );
64  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
65  while(!(*theElastic[i]).Register(theFS)) ;
66  }
67  delete theFS;
68 */
69  SetMinEnergy(0.*eV);
70  SetMaxEnergy(20.*MeV);
71  }
72 
74  {
75  /*
76  //delete [] theElastic;
77  for ( std::vector<G4NeutronHPChannel*>::iterator
78  it = theElastic.begin() ; it != theElastic.end() ; it++ )
79  {
80  delete *it;
81  }
82  theElastic.clear();
83  */
84  }
85 
87 
89  {
90 
91 //G4cout << "TKDB::G4NeutronHPElastic::ApplyYourself this " << this << " numEle " << numEle << " &numEle " << &numEle << " G4Element::GetNumberOfElements() " <<G4Element::GetNumberOfElements() << G4endl;
92 
93  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
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  G4int i;
102  xSec = new G4double[n];
103  G4double sum=0;
104  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
105  G4double rWeight;
106  G4NeutronHPThermalBoost aThermalE;
107  for (i=0; i<n; i++)
108  {
109  index = theMaterial->GetElement(i)->GetIndex();
110  rWeight = NumAtomsPerVolume[i];
111  //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
112  xSec[i] = (*(*theElastic)[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
113  theMaterial->GetElement(i),
114  theMaterial->GetTemperature()));
115  xSec[i] *= rWeight;
116  sum+=xSec[i];
117  }
118  G4double random = G4UniformRand();
119  G4double running = 0;
120  for (i=0; i<n; i++)
121  {
122  running += xSec[i];
123  index = theMaterial->GetElement(i)->GetIndex();
124  //if(random<=running/sum) break;
125  if( sum == 0 || random <= running/sum ) break;
126  }
127  delete [] xSec;
128  // it is element-wise initialised.
129  }
130  //G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
131  G4HadFinalState* finalState = (*(*theElastic)[index]).ApplyYourself(aTrack);
132  if (overrideSuspension) finalState->SetStatusChange(isAlive);
133 
134  //Overwrite target parameters
135  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
136  const G4Element* target_element = (*G4Element::GetElementTable())[index];
137  const G4Isotope* target_isotope=NULL;
138  G4int iele = target_element->GetNumberOfIsotopes();
139  for ( G4int j = 0 ; j != iele ; j++ ) {
140  target_isotope=target_element->GetIsotope( j );
141  if ( target_isotope->GetN() == G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
142  }
143  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
144  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
145  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
146  aNucleus.SetIsotope( target_isotope );
147 
149  return finalState;
150  }
151 
152 const std::pair<G4double, G4double> G4NeutronHPElastic::GetFatalEnergyCheckLevels() const
153 {
154  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
155  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
156 }
157 
158 /*
159 void G4NeutronHPElastic::addChannelForNewElement()
160 {
161  G4NeutronHPElasticFS* theFS = new G4NeutronHPElasticFS;
162  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
163  {
164  G4cout << "G4NeutronHPElastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
165  theElastic.push_back( new G4NeutronHPChannel );
166  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
167  while(!(*theElastic[i]).Register(theFS)) ;
168  }
169  delete theFS;
170  numEle = (G4int)G4Element::GetNumberOfElements();
171 }
172 */
173 
175 {
177 }
179 {
181 }
182 
184 {
185 
187 
188  theElastic = hpmanager->GetElasticFinalStates();
189 
190  if ( !G4Threading::IsWorkerThread() ) {
191 
192  if ( theElastic == NULL ) theElastic = new std::vector<G4NeutronHPChannel*>;
193 
194  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
195 
196  if ( theElastic->size() == G4Element::GetNumberOfElements() ) {
198  return;
199  }
200 
202  if(!getenv("G4NEUTRONHPDATA"))
203  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
204  dirName = getenv("G4NEUTRONHPDATA");
205  G4String tString = "/Elastic";
206  dirName = dirName + tString;
207  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
208  (*theElastic).push_back( new G4NeutronHPChannel );
209  (*(*theElastic)[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
210  while(!(*(*theElastic)[i]).Register(theFS)) ;
211  }
212  delete theFS;
213 
215 
216  }
217 
219 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
static const double MeV
Definition: G4SIunits.hh:193
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4int GetVerboseLevel() const
std::vector< G4NeutronHPChannel * > * GetElasticFinalStates()
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
std::vector< G4NeutronHPChannel * > * theElastic
void SetStatusChange(G4HadFinalStateStatus aS)
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:93
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
G4bool IsWorkerThread()
Definition: G4Threading.cc:128
void SetVerboseLevel(G4int i)
static const double eV
Definition: G4SIunits.hh:194
void RegisterElasticFinalStates(std::vector< G4NeutronHPChannel * > *val)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4double GetTemperature() const
Definition: G4Material.hh:182
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
void BuildPhysicsTable(const G4ParticleDefinition &)