Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronElasticXS.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 // $Id$
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4NeutronElasticXS
34 //
35 // Author Ivantchenko, Geant4, 3-Aug-09
36 //
37 // Modifications:
38 //
39 
40 #include "G4HadronicException.hh"
41 #include "G4NeutronElasticXS.hh"
42 #include "G4Neutron.hh"
43 #include "G4DynamicParticle.hh"
44 #include "G4Element.hh"
45 #include "G4ElementTable.hh"
46 #include "G4PhysicsLogVector.hh"
47 #include "G4PhysicsVector.hh"
49 #include "G4HadronNucleonXsc.hh"
50 #include "G4NistManager.hh"
51 #include "G4Proton.hh"
52 
53 #include <iostream>
54 #include <fstream>
55 #include <sstream>
56 
57 using namespace std;
58 
60  : G4VCrossSectionDataSet("G4NeutronElasticXS"),
61  proton(G4Proton::Proton()), maxZ(92)
62 {
63  // verboseLevel = 0;
64  if(verboseLevel > 0){
65  G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
66  << maxZ + 1 << G4endl;
67  }
68  data.resize(maxZ+1, 0);
69  coeff.resize(maxZ+1, 1.0);
70  ggXsection = new G4GlauberGribovCrossSection();
71  fNucleon = new G4HadronNucleonXsc();
72  isInitialized = false;
73 }
74 
76 {
77  delete fNucleon;
78  for(G4int i=0; i<=maxZ; ++i) {
79  delete data[i];
80  }
81 }
82 
84 {
85  outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
86  << "cross section on nuclei using data from the high precision\n"
87  << "neutron database. These data are simplified and smoothed over\n"
88  << "the resonance region in order to reduce CPU time.\n"
89  << "G4NeutronElasticXS is valid for energies up to 20 MeV, for all\n"
90  << "targets through U.\n";
91 }
92 
93 G4bool
95  G4int, const G4Material*)
96 {
97  return true;
98 }
99 
100 G4double
102  G4int Z, const G4Material*)
103 {
104  G4double xs = 0.0;
105  G4double ekin = aParticle->GetKineticEnergy();
106 
107  if(Z < 1 || Z > maxZ) { return xs; }
108 
109  G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
110  G4PhysicsVector* pv = data[Z];
111  // G4cout << "G4NeutronElasticXS::GetCrossSection e= " << ekin << " Z= " << Z << G4endl;
112 
113  // element was not initialised
114  if(!pv) {
115  Initialise(Z);
116  pv = data[Z];
117  if(!pv) { return xs; }
118  }
119 
120  G4double e1 = pv->Energy(0);
121  if(ekin <= e1) { return (*pv)[0]; }
122 
123  G4int n = pv->GetVectorLength() - 1;
124  G4double e2 = pv->Energy(n);
125 
126  if(ekin <= e2) {
127  xs = pv->Value(ekin);
128  } else if(1 == Z) {
129  fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
130  xs = coeff[1]*fNucleon->GetElasticHadronNucleonXsc();
131  } else {
132  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
133  xs = coeff[Z]*ggXsection->GetElasticGlauberGribovXsc();
134  }
135 
136  if(verboseLevel > 0){
137  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
138  }
139  return xs;
140 }
141 
142 void
144 {
145  if(isInitialized) { return; }
146  if(verboseLevel > 0){
147  G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
148  << p.GetParticleName() << G4endl;
149  }
150  if(p.GetParticleName() != "neutron") {
151  throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
152  return;
153  }
154  isInitialized = true;
155 
156  // check environment variable
157  // Build the complete string identifying the file with the data set
158  char* path = getenv("G4NEUTRONXSDATA");
159  if (!path){
160  throw G4HadronicException(__FILE__, __LINE__,
161  "G4NEUTRONXSDATA environment variable not defined");
162  return;
163  }
164 
165  G4DynamicParticle* dynParticle =
167 
168  // Access to elements
169  const G4ElementTable* theElmTable = G4Element::GetElementTable();
170  size_t numOfElm = G4Element::GetNumberOfElements();
171  if(numOfElm > 0) {
172  for(size_t i=0; i<numOfElm; ++i) {
173  G4int Z = G4int(((*theElmTable)[i])->GetZ());
174  if(Z < 1) { Z = 1; }
175  else if(Z > maxZ) { Z = maxZ; }
176  //G4cout << "Z= " << Z << G4endl;
177  // Initialisation
178  if(!data[Z]) { Initialise(Z, dynParticle, path); }
179  }
180  }
181  delete dynParticle;
182 }
183 
184 void
185 G4NeutronElasticXS::Initialise(G4int Z, G4DynamicParticle* dp,
186  const char* p)
187 {
188  if(data[Z]) { return; }
189  const char* path = p;
190  if(!p) {
191  // check environment variable
192  // Build the complete string identifying the file with the data set
193  path = getenv("G4NEUTRONXSDATA");
194  if (!path) {
195  throw G4HadronicException(__FILE__, __LINE__,
196  "G4NEUTRONXSDATA environment variable not defined");
197  return;
198  }
199  }
200  G4DynamicParticle* dynParticle = dp;
201  if(!dp) {
202  dynParticle =
204  }
205 
206  G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
207 
208  // upload data from file
209  data[Z] = new G4PhysicsLogVector();
210 
211  std::ostringstream ost;
212  ost << path << "/elast" << Z ;
213  std::ifstream filein(ost.str().c_str());
214  if (!(filein)) {
215  G4cout << ost.str()
216  << " is not opened by G4NeutronElasticXS" << G4endl;
217  throw G4HadronicException(__FILE__, __LINE__,"NO data sets opened");
218  return;
219  }else{
220  if(verboseLevel > 1) {
221  G4cout << "file " << ost.str()
222  << " is opened by G4NeutronElasticXS" << G4endl;
223  }
224 
225  // retrieve data from DB
226  data[Z]->Retrieve(filein, true);
227 
228  // smooth transition
229  size_t n = data[Z]->GetVectorLength() - 1;
230  G4double emax = data[Z]->Energy(n);
231  G4double sig1 = (*data[Z])[n];
232  dynParticle->SetKineticEnergy(emax);
233  G4double sig2 = 0.0;
234  if(1 == Z) {
235  fNucleon->GetHadronNucleonXscPDG(dynParticle, proton);
236  sig2 = fNucleon->GetElasticHadronNucleonXsc();
237  } else {
238  ggXsection->GetIsoCrossSection(dynParticle, Z, Amean);
239  sig2 = ggXsection->GetElasticGlauberGribovXsc();
240  }
241  if(sig2 > 0.) { coeff[Z] = sig1/sig2; }
242  }
243  if(!dp) { delete dynParticle; }
244 }