Geant4  10.02.p03
G4NeutronInelasticXS Class Reference

#include <G4NeutronInelasticXS.hh>

Inheritance diagram for G4NeutronInelasticXS:
Collaboration diagram for G4NeutronInelasticXS:

Public Member Functions

 G4NeutronInelasticXS ()
 
virtual ~G4NeutronInelasticXS ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *)
 
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 G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Static Public Member Functions

static const char * Default_Name ()
 

Private Member Functions

void Initialise (G4int Z, G4DynamicParticle *dp=0, const char *=0)
 
G4PhysicsVectorRetrieveVector (std::ostringstream &in, G4bool warn)
 
G4double IsoCrossSection (G4double ekin, G4int Z, G4int A)
 
G4NeutronInelasticXSoperator= (const G4NeutronInelasticXS &right)
 
 G4NeutronInelasticXS (const G4NeutronInelasticXS &)
 

Private Attributes

G4ComponentGGHadronNucleusXscggXsection
 
G4HadronNucleonXscfNucleon
 
const G4ParticleDefinitionproton
 
G4bool isMaster
 
std::vector< G4doubletemp
 

Static Private Attributes

static G4ElementDatadata = 0
 
static G4double coeff [MAXZINEL] = {1.0}
 
static const G4int amin [MAXZINEL]
 
static const G4int amax [MAXZINEL]
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 62 of file G4NeutronInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronInelasticXS() [1/2]

G4NeutronInelasticXS::G4NeutronInelasticXS ( )

Definition at line 93 of file G4NeutronInelasticXS.cc.

96 {
97  // verboseLevel = 0;
98  if(verboseLevel > 0){
99  G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
100  << MAXZINEL << G4endl;
101  }
104  isMaster = false;
105 }
G4VCrossSectionDataSet(const G4String &nam="")
const G4int MAXZINEL
const G4ParticleDefinition * proton
static const char * Default_Name()
G4GLOB_DLL std::ostream G4cout
G4ComponentGGHadronNucleusXsc * ggXsection
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4HadronNucleonXsc * fNucleon
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ ~G4NeutronInelasticXS()

G4NeutronInelasticXS::~G4NeutronInelasticXS ( )
virtual

Definition at line 107 of file G4NeutronInelasticXS.cc.

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 }
G4ComponentGGHadronNucleusXsc * ggXsection
G4HadronNucleonXsc * fNucleon
static G4ElementData * data

◆ G4NeutronInelasticXS() [2/2]

G4NeutronInelasticXS::G4NeutronInelasticXS ( const G4NeutronInelasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 279 of file G4NeutronInelasticXS.cc.

280 {
281  if(verboseLevel > 0){
282  G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
283  << p.GetParticleName() << G4endl;
284  }
285  if(p.GetParticleName() != "neutron") {
287  ed << p.GetParticleName() << " is a wrong particle type -"
288  << " only neutron is allowed";
289  G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
290  FatalException, ed, "");
291  return;
292  }
293 
294  if(!data) {
295  isMaster = true;
296  data = new G4ElementData();
297  data->SetName("NeutronInelastic");
298  temp.resize(13,0.0);
299  }
300 
301  // it is possible re-initialisation for the new run
302  if(isMaster) {
303 
304  // check environment variable
305  // Build the complete string identifying the file with the data set
306  char* path = getenv("G4NEUTRONXSDATA");
307 
308  G4DynamicParticle* dynParticle =
310 
311  // Access to elements
312  const G4ElementTable* theElmTable = G4Element::GetElementTable();
313  size_t numOfElm = G4Element::GetNumberOfElements();
314  if(numOfElm > 0) {
315  for(size_t i=0; i<numOfElm; ++i) {
316  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
317  if(Z < 1) { Z = 1; }
318  else if(Z >= MAXZINEL) { Z = MAXZINEL-1; }
319  //G4cout << "Z= " << Z << G4endl;
320  // Initialisation
321  if(!(data->GetElementData(Z))) {
322  Initialise(Z, dynParticle, path);
323  }
324  }
325  }
326  delete dynParticle;
327  }
328 }
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
G4PhysicsVector * GetElementData(G4int Z)
const G4int MAXZINEL
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
Float_t Z
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
int G4lrint(double ad)
Definition: templates.hh:163
std::vector< G4double > temp
#define G4endl
Definition: G4ios.hh:61
void SetName(const G4String &nam)
std::vector< G4Element * > G4ElementTable
static G4ElementData * data
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionDescription()

