Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4IonICRU73Data Class Reference

#include <G4IonICRU73Data.hh>

Public Member Functions

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

Private Member Functions

G4PhysicsLogVectorFindOrBuildElementData (const G4int Z, const G4int Z1, G4bool useICRU90)
 
void ReadElementData (const G4Material *mat, G4bool type)
 
void ReadMaterialData (const G4String &name, G4bool type)
 
G4PhysicsLogVectorRetrieveVector (std::ostringstream &in, G4bool warn)
 

Private Attributes

G4String fDataDirectory = ""
 
G4PhysicsLogVectorfElmData [81][81] = {{nullptr}}
 
G4double fEmax
 
G4double fEmin
 
std::vector< G4PhysicsLogVector * > * fMatData [81] = {nullptr}
 
std::vector< G4intfMatIndex
 
G4int fNbins = 0
 
G4int fNbinsPerDecade = 10
 
G4int fNmat = 0
 
G4bool fSpline = false
 
G4PhysicsFreeVectorfVector = nullptr
 
G4int fVerbose = 0
 

Detailed Description

Definition at line 55 of file G4IonICRU73Data.hh.

Constructor & Destructor Documentation

◆ G4IonICRU73Data() [1/2]

G4IonICRU73Data::G4IonICRU73Data ( )
explicit

Definition at line 90 of file G4IonICRU73Data.cc.

91{
92 fEmin = 0.025*CLHEP::MeV;
93 fEmax = 2.5*CLHEP::MeV;
94 fNbins = fNbinsPerDecade*G4lrint(std::log10(fEmax/fEmin) + 0.000001);
96 for(G4int i=0; i<81; ++i) {
97 fMatData[i] = new std::vector<G4PhysicsLogVector*>;
98 }
99}
int G4int
Definition: G4Types.hh:85
G4PhysicsFreeVector * fVector
std::vector< G4PhysicsLogVector * > * fMatData[81]
static constexpr double MeV
int G4lrint(double ad)
Definition: templates.hh:134

References fEmax, fEmin, fMatData, fNbins, fNbinsPerDecade, fSpline, fVector, G4lrint(), and CLHEP::MeV.

◆ ~G4IonICRU73Data()

G4IonICRU73Data::~G4IonICRU73Data ( )

Definition at line 103 of file G4IonICRU73Data.cc.

104{
105 delete fVector;
106 for(G4int i=0; i<81; ++i) {
107 auto v = fMatData[i];
108 for(G4int j=0; j<fNmat; ++j) {
109 delete (*v)[j];
110 }
111 delete v;
112 for(G4int j=0; j<81; ++j) { delete fElmData[i][j]; }
113 }
114}
G4PhysicsLogVector * fElmData[81][81]

References fElmData, fMatData, fNmat, and fVector.

◆ G4IonICRU73Data() [2/2]

G4IonICRU73Data::G4IonICRU73Data ( const G4IonICRU73Data )
delete

Member Function Documentation

◆ FindOrBuildElementData()

G4PhysicsLogVector * G4IonICRU73Data::FindOrBuildElementData ( const G4int  Z,
const G4int  Z1,
G4bool  useICRU90 
)
private

Definition at line 262 of file G4IonICRU73Data.cc.

264{
265 G4PhysicsLogVector* v = nullptr;
266 if(Z <= 80 && Z1 <= 80) {
267 v = fElmData[Z][Z1];
268 if(nullptr == v) {
269 G4bool isICRU90 = (useICRU90 && Z <= 18 &&
270 (Z1 == 1 || Z1 == 6 || Z1 == 7 || Z1 == 8));
271 std::ostringstream ost;
272 ost << fDataDirectory << "icru";
273 if(isICRU90) { ost << "90"; }
274 else { ost << "73"; }
275 ost << "/z" << Z << "_" << Z1 << ".dat";
276 v = RetrieveVector(ost, false);
277 fElmData[Z][Z1] = v;
278 }
279 }
280 return v;
281}
bool G4bool
Definition: G4Types.hh:86
const G4int Z[17]
G4PhysicsLogVector * RetrieveVector(std::ostringstream &in, G4bool warn)
static const G4double Z1[5]
Definition: paraMaker.cc:41

