Geant4 11.4.0
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
69G4bool G4LossTableBuilder::baseMatFlag = false;
70std::vector<G4double>* G4LossTableBuilder::theDensityFactor = nullptr;
71std::vector<G4int>* G4LossTableBuilder::theDensityIdx = nullptr;
72std::vector<G4bool>* G4LossTableBuilder::theFlag = nullptr;
73std::vector<G4bool>* G4LossTableBuilder::theFluct = nullptr;
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
78 : isInitializer(master)
79{
80 theParameters = G4EmParameters::Instance();
81 if (nullptr == theFlag) {
82 theDensityFactor = new std::vector<G4double>;
83 theDensityIdx = new std::vector<G4int>;
84 theFlag = new std::vector<G4bool>;
85 theFluct = new std::vector<G4bool>;
86 }
87}
88
89//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90
92{
93 if (isInitializer) {
94 delete theDensityFactor;
95 delete theDensityIdx;
96 delete theFlag;
97 delete theFluct;
98 theDensityFactor = nullptr;
99 theDensityIdx = nullptr;
100 theFlag = nullptr;
101 theFluct = nullptr;
102 }
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
107const std::vector<G4int>* G4LossTableBuilder::GetCoupleIndexes()
108{
109 return theDensityIdx;
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
114const std::vector<G4double>* G4LossTableBuilder::GetDensityFactors()
115{
116 return theDensityFactor;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
120
121const std::vector<G4bool>* G4LossTableBuilder::GetFluctuationFlags()
122{
123 return theFluct;
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127
129{
130 return (idx < theFlag->size()) ? (*theFlag)[idx] : false;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
134
136{
137 return baseMatFlag;
138}
139
140//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141
142void
144 const std::vector<G4PhysicsTable*>& list)
145{
146 InitialiseBaseMaterials(dedxTable);
147 std::size_t n_processes = list.size();
148 if(1 >= n_processes) { return; }
149
150 std::size_t nCouples = dedxTable->size();
151 //G4cout << "Nproc= " << n_processes << " nCouples=" << nCouples << " Nv= "
152 // << dedxTable->size() << G4endl;
153 if(0 >= nCouples) { return; }
154
155 for (std::size_t i=0; i<nCouples; ++i) {
156 auto pv0 = static_cast<G4PhysicsLogVector*>((*(list[0]))[i]);
157 //if (0 == i) G4cout << i << ". pv0=" << pv0 << " t:" << list[0] << G4endl;
158 if(pv0 == nullptr) { continue; }
159 std::size_t npoints = pv0->GetVectorLength();
160 auto pv = new G4PhysicsLogVector(*pv0);
161 for (std::size_t j=0; j<npoints; ++j) {
162 G4double dedx = 0.0;
163 for (std::size_t k=0; k<n_processes; ++k) {
164 const G4PhysicsVector* pv1 = (*(list[k]))[i];
165 //if (0 == i) G4cout << " " << k << ". pv1=" << pv1 << " t:" << list[k] << G4endl;
166 dedx += (*pv1)[j];
167 }
168 pv->PutValue(j, dedx);
169 }
170 if(splineFlag) { pv->FillSecondDerivatives(); }
172 }
173 //G4cout << "### G4LossTableBuilder::BuildDEDXTable " << G4endl;
174 //G4cout << *dedxTable << G4endl;
175}
176
177//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178
180 G4PhysicsTable* rangeTable)
181// Build range table from the energy loss table
182{
183 //G4cout << "### G4LossTableBuilder::BuildRangeTable: DEDX table" << G4endl;
184 //G4cout << *const_cast<G4PhysicsTable*>(dedxTable) << G4endl;
185 const std::size_t nCouples = dedxTable->size();
186 if(0 >= nCouples) { return; }
187
188 const std::size_t n = 100;
189 const G4double del = 1.0/(G4double)n;
190
191 for (std::size_t i=0; i<nCouples; ++i) {
192 auto pv = static_cast<G4PhysicsLogVector*>((*dedxTable)[i]);
193 if((pv == nullptr) || (isBaseMatActive && !(*theFlag)[i])) { continue; }
194 std::size_t npoints = pv->GetVectorLength();
195 std::size_t bin0 = 0;
196 G4double elow = pv->Energy(0);
197 G4double ehigh = pv->Energy(npoints-1);
198 G4double dedx1 = (*pv)[0];
199
200 // protection for specific cases dedx=0
201 if(dedx1 == 0.0) {
202 for (std::size_t k=1; k<npoints; ++k) {
203 ++bin0;
204 elow = pv->Energy(k);
205 dedx1 = (*pv)[k];
206 if(dedx1 > 0.0) { break; }
207 }
208 npoints -= bin0;
209 }
210
211 // initialisation of a new vector
212 if(npoints < 3) { npoints = 3; }
213
214 delete (*rangeTable)[i];
216 if(0 == bin0) { v = new G4PhysicsLogVector(*pv); }
217 else { v = new G4PhysicsLogVector(elow, ehigh, npoints-1, splineFlag); }
218
219 // assumed dedx proportional to beta
220 G4double energy1 = v->Energy(0);
221 G4double range = 2.*energy1/dedx1;
222 /*
223 G4cout << "New Range vector Npoints=" << v->GetVectorLength()
224 << " coupleIdx=" << i << " spline=" << v->GetSpline()
225 << " Elow=" << v->GetMinEnergy() <<" Ehigh=" << v->GetMinEnergy()
226 << " DEDX(Elow)=" << dedx1 << " R(Elow)=" << range << G4endl;
227 */
228 v->PutValue(0,range);
229
230 for (std::size_t j=1; j<npoints; ++j) {
231
232 G4double energy2 = v->Energy(j);
233 G4double de = (energy2 - energy1) * del;
234 G4double energy = energy2 + de*0.5;
235 G4double sum = 0.0;
236 std::size_t idx = j - 1;
237 for (std::size_t k=0; k<n; ++k) {
238 energy -= de;
239 dedx1 = pv->Value(energy, idx);
240 if(dedx1 > 0.0) { sum += de/dedx1; }
241 }
242 range += sum;
243 /*
244 if(energy < 10.)
245 G4cout << "j= " << j << " e1= " << energy1 << " e2= " << energy2
246 << " n= " << n << " range=" << range<< G4endl;
247 */
248 v->PutValue(j,range);
249 energy1 = energy2;
250 }
251 if(splineFlag) { v->FillSecondDerivatives(); }
253 }
254 //G4cout << "### Range table" << G4endl;
255 //G4cout << *rangeTable << G4endl;
256}
257
258//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
259
260void
262 G4PhysicsTable* invRangeTable)
263// Build inverse range table from the energy loss table
264{
265 std::size_t nCouples = rangeTable->size();
266 if(0 >= nCouples) { return; }
267
268 for (std::size_t i=0; i<nCouples; ++i) {
269 G4PhysicsVector* pv = (*rangeTable)[i];
270 if((pv == nullptr) || (isBaseMatActive && !(*theFlag)[i])) { continue; }
271 std::size_t npoints = pv->GetVectorLength();
272
273 delete (*invRangeTable)[i];
274 auto v = new G4PhysicsFreeVector(npoints, splineFlag);
275
276 for (std::size_t j=0; j<npoints; ++j) {
277 G4double e = pv->Energy(j);
278 G4double r = (*pv)[j];
279 v->PutValues(j,r,e);
280 }
281 if (splineFlag) { v->FillSecondDerivatives(); }
282 v->EnableLogBinSearch(theParameters->NumberForFreeVector());
283
284 G4PhysicsTableHelper::SetPhysicsVector(invRangeTable, i, v);
285 }
286 //G4cout << "### Inverse range table" << G4endl;
287 //G4cout << *invRangeTable << G4endl;
288}
289
290//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
291
293{
294 if(!isInitializer) { return; }
295 const G4ProductionCutsTable* theCoupleTable=
297 std::size_t nCouples = theCoupleTable->GetTableSize();
298 std::size_t nFlags = theFlag->size();
299 /*
300 G4cout << "### InitialiseBaseMaterials: nCouples=" << nCouples
301 << " nFlags=" << nFlags << " isInit:" << isInitialized
302 << " baseMat:" << baseMatFlag << G4endl;
303 */
304 // define base material flag
305 if(isBaseMatActive && !baseMatFlag) {
306 for(G4int i=0; i<(G4int)nCouples; ++i) {
307 if(nullptr != theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetBaseMaterial()) {
308 baseMatFlag = true;
309 isInitialized = false;
310 break;
311 }
312 }
313 }
314
315 if(nFlags != nCouples) { isInitialized = false; }
316 if(isInitialized) { return; }
317
318 // reserve and fill memory
319 theFlag->resize(nCouples, true);
320 theFluct->resize(nCouples, theParameters->LossFluctuation());
321 theParameters->DefineFluctuationFlags(theFluct);
322
323 theDensityFactor->resize(nCouples,1.0);
324 theDensityIdx->resize(nCouples, 0);
325
326 // define default flag and index of used material cut couple
327 for (G4int i=0; i<(G4int)nCouples; ++i) {
328 (*theFlag)[i] = (nullptr == table) ? true : table->GetFlag(i);
329 (*theDensityIdx)[i] = i;
330 }
331 isInitialized = true;
332 if (!baseMatFlag) { return; }
333
334 // use base materials
335 for (G4int i=0; i<(G4int)nCouples; ++i) {
336 // base material is needed only for a couple which is not
337 // initialised and for which tables will be computed
338 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
339 auto pcuts = couple->GetProductionCuts();
340 auto mat = couple->GetMaterial();
341 auto bmat = mat->GetBaseMaterial();
342
343 // base material exists - find it and check if it can be reused
344 if(nullptr != bmat) {
345 for(G4int j=0; j<(G4int)nCouples; ++j) {
346 if(j == i) { continue; }
347 auto bcouple = theCoupleTable->GetMaterialCutsCouple(j);
348
349 if(bcouple->GetMaterial() == bmat &&
350 bcouple->GetProductionCuts() == pcuts) {
351
352 // based couple exist in the same region
353 (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
354 (*theDensityIdx)[i] = j;
355 (*theFlag)[i] = false;
356
357 // ensure that there will no double initialisation
358 (*theDensityFactor)[j] = 1.0;
359 (*theDensityIdx)[j] = j;
360 (*theFlag)[j] = true;
361 break;
362 }
363 }
364 }
365 }
366}
367
368//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
369
372 G4VEmModel* model,
373 const G4ParticleDefinition* part,
374 G4double emin, G4double emax,
375 G4bool spline)
376{
377 // check input
379 if (nullptr == table) { return table; }
380 if (aTable != nullptr && aTable != table) {
381 aTable->clearAndDestroy();
382 delete aTable;
383 }
384
386 G4int nbins = theParameters->NumberOfBinsPerDecade();
387
388 // Access to materials
389 const G4ProductionCutsTable* theCoupleTable=
391 std::size_t numOfCouples = theCoupleTable->GetTableSize();
392 /*
393 G4cout << " G4LossTableBuilder::BuildTableForModel Ncouple=" << numOfCouples
394 << " isMaster=" << isInitializer << " model:" << model->GetName()
395 << " " << part->GetParticleName() << G4endl;
396 */
397 G4PhysicsLogVector* aVector = nullptr;
398
399 for(G4int i=0; i<(G4int)numOfCouples; ++i) {
400 //G4cout << i << ". " << (*theFlag)[i] << " " << table->GetFlag(i) << G4endl;
401 if (table->GetFlag(i)) {
402
403 // create physics vector and fill it
404 auto couple = theCoupleTable->GetMaterialCutsCouple(i);
405 delete (*table)[i];
406
407 // if start from zero then change the scale
408 const G4Material* mat = couple->GetMaterial();
409
410 G4double tmin = std::max(emin, model->MinPrimaryEnergy(mat,part));
411 if(0.0 >= tmin) { tmin = CLHEP::eV; }
412 G4int n = nbins;
413
414 if(tmin >= emax) {
415 aVector = nullptr;
416 } else {
417 n *= G4lrint(std::log10(emax/tmin));
418 n = std::max(n, 3);
419 aVector = new G4PhysicsLogVector(tmin, emax, n, spline);
420 }
421
422 if(nullptr != aVector) {
423 //G4cout << part->GetParticleName() << " in " << mat->GetName()
424 // << " emin= " << tmin << " emax=" << emax << " n=" << n << G4endl;
425 for(G4int j=0; j<=n; ++j) {
426 G4double e = aVector->Energy(j);
427 G4double y = model->Value(couple, part, e);
428 //G4cout << " " << j << ") E=" << e << " y=" << y << G4endl;
429 aVector->PutValue(j, y);
430 }
431 if(spline) { aVector->FillSecondDerivatives(); }
432 }
434 }
435 }
436 /*
437 G4cout << "G4LossTableBuilder::BuildTableForModel done for "
438 << part->GetParticleName() << " and "<< model->GetName()
439 << " " << table << G4endl;
440 */
441 //G4cout << *table << G4endl;
442 return table;
443}
444
445//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static G4EmParameters * Instance()
static G4bool GetBaseMaterialFlag()
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable)
static const std::vector< G4double > * GetDensityFactors()
void BuildDEDXTable(G4PhysicsTable *dedxTable, const std::vector< G4PhysicsTable * > &)
static const std::vector< G4bool > * GetFluctuationFlags()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
static G4bool GetFlag(std::size_t idx)
static const std::vector< G4int > * GetCoupleIndexes()
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable)
G4LossTableBuilder(G4bool master)
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()
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