void G4NeutronInelasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 116 of file G4NeutronInelasticXS.cc.

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 }
Here is the caller graph for this function:

◆ Default_Name()

static const char* G4NeutronInelasticXS::Default_Name ( )
inlinestatic

Definition at line 70 of file G4NeutronInelasticXS.hh.

70 {return "G4NeutronInelasticXS";}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetElementCrossSection()

G4double G4NeutronInelasticXS::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 141 of file G4NeutronInelasticXS.cc.

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 }
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
G4PhysicsVector * GetElementData(G4int Z)
static G4double coeff[MAXZINEL]
static const G4double e2
const G4int MAXZINEL
const G4ParticleDefinition * proton
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
Float_t Z
G4ComponentGGHadronNucleusXsc * ggXsection
G4double Value(G4double theEnergy, size_t &lastidx) const
static const G4double e1
G4double GetMaxEnergy() const
int G4lrint(double ad)
Definition: templates.hh:163
G4HadronNucleonXsc * fNucleon
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
static G4ElementData * data
G4double GetInelasticHadronNucleonXsc()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetIsoCrossSection()

G4double G4NeutronInelasticXS::GetIsoCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A,
const G4Isotope iso,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 185 of file G4NeutronInelasticXS.cc.

190 {
191  if(Z >= MAXZINEL) { Z = MAXZINEL - 1; }
192  else if(Z < 1) { Z = 1; }
193 
194  return IsoCrossSection(aParticle->GetKineticEnergy(), Z, A);
195 }
const G4int MAXZINEL
G4double GetKineticEnergy() const
double A(double temperature)
Float_t Z
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

void G4NeutronInelasticXS::Initialise ( G4int  Z,
G4DynamicParticle dp = 0,
const char *  p = 0 
)
private

Definition at line 331 of file G4NeutronInelasticXS.cc.

333 {
334  if(data->GetElementData(Z) || Z < 1 || Z >= MAXZINEL) { return; }
335  const char* path = p;
336  if(!p) {
337  // check environment variable
338  // Build the complete string identifying the file with the data set
339  path = getenv("G4NEUTRONXSDATA");
340  if (!path) {
341  G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
343  "Environment variable G4NEUTRONXSDATA is not defined");
344  return;
345  }
346  }
347  G4DynamicParticle* dynParticle = dp;
348  if(!dp) {
349  dynParticle =
351  }
352 
353  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
354 
355  // upload element data
356  std::ostringstream ost;
357  ost << path << "/inelast" << Z ;
358  G4PhysicsVector* v = RetrieveVector(ost, true);
359  data->InitialiseForElement(Z, v);
360  /*
361  G4cout << "G4NeutronInelasticXS::Initialise for Z= " << Z
362  << " A= " << Amean << " Amin= " << amin[Z]
363  << " Amax= " << amax[Z] << G4endl;
364  */
365  // upload isotope data
366  if(amin[Z] > 0) {
367  size_t nmax = (size_t)(amax[Z]-amin[Z]+1);
368  data->InitialiseForComponent(Z, nmax);
369 
370  for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
371  std::ostringstream ost1;
372  ost1 << path << "/inelast" << Z << "_" << A;
373  G4PhysicsVector* v1 = RetrieveVector(ost1, false);
374  data->AddComponent(Z, A, v1);
375  }
376  }
377 
378  // smooth transition
379  G4double emax = v->GetMaxEnergy();
380  G4double sig1 = (*v)[v->GetVectorLength() - 1];
381  dynParticle->SetKineticEnergy(emax);
382  G4double sig2 = 0.0;
383  if(1 == Z) {
384  fNucleon->GetHadronNucleonXscPDG(dynParticle, proton);
386  } else {
387  ggXsection->GetIsoCrossSection(dynParticle, Z, Amean);
389  }
390  if(sig2 > 0.) { coeff[Z] = sig1/sig2; }
391  if(!dp) { delete dynParticle; }
392 }
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
CLHEP::Hep3Vector G4ThreeVector
G4PhysicsVector * GetElementData(G4int Z)
static G4double coeff[MAXZINEL]
const G4int MAXZINEL
static const G4int amax[MAXZINEL]
const G4ParticleDefinition * proton
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
double A(double temperature)
Float_t Z
G4ComponentGGHadronNucleusXsc * ggXsection
static const G4int amin[MAXZINEL]
const G4int nmax
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
size_t GetVectorLength() const
void SetKineticEnergy(G4double aEnergy)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
G4double GetMaxEnergy() const
int G4lrint(double ad)
Definition: templates.hh:163
G4HadronNucleonXsc * fNucleon
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
double G4double
Definition: G4Types.hh:76
static G4ElementData * data
void InitialiseForComponent(G4int Z, G4int nComponents=0)
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool warn)
G4double GetInelasticHadronNucleonXsc()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