References fDataDirectory, fElmData, RetrieveVector(), Z, and anonymous_namespace{paraMaker.cc}::Z1.

Referenced by ReadElementData().

◆ GetDEDX()

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

Definition at line 118 of file G4IonICRU73Data.cc.

120{
121 G4PhysicsLogVector* v = nullptr;
122 if(1 == mat->GetNumberOfElements()) {
123 G4int Z1 = (*(mat->GetElementVector()))[0]->GetZasInt();
124 if(Z <= 80 && Z1 <= 80) {
125 v = fElmData[Z][Z1];
126 }
127 } else {
128 G4int idx = fMatIndex[mat->GetIndex()];
129 if(idx < fNmat && Z <= 80) {
130 v = (*(fMatData[Z]))[idx];
131 }
132 }
133 G4double res = (nullptr != v) ? v->LogVectorValue(e, loge) : 0.0;
134 return res;
135}
double G4double
Definition: G4Types.hh:83
std::vector< G4int > fMatIndex
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
size_t GetIndex() const
Definition: G4Material.hh:256
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

References fElmData, fMatData, fMatIndex, fNmat, G4Material::GetElementVector(), G4Material::GetIndex(), G4Material::GetNumberOfElements(), G4PhysicsVector::LogVectorValue(), Z, and anonymous_namespace{paraMaker.cc}::Z1.

Referenced by G4LindhardSorensenIonModel::CorrectionsAlongStep().

◆ Initialise()

void G4IonICRU73Data::Initialise ( )

Definition at line 139 of file G4IonICRU73Data.cc.

140{
141 // fill directory path
142 if(fDataDirectory.empty()) {
143 char* path = std::getenv("G4LEDATA");
144 if (nullptr != path) {
145 std::ostringstream ost;
146 ost << path << "/ion_stopping_data/";
147 fDataDirectory = ost.str();
148 } else {
149 G4Exception("G4IonICRU73Data::Initialise(..)","em013",
151 "Environment variable G4LEDATA is not defined");
152 }
153 }
154
156 size_t nmat = G4Material::GetNumberOfMaterials();
157 fMatIndex.resize(nmat, -1);
158
159 const G4ProductionCutsTable* theCoupleTable=
161 size_t numOfCouples = theCoupleTable->GetTableSize();
162
163 for(size_t i=0; i<numOfCouples; ++i) {
164 const G4MaterialCutsCouple* couple =
165 theCoupleTable->GetMaterialCutsCouple(i);
166 const G4Material* mat = couple->GetMaterial();
167 G4int idx = mat->GetIndex();
168 if(fMatIndex[idx] == -1) {
169 G4String matname = mat->GetName();
170 G4bool isOK = false;
171 if(1 == mat->GetNumberOfElements()) {
172 ReadElementData(mat, useICRU90);
173 isOK = true;
174 }
175 if(!isOK && useICRU90) {
176 for(G4int j=0; j<3; ++j) {
177 if(matname == namesICRU90[j]) {
178 ReadMaterialData(matname, true);
179 isOK = true;
180 fMatIndex[idx] = fNmat;
181 ++fNmat;
182 break;
183 }
184 }
185 }
186 if(!isOK) {
187 for(G4int j=0; j<31; ++j) {
188 if(matname == namesICRU73[j]) {
189 ReadMaterialData(matname, false);
190 isOK = true;
191 fMatIndex[idx] = fNmat;
192 ++fNmat;
193 break;
194 }
195 }
196 }
197 if(!isOK) {
198 ReadElementData(mat, useICRU90);
199 fMatIndex[idx] = fNmat;
200 ++fNmat;
201 }
202 }
203 }
204}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
const G4String namesICRU73[31]
const G4String namesICRU90[3]
static G4EmParameters * Instance()
G4bool UseICRU90Data() const
void ReadElementData(const G4Material *mat, G4bool type)
void ReadMaterialData(const G4String &name, G4bool type)
const G4Material * GetMaterial() const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:679
const G4String & GetName() const
Definition: G4Material.hh:173
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()

