Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4LossTableBuilder.cc 96088 2016-03-14 16:03:38Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4LossTableBuilder
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 //
41 // 23-01-03 V.Ivanchenko Cut per region
42 // 21-07-04 V.Ivanchenko Fix problem of range for dedx=0
43 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
44 // 07-12-04 Fix of BuildDEDX table (V.Ivanchenko)
45 // 27-03-06 Add bool options isIonisation (V.Ivanchenko)
46 // 16-01-07 Fill new (not old) DEDX table (V.Ivanchenko)
47 // 12-02-07 Use G4LPhysicsFreeVector for the inverse range table (V.Ivanchenko)
48 // 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
49 //
50 // Class Description:
51 //
52 // -------------------------------------------------------------------
53 //
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56 
57 #include "G4LossTableBuilder.hh"
58 #include "G4SystemOfUnits.hh"
59 #include "G4PhysicsTable.hh"
60 #include "G4PhysicsLogVector.hh"
61 #include "G4PhysicsTableHelper.hh"
62 #include "G4LPhysicsFreeVector.hh"
63 #include "G4ProductionCutsTable.hh"
64 #include "G4MaterialCutsCouple.hh"
65 #include "G4Material.hh"
66 #include "G4VEmModel.hh"
67 #include "G4ParticleDefinition.hh"
68 #include "G4LossTableManager.hh"
69 #include "G4EmParameters.hh"
70 
71 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
72 
74 {
75  theParameters = G4EmParameters::Instance();
76  splineFlag = true;
77  isInitialized = false;
78 
79  theDensityFactor = new std::vector<G4double>;
80  theDensityIdx = new std::vector<G4int>;
81  theFlag = new std::vector<G4bool>;
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {
88  delete theDensityFactor;
89  delete theDensityIdx;
90  delete theFlag;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 
95 void
97  const std::vector<G4PhysicsTable*>& list)
98 {
99  size_t n_processes = list.size();
100  //G4cout << "Nproc= " << n_processes << " Ncoup= "
101  //<< dedxTable->size() << G4endl;
102  if(1 >= n_processes) { return; }
103 
104  size_t nCouples = dedxTable->size();
105  if(0 >= nCouples) { return; }
106 
107  for (size_t i=0; i<nCouples; ++i) {
108  // if ((*theFlag)[i]) {
109  G4PhysicsLogVector* pv0 =
110  static_cast<G4PhysicsLogVector*>((*(list[0]))[i]);
111  if(pv0) {
112  size_t npoints = pv0->GetVectorLength();
113  G4PhysicsLogVector* pv = new G4PhysicsLogVector(*pv0);
114  pv->SetSpline(splineFlag);
115  for (size_t j=0; j<npoints; ++j) {
116  G4double dedx = 0.0;
117  for (size_t k=0; k<n_processes; ++k) {
118  G4PhysicsVector* pv1 = (*(list[k]))[i];
119  dedx += (*pv1)[j];
120  }
121  pv->PutValue(j, dedx);
122  }
123  if(splineFlag) { pv->FillSecondDerivatives(); }
124  G4PhysicsTableHelper::SetPhysicsVector(dedxTable, i, pv);
125  }
126  }
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
132  G4PhysicsTable* rangeTable,
133  G4bool isIonisation)
134 // Build range table from the energy loss table
135 {
136  size_t nCouples = dedxTable->size();
137  if(0 >= nCouples) { return; }
138 
139  size_t n = 100;
140  G4double del = 1.0/(G4double)n;
141 
142  for (size_t i=0; i<nCouples; ++i) {
143  if(isIonisation) {
144  if( !(*theFlag)[i] ) { continue; }
145  }
146  G4PhysicsLogVector* pv = static_cast<G4PhysicsLogVector*>((*dedxTable)[i]);
147  size_t npoints = pv->GetVectorLength();
148  size_t bin0 = 0;
149  G4double elow = pv->Energy(0);
150  G4double ehigh = pv->Energy(npoints-1);
151  G4double dedx1 = (*pv)[0];
152 
153  //G4cout << "i= " << i << "npoints= " << npoints << " dedx1= "
154  //<< dedx1 << G4endl;
155 
156  // protection for specific cases dedx=0
157  if(dedx1 == 0.0) {
158  for (size_t k=1; k<npoints; ++k) {
159  bin0++;
160  elow = pv->Energy(k);
161  dedx1 = (*pv)[k];
162  if(dedx1 > 0.0) { break; }
163  }
164  npoints -= bin0;
165  }
166  //G4cout<<"New Range vector" << G4endl;
167  //G4cout<<"nbins= "<<npoints-1<<" elow= "<<elow<<" ehigh= "<<ehigh
168  // <<" bin0= " << bin0 <<G4endl;
169 
170  // initialisation of a new vector
171  if(npoints < 2) { npoints = 2; }
172 
173  delete (*rangeTable)[i];
175  if(0 == bin0) { v = new G4PhysicsLogVector(*pv); }
176  else { v = new G4PhysicsLogVector(elow, ehigh, npoints-1); }
177 
178  // dedx is exact zero cannot build range table
179  if(2 == npoints) {
180  v->PutValue(0,1000.);
181  v->PutValue(1,2000.);
182  G4PhysicsTableHelper::SetPhysicsVector(rangeTable, i, v);
183  return;
184  }
185  v->SetSpline(splineFlag);
186 
187  // assumed dedx proportional to beta
188  G4double energy1 = v->Energy(0);
189  G4double range = 2.*energy1/dedx1;
190  //G4cout << "range0= " << range << G4endl;
191  v->PutValue(0,range);
192 
193  for (size_t j=1; j<npoints; ++j) {
194 
195  G4double energy2 = v->Energy(j);
196  G4double de = (energy2 - energy1) * del;
197  G4double energy = energy2 + de*0.5;
198  G4double sum = 0.0;
199  //G4cout << "j= " << j << " e1= " << energy1 << " e2= " << energy2
200  // << " n= " << n << G4endl;
201  for (size_t k=0; k<n; ++k) {
202  energy -= de;
203  dedx1 = pv->Value(energy);
204  if(dedx1 > 0.0) { sum += de/dedx1; }
205  }
206  range += sum;
207  v->PutValue(j,range);
208  energy1 = energy2;
209  }
210  if(splineFlag) { v->FillSecondDerivatives(); }
211  G4PhysicsTableHelper::SetPhysicsVector(rangeTable, i, v);
212  }
213 }
214 
215 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
216 
217 void
219  G4PhysicsTable* invRangeTable,
220  G4bool isIonisation)
221 // Build inverse range table from the energy loss table
222 {
223  size_t nCouples = rangeTable->size();
224  if(0 >= nCouples) { return; }
225 
226  for (size_t i=0; i<nCouples; ++i) {
227 
228  if(isIonisation) {
229  if( !(*theFlag)[i] ) { continue; }
230  }
231  G4PhysicsVector* pv = (*rangeTable)[i];
232  size_t npoints = pv->GetVectorLength();
233  G4double rlow = (*pv)[0];
234  G4double rhigh = (*pv)[npoints-1];
235 
236  delete (*invRangeTable)[i];
237  G4LPhysicsFreeVector* v = new G4LPhysicsFreeVector(npoints,rlow,rhigh);
238  v->SetSpline(splineFlag);
239 
240  for (size_t j=0; j<npoints; ++j) {
241  G4double e = pv->Energy(j);
242  G4double r = (*pv)[j];
243  v->PutValues(j,r,e);
244  }
245  if(splineFlag) { v->FillSecondDerivatives(); }
246 
247  G4PhysicsTableHelper::SetPhysicsVector(invRangeTable, i, v);
248  }
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252 
253 void
255 {
256  size_t nCouples = table->size();
257  size_t nFlags = theFlag->size();
258 
259  if(nCouples == nFlags && isInitialized) { return; }
260 
261  isInitialized = true;
262 
263  //G4cout << "%%%%%% G4LossTableBuilder::InitialiseBaseMaterials Ncouples= "
264  // << nCouples << " FlagSize= " << nFlags << G4endl;
265 
266  // variable density check
267  const G4ProductionCutsTable* theCoupleTable=
269 
270  /*
271  for(size_t i=0; i<nFlags; ++i) {
272  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i]
273  << " tableFlag= " << table->GetFlag(i) << " "
274  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
275  << G4endl;
276  }
277  */
278 
279  // expand vectors
280  if(nFlags < nCouples) {
281  for(size_t i=nFlags; i<nCouples; ++i) {
282  theDensityFactor->push_back(1.0);
283  }
284  for(size_t i=nFlags; i<nCouples; ++i) { theDensityIdx->push_back(-1); }
285  for(size_t i=nFlags; i<nCouples; ++i) { theFlag->push_back(true); }
286  }
287  for(size_t i=0; i<nCouples; ++i) {
288 
289  // base material is needed only for a couple which is not
290  // initialised and for which tables will be computed
291  (*theFlag)[i] = table->GetFlag(i);
292  if ((*theDensityIdx)[i] < 0) {
293  (*theDensityIdx)[i] = i;
294  const G4MaterialCutsCouple* couple =
295  theCoupleTable->GetMaterialCutsCouple(i);
296  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
297  const G4Material* mat = couple->GetMaterial();
298  const G4Material* bmat = mat->GetBaseMaterial();
299 
300  // base material exists - find it and check if it can be reused
301  if(bmat) {
302  for(size_t j=0; j<nCouples; ++j) {
303 
304  if(j == i) { continue; }
305  const G4MaterialCutsCouple* bcouple =
306  theCoupleTable->GetMaterialCutsCouple(j);
307 
308  if(bcouple->GetMaterial() == bmat &&
309  bcouple->GetProductionCuts() == pcuts) {
310 
311  // based couple exist in the same region
312  (*theDensityIdx)[i] = j;
313  (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
314  (*theFlag)[i] = false;
315 
316  // ensure that there will no double initialisation
317  (*theDensityIdx)[j] = j;
318  (*theDensityFactor)[j] = 1.0;
319  (*theFlag)[j] = true;
320  break;
321  }
322  }
323  }
324  }
325  }
326  /*
327  for(size_t i=0; i<nCouples; ++i) {
328  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i]
329  << " TableFlag= " << table->GetFlag(i) << " "
330  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
331  << G4endl;
332  }
333  G4cout << "%%%%%% G4LossTableBuilder::InitialiseBaseMaterials end"
334  << G4endl;
335  */
336 }
337 
338 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
339 
340 void G4LossTableBuilder::InitialiseCouples()
341 {
342  isInitialized = true;
343 
344  // variable density initialisation for the cas without tables
345  const G4ProductionCutsTable* theCoupleTable=
347  size_t nCouples = theCoupleTable->GetTableSize();
348  //G4cout << "%%%%%% G4LossTableBuilder::InitialiseCouples() nCouples= "
349  // << nCouples << G4endl;
350 
351  theDensityFactor->resize(nCouples, 1.0);
352  theDensityIdx->resize(nCouples, -1);
353  theFlag->resize(nCouples, true);
354 
355  for(size_t i=0; i<nCouples; ++i) {
356 
357  // base material is needed only for a couple which is not
358  // initialised and for which tables will be computed
359  if ((*theDensityIdx)[i] < 0) {
360  (*theDensityIdx)[i] = i;
361  const G4MaterialCutsCouple* couple =
362  theCoupleTable->GetMaterialCutsCouple(i);
363  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
364  const G4Material* mat = couple->GetMaterial();
365  const G4Material* bmat = mat->GetBaseMaterial();
366 
367  // base material exists - find it and check if it can be reused
368  if(bmat) {
369  for(size_t j=0; j<nCouples; ++j) {
370 
371  if(j == i) { continue; }
372  const G4MaterialCutsCouple* bcouple =
373  theCoupleTable->GetMaterialCutsCouple(j);
374 
375  if(bcouple->GetMaterial() == bmat &&
376  bcouple->GetProductionCuts() == pcuts) {
377 
378  // based couple exist in the same region
379  (*theDensityIdx)[i] = j;
380  (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
381  (*theFlag)[i] = false;
382 
383  // ensure that there will no double initialisation
384  (*theDensityIdx)[j] = j;
385  (*theDensityFactor)[j] = 1.0;
386  (*theFlag)[j] = true;
387  break;
388  }
389  }
390  }
391  }
392  }
393  /*
394  for(size_t i=0; i<nCouples; ++i) {
395  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i] << " "
396  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
397  << " DensityFactor= " << (*theDensityFactor)[i]
398  << G4endl;
399  }
400  G4cout << "%%%%%% G4LossTableBuilder::InitialiseCouples() end" << G4endl;
401  */
402 }
403 
404 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
405 
408  G4VEmModel* model,
409  const G4ParticleDefinition* part,
410  G4double emin, G4double emax,
411  G4bool spline)
412 {
413  // check input
415  if(!table) { return table; }
416  if(emin >= emax) {
417  table->clearAndDestroy();
418  delete table;
419  table = 0;
420  return table;
421  }
423 
424  G4int nbins = theParameters->NumberOfBinsPerDecade();
425 
426  // Access to materials
427  const G4ProductionCutsTable* theCoupleTable=
429  size_t numOfCouples = theCoupleTable->GetTableSize();
430 
431  G4PhysicsLogVector* aVector = 0;
432 
433  for(size_t i=0; i<numOfCouples; ++i) {
434 
435  //G4cout<< "i= " << i << " Flag= " << GetFlag(i) << G4endl;
436 
437  if (GetFlag(i)) {
438 
439  // create physics vector and fill it
440  const G4MaterialCutsCouple* couple =
441  theCoupleTable->GetMaterialCutsCouple(i);
442  delete (*table)[i];
443 
444  // if start from zero then change the scale
445 
446  const G4Material* mat = couple->GetMaterial();
447 
448  G4double tmin = std::max(emin,model->MinPrimaryEnergy(mat,part));
449  if(0.0 >= tmin) { tmin = eV; }
450  G4int n = nbins;
451 
452  if(tmin >= emax) {
453  aVector = 0;
454  } else {
455  n *= G4int(std::log10(emax/tmin) + 0.5);
456  if(n < 3) { n = 3; }
457  aVector = new G4PhysicsLogVector(tmin, emax, n);
458  }
459 
460  if(aVector) {
461  aVector->SetSpline(spline);
462  for(G4int j=0; j<=n; ++j) {
463  aVector->PutValue(j, model->Value(couple, part,
464  aVector->Energy(j)));
465  }
466  if(spline) { aVector->FillSecondDerivatives(); }
467  }
468  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
469  }
470  }
471  /*
472  G4cout << "G4LossTableBuilder::BuildTableForModel done for "
473  << part->GetParticleName() << " and "<< model->GetName()
474  << " " << table << G4endl;
475  */
476  return table;
477 }
478 
479 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
480 
G4int NumberOfBinsPerDecade() const
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable, G4bool isIonisation=false)
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
virtual G4double MinPrimaryEnergy(const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
Definition: G4VEmModel.cc:386
G4double GetDensity() const
Definition: G4Material.hh:180
G4bool GetFlag(size_t idx) const
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void SetSpline(G4bool)
const G4ParticleDefinition const G4Material *G4double range
bool G4bool
Definition: G4Types.hh:79
void PutValue(size_t index, G4double theValue)
static constexpr double eV
Definition: G4SIunits.hh:215
const G4int n
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable, G4bool isIonisation=false)
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
tuple v
Definition: test.py:18
static const G4double emax
void BuildDEDXTable(G4PhysicsTable *dedxTable, const std::vector< G4PhysicsTable * > &)
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void InitialiseBaseMaterials(G4PhysicsTable *table)
G4double energy(const ThreeVector &p, const G4double m)
static G4EmParameters * Instance()
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:233
void PutValues(size_t index, G4double e, G4double dataValue)
G4bool GetFlag(size_t i) const
double G4double
Definition: G4Types.hh:76
G4ProductionCuts * GetProductionCuts() const
const XML_Char XML_Content * model
Definition: expat.h:151
void clearAndDestroy()
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:377
const G4Material * GetMaterial() const