Geant4  9.6.p02
 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$
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  // pv = new G4PhysicsLogVector(elow, ehigh, npoints-1);
111  pv->SetSpline(splineFlag);
112  for (size_t j=0; j<npoints; ++j) {
113  G4double dedx = 0.0;
114  for (size_t k=0; k<n_processes; ++k) {
115  G4PhysicsVector* pv1 = (*(list[k]))[i];
116  dedx += (*pv1)[j];
117  }
118  pv->PutValue(j, dedx);
119  }
120  if(splineFlag) { pv->FillSecondDerivatives(); }
121  G4PhysicsTableHelper::SetPhysicsVector(dedxTable, i, pv);
122  }
123  }
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127 
129  G4PhysicsTable* rangeTable,
130  G4bool isIonisation)
131 // Build range table from the energy loss table
132 {
133  size_t nCouples = dedxTable->size();
134  if(0 >= nCouples) { return; }
135 
136  size_t n = 100;
137  G4double del = 1.0/(G4double)n;
138 
139  for (size_t i=0; i<nCouples; ++i) {
140  if(isIonisation) {
141  if( !(*theFlag)[i] ) { continue; }
142  }
143  G4PhysicsLogVector* pv = static_cast<G4PhysicsLogVector*>((*dedxTable)[i]);
144  size_t npoints = pv->GetVectorLength();
145  size_t bin0 = 0;
146  G4double elow = pv->Energy(0);
147  G4double ehigh = pv->Energy(npoints-1);
148  G4double dedx1 = (*pv)[0];
149 
150  //G4cout << "i= " << i << "npoints= " << npoints << " dedx1= " << dedx1 << G4endl;
151 
152  // protection for specific cases dedx=0
153  if(dedx1 == 0.0) {
154  for (size_t k=1; k<npoints; ++k) {
155  bin0++;
156  elow = pv->Energy(k);
157  dedx1 = (*pv)[k];
158  if(dedx1 > 0.0) { break; }
159  }
160  npoints -= bin0;
161  }
162  //G4cout<<"New Range vector" << G4endl;
163  //G4cout<<"nbins= "<<npoints-1<<" elow= "<<elow<<" ehigh= "<<ehigh
164  // <<" bin0= " << bin0 <<G4endl;
165 
166  // initialisation of a new vector
167  if(npoints < 2) { npoints = 2; }
168 
169  delete (*rangeTable)[i];
171  if(0 == bin0) { v = new G4PhysicsLogVector(*pv); }
172  else { v = new G4PhysicsLogVector(elow, ehigh, npoints-1); }
173 
174  // dedx is exact zero cannot build range table
175  if(2 == npoints) {
176  v->PutValue(0,1000.);
177  v->PutValue(1,2000.);
178  G4PhysicsTableHelper::SetPhysicsVector(rangeTable, i, v);
179  return;
180  }
181  v->SetSpline(splineFlag);
182 
183  // assumed dedx proportional to beta
184  G4double energy1 = v->Energy(0);
185  G4double range = 2.*energy1/dedx1;
186  //G4cout << "range0= " << range << G4endl;
187  v->PutValue(0,range);
188 
189  for (size_t j=1; j<npoints; ++j) {
190 
191  G4double energy2 = v->Energy(j);
192  G4double de = (energy2 - energy1) * del;
193  G4double energy = energy2 + de*0.5;
194  G4double sum = 0.0;
195  //G4cout << "j= " << j << " e1= " << energy1 << " e2= " << energy2
196  // << " n= " << n << G4endl;
197  for (size_t k=0; k<n; ++k) {
198  energy -= de;
199  dedx1 = pv->Value(energy);
200  if(dedx1 > 0.0) { sum += de/dedx1; }
201  }
202  range += sum;
203  v->PutValue(j,range);
204  energy1 = energy2;
205  }
206  if(splineFlag) { v->FillSecondDerivatives(); }
207  G4PhysicsTableHelper::SetPhysicsVector(rangeTable, i, v);
208  }
209 }
210 
211 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
212 
214  G4PhysicsTable* invRangeTable,
215  G4bool isIonisation)
216 // Build inverse range table from the energy loss table
217 {
218  size_t nCouples = rangeTable->size();
219  if(0 >= nCouples) { return; }
220 
221  for (size_t i=0; i<nCouples; ++i) {
222 
223  if(isIonisation) {
224  if( !(*theFlag)[i] ) { continue; }
225  }
226  G4PhysicsVector* pv = (*rangeTable)[i];
227  size_t npoints = pv->GetVectorLength();
228  G4double rlow = (*pv)[0];
229  G4double rhigh = (*pv)[npoints-1];
230 
231  delete (*invRangeTable)[i];
232  G4LPhysicsFreeVector* v = new G4LPhysicsFreeVector(npoints,rlow,rhigh);
233  v->SetSpline(splineFlag);
234 
235  for (size_t j=0; j<npoints; ++j) {
236  G4double e = pv->Energy(j);
237  G4double r = (*pv)[j];
238  v->PutValues(j,r,e);
239  }
240  if(splineFlag) { v->FillSecondDerivatives(); }
241 
242  G4PhysicsTableHelper::SetPhysicsVector(invRangeTable, i, v);
243  }
244 }
245 
246 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
247 
248 void
250 {
251  size_t nCouples = table->size();
252  size_t nFlags = theFlag->size();
253 
254  if(nCouples == nFlags && isInitialized) { return; }
255 
256  isInitialized = true;
257 
258  //G4cout << "%%%%%% G4LossTableBuilder::InitialiseBaseMaterials Ncouples= "
259  // << nCouples << " FlagSize= " << nFlags << G4endl;
260 
261  // variable density check
262  const G4ProductionCutsTable* theCoupleTable=
264 
265  /*
266  for(size_t i=0; i<nFlags; ++i) {
267  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i]
268  << " tableFlag= " << table->GetFlag(i) << " "
269  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
270  << G4endl;
271  }
272  */
273 
274  // expand vectors
275  if(nFlags < nCouples) {
276  for(size_t i=nFlags; i<nCouples; ++i) { theDensityFactor->push_back(1.0); }
277  for(size_t i=nFlags; i<nCouples; ++i) { theDensityIdx->push_back(-1); }
278  for(size_t i=nFlags; i<nCouples; ++i) { theFlag->push_back(true); }
279  }
280  for(size_t i=0; i<nCouples; ++i) {
281 
282  // base material is needed only for a couple which is not
283  // initialised and for which tables will be computed
284  (*theFlag)[i] = table->GetFlag(i);
285  if ((*theDensityIdx)[i] < 0) {
286  (*theDensityIdx)[i] = i;
287  const G4MaterialCutsCouple* couple =
288  theCoupleTable->GetMaterialCutsCouple(i);
289  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
290  const G4Material* mat = couple->GetMaterial();
291  const G4Material* bmat = mat->GetBaseMaterial();
292 
293  // base material exists - find it and check if it can be reused
294  if(bmat) {
295  for(size_t j=0; j<nCouples; ++j) {
296 
297  if(j == i) { continue; }
298  const G4MaterialCutsCouple* bcouple =
299  theCoupleTable->GetMaterialCutsCouple(j);
300 
301  if(bcouple->GetMaterial() == bmat &&
302  bcouple->GetProductionCuts() == pcuts) {
303 
304  // based couple exist in the same region
305  (*theDensityIdx)[i] = j;
306  (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
307  (*theFlag)[i] = false;
308 
309  // ensure that there will no double initialisation
310  (*theDensityIdx)[j] = j;
311  (*theDensityFactor)[j] = 1.0;
312  (*theFlag)[j] = true;
313  break;
314  }
315  }
316  }
317  }
318  }
319  /*
320  for(size_t i=0; i<nCouples; ++i) {
321  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i]
322  << " TableFlag= " << table->GetFlag(i) << " "
323  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
324  << G4endl;
325  }
326  G4cout << "%%%%%% G4LossTableBuilder::InitialiseBaseMaterials end"
327  << G4endl;
328  */
329 }
330 
331 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
332 
333 void G4LossTableBuilder::InitialiseCouples()
334 {
335  isInitialized = true;
336 
337  // variable density initialisation for the cas without tables
338  const G4ProductionCutsTable* theCoupleTable=
340  size_t nCouples = theCoupleTable->GetTableSize();
341  //G4cout << "%%%%%% G4LossTableBuilder::InitialiseCouples() nCouples= "
342  // << nCouples << G4endl;
343 
344  theDensityFactor->resize(nCouples, 1.0);
345  theDensityIdx->resize(nCouples, -1);
346  theFlag->resize(nCouples, true);
347 
348  for(size_t i=0; i<nCouples; ++i) {
349 
350  // base material is needed only for a couple which is not
351  // initialised and for which tables will be computed
352  if ((*theDensityIdx)[i] < 0) {
353  (*theDensityIdx)[i] = i;
354  const G4MaterialCutsCouple* couple =
355  theCoupleTable->GetMaterialCutsCouple(i);
356  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
357  const G4Material* mat = couple->GetMaterial();
358  const G4Material* bmat = mat->GetBaseMaterial();
359 
360  // base material exists - find it and check if it can be reused
361  if(bmat) {
362  for(size_t j=0; j<nCouples; ++j) {
363 
364  if(j == i) { continue; }
365  const G4MaterialCutsCouple* bcouple =
366  theCoupleTable->GetMaterialCutsCouple(j);
367 
368  if(bcouple->GetMaterial() == bmat &&
369  bcouple->GetProductionCuts() == pcuts) {
370 
371  // based couple exist in the same region
372  (*theDensityIdx)[i] = j;
373  (*theDensityFactor)[i] = mat->GetDensity()/bmat->GetDensity();
374  (*theFlag)[i] = false;
375 
376  // ensure that there will no double initialisation
377  (*theDensityIdx)[j] = j;
378  (*theDensityFactor)[j] = 1.0;
379  (*theFlag)[j] = true;
380  break;
381  }
382  }
383  }
384  }
385  }
386  /*
387  for(size_t i=0; i<nCouples; ++i) {
388  G4cout << "CoupleIdx= " << i << " Flag= " << (*theFlag)[i] << " "
389  << theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial()->GetName()
390  << " DensityFactor= " << (*theDensityFactor)[i]
391  << G4endl;
392  }
393  G4cout << "%%%%%% G4LossTableBuilder::InitialiseCouples() end" << G4endl;
394  */
395 }
396 
397 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
398 
401  G4VEmModel* model,
402  const G4ParticleDefinition* part,
403  G4double emin, G4double emax,
404  G4bool spline)
405 {
406  // check input
408  if(!table) { return table; }
409  if(emin >= emax) {
410  table->clearAndDestroy();
411  delete table;
412  table = 0;
413  return table;
414  }
416 
417  G4int nbins = G4int(std::log10(emax/emin) + 0.5)
419  if(nbins < 3) { nbins = 3; }
420 
421  // Access to materials
422  const G4ProductionCutsTable* theCoupleTable=
424  size_t numOfCouples = theCoupleTable->GetTableSize();
425 
426  G4PhysicsLogVector* aVector = 0;
427  G4PhysicsLogVector* bVector = 0;
428 
429  for(size_t i=0; i<numOfCouples; ++i) {
430 
431  //G4cout<< "i= " << i << " Flag= " << GetFlag(i) << G4endl;
432 
433  if (GetFlag(i)) {
434 
435  // create physics vector and fill it
436  const G4MaterialCutsCouple* couple =
437  theCoupleTable->GetMaterialCutsCouple(i);
438  delete (*table)[i];
439 
440  // if start from zero then change the scale
441 
442  const G4Material* mat = couple->GetMaterial();
443 
444  G4double tmin = std::max(emin,model->MinPrimaryEnergy(mat,part));
445  if(0.0 >= tmin) { tmin = eV; }
446  G4int n = nbins + 1;
447 
448  if(tmin >= emax) {
449  aVector = 0;
450  } else if(tmin > emin) {
451  G4int bin = nbins*G4int(std::log10(emax/tmin) + 0.5);
452  if(bin < 3) { bin = 3; }
453  n = bin + 1;
454  aVector = new G4PhysicsLogVector(tmin, emax, bin);
455 
456  } else if(!bVector) {
457  aVector = new G4PhysicsLogVector(emin, emax, nbins);
458  bVector = aVector;
459 
460  } else {
461  aVector = new G4PhysicsLogVector(*bVector);
462  }
463 
464  if(aVector) {
465  aVector->SetSpline(spline);
466  for(G4int j=0; j<n; ++j) {
467  aVector->PutValue(j, model->Value(couple, part, aVector->Energy(j)));
468  }
469  if(spline) { aVector->FillSecondDerivatives(); }
470  }
471  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
472  }
473  }
474  /*
475  G4cout << "G4LossTableBuilder::BuildTableForModel done for "
476  << part->GetParticleName() << " and "<< model->GetName()
477  << " " << table << G4endl;
478  */
479  return table;
480 }
481 
482 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
483