Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 ()
 

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

Definition at line 66 of file G4NeutronElasticXS.cc.

68  proton(G4Proton::Proton())
69 {
70  // verboseLevel = 0;
71  if(verboseLevel > 0){
72  G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
73  << MAXZEL << G4endl;
74  }
75  ggXsection = new G4ComponentGGHadronNucleusXsc();
76  fNucleon = new G4HadronNucleonXsc();
77  isMaster = false;
78 }
G4VCrossSectionDataSet(const G4String &nam="")
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const char * Default_Name()
#define G4endl
Definition: G4ios.hh:61
const G4int MAXZEL
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 }
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268
const G4int MAXZEL

Member Function Documentation

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
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
const XML_Char const XML_Char * data
Definition: expat.h:268
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
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:398
const G4int MAXZEL

Here is the call graph for this function:

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 }
static const char* G4NeutronElasticXS::Default_Name ( )
inlinestatic

Definition at line 69 of file G4NeutronElasticXS.hh.

69 {return "G4NeutronElasticXS";}

Here is the caller graph for this function:

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) {
141  fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
142  xs = coeff[1]*fNucleon->GetElasticHadronNucleonXsc();
143  } else {
144  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
145  xs = coeff[Z]*ggXsection->GetElasticGlauberGribovXsc();
146  }
147 
148  if(verboseLevel > 0){
149  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
150  }
151  return xs;
152 }
G4double GetElasticHadronNucleonXsc()
G4double GetKineticEnergy() const
size_t GetVectorLength() const
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)
G4GLOB_DLL std::ostream G4cout
const G4int n
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 G4int MAXZEL

Here is the call graph for this function:

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 }

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