Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CrossSectionHandler.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// Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
29//
30// History:
31// -----------
32// 1 Aug 2001 MGP Created
33// 19 Jul 2002 VI Create composite data set for material
34// 24 Apr 2003 VI Cut per region mfpt
35//
36// 15 Jul 2009 Nicolas A. Karakatsanis
37//
38// - BuildCrossSectionForMaterials method was revised in order to calculate the
39// logarithmic values of the loaded data.
40// It retrieves the data values from the G4EMLOW data files but, then, calculates the
41// respective log values and loads them to seperate data structures.
42// The EM data sets, initialized this way, contain both non-log and log values.
43// These initialized data sets can enhance the computing performance of data interpolation
44// operations
45//
46// -------------------------------------------------------------------
47
50#include "G4VEMDataSet.hh"
51#include "G4EMDataSet.hh"
53#include "G4ShellEMDataSet.hh"
55#include "G4Material.hh"
56#include "G4Element.hh"
57#include "Randomize.hh"
58#include <map>
59#include <vector>
60
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
66{ }
67
68//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69
71{ }
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
75std::vector<G4VEMDataSet*>*
77 const G4DataVector*)
78{
79 G4DataVector* energies;
80 G4DataVector* data;
81
82 G4DataVector* log_energies;
83 G4DataVector* log_data;
84
85 std::vector<G4VEMDataSet*>* matCrossSections = new std::vector<G4VEMDataSet*>;
86
87 const G4ProductionCutsTable* theCoupleTable=
89 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
90
91 std::size_t nOfBins = energyVector.size();
92 const G4VDataSetAlgorithm* interpolationAlgo = CreateInterpolation();
93
94 for (G4int mLocal=0; mLocal<numOfCouples; ++mLocal)
95 {
96 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(mLocal);
97 const G4Material* material= couple->GetMaterial();
98 G4int nElements = (G4int)material->GetNumberOfElements();
99 const G4ElementVector* elementVector = material->GetElementVector();
100 const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector();
101
102 G4VDataSetAlgorithm* algo = interpolationAlgo->Clone();
103
104 G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.);
105
106 for (G4int i=0; i<nElements; ++i) {
107
108 G4int Z = (G4int) (*elementVector)[i]->GetZ();
109 G4double density = nAtomsPerVolume[i];
110
111 energies = new G4DataVector;
112 data = new G4DataVector;
113
114 log_energies = new G4DataVector;
115 log_data = new G4DataVector;
116
117 for (std::size_t bin=0; bin<nOfBins; ++bin)
118 {
119 G4double e = energyVector[bin];
120 energies->push_back(e);
121 if (e==0.) e=1e-300;
122 log_energies->push_back(std::log10(e));
123 G4double cross = density*FindValue(Z,e);
124 data->push_back(cross);
125 if (cross==0.) cross=1e-300;
126 log_data->push_back(std::log10(cross));
127 }
128
129 G4VDataSetAlgorithm* algo1 = interpolationAlgo->Clone();
130 G4VEMDataSet* elSet = new G4EMDataSet(i,energies,data,log_energies,log_data,algo1,1.,1.);
131
132 setForMat->AddComponent(elSet);
133 }
134
135 matCrossSections->push_back(setForMat);
136 }
137 delete interpolationAlgo;
138 return matCrossSections;
139}
140
std::vector< const G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
std::vector< G4VEMDataSet * > * BuildCrossSectionsForMaterials(const G4DataVector &energyVector, const G4DataVector *energyCuts=0) override
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:211
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double FindValue(G4int Z, G4double e) const
virtual G4VDataSetAlgorithm * CreateInterpolation()
virtual G4VDataSetAlgorithm * Clone() const =0
virtual void AddComponent(G4VEMDataSet *dataSet)=0