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

#include <G4IonICRU73Data.hh>

Public Member Functions

 G4IonICRU73Data ()
 
 ~G4IonICRU73Data ()
 
void Initialise ()
 
G4double GetDEDX (const G4Material *, const G4int Z, const G4double e, const G4double loge) const
 
G4IonICRU73Dataoperator= (const G4IonICRU73Data &right)=delete
 
 G4IonICRU73Data (const G4IonICRU73Data &)=delete
 

Detailed Description

Definition at line 55 of file G4IonICRU73Data.hh.

Constructor & Destructor Documentation

◆ G4IonICRU73Data() [1/2]

G4IonICRU73Data::G4IonICRU73Data ( )
explicit

Definition at line 98 of file G4IonICRU73Data.cc.

99{
100 fEmin = 0.025*CLHEP::MeV;
101 fEmax = 2.5*CLHEP::MeV;
102 fNbins = fNbinsPerDecade*G4lrint(std::log10(fEmax/fEmin));
103 fVector = new G4PhysicsFreeVector(fSpline);
104 for(G4int i=0; i<81; ++i) {
105 fMatData[i] = new std::vector<G4PhysicsLogVector*>;
106 }
107}
int G4int
Definition: G4Types.hh:85
int G4lrint(double ad)
Definition: templates.hh:134

◆ ~G4IonICRU73Data()

G4IonICRU73Data::~G4IonICRU73Data ( )

Definition at line 111 of file G4IonICRU73Data.cc.

112{
113 delete fVector;
114 for(G4int i=0; i<81; ++i) {
115 auto v = fMatData[i];
116 for(G4int j=0; j<fNmat; ++j) {
117 delete (*v)[j];
118 }
119 delete v;
120 for(G4int j=0; j<93; ++j) { delete fElmData[i][j]; }
121 }
122}

◆ G4IonICRU73Data() [2/2]

G4IonICRU73Data::G4IonICRU73Data ( const G4IonICRU73Data )
delete

Member Function Documentation

◆ GetDEDX()

G4double G4IonICRU73Data::GetDEDX ( const G4Material mat,
const G4int  Z,
const G4double  e,
const G4double  loge 
) const

Definition at line 126 of file G4IonICRU73Data.cc.

128{
129 G4PhysicsLogVector* v = nullptr;
130 G4int Z2 = std::min(Z, 80);
131 if(1 == mat->GetNumberOfElements()) {
132 G4int Z1 = std::min((*(mat->GetElementVector()))[0]->GetZasInt(), 80);
133 v = fElmData[Z2][Z1];
134 } else {
135 G4int idx = fMatIndex[mat->GetIndex()];
136 if(idx < fNmat) {
137 v = (*(fMatData[Z2]))[idx];
138 }
139 }
140 if(nullptr == v) { return 0.0; }
141 G4double res = (e > fEmin) ? v->LogVectorValue(e, loge)
142 : (*v)[0]*std::sqrt(e/fEmin);
143 return res;
144}
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
size_t GetIndex() const
Definition: G4Material.hh:255
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

Referenced by G4LindhardSorensenIonModel::CorrectionsAlongStep().

◆ Initialise()

void G4IonICRU73Data::Initialise ( )

Definition at line 148 of file G4IonICRU73Data.cc.

149{
150 // fill directory path
151 if(fDataDirectory.empty()) {
152 const char* path = G4FindDataDir("G4LEDATA");
153 if (nullptr != path) {
154 std::ostringstream ost;
155 ost << path << "/ion_stopping_data/";
156 fDataDirectory = ost.str();
157 } else {
158 G4Exception("G4IonICRU73Data::Initialise(..)","em013",
160 "Environment variable G4LEDATA is not defined");
161 }
162 }
163
164 std::size_t nmat = G4Material::GetNumberOfMaterials();
165 if(nmat == fMatIndex.size()) { return; }
166
167 if(0 < fVerbose) {
168 G4cout << "### G4IonICRU73Data::Initialise() for " << nmat
169 << " materials" << G4endl;
170 }
171 fMatIndex.resize(nmat, -1);
172 for(G4int j=0; j<81; ++j) {
173 fMatData[j]->resize(nmat, nullptr);
174 }
176 auto mtable = G4Material::GetMaterialTable();
177
178 for(G4int i=0; i<(G4int)nmat; ++i) {
179 fNmat = i;
180 const G4Material* mat = (*mtable)[i];
181 G4int idx = (G4int)mat->GetIndex();
182 if(1 < fVerbose) {
183 G4cout << i << ". material:" << mat->GetName()
184 << " idx=" << idx << G4endl;
185 }
186 if(fMatIndex[idx] == -1) {
187 fMatIndex[idx] = i;
188 G4String matname = mat->GetName();
189 G4bool isOK = false;
190 if(1 == mat->GetNumberOfElements()) {
191 ReadElementData(mat, useICRU90);
192 isOK = true;
193 if(1 < fVerbose) {
194 G4cout << "Material from single element" << G4endl;
195 }
196 }
197 if(!isOK && useICRU90) {
198 for(G4int j=0; j<3; ++j) {
199 if(matname == namesICRU90[j]) {
200 ReadMaterialData(mat, densityCoef[j], true);
201 isOK = true;
202 if(1 < fVerbose) {
203 G4cout << "ICRU90 material" << G4endl;
204 }
205 break;
206 }
207 }
208 }
209 if(!isOK) {
210 for(G4int j=0; j<31; ++j) {
211 if(matname == namesICRU73[j]) {
212 ReadMaterialData(mat, 1.0, false);
213 isOK = true;
214 if(1 < fVerbose) {
215 G4cout << "ICRU73 material" << G4endl;
216 }
217 break;
218 }
219 }
220 }
221 if(!isOK) {
222 ReadElementData(mat, useICRU90);
223 if(1 < fVerbose) {
224 G4cout << "Read element data" << G4endl;
225 }
226 }
227 }
228 }
229}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4EmParameters * Instance()
G4bool UseICRU90Data() const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172

Referenced by G4LindhardSorensenIonModel::Initialise().

◆ operator=()

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

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