Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronInelasticXS.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: G4NeutronInelasticXS
34//
35// Author Ivantchenko, Geant4, 3-Aug-09
36//
37// Modifications:
38//
39
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#include "Randomize.hh"
53
54#include <iostream>
55#include <fstream>
56#include <sstream>
57
58using namespace std;
59
61 : G4VCrossSectionDataSet("G4NeutronInelasticXS"),
62 proton(G4Proton::Proton()), maxZ(92)
63{
64 // verboseLevel = 0;
65 if(verboseLevel > 0){
66 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
67 << maxZ + 1 << G4endl;
68 }
69 data.resize(maxZ+1, 0);
70 coeff.resize(maxZ+1, 1.0);
71 ggXsection = new G4GlauberGribovCrossSection();
72 fNucleon = new G4HadronNucleonXsc();
73 isInitialized = false;
74}
75
77{
78 delete fNucleon;
79 for(G4int i=0; i<=maxZ; ++i) {
80 delete data[i];
81 }
82}
83
84void G4NeutronInelasticXS::CrossSectionDescription(std::ostream& outFile) const
85{
86 outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
87 << "cross section on nuclei using data from the high precision\n"
88 << "neutron database. These data are simplified and smoothed over\n"
89 << "the resonance region in order to reduce CPU time.\n"
90 << "G4NeutronInelasticXS is valid for energies up to 20 MeV, for\n"
91 << "nuclei through U.\n";
92}
93
94G4bool
96 G4int, const G4Material*)
97{
98 return true;
99}
100
101G4bool
103 G4int /*ZZ*/, G4int /*AA*/,
104 const G4Element*, const G4Material*)
105{
106 return false;
107}
108
111 G4int Z, const G4Material*)
112{
113 G4double xs = 0.0;
114 G4double ekin = aParticle->GetKineticEnergy();
115
116 if(Z < 1 || Z > maxZ) { return xs; }
117 G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
118 G4PhysicsVector* pv = data[Z];
119 // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin << " Z= " << Z << G4endl;
120
121 // element was not initialised
122 if(!pv) {
123 Initialise(Z);
124 pv = data[Z];
125 if(!pv) { return xs; }
126 }
127
128 G4double e1 = pv->Energy(0);
129 if(ekin <= e1) { return xs; }
130
131 G4int n = pv->GetVectorLength() - 1;
132 G4double e2 = pv->Energy(n);
133 if(ekin <= e2) {
134 xs = pv->Value(ekin);
135 } else if(1 == Z) {
136 fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
137 xs = coeff[1]*fNucleon->GetInelasticHadronNucleonXsc();
138 } else {
139 ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
140 xs = coeff[Z]*ggXsection->GetInelasticGlauberGribovXsc();
141 }
142
143 if(verboseLevel > 0) {
144 G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
145 }
146 return xs;
147}
148
149/*
150G4double
151G4NeutronInelasticXS::GetIsoCrossSection(const G4DynamicParticle* aParticle,
152 G4int Z, G4int A,
153 const G4Isotope*, const G4Element*,
154 const G4Material*)
155{
156 G4double xs = 0.0;
157 G4double ekin = aParticle->GetKineticEnergy();
158 if(ekin > emax || Z < 1 || Z > maxZ) { return xs; }
159 const G4double elimit = 1.0e-10*eV;
160 if(ekin < elimit) { ekin = elimit; }
161
162 // G4PhysicsVector* pv = data[Z];
163 G4PhysicsVector* pv = data.GetElementData(Z);
164
165 // element was not initialised
166 if(!pv) {
167 Initialise(Z);
168 // pv = data[Z];
169 pv = data.GetElementData(Z);
170 if(!pv) { return xs; }
171 }
172 pv = data.GetComponentDataByID(Z, A);
173 if(!pv) { return xs; }
174
175 xs = pv->Value(ekin);
176
177 if(verboseLevel > 0){
178 G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
179 }
180 return xs;
181}
182
183G4Isotope* G4NeutronInelasticXS::SelectIsotope(const G4Element* anElement,
184 G4double kinEnergy)
185{
186 G4int nIso = anElement->GetNumberOfIsotopes();
187 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
188 G4Isotope* iso = (*isoVector)[0];
189
190 // more than 1 isotope
191 if(1 < nIso) {
192 G4int Z = G4lrint(anElement->GetZ());
193 if(Z > maxZ) { Z = maxZ; }
194 G4double* abundVector = anElement->GetRelativeAbundanceVector();
195 G4double q = G4UniformRand();
196 G4double sum = 0.0;
197
198 // is there isotope wise cross section?
199 if(0 == amin[Z]) {
200 for (G4int j = 0; j<nIso; ++j) {
201 sum += abundVector[j];
202 if(q <= sum) {
203 iso = (*isoVector)[j];
204 break;
205 }
206 }
207 } else {
208 size_t nmax = data.GetNumberOfComponents(Z);
209 if(temp.size() < nmax) { temp.resize(nmax,0.0); }
210 for (size_t i=0; i<nmax; ++i) {
211 G4int A = (*isoVector)[i]->GetN();
212 G4PhysicsVector* v = data.GetComponentDataByID(Z, A);
213 if(v) { sum += abundVector[i]*v->Value(kinEnergy); }
214 temp[i] = sum;
215 }
216 sum *= q;
217 for (size_t j = 0; j<nmax; ++j) {
218 if(temp[j] >= sum) {
219 iso = (*isoVector)[j];
220 break;
221 }
222 }
223 }
224 }
225 return iso;
226}
227*/
228void
230{
231 if(isInitialized) { return; }
232 if(verboseLevel > 0){
233 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
234 << p.GetParticleName() << G4endl;
235 }
236 if(p.GetParticleName() != "neutron") {
237 throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
238 return;
239 }
240 isInitialized = true;
241
242 // check environment variable
243 // Build the complete string identifying the file with the data set
244 char* path = getenv("G4NEUTRONXSDATA");
245 if (!path){
246 throw G4HadronicException(__FILE__, __LINE__,
247 "G4NEUTRONXSDATA environment variable not defined");
248 return;
249 }
250
251 G4DynamicParticle* dynParticle =
253
254 // Access to elements
255 const G4ElementTable* theElmTable = G4Element::GetElementTable();
256 size_t numOfElm = G4Element::GetNumberOfElements();
257 if(numOfElm > 0) {
258 for(size_t i=0; i<numOfElm; ++i) {
259 G4int Z = G4int(((*theElmTable)[i])->GetZ());
260 if(Z < 1) { Z = 1; }
261 else if(Z > maxZ) { Z = maxZ; }
262 //G4cout << "Z= " << Z << G4endl;
263 // Initialisation
264 if(!data[Z]) { Initialise(Z, dynParticle, path); }
265 }
266 }
267 delete dynParticle;
268}
269
270void
271G4NeutronInelasticXS::Initialise(G4int Z, G4DynamicParticle* dp,
272 const char* p)
273{
274 if(data[Z]) { return; }
275 const char* path = p;
276 if(!p) {
277 // check environment variable
278 // Build the complete string identifying the file with the data set
279 path = getenv("G4NEUTRONXSDATA");
280 if (!path) {
281 throw G4HadronicException(__FILE__, __LINE__,
282 "G4NEUTRONXSDATA environment variable not defined");
283 return;
284 }
285 }
286 G4DynamicParticle* dynParticle = dp;
287 if(!dp) {
288 dynParticle =
290 }
291
292 G4int Amean = G4int(G4NistManager::Instance()->GetAtomicMassAmu(Z)+0.5);
293
294 // upload data from file
295 data[Z] = new G4PhysicsLogVector();
296
297 std::ostringstream ost;
298 ost << path << "/inelast" << Z ;
299 std::ifstream filein(ost.str().c_str());
300
301 if (!(filein)) {
302 throw G4HadronicException(__FILE__, __LINE__,"NO data sets opened");
303 return;
304 }else{
305 if(verboseLevel > 1) {
306 G4cout << "file " << ost.str()
307 << " is opened by G4NeutronInelasticXS" << G4endl;
308 }
309
310 // retrieve data from DB
311 data[Z]->Retrieve(filein, true);
312
313 // smooth transition
314 size_t n = data[Z]->GetVectorLength() - 1;
315 G4double emax = data[Z]->Energy(n);
316 G4double sig1 = (*data[Z])[n];
317 dynParticle->SetKineticEnergy(emax);
318 G4double sig2 = 0.0;
319 if(1 == Z) {
320 fNucleon->GetHadronNucleonXscPDG(dynParticle, proton);
321 sig2 = fNucleon->GetInelasticHadronNucleonXsc();
322 } else {
323 ggXsection->GetIsoCrossSection(dynParticle, Z, Amean);
324 sig2 = ggXsection->GetInelasticGlauberGribovXsc();
325 }
326 if(sig2 > 0.) { coeff[Z] = sig1/sig2; }
327 }
328 if(!dp) { delete dynParticle; }
329}
std::vector< G4Element * > G4ElementTable
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticHadronNucleonXsc()
virtual void CrossSectionDescription(std::ostream &) const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4double Value(G4double theEnergy)
size_t GetVectorLength() const
G4double Energy(size_t index) const