Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ICRU90StoppingData Class Reference

#include <G4ICRU90StoppingData.hh>

Public Member Functions

 G4ICRU90StoppingData ()
 
 ~G4ICRU90StoppingData ()
 
G4ICRU90StoppingDataoperator= (const G4ICRU90StoppingData &right)=delete
 
 G4ICRU90StoppingData (const G4ICRU90StoppingData &)=delete
 
void Initialise ()
 
G4double GetElectronicDEDXforProton (const G4Material *, G4double kinEnergy) const
 
G4double GetElectronicDEDXforAlpha (const G4Material *, G4double scaledKinEnergy) const
 
G4int GetIndex (const G4Material *) const
 
G4int GetIndex (const G4String &) const
 
G4double GetElectronicDEDXforProton (G4int idx, G4double kinEnergy) const
 
G4double GetElectronicDEDXforAlpha (G4int idx, G4double scaledKinEnergy) const
 
G4bool IsApplicable (const G4Material *) const
 

Detailed Description

Definition at line 54 of file G4ICRU90StoppingData.hh.

Constructor & Destructor Documentation

◆ G4ICRU90StoppingData() [1/2]

G4ICRU90StoppingData::G4ICRU90StoppingData ( )

Definition at line 49 of file G4ICRU90StoppingData.cc.

50{
51 // 1st initialisation
52 for (std::size_t i = 0; i < nvectors; ++i) {
53 materials[i] = nullptr;
54 sdata_proton[i] = nullptr;
55 sdata_alpha[i] = nullptr;
56 }
57 FillData();
58
59 Initialise();
60}

◆ ~G4ICRU90StoppingData()

G4ICRU90StoppingData::~G4ICRU90StoppingData ( )

Definition at line 64 of file G4ICRU90StoppingData.cc.

65{
66 for (std::size_t i = 0; i < nvectors; ++i) {
67 delete sdata_proton[i];
68 delete sdata_alpha[i];
69 }
70}

◆ G4ICRU90StoppingData() [2/2]

G4ICRU90StoppingData::G4ICRU90StoppingData ( const G4ICRU90StoppingData & )
delete

Member Function Documentation

◆ GetElectronicDEDXforAlpha() [1/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforAlpha ( const G4Material * mat,
G4double scaledKinEnergy ) const

Definition at line 127 of file G4ICRU90StoppingData.cc.

129{
130 const G4int idx = GetIndex(mat);
131 return (idx >= 0) ? GetDEDX(sdata_alpha[idx], scaledKinEnergy) : 0.0;
132}
int G4int
Definition G4Types.hh:85
G4int GetIndex(const G4Material *) const

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume().

◆ GetElectronicDEDXforAlpha() [2/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforAlpha ( G4int idx,
G4double scaledKinEnergy ) const
inline

Definition at line 157 of file G4ICRU90StoppingData.hh.

159{
160 return (idx >= 0 && idx < nvectors) ? GetDEDX(sdata_alpha[idx], scaledKinEnergy) : 0.0;
161}

◆ GetElectronicDEDXforProton() [1/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforProton ( const G4Material * mat,
G4double kinEnergy ) const

Definition at line 118 of file G4ICRU90StoppingData.cc.

120{
121 const G4int idx = GetIndex(mat);
122 return (idx >= 0) ? GetDEDX(sdata_proton[idx], kinEnergy) : 0.0;
123}

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume().

◆ GetElectronicDEDXforProton() [2/2]

G4double G4ICRU90StoppingData::GetElectronicDEDXforProton ( G4int idx,
G4double kinEnergy ) const
inline

Definition at line 149 of file G4ICRU90StoppingData.hh.

151{
152 return (idx >= 0 && idx < nvectors) ? GetDEDX(sdata_proton[idx], kinEnergy) : 0.0;
153}

◆ GetIndex() [1/2]

G4int G4ICRU90StoppingData::GetIndex ( const G4Material * mat) const
inline

Definition at line 106 of file G4ICRU90StoppingData.hh.

107{
108 G4int idx = -1;
109 if (mat == materials[1]) {
110 idx = 1;
111 }
112 else if (mat == materials[0]) {
113 idx = 0;
114 }
115 else if (mat == materials[2]) {
116 idx = 2;
117 }
118 return idx;
119}

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume(), GetElectronicDEDXforAlpha(), and GetElectronicDEDXforProton().

◆ GetIndex() [2/2]

G4int G4ICRU90StoppingData::GetIndex ( const G4String & nam) const
inline

Definition at line 123 of file G4ICRU90StoppingData.hh.

124{
125 G4int idx = -1;
126 if (nam == materials[1]->GetName()) {
127 idx = 1;
128 }
129 else if (nam == materials[0]->GetName()) {
130 idx = 0;
131 }
132 else if (nam == materials[2]->GetName()) {
133 idx = 2;
134 }
135 return idx;
136}

◆ Initialise()

void G4ICRU90StoppingData::Initialise ( )

Definition at line 74 of file G4ICRU90StoppingData.cc.

75{
76 if (isInitialized) {
77 return;
78 }
79 // this method may be called several times during initialisation
81 if (nmat == (G4int)nvectors) {
82 return;
83 }
84
85 static const G4String nameNIST_ICRU90[3] = {"G4_AIR", "G4_WATER", "G4_GRAPHITE"};
86
87 // loop via material list to add extra data
88 for (G4int i = 0; i < nmat; ++i) {
89 const G4Material* mat = (*(G4Material::GetMaterialTable()))[i];
90
91 G4bool isThere = false;
92 for (auto const & material : materials) {
93 if (mat == material) {
94 isThere = true;
95 break;
96 }
97 }
98 if (! isThere) {
99 // check list of NIST materials
100 G4String mname = mat->GetName();
101 for (G4int j = 0; j < nvectors; ++j) {
102 if (mname == nameNIST_ICRU90[j]) {
103 materials[j] = mat;
104 break;
105 }
106 }
107 }
108 isInitialized =
109 ((materials[0] != nullptr) && (materials[1] != nullptr) && (materials[2] != nullptr));
110 if (isInitialized) {
111 return;
112 }
113 }
114}
bool G4bool
Definition G4Types.hh:86
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const

Referenced by G4ICRU90StoppingData(), G4BetheBlochModel::Initialise(), and G4BraggModel::Initialise().

◆ IsApplicable()

G4bool G4ICRU90StoppingData::IsApplicable ( const G4Material * mat) const
inline

Definition at line 99 of file G4ICRU90StoppingData.hh.

100{
101 return (mat == materials[1] || mat == materials[0] || mat == materials[2]);
102}

◆ operator=()

G4ICRU90StoppingData & G4ICRU90StoppingData::operator= ( const G4ICRU90StoppingData & right)
delete

The documentation for this class was generated from the following files: