Geant4  10.01.p02
G4ParticleHPElastic.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 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
34 //
35 #include "G4ParticleHPElastic.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "G4ParticleHPElasticFS.hh"
38 #include "G4ParticleHPManager.hh"
39 
41  :G4HadronicInteraction("NeutronHPElastic")
42  {
43  overrideSuspension = false;
45  if(!getenv("G4NEUTRONHPDATA"))
46  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
47  dirName = getenv("G4NEUTRONHPDATA");
48  G4String tString = "/Elastic";
49  dirName = dirName + tString;
50 // G4cout <<"G4ParticleHPElastic::G4ParticleHPElastic testit "<<dirName<<G4endl;
52  //theElastic = new G4ParticleHPChannel[numEle];
53  //for (G4int i=0; i<numEle; i++)
54  //{
55  // theElastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
56  // while(!theElastic[i].Register(theFS)) ;
57  //}
58  for ( G4int i = 0 ; i < numEle ; i++ )
59  {
60  theElastic.push_back( new G4ParticleHPChannel );
62  while(!(*theElastic[i]).Register(theFS)) ;
63  }
64  delete theFS;
65  SetMinEnergy(0.*eV);
66  SetMaxEnergy(20.*MeV);
67  }
68 
70  {
71  //delete [] theElastic;
72  for ( std::vector<G4ParticleHPChannel*>::iterator
73  it = theElastic.begin() ; it != theElastic.end() ; it++ )
74  {
75  delete *it;
76  }
77  theElastic.clear();
78  }
79 
81 
83  {
84 
86 
88  const G4Material * theMaterial = aTrack.GetMaterial();
89  G4int n = theMaterial->GetNumberOfElements();
90  G4int index = theMaterial->GetElement(0)->GetIndex();
91  if(n!=1)
92  {
93  G4int i;
94  xSec = new G4double[n];
95  G4double sum=0;
96  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
97  G4double rWeight;
98  G4ParticleHPThermalBoost aThermalE;
99  for (i=0; i<n; i++)
100  {
101  index = theMaterial->GetElement(i)->GetIndex();
102  rWeight = NumAtomsPerVolume[i];
103  //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
104  xSec[i] = (*theElastic[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
105  theMaterial->GetElement(i),
106  theMaterial->GetTemperature()));
107  xSec[i] *= rWeight;
108  sum+=xSec[i];
109  }
110  G4double random = G4UniformRand();
111  G4double running = 0;
112  for (i=0; i<n; i++)
113  {
114  running += xSec[i];
115  index = theMaterial->GetElement(i)->GetIndex();
116  //if(random<=running/sum) break;
117  if( sum == 0 || random <= running/sum ) break;
118  }
119  delete [] xSec;
120  // it is element-wise initialised.
121  }
122  //G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
123  G4HadFinalState* finalState = (*theElastic[index]).ApplyYourself(aTrack);
124  if (overrideSuspension) finalState->SetStatusChange(isAlive);
125 
126  //Overwrite target parameters
127  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
128  const G4Element* target_element = (*G4Element::GetElementTable())[index];
129  const G4Isotope* target_isotope=NULL;
130  G4int iele = target_element->GetNumberOfIsotopes();
131  for ( G4int j = 0 ; j != iele ; j++ ) {
132  target_isotope=target_element->GetIsotope( j );
133  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
134  }
135  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
136  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
137  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
138  aNucleus.SetIsotope( target_isotope );
139 
141  return finalState;
142  }
143 
144 const std::pair<G4double, G4double> G4ParticleHPElastic::GetFatalEnergyCheckLevels() const
145 {
146  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
147  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
148 }
149 
151 {
153  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
154  {
155  G4cout << "G4ParticleHPElastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
156  theElastic.push_back( new G4ParticleHPChannel );
158  while(!(*theElastic[i]).Register(theFS)) ;
159  }
160  delete theFS;
162 }
163 
165 {
167 }
169 {
171 }
static G4ParticleHPManager * GetInstance()
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)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
std::vector< G4ParticleHPChannel * > 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
G4GLOB_DLL std::ostream G4cout
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
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
static const double eV
Definition: G4SIunits.hh:194
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4int GetVerboseLevel() const
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198