Geant4  10.02.p01
G4VCrossSectionDataSet.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4VCrossSectionDataSet.cc 89024 2015-03-18 08:17:25Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4VCrossSectionDataSet
34 //
35 // Author F.W. Jones, TRIUMF, 20-JAN-97
36 //
37 // Modifications:
38 // 23.01.2009 V.Ivanchenko move constructor and destructor to source
39 // 12.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class
40 //
41 
43 #include "G4SystemOfUnits.hh"
45 #include "G4DynamicParticle.hh"
46 #include "G4Material.hh"
47 #include "G4Element.hh"
48 #include "G4Isotope.hh"
49 #include "G4NistManager.hh"
50 #include "G4HadronicException.hh"
51 #include "G4HadTmpUtil.hh"
52 #include "Randomize.hh"
53 
55  verboseLevel(0),minKinEnergy(0.0),maxKinEnergy(100*TeV),name(nam)
56 {
58  registry->Register(this);
59 }
60 
62 {
63  registry->DeRegister(this);
64 }
65 
66 G4bool
68  G4int,
69  const G4Material*)
70 {
71  return false;
72 }
73 
74 G4bool
76  G4int, G4int,
77  const G4Element*,
78  const G4Material*)
79 {
80  return false;
81 }
82 
83 G4double
85  const G4Element* elm,
86  const G4Material* mat)
87 {
88  G4int Z = G4lrint(elm->GetZ());
89 
90  if (IsElementApplicable(part, Z, mat)) {
91  return GetElementCrossSection(part, Z, mat);
92  }
93 
94  // isotope-wise cross section making sum over available
95  // isotope cross sections, which may be incomplete, so
96  // the result is corrected
97  G4int nIso = elm->GetNumberOfIsotopes();
98  G4double fact = 0.0;
99  G4double xsec = 0.0;
100  G4Isotope* iso = 0;
101 
102  if (0 < nIso) {
103 
104  // user-defined isotope abundances
105  G4IsotopeVector* isoVector = elm->GetIsotopeVector();
106  G4double* abundVector = elm->GetRelativeAbundanceVector();
107 
108  for (G4int j = 0; j<nIso; ++j) {
109  iso = (*isoVector)[j];
110  G4int A = iso->GetN();
111  if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
112  fact += abundVector[j];
113  xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
114  }
115  }
116 
117  } else {
118 
119  // natural isotope abundances
121  G4int n0 = nist->GetNistFirstIsotopeN(Z);
122  G4int nn = nist->GetNumberOfNistIsotopes(Z);
123  for (G4int A = n0; A < n0+nn; ++A) {
124  G4double abund = nist->GetIsotopeAbundance(Z, A);
125  if(abund > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
126  fact += abund;
127  xsec += abund*GetIsoCrossSection(part, Z, A, iso, elm, mat);
128  }
129  }
130  }
131  if(fact > 0.0) { xsec /= fact; }
132  return xsec;
133 }
134 
135 G4double
137  G4int Z,
138  const G4Material* mat)
139 {
140  G4cout << "G4VCrossSectionDataSet::GetCrossSection per element ERROR: "
141  << " there is no cross section for "
142  << dynPart->GetDefinition()->GetParticleName()
143  << " E(MeV)= " << dynPart->GetKineticEnergy()/MeV;
144  if(mat) { G4cout << " inside " << mat->GetName(); }
145  G4cout << " for Z= " << Z << G4endl;
146  throw G4HadronicException(__FILE__, __LINE__,
147  "G4VCrossSectionDataSet::GetElementCrossSection is absent");
148  return 0.0;
149 }
150 
151 G4double
153  G4int Z, G4int A,
154  const G4Isotope*,
155  const G4Element* elm,
156  const G4Material* mat)
157 {
158  G4cout << "G4VCrossSectionDataSet::GetCrossSection per isotope ERROR: "
159  << " there is no cross section for "
160  << dynPart->GetDefinition()->GetParticleName()
161  << " E(MeV)= " << dynPart->GetKineticEnergy()/MeV;
162  if(mat) { G4cout << " inside " << mat->GetName(); }
163  if(elm) { G4cout << " for " << elm->GetName(); }
164  G4cout << " Z= " << Z << " A= " << A << G4endl;
165  throw G4HadronicException(__FILE__, __LINE__,
166  "G4VCrossSectionDataSet::GetIsoCrossSection is absent");
167  return 0.0;
168 }
169 
170 G4Isotope*
172 {
173  G4int nIso = anElement->GetNumberOfIsotopes();
174  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
175  G4Isotope* iso = (*isoVector)[0];
176 
177  // more than 1 isotope
178  if(1 < nIso) {
179  G4double* abundVector = anElement->GetRelativeAbundanceVector();
180  G4double sum = 0.0;
181  G4double q = G4UniformRand();
182  for (G4int j = 0; j<nIso; ++j) {
183  sum += abundVector[j];
184  if(q <= sum) {
185  iso = (*isoVector)[j];
186  break;
187  }
188  }
189  }
190  return iso;
191 }
192 
194 {}
195 
197 {}
198 
199 void G4VCrossSectionDataSet::CrossSectionDescription(std::ostream& outFile) const
200 {
201  outFile << "The description for this cross section data set has not been written yet.\n";
202 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
static const double MeV
Definition: G4SIunits.hh:211
std::vector< G4Isotope * > G4IsotopeVector
const G4double fact
G4double GetKineticEnergy() const
void DeRegister(G4VCrossSectionDataSet *)
virtual void CrossSectionDescription(std::ostream &) const
G4VCrossSectionDataSet(const G4String &nam="")
G4String name
Definition: TRTMaterials.hh:40
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const G4String & GetName() const
Definition: G4Material.hh:178
G4double GetZ() const
Definition: G4Element.hh:131
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4int GetNistFirstIsotopeN(G4int Z) const
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
G4int GetN() const
Definition: G4Isotope.hh:94
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
static G4CrossSectionDataSetRegistry * Instance()
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetIsotopeAbundance(G4int Z, G4int N) const
G4int GetNumberOfNistIsotopes(G4int Z) const
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
void Register(G4VCrossSectionDataSet *)
#define G4endl
Definition: G4ios.hh:61
G4CrossSectionDataSetRegistry * registry
static const double TeV
Definition: G4SIunits.hh:215
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
virtual G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy)