Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eIonisationCrossSectionHandler.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: G4eIonisationCrossSectionHandler.cc 66241 2012-12-13 18:34:42Z gunter $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eIonisationCrossSectionHandler
34 //
35 // Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
36 //
37 // Creation date: 25 Sept 2001
38 //
39 // Modifications:
40 // 10 Oct 2001 M.G. Pia Revision to improve code quality and consistency with design
41 // 19 Jul 2002 VI Create composite data set for material
42 // 21 Jan 2003 V.Ivanchenko Cut per region
43 // 28 Jan 2009 L.Pandola Added public method to make a easier migration of
44 // G4LowEnergyIonisation to G4LivermoreIonisationModel
45 // 15 Jul 2009 Nicolas A. Karakatsanis
46 //
47 // - BuildCrossSectionForMaterials method was revised in order to calculate the
48 // logarithmic values of the loaded data.
49 // It retrieves the data values from the G4EMLOW data files but, then, calculates the
50 // respective log values and loads them to seperate data structures.
51 // The EM data sets, initialized this way, contain both non-log and log values.
52 // These initialized data sets can enhance the computing performance of data interpolation
53 // operations
54 //
55 //
56 //
57 // -------------------------------------------------------------------
58 
60 #include "G4SystemOfUnits.hh"
61 #include "G4VEnergySpectrum.hh"
62 #include "G4DataVector.hh"
63 #include "G4CompositeEMDataSet.hh"
64 #include "G4VDataSetAlgorithm.hh"
67 #include "G4VEMDataSet.hh"
68 #include "G4EMDataSet.hh"
69 #include "G4Material.hh"
70 #include "G4ProductionCutsTable.hh"
71 
72 
74  const G4VEnergySpectrum* spec, G4VDataSetAlgorithm* alg,
75  G4double emin, G4double emax, G4int nbin)
77  theParam(spec),verbose(0)
78 {
79  G4VCrossSectionHandler::Initialise(alg, emin, emax, nbin);
80  interp = new G4LinLogLogInterpolation();
81 }
82 
83 
85 {
86  delete interp;
87 }
88 
89 
91  const G4DataVector& energyVector,
92  const G4DataVector* energyCuts)
93 {
94  std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>;
95 
96  G4DataVector* energies;
97  G4DataVector* cs;
98 
99  G4DataVector* log_energies;
100  G4DataVector* log_cs;
101 
102  G4int nOfBins = energyVector.size();
103 
104  const G4ProductionCutsTable* theCoupleTable=
106  size_t numOfCouples = theCoupleTable->GetTableSize();
107 
108  for (size_t mLocal=0; mLocal<numOfCouples; mLocal++) {
109 
110  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(mLocal);
111  const G4Material* material= couple->GetMaterial();
112  const G4ElementVector* elementVector = material->GetElementVector();
113  const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
114  G4int nElements = material->GetNumberOfElements();
115 
116  if(verbose > 0)
117  {
118  G4cout << "eIonisation CS for " << mLocal << "th material "
119  << material->GetName()
120  << " eEl= " << nElements << G4endl;
121  }
122 
123  G4double tcut = (*energyCuts)[mLocal];
124 
125  G4VDataSetAlgorithm* algo = interp->Clone();
126  G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
127 
128  for (G4int i=0; i<nElements; i++) {
129 
130  G4int Z = (G4int) (*elementVector)[i]->GetZ();
131  G4int nShells = NumberOfComponents(Z);
132 
133  energies = new G4DataVector;
134  cs = new G4DataVector;
135 
136  log_energies = new G4DataVector;
137  log_cs = new G4DataVector;
138 
139  G4double density = nAtomsPerVolume[i];
140 
141  for (G4int bin=0; bin<nOfBins; bin++) {
142 
143  G4double e = energyVector[bin];
144  energies->push_back(e);
145  log_energies->push_back(std::log10(e));
146  G4double value = 0.0;
147  G4double log_value = -300;
148 
149  if(e > tcut) {
150  for (G4int n=0; n<nShells; n++) {
151  G4double cross = FindValue(Z, e, n);
152  G4double p = theParam->Probability(Z, tcut, e, e, n);
153  value += cross * p * density;
154 
155  if(verbose>0 && mLocal == 0 && e>=1. && e<=0.)
156  {
157  G4cout << "G4eIonCrossSH: e(MeV)= " << e/MeV
158  << " n= " << n
159  << " cross= " << cross
160  << " p= " << p
161  << " value= " << value
162  << " tcut(MeV)= " << tcut/MeV
163  << " rho= " << density
164  << " Z= " << Z
165  << G4endl;
166  }
167 
168  }
169  if (value == 0.) value = 1e-300;
170  log_value = std::log10(value);
171  }
172  cs->push_back(value);
173  log_cs->push_back(log_value);
174  }
175  G4VDataSetAlgorithm* algoLocal = interp->Clone();
176 
177  //G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algoLocal,1.,1.);
178 
179  G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,log_energies,log_cs,algoLocal,1.,1.);
180 
181  setForMat->AddComponent(elSet);
182  }
183  set->push_back(setForMat);
184  }
185 
186  return set;
187 }
188 
190  G4double cutEnergy,
191  G4int Z)
192 {
193  G4int nShells = NumberOfComponents(Z);
194  G4double value = 0.;
195  if(energy > cutEnergy)
196  {
197  for (G4int n=0; n<nShells; n++) {
198  G4double cross = FindValue(Z, energy, n);
199  G4double p = theParam->Probability(Z, cutEnergy, energy, energy, n);
200  value += cross * p;
201  }
202  }
203  return value;
204 }
std::vector< G4Element * > G4ElementVector
tuple bin
Definition: plottest35.py:22
G4double GetCrossSectionAboveThresholdForElement(G4double energy, G4double cutEnergy, G4int Z)
const char * p
Definition: xmltok.h:285
const G4String & GetName() const
Definition: G4Material.hh:178
std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
virtual G4VDataSetAlgorithm * Clone() const =0
string material
Definition: eplot.py:19
G4int NumberOfComponents(G4int Z) const
G4double FindValue(G4int Z, G4double e) const
G4GLOB_DLL std::ostream G4cout
const XML_Char int const XML_Char * value
Definition: expat.h:331
virtual void AddComponent(G4VEMDataSet *dataSet)=0
const G4int n
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:216
virtual G4double Probability(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
static const G4double emax
G4eIonisationCrossSectionHandler(const G4VEnergySpectrum *spec, G4VDataSetAlgorithm *alg, G4double emin, G4double emax, G4int nbin)
static G4ProductionCutsTable * GetProductionCutsTable()
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4double energy(const ThreeVector &p, const G4double m)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const