Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VCrossSectionDataSet.cc
이 파일의 문서화 페이지로 가기
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 110787 2018-06-14 06:43:31Z 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),
56  isForAllAtomsAndEnergies(false),name(nam)
57 {
59  registry->Register(this);
60 }
61 
63 {
64  registry->DeRegister(this);
65 }
66 
67 G4bool
69  G4int,
70  const G4Material*)
71 {
72  return false;
73 }
74 
75 G4bool
77  G4int, G4int,
78  const G4Element*,
79  const G4Material*)
80 {
81  return false;
82 }
83 
84 G4double
86  const G4Element* elm,
87  const G4Material* mat)
88 {
89  G4int Z = elm->GetZasInt();
90 
91  if (IsElementApplicable(part, Z, mat)) {
92  return GetElementCrossSection(part, Z, mat);
93  }
94 
95  // isotope-wise cross section making sum over available
96  // isotope cross sections, which may be incomplete, so
97  // the result is corrected
98  G4int nIso = elm->GetNumberOfIsotopes();
99  G4double fact = 0.0;
100  G4double xsec = 0.0;
101  const G4Isotope* iso = nullptr;
102 
103  // user-defined isotope abundances
104  const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
105  const 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  if(fact > 0.0) { xsec /= fact; }
116  return xsec;
117 }
118 
119 G4double
121  G4int Z,
122  const G4Material* mat)
123 {
124  G4cout << "G4VCrossSectionDataSet::GetCrossSection per element ERROR: "
125  << " there is no cross section for "
126  << dynPart->GetDefinition()->GetParticleName()
127  << " E(MeV)= " << dynPart->GetKineticEnergy()/MeV;
128  if(mat) { G4cout << " inside " << mat->GetName(); }
129  G4cout << " for Z= " << Z << G4endl;
130  throw G4HadronicException(__FILE__, __LINE__,
131  "G4VCrossSectionDataSet::GetElementCrossSection is absent");
132 }
133 
134 G4double
136  G4int Z, G4int A,
137  const G4Isotope*,
138  const G4Element* elm,
139  const G4Material* mat)
140 {
141  G4cout << "G4VCrossSectionDataSet::GetCrossSection per isotope ERROR: "
142  << " there is no cross section for "
143  << dynPart->GetDefinition()->GetParticleName()
144  << " E(MeV)= " << dynPart->GetKineticEnergy()/MeV;
145  if(mat) { G4cout << " inside " << mat->GetName(); }
146  if(elm) { G4cout << " for " << elm->GetName(); }
147  G4cout << " Z= " << Z << " A= " << A << G4endl;
148  throw G4HadronicException(__FILE__, __LINE__,
149  "G4VCrossSectionDataSet::GetIsoCrossSection is absent");
150 }
151 
152 const G4Isotope*
154 {
155  G4int nIso = anElement->GetNumberOfIsotopes();
156  const G4Isotope* iso = anElement->GetIsotope(0);
157 
158  // more than 1 isotope
159  if(1 < nIso) {
160  const G4double* abundVector = anElement->GetRelativeAbundanceVector();
161  G4double sum = 0.0;
162  G4double q = G4UniformRand();
163  for (G4int j = 0; j<nIso; ++j) {
164  sum += abundVector[j];
165  if(q <= sum) {
166  iso = anElement->GetIsotope(j);
167  break;
168  }
169  }
170  }
171  return iso;
172 }
173 
175 {}
176 
178 {}
179 
180 void G4VCrossSectionDataSet::CrossSectionDescription(std::ostream& outFile) const
181 {
182  outFile << "The description for this cross section data set has not been written yet.\n";
183 }
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
const XML_Char * name
Definition: expat.h:151
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4CrossSectionDataSetRegistry * registry
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4String & GetName() const
Definition: G4Element.hh:127
#define G4endl
Definition: G4ios.hh:61
const G4String & GetParticleName() const
static G4CrossSectionDataSetRegistry * Instance()
const G4String & GetName() const
Definition: G4Material.hh:179
Float_t Z
static constexpr double TeV
Definition: G4SIunits.hh:218
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
virtual const G4Isotope * SelectIsotope(const G4Element *, G4double kinEnergy)
TString part[npart]
Definition: Style.C:32
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
G4ParticleDefinition * GetDefinition() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
std::vector< G4Isotope * > G4IsotopeVector
G4VCrossSectionDataSet(const G4String &nam="")
#define G4UniformRand()
Definition: Randomize.hh:53
double A(double temperature)
void DeRegister(G4VCrossSectionDataSet *)
Float_t mat
G4int GetZasInt() const
Definition: G4Element.hh:132
G4int GetN() const
Definition: G4Isotope.hh:94
int G4int
Definition: G4Types.hh:78
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
G4double GetKineticEnergy() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
void Register(G4VCrossSectionDataSet *)
G4GLOB_DLL std::ostream G4cout
virtual void CrossSectionDescription(std::ostream &) const
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159