Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 66241 2012-12-13 18:34:42Z gunter $
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 }
59 
61 {
63 }
64 
65 G4bool
67  G4int,
68  const G4Material*)
69 {
70  return false;
71 }
72 
73 G4bool
75  G4int, G4int,
76  const G4Element*,
77  const G4Material*)
78 {
79  return false;
80 }
81 
82 G4double
84  const G4Element* elm,
85  const G4Material* mat)
86 {
87  G4int Z = G4lrint(elm->GetZ());
88 
89  if (IsElementApplicable(part, Z, mat)) {
90  return GetElementCrossSection(part, Z, mat);
91  }
92 
93  // isotope-wise cross section making sum over available
94  // isotope cross sections, which may be incomplete, so
95  // the result is corrected
96  G4int nIso = elm->GetNumberOfIsotopes();
97  G4double fact = 0.0;
98  G4double xsec = 0.0;
99  G4Isotope* iso = 0;
100 
101  if (0 < nIso) {
102 
103  // user-defined isotope abundances
104  G4IsotopeVector* isoVector = elm->GetIsotopeVector();
105  G4double* abundVector = elm->GetRelativeAbundanceVector();
106 
107  for (G4int j = 0; j<nIso; ++j) {
108  iso = (*isoVector)[j];
109  G4int A = iso->GetN();
110  if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
111  fact += abundVector[j];
112  xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
113  }
114  }
115 
116  } else {
117 
118  // natural isotope abundances
120  G4int n0 = nist->GetNistFirstIsotopeN(Z);
121  G4int nn = nist->GetNumberOfNistIsotopes(Z);
122  for (G4int A = n0; A < n0+nn; ++A) {
123  G4double abund = nist->GetIsotopeAbundance(Z, A);
124  if(abund > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
125  fact += abund;
126  xsec += abund*GetIsoCrossSection(part, Z, A, iso, elm, mat);
127  }
128  }
129  }
130  if(fact > 0.0) { xsec /= fact; }
131  return xsec;
132 }
133 
134 G4double
136  G4int Z,
137  const G4Material* mat)
138 {
139  G4cout << "G4VCrossSectionDataSet::GetCrossSection per element ERROR: "
140  << " there is no cross section for "
141  << dynPart->GetDefinition()->GetParticleName()
142  << " E(MeV)= " << dynPart->GetKineticEnergy()/MeV;
143  if(mat) { G4cout << " inside " << mat->GetName(); }
144  G4cout << " for Z= " << Z << G4endl;
145  throw G4HadronicException(__FILE__, __LINE__,
146  "G4VCrossSectionDataSet::GetElementCrossSection is absent");
147  return 0.0;
148 }
149 
150 G4double
152  G4int Z, G4int A,
153  const G4Isotope*,
154  const G4Element* elm,
155  const G4Material* mat)
156 {
157  G4cout << "G4VCrossSectionDataSet::GetCrossSection per isotope ERROR: "
158  << " there is no cross section for "
159  << dynPart->GetDefinition()->GetParticleName()
160  << " E(MeV)= " << dynPart->GetKineticEnergy()/MeV;
161  if(mat) { G4cout << " inside " << mat->GetName(); }
162  if(elm) { G4cout << " for " << elm->GetName(); }
163  G4cout << " Z= " << Z << " A= " << A << G4endl;
164  throw G4HadronicException(__FILE__, __LINE__,
165  "G4VCrossSectionDataSet::GetIsoCrossSection is absent");
166  return 0.0;
167 }
168 
169 G4Isotope*
171 {
172  G4int nIso = anElement->GetNumberOfIsotopes();
173  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
174  G4Isotope* iso = (*isoVector)[0];
175 
176  // more than 1 isotope
177  if(1 < nIso) {
178  G4double* abundVector = anElement->GetRelativeAbundanceVector();
179  G4double sum = 0.0;
180  G4double q = G4UniformRand();
181  for (G4int j = 0; j<nIso; ++j) {
182  sum += abundVector[j];
183  if(q <= sum) {
184  iso = (*isoVector)[j];
185  break;
186  }
187  }
188  }
189  return iso;
190 }
191 
193 {}
194 
196 {}
197 
199 {
200  outFile << "The description for this cross section data set has not been written yet.\n";
201 }
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)
std::vector< G4Isotope * > G4IsotopeVector
G4double GetKineticEnergy() const
std::ofstream outFile
Definition: GammaRayTel.cc:68
void DeRegister(G4VCrossSectionDataSet *)
virtual void CrossSectionDescription(std::ostream &) const
G4VCrossSectionDataSet(const G4String &nam="")
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const G4String & GetName() const
Definition: G4Material.hh:176
G4double GetZ() const
Definition: G4Element.hh:131
const XML_Char * name
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:87
G4GLOB_DLL std::ostream G4cout
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
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)