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

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 93 of file G4NeutronInelasticXS.cc.

95  proton(G4Proton::Proton())
96 {
97  // verboseLevel = 0;
98  if(verboseLevel > 0){
99  G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
100  << MAXZINEL << G4endl;
101  }
102  ggXsection = new G4ComponentGGHadronNucleusXsc();
103  fNucleon = new G4HadronNucleonXsc();
104  isMaster = false;
105 }
G4VCrossSectionDataSet(const G4String &nam="")
const G4int MAXZINEL
static const char * Default_Name()
G4GLOB_DLL std::ostream G4cout
static G4Proton * Proton()
Definition: G4Proton.cc:93
#define G4endl
Definition: G4ios.hh:61
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 }
const XML_Char const XML_Char * data
Definition: expat.h:268

Member Function Documentation

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 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
const G4int MAXZINEL
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
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398

Here is the call graph for this function:

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

Definition at line 70 of file G4NeutronInelasticXS.hh.

70 {return "G4NeutronInelasticXS";}

Here is the caller graph for this function:

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 
153  G4PhysicsVector* pv = data->GetElementData(Z);
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) {
172  fNucleon->GetHadronNucleonXscPDG(aParticle, proton);
173  xs = coeff[1]*fNucleon->GetInelasticHadronNucleonXsc();
174  } else {
175  ggXsection->GetIsoCrossSection(aParticle, Z, Amean);
176  xs = coeff[Z]*ggXsection->GetInelasticGlauberGribovXsc();
177  }
178 
179  if(verboseLevel > 0) {
180  G4cout << "ekin= " << ekin << ", XSinel= " << xs << G4endl;
181  }
182  return xs;
183 }
G4double GetMaxEnergy() const
G4double GetKineticEnergy() const
const G4int MAXZINEL
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)
const XML_Char const XML_Char * data
Definition: expat.h:268
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
G4double GetInelasticHadronNucleonXsc()

Here is the call graph for this function:

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 }
G4double GetKineticEnergy() const
const G4int MAXZINEL
double A(double temperature)

Here is the call graph for this function:

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 }
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 }
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 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
std::vector< G4Isotope * > G4IsotopeVector
G4double GetZ() const
Definition: G4Element.hh:131
const G4int MAXZINEL
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268
#define G4UniformRand()
Definition: Randomize.hh:97
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
int G4lrint(double ad)
Definition: templates.hh:163
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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