Geant4.10
 All Data Structures 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 71484 2013-06-17 08:12:51Z 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 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71 
73 {
74  splineFlag = true;
75  isInitialized = false;
76 
77  theDensityFactor = new std::vector<G4double>;
78  theDensityIdx = new std::vector<G4int>;
79  theFlag = new std::vector<G4bool>;
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
85 {
86  delete theDensityFactor;
87  delete theDensityIdx;
88  delete theFlag;
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
93 void
95  const std::vector<G4PhysicsTable*>& list)
96 {
97  size_t n_processes = list.size();
98  //G4cout << "Nproc= " << n_processes << " Ncoup= " << dedxTable->size() << G4endl;
99  if(1 >= n_processes) { return; }
100 
101  size_t nCouples = dedxTable->size();
102  if(0 >= nCouples) { return; }
103 
104  for (size_t i=0; i<nCouples; ++i) {
105  // if ((*theFlag)[i]) {
106  G4PhysicsLogVector* pv0 = static_cast<G4PhysicsLogVector*>((*(list[0]))[i]);
107  if(pv0) {
108  size_t npoints = pv0->GetVectorLength();
109  G4PhysicsLogVector* pv = new G4PhysicsLogVector(*pv0);
110  /*
111  G4PhysicsLogVector* pv = new G4PhysicsLogVector(pv0->Energy(0),
112  pv0->GetMaxEnergy(),
113  npoints-1);
114  */
115  pv->SetSpline(splineFlag);
116  for (size_t j=0; j<npoints; ++j) {
117  G4double dedx = 0.0;
118  for (size_t k=0; k<n_processes; ++k) {
119  G4PhysicsVector* pv1 = (*(list[k]))[i];
120  dedx += (*pv1)[j];
121  }
122  pv->PutValue(j, dedx);
123  }
124  if(splineFlag) { pv->FillSecondDerivatives(); }
125  G4PhysicsTableHelper::SetPhysicsVector(dedxTable, i, pv);
126  }
127  }
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
133  G4PhysicsTable* rangeTable,
134  G4bool isIonisation)
135 // Build range table from the energy loss table
136 {
137  size_t nCouples = dedxTable->size();
138  if(0 >= nCouples) { return; }
139 
140  size_t n = 100;
141  G4double del = 1.0/(G4double)n;
142 
143  for (size_t i=0; i<nCouples; ++i) {
144  if(isIonisation) {
145  if( !(*theFlag)[i] ) { continue; }
146  }
147  G4PhysicsLogVector* pv = static_cast<G4PhysicsLogVector*>((*dedxTable)[i]);
148  size_t npoints = pv->GetVectorLength();
149  size_t bin0 = 0;
150  G4double elow = pv->Energy(0);
151  G4double ehigh = pv->Energy(npoints-1);
152  G4double dedx1 = (*pv)[0];
153 
154  //G4cout << "i= " << i << "npoints= " << npoints << " dedx1= " << 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 
218  G4PhysicsTable* invRangeTable,
219  G4bool isIonisation)
220 // Build inverse range table from the energy loss table
221 {
222  size_t nCouples = rangeTable->size();
223  if(0 >= nCouples) { return; }
224 
225  for (size_t i=0; i<nCouples; ++i) {
226 
227  if(isIonisation) {
228  if( !(*theFlag)[i] ) { continue; }
229  }
230  G4PhysicsVector* pv = (*rangeTable)[i];
231  size_t npoints = pv->GetVectorLength();
232  G4double rlow = (*pv)[0];
233  G4double rhigh = (*pv)[npoints-1];
234 
235  delete (*invRangeTable)[i];
236  G4LPhysicsFreeVector* v = new G4LPhysicsFreeVector(npoints,rlow,rhigh);
237  v->SetSpline(splineFlag);
238 
239  for (size_t j=0; j<npoints; ++j) {
240  G4double e = pv->Energy(j);
241  G4double r = (*pv)[j];
242  v->PutValues(j,r,e);
243  }
244  if(splineFlag) { v->FillSecondDerivatives(); }
245 
246  G4PhysicsTableHelper::SetPhysicsVector(invRangeTable, i, v);
247  }
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
251 
252 void
254 {
255  size_t nCouples = table->size();
256  size_t nFlags = theFlag->size();
257 
258  if(nCouples == nFlags && isInitialized) { return; }
259 
260  isInitialized = true;
261 
262  //G4cout << "%%%%%% G4LossTableBuilder::InitialiseBaseMaterials Ncouples= "
263  // << nCouples << " FlagSize= " << nFlags << G4endl;
264 
265  // variable density check
266  const G4ProductionCutsTable* theCoupleTable=
268 
269  /*
270  for(size_t i=0; i<nFlags; ++i) {
271  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i]
272  << " tableFlag= " << table->GetFlag(i) << " "
273  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
274  << G4endl;
275  }
276  */
277 
278  // expand vectors
279  if(nFlags < nCouples) {
280  for(size_t i=nFlags; i<nCouples; ++i) { theDensityFactor->push_back(1.0); }
281  for(size_t i=nFlags; i<nCouples; ++i) { theDensityIdx->push_back(-1); }
282  for(size_t i=nFlags; i<nCouples; ++i) { theFlag->push_back(true); }
283  }
284  for(size_t i=0; i<nCouples; ++i) {
285 
286  // base material is needed only for a couple which is not
287  // initialised and for which tables will be computed
288  (*theFlag)[i] = table->GetFlag(i);
289  if ((*theDensityIdx)[i] < 0) {
290  (*theDensityIdx)[i] = i;
291  const G4MaterialCutsCouple* couple =
292  theCoupleTable->GetMaterialCutsCouple(i);
293  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
294  const G4Material* mat = couple->GetMaterial();
295  const G4Material* bmat = mat->GetBaseMaterial();
296 
297  // base material exists - find it and check if it can be reused
298  if(bmat) {
299  for(size_t j=0; j<nCouples; ++j) {
300 
301  if(j == i) { continue; }
302  const G4MaterialCutsCouple* bcouple =
303  theCoupleTable->GetMaterialCutsCouple(j);
304 
305  if(bcouple->GetMaterial() == bmat &&
306  bcouple->GetProductionCuts() == pcuts) {
307 
308  // based couple exist in the same region
309  (*theDensityIdx)[i] = j;
310  (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
311  (*theFlag)[i] = false;
312 
313  // ensure that there will no double initialisation
314  (*theDensityIdx)[j] = j;
315  (*theDensityFactor)[j] = 1.0;
316  (*theFlag)[j] = true;
317  break;
318  }
319  }
320  }
321  }
322  }
323  /*
324  for(size_t i=0; i<nCouples; ++i) {
325  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i]
326  << " TableFlag= " << table->GetFlag(i) << " "
327  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
328  << G4endl;
329  }
330  G4cout << "%%%%%% G4LossTableBuilder::InitialiseBaseMaterials end"
331  << G4endl;
332  */
333 }
334 
335 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
336 
337 void G4LossTableBuilder::InitialiseCouples()
338 {
339  isInitialized = true;
340 
341  // variable density initialisation for the cas without tables
342  const G4ProductionCutsTable* theCoupleTable=
344  size_t nCouples = theCoupleTable->GetTableSize();
345  //G4cout << "%%%%%% G4LossTableBuilder::InitialiseCouples() nCouples= "
346  // << nCouples << G4endl;
347 
348  theDensityFactor->resize(nCouples, 1.0);
349  theDensityIdx->resize(nCouples, -1);
350  theFlag->resize(nCouples, true);
351 
352  for(size_t i=0; i<nCouples; ++i) {
353 
354  // base material is needed only for a couple which is not
355  // initialised and for which tables will be computed
356  if ((*theDensityIdx)[i] < 0) {
357  (*theDensityIdx)[i] = i;
358  const G4MaterialCutsCouple* couple =
359  theCoupleTable->GetMaterialCutsCouple(i);
360  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
361  const G4Material* mat = couple->GetMaterial();
362  const G4Material* bmat = mat->GetBaseMaterial();
363 
364  // base material exists - find it and check if it can be reused
365  if(bmat) {
366  for(size_t j=0; j<nCouples; ++j) {
367 
368  if(j == i) { continue; }
369  const G4MaterialCutsCouple* bcouple =
370  theCoupleTable->GetMaterialCutsCouple(j);
371 
372  if(bcouple->GetMaterial() == bmat &&
373  bcouple->GetProductionCuts() == pcuts) {
374 
375  // based couple exist in the same region
376  (*theDensityIdx)[i] = j;
377  (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
378  (*theFlag)[i] = false;
379 
380  // ensure that there will no double initialisation
381  (*theDensityIdx)[j] = j;
382  (*theDensityFactor)[j] = 1.0;
383  (*theFlag)[j] = true;
384  break;
385  }
386  }
387  }
388  }
389  }
390  /*
391  for(size_t i=0; i<nCouples; ++i) {
392  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i] << " "
393  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
394  << " DensityFactor= " << (*theDensityFactor)[i]
395  << G4endl;
396  }
397  G4cout << "%%%%%% G4LossTableBuilder::InitialiseCouples() end" << G4endl;
398  */
399 }
400 
401 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
402 
405  G4VEmModel* model,
406  const G4ParticleDefinition* part,
407  G4double emin, G4double emax,
408  G4bool spline)
409 {
410  // check input
412  if(!table) { return table; }
413  if(emin >= emax) {
414  table->clearAndDestroy();
415  delete table;
416  table = 0;
417  return table;
418  }
420 
422 
423  // Access to materials
424  const G4ProductionCutsTable* theCoupleTable=
426  size_t numOfCouples = theCoupleTable->GetTableSize();
427 
428  G4PhysicsLogVector* aVector = 0;
429  //G4PhysicsLogVector* bVector = 0;
430 
431  for(size_t i=0; i<numOfCouples; ++i) {
432 
433  //G4cout<< "i= " << i << " Flag= " << GetFlag(i) << G4endl;
434 
435  if (GetFlag(i)) {
436 
437  // create physics vector and fill it
438  const G4MaterialCutsCouple* couple =
439  theCoupleTable->GetMaterialCutsCouple(i);
440  delete (*table)[i];
441 
442  // if start from zero then change the scale
443 
444  const G4Material* mat = couple->GetMaterial();
445 
446  G4double tmin = std::max(emin,model->MinPrimaryEnergy(mat,part));
447  if(0.0 >= tmin) { tmin = eV; }
448  G4int n = nbins;
449 
450  if(tmin >= emax) {
451  aVector = 0;
452  } else {
453  n *= G4int(std::log10(emax/tmin) + 0.5);
454  if(n < 3) { n = 3; }
455  aVector = new G4PhysicsLogVector(tmin, emax, n);
456  }
457 
458  if(aVector) {
459  aVector->SetSpline(spline);
460  for(G4int j=0; j<=n; ++j) {
461  aVector->PutValue(j, model->Value(couple, part, aVector->Energy(j)));
462  }
463  if(spline) { aVector->FillSecondDerivatives(); }
464  }
465  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
466  }
467  }
468  /*
469  G4cout << "G4LossTableBuilder::BuildTableForModel done for "
470  << part->GetParticleName() << " and "<< model->GetName()
471  << " " << table << G4endl;
472  */
473  return table;
474 }
475 
476 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
477 
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable, G4bool isIonisation=false)
void PutValues(size_t binNumber, G4double binValue, G4double dataValue)
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
static G4LossTableManager * Instance()
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:355
G4double GetDensity() const
Definition: G4Material.hh:178
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)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
bool G4bool
Definition: G4Types.hh:79
void PutValue(size_t index, G4double theValue)
const XML_Char XML_Content * model
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
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)
G4int GetNumberOfBinsPerDecade() const
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:231
G4bool GetFlag(size_t i) const
double G4double
Definition: G4Types.hh:76
G4ProductionCuts * GetProductionCuts() const
void clearAndDestroy()
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:346
const G4Material * GetMaterial() const