Geant4  10.02
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 #include "G4Threading.hh"
40 
42  :G4HadronicInteraction("NeutronHPElastic")
43  ,theElastic(NULL)
44  ,numEle(0)
45  {
46  overrideSuspension = false;
47 /*
48  G4ParticleHPElasticFS * theFS = new G4ParticleHPElasticFS;
49  if(!getenv("G4NEUTRONHPDATA"))
50  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
51  dirName = getenv("G4NEUTRONHPDATA");
52  G4String tString = "/Elastic";
53  dirName = dirName + tString;
54 // G4cout <<"G4ParticleHPElastic::G4ParticleHPElastic testit "<<dirName<<G4endl;
55  numEle = G4Element::GetNumberOfElements();
56  //theElastic = new G4ParticleHPChannel[numEle];
57  //for (G4int i=0; i<numEle; i++)
58  //{
59  // theElastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
60  // while(!theElastic[i].Register(theFS)) ;
61  //}
62  for ( G4int i = 0 ; i < numEle ; i++ )
63  {
64  theElastic.push_back( new G4ParticleHPChannel );
65  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
66  while(!(*theElastic[i]).Register(theFS)) ;
67  }
68  delete theFS;
69 */
70  SetMinEnergy(0.*eV);
71  SetMaxEnergy(20.*MeV);
72  }
73 
75  {
76  //delete [] theElastic;
77  /*
78  for ( std::vector<G4ParticleHPChannel*>::iterator
79  it = theElastic.begin() ; it != theElastic.end() ; it++ )
80  {
81  delete *it;
82  }
83  */
84  theElastic->clear();
85  }
86 
88 
90  {
91 
92  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
93 
95  const G4Material * theMaterial = aTrack.GetMaterial();
96  G4int n = theMaterial->GetNumberOfElements();
97  G4int index = theMaterial->GetElement(0)->GetIndex();
98  if(n!=1)
99  {
100  G4int i;
101  xSec = new G4double[n];
102  G4double sum=0;
103  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
104  G4double rWeight;
105  G4ParticleHPThermalBoost aThermalE;
106  for (i=0; i<n; i++)
107  {
108  index = theMaterial->GetElement(i)->GetIndex();
109  rWeight = NumAtomsPerVolume[i];
110  //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
111  xSec[i] = ((*theElastic)[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  delete [] xSec;
127  // it is element-wise initialised.
128  }
129  //G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
130  G4HadFinalState* finalState = ((*theElastic)[index])->ApplyYourself(aTrack);
131  if (overrideSuspension) finalState->SetStatusChange(isAlive);
132 
133  //Overwrite target parameters
134  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
135  const G4Element* target_element = (*G4Element::GetElementTable())[index];
136  const G4Isotope* target_isotope=NULL;
137  G4int iele = target_element->GetNumberOfIsotopes();
138  for ( G4int j = 0 ; j != iele ; j++ ) {
139  target_isotope=target_element->GetIsotope( j );
140  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
141  }
142  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
143  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
144  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
145  aNucleus.SetIsotope( target_isotope );
146 
148  return finalState;
149  }
150 
151 const std::pair<G4double, G4double> G4ParticleHPElastic::GetFatalEnergyCheckLevels() const
152 {
153  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
154  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
155 }
156 
157 /*
158 void G4ParticleHPElastic::addChannelForNewElement()
159 {
160  G4ParticleHPElasticFS* theFS = new G4ParticleHPElasticFS;
161  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
162  {
163  G4cout << "G4ParticleHPElastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
164  theElastic.push_back( new G4ParticleHPChannel );
165  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
166  while(!(*theElastic[i]).Register(theFS)) ;
167  }
168  delete theFS;
169  numEle = (G4int)G4Element::GetNumberOfElements();
170 }
171 */
172 
174 {
176 }
178 {
180 }
182 {
183 
185 
186  theElastic = hpmanager->GetElasticFinalStates();
187 
188  if ( G4Threading::IsMasterThread() ) {
189 
190  if ( theElastic == NULL ) theElastic = new std::vector<G4ParticleHPChannel*>;
191 
192  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
193 
194  if ( theElastic->size() == G4Element::GetNumberOfElements() ) {
196  return;
197  }
198 
200  if(!getenv("G4NEUTRONHPDATA"))
201  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
202  dirName = getenv("G4NEUTRONHPDATA");
203  G4String tString = "/Elastic";
204  dirName = dirName + tString;
205  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
206  theElastic->push_back( new G4ParticleHPChannel );
207  ((*theElastic)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
208  //while(!((*theElastic)[i])->Register(theFS)) ;
209  ((*theElastic)[i])->Register(theFS) ;
210  }
211  delete theFS;
213 
214  }
216 }
217 void G4ParticleHPElastic::ModelDescription(std::ostream& outFile) const
218 {
219  outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for inelastic reaction of neutrons below 20MeV\n";
220 }
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void RegisterElasticFinalStates(std::vector< G4ParticleHPChannel * > *val)
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
static const double MeV
Definition: G4SIunits.hh:211
void BuildPhysicsTable(const G4ParticleDefinition &)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
void SetStatusChange(G4HadFinalStateStatus aS)
void SetMinEnergy(G4double anEnergy)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
std::vector< G4ParticleHPChannel * > * theElastic
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:329
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:212
std::vector< G4ParticleHPChannel * > * GetElasticFinalStates()
virtual void ModelDescription(std::ostream &outFile) const
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
G4bool IsMasterThread()
Definition: G4Threading.cc:130
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4int GetVerboseLevel() const
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:211