Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronXSDataTable.hh
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#ifndef HadronXSDataTable_h
41#define HadronXSDataTable_h 1
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
45#include "globals.hh"
46#include "G4Element.hh"
47#include "G4ElementVector.hh"
48#include "G4PhysicsVector.hh"
49#include "Randomize.hh"
50#include <vector>
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
54class G4Material;
57
59{
60public:
61
63 const G4Material*, G4int bins,
64 G4double emin, G4double emax, G4bool spline);
65
67
68 void Dump();
69
70 inline const G4Element* SelectRandomAtom(G4double e) const
71 {
72 const G4Element* element = (*theElementVector)[nElmMinusOne];
73 if (nElmMinusOne > 0) {
75 for(G4int i=0; i<nElmMinusOne; ++i) {
76 if (x <= xSections[i]->Value(e)) {
77 element = (*theElementVector)[i];
78 break;
79 }
80 }
81 }
82 return element;
83 }
84
85private:
86
88 G4HadElementSelector& operator=(const G4HadElementSelector &right) = delete;
89
90 G4int nElmMinusOne;
91 const G4ElementVector* theElementVector;
92 std::vector<G4PhysicsVector*> xSections;
93};
94
96{
97
98public:
99
100 explicit G4HadronXSDataTable();
101
103
105 G4int bins, G4double emin, G4double emax,
106 G4bool spline);
107
108 inline const G4PhysicsVector* HasData(size_t idx) const
109 {
110 return xsData[idx];
111 };
112
113 inline G4double GetCrossSection(G4double e, size_t idx) const
114 {
115 return xsData[idx]->Value(e);
116 };
117
118 inline const G4Element* SelectRandomAtom(G4double e, size_t idx) const
119 {
120 return elmSelectors[idx]->SelectRandomAtom(e);
121 };
122
123 void Dump();
124
125private:
126
127 // Assignment operator and copy constructor
128 G4HadronXSDataTable & operator=
129 (const G4HadronXSDataTable &right) = delete;
131
132 std::vector<G4PhysicsVector*> xsData;
133 std::vector<G4HadElementSelector*> elmSelectors;
134
135 size_t nMaterials;
136};
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
140#endif
141
std::vector< const G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Element * SelectRandomAtom(G4double e) const
void Initialise(G4DynamicParticle *, G4CrossSectionDataStore *, G4int bins, G4double emin, G4double emax, G4bool spline)
G4double GetCrossSection(G4double e, size_t idx) const
const G4PhysicsVector * HasData(size_t idx) const
const G4Element * SelectRandomAtom(G4double e, size_t idx) const