Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronXSDataTable.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// Description: Data structure for cross sections per materials
32//
33// Author: V.Ivanchenko 31.05.2018
34//
35// Modifications:
36//
37//----------------------------------------------------------------------------
38//
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
43#include "G4PhysicsLogVector.hh"
44#include "G4Material.hh"
45#include "G4MaterialTable.hh"
46#include "G4DynamicParticle.hh"
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50
53 const G4Material* mat,
54 G4int bins, G4double emin,
55 G4double emax, G4bool)
56{
57 std::size_t n = mat->GetNumberOfElements();
58 nElmMinusOne = G4int(n - 1);
59 theElementVector = mat->GetElementVector();
60 if(nElmMinusOne > 0) {
61 G4PhysicsVector* first = nullptr;
62 xSections.resize(n, first);
63 first = new G4PhysicsLogVector(emin,emax,bins,false);
64 xSections[0] = first;
65 for(std::size_t i=1; i<n; ++i) {
66 xSections[i] = new G4PhysicsVector(*first);
67 }
68 std::vector<G4double> temp;
69 temp.resize(n, 0.0);
70 for(G4int j=0; j<=bins; ++j) {
71 G4double cross = 0.0;
72 G4double e = first->Energy(j);
73 dp->SetKineticEnergy(e);
74 for(std::size_t i=0; i<n; ++i) {
75 cross += xs->GetCrossSection(dp, (*theElementVector)[i], mat);
76 temp[i] = cross;
77 }
78 G4double fact = (cross > 0.0) ? 1.0/cross : 0.0;
79 for(std::size_t i=0; i<n; ++i) {
80 G4double y = (i<n-1) ? temp[i]*fact : 1.0;
81 xSections[i]->PutValue(j, y);
82 }
83 }
84 }
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90{
91 if(nElmMinusOne > 0) {
92 for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
93 }
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
99{}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104{}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
110 G4int bins, G4double emin, G4double emax,
111 G4bool spline)
112{
113 std::size_t nn = G4Material::GetNumberOfMaterials();
114 if(nn > nMaterials) {
115 if(0 == nMaterials) {
116 xsData.reserve(nn);
117 elmSelectors.reserve(nn);
118 }
119 G4PhysicsLogVector* first = nullptr;
120 G4int sbins = std::max(10, bins/5);
122 for(std::size_t i=nMaterials; i<nn; ++i) {
123 const G4Material* mat = (*mtable)[i];
124 G4PhysicsVector* v = nullptr;
125 G4HadElementSelector* es = nullptr;
126 // create real vector only for complex materials
127 if(mat->GetNumberOfElements() > 1) {
128 if(nullptr == first) {
129 first = new G4PhysicsLogVector(emin, emax, bins, spline);
130 v = first;
131 } else {
132 v = new G4PhysicsVector(*first);
133 }
134 for(G4int j=0; j<=bins; ++j) {
135 G4double e = first->Energy(j);
136 dp->SetKineticEnergy(e);
137 G4double cros = xs->ComputeCrossSection(dp, mat);
138 v->PutValue(j, cros);
139 }
140 if(spline) v->FillSecondDerivatives();
141 elmSelectors[i] = new G4HadElementSelector(dp, xs, mat, sbins, emin, emax, spline);
142 }
143 xsData.push_back(v);
144 elmSelectors.push_back(es);
145 }
146 nMaterials = nn;
147 }
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151
153{
154 for(std::size_t i=0; i<nMaterials; ++i) {
155 delete xsData[i];
156 delete elmSelectors[i];
157 }
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161
163{}
164
165//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void SetKineticEnergy(G4double aEnergy)
G4HadElementSelector(G4DynamicParticle *, G4CrossSectionDataStore *, const G4Material *, G4int bins, G4double emin, G4double emax, G4bool spline)
void Initialise(G4DynamicParticle *, G4CrossSectionDataStore *, G4int bins, G4double emin, G4double emax, G4bool spline)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)