Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 // $Id: G4CrossSectionDataStore.cc 79192 2014-02-20 10:07:18Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4CrossSectionDataStore
34 //
35 // Modifications:
36 // 23.01.2009 V.Ivanchenko add destruction of data sets
37 // 29.04.2010 G.Folger modifictaions for integer A & Z
38 // 14.03.2011 V.Ivanchenko fixed DumpPhysicsTable
39 // 15.08.2011 G.Folger, V.Ivanchenko, T.Koi, D.Wright redesign the class
40 // 07.03.2013 M.Maire cosmetic in DumpPhysicsTable
41 //
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
44 
46 #include "G4SystemOfUnits.hh"
47 #include "G4UnitsTable.hh"
48 #include "G4HadronicException.hh"
49 #include "G4HadTmpUtil.hh"
50 #include "Randomize.hh"
51 #include "G4Nucleus.hh"
52 
53 #include "G4DynamicParticle.hh"
54 #include "G4Isotope.hh"
55 #include "G4Element.hh"
56 #include "G4Material.hh"
57 #include "G4NistManager.hh"
58 #include <iostream>
59 
60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
61 
63  nDataSetList(0), verboseLevel(0)
64 {
65  nist = G4NistManager::Instance();
66  currentMaterial = elmMaterial = 0;
67  currentElement = 0; //ALB 14-Aug-2012 Coverity fix.
68  matParticle = elmParticle = 0;
69  matKinEnergy = elmKinEnergy = matCrossSection = elmCrossSection = 0.0;
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
73 
75 {}
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
78 
81  const G4Material* mat)
82 {
83  if(mat == currentMaterial && part->GetDefinition() == matParticle
84  && part->GetKineticEnergy() == matKinEnergy)
85  { return matCrossSection; }
86 
87  currentMaterial = mat;
88  matParticle = part->GetDefinition();
89  matKinEnergy = part->GetKineticEnergy();
90  matCrossSection = 0;
91 
92  G4int nElements = mat->GetNumberOfElements();
93  const G4double* nAtomsPerVolume = mat->GetVecNbOfAtomsPerVolume();
94 
95  if(G4int(xsecelm.size()) < nElements) { xsecelm.resize(nElements); }
96 
97  for(G4int i=0; i<nElements; ++i) {
98  matCrossSection += nAtomsPerVolume[i] *
99  GetCrossSection(part, (*mat->GetElementVector())[i], mat);
100  xsecelm[i] = matCrossSection;
101  }
102  return matCrossSection;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
106 
107 G4double
109  const G4Element* elm,
110  const G4Material* mat)
111 {
112  if(mat == elmMaterial && elm == currentElement &&
113  part->GetDefinition() == elmParticle &&
114  part->GetKineticEnergy() == elmKinEnergy)
115  { return elmCrossSection; }
116 
117  elmMaterial = mat;
118  currentElement = elm;
119  elmParticle = part->GetDefinition();
120  elmKinEnergy = part->GetKineticEnergy();
121  elmCrossSection = 0.0;
122 
123  G4int i = nDataSetList-1;
124  G4int Z = G4lrint(elm->GetZ());
125  if (elm->GetNaturalAbundanceFlag() &&
126  dataSetList[i]->IsElementApplicable(part, Z, mat)) {
127 
128  // element wise cross section
129  elmCrossSection = dataSetList[i]->GetElementCrossSection(part, Z, mat);
130 
131  //G4cout << "Element wise " << elmParticle->GetParticleName()
132  // << " xsec(barn)= " << elmCrossSection/barn
133  // << " E(MeV)= " << elmKinEnergy/MeV
134  // << " Z= " << Z << " AbundFlag= " << elm->GetNaturalAbandancesFlag()
135  // <<G4endl;
136 
137  } else {
138  // isotope wise cross section
139  G4int nIso = elm->GetNumberOfIsotopes();
140  G4Isotope* iso = 0;
141 
142  // user-defined isotope abundances
143  G4IsotopeVector* isoVector = elm->GetIsotopeVector();
144  G4double* abundVector = elm->GetRelativeAbundanceVector();
145 
146  for (G4int j = 0; j<nIso; ++j) {
147  if(abundVector[j] > 0.0) {
148  iso = (*isoVector)[j];
149  elmCrossSection += abundVector[j]*
150  GetIsoCrossSection(part, Z, iso->GetN(), iso, elm, mat, i);
151  //G4cout << "Isotope wise " << elmParticle->GetParticleName()
152  // << " xsec(barn)= " << elmCrossSection/barn
153  // << " E(MeV)= " << elmKinEnergy/MeV
154  // << " Z= " << Z << " A= " << iso->GetN() << " j= " << j << G4endl;
155  }
156  }
157  }
158  //G4cout << " E(MeV)= " << elmKinEnergy/MeV
159  // << "xsec(barn)= " << elmCrossSection/barn <<G4endl;
160  return elmCrossSection;
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
164 
165 G4double
166 G4CrossSectionDataStore::GetIsoCrossSection(const G4DynamicParticle* part,
167  G4int Z, G4int A,
168  const G4Isotope* iso,
169  const G4Element* elm,
170  const G4Material* mat,
171  G4int idx)
172 {
173  // this methods is called after the check that dataSetList[idx]
174  // depend on isotopes, so for this DataSet only isotopes are checked
175 
176  // isotope-wise cross section does exist
177  if(dataSetList[idx]->IsIsoApplicable(part, Z, A, elm, mat) ) {
178  return dataSetList[idx]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
179 
180  } else {
181  // seach for other dataSet
182  for (G4int j = nDataSetList-1; j >= 0; --j) {
183  if (dataSetList[j]->IsElementApplicable(part, Z, mat)) {
184  return dataSetList[j]->GetElementCrossSection(part, Z, mat);
185  } else if (dataSetList[j]->IsIsoApplicable(part, Z, A, elm, mat)) {
186  return dataSetList[j]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
187  }
188  }
189  }
190  G4cout << "G4CrossSectionDataStore::GetCrossSection ERROR: "
191  << " no isotope cross section found"
192  << G4endl;
193  G4cout << " for " << part->GetDefinition()->GetParticleName()
194  << " off Element " << elm->GetName()
195  << " in " << mat->GetName()
196  << " Z= " << Z << " A= " << A
197  << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
198  throw G4HadronicException(__FILE__, __LINE__,
199  " no applicable data set found for the isotope");
200  return 0.0;
201  //return dataSetList[idx]->ComputeCrossSection(part, elm, mat);
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
205 
206 G4double
208  G4int Z, G4int A,
209  const G4Isotope* iso,
210  const G4Element* elm,
211  const G4Material* mat)
212 {
213  for (G4int i = nDataSetList-1; i >= 0; --i) {
214  if (dataSetList[i]->IsIsoApplicable(part, Z, A, elm, mat) ) {
215  return dataSetList[i]->GetIsoCrossSection(part, Z, A, iso, elm, mat);
216  }
217  }
218  G4cout << "G4CrossSectionDataStore::GetCrossSection ERROR: "
219  << " no isotope cross section found"
220  << G4endl;
221  G4cout << " for " << part->GetDefinition()->GetParticleName()
222  << " off Element " << elm->GetName()
223  << " in " << mat->GetName()
224  << " Z= " << Z << " A= " << A
225  << " E(MeV)= " << part->GetKineticEnergy()/MeV << G4endl;
226  throw G4HadronicException(__FILE__, __LINE__,
227  " no applicable data set found for the isotope");
228  return 0.0;
229 }
230 
231 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
232 
233 G4Element*
235  const G4Material* mat,
236  G4Nucleus& target)
237 {
238  G4int nElements = mat->GetNumberOfElements();
239  const G4ElementVector* theElementVector = mat->GetElementVector();
240  G4Element* anElement = (*theElementVector)[0];
241 
242  G4double cross = GetCrossSection(part, mat);
243 
244  // select element from a compound
245  if(1 < nElements) {
246  cross *= G4UniformRand();
247  for(G4int i=0; i<nElements; ++i) {
248  if(cross <= xsecelm[i]) {
249  anElement = (*theElementVector)[i];
250  break;
251  }
252  }
253  }
254 
255  G4int Z = G4lrint(anElement->GetZ());
256  G4Isotope* iso = 0;
257 
258  G4int i = nDataSetList-1;
259  if (dataSetList[i]->IsElementApplicable(part, Z, mat)) {
260 
261  //----------------------------------------------------------------
262  // element-wise cross section
263  // isotope cross section is not computed
264  //----------------------------------------------------------------
265  G4int nIso = anElement->GetNumberOfIsotopes();
266  if (0 >= nIso) {
267  G4cout << " Element " << anElement->GetName() << " Z= " << Z
268  << " has no isotopes " << G4endl;
269  throw G4HadronicException(__FILE__, __LINE__,
270  " Isotope vector is not defined");
271  return anElement;
272  }
273  // isotope abundances
274  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
275  iso = (*isoVector)[0];
276 
277  // more than 1 isotope
278  if(1 < nIso) {
279  iso = dataSetList[i]->SelectIsotope(anElement, part->GetKineticEnergy());
280  }
281 
282  } else {
283 
284  //----------------------------------------------------------------
285  // isotope-wise cross section
286  // isotope cross section is computed
287  //----------------------------------------------------------------
288  G4int nIso = anElement->GetNumberOfIsotopes();
289  cross = 0.0;
290 
291  if (0 >= nIso) {
292  G4cout << " Element " << anElement->GetName() << " Z= " << Z
293  << " has no isotopes " << G4endl;
294  throw G4HadronicException(__FILE__, __LINE__,
295  " Isotope vector is not defined");
296  return anElement;
297  }
298 
299  // user-defined isotope abundances
300  G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
301  iso = (*isoVector)[0];
302 
303  // more than 1 isotope
304  if(1 < nIso) {
305  G4double* abundVector = anElement->GetRelativeAbundanceVector();
306  if(G4int(xseciso.size()) < nIso) { xseciso.resize(nIso); }
307 
308  for (G4int j = 0; j<nIso; ++j) {
309  G4double xsec = 0.0;
310  if(abundVector[j] > 0.0) {
311  iso = (*isoVector)[j];
312  xsec = abundVector[j]*
313  GetIsoCrossSection(part, Z, iso->GetN(), iso, anElement, mat, i);
314  }
315  cross += xsec;
316  xseciso[j] = cross;
317  }
318  cross *= G4UniformRand();
319  for (G4int j = 0; j<nIso; ++j) {
320  if(cross <= xseciso[j]) {
321  iso = (*isoVector)[j];
322  break;
323  }
324  }
325  }
326  }
327  target.SetIsotope(iso);
328  return anElement;
329 }
330 
331 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
332 
333 void
335 {
336  if (nDataSetList == 0)
337  {
338  throw G4HadronicException(__FILE__, __LINE__,
339  "G4CrossSectionDataStore: no data sets registered");
340  return;
341  }
342  for (G4int i=0; i<nDataSetList; ++i) {
343  dataSetList[i]->BuildPhysicsTable(aParticleType);
344  }
345 }
346 
347 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
348 
349 void
351 {
352  // Print out all cross section data sets used and the energies at
353  // which they apply
354 
355  if (nDataSetList == 0) {
356  G4cout << "WARNING - G4CrossSectionDataStore::DumpPhysicsTable: "
357  << " no data sets registered" << G4endl;
358  return;
359  }
360 
361  for (G4int i = nDataSetList-1; i >= 0; --i) {
362  G4double e1 = dataSetList[i]->GetMinKinEnergy();
363  G4double e2 = dataSetList[i]->GetMaxKinEnergy();
364  G4cout
365  << " Cr_sctns: " << std::setw(25) << dataSetList[i]->GetName() << ": "
366  << G4BestUnit(e1, "Energy")
367  << " ---> "
368  << G4BestUnit(e2, "Energy") << "\n";
369  if (dataSetList[i]->GetName() == "G4CrossSectionPairGG") {
370  dataSetList[i]->DumpPhysicsTable(aParticleType);
371  }
372  }
373 }
374 
375 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
376 
378  std::ofstream& outFile)
379 {
380  // Write cross section data set info to html physics list
381  // documentation page
382 
383  G4double ehi = 0;
384  G4double elo = 0;
385  for (G4int i = nDataSetList-1; i > 0; i--) {
386  elo = dataSetList[i]->GetMinKinEnergy()/GeV;
387  ehi = dataSetList[i]->GetMaxKinEnergy()/GeV;
388  outFile << " <li><b><a href=\"" << dataSetList[i]->GetName() << ".html\"> "
389  << dataSetList[i]->GetName() << "</a> from "
390  << elo << " GeV to " << ehi << " GeV </b></li>\n";
391  }
392 
393  G4double defaultHi = dataSetList[0]->GetMaxKinEnergy()/GeV;
394  if (ehi < defaultHi) {
395  outFile << " <li><b><a href=\"" << dataSetList[0]->GetName() << ".html\"> "
396  << dataSetList[0]->GetName() << "</a> from "
397  << ehi << " GeV to " << defaultHi << " GeV </b></li>\n";
398  }
399 }
400 
401 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
G4bool GetNaturalAbundanceFlag() const
Definition: G4Element.hh:263
std::vector< G4Isotope * > G4IsotopeVector
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
std::ofstream outFile
Definition: GammaRayTel.cc:68
const G4String & GetName() const
Definition: G4Material.hh:176
G4double GetZ() const
Definition: G4Element.hh:131
G4ParticleDefinition * GetDefinition() const
void DumpHtml(const G4ParticleDefinition &, std::ofstream &)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const XML_Char * target
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
void BuildPhysicsTable(const G4ParticleDefinition &)
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:162
void DumpPhysicsTable(const G4ParticleDefinition &)
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127