References FatalException, fDataDirectory, fMatIndex, fNmat, G4Exception(), G4Material::GetIndex(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetName(), G4Material::GetNumberOfElements(), G4Material::GetNumberOfMaterials(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4EmParameters::Instance(), namesICRU73, namesICRU90, ReadElementData(), ReadMaterialData(), and G4EmParameters::UseICRU90Data().

Referenced by G4LindhardSorensenIonModel::Initialise().

◆ operator=()

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

◆ ReadElementData()

void G4IonICRU73Data::ReadElementData ( const G4Material mat,
G4bool  type 
)
private

Definition at line 223 of file G4IonICRU73Data.cc.

224{
225 const G4ElementVector* elmv = mat->GetElementVector();
226 const G4double* dens = mat->GetFractionVector();
227 const G4int nelm = mat->GetNumberOfElements();
228 for(G4int Z=3; Z<81; ++Z) {
229 if(1 == nelm) {
230 FindOrBuildElementData(Z, (*elmv)[0]->GetZasInt(), useICRU90);
231 } else {
232 G4PhysicsLogVector* v2 = nullptr;
233 for(G4int j=0; j<nelm; ++j) {
234 v2 = FindOrBuildElementData(Z, (*elmv)[j]->GetZasInt(), useICRU90);
235 if(nullptr == v2) { break; }
236 }
237 // for one or more elements there is no data
238 if(nullptr == v2) {
239 fMatData[Z]->push_back(v2);
240 continue;
241 }
244 for(G4int i=0; i<=fNbins; ++i) {
245 G4double dedx = 0;
246 for(G4int j=0; j<nelm; ++j) {
247 G4PhysicsLogVector* v1 =
248 FindOrBuildElementData(Z, (*elmv)[j]->GetZasInt(), useICRU90);
249 dedx += (*v1)[i]*dens[j];
250 }
251 v->PutValue(i, dedx);
252 }
253 if(fSpline) { v->FillSecondDerivatives(); }
254 fMatData[Z]->push_back(v);
255 }
256 }
257}
std::vector< const G4Element * > G4ElementVector
G4PhysicsLogVector * FindOrBuildElementData(const G4int Z, const G4int Z1, G4bool useICRU90)
const G4double * GetFractionVector() const
Definition: G4Material.hh:190
void PutValue(const std::size_t index, const G4double value)
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)

References fEmax, fEmin, G4PhysicsVector::FillSecondDerivatives(), FindOrBuildElementData(), fMatData, fNbins, fSpline, G4Material::GetElementVector(), G4Material::GetFractionVector(), G4Material::GetNumberOfElements(), G4PhysicsVector::PutValue(), and Z.

Referenced by Initialise().

◆ ReadMaterialData()

void G4IonICRU73Data::ReadMaterialData ( const G4String name,
G4bool  type 
)
private

Definition at line 208 of file G4IonICRU73Data.cc.

209{
210 for(G4int Z=3; Z<81; ++Z) {
211 std::ostringstream ost;
212 ost << fDataDirectory << "icru";
213 if(isICRU90 && Z <= 18) { ost << "90"; }
214 else { ost << "73"; }
215 ost << "/z" << Z << "_" << name << ".dat";
216 G4PhysicsLogVector* v = RetrieveVector(ost, false);
217 (fMatData[Z])->push_back(v);
218 }
219}
const char * name(G4int ptype)

References fDataDirectory, fMatData, G4InuclParticleNames::name(), RetrieveVector(), and Z.

Referenced by Initialise().

◆ RetrieveVector()

G4PhysicsLogVector * G4IonICRU73Data::RetrieveVector ( std::ostringstream &  in,
G4bool  warn 
)
private

Definition at line 286 of file G4IonICRU73Data.cc.

287{
288 G4PhysicsLogVector* v = nullptr;
289 std::ifstream filein(ost.str().c_str());
290 if (!filein.is_open()) {
291 if(warn) {
293 ed << "Data file <" << ost.str().c_str()
294 << "> is not opened";
295 G4Exception("G4IonICRU73Data::RetrieveVector(..)","em013",
296 FatalException, ed, "Check G4LEDATA");
297 }
298 } else {
299 if(fVerbose > 0) {
300 G4cout << "File " << ost.str()
301 << " is opened by G4IonICRU73Data" << G4endl;
302 }
303 // retrieve data from DB
304 if(!fVector->Retrieve(filein, true)) {
306 ed << "Data file <" << ost.str().c_str()
307 << "> is not retrieved!";
308 G4Exception("G4IonICRU73Data::RetrieveVector(..)","had015",
309 FatalException, ed, "Check G4LEDATA");
310 } else {
313 for(G4int i=0; i<=fNbins; ++i) {
314 G4double e = v->Energy(i);
315 G4double dedx = fVector->Value(e);
316 v->PutValue(i, dedx);
317 }
318 const G4double fact = CLHEP::MeV * CLHEP::cm2 /( 0.001 * CLHEP::g);
319 v->ScaleVector(CLHEP::MeV, fact);
320 if(fSpline) { v->FillSecondDerivatives(); }
321 if(fVerbose > 1) { G4cout << *v << G4endl; }
322 }
323 }
324 return v;
325}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void ScaleVector(const G4double factorE, const G4double factorV)
G4double Energy(const std::size_t index) const
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)
G4double Value(const G4double energy, std::size_t &lastidx) const
static constexpr double g
static constexpr double cm2

