Geant4 11.2.2
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=3; i<=ZPROJMAX; ++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=3; i<=ZPROJMAX; ++i) {
115 auto v = fMatData[i];
116 if(nullptr != v) {
117 for(auto & dat : *v) {
118 delete dat;
119 }
120 delete v;
121 }
122 for(G4int j=1; j<=ZTARGMAX; ++j) {
123 delete fElmData[i][j];
124 }
125 }
126}

◆ 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 130 of file G4IonICRU73Data.cc.

132{
133 // if data does not exist the result is zero
134 if (Z > ZPROJMAX) { return 0.0; }
135
136 G4PhysicsLogVector* v = nullptr;
137 if (1 == mat->GetNumberOfElements()) {
138 G4int Z1 = (*(mat->GetElementVector()))[0]->GetZasInt();
139 if(Z1 > ZTARGMAX) { return 0.0; }
140 v = fElmData[Z][Z1];
141 }
142 else {
143 G4int idx = fMatIndex[mat->GetIndex()];
144 if (idx < 0) { return 0.0; }
145 v = (*(fMatData[Z]))[idx];
146 }
147 if(nullptr == v) { return 0.0; }
148 G4double res = (e > fEmin) ? v->LogVectorValue(e, loge)
149 : (*v)[0]*std::sqrt(e/fEmin);
150 return res;
151}
double G4double
Definition G4Types.hh:83
const G4ElementVector * GetElementVector() const
std::size_t GetIndex() const
std::size_t GetNumberOfElements() const
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

Referenced by G4LindhardSorensenIonModel::CorrectionsAlongStep().

◆ Initialise()

void G4IonICRU73Data::Initialise ( )

Definition at line 155 of file G4IonICRU73Data.cc.

156{
157 // fill directory path
158 if(fDataDirectory.empty()) {
159 std::ostringstream ost;
160 ost << G4EmParameters::Instance()->GetDirLEDATA() << "/ion_stopping_data/";
161 fDataDirectory = ost.str();
162 }
163
164 const std::size_t nmat = G4Material::GetNumberOfMaterials();
165
166 // do nothing if for a new run list of materials is not changed
167 if(nmat == fMatIndex.size()) { return; }
168
169 // look over all materials and recompute all materials
170 // do not recompute element data if exist
171 if(1 < fVerbose) {
172 G4cout << "### G4IonICRU73Data::Initialise() for " << nmat
173 << " materials" << G4endl;
174 }
175 fMatIndex.resize(nmat, -1);
176 for(G4int j=3; j<=ZPROJMAX; ++j) {
177 fMatData[j]->resize(nmat, nullptr);
178 }
180 auto mtable = G4Material::GetMaterialTable();
181
182 for(G4int i=0; i<(G4int)nmat; ++i) {
183 const G4Material* mat = (*mtable)[i];
184 const G4String matname = mat->GetName();
185 G4int idx = (G4int)mat->GetIndex();
186 if(1 < fVerbose) {
187 G4cout << i << ". material: " << matname
188 << " idx=" << idx << " matIdx=" << fMatIndex[idx]
189 << G4endl;
190 }
191 if(fMatIndex[idx] == -1) {
192 G4bool isOK = false;
193
194 // for material in ICRU90
195 if(useICRU90) {
196 for(G4int j=0; j<3; ++j) {
197 if(matname == namesICRU90[j]) {
198 ReadMaterialData(mat, densityCoef[j], true);
199 isOK = true;
200 if(1 < fVerbose) {
201 G4cout << "ICRU90 material " << matname << G4endl;
202 }
203 break;
204 }
205 }
206 }
207 // for material in ICRU73
208 if(!isOK) {
209 for(G4int j=0; j<31; ++j) {
210 if(matname == namesICRU73[j]) {
211 ReadMaterialData(mat, 1.0, false);
212 isOK = true;
213 if(1 < fVerbose) {
214 G4cout << "ICRU73 material " << matname << G4endl;
215 }
216 break;
217 }
218 }
219 }
220 // dEdx is combined from element stopping power
221 // Target element Z <= ZTARGMAX
222 if(!isOK) {
223 const std::size_t nelm = mat->GetNumberOfElements();
224 auto elmv = mat->GetElementVector();
225 G4bool isOut = false;
226 for (std::size_t j=0; j<nelm; ++j) {
227 if ((*elmv)[j]->GetZasInt() > ZTARGMAX) {
228 isOut = true;
229 break;
230 }
231 }
232 if (!isOut) {
233 ReadElementData(mat, useICRU90);
234 isOK = true;
235 if(1 < fVerbose) {
236 G4cout << "Data via elements for " << matname << G4endl;
237 }
238 }
239 }
240 if(isOK) { fMatIndex[idx] = i; }
241 }
242 if(1 < fVerbose) {
243 G4cout << " matData: " << fMatData[i] << G4endl;
244 }
245 }
246}
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4EmParameters * Instance()
G4bool UseICRU90Data() const
const G4String & GetDirLEDATA() const
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const

Referenced by G4LindhardSorensenIonModel::Initialise().

◆ operator=()

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

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