Geant4  10.01.p03
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: G4NeutronInelasticXS.cc 95293 2016-02-04 08:27:12Z gcosmo $
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 
40 #include "G4NeutronInelasticXS.hh"
41 #include "G4Neutron.hh"
42 #include "G4DynamicParticle.hh"
43 #include "G4Element.hh"
44 #include "G4ElementTable.hh"
45 #include "G4PhysicsLogVector.hh"
46 #include "G4PhysicsVector.hh"
48 #include "G4HadronNucleonXsc.hh"
49 #include "G4NistManager.hh"
50 #include "G4Proton.hh"
51 #include "Randomize.hh"
52 
53 #include <iostream>
54 #include <fstream>
55 #include <sstream>
56 
57 // factory
58 #include "G4CrossSectionFactory.hh"
59 //
61 
62 using namespace std;
63 
65  0,
66  0, 0, 6, 0,10,12,14,16, 0, 0, //1-10
67  0, 0, 0,28, 0, 0, 0,36, 0,40, //11-20
68  0, 0, 0, 0, 0,54, 0,58,63,64, //21-30
69  0,70, 0, 0, 0, 0, 0, 0, 0,90, //31-40
70  0, 0, 0, 0, 0, 0,107,106, 0,112, //41-50
71  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //51-60
72  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //61-70
73  0, 0, 0,180, 0, 0, 0, 0, 0, 0, //71-80
74  0,204, 0, 0, 0, 0, 0, 0, 0, 0, //81-90
75  0,235};
77  0,
78  0, 0, 7, 0,11,13,15,18, 0, 0, //1-10
79  0, 0, 0,30, 0, 0, 0,40, 0,48, //11-20
80  0, 0, 0, 0, 0,58, 0,64,65,70, //21-30
81  0,76, 0, 0, 0, 0, 0, 0, 0,96, //31-40
82  0, 0, 0, 0, 0, 0,109,116, 0,124, //41-50
83  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //51-60
84  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, //61-70
85  0, 0, 0,186, 0, 0, 0, 0, 0, 0, //71-80
86  0,208, 0, 0, 0, 0, 0, 0, 0, 0, //81-90
87  0,238};
88 
90 
92 
94  : G4VCrossSectionDataSet(Default_Name()),
96 {
97  // verboseLevel = 0;
98  if(verboseLevel > 0){
99  G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
100  << MAXZINEL << G4endl;
101  }
104  isMaster = false;
105 }
106 
108 {
109  //G4cout << "G4NeutronInelasticXS::~G4NeutronInelasticXS() "
110  // << " isMaster= " << isMaster << " data: " << data << G4endl;
111  delete fNucleon;
112  delete ggXsection;
113  if(isMaster) { delete data; data = 0; }
114 }
115 
116 void G4NeutronInelasticXS::CrossSectionDescription(std::ostream& outFile) const
117 {
118  outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
119  << "cross section on nuclei using data from the high precision\n"
120  << "neutron database. These data are simplified and smoothed over\n"
121  << "the resonance region in order to reduce CPU time.\n"
122  << "G4NeutronInelasticXS is valid for energies up to 20 MeV, for\n"
123  << "nuclei through U.\n";
124 }
125 
126 G4bool
128  G4int, const G4Material*)
129 {
130  return true;
131 }
132 
133 G4bool
135  G4int /*ZZ*/, G4int /*AA*/,
136  const G4Element*, const G4Material*)
137 {
138  return true;
139 }
140 
142  const G4DynamicParticle* aParticle,
143  G4int Z, const G4Material*)
144 {
145  G4double xs = 0.0;
146  G4double ekin = aParticle->GetKineticEnergy();
147 
148  if(Z >= MAXZINEL) { Z = MAXZINEL - 1; }
149  else if(Z < 1) { Z = 1; }
150 
151  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
152 
154  // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin
155  // << " Z= " << Z << G4endl;
156 
157  // element was not initialised
158  if(!pv) {
159  Initialise(Z);
160  pv = data->GetElementData(Z);
161  if(!pv) { return xs; }
162  }
163 
164  G4double e1 = pv->Energy(0);
165  if(ekin <= e1) { return xs; }
166 
167  G4double e2 = pv->GetMaxEnergy();
168 
169  if(ekin <= e2) {
170  xs = pv->Value(ekin);
171  } else if(1 == Z) {
174  } else {
175  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
177  }
178 
179  if(verboseLevel > 0) {
180  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
181  }
182  return xs;
183 }
184 
186  const G4DynamicParticle* aParticle,
187  G4int Z, G4int A,
188  const G4Isotope*, const G4Element*,
189  const G4Material*)
190 {
191  if(Z >= MAXZINEL) { Z = MAXZINEL - 1; }
192  else if(Z < 1) { Z = 1; }
193 
194  return IsoCrossSection(aParticle->GetKineticEnergy(), Z, A);
195 }
196 
197 G4double
199 {
200  G4double xs = 0.0;
201 
202  // element was not initialised
204  if(!pv) {
205  Initialise(Z);
206  pv = data->GetElementData(Z);
207  }
208 
209  // isotope cross section exist
210  if(pv && amin[Z] > 0 && A >= amin[Z] && A <= amax[Z]) {
211  pv = data->GetComponentDataByIndex(Z, A - amin[Z]);
212  if(pv && ekin > pv->Energy(0)) { xs = pv->Value(ekin); }
213  }
214  if(verboseLevel > 0){
215  G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
216  }
217  return xs;
218 }
219 
221  const G4Element* anElement, G4double kinEnergy)
222 {
223  size_t nIso = anElement->GetNumberOfIsotopes();
224  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
225  G4Isotope* iso = (*isoVector)[0];
226 
227  // more than 1 isotope
228  if(1 < nIso) {
229  G4int Z = G4lrint(anElement->GetZ());
230 
231  G4double* abundVector = anElement->GetRelativeAbundanceVector();
232  G4double q = G4UniformRand();
233  G4double sum = 0.0;
234 
235  // is there isotope wise cross section?
236  size_t j;
237  if(0 == amin[Z] || Z >= MAXZINEL) {
238  for (j = 0; j<nIso; ++j) {
239  sum += abundVector[j];
240  if(q <= sum) {
241  iso = (*isoVector)[j];
242  break;
243  }
244  }
245  } else {
246 
247  // element may be not initialised in unit test
248  if(!data->GetElementData(Z)) { Initialise(Z); }
249  size_t nn = temp.size();
250  if(nn < nIso) { temp.resize(nIso, 0.); }
251 
252  for (j=0; j<nIso; ++j) {
253  sum += abundVector[j]*IsoCrossSection(kinEnergy, Z,
254  (*isoVector)[j]->GetN());
255  temp[j] = sum;
256  }
257  sum *= q;
258  for (j = 0; j<nIso; ++j) {
259  if(temp[j] >= sum) {
260  iso = (*isoVector)[j];
261  break;
262  }
263  }
264  }
265  }
266  return iso;
267 }
268 
269 void
271 {
272  if(verboseLevel > 0){
273  G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
274  << p.GetParticleName() << G4endl;
275  }
276  if(p.GetParticleName() != "neutron") {
278  ed << p.GetParticleName() << " is a wrong particle type -"
279  << " only neutron is allowed";
280  G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
281  FatalException, ed, "");
282  return;
283  }
284 
285  if(!data) {
286  isMaster = true;
287  data = new G4ElementData();
288  data->SetName("NeutronInelastic");
289  temp.resize(13,0.0);
290  }
291 
292  // it is possible re-initialisation for the new run
293  if(isMaster) {
294 
295  // check environment variable
296  // Build the complete string identifying the file with the data set
297  char* path = getenv("G4NEUTRONXSDATA");
298 
299  G4DynamicParticle* dynParticle =
301 
302  // Access to elements
303  const G4ElementTable* theElmTable = G4Element::GetElementTable();
304  size_t numOfElm = G4Element::GetNumberOfElements();
305  if(numOfElm > 0) {
306  for(size_t i=0; i<numOfElm; ++i) {
307  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
308  if(Z < 1) { Z = 1; }
309  else if(Z >= MAXZINEL) { Z = MAXZINEL-1; }
310  //G4cout << "Z= " << Z << G4endl;
311  // Initialisation
312  if(!(data->GetElementData(Z))) {
313  Initialise(Z, dynParticle, path);
314  }
315  }
316  }
317  delete dynParticle;
318  }
319 }
320 
321 void
323  const char* p)
324 {
325  if(data->GetElementData(Z) || Z < 1 || Z >= MAXZINEL) { return; }
326  const char* path = p;
327  if(!p) {
328  // check environment variable
329  // Build the complete string identifying the file with the data set
330  path = getenv("G4NEUTRONXSDATA");
331  if (!path) {
332  G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
334  "Environment variable G4NEUTRONXSDATA is not defined");
335  return;
336  }
337  }
338  G4DynamicParticle* dynParticle = dp;
339  if(!dp) {
340  dynParticle =
342  }
343 
344  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
345 
346  // upload element data
347  std::ostringstream ost;
348  ost << path << "/inelast" << Z ;
349  G4PhysicsVector* v = RetrieveVector(ost, true);
350  data->InitialiseForElement(Z, v);
351  /*
352  G4cout << "G4NeutronInelasticXS::Initialise for Z= " << Z
353  << " A= " << Amean << " Amin= " << amin[Z]
354  << " Amax= " << amax[Z] << G4endl;
355  */
356  // upload isotope data
357  if(amin[Z] > 0) {
358  size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
359  data->InitialiseForComponent(Z, nmax);
360 
361  for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
362  std::ostringstream ost1;
363  ost1 << path << "/inelast" << Z << "_" << A;
364  G4PhysicsVector* v1 = RetrieveVector(ost1, false);
365  data->AddComponent(Z, A, v1);
366  }
367  }
368 
369  // smooth transition
370  G4double emax = v->GetMaxEnergy();
371  G4double sig1 = (*v)[v->GetVectorLength() - 1];
372  dynParticle->SetKineticEnergy(emax);
373  G4double sig2 = 0.0;
374  if(1 == Z) {
375  fNucleon->GetHadronNucleonXscPDG(dynParticle, proton);
377  } else {
378  ggXsection->GetIsoCrossSection(dynParticle, Z, Amean);
380  }
381  if(sig2 > 0.) { coeff[Z] = sig1/sig2; }
382  if(!dp) { delete dynParticle; }
383 }
384 
386 G4NeutronInelasticXS::RetrieveVector(std::ostringstream& ost, G4bool warn)
387 {
388  G4PhysicsLogVector* v = 0;
389  std::ifstream filein(ost.str().c_str());
390  if (!(filein)) {
391  if(warn) {
393  ed << "Data file <" << ost.str().c_str()
394  << "> is not opened!";
395  G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
396  FatalException, ed, "Check G4NEUTRONXSDATA");
397  }
398  } else {
399  if(verboseLevel > 1) {
400  G4cout << "File " << ost.str()
401  << " is opened by G4NeutronInelasticXS" << G4endl;
402  }
403  // retrieve data from DB
404  v = new G4PhysicsLogVector();
405  if(!v->Retrieve(filein, true)) {
407  ed << "Data file <" << ost.str().c_str()
408  << "> is not retrieved!";
409  G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
410  FatalException, ed, "Check G4NEUTRONXSDATA");
411  }
412  }
413  return v;
414 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
G4PhysicsVector * GetComponentDataByIndex(G4int Z, size_t idx)
G4_DECLARE_XS_FACTORY(G4NeutronInelasticXS)
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii)
G4double GetMaxEnergy() const
std::vector< G4Isotope * > G4IsotopeVector
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4PhysicsVector * GetElementData(G4int Z)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetZ() const
Definition: G4Element.hh:131
static G4double coeff[MAXZINEL]
static const G4double e2
const G4int MAXZINEL
static const G4int amax[MAXZINEL]
const G4ParticleDefinition * proton
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
#define G4UniformRand()
Definition: Randomize.hh:93
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
bool G4bool
Definition: G4Types.hh:79
static const G4int amin[MAXZINEL]
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4int nmax
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
static const G4double A[nN]
void SetKineticEnergy(G4double aEnergy)
static const G4double e1
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
int G4lrint(double ad)
Definition: templates.hh:163
std::vector< G4double > temp
virtual void CrossSectionDescription(std::ostream &) const
G4HadronNucleonXsc * fNucleon
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)
virtual G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void SetName(const G4String &nam)
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementData * data
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
void InitialiseForComponent(G4int Z, G4int nComponents=0)
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool warn)
G4double GetInelasticHadronNucleonXsc()
G4GlauberGribovCrossSection * ggXsection