References CLHEP::cm2, G4PhysicsVector::Energy(), FatalException, fEmax, fEmin, G4PhysicsVector::FillSecondDerivatives(), fNbins, fSpline, fVector, fVerbose, CLHEP::g, G4cout, G4endl, G4Exception(), CLHEP::MeV, G4PhysicsVector::PutValue(), G4PhysicsVector::Retrieve(), G4PhysicsVector::ScaleVector(), and G4PhysicsVector::Value().

Referenced by FindOrBuildElementData(), and ReadMaterialData().

Field Documentation

◆ fDataDirectory

G4String G4IonICRU73Data::fDataDirectory = ""
private

Definition at line 101 of file G4IonICRU73Data.hh.

Referenced by FindOrBuildElementData(), Initialise(), and ReadMaterialData().

◆ fElmData

G4PhysicsLogVector* G4IonICRU73Data::fElmData[81][81] = {{nullptr}}
private

Definition at line 93 of file G4IonICRU73Data.hh.

Referenced by FindOrBuildElementData(), GetDEDX(), and ~G4IonICRU73Data().

◆ fEmax

G4double G4IonICRU73Data::fEmax
private

Definition at line 89 of file G4IonICRU73Data.hh.

Referenced by G4IonICRU73Data(), ReadElementData(), and RetrieveVector().

◆ fEmin

G4double G4IonICRU73Data::fEmin
private

Definition at line 88 of file G4IonICRU73Data.hh.

Referenced by G4IonICRU73Data(), ReadElementData(), and RetrieveVector().

◆ fMatData

std::vector<G4PhysicsLogVector*>* G4IonICRU73Data::fMatData[81] = {nullptr}
private

◆ fMatIndex

std::vector<G4int> G4IonICRU73Data::fMatIndex
private

Definition at line 91 of file G4IonICRU73Data.hh.

Referenced by GetDEDX(), and Initialise().

◆ fNbins

G4int G4IonICRU73Data::fNbins = 0
private

Definition at line 97 of file G4IonICRU73Data.hh.

Referenced by G4IonICRU73Data(), ReadElementData(), and RetrieveVector().

◆ fNbinsPerDecade

G4int G4IonICRU73Data::fNbinsPerDecade = 10
private

Definition at line 98 of file G4IonICRU73Data.hh.

Referenced by G4IonICRU73Data().

◆ fNmat

G4int G4IonICRU73Data::fNmat = 0
private

Definition at line 96 of file G4IonICRU73Data.hh.

Referenced by GetDEDX(), Initialise(), and ~G4IonICRU73Data().

◆ fSpline

G4bool G4IonICRU73Data::fSpline = false
private

Definition at line 100 of file G4IonICRU73Data.hh.

Referenced by G4IonICRU73Data(), ReadElementData(), and RetrieveVector().

◆ fVector

G4PhysicsFreeVector* G4IonICRU73Data::fVector = nullptr
private

Definition at line 94 of file G4IonICRU73Data.hh.

Referenced by G4IonICRU73Data(), RetrieveVector(), and ~G4IonICRU73Data().

◆ fVerbose

G4int G4IonICRU73Data::fVerbose = 0
private

Definition at line 99 of file G4IonICRU73Data.hh.

Referenced by RetrieveVector().


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