G4bool G4NeutronInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 127 of file G4NeutronInelasticXS.cc.

129 {
130  return true;
131 }
Here is the caller graph for this function:

◆ IsIsoApplicable()

G4bool G4NeutronInelasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element ,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 134 of file G4NeutronInelasticXS.cc.

137 {
138  return true;
139 }
Here is the caller graph for this function:

◆ IsoCrossSection()

G4double G4NeutronInelasticXS::IsoCrossSection ( G4double  ekin,
G4int  Z,
G4int  A 
)
private

Definition at line 198 of file G4NeutronInelasticXS.cc.

199 {
200  G4double xs = 0.0;
201  /*
202  G4cout << "IsoCrossSection Z= " << Z << " A= " << A
203  << " Amin= " << amin[Z] << " Amax= " << amax[Z]
204  << " E(MeV)= " << ekin << G4endl;
205  */
206  // element was not initialised
208  if(!pv) {
209  Initialise(Z);
210  pv = data->GetElementData(Z);
211  }
212 
213  // isotope cross section exist
214  if(pv && amin[Z] > 0 && A >= amin[Z] && A <= amax[Z]) {
215  pv = data->GetComponentDataByIndex(Z, A - amin[Z]);
216  if(pv && ekin > pv->Energy(0)) { xs = pv->Value(ekin); }
217  }
218  if(verboseLevel > 0){
219  G4cout << "ekin= " << ekin << ", xs= " << xs << G4endl;
220  }
221  return xs;
222 }
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
G4PhysicsVector * GetComponentDataByIndex(G4int Z, size_t idx)
G4PhysicsVector * GetElementData(G4int Z)
static const G4int amax[MAXZINEL]
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
static const G4int amin[MAXZINEL]
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
static G4ElementData * data
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator=()

G4NeutronInelasticXS& G4NeutronInelasticXS::operator= ( const G4NeutronInelasticXS right)
private
Here is the caller graph for this function:

◆ RetrieveVector()

G4PhysicsVector * G4NeutronInelasticXS::RetrieveVector ( std::ostringstream &  in,
G4bool  warn 
)
private

Definition at line 395 of file G4NeutronInelasticXS.cc.

396 {
397  G4PhysicsLogVector* v = 0;
398  std::ifstream filein(ost.str().c_str());
399  if (!(filein)) {
400  if(warn) {
402  ed << "Data file <" << ost.str().c_str()
403  << "> is not opened!";
404  G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
405  FatalException, ed, "Check G4NEUTRONXSDATA");
406  }
407  } else {
408  if(verboseLevel > 1) {
409  G4cout << "File " << ost.str()
410  << " is opened by G4NeutronInelasticXS" << G4endl;
411  }
412  // retrieve data from DB
413  v = new G4PhysicsLogVector();
414  if(!v->Retrieve(filein, true)) {
416  ed << "Data file <" << ost.str().c_str()
417  << "> is not retrieved!";
418  G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
419  FatalException, ed, "Check G4NEUTRONXSDATA");
420  }
421  }
422  return v;
423 }
virtual G4bool Retrieve(std::ifstream &fIn, G4bool ascii)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SelectIsotope()

G4Isotope * G4NeutronInelasticXS::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 224 of file G4NeutronInelasticXS.cc.

