Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LossTableBuilder.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//
31// File name: G4LossTableBuilder
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 03.01.2002
36//
37// Modifications:
38//
39// 23-01-03 V.Ivanchenko Cut per region
40// 21-07-04 V.Ivanchenko Fix problem of range for dedx=0
41// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
42// 07-12-04 Fix of BuildDEDX table (V.Ivanchenko)
43// 27-03-06 Add bool options isIonisation (V.Ivanchenko)
44// 16-01-07 Fill new (not old) DEDX table (V.Ivanchenko)
45// 12-02-07 Use G4LPhysicsFreeVector for the inverse range table (V.Ivanchenko)
46// 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
47//
48// Class Description:
49//
50// -------------------------------------------------------------------
51//
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
55#include "G4LossTableBuilder.hh"
56#include "G4SystemOfUnits.hh"
57#include "G4PhysicsTable.hh"
58#include "G4PhysicsLogVector.hh"
63#include "G4Material.hh"
64#include "G4VEmModel.hh"
66#include "G4LossTableManager.hh"
67#include "G4EmParameters.hh"
68
69std::vector<G4double>* G4LossTableBuilder::theDensityFactor = nullptr;
70std::vector<G4int>* G4LossTableBuilder::theDensityIdx = nullptr;
71std::vector<G4bool>* G4LossTableBuilder::theFlag = nullptr;
72
73//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74
76{
77 theParameters = G4EmParameters::Instance();
78 if (nullptr == theFlag) {
79 if (!master) {
81 ed << "The table builder is instantiated in a worker thread ";
82 G4Exception("G4LossTableBuilder::G4LossTableBuilder ", "em0001",
83 JustWarning, ed);
84 }
85 theDensityFactor = new std::vector<G4double>;
86 theDensityIdx = new std::vector<G4int>;
87 theFlag = new std::vector<G4bool>;
88 isInitializer = true;
89 }
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93
95{
96 if (isInitializer) {
97 delete theDensityFactor;
98 delete theDensityIdx;
99 delete theFlag;
100 theDensityFactor = nullptr;
101 theDensityIdx = nullptr;
102 theFlag = nullptr;
103 }
104}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107
108const std::vector<G4int>* G4LossTableBuilder::GetCoupleIndexes() const
109{
110 return theDensityIdx;
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114
115const std::vector<G4double>* G4LossTableBuilder::GetDensityFactors() const
116{
117 return theDensityFactor;
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
123{
124 if (theFlag->empty()) { InitialiseBaseMaterials(); }
125 return (idx < theFlag->size()) ? (*theFlag)[idx] : false;
126}
127
128//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
129
131{
132 if (theFlag->empty()) { InitialiseBaseMaterials(); }
133 return baseMatFlag;
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
137
138void
140 const std::vector<G4PhysicsTable*>& list)
141{
142 InitialiseBaseMaterials(dedxTable);
143 std::size_t n_processes = list.size();
144 if(1 >= n_processes) { return; }
145
146 std::size_t nCouples = dedxTable->size();
147 //G4cout << "Nproc= " << n_processes << " nCouples=" << nCouples << " Nv= "
148 // << dedxTable->size() << G4endl;
149 if(0 >= nCouples) { return; }
150
151 for (std::size_t i=0; i<nCouples; ++i) {
152 auto pv0 = static_cast<G4PhysicsLogVector*>((*(list[0]))[i]);
153 if(pv0 == nullptr) { continue; }
154 std::size_t npoints = pv0->GetVectorLength();
155 auto pv = new G4PhysicsLogVector(*pv0);
156 for (std::size_t j=0; j<npoints; ++j) {
157 G4double dedx = 0.0;
158 for (std::size_t k=0; k<n_processes; ++k) {
159 const G4PhysicsVector* pv1 = (*(list[k]))[i];
160 dedx += (*pv1)[j];
161 }
162 pv->PutValue(j, dedx);
163 }
164 if(splineFlag) { pv->FillSecondDerivatives(); }
166 }
167 //G4cout << "### G4LossTableBuilder::BuildDEDXTable " << G4endl;
168 //G4cout << *dedxTable << G4endl;
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172
174 G4PhysicsTable* rangeTable)
175// Build range table from the energy loss table
176{
177 //G4cout << "### G4LossTableBuilder::BuildRangeTable: DEDX table" << G4endl;
178 //G4cout << *const_cast<G4PhysicsTable*>(dedxTable) << G4endl;
179 const std::size_t nCouples = dedxTable->size();
180 if(0 >= nCouples) { return; }
181
182 const std::size_t n = 100;
183 const G4double del = 1.0/(G4double)n;
184
185 for (std::size_t i=0; i<nCouples; ++i) {
186 auto pv = static_cast<G4PhysicsLogVector*>((*dedxTable)[i]);
187 if((pv == nullptr) || (isBaseMatActive && !(*theFlag)[i])) { continue; }
188 std::size_t npoints = pv->GetVectorLength();
189 std::size_t bin0 = 0;
190 G4double elow = pv->Energy(0);
191 G4double ehigh = pv->Energy(npoints-1);
192 G4double dedx1 = (*pv)[0];
193
194 // protection for specific cases dedx=0
195 if(dedx1 == 0.0) {
196 for (std::size_t k=1; k<npoints; ++k) {
197 ++bin0;
198 elow = pv->Energy(k);
199 dedx1 = (*pv)[k];
200 if(dedx1 > 0.0) { break; }
201 }
202 npoints -= bin0;
203 }
204
205 // initialisation of a new vector
206 if(npoints < 3) { npoints = 3; }
207
208 delete (*rangeTable)[i];
210 if(0 == bin0) { v = new G4PhysicsLogVector(*pv); }
211 else { v = new G4PhysicsLogVector(elow, ehigh, npoints-1, splineFlag); }
212
213 // assumed dedx proportional to beta
214 G4double energy1 = v->Energy(0);
215 G4double range = 2.*energy1/dedx1;
216 /*
217 G4cout << "New Range vector Npoints=" << v->GetVectorLength()
218 << " coupleIdx=" << i << " spline=" << v->GetSpline()
219 << " Elow=" << v->GetMinEnergy() <<" Ehigh=" << v->GetMinEnergy()
220 << " DEDX(Elow)=" << dedx1 << " R(Elow)=" << range << G4endl;
221 */
222 v->PutValue(0,range);
223
224 for (std::size_t j=1; j<npoints; ++j) {
225
226 G4double energy2 = v->Energy(j);
227 G4double de = (energy2 - energy1) * del;
228 G4double energy = energy2 + de*0.5;
229 G4double sum = 0.0;
230 std::size_t idx = j - 1;
231 for (std::size_t k=0; k<n; ++k) {
232 energy -= de;
233 dedx1 = pv->Value(energy, idx);
234 if(dedx1 > 0.0) { sum += de/dedx1; }
235 }
236 range += sum;
237 /*
238 if(energy < 10.)
239 G4cout << "j= " << j << " e1= " << energy1 << " e2= " << energy2
240 << " n= " << n << " range=" << range<< G4endl;
241 */
242 v->PutValue(j,range);
243 energy1 = energy2;
244 }
245 if(splineFlag) { v->FillSecondDerivatives(); }
247 }
248 //G4cout << "### Range table" << G4endl;
249 //G4cout << *rangeTable << G4endl;
250}
251
252//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
253
254void
256 G4PhysicsTable* invRangeTable)
257// Build inverse range table from the energy loss table
258{
259 std::size_t nCouples = rangeTable->size();
260 if(0 >= nCouples) { return; }
261
262 for (std::size_t i=0; i<nCouples; ++i) {
263 G4PhysicsVector* pv = (*rangeTable)[i];
264 if((pv == nullptr) || (isBaseMatActive && !(*theFlag)[i])) { continue; }
265 std::size_t npoints = pv->GetVectorLength();
266
267 delete (*invRangeTable)[i];
268 auto v = new G4PhysicsFreeVector(npoints, splineFlag);
269
270 for (std::size_t j=0; j<npoints; ++j) {
271 G4double e = pv->Energy(j);
272 G4double r = (*pv)[j];
273 v->PutValues(j,r,e);
274 }
275 if (splineFlag) { v->FillSecondDerivatives(); }
276 v->EnableLogBinSearch(theParameters->NumberForFreeVector());
277
278 G4PhysicsTableHelper::SetPhysicsVector(invRangeTable, i, v);
279 }
280 //G4cout << "### Inverse range table" << G4endl;
281 //G4cout << *invRangeTable << G4endl;
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
285
287{
288 if(!isInitializer) { return; }
289 const G4ProductionCutsTable* theCoupleTable=
291 std::size_t nCouples = theCoupleTable->GetTableSize();
292 std::size_t nFlags = theFlag->size();
293 /*
294 G4cout << "### InitialiseBaseMaterials: nCouples=" << nCouples
295 << " nFlags=" << nFlags << " isInit:" << isInitialized
296 << " baseMat:" << baseMatFlag << G4endl;
297 */
298 // define base material flag
299 if(isBaseMatActive && !baseMatFlag) {
300 for(G4int i=0; i<(G4int)nCouples; ++i) {
301 if(nullptr != theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetBaseMaterial()) {
302 baseMatFlag = true;
303 isInitialized = false;
304 break;
305 }
306 }
307 }
308
309 if(nFlags != nCouples) { isInitialized = false; }
310 if(isInitialized) { return; }
311
312 // reserve memory
313 theFlag->resize(nCouples, true);
314 if(nullptr == table) { return; }
315
316 if(baseMatFlag) {
317 theDensityFactor->resize(nCouples,1.0);
318 theDensityIdx->resize(nCouples);
319 }
320
321 // define default flag and index of used material cut couple
322 for(G4int i=0; i<(G4int)nCouples; ++i) {
323 (*theFlag)[i] = table->GetFlag(i);
324 if(baseMatFlag) { (*theDensityIdx)[i] = i; }
325 }
326 isInitialized = true;
327 if(baseMatFlag) {
328 // use base materials
329 for(G4int i=0; i<(G4int)nCouples; ++i) {
330 // base material is needed only for a couple which is not
331 // initialised and for which tables will be computed
332 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
333 auto pcuts = couple->GetProductionCuts();
334 auto mat = couple->GetMaterial();
335 auto bmat = mat->GetBaseMaterial();
336
337 // base material exists - find it and check if it can be reused
338 if(nullptr != bmat) {
339 for(G4int j=0; j<(G4int)nCouples; ++j) {
340 if(j == i) { continue; }
341 auto bcouple = theCoupleTable->GetMaterialCutsCouple(j);
342
343 if(bcouple->GetMaterial() == bmat &&
344 bcouple->GetProductionCuts() == pcuts) {
345
346 // based couple exist in the same region
347 (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
348 (*theDensityIdx)[i] = j;
349 (*theFlag)[i] = false;
350
351 // ensure that there will no double initialisation
352 (*theDensityFactor)[j] = 1.0;
353 (*theDensityIdx)[j] = j;
354 (*theFlag)[j] = true;
355 break;
356 }
357 }
358 }
359 }
360 }
361 /*
362 G4cout << "### G4LossTableBuilder::InitialiseBaseMaterials: flag="
363 << baseMatFlag << G4endl;
364 for(std::size_t i=0; i<nCouples; ++i) {
365 G4cout << "CoupleIdx=" << i << " Flag= " << (*theFlag)[i] << " "
366 << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
367 << " TableFlag= " << table->GetFlag(i)
368 << " " << (*table)[i]
369 << G4endl;
370 }
371 */
372}
373
374//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
375
378 G4VEmModel* model,
379 const G4ParticleDefinition* part,
380 G4double emin, G4double emax,
381 G4bool spline)
382{
383 // check input
385 if(nullptr == table) { return table; }
386 if(emin >= emax) {
387 table->clearAndDestroy();
388 delete table;
389 table = nullptr;
390 return table;
391 }
393 G4int nbins = theParameters->NumberOfBinsPerDecade();
394
395 // Access to materials
396 const G4ProductionCutsTable* theCoupleTable=
398 std::size_t numOfCouples = theCoupleTable->GetTableSize();
399
400 G4PhysicsLogVector* aVector = nullptr;
401
402 for(G4int i=0; i<(G4int)numOfCouples; ++i) {
403 if ((*theFlag)[i]) {
404
405 // create physics vector and fill it
406 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
407 delete (*table)[i];
408
409 // if start from zero then change the scale
410
411 const G4Material* mat = couple->GetMaterial();
412
413 G4double tmin = std::max(emin,model->MinPrimaryEnergy(mat,part));
414 if(0.0 >= tmin) { tmin = CLHEP::eV; }
415 G4int n = nbins;
416
417 if(tmin >= emax) {
418 aVector = nullptr;
419 } else {
420 n *= G4lrint(std::log10(emax/tmin));
421 n = std::max(n, 3);
422 aVector = new G4PhysicsLogVector(tmin, emax, n, spline);
423 }
424
425 if(nullptr != aVector) {
426 //G4cout << part->GetParticleName() << " in " << mat->GetName()
427 // << " tmin= " << tmin << G4endl;
428 for(G4int j=0; j<=n; ++j) {
429 aVector->PutValue(j, model->Value(couple, part,
430 aVector->Energy(j)));
431 }
432 if(spline) { aVector->FillSecondDerivatives(); }
433 }
435 }
436 }
437 /*
438 G4cout << "G4LossTableBuilder::BuildTableForModel done for "
439 << part->GetParticleName() << " and "<< model->GetName()
440 << " " << table << G4endl;
441 */
442 //G4cout << *table << G4endl;
443 return table;
444}
445
446//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static G4EmParameters * Instance()
G4int NumberOfBinsPerDecade() const
G4int NumberForFreeVector() const
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable)
void BuildDEDXTable(G4PhysicsTable *dedxTable, const std::vector< G4PhysicsTable * > &)
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
G4bool GetFlag(size_t idx)
G4LossTableBuilder(G4bool master=true)
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable)
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4Material * GetBaseMaterial() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
void clearAndDestroy()
G4bool GetFlag(std::size_t i) const
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
int G4lrint(double ad)
Definition templates.hh:134