Geant4 10.7.0
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 spline)
56{
57 G4int n = mat->GetNumberOfElements();
58 nElmMinusOne = 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);
64 first->SetSpline(spline);
65 xSections[0] = first;
66 for(G4int i=1; i<n; ++i) {
67 xSections[i] = new G4PhysicsVector(*first);
68 }
69 std::vector<double> temp;
70 temp.resize(n, 0.0);
71 for(G4int j=0; j<=bins; ++j) {
72 G4double cross = 0.0;
73 G4double e = first->Energy(j);
74 dp->SetKineticEnergy(e);
75 for(G4int i=0; i<n; ++i) {
76 cross += xs->GetCrossSection(dp, (*theElementVector)[i], mat);
77 temp[i] = cross;
78 }
79 G4double fact = (cross > 0.0) ? 1.0/cross : 0.0;
80 for(G4int i=0; i<n; ++i) {
81 G4double y = (i<n-1) ? temp[i]*fact : 1.0;
82 xSections[i]->PutValue(j, y);
83 }
84 }
85 }
86}
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
91{
92 if(nElmMinusOne > 0) {
93 for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
94 }
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98
100{}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
105{}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
111 G4int bins, G4double emin, G4double emax,
112 G4bool spline)
113{
115 if(nn > nMaterials) {
116 G4PhysicsLogVector* first = nullptr;
117 G4int sbins = std::max(10, bins/5);
119 for(size_t i=nMaterials; i<nn; ++i) {
120 const G4Material* mat = (*mtable)[i];
121 G4PhysicsVector* v = nullptr;
122 G4HadElementSelector* es = nullptr;
123 // create real vector only for complex materials
124 if(mat->GetNumberOfElements() > 1) {
125 if(!first) {
126 first = new G4PhysicsLogVector(emin, emax, bins);
127 first->SetSpline(spline);
128 v = first;
129 } else {
130 v = new G4PhysicsVector(*first);
131 }
132 for(G4int j=0; j<=bins; ++j) {
133 G4double e = first->Energy(j);
134 dp->SetKineticEnergy(e);
135 G4double cros = xs->ComputeCrossSection(dp, mat);
136 v->PutValue(j, cros);
137 }
138 elmSelectors[i] = new G4HadElementSelector(dp, xs, mat, sbins, emin, emax, spline);
139 }
140 xsData.push_back(v);
141 elmSelectors.push_back(es);
142 }
143 nMaterials = nn;
144 }
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
150{
151 for(size_t i=0; i<nMaterials; ++i) {
152 delete xsData[i];
153 delete elmSelectors[i];
154 }
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
158
160{}
161
162//....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:188
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
G4double Energy(std::size_t index) const
void PutValue(std::size_t index, G4double theValue)
void SetSpline(G4bool)