Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4NeutronCaptureXS.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: G4NeutronCaptureXS
34//
35// Author Ivantchenko, Geant4, 3-Aug-09
36//
37// Modifications:
38//
39
40#include <fstream>
41#include <sstream>
42
44#include "G4SystemOfUnits.hh"
45#include "G4NeutronCaptureXS.hh"
46#include "G4Element.hh"
47#include "G4ElementTable.hh"
48#include "G4PhysicsLogVector.hh"
49#include "G4PhysicsVector.hh"
50#include "G4DynamicParticle.hh"
51#include "Randomize.hh"
52
53using namespace std;
54
55const G4int G4NeutronCaptureXS::amin[] = {0,
56 0, 0, 6, 0,10,12,14,16, 0, 0, //1-10
57 0, 0, 0,28, 0, 0, 0,36, 0,40, //11-20
58 0, 0, 0, 0, 0,54, 0,58,63,64, //21-30
59 0,70, 0, 0, 0, 0, 0, 0, 0,90, //31-40
60 0, 0, 0, 0, 0, 0,107,106, 0,112, //41-50
61 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //51-60
62 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //61-70
63 0, 0, 0,180, 0, 0, 0, 0, 0, 0, //71-80
64 0,204, 0, 0, 0, 0, 0, 0, 0, 0, //81-90
65 0,235};
66const G4int G4NeutronCaptureXS::amax[] = {0,
67 0, 0, 7, 0,11,13,15,18, 0, 0, //1-10
68 0, 0, 0,30, 0, 0, 0,40, 0,48, //11-20
69 0, 0, 0, 0, 0,58, 0,64,65,70, //21-30
70 0,76, 0, 0, 0, 0, 0, 0, 0,96, //31-40
71 0, 0, 0, 0, 0, 0,109,116, 0,124, //41-50
72 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //51-60
73 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //61-70
74 0, 0, 0,186, 0, 0, 0, 0, 0, 0, //71-80
75 0,208, 0, 0, 0, 0, 0, 0, 0, 0, //81-90
76 0,238};
77
79 : G4VCrossSectionDataSet("G4NeutronCaptureXS"),
80 emax(20*MeV),maxZ(92)
81{
82 // verboseLevel = 0;
83 if(verboseLevel > 0){
84 G4cout << "G4NeutronCaptureXS::G4NeutronCaptureXS: Initialise for Z < "
85 << maxZ + 1 << G4endl;
86 }
87 //data.resize(maxZ+1, 0);
88 data.SetName("NeutronCapture");
89 work.resize(13,0);
90 temp.resize(13,0.0);
91 isInitialized = false;
92}
93
95{
96 /*
97 for(G4int i=0; i<=maxZ; ++i) {
98 delete data[i];
99 }
100 */
101}
102
103void G4NeutronCaptureXS::CrossSectionDescription(std::ostream& outFile) const
104{
105 outFile << "G4NeutronCaptureXS calculates the neutron capture cross sections\n"
106 << "on nuclei using data from the high precision neutron database.\n"
107 << "These data are simplified and smoothed over the resonance region\n"
108 << "in order to reduce CPU time. G4NeutronCaptureXS is valid up to\n"
109 << "20 MeV for all targets through U.\n";
110}
111
112G4bool
114 G4int, const G4Material*)
115{
116 return true;
117}
118
119G4bool
121 G4int /*ZZ*/, G4int /*AA*/,
122 const G4Element*, const G4Material*)
123{
124 return false;
125}
126
129 G4int Z, const G4Material*)
130{
131 G4double xs = 0.0;
132 G4double ekin = aParticle->GetKineticEnergy();
133 if(ekin > emax || Z < 1 || Z > maxZ) { return xs; }
134 const G4double elimit = 1.0e-10*eV;
135 if(ekin < elimit) { ekin = elimit; }
136
137 // G4PhysicsVector* pv = data[Z];
138 G4PhysicsVector* pv = data.GetElementData(Z);
139
140 // element was not initialised
141 if(!pv) {
142 Initialise(Z);
143 // pv = data[Z];
144 pv = data.GetElementData(Z);
145 if(!pv) { return xs; }
146 }
147
148 G4double e1 = pv->Energy(0);
149 if(ekin < e1) { xs = (*pv)[0]*std::sqrt(e1/ekin); }
150 else { xs = pv->Value(ekin); }
151
152 if(verboseLevel > 0){
153 G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
154 }
155 return xs;
156}
157
160 G4int Z, G4int A,
161 const G4Isotope*, const G4Element*,
162 const G4Material*)
163{
164 G4double xs = 0.0;
165 G4double ekin = aParticle->GetKineticEnergy();
166 if(ekin > emax || Z < 1 || Z > maxZ) { return xs; }
167 const G4double elimit = 1.0e-10*eV;
168 if(ekin < elimit) { ekin = elimit; }
169
170 // G4PhysicsVector* pv = data[Z];
171 G4PhysicsVector* pv = data.GetElementData(Z);
172
173 // element was not initialised
174 if(!pv) {
175 Initialise(Z);
176 // pv = data[Z];
177 pv = data.GetElementData(Z);
178 if(!pv) { return xs; }
179 }
180 pv = data.GetComponentDataByID(Z, A);
181 if(!pv) { return xs; }
182
183 G4double e1 = pv->Energy(0);
184 if(ekin < e1) { xs = (*pv)[0]*std::sqrt(e1/ekin); }
185 else { xs = pv->Value(ekin); }
186
187 if(verboseLevel > 0){
188 G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
189 }
190 return xs;
191}
192
194 G4double kinEnergy)
195{
196 G4int nIso = anElement->GetNumberOfIsotopes();
197 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
198 G4Isotope* iso = (*isoVector)[0];
199
200 // more than 1 isotope
201 if(1 < nIso) {
202 G4int Z = G4lrint(anElement->GetZ());
203 if(Z > maxZ) { Z = maxZ; }
204 G4double* abundVector = anElement->GetRelativeAbundanceVector();
206 G4double sum = 0.0;
207
208 // is there isotope wise cross section?
209 if(0 == amin[Z]) {
210 for (G4int j = 0; j<nIso; ++j) {
211 sum += abundVector[j];
212 if(q <= sum) {
213 iso = (*isoVector)[j];
214 break;
215 }
216 }
217 } else {
218 size_t nmax = data.GetNumberOfComponents(Z);
219 if(temp.size() < nmax) { temp.resize(nmax,0.0); }
220 for (size_t i=0; i<nmax; ++i) {
221 G4int A = (*isoVector)[i]->GetN();
223 if(v) { sum += abundVector[i]*v->Value(kinEnergy); }
224 temp[i] = sum;
225 }
226 sum *= q;
227 for (size_t j = 0; j<nmax; ++j) {
228 if(temp[j] >= sum) {
229 iso = (*isoVector)[j];
230 break;
231 }
232 }
233 }
234 }
235 return iso;
236}
237
238void
240{
241 if(isInitialized) { return; }
242 if(verboseLevel > 0){
243 G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
244 << p.GetParticleName() << G4endl;
245 }
246 if(p.GetParticleName() != "neutron") {
247 throw G4HadronicException(__FILE__, __LINE__,"Wrong particle type");
248 return;
249 }
250 isInitialized = true;
251
252 // check environment variable
253 // Build the complete string identifying the file with the data set
254 char* path = getenv("G4NEUTRONXSDATA");
255 if (!path){
256 throw G4HadronicException(__FILE__, __LINE__,
257 "G4NEUTRONXSDATA environment variable not defined");
258 return;
259 }
260
261 // Access to elements
262 const G4ElementTable* theElmTable = G4Element::GetElementTable();
263 size_t numOfElm = G4Element::GetNumberOfElements();
264 if(numOfElm > 0) {
265 for(size_t i=0; i<numOfElm; ++i) {
266 G4int Z = G4int(((*theElmTable)[i])->GetZ());
267 if(Z < 1) { Z = 1; }
268 else if(Z > maxZ) { Z = maxZ; }
269 //G4cout << "Z= " << Z << G4endl;
270 // Initialisation
271 // if(!data[Z]) { Initialise(Z, path); }
272 if(!data.GetElementData(Z)) { Initialise(Z, path); }
273 }
274 }
275}
276
277void
278G4NeutronCaptureXS::Initialise(G4int Z, const char* p)
279{
280 if(data.GetElementData(Z)) { return; }
281 const char* path = p;
282
283 // check environment variable
284 if(!p) {
285 path = getenv("G4NEUTRONXSDATA");
286 if (!path) {
287 throw G4HadronicException(__FILE__, __LINE__,
288 "G4NEUTRONXSDATA environment variable not defined");
289 return;
290 }
291 }
292
293 // upload element data
294 std::ostringstream ost;
295 ost << path << "/cap" << Z ;
296 G4PhysicsVector* v = RetrieveVector(ost, true);
297 data.InitialiseForElement(Z, v);
298
299 // upload isotope data
300 if(amin[Z] > 0) {
301 size_t n = 0;
302 size_t i = 0;
303 size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
304 if(work.size() < nmax) { work.resize(nmax,0); }
305 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
306 std::ostringstream ost1;
307 ost1 << path << "/cap" << Z << "_" << A;
308 v = RetrieveVector(ost1, false);
309 if(v) { ++n; }
310 work[i] = v;
311 ++i;
312 }
313 data.InitialiseForComponent(Z, n);
314 for(size_t j=0; j<i; ++j) {
315 if(work[j]) { data.AddComponent(Z, amin[Z]+j, work[j]); }
316 }
317 }
318}
319
321G4NeutronCaptureXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
322{
323 G4PhysicsLogVector* v = 0;
324 std::ifstream filein(ost.str().c_str());
325 if (!(filein)) {
326 if(!warn) { return v; }
327 G4cout << ost.str() << " is not opened by G4NeutronCaptureXS" << G4endl;
328 throw G4HadronicException(__FILE__, __LINE__,"NO data sets opened");
329 }else{
330 if(verboseLevel > 1) {
331 G4cout << "File " << ost.str()
332 << " is opened by G4NeutronCaptureXS" << G4endl;
333 }
334 // retrieve data from DB
335 v = new G4PhysicsLogVector();
336 if(!v->Retrieve(filein, true)) {
337 G4cout << ost.str() << " is not retrieved in G4NeutronCaptureXS" << G4endl;
338 throw G4HadronicException(__FILE__, __LINE__,"ERROR: retrieve data fail");
339 }
340 }
341 return v;
342}
std::vector< G4Element * > G4ElementTable
std::vector< G4Isotope * > G4IsotopeVector
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
#define G4UniformRand()
Definition: Randomize.hh:53
G4double GetKineticEnergy() const
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
G4PhysicsVector * GetElementData(G4int Z)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
void SetName(const G4String &nam)
G4PhysicsVector * GetComponentDataByID(G4int Z, G4int id)
size_t GetNumberOfComponents(G4int Z)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:131
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
virtual void CrossSectionDescription(std::ostream &) const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
virtual G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *)
const G4String & GetParticleName() const
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii)
G4double Value(G4double theEnergy)
G4double Energy(size_t index) const
int G4lrint(double ad)
Definition: templates.hh:163