Geant4  10.02.p02
G4PenelopeBremsstrahlungAngular.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: G4PenelopeBremsstrahlungAngular.cc 97613 2016-06-06 12:24:51Z gcosmo $
27 //
28 // --------------------------------------------------------------
29 //
30 // File name: G4PenelopeBremsstrahlungAngular
31 //
32 // Author: Luciano Pandola
33 //
34 // Creation date: November 2010
35 //
36 // History:
37 // -----------
38 // 23 Nov 2010 L. Pandola 1st implementation
39 // 24 May 2011 L. Pandola Renamed (make v2008 as default Penelope)
40 // 13 Mar 2012 L. Pandola Made a derived class of G4VEmAngularDistribution
41 // and update the interface accordingly
42 // 18 Jul 2012 L. Pandola Migrated to the new basic interface of G4VEmAngularDistribution
43 // Now returns a G4ThreeVector and takes care of the rotation
44 // 03 Oct 2013 L. Pandola Migrated to MT: only the master model handles tables
45 // 17 Oct 2013 L. Pandola Partially revert MT migration. The angular generator is kept as
46 // thread-local, and each worker has full access to it.
47 //
48 //----------------------------------------------------------------
49 
51 
52 #include "globals.hh"
53 #include "G4PhysicalConstants.hh"
54 #include "G4SystemOfUnits.hh"
55 #include "G4PhysicsFreeVector.hh"
56 #include "G4PhysicsTable.hh"
57 #include "G4Material.hh"
58 #include "Randomize.hh"
59 #include "G4Exp.hh"
60 
62  G4VEmAngularDistribution("Penelope"), theEffectiveZSq(0),
63  theLorentzTables1(0),theLorentzTables2(0)
64 
65 {
66  dataRead = false;
67  verbosityLevel = 0;
68 }
69 
70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71 
73 {
74  ClearTables();
75 }
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
80 {
81  ClearTables();
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {
88  std::map<G4double,G4PhysicsTable*>::iterator j;
89 
91  {
92  for (j=theLorentzTables1->begin(); j != theLorentzTables1->end(); j++)
93  {
94  G4PhysicsTable* tab = j->second;
95  //tab->clearAndDestroy();
96  delete tab;
97  }
98  delete theLorentzTables1;
100  }
101 
102  if (theLorentzTables2)
103  {
104  for (j=theLorentzTables2->begin(); j != theLorentzTables2->end(); j++)
105  {
106  G4PhysicsTable* tab = j->second;
107  //tab->clearAndDestroy();
108  delete tab;
109  }
110  delete theLorentzTables2;
111  theLorentzTables2 = 0;
112  }
113  if (theEffectiveZSq)
114  {
115  delete theEffectiveZSq;
116  theEffectiveZSq = 0;
117  }
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
123 {
124  //Read information from DataBase file
125  char* path = getenv("G4LEDATA");
126  if (!path)
127  {
128  G4String excep =
129  "G4PenelopeBremsstrahlungAngular - G4LEDATA environment variable not set!";
130  G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
131  "em0006",FatalException,excep);
132  return;
133  }
134  G4String pathString(path);
135  G4String pathFile = pathString + "/penelope/bremsstrahlung/pdbrang.p08";
136  std::ifstream file(pathFile);
137 
138  if (!file.is_open())
139  {
140  G4String excep = "G4PenelopeBremsstrahlungAngular - data file " + pathFile + " not found!";
141  G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
142  "em0003",FatalException,excep);
143  return;
144  }
145  G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
146 
147  for (k=0;k<NumberofKPoints;k++)
148  for (i=0;i<NumberofZPoints;i++)
149  for (j=0;j<NumberofEPoints;j++)
150  {
151  G4double a1,a2;
152  G4int ik1,iz1,ie1;
153  G4double zr,er,kr;
154  file >> iz1 >> ie1 >> ik1 >> zr >> er >> kr >> a1 >> a2;
155  //check the data are correct
156  if ((iz1-1 == i) && (ik1-1 == k) && (ie1-1 == j))
157  {
158  QQ1[i][j][k]=a1;
159  QQ2[i][j][k]=a2;
160  }
161  else
162  {
164  ed << "Corrupted data file " << pathFile << "?" << G4endl;
165  G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
166  "em0005",FatalException,ed);
167  }
168  }
169  file.close();
170  dataRead = true;
171 }
172 
173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
174 
176 {
177  //Unused at the moment: the G4PenelopeBremsstrahlungAngular is thread-local, so each worker
178  //builds its own version of the tables.
179  /*
180  if (!isMaster)
181  //Should not be here!
182  G4Exception("G4PenelopeBremsstrahlungAngular::PrepareTables()",
183  "em0100",FatalException,"Worker thread in this method");
184  */
185 
186  //Check if data file has already been read
187  if (!dataRead)
188  {
189  ReadDataFile();
190  if (!dataRead)
191  G4Exception("G4PenelopeBremsstrahlungAngular::PrepareInterpolationTables()",
192  "em2001",FatalException,"Unable to build interpolation table");
193  }
194 
195  if (!theLorentzTables1)
196  theLorentzTables1 = new std::map<G4double,G4PhysicsTable*>;
197  if (!theLorentzTables2)
198  theLorentzTables2 = new std::map<G4double,G4PhysicsTable*>;
199 
200  G4double Zmat = CalculateEffectiveZ(material);
201 
202  const G4int reducedEnergyGrid=21;
203  //Support arrays.
204  G4double betas[NumberofEPoints]; //betas for interpolation
205  //tables for interpolation
208  //expanded tables for interpolation
209  G4double Q1E[NumberofEPoints][reducedEnergyGrid];
210  G4double Q2E[NumberofEPoints][reducedEnergyGrid];
211  G4double pZ[NumberofZPoints] = {2.0,8.0,13.0,47.0,79.0,92.0};
212 
213  G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
214  //Interpolation in Z
215  for (i=0;i<NumberofEPoints;i++)
216  {
217  for (j=0;j<NumberofKPoints;j++)
218  {
221 
222  //fill vectors
223  for (k=0;k<NumberofZPoints;k++)
224  {
225  QQ1vector->PutValue(k,pZ[k],std::log(QQ1[k][i][j]));
226  QQ2vector->PutValue(k,pZ[k],QQ2[k][i][j]);
227  }
228 
229  QQ1vector->SetSpline(true);
230  QQ2vector->SetSpline(true);
231 
232  Q1[i][j]= G4Exp(QQ1vector->Value(Zmat));
233  Q2[i][j]=QQ2vector->Value(Zmat);
234  delete QQ1vector;
235  delete QQ2vector;
236  }
237  }
238  G4double pE[NumberofEPoints] = {1.0e-03*MeV,5.0e-03*MeV,1.0e-02*MeV,5.0e-02*MeV,
239  1.0e-01*MeV,5.0e-01*MeV};
240  G4double pK[NumberofKPoints] = {0.0,0.6,0.8,0.95};
241  G4double ppK[reducedEnergyGrid];
242 
243  for(i=0;i<reducedEnergyGrid;i++)
244  ppK[i]=((G4double) i) * 0.05;
245 
246 
247  for(i=0;i<NumberofEPoints;i++)
248  betas[i]=std::sqrt(pE[i]*(pE[i]+2*electron_mass_c2))/(pE[i]+electron_mass_c2);
249 
250 
251  for (i=0;i<NumberofEPoints;i++)
252  {
253  for (j=0;j<NumberofKPoints;j++)
254  Q1[i][j]=Q1[i][j]/Zmat;
255  }
256 
257  //Expanded table of distribution parameters
258  for (i=0;i<NumberofEPoints;i++)
259  {
262 
263  for (j=0;j<NumberofKPoints;j++)
264  {
265  Q1vector->PutValue(j,pK[j],std::log(Q1[i][j])); //logarithmic
266  Q2vector->PutValue(j,pK[j],Q2[i][j]);
267  }
268 
269  for (j=0;j<reducedEnergyGrid;j++)
270  {
271  Q1E[i][j]=Q1vector->Value(ppK[j]);
272  Q2E[i][j]=Q2vector->Value(ppK[j]);
273  }
274  delete Q1vector;
275  delete Q2vector;
276  }
277  //
278  //TABLES to be stored
279  //
280  G4PhysicsTable* theTable1 = new G4PhysicsTable();
281  G4PhysicsTable* theTable2 = new G4PhysicsTable();
282  // the table will contain reducedEnergyGrid G4PhysicsFreeVectors with different
283  // values of k,
284  // Each of the G4PhysicsFreeVectors has a profile of
285  // y vs. E
286  //
287  //reserve space of the vectors.
288  for (j=0;j<reducedEnergyGrid;j++)
289  {
290  G4PhysicsFreeVector* thevec = new G4PhysicsFreeVector(NumberofEPoints);
291  theTable1->push_back(thevec);
292  G4PhysicsFreeVector* thevec2 = new G4PhysicsFreeVector(NumberofEPoints);
293  theTable2->push_back(thevec2);
294  }
295 
296  for (j=0;j<reducedEnergyGrid;j++)
297  {
298  G4PhysicsFreeVector* thevec = (G4PhysicsFreeVector*) (*theTable1)[j];
299  G4PhysicsFreeVector* thevec2 = (G4PhysicsFreeVector*) (*theTable2)[j];
300  for (i=0;i<NumberofEPoints;i++)
301  {
302  thevec->PutValue(i,betas[i],Q1E[i][j]);
303  thevec2->PutValue(i,betas[i],Q2E[i][j]);
304  }
305  thevec->SetSpline(true);
306  thevec2->SetSpline(true);
307  }
308 
310  {
311  theLorentzTables1->insert(std::make_pair(Zmat,theTable1));
312  theLorentzTables2->insert(std::make_pair(Zmat,theTable2));
313  }
314  else
315  {
317  ed << "Unable to create tables of Lorentz coefficients for " << G4endl;
318  ed << "<Z>= " << Zmat << " in G4PenelopeBremsstrahlungAngular" << G4endl;
319  delete theTable1;
320  delete theTable2;
321  G4Exception("G4PenelopeBremsstrahlungAngular::PrepareInterpolationTables()",
322  "em2005",FatalException,ed);
323  }
324 
325  return;
326 }
327 
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
329 
331  G4double eGamma,
332  G4int,
333  const G4Material* material)
334 {
335  if (!material)
336  {
337  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
338  "em2040",FatalException,"The pointer to G4Material* is NULL");
339  return fLocalDirection;
340  }
341 
342  //Retrieve the effective Z
343  G4double Zmat = 0;
344 
345  if (!theEffectiveZSq)
346  {
347  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
348  "em2040",FatalException,"EffectiveZ table not available");
349  return fLocalDirection;
350  }
351 
352  //found in the table: return it
353  if (theEffectiveZSq->count(material))
354  Zmat = theEffectiveZSq->find(material)->second;
355  else
356  {
357  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
358  "em2040",FatalException,"Material not found in the effectiveZ table");
359  return fLocalDirection;
360  }
361 
362  if (verbosityLevel > 0)
363  {
364  G4cout << "Effective <Z> for material : " << material->GetName() <<
365  " = " << Zmat << G4endl;
366  }
367 
368  G4double ePrimary = dp->GetKineticEnergy();
369 
370  G4double beta = std::sqrt(ePrimary*(ePrimary+2*electron_mass_c2))/
371  (ePrimary+electron_mass_c2);
372  G4double cdt = 0;
373  G4double sinTheta = 0;
374  G4double phi = 0;
375 
376  //Use a pure dipole distribution for energy above 500 keV
377  if (ePrimary > 500*keV)
378  {
379  cdt = 2.0*G4UniformRand() - 1.0;
380  if (G4UniformRand() > 0.75)
381  {
382  if (cdt<0)
383  cdt = -1.0*std::pow(-cdt,1./3.);
384  else
385  cdt = std::pow(cdt,1./3.);
386  }
387  cdt = (cdt+beta)/(1.0+beta*cdt);
388  //Get primary kinematics
389  sinTheta = std::sqrt(1. - cdt*cdt);
390  phi = twopi * G4UniformRand();
391  fLocalDirection.set(sinTheta* std::cos(phi),
392  sinTheta* std::sin(phi),
393  cdt);
394  //rotate
395  fLocalDirection.rotateUz(dp->GetMomentumDirection());
396  //return
397  return fLocalDirection;
398  }
399 
400  if (!(theLorentzTables1->count(Zmat)) || !(theLorentzTables2->count(Zmat)))
401  {
403  ed << "Unable to retrieve Lorentz tables for Z= " << Zmat << G4endl;
404  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
405  "em2006",FatalException,ed);
406  }
407 
408  //retrieve actual tables
409  const G4PhysicsTable* theTable1 = theLorentzTables1->find(Zmat)->second;
410  const G4PhysicsTable* theTable2 = theLorentzTables2->find(Zmat)->second;
411 
412  G4double RK=20.0*eGamma/ePrimary;
413  G4int ik=std::min((G4int) RK,19);
414 
415  G4double P10=0,P11=0,P1=0;
416  G4double P20=0,P21=0,P2=0;
417 
418  //First coefficient
419  const G4PhysicsFreeVector* v1 = (G4PhysicsFreeVector*) (*theTable1)[ik];
420  const G4PhysicsFreeVector* v2 = (G4PhysicsFreeVector*) (*theTable1)[ik+1];
421  P10 = v1->Value(beta);
422  P11 = v2->Value(beta);
423  P1=P10+(RK-(G4double) ik)*(P11-P10);
424 
425  //Second coefficient
426  const G4PhysicsFreeVector* v3 = (G4PhysicsFreeVector*) (*theTable2)[ik];
427  const G4PhysicsFreeVector* v4 = (G4PhysicsFreeVector*) (*theTable2)[ik+1];
428  P20=v3->Value(beta);
429  P21=v4->Value(beta);
430  P2=P20+(RK-(G4double) ik)*(P21-P20);
431 
432  //Sampling from the Lorenz-trasformed dipole distributions
433  P1=std::min(G4Exp(P1)/beta,1.0);
434  G4double betap = std::min(std::max(beta*(1.0+P2/beta),0.0),0.9999);
435 
436  G4double testf=0;
437 
438  if (G4UniformRand() < P1)
439  {
440  do{
441  cdt = 2.0*G4UniformRand()-1.0;
442  testf=2.0*G4UniformRand()-(1.0+cdt*cdt);
443  }while(testf>0);
444  }
445  else
446  {
447  do{
448  cdt = 2.0*G4UniformRand()-1.0;
449  testf=G4UniformRand()-(1.0-cdt*cdt);
450  }while(testf>0);
451  }
452  cdt = (cdt+betap)/(1.0+betap*cdt);
453 
454  //Get primary kinematics
455  sinTheta = std::sqrt(1. - cdt*cdt);
456  phi = twopi * G4UniformRand();
457  fLocalDirection.set(sinTheta* std::cos(phi),
458  sinTheta* std::sin(phi),
459  cdt);
460  //rotate
461  fLocalDirection.rotateUz(dp->GetMomentumDirection());
462  //return
463  return fLocalDirection;
464 
465 }
466 
467 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
468 
470  const G4double ,
471  const G4int )
472 {
473  G4cout << "WARNING: G4PenelopeBremsstrahlungAngular() does NOT support PolarAngle()" << G4endl;
474  G4cout << "Please use the alternative interface SampleDirection()" << G4endl;
475  G4Exception("G4PenelopeBremsstrahlungAngular::PolarAngle()",
476  "em0005",FatalException,"Unsupported interface");
477  return 0;
478 }
479 
480 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
481 
483 {
484  if (!theEffectiveZSq)
485  theEffectiveZSq = new std::map<const G4Material*,G4double>;
486 
487  //Already exists: return it
488  if (theEffectiveZSq->count(material))
489  return theEffectiveZSq->find(material)->second;
490 
491  //Helper for the calculation
492  std::vector<G4double> *StechiometricFactors = new std::vector<G4double>;
493  G4int nElements = material->GetNumberOfElements();
494  const G4ElementVector* elementVector = material->GetElementVector();
495  const G4double* fractionVector = material->GetFractionVector();
496  for (G4int i=0;i<nElements;i++)
497  {
498  G4double fraction = fractionVector[i];
499  G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
500  StechiometricFactors->push_back(fraction/atomicWeigth);
501  }
502  //Find max
503  G4double MaxStechiometricFactor = 0.;
504  for (G4int i=0;i<nElements;i++)
505  {
506  if ((*StechiometricFactors)[i] > MaxStechiometricFactor)
507  MaxStechiometricFactor = (*StechiometricFactors)[i];
508  }
509  //Normalize
510  for (G4int i=0;i<nElements;i++)
511  (*StechiometricFactors)[i] /= MaxStechiometricFactor;
512 
513  G4double sumz2 = 0;
514  G4double sums = 0;
515  for (G4int i=0;i<nElements;i++)
516  {
517  G4double Z = (*elementVector)[i]->GetZ();
518  sumz2 += (*StechiometricFactors)[i]*Z*Z;
519  sums += (*StechiometricFactors)[i];
520  }
521  delete StechiometricFactors;
522 
523  G4double ZBR = std::sqrt(sumz2/sums);
524  theEffectiveZSq->insert(std::make_pair(material,ZBR));
525 
526  return ZBR;
527 }
static const double MeV
Definition: G4SIunits.hh:211
std::map< G4double, G4PhysicsTable * > * theLorentzTables2
static const G4double * P1[nN]
std::vector< G4Element * > G4ElementVector
void PutValue(size_t binNumber, G4double binValue, G4double dataValue)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
static const G4double a1
const G4String & GetName() const
Definition: G4Material.hh:178
void push_back(G4PhysicsVector *)
static const G4double P20[nE]
static const G4double P10[nE]
static const G4double P11[nE]
void PrepareTables(const G4Material *material, G4bool isMaster)
Reserved for Master Model.
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
Samples the direction of the outgoing photon (in global coordinates).
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
void SetSpline(G4bool)
std::map< G4double, G4PhysicsTable * > * theLorentzTables1
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
std::ostream & tab(std::ostream &)
Definition: CCalutils.cc:89
const G4ThreeVector & GetMomentumDirection() const
G4double QQ1[NumberofZPoints][NumberofEPoints][NumberofKPoints]
static const double twopi
Definition: G4SIunits.hh:75
G4double Value(G4double theEnergy, size_t &lastidx) const
void Initialize()
Reserved for Master Model The Initialize() method forces the cleaning of tables.
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double CalculateEffectiveZ(const G4Material *material)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static const double g
Definition: G4SIunits.hh:180
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4double * P2[nN]
G4double QQ2[NumberofZPoints][NumberofEPoints][NumberofKPoints]
static const double mole
Definition: G4SIunits.hh:283
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
G4double PolarAngle(const G4double initial_energy, const G4double final_energy, const G4int Z)
Old interface, backwards compatibility. Will not work in this case it will produce a G4Exception()...
std::map< const G4Material *, G4double > * theEffectiveZSq
const G4double * GetFractionVector() const
Definition: G4Material.hh:194
static const G4double P21[nE]
static const G4double a2