226 {
227  size_t nIso = anElement->GetNumberOfIsotopes();
228  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
229  G4Isotope* iso = (*isoVector)[0];
230 
231  //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
232 
233  // more than 1 isotope
234  if(1 < nIso) {
235  G4int Z = G4lrint(anElement->GetZ());
236  //G4cout << "SelectIsotope Z= " << Z << G4endl;
237 
238  G4double* abundVector = anElement->GetRelativeAbundanceVector();
239  G4double q = G4UniformRand();
240  G4double sum = 0.0;
241 
242  // is there isotope wise cross section?
243  size_t j;
244  if(0 == amin[Z] || Z >= MAXZINEL) {
245  for (j = 0; j<nIso; ++j) {
246  sum += abundVector[j];
247  if(q <= sum) {
248  iso = (*isoVector)[j];
249  break;
250  }
251  }
252  } else {
253 
254  // element may be not initialised in unit test
255  if(!data->GetElementData(Z)) { Initialise(Z); }
256  size_t nn = temp.size();
257  if(nn < nIso) { temp.resize(nIso, 0.); }
258 
259  for (j=0; j<nIso; ++j) {
260  //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
261  // << " abund= " << abundVector[j] << G4endl;
262  sum += abundVector[j]*IsoCrossSection(kinEnergy, Z,
263  (*isoVector)[j]->GetN());
264  temp[j] = sum;
265  }
266  sum *= q;
267  for (j = 0; j<nIso; ++j) {
268  if(temp[j] >= sum) {
269  iso = (*isoVector)[j];
270  break;
271  }
272  }
273  }
274  }
275  return iso;
276 }
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
std::vector< G4Isotope * > G4IsotopeVector
G4PhysicsVector * GetElementData(G4int Z)
const G4int MAXZINEL
int G4int
Definition: G4Types.hh:78
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
#define G4UniformRand()
Definition: Randomize.hh:97
Float_t Z
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
static const G4int amin[MAXZINEL]
G4double IsoCrossSection(G4double ekin, G4int Z, G4int A)
int G4lrint(double ad)
Definition: templates.hh:163
std::vector< G4double > temp
double G4double
Definition: G4Types.hh:76
static G4ElementData * data
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ amax

const G4int G4NeutronInelasticXS::amax
staticprivate
Initial value:
= {
0,
0, 0, 7, 0,11,13,15,18, 0, 0,
0, 0, 0,30, 0, 0, 0,40, 0,48,
0, 0, 0, 0, 0,58, 0,64,65,70,
0,76, 0, 0, 0, 0, 0, 0, 0,96,
0, 0, 0, 0, 0, 0,109,116, 0,124,
0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
0, 0, 0,186, 0, 0, 0, 0, 0, 0,
0,208, 0, 0, 0, 0, 0, 0, 0, 0,
0,238}

Definition at line 121 of file G4NeutronInelasticXS.hh.

◆ amin

const G4int G4NeutronInelasticXS::amin
staticprivate
Initial value:
= {
0,
0, 0, 6, 0,10,12,14,16, 0, 0,
0, 0, 0,28, 0, 0, 0,36, 0,40,
0, 0, 0, 0, 0,54, 0,58,63,64,
0,70, 0, 0, 0, 0, 0, 0, 0,90,
0, 0, 0, 0, 0, 0,107,106, 0,112,
0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
0, 0, 0,180, 0, 0, 0, 0, 0, 0,
0,204, 0, 0, 0, 0, 0, 0, 0, 0,
0,235}

Definition at line 120 of file G4NeutronInelasticXS.hh.

◆ coeff

G4double G4NeutronInelasticXS::coeff = {1.0}
staticprivate

Definition at line 118 of file G4NeutronInelasticXS.hh.

◆ data

G4ElementData * G4NeutronInelasticXS::data = 0
staticprivate

Definition at line 115 of file G4NeutronInelasticXS.hh.

◆ fNucleon

G4HadronNucleonXsc* G4NeutronInelasticXS::fNucleon
private

Definition at line 109 of file G4NeutronInelasticXS.hh.

◆ ggXsection

G4ComponentGGHadronNucleusXsc* G4NeutronInelasticXS::ggXsection
private

Definition at line 108 of file G4NeutronInelasticXS.hh.

◆ isMaster

G4bool G4NeutronInelasticXS::isMaster
private

Definition at line 113 of file G4NeutronInelasticXS.hh.

◆ proton

const G4ParticleDefinition* G4NeutronInelasticXS::proton
private

Definition at line 111 of file G4NeutronInelasticXS.hh.

◆ temp

std::vector<G4double> G4NeutronInelasticXS::temp
private

Definition at line 116 of file G4NeutronInelasticXS.hh.


The documentation for this class was generated from the following files: