Geant4  10.02.p03
G4NeutronElasticXS Class Reference

#include <G4NeutronElasticXS.hh>

Inheritance diagram for G4NeutronElasticXS:
Collaboration diagram for G4NeutronElasticXS:

Public Member Functions

 G4NeutronElasticXS ()
 
virtual ~G4NeutronElasticXS ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
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)
 
G4NeutronElasticXSoperator= (const G4NeutronElasticXS &right)
 
 G4NeutronElasticXS (const G4NeutronElasticXS &)
 

Private Attributes

G4ComponentGGHadronNucleusXscggXsection
 
G4HadronNucleonXscfNucleon
 
const G4ParticleDefinitionproton
 
G4bool isMaster
 

Static Private Attributes

static std::vector< G4PhysicsVector * > * data = 0
 
static G4double coeff [MAXZEL] = {0.0}
 

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 61 of file G4NeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronElasticXS() [1/2]

G4NeutronElasticXS::G4NeutronElasticXS ( )

Definition at line 66 of file G4NeutronElasticXS.cc.

69 {
70  // verboseLevel = 0;
71  if(verboseLevel > 0){
72  G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
73  << MAXZEL << G4endl;
74  }
77  isMaster = false;
78 }
const G4ParticleDefinition * proton
G4HadronNucleonXsc * fNucleon
G4VCrossSectionDataSet(const G4String &nam="")
G4GLOB_DLL std::ostream G4cout
G4ComponentGGHadronNucleusXsc * ggXsection
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const char * Default_Name()
#define G4endl
Definition: G4ios.hh:61
const G4int MAXZEL
Here is the caller graph for this function:

◆ ~G4NeutronElasticXS()

G4NeutronElasticXS::~G4NeutronElasticXS ( )
virtual

Definition at line 80 of file G4NeutronElasticXS.cc.

81 {
82  delete fNucleon;
83  delete ggXsection;
84  if(isMaster) {
85  for(G4int i=0; i<MAXZEL; ++i) {
86  delete (*data)[i];
87  (*data)[i] = 0;
88  }
89  delete data;
90  data = 0;
91  }
92 }
G4HadronNucleonXsc * fNucleon
int G4int
Definition: G4Types.hh:78
G4ComponentGGHadronNucleusXsc * ggXsection
static std::vector< G4PhysicsVector * > * data
const G4int MAXZEL

◆ G4NeutronElasticXS() [2/2]

