Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPorLElastic.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 //
27 // 05-11-21 NeutronHP or Low Energy Parameterization Models
28 // Implemented by T. Koi (SLAC/SCCS)
29 // If NeutronHP data do not available for an element, then Low Energy
30 // Parameterization models handle the interactions of the element.
31 // 080319 Compilation warnings - gcc-4.3.0 fix by T. Koi
32 //
33 
34 // neutron_hp -- source file
35 // J.P. Wellisch, Nov-1996
36 // A prototype of the low energy neutron transport model.
37 //
38 #include "G4NeutronHPorLElastic.hh"
39 #include "G4SystemOfUnits.hh"
40 #include "G4NeutronHPElasticFS.hh"
41 
43  :G4HadronicInteraction("NeutronHPorLElastic")
44 {
45  overrideSuspension = false;
47  if(!getenv("G4NEUTRONHPDATA"))
48  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
49  dirName = getenv("G4NEUTRONHPDATA");
50  G4String tString = "/Elastic/";
51  dirName = dirName + tString;
52 // G4cout <<"G4NeutronHPorLElastic::G4NeutronHPorLElastic testit "<<dirName<<G4endl;
54  theElastic = new G4NeutronHPChannel[numEle];
55  unavailable_elements.clear();
56  for (G4int i=0; i<numEle; i++)
57  {
58  theElastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
59  //while(!theElastic[i].Register(theFS));
60  try { while(!theElastic[i].Register(theFS)) ; }
61  catch ( G4HadronicException )
62  {
63  unavailable_elements.insert ( (*(G4Element::GetElementTable()))[i]->GetName() );
64  }
65  }
66  delete theFS;
67  SetMinEnergy(0.*eV);
68  SetMaxEnergy(20.*MeV);
69  if ( unavailable_elements.size() > 0 )
70  {
71  std::set< G4String>::iterator it;
72  G4cout << "HP Elastic data are not available for thess elements "<< G4endl;
73  for ( it = unavailable_elements.begin() ; it != unavailable_elements.end() ; it++ )
74  G4cout << *it << G4endl;
75  G4cout << "Low Energy Parameterization Models will be used."<< G4endl;
76  }
77 
78  createXSectionDataSet();
79 }
80 
82 {
83  delete [] theElastic;
84  delete theDataSet;
85 }
86 
88 
90 {
91  const G4Material * theMaterial = aTrack.GetMaterial();
92  G4int n = theMaterial->GetNumberOfElements();
93  G4int index = theMaterial->GetElement(0)->GetIndex();
94  if(n!=1)
95  {
96  G4int i;
97  xSec = new G4double[n];
98  G4double sum=0;
99  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
100  G4double rWeight;
101  G4NeutronHPThermalBoost aThermalE;
102  for (i=0; i<n; i++)
103  {
104  index = theMaterial->GetElement(i)->GetIndex();
105  rWeight = NumAtomsPerVolume[i];
106  G4double x = aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i), theMaterial->GetTemperature());
107 
108  //xSec[i] = theElastic[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
109  // theMaterial->GetElement(i),
110  // theMaterial->GetTemperature()));
111  xSec[i] = theElastic[index].GetXsec(x);
112 
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  }
124  delete [] xSec;
125  // it is element-wise initialised.
126  }
127  G4HadFinalState* finalState = theElastic[index].ApplyYourself(aTrack);
128  if (overrideSuspension) finalState->SetStatusChange(isAlive);
129  return finalState;
130 }
131 
132 
133 
135 {
136  if ( unavailable_elements.find( name ) == unavailable_elements.end() )
137  return TRUE;
138  else
139  return FALSE;
140 }
141 
142 
143 
144 void G4NeutronHPorLElastic::createXSectionDataSet()
145 {
146  theDataSet = new G4NeutronHPorLElasticData ( theElastic , &unavailable_elements );
147 }
148 
149 const std::pair<G4double, G4double> G4NeutronHPorLElastic::GetFatalEnergyCheckLevels() const
150 {
151  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
152  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
153 }