Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionDataStore.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4CrossSectionDataStore
33//
34// Modifications:
35// 23.01.2009 V.Ivanchenko add destruction of data sets
36// 29.04.2010 G.Folger modifictaions for integer A & Z
37// 14.03.2011 V.Ivanchenko fixed DumpPhysicsTable
38// 15.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class
39// 07.03.2013 M.Maire cosmetic in DumpPhysicsTable
40//
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
43
45#include "G4SystemOfUnits.hh"
46#include "G4UnitsTable.hh"
47#include "Randomize.hh"
48#include "G4Nucleus.hh"
49
50#include "G4DynamicParticle.hh"
51#include "G4Isotope.hh"
52#include "G4Element.hh"
53#include "G4Material.hh"
54#include "G4MaterialTable.hh"
55#include "G4NistManager.hh"
56#include <algorithm>
57#include <typeinfo>
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
60
62 : nist(G4NistManager::Instance())
63{}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
66
69 const G4Material* mat)
70{
71 currentMaterial = mat;
72 matParticle = dp->GetDefinition();
73 matKinEnergy = dp->GetKineticEnergy();
74 matCrossSection = 0.0;
75
76 std::size_t nElements = mat->GetNumberOfElements();
77 const G4double* nAtomsPerVolume = mat->GetVecNbOfAtomsPerVolume();
78
79 if(xsecelm.size() < nElements) { xsecelm.resize(nElements); }
80
81 for(G4int i=0; i<(G4int)nElements; ++i) {
82 G4double xs =
83 nAtomsPerVolume[i]*GetCrossSection(dp, mat->GetElement(i), mat);
84 matCrossSection += std::max(xs, 0.0);
85 xsecelm[i] = matCrossSection;
86 }
87 return matCrossSection;
88}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
91
93 const G4Element* elm,
94 const G4Material* mat)
95{
96 // first check the most last cross section
97 G4int i = nDataSetList-1;
98 G4int Z = elm->GetZasInt();
99
100 if(elm->GetNaturalAbundanceFlag() &&
101 dataSetList[i]->IsElementApplicable(dp, Z, mat))
102 {
103 // element wise cross section
104 return dataSetList[i]->GetElementCrossSection(dp, Z, mat);
105 }
106
107 // isotope wise cross section
108 G4int nIso = (G4int)elm->GetNumberOfIsotopes();
109
110 // user-defined isotope abundances
111 const G4double* abundVector = elm->GetRelativeAbundanceVector();
112
113 G4double sigma = 0.0;
114
115 // isotope and element wise cross sections
116 for(G4int j = 0; j < nIso; ++j)
117 {
118 const G4Isotope* iso = elm->GetIsotope(j);
119 sigma += abundVector[j] *
120 GetIsoCrossSection(dp, Z, iso->GetN(), iso, elm, mat, i);
121 }
122 return sigma;
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
126
128G4CrossSectionDataStore::GetIsoCrossSection(const G4DynamicParticle* dp,
129 G4int Z, G4int A,
130 const G4Isotope* iso,
131 const G4Element* elm,
132 const G4Material* mat,
133 G4int idx)
134{
135 // this methods is called after the check that dataSetList[idx]
136 // depend on isotopes, so first isotopes are checked
137 if(dataSetList[idx]->IsIsoApplicable(dp, Z, A, elm, mat) ) {
138 return dataSetList[idx]->GetIsoCrossSection(dp, Z, A, iso, elm, mat);
139 }
140
141 // no isotope wise cross section - check other datasets
142 for (G4int j = nDataSetList-1; j >= 0; --j) {
143 if(dataSetList[j]->IsElementApplicable(dp, Z, mat)) {
144 return dataSetList[j]->GetElementCrossSection(dp, Z, mat);
145 } else if (dataSetList[j]->IsIsoApplicable(dp, Z, A, elm, mat)) {
146 return dataSetList[j]->GetIsoCrossSection(dp, Z, A, iso, elm, mat);
147 }
148 }
150 ed << "No isotope cross section found for "
152 << " off target Element " << elm->GetName()
153 << " Z= " << Z << " A= " << A;
154 if(nullptr != mat) ed << " from " << mat->GetName();
155 ed << " E(MeV)=" << dp->GetKineticEnergy()/MeV << G4endl;
156 G4Exception("G4CrossSectionDataStore::GetIsoCrossSection", "had001",
157 FatalException, ed);
158 return 0.0;
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
162
165 G4int Z, G4int A,
166 const G4Isotope* iso,
167 const G4Element* elm,
168 const G4Material* mat)
169{
170 for (G4int i = nDataSetList-1; i >= 0; --i) {
171 if (dataSetList[i]->IsIsoApplicable(dp, Z, A, elm, mat) ) {
172 return dataSetList[i]->GetIsoCrossSection(dp, Z, A, iso, elm, mat);
173 } else if(dataSetList[i]->IsElementApplicable(dp, Z, mat)) {
174 return dataSetList[i]->GetElementCrossSection(dp, Z, mat);
175 }
176 }
178 ed << "No isotope cross section found for "
180 << " off target Element " << elm->GetName()
181 << " Z= " << Z << " A= " << A;
182 if(nullptr != mat) ed << " from " << mat->GetName();
183 ed << " E(MeV)=" << dp->GetKineticEnergy()/MeV << G4endl;
184 G4Exception("G4CrossSectionDataStore::GetCrossSection", "had001",
185 FatalException, ed);
186 return 0.0;
187}
188
189//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
190
191const G4Element*
193 const G4Material* mat,
194 G4Nucleus& target)
195{
196 if(nullptr != forcedElement) { return forcedElement; }
197 std::size_t nElements = mat->GetNumberOfElements();
198 const G4Element* anElement = mat->GetElement(0);
199
200 // select element from a compound
201 if(1 < nElements) {
202 G4double cross = matCrossSection*G4UniformRand();
203 for(G4int i=0; i<(G4int)nElements; ++i) {
204 if(cross <= xsecelm[i]) {
205 anElement = mat->GetElement(i);
206 break;
207 }
208 }
209 }
210
211 G4int Z = anElement->GetZasInt();
212 const G4Isotope* iso = nullptr;
213
214 G4int i = nDataSetList-1;
215 if (dataSetList[i]->IsElementApplicable(dp, Z, mat)) {
216
217 //----------------------------------------------------------------
218 // element-wise cross section
219 // isotope cross section is not computed
220 //----------------------------------------------------------------
221 std::size_t nIso = anElement->GetNumberOfIsotopes();
222 iso = anElement->GetIsotope(0);
223
224 // more than 1 isotope
225 if(1 < nIso) {
226 iso = dataSetList[i]->SelectIsotope(anElement,
227 dp->GetKineticEnergy(),
228 dp->GetLogKineticEnergy());
229 }
230 } else {
231
232 //----------------------------------------------------------------
233 // isotope-wise cross section
234 // isotope cross section is computed
235 //----------------------------------------------------------------
236 std::size_t nIso = anElement->GetNumberOfIsotopes();
237 iso = anElement->GetIsotope(0);
238
239 // more than 1 isotope
240 if(1 < nIso) {
241 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
242 if(xseciso.size() < nIso) { xseciso.resize(nIso); }
243
244 G4double cross = 0.0;
245 G4int j;
246 for (j = 0; j<(G4int)nIso; ++j) {
247 G4double xsec = 0.0;
248 if(abundVector[j] > 0.0) {
249 iso = anElement->GetIsotope(j);
250 xsec = abundVector[j]*
251 GetIsoCrossSection(dp, Z, iso->GetN(), iso, anElement, mat, i);
252 }
253 cross += xsec;
254 xseciso[j] = cross;
255 }
256 cross *= G4UniformRand();
257 for (j = 0; j<(G4int)nIso; ++j) {
258 if(cross <= xseciso[j]) {
259 iso = anElement->GetIsotope(j);
260 break;
261 }
262 }
263 }
264 }
265 target.SetIsotope(iso);
266 return anElement;
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
270
271void
273{
274 if (nDataSetList == 0) {
276 ed << "No cross section is registered for "
277 << part.GetParticleName() << G4endl;
278 G4Exception("G4CrossSectionDataStore::BuildPhysicsTable", "had001",
279 FatalException, ed);
280 return;
281 }
282 matParticle = &part;
283 for (G4int i=0; i<nDataSetList; ++i) {
284 dataSetList[i]->BuildPhysicsTable(part);
285 }
286 const G4MaterialTable* theMatTable = G4Material::GetMaterialTable();
287 std::size_t nelm = 0;
288 std::size_t niso = 0;
289 for(auto mat : *theMatTable) {
290 std::size_t nElements = mat->GetNumberOfElements();
291 nelm = std::max(nelm, nElements);
292 for(G4int j=0; j<(G4int)nElements; ++j) {
293 niso = std::max(niso, mat->GetElement(j)->GetNumberOfIsotopes());
294 }
295 }
296 // define vectors for a run
297 xsecelm.resize(nelm, 0.0);
298 xseciso.resize(niso, 0.0);
299}
300
301//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
302
303void
305{
306 // Print out all cross section data sets used and the energies at
307 // which they apply
308
309 if (nDataSetList == 0) {
310 G4cout << "WARNING - G4CrossSectionDataStore::DumpPhysicsTable: "
311 << " no data sets registered" << G4endl;
312 return;
313 }
314
315 for (G4int i = nDataSetList-1; i >= 0; --i) {
316 G4double e1 = dataSetList[i]->GetMinKinEnergy();
317 G4double e2 = dataSetList[i]->GetMaxKinEnergy();
318 G4cout
319 << " Cr_sctns: " << std::setw(25) << dataSetList[i]->GetName() << ": "
320 << G4BestUnit(e1, "Energy") << " ---> "
321 << G4BestUnit(e2, "Energy") << "\n";
322 if (dataSetList[i]->GetName() == "G4CrossSectionPairGG") {
323 dataSetList[i]->DumpPhysicsTable(part);
324 }
325 G4cout << G4endl;
326 }
327}
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
330
332 std::ofstream& outFile) const
333{
334 // Write cross section data set info to html physics list
335 // documentation page
336
337 G4double ehi = 0;
338 G4double elo = 0;
339 G4String physListName(std::getenv("G4PhysListName"));
340 for (G4int i = nDataSetList-1; i > 0; i--) {
341 elo = dataSetList[i]->GetMinKinEnergy()/GeV;
342 ehi = dataSetList[i]->GetMaxKinEnergy()/GeV;
343 outFile << " <li><b><a href=\"" << physListName << "_"
344 << dataSetList[i]->GetName() << ".html\"> "
345 << dataSetList[i]->GetName() << "</a> from "
346 << elo << " GeV to " << ehi << " GeV </b></li>\n";
347 PrintCrossSectionHtml(dataSetList[i]);
348 }
349
350 G4double defaultHi = dataSetList[0]->GetMaxKinEnergy()/GeV;
351 if (ehi < defaultHi) {
352 outFile << " <li><b><a href=\"" << dataSetList[0]->GetName()
353 << ".html\"> "
354 << dataSetList[0]->GetName() << "</a> from "
355 << ehi << " GeV to " << defaultHi << " GeV </b></li>\n";
356 PrintCrossSectionHtml(dataSetList[0]);
357 }
358}
359
360//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
361
363{
364 G4String dirName(std::getenv("G4PhysListDocDir"));
365 G4String physListName(std::getenv("G4PhysListName"));
366
367 G4String pathName = dirName + "/" + physListName + "_" + HtmlFileName(cs->GetName());
368 std::ofstream outCS;
369 outCS.open(pathName);
370 outCS << "<html>\n";
371 outCS << "<head>\n";
372 outCS << "<title>Description of " << cs->GetName()
373 << "</title>\n";
374 outCS << "</head>\n";
375 outCS << "<body>\n";
376
377 cs->CrossSectionDescription(outCS);
378
379 outCS << "</body>\n";
380 outCS << "</html>\n";
381}
382
383//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
384
385G4String G4CrossSectionDataStore::HtmlFileName(const G4String & in) const
386{
387 G4String str(in);
388 // replace blanks by _ C++11 version:
389 std::transform(str.begin(), str.end(), str.begin(), [](char ch) {
390 return ch == ' ' ? '_' : ch;
391 });
392 str=str + ".html";
393 return str;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
397
399{
400 if(p->ForAllAtomsAndEnergies()) {
401 dataSetList.clear();
402 nDataSetList = 0;
403 }
404 dataSetList.push_back(p);
405 ++nDataSetList;
406}
407
408//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
409
411{
412 if(p->ForAllAtomsAndEnergies()) {
413 dataSetList.clear();
414 dataSetList.push_back(p);
415 nDataSetList = 1;
416 } else if ( i >= dataSetList.size() ) {
417 dataSetList.push_back(p);
418 ++nDataSetList;
419 } else {
420 std::vector< G4VCrossSectionDataSet* >::iterator it = dataSetList.end() - i;
421 dataSetList.insert(it , p);
422 ++nDataSetList;
423 }
424}
425
426//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::vector< G4Material * > G4MaterialTable
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void DumpHtml(const G4ParticleDefinition &, std::ofstream &) const
void BuildPhysicsTable(const G4ParticleDefinition &)
void AddDataSet(G4VCrossSectionDataSet *)
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
G4double GetLogKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4bool GetNaturalAbundanceFlag() const
Definition: G4Element.hh:262
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4String & GetName() const
Definition: G4Element.hh:127
G4int GetZasInt() const
Definition: G4Element.hh:132
G4int GetN() const
Definition: G4Isotope.hh:93
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:114
const G4String & GetParticleName() const
const G4String & GetName() const
virtual void CrossSectionDescription(std::ostream &) const