Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4NeutronInelasticXS Class Reference

#include <G4NeutronInelasticXS.hh>

Inheritance diagram for G4NeutronInelasticXS:
G4VCrossSectionDataSet

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
 

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::G4NeutronInelasticXS ( )

Definition at line 82 of file G4NeutronInelasticXS.cc.

References G4cout, G4endl, MAXZINEL, and G4VCrossSectionDataSet::verboseLevel.

83  : G4VCrossSectionDataSet("G4NeutronInelasticXS"),
84  proton(G4Proton::Proton())
85 {
86  // verboseLevel = 0;
87  if(verboseLevel > 0){
88  G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
89  << MAXZINEL << G4endl;
90  }
91  data.SetName("NeutronInelastic");
92  work.resize(13,0);
93  temp.resize(13,0.0);
94  coeff.resize(MAXZINEL, 1.0);
95  ggXsection = new G4GlauberGribovCrossSection();
96  fNucleon = new G4HadronNucleonXsc();
97  isInitialized = false;
98 }
G4VCrossSectionDataSet(const G4String &nam="")
const G4int MAXZINEL
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
#define G4endl
Definition: G4ios.hh:61
const XML_Char const XML_Char * data
G4NeutronInelasticXS::~G4NeutronInelasticXS ( )
virtual

Definition at line 100 of file G4NeutronInelasticXS.cc.

101 {
102  delete fNucleon;
103 }

Member Function Documentation

void G4NeutronInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 253 of file G4NeutronInelasticXS.cc.

References FatalException, G4cout, G4endl, G4Exception(), G4lrint(), G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4ParticleDefinition::GetParticleName(), MAXZINEL, G4Neutron::Neutron(), and G4VCrossSectionDataSet::verboseLevel.

254 {
255  if(isInitialized) { return; }
256  if(verboseLevel > 0){
257  G4cout << "G4NeutronCaptureXS::BuildPhysicsTable for "
258  << p.GetParticleName() << G4endl;
259  }
260  if(p.GetParticleName() != "neutron") {
262  ed << p.GetParticleName() << " is a wrong particle type -"
263  << " only neutron is allowed";
264  G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
265  FatalException, ed, "");
266  return;
267  }
268  isInitialized = true;
269 
270  // check environment variable
271  // Build the complete string identifying the file with the data set
272  char* path = getenv("G4NEUTRONXSDATA");
273 
274  G4DynamicParticle* dynParticle =
276 
277  // Access to elements
278  const G4ElementTable* theElmTable = G4Element::GetElementTable();
279  size_t numOfElm = G4Element::GetNumberOfElements();
280  if(numOfElm > 0) {
281  for(size_t i=0; i<numOfElm; ++i) {
282  G4int Z = G4lrint(((*theElmTable)[i])->GetZ());
283  if(Z < 1) { Z = 1; }
284  else if(Z >= MAXZINEL) { Z = MAXZINEL-1; }
285  //G4cout << "Z= " << Z << G4endl;
286  // Initialisation
287  if(!data.GetElementData(Z)) { Initialise(Z, dynParticle, path); }
288  }
289  }
290  delete dynParticle;
291 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
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
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
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
const XML_Char const XML_Char * data
void G4NeutronInelasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 105 of file G4NeutronInelasticXS.cc.

