Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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  if ( theElastic != NULL ) {
78  for ( std::vector<G4ParticleHPChannel*>::iterator
79  it = theElastic->begin() ; it != theElastic->end() ; it++ ) {
80  delete *it;
81  }
82  theElastic->clear();
83  }
84  }
85 
87 
89  {
90 
91  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
92 
94  const G4Material * theMaterial = aTrack.GetMaterial();
95  G4int n = theMaterial->GetNumberOfElements();
96  G4int index = theMaterial->GetElement(0)->GetIndex();
97  if(n!=1)
98  {
99  G4int i;
100  G4double* xSec = new G4double[n];
101  G4double sum=0;
102  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
103  G4double rWeight;
104  G4ParticleHPThermalBoost aThermalE;
105  for (i=0; i<n; i++)
106  {
107  index = theMaterial->GetElement(i)->GetIndex();
108  rWeight = NumAtomsPerVolume[i];
109  //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
110  xSec[i] = ((*theElastic)[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  delete [] xSec;
126  // it is element-wise initialised.
127  }
128  //G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
129  G4HadFinalState* finalState = ((*theElastic)[index])->ApplyYourself(aTrack);
130  if (overrideSuspension) finalState->SetStatusChange(isAlive);
131 
132  //Overwrite target parameters
133  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
134  const G4Element* target_element = (*G4Element::GetElementTable())[index];
135  const G4Isotope* target_isotope=NULL;
136  G4int iele = target_element->GetNumberOfIsotopes();
137  for ( G4int j = 0 ; j != iele ; j++ ) {
138  target_isotope=target_element->GetIsotope( j );
139  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
140  }
141  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
142  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
143  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
144  aNucleus.SetIsotope( target_isotope );
145 
147  return finalState;
148  }
149 
150 const std::pair<G4double, G4double> G4ParticleHPElastic::GetFatalEnergyCheckLevels() const
151 {
152  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
153  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
154 }
155 
156 /*
157 void G4ParticleHPElastic::addChannelForNewElement()
158 {
159  G4ParticleHPElasticFS* theFS = new G4ParticleHPElasticFS;
160  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
161  {
162  G4cout << "G4ParticleHPElastic Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
163  theElastic.push_back( new G4ParticleHPChannel );
164  (*theElastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
165  while(!(*theElastic[i]).Register(theFS)) ;
166  }
167  delete theFS;
168  numEle = (G4int)G4Element::GetNumberOfElements();
169 }
170 */
171 
173 {
175 }
177 {
179 }
181 {
182 
184 
185  theElastic = hpmanager->GetElasticFinalStates();
186 
187  if ( G4Threading::IsMasterThread() ) {
188 
189  if ( theElastic == NULL ) theElastic = new std::vector<G4ParticleHPChannel*>;
190 
191  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
192 
193  if ( theElastic->size() == G4Element::GetNumberOfElements() ) {
195  return;
196  }
197 
199  if(!getenv("G4NEUTRONHPDATA"))
200  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
201  dirName = getenv("G4NEUTRONHPDATA");
202  G4String tString = "/Elastic";
203  dirName = dirName + tString;
204  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ ) {
205  theElastic->push_back( new G4ParticleHPChannel );
206  ((*theElastic)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
207  //while(!((*theElastic)[i])->Register(theFS)) ;
208  ((*theElastic)[i])->Register(theFS) ;
209  }
210  delete theFS;
211  hpmanager->RegisterElasticFinalStates( theElastic );
212 
213  }
215 }
216 void G4ParticleHPElastic::ModelDescription(std::ostream& outFile) const
217 {
218  outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for inelastic reaction of neutrons below 20MeV\n";
219 }
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
void RegisterElasticFinalStates(std::vector< G4ParticleHPChannel * > *val)
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:90
void BuildPhysicsTable(const G4ParticleDefinition &)
static constexpr double perCent
Definition: G4SIunits.hh:332
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:405
size_t GetIndex() const
Definition: G4Element.hh:182
static constexpr double eV
Definition: G4SIunits.hh:215
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
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:170
G4double GetTemperature() const
Definition: G4Material.hh:182
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
G4bool IsMasterThread()
Definition: G4Threading.cc:146
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4int GetVerboseLevel() const
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:212