Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4IonICRU73Data.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// GEANT4 Class file
29//
30// Description: Data on stopping power
31//
32// Author: Vladimir Ivanchenko
33//
34// Creation date: 23.10.2021
35//
36//----------------------------------------------------------------------------
37//
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40
41#include <fstream>
42#include <sstream>
43
44#include "G4IonICRU73Data.hh"
45#include "G4SystemOfUnits.hh"
46#include "G4PhysicsLogVector.hh"
48#include "G4EmParameters.hh"
51#include "G4Material.hh"
52#include "G4Element.hh"
53
54namespace {
55
56const G4String namesICRU73[31] = {
57"G4_A-150_TISSUE"
58"G4_ADIPOSE_TISSUE_ICRP"
59"G4_AIR"
60"G4_ALUMINUM_OXIDE"
61"G4_BONE_COMPACT_ICRU"
62"G4_BONE_CORTICAL_ICRP"
63"G4_C-552"
64"G4_CALCIUM_FLUORIDE"
65"G4_CARBON_DIOXIDE"
66"G4_KAPTON"
67"G4_LITHIUM_FLUORIDE"
68"G4_LITHIUM_TETRABORATE"
69"G4_LUCITE"
70"G4_METHANE"
71"G4_MUSCLE_STRIATED_ICRU"
72"G4_MYLAR"
73"G4_NYLON-6-6"
74"G4_PHOTO_EMULSION"
75"G4_PLASTIC_SC_VINYLTOLUENE"
76"G4_POLYCARBONATE"
77"G4_POLYETHYLENE"
78"G4_POLYSTYRENE"
79"G4_PROPANE"
80"G4_Pyrex_Glass"
81"G4_SILICON_DIOXIDE"
82"G4_SODIUM_IODIDE"
83"G4_TEFLON"
84"G4_TISSUE-METHANE"
85"G4_TISSUE-PROPANE"
86"G4_WATER"
87"G4_WATER_VAPOR"};
88 const G4String namesICRU90[3] = {"G4_AIR","G4_GRAPHITE","G4_WATER"};
89 const G4double densityCoef[3] = {0.996, 1.025, 0.998};
90 const G4int NZ = 27;
91 const G4int zdat[NZ] = { 5, 6, 7, 8, 13, 14, 15, 16, 22, 26,
92 28, 29, 30, 31, 32, 33, 34, 40, 47, 48,
93 49, 51, 52, 72, 73, 74, 79 };
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
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}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
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}
123
124//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
125
127 const G4double e, const G4double loge) const
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}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147
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}
230
231//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
232
233void G4IonICRU73Data::ReadMaterialData(const G4Material* mat,
234 const G4double coeff,
235 const G4bool useICRU90)
236{
237 G4String name = mat->GetName();
238 for(G4int Z=3; Z<81; ++Z) {
239 std::ostringstream ost;
240 ost << fDataDirectory << "icru";
241 G4int Z1 = Z;
242 G4double scale = 1.0;
243 if(useICRU90 && Z <= 18) {
244 ost << "90";
245 } else {
246 ost << "73";
247 for(G4int i=0; i<NZ; ++i) {
248 if(Z == zdat[i]) {
249 break;
250 } else if(i == NZ-1) {
251 Z1 = zdat[NZ - 1];
252 scale = (G4double)(Z*Z)/(G4double)(Z1*Z1);
253 } else if(Z > zdat[i] && Z < zdat[i+1]) {
254 if(Z - zdat[i] <= zdat[i + 1] - Z) {
255 Z1 = zdat[i];
256 } else {
257 Z1 = zdat[i+1];
258 }
259 scale = (G4double)(Z*Z)/(G4double)(Z1*Z1);
260 break;
261 }
262 }
263 }
264 if(nullptr == (*(fMatData[Z1]))[fNmat]) {
265 ost << "/z" << Z1 << "_" << name << ".dat";
266 G4PhysicsLogVector* v = RetrieveVector(ost, false);
267 if(nullptr != v) {
268 const G4double fact = coeff *
269 mat->GetDensity() * CLHEP::MeV * 1000 * CLHEP::cm2 / CLHEP::g;
270 v->ScaleVector(CLHEP::MeV, fact);
271 if(2 < fVerbose) {
272 G4cout << "### Data for " << name
273 << " and projectile Z=" << Z1 << G4endl;
274 G4cout << *v << G4endl;
275 }
276 (*(fMatData[Z1]))[fNmat] = v;
277 }
278 }
279 if(Z != Z1) {
280 auto v2 = (*(fMatData[Z1]))[fNmat];
281 if(nullptr != v2) {
282 auto v1 = new G4PhysicsLogVector(*v2);
283 (*(fMatData[Z]))[fNmat] = v1;
284 v1->ScaleVector(1.0, scale);
285 }
286 }
287 }
288}
289
290//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
291
292void G4IonICRU73Data::ReadElementData(const G4Material* mat, G4bool useICRU90)
293{
294 const G4ElementVector* elmv = mat->GetElementVector();
295 const G4double* dens = mat->GetFractionVector();
296 const G4int nelm = (G4int)mat->GetNumberOfElements();
297 for(G4int Z=3; Z<81; ++Z) {
298 G4PhysicsLogVector* v = nullptr;
299 if(1 == nelm) {
300 v = FindOrBuildElementData(Z, (*elmv)[0]->GetZasInt(), useICRU90);
301 } else {
302 v = new G4PhysicsLogVector(fEmin, fEmax, fNbins, fSpline);
303 for(G4int i=0; i<=fNbins; ++i) {
304 G4double dedx = 0.0;
305 for(G4int j=0; j<nelm; ++j) {
306 G4PhysicsLogVector* v1 =
307 FindOrBuildElementData(Z, (*elmv)[j]->GetZasInt(), useICRU90);
308 dedx += (*v1)[i]*dens[j];
309 }
310 v->PutValue(i, dedx);
311 }
312 if(fSpline) { v->FillSecondDerivatives(); }
313 }
314 (*(fMatData[Z]))[fNmat] = v;
315 // scale data for correct units
316 if(nullptr != v) {
317 const G4double fact =
318 mat->GetDensity() * CLHEP::MeV * 1000 * CLHEP::cm2 / CLHEP::g;
319 v->ScaleVector(CLHEP::MeV, fact);
320 if(2 < fVerbose) {
321 G4cout << "### Data for "<< mat->GetName()
322 << " for projectile Z=" << Z << G4endl;
323 G4cout << *v << G4endl;
324 }
325 }
326 }
327}
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
330
332G4IonICRU73Data::FindOrBuildElementData(const G4int Z, const G4int Z1,
333 G4bool useICRU90)
334{
335 G4PhysicsLogVector* v = nullptr;
336 if(Z <= 80 && Z1 <= 92) {
337 v = fElmData[Z][Z1];
338 if(nullptr == v) {
339 G4int Z2 = Z1;
340 G4bool isICRU90 = (useICRU90 && Z <= 18 &&
341 (Z1 == 1 || Z1 == 6 || Z1 == 7 || Z1 == 8));
342
343 G4double scale = 1.0;
344 if(!isICRU90) {
345 for(G4int i=0; i<NZ; ++i) {
346 if(Z1 == zdat[i]) {
347 break;
348 } else if(i == NZ-1) {
349 Z2 = zdat[NZ - 1];
350 scale = (G4double)Z1/(G4double)Z2;
351 } else if(Z1 > zdat[i] && Z1 < zdat[i+1]) {
352 if(Z1 - zdat[i] <= zdat[i + 1] - Z1) {
353 Z2 = zdat[i];
354 } else {
355 Z2 = zdat[i+1];
356 }
357 scale = (G4double)Z1/(G4double)Z2;
358 break;
359 }
360 }
361 }
362
363 std::ostringstream ost;
364 ost << fDataDirectory << "icru";
365 if(isICRU90) { ost << "90"; }
366 else { ost << "73"; }
367 ost << "/z" << Z << "_" << Z2 << ".dat";
368 v = RetrieveVector(ost, false);
369 fElmData[Z][Z2] = v;
370 if(Z1 != Z2 && nullptr != v) {
371 fElmData[Z][Z1] = new G4PhysicsLogVector(*v);
372 fElmData[Z][Z1]->ScaleVector(1.0, scale);
373 }
374 }
375 }
376 return v;
377}
378
379//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
380
382G4IonICRU73Data::RetrieveVector(std::ostringstream& ost, G4bool warn)
383{
384 G4PhysicsLogVector* v = nullptr;
385 std::ifstream filein(ost.str().c_str());
386 if (!filein.is_open()) {
387 if(warn) {
389 ed << "Data file <" << ost.str().c_str()
390 << "> is not opened";
391 G4Exception("G4IonICRU73Data::RetrieveVector(..)","em013",
392 FatalException, ed, "Check G4LEDATA");
393 }
394 } else {
395 if(fVerbose > 0) {
396 G4cout << "File " << ost.str()
397 << " is opened by G4IonICRU73Data" << G4endl;
398 }
399
400 // retrieve data from DB
401 if(!fVector->Retrieve(filein, true)) {
403 ed << "Data file <" << ost.str().c_str()
404 << "> is not retrieved!";
405 G4Exception("G4IonICRU73Data::RetrieveVector(..)","had015",
406 FatalException, ed, "Check G4LEDATA");
407 } else {
408 if(fSpline) { fVector->FillSecondDerivatives(); }
409 v = new G4PhysicsLogVector(fEmin, fEmax, fNbins, fSpline);
410 for(G4int i=0; i<=fNbins; ++i) {
411 G4double e = v->Energy(i);
412 G4double dedx = fVector->Value(e);
413 v->PutValue(i, dedx);
414 }
415 if(fSpline) { v->FillSecondDerivatives(); }
416 if(fVerbose > 1) { G4cout << *v << G4endl; }
417 }
418 }
419 return v;
420}
421
422//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4EmParameters * Instance()
G4bool UseICRU90Data() const
G4double GetDEDX(const G4Material *, const G4int Z, const G4double e, const G4double loge) const
G4double GetDensity() const
Definition: G4Material.hh:175
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
const G4double * GetFractionVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
size_t GetIndex() const
Definition: G4Material.hh:255
void PutValue(const std::size_t index, const G4double value)
void ScaleVector(const G4double factorE, const G4double factorV)
G4double Energy(const std::size_t index) const
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
const char * name(G4int ptype)
int G4lrint(double ad)
Definition: templates.hh:134