Geant4 10.7.0
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 ()
 
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 56 of file G4ICRU90StoppingData.hh.

Constructor & Destructor Documentation

◆ G4ICRU90StoppingData()

G4ICRU90StoppingData::G4ICRU90StoppingData ( )
explicit

Definition at line 53 of file G4ICRU90StoppingData.cc.

53 : isInitialized(false)
54{
55 // 1st initialisation
56 for(size_t i=0; i<nvectors; ++i) {
57 materials[i] = nullptr;
58 sdata_proton[i] = nullptr;
59 sdata_alpha[i] = nullptr;
60 }
61 FillData();
62
63 Initialise();
64}

◆ ~G4ICRU90StoppingData()

G4ICRU90StoppingData::~G4ICRU90StoppingData ( )

Definition at line 68 of file G4ICRU90StoppingData.cc.

69{
70 for(size_t i=0; i<nvectors; ++i) {
71 delete sdata_proton[i];
72 delete sdata_alpha[i];
73 }
74}

Member Function Documentation

◆ GetElectronicDEDXforAlpha() [1/2]

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

Definition at line 125 of file G4ICRU90StoppingData.cc.

127{
128 G4int idx = GetIndex(mat);
129 return (idx < 0) ? 0.0 : GetDEDX(sdata_alpha[idx], scaledKinEnergy);
130}
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 151 of file G4ICRU90StoppingData.hh.

153{
154 return (idx < 0 || idx >= nvectors) ? 0.0
155 : GetDEDX(sdata_alpha[idx], scaledKinEnergy);
156}

◆ GetElectronicDEDXforProton() [1/2]

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

Definition at line 116 of file G4ICRU90StoppingData.cc.

118{
119 G4int idx = GetIndex(mat);
120 return (idx < 0) ? 0.0 : GetDEDX(sdata_proton[idx], kinEnergy);
121}

Referenced by G4BetheBlochModel::ComputeDEDXPerVolume().

◆ GetElectronicDEDXforProton() [2/2]

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

Definition at line 142 of file G4ICRU90StoppingData.hh.

144{
145 return (idx < 0 || idx >= nvectors) ? 0.0
146 : GetDEDX(sdata_proton[idx], kinEnergy);
147}

◆ GetIndex() [1/2]

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

Definition at line 111 of file G4ICRU90StoppingData.hh.

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

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

◆ GetIndex() [2/2]

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

Definition at line 122 of file G4ICRU90StoppingData.hh.

123{
124 G4int idx = -1;
125 if (nam == materials[0]->GetName()) { idx = 0; }
126 else if(nam == materials[1]->GetName()) { idx = 1; }
127 else if(nam == materials[2]->GetName()) { idx = 2; }
128 return idx;
129}

◆ Initialise()

void G4ICRU90StoppingData::Initialise ( )

Definition at line 78 of file G4ICRU90StoppingData.cc.

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

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

◆ IsApplicable()

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

Definition at line 106 of file G4ICRU90StoppingData.hh.

107{
108 return (mat == materials[0] || mat == materials[1] || mat == materials[2]);
109}

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