106 {
107  outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
108  << "cross section on nuclei using data from the high precision\n"
109  << "neutron database. These data are simplified and smoothed over\n"
110  << "the resonance region in order to reduce CPU time.\n"
111  << "G4NeutronInelasticXS is valid for energies up to 20 MeV, for\n"
112  << "nuclei through U.\n";
113 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4NeutronInelasticXS::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 131 of file G4NeutronInelasticXS.cc.

References G4PhysicsVector::Energy(), G4cout, G4endl, G4lrint(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4GlauberGribovCrossSection::GetInelasticGlauberGribovXsc(), G4HadronNucleonXsc::GetInelasticHadronNucleonXsc(), G4GlauberGribovCrossSection::GetIsoCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4PhysicsVector::GetMaxEnergy(), G4NistManager::Instance(), MAXZINEL, G4PhysicsVector::Value(), and G4VCrossSectionDataSet::verboseLevel.

133 {
134  G4double xs = 0.0;
135  G4double ekin = aParticle->GetKineticEnergy();
136 
137  if(Z < 1 || Z >= MAXZINEL) { return xs; }
138  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
139 
140  G4PhysicsVector* pv = data.GetElementData(Z);
141  // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin << " Z= "
142  // << Z << G4endl;
143 
144  // element was not initialised
145  if(!pv) {
146  Initialise(Z);
147  pv = data.GetElementData(Z);
148  if(!pv) { return xs; }
149  }
150 
151  G4double e1 = pv->Energy(0);
152  if(ekin <= e1) { return xs; }
153 
154  G4double e2 = pv->GetMaxEnergy();
155 
156  if(ekin <= e2) {
157  xs = pv->Value(ekin);
158  } else if(1 == Z) {
159  fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
160  xs = coeff[1]*fNucleon->GetInelasticHadronNucleonXsc();
161  } else {
162  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
163  xs = coeff[Z]*ggXsection->GetInelasticGlauberGribovXsc();
164  }
165 
166  if(verboseLevel > 0) {
167  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
168  }
169  return xs;
170 }
G4double GetMaxEnergy() const
G4double GetKineticEnergy() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
const G4int MAXZINEL
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4GLOB_DLL std::ostream G4cout
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data
G4double GetInelasticHadronNucleonXsc()
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 173 of file G4NeutronInelasticXS.cc.

References G4DynamicParticle::GetKineticEnergy(), and MAXZINEL.

177 {
178  G4double xs = 0.0;
179  G4double ekin = aParticle->GetKineticEnergy();
180  if(Z > 0 && Z < MAXZINEL) { xs = IsoCrossSection(ekin, Z, A); }
181  return xs;
182 }
G4double GetKineticEnergy() const
const G4int MAXZINEL
double G4double
Definition: G4Types.hh:76
G4bool G4NeutronInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 116 of file G4NeutronInelasticXS.cc.

118 {
119  return true;
120 }
G4bool G4NeutronInelasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element ,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 123 of file G4NeutronInelasticXS.cc.

126 {
127  return true;
128 }
G4Isotope * G4NeutronInelasticXS::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 208 of file G4NeutronInelasticXS.cc.

References G4lrint(), G4UniformRand, G4Element::GetIsotopeVector(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZ(), MAXZINEL, and nmax.

210 {
211  G4int nIso = anElement->GetNumberOfIsotopes();
212  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
213  G4Isotope* iso = (*isoVector)[0];
214 
215  // more than 1 isotope
216  if(1 < nIso) {
217  G4int Z = G4lrint(anElement->GetZ());
218  if(Z >= MAXZINEL) { Z = MAXZINEL - 1; }
219  G4double* abundVector = anElement->GetRelativeAbundanceVector();
220  G4double q = G4UniformRand();
221  G4double sum = 0.0;
222 
223  // is there isotope wise cross section?
224  if(0 == amin[Z]) {
225  for (G4int j = 0; j<nIso; ++j) {
226  sum += abundVector[j];
227  if(q <= sum) {
228  iso = (*isoVector)[j];
229  break;
230  }
231  }
232  } else {
233  size_t nmax = data.GetNumberOfComponents(Z);
234  if(temp.size() < nmax) { temp.resize(nmax,0.0); }
235  for (size_t i=0; i<nmax; ++i) {
236  G4int A = (*isoVector)[i]->GetN();
237  sum += abundVector[i]*IsoCrossSection(kinEnergy, Z, A);
238  temp[i] = sum;
239  }
240  sum *= q;
241  for (size_t j = 0; j<nmax; ++j) {
242  if(temp[j] >= sum) {
243  iso = (*isoVector)[j];
244  break;
245  }
246  }
247  }
248  }
249  return iso;
250 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
std::vector< G4Isotope * > G4IsotopeVector
G4double GetZ() const
Definition: G4Element.hh:131
const G4int MAXZINEL
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4int nmax
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data

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