Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "G4NistManager.hh"
55#include <algorithm>
56
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
59
61 : nist(G4NistManager::Instance())
62 , currentMaterial(nullptr)
63 , matParticle(nullptr)
64 , matKinEnergy(0.0)
65 , matCrossSection(0.0)
66 , nDataSetList(0)
67 , verboseLevel(0)
68{}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
71
73{}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
76
79 const G4Material* mat)
80{
81 if(mat == currentMaterial && part->GetDefinition() == matParticle
82 && part->GetKineticEnergy() == matKinEnergy)
83 { return matCrossSection; }
84
85 currentMaterial = mat;
86 matParticle = part->GetDefinition();
89
90 G4int nElements = mat->GetNumberOfElements();
91 const G4double* nAtomsPerVolume = mat->GetVecNbOfAtomsPerVolume();
92
93 if(G4int(xsecelm.size()) < nElements) { xsecelm.resize(nElements); }
94
95 for(G4int i=0; i<nElements; ++i) {
96 matCrossSection += nAtomsPerVolume[i] *
97 GetCrossSection(part, (*mat->GetElementVector())[i], mat);
99 }
100 return matCrossSection;
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
104
107 const G4Material* mat)
108{
109 if(part->GetKineticEnergy() == matKinEnergy && mat == currentMaterial &&
110 part->GetDefinition() == matParticle)
111 return matCrossSection;
112
113 currentMaterial = mat;
114 matParticle = part->GetDefinition();
116 matCrossSection = 0.0;
117
118 size_t nElements = mat->GetNumberOfElements();
119 const G4double* nAtomsPerVolume = mat->GetVecNbOfAtomsPerVolume();
120
121 if(xsecelm.size() < nElements) { xsecelm.resize(nElements); }
122
123 for(size_t i=0; i<nElements; ++i) {
124 matCrossSection += nAtomsPerVolume[i] *
125 GetCrossSection(part, mat->GetElement(i), mat);
127 }
128 return matCrossSection;
129}
130
131//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
132
134 const G4Element* elm,
135 const G4Material* mat)
136{
137 G4int i = nDataSetList - 1;
138 G4int Z = elm->GetZasInt();
139
140 if(elm->GetNaturalAbundanceFlag() &&
141 dataSetList[i]->IsElementApplicable(part, Z, mat))
142 {
143 // element wise cross section
144 return dataSetList[i]->GetElementCrossSection(part, Z, mat);
145 }
146
147 // isotope wise cross section
148 size_t nIso = elm->GetNumberOfIsotopes();
149
150 // user-defined isotope abundances
151 const G4double* abundVector = elm->GetRelativeAbundanceVector();
152
153 G4double sigma = 0.0;
154
155 for(size_t j = 0; j < nIso; ++j)
156 {
157 const G4Isotope* iso = elm->GetIsotope(j);
158 sigma += abundVector[j] *
159 GetIsoCrossSection(part, Z, iso->GetN(), iso, elm, mat, i);
160 }
161
162 return sigma;
163}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
166
169 G4int Z, G4int A,
170 const G4Isotope* iso,
171 const G4Element* elm,
172 const G4Material* mat,
173 G4int idx)
174{
175 // this methods is called after the check that dataSetList[idx]
176 // depend on isotopes, so for this DataSet only isotopes are checked
177
178 // isotope-wise cross section does exist
179 if(dataSetList[idx]->IsIsoApplicable(part, Z, A, elm, mat) ) {
180 return dataSetList[idx]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
181
182 } else {
183 // seach for other dataSet
184 for (G4int j = nDataSetList-1; j >= 0; --j) {
185 if (dataSetList[j]->IsElementApplicable(part, Z, mat)) {
186 return dataSetList[j]->GetElementCrossSection(part, Z, mat);
187 } else if (dataSetList[j]->IsIsoApplicable(part, Z, A, elm, mat)) {
188 return dataSetList[j]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
189 }
190 }
191 }
193 ed << "No isotope cross section found for "
194 << part->GetDefinition()->GetParticleName()
195 << " off Element " << elm->GetName()
196 << " in " << mat->GetName() << " Z= " << Z << " A= " << A
197 << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
198 G4Exception("G4CrossSectionDataStore::GetIsoCrossSection", "had001",
199 FatalException, ed);
200 return 0.0;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
204
207 G4int Z, G4int A,
208 const G4Isotope* iso,
209 const G4Element* elm,
210 const G4Material* mat)
211{
212 for (G4int i = nDataSetList-1; i >= 0; --i) {
213 if (dataSetList[i]->IsIsoApplicable(part, Z, A, elm, mat) ) {
214 return dataSetList[i]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
215 }
216 }
218 ed << "No isotope cross section found for "
219 << part->GetDefinition()->GetParticleName()
220 << " off Element " << elm->GetName()
221 << " in " << mat->GetName() << " Z= " << Z << " A= " << A
222 << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
223 G4Exception("G4CrossSectionDataStore::GetCrossSection", "had001",
224 FatalException, ed);
225 return 0.0;
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
229
230const G4Element*
232 const G4Material* mat,
233 G4Nucleus& target)
234{
235 size_t nElements = mat->GetNumberOfElements();
236 const G4Element* anElement = mat->GetElement(0);
237
238 // select element from a compound
239 if(1 < nElements) {
241 for(size_t i=0; i<nElements; ++i) {
242 if(cross <= xsecelm[i]) {
243 anElement = mat->GetElement(i);
244 break;
245 }
246 }
247 }
248
249 G4int Z = anElement->GetZasInt();
250 const G4Isotope* iso = nullptr;
251
252 G4int i = nDataSetList-1;
253 if (dataSetList[i]->IsElementApplicable(part, Z, mat)) {
254
255 //----------------------------------------------------------------
256 // element-wise cross section
257 // isotope cross section is not computed
258 //----------------------------------------------------------------
259 size_t nIso = anElement->GetNumberOfIsotopes();
260 iso = anElement->GetIsotope(0);
261
262 // more than 1 isotope
263 if(1 < nIso) {
264 iso = dataSetList[i]->SelectIsotope(anElement,
265 part->GetKineticEnergy(),
266 part->GetLogKineticEnergy());
267 }
268 } else {
269
270 //----------------------------------------------------------------
271 // isotope-wise cross section
272 // isotope cross section is computed
273 //----------------------------------------------------------------
274 size_t nIso = anElement->GetNumberOfIsotopes();
275 iso = anElement->GetIsotope(0);
276
277 // more than 1 isotope
278 if(1 < nIso) {
279 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
280 if(xseciso.size() < nIso) { xseciso.resize(nIso); }
281
282 G4double cross = 0.0;
283 size_t j;
284 for (j = 0; j<nIso; ++j) {
285 G4double xsec = 0.0;
286 if(abundVector[j] > 0.0) {
287 iso = anElement->GetIsotope(j);
288 xsec = abundVector[j]*
289 GetIsoCrossSection(part, Z, iso->GetN(), iso, anElement, mat, i);
290 }
291 cross += xsec;
292 xseciso[j] = cross;
293 }
294 cross *= G4UniformRand();
295 for (j = 0; j<nIso; ++j) {
296 if(cross <= xseciso[j]) {
297 iso = anElement->GetIsotope(j);
298 break;
299 }
300 }
301 }
302 }
303 target.SetIsotope(iso);
304 return anElement;
305}
306
307//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
308
309void
311{
312 if (nDataSetList == 0) {
314 ed << "No cross section is registered for "
315 << aParticleType.GetParticleName() << G4endl;
316 G4Exception("G4CrossSectionDataStore::BuildPhysicsTable", "had001",
317 FatalException, ed);
318 return;
319 }
320 for (G4int i=0; i<nDataSetList; ++i) {
321 dataSetList[i]->BuildPhysicsTable(aParticleType);
322 }
323}
324
325//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
326
327void
329{
330 // Print out all cross section data sets used and the energies at
331 // which they apply
332
333 if (nDataSetList == 0) {
334 G4cout << "WARNING - G4CrossSectionDataStore::DumpPhysicsTable: "
335 << " no data sets registered" << G4endl;
336 return;
337 }
338
339 for (G4int i = nDataSetList-1; i >= 0; --i) {
340 G4double e1 = dataSetList[i]->GetMinKinEnergy();
341 G4double e2 = dataSetList[i]->GetMaxKinEnergy();
342 G4cout
343 << " Cr_sctns: " << std::setw(25) << dataSetList[i]->GetName() << ": "
344 << G4BestUnit(e1, "Energy")
345 << " ---> "
346 << G4BestUnit(e2, "Energy") << "\n";
347 if (dataSetList[i]->GetName() == "G4CrossSectionPairGG") {
348 dataSetList[i]->DumpPhysicsTable(aParticleType);
349 }
350 }
351}
352
353//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
354#include <typeinfo>
356 std::ofstream& outFile) const
357{
358 // Write cross section data set info to html physics list
359 // documentation page
360
361 G4double ehi = 0;
362 G4double elo = 0;
363 G4String physListName(std::getenv("G4PhysListName"));
364 for (G4int i = nDataSetList-1; i > 0; i--) {
365 elo = dataSetList[i]->GetMinKinEnergy()/GeV;
366 ehi = dataSetList[i]->GetMaxKinEnergy()/GeV;
367 outFile << " <li><b><a href=\"" << physListName << "_"
368 << dataSetList[i]->GetName() << ".html\"> "
369 << dataSetList[i]->GetName() << "</a> from "
370 << elo << " GeV to " << ehi << " GeV </b></li>\n";
371 //G4cerr << i << ": XS for " << pD.GetParticleName() << " : " << dataSetList[i]->GetName()
372 // << " typeid : " << typeid(dataSetList[i]).name()<< G4endl;
374 }
375
376 G4double defaultHi = dataSetList[0]->GetMaxKinEnergy()/GeV;
377 if (ehi < defaultHi) {
378 outFile << " <li><b><a href=\"" << dataSetList[0]->GetName() << ".html\"> "
379 << dataSetList[0]->GetName() << "</a> from "
380 << ehi << " GeV to " << defaultHi << " GeV </b></li>\n";
382 }
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
386
388{
389 G4String dirName(std::getenv("G4PhysListDocDir"));
390 G4String physListName(std::getenv("G4PhysListName"));
391
392 G4String pathName = dirName + "/" + physListName + "_" + HtmlFileName(cs->GetName());
393 std::ofstream outCS;
394 outCS.open(pathName);
395 outCS << "<html>\n";
396 outCS << "<head>\n";
397 outCS << "<title>Description of " << cs->GetName()
398 << "</title>\n";
399 outCS << "</head>\n";
400 outCS << "<body>\n";
401
402 cs->CrossSectionDescription(outCS);
403
404 outCS << "</body>\n";
405 outCS << "</html>\n";
406
407}
408//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
409
411{
412 G4String str(in);
413 // replace blanks by _ C++11 version:
414 std::transform(str.begin(), str.end(), str.begin(), [](char ch) {
415 return ch == ' ' ? '_' : ch;
416 });
417 str=str + ".html";
418 return str;
419}
420
421//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
422
424{
425 if(p->ForAllAtomsAndEnergies()) {
426 dataSetList.clear();
427 nDataSetList = 0;
428 }
429 dataSetList.push_back(p);
430 ++nDataSetList;
431}
432
433
434
435//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
436
438{
439 if(p->ForAllAtomsAndEnergies()) {
440 dataSetList.clear();
441 dataSetList.push_back(p);
442 nDataSetList = 1;
443 } else {
444 if ( i > dataSetList.size() ) i = dataSetList.size();
445 std::vector< G4VCrossSectionDataSet* >::iterator it = dataSetList.end() - i;
446 dataSetList.insert(it , p);
447 ++nDataSetList;
448 }
449}
450//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
static const G4double e1[44]
static const G4double e2[44]
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
#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
G4String HtmlFileName(const G4String &in) const
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *aMaterial, G4int index)
void BuildPhysicsTable(const G4ParticleDefinition &)
void AddDataSet(G4VCrossSectionDataSet *)
std::vector< G4VCrossSectionDataSet * > dataSetList
void PrintCrossSectionHtml(const G4VCrossSectionDataSet *cs) const
void DumpPhysicsTable(const G4ParticleDefinition &)
const G4ParticleDefinition * matParticle
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
std::vector< G4double > xseciso
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
std::vector< G4double > xsecelm
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:263
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 G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
const G4String & GetName() const
Definition: G4Material.hh:173
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:114
const G4String & GetParticleName() const
const G4String & GetName() const
virtual void CrossSectionDescription(std::ostream &) const
G4bool transform(G4String &input, const G4String &type)