G4NeutronElasticXS::G4NeutronElasticXS ( const G4NeutronElasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 155 of file G4NeutronElasticXS.cc.

156 {
157  if(verboseLevel > 0){
158  G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
159  << p.GetParticleName() << G4endl;
160  }
161  if(p.GetParticleName() != "neutron") {
163  ed << p.GetParticleName() << " is a wrong particle type -"
164  << " only neutron is allowed";
165  G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
166  FatalException, ed, "");
167  return;
168  }
169 
170  if(0 == coeff[0]) {
171  isMaster = true;
172  for(G4int i=0; i<MAXZEL; ++i) { coeff[i] = 1.0; }
173  data = new std::vector<G4PhysicsVector*>;
174  data->resize(MAXZEL, 0);
175  }
176 
177  // it is possible re-initialisation for the second run
178  if(isMaster) {
179 
180  // check environment variable
181  // Build the complete string identifying the file with the data set
182  char* path = getenv("G4NEUTRONXSDATA");
183 
184  G4DynamicParticle* dynParticle =
186 
187  // Access to elements
188  const G4ElementTable* theElmTable = G4Element::GetElementTable();
189  size_t numOfElm = G4Element::GetNumberOfElements();
190  if(numOfElm > 0) {
191  for(size_t i=0; i<numOfElm; ++i) {
192  G4int Z = G4int(((*theElmTable)[i])->GetZ());
193  if(Z < 1) { Z = 1; }
194  else if(Z >= MAXZEL) { Z = MAXZEL-1; }
195  //G4cout << "Z= " << Z << G4endl;
196  // Initialisation
197  if(!(*data)[Z]) { Initialise(Z, dynParticle, path); }
198  }
199  }
200  delete dynParticle;
201  }
202 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
static G4double coeff[MAXZEL]
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
static std::vector< G4PhysicsVector * > * data
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
const G4int MAXZEL
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4NeutronElasticXS.cc.

95 {
96  outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
97  << "cross section on nuclei using data from the high precision\n"
98  << "neutron database. These data are simplified and smoothed over\n"
99  << "the resonance region in order to reduce CPU time.\n"
100  << "G4NeutronElasticXS is valid for energies up to 20 MeV, for all\n"
101  << "targets through U.\n";
102 }
Here is the caller graph for this function:

◆ Default_Name()

static const char* G4NeutronElasticXS::Default_Name ( )
inlinestatic

Definition at line 69 of file G4NeutronElasticXS.hh.

69 {return "G4NeutronElasticXS";}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 112 of file G4NeutronElasticXS.cc.

114 {
115  G4double xs = 0.0;
116  G4double ekin = aParticle->GetKineticEnergy();
117 
118  if(Z < 1 || Z >=MAXZEL) { return xs; }
119 
120  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
121  G4PhysicsVector* pv = (*data)[Z];
122  // G4cout << "G4NeutronElasticXS::GetCrossSection e= " << ekin
123  // << " Z= " << Z << G4endl;
124 
125  // element was not initialised
126  if(!pv) {
127  Initialise(Z);
128  pv = (*data)[Z];
129  if(!pv) { return xs; }
130  }
131 
132  G4double e1 = pv->Energy(0);
133  if(ekin <= e1) { return (*pv)[0]; }
134 
135  G4int n = pv->GetVectorLength() - 1;
136  G4double e2 = pv->Energy(n);
137 
138  if(ekin <= e2) {
139  xs = pv->Value(ekin);
140  } else if(1 == Z) {
143  } else {
144  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
146  }
147 
148  if(verboseLevel > 0){
149  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
150  }
151  return xs;
152 }
G4double GetElasticHadronNucleonXsc()
const G4ParticleDefinition * proton
G4HadronNucleonXsc * fNucleon
static G4double coeff[MAXZEL]
static const G4double e2
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)
Char_t n[5]
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ComponentGGHadronNucleusXsc * ggXsection
Float_t Z
size_t GetVectorLength() const
G4double Value(G4double theEnergy, size_t &lastidx) const
static const G4double e1
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
const G4int MAXZEL
void Initialise(G4int Z, G4DynamicParticle *dp=0, const char *=0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

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

Definition at line 205 of file G4NeutronElasticXS.cc.

207 {
208  if((*data)[Z]) { return; }
209  const char* path = p;
210  if(!p) {
211  // check environment variable
212  // Build the complete string identifying the file with the data set
213  path = getenv("G4NEUTRONXSDATA");
214  if (!path) {
215  G4Exception("G4NeutronElasticXS::Initialise(..)","had013",
217  "Environment variable G4NEUTRONXSDATA is not defined");
218  return;
219  }
220  }
221  G4DynamicParticle* dynParticle = dp;
222  if(!dp) {
223  dynParticle =
225  }
226 
227  G4int Amean = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
228 
229  // upload data from file
230  (*data)[Z] = new G4PhysicsLogVector();
231 
232  std::ostringstream ost;
233  ost << path << "/elast" << Z ;
234  std::ifstream filein(ost.str().c_str());
235  if (!(filein)) {
237  ed << "Data file <" << ost.str().c_str()
238  << "> is not opened!";
239  G4Exception("G4NeutronElasticXS::Initialise(..)","had014",
240  FatalException, ed, "Check G4NEUTRONXSDATA");
241  return;
242  }else{
243  if(verboseLevel > 1) {
244  G4cout << "file " << ost.str()
245  << " is opened by G4NeutronElasticXS" << G4endl;
246  }
247 
248  // retrieve data from DB
249  if(!((*data)[Z]->Retrieve(filein, true))) {
251  ed << "Data file <" << ost.str().c_str()
252  << "> is not retrieved!";
253  G4Exception("G4NeutronElasticXS::Initialise(..)","had015",
254  FatalException, ed, "Check G4NEUTRONXSDATA");
255  return;
256  }
257 
258  // smooth transition
259  size_t n = (*data)[Z]->GetVectorLength() - 1;
260  G4double emax = (*data)[Z]->Energy(n);
261  G4double sig1 = (*(*data)[Z])[n];
262  dynParticle->SetKineticEnergy(emax);
263  G4double sig2 = 0.0;
264  if(1 == Z) {
265  fNucleon->GetHadronNucleonXscPDG(dynParticle, proton);
267  } else {
268  ggXsection->GetIsoCrossSection(dynParticle, Z, Amean);
270  }
271  if(sig2 > 0.) { coeff[Z] = sig1/sig2; }
272  }
273  if(!dp) { delete dynParticle; }
274 }
G4double GetElasticHadronNucleonXsc()
const G4ParticleDefinition * proton
G4HadronNucleonXsc * fNucleon
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
static G4double coeff[MAXZEL]
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)
Char_t n[5]
G4GLOB_DLL std::ostream G4cout
G4ComponentGGHadronNucleusXsc * ggXsection
Float_t Z
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void SetKineticEnergy(G4double aEnergy)
static std::vector< G4PhysicsVector * > * data
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const G4double emax
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 105 of file G4NeutronElasticXS.cc.

107 {
108  return true;
109 }
Here is the caller graph for this function:

◆ operator=()

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

Member Data Documentation

◆ coeff

G4double G4NeutronElasticXS::coeff = {0.0}
staticprivate

Definition at line 97 of file G4NeutronElasticXS.hh.

◆ data

std::vector< G4PhysicsVector * > * G4NeutronElasticXS::data = 0
staticprivate

Definition at line 96 of file G4NeutronElasticXS.hh.

◆ fNucleon

G4HadronNucleonXsc* G4NeutronElasticXS::fNucleon
private

Definition at line 92 of file G4NeutronElasticXS.hh.

◆ ggXsection

G4ComponentGGHadronNucleusXsc* G4NeutronElasticXS::ggXsection
private

Definition at line 91 of file G4NeutronElasticXS.hh.

◆ isMaster

G4bool G4NeutronElasticXS::isMaster
private

Definition at line 99 of file G4NeutronElasticXS.hh.

◆ proton

const G4ParticleDefinition* G4NeutronElasticXS::proton
private

Definition at line 94 of file G4NeutronElasticXS.hh.


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