Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4PenelopeBremsstrahlungAngular.cc
이 파일의 문서화 페이지로 가기
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 99415 2016-09-21 09:05:43Z 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 {
89  {
90  for (auto j = theLorentzTables1->begin(); j != theLorentzTables1->end(); j++)
91  {
92  G4PhysicsTable* tab = j->second;
93  //tab->clearAndDestroy();
94  delete tab;
95  }
96  delete theLorentzTables1;
97  theLorentzTables1 = nullptr;
98  }
99 
100  if (theLorentzTables2)
101  {
102  for (auto j=theLorentzTables2->begin(); j != theLorentzTables2->end(); j++)
103  {
104  G4PhysicsTable* tab = j->second;
105  //tab->clearAndDestroy();
106  delete tab;
107  }
108  delete theLorentzTables2;
109  theLorentzTables2 = nullptr;
110  }
111  if (theEffectiveZSq)
112  {
113  delete theEffectiveZSq;
114  theEffectiveZSq = nullptr;
115  }
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119 
121 {
122  //Read information from DataBase file
123  char* path = getenv("G4LEDATA");
124  if (!path)
125  {
126  G4String excep =
127  "G4PenelopeBremsstrahlungAngular - G4LEDATA environment variable not set!";
128  G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
129  "em0006",FatalException,excep);
130  return;
131  }
132  G4String pathString(path);
133  G4String pathFile = pathString + "/penelope/bremsstrahlung/pdbrang.p08";
134  std::ifstream file(pathFile);
135 
136  if (!file.is_open())
137  {
138  G4String excep = "G4PenelopeBremsstrahlungAngular - data file " + pathFile + " not found!";
139  G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
140  "em0003",FatalException,excep);
141  return;
142  }
143  G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
144 
145  for (k=0;k<NumberofKPoints;k++)
146  for (i=0;i<NumberofZPoints;i++)
147  for (j=0;j<NumberofEPoints;j++)
148  {
149  G4double a1,a2;
150  G4int ik1,iz1,ie1;
151  G4double zr,er,kr;
152  file >> iz1 >> ie1 >> ik1 >> zr >> er >> kr >> a1 >> a2;
153  //check the data are correct
154  if ((iz1-1 == i) && (ik1-1 == k) && (ie1-1 == j))
155  {
156  QQ1[i][j][k]=a1;
157  QQ2[i][j][k]=a2;
158  }
159  else
160  {
162  ed << "Corrupted data file " << pathFile << "?" << G4endl;
163  G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
164  "em0005",FatalException,ed);
165  }
166  }
167  file.close();
168  dataRead = true;
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
174 {
175  //Unused at the moment: the G4PenelopeBremsstrahlungAngular is thread-local, so each worker
176  //builds its own version of the tables.
177  /*
178  if (!isMaster)
179  //Should not be here!
180  G4Exception("G4PenelopeBremsstrahlungAngular::PrepareTables()",
181  "em0100",FatalException,"Worker thread in this method");
182  */
183 
184  //Check if data file has already been read
185  if (!dataRead)
186  {
187  ReadDataFile();
188  if (!dataRead)
189  G4Exception("G4PenelopeBremsstrahlungAngular::PrepareInterpolationTables()",
190  "em2001",FatalException,"Unable to build interpolation table");
191  }
192 
193  if (!theLorentzTables1)
194  theLorentzTables1 = new std::map<G4double,G4PhysicsTable*>;
195  if (!theLorentzTables2)
196  theLorentzTables2 = new std::map<G4double,G4PhysicsTable*>;
197 
198  G4double Zmat = CalculateEffectiveZ(material);
199 
200  const G4int reducedEnergyGrid=21;
201  //Support arrays.
202  G4double betas[NumberofEPoints]; //betas for interpolation
203  //tables for interpolation
206  //expanded tables for interpolation
207  G4double Q1E[NumberofEPoints][reducedEnergyGrid];
208  G4double Q2E[NumberofEPoints][reducedEnergyGrid];
209  G4double pZ[NumberofZPoints] = {2.0,8.0,13.0,47.0,79.0,92.0};
210 
211  G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
212  //Interpolation in Z
213  for (i=0;i<NumberofEPoints;i++)
214  {
215  for (j=0;j<NumberofKPoints;j++)
216  {
219 
220  //fill vectors
221  for (k=0;k<NumberofZPoints;k++)
222  {
223  QQ1vector->PutValue(k,pZ[k],std::log(QQ1[k][i][j]));
224  QQ2vector->PutValue(k,pZ[k],QQ2[k][i][j]);
225  }
226 
227  QQ1vector->SetSpline(true);
228  QQ2vector->SetSpline(true);
229 
230  Q1[i][j]= G4Exp(QQ1vector->Value(Zmat));
231  Q2[i][j]=QQ2vector->Value(Zmat);
232  delete QQ1vector;
233  delete QQ2vector;
234  }
235  }
236  G4double pE[NumberofEPoints] = {1.0e-03*MeV,5.0e-03*MeV,1.0e-02*MeV,5.0e-02*MeV,
237  1.0e-01*MeV,5.0e-01*MeV};
238  G4double pK[NumberofKPoints] = {0.0,0.6,0.8,0.95};
239  G4double ppK[reducedEnergyGrid];
240 
241  for(i=0;i<reducedEnergyGrid;i++)
242  ppK[i]=((G4double) i) * 0.05;
243 
244 
245  for(i=0;i<NumberofEPoints;i++)
246  betas[i]=std::sqrt(pE[i]*(pE[i]+2*electron_mass_c2))/(pE[i]+electron_mass_c2);
247 
248 
249  for (i=0;i<NumberofEPoints;i++)
250  {
251  for (j=0;j<NumberofKPoints;j++)
252  Q1[i][j]=Q1[i][j]/Zmat;
253  }
254 
255  //Expanded table of distribution parameters
256  for (i=0;i<NumberofEPoints;i++)
257  {
260 
261  for (j=0;j<NumberofKPoints;j++)
262  {
263  Q1vector->PutValue(j,pK[j],std::log(Q1[i][j])); //logarithmic
264  Q2vector->PutValue(j,pK[j],Q2[i][j]);
265  }
266 
267  for (j=0;j<reducedEnergyGrid;j++)
268  {
269  Q1E[i][j]=Q1vector->Value(ppK[j]);
270  Q2E[i][j]=Q2vector->Value(ppK[j]);
271  }
272  delete Q1vector;
273  delete Q2vector;
274  }
275  //
276  //TABLES to be stored
277  //
278  G4PhysicsTable* theTable1 = new G4PhysicsTable();
279  G4PhysicsTable* theTable2 = new G4PhysicsTable();
280  // the table will contain reducedEnergyGrid G4PhysicsFreeVectors with different
281  // values of k,
282  // Each of the G4PhysicsFreeVectors has a profile of
283  // y vs. E
284  //
285  //reserve space of the vectors.
286  for (j=0;j<reducedEnergyGrid;j++)
287  {
288  G4PhysicsFreeVector* thevec = new G4PhysicsFreeVector(NumberofEPoints);
289  theTable1->push_back(thevec);
290  G4PhysicsFreeVector* thevec2 = new G4PhysicsFreeVector(NumberofEPoints);
291  theTable2->push_back(thevec2);
292  }
293 
294  for (j=0;j<reducedEnergyGrid;j++)
295  {
296  G4PhysicsFreeVector* thevec = (G4PhysicsFreeVector*) (*theTable1)[j];
297  G4PhysicsFreeVector* thevec2 = (G4PhysicsFreeVector*) (*theTable2)[j];
298  for (i=0;i<NumberofEPoints;i++)
299  {
300  thevec->PutValue(i,betas[i],Q1E[i][j]);
301  thevec2->PutValue(i,betas[i],Q2E[i][j]);
302  }
303  thevec->SetSpline(true);
304  thevec2->SetSpline(true);
305  }
306 
308  {
309  theLorentzTables1->insert(std::make_pair(Zmat,theTable1));
310  theLorentzTables2->insert(std::make_pair(Zmat,theTable2));
311  }
312  else
313  {
315  ed << "Unable to create tables of Lorentz coefficients for " << G4endl;
316  ed << "<Z>= " << Zmat << " in G4PenelopeBremsstrahlungAngular" << G4endl;
317  delete theTable1;
318  delete theTable2;
319  G4Exception("G4PenelopeBremsstrahlungAngular::PrepareInterpolationTables()",
320  "em2005",FatalException,ed);
321  }
322 
323  return;
324 }
325 
326 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
327 
329  G4double eGamma,
330  G4int,
331  const G4Material* material)
332 {
333  if (!material)
334  {
335  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
336  "em2040",FatalException,"The pointer to G4Material* is nullptr");
337  return fLocalDirection;
338  }
339 
340  //Retrieve the effective Z
341  G4double Zmat = 0;
342 
343  if (!theEffectiveZSq)
344  {
345  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
346  "em2040",FatalException,"EffectiveZ table not available");
347  return fLocalDirection;
348  }
349 
350  //found in the table: return it
351  if (theEffectiveZSq->count(material))
352  Zmat = theEffectiveZSq->find(material)->second;
353  else
354  {
355  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
356  "em2040",FatalException,"Material not found in the effectiveZ table");
357  return fLocalDirection;
358  }
359 
360  if (verbosityLevel > 0)
361  {
362  G4cout << "Effective <Z> for material : " << material->GetName() <<
363  " = " << Zmat << G4endl;
364  }
365 
366  G4double ePrimary = dp->GetKineticEnergy();
367 
368  G4double beta = std::sqrt(ePrimary*(ePrimary+2*electron_mass_c2))/
369  (ePrimary+electron_mass_c2);
370  G4double cdt = 0;
371  G4double sinTheta = 0;
372  G4double phi = 0;
373 
374  //Use a pure dipole distribution for energy above 500 keV
375  if (ePrimary > 500*keV)
376  {
377  cdt = 2.0*G4UniformRand() - 1.0;
378  if (G4UniformRand() > 0.75)
379  {
380  if (cdt<0)
381  cdt = -1.0*std::pow(-cdt,1./3.);
382  else
383  cdt = std::pow(cdt,1./3.);
384  }
385  cdt = (cdt+beta)/(1.0+beta*cdt);
386  //Get primary kinematics
387  sinTheta = std::sqrt(1. - cdt*cdt);
388  phi = twopi * G4UniformRand();
389  fLocalDirection.set(sinTheta* std::cos(phi),
390  sinTheta* std::sin(phi),
391  cdt);
392  //rotate
394  //return
395  return fLocalDirection;
396  }
397 
398  if (!(theLorentzTables1->count(Zmat)) || !(theLorentzTables2->count(Zmat)))
399  {
401  ed << "Unable to retrieve Lorentz tables for Z= " << Zmat << G4endl;
402  G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
403  "em2006",FatalException,ed);
404  }
405 
406  //retrieve actual tables
407  const G4PhysicsTable* theTable1 = theLorentzTables1->find(Zmat)->second;
408  const G4PhysicsTable* theTable2 = theLorentzTables2->find(Zmat)->second;
409 
410  G4double RK=20.0*eGamma/ePrimary;
411  G4int ik=std::min((G4int) RK,19);
412 
413  G4double P10=0,P11=0,P1=0;
414  G4double P20=0,P21=0,P2=0;
415 
416  //First coefficient
417  const G4PhysicsFreeVector* v1 = (G4PhysicsFreeVector*) (*theTable1)[ik];
418  const G4PhysicsFreeVector* v2 = (G4PhysicsFreeVector*) (*theTable1)[ik+1];
419  P10 = v1->Value(beta);
420  P11 = v2->Value(beta);
421  P1=P10+(RK-(G4double) ik)*(P11-P10);
422 
423  //Second coefficient
424  const G4PhysicsFreeVector* v3 = (G4PhysicsFreeVector*) (*theTable2)[ik];
425  const G4PhysicsFreeVector* v4 = (G4PhysicsFreeVector*) (*theTable2)[ik+1];
426  P20=v3->Value(beta);
427  P21=v4->Value(beta);
428  P2=P20+(RK-(G4double) ik)*(P21-P20);
429 
430  //Sampling from the Lorenz-trasformed dipole distributions
431  P1=std::min(G4Exp(P1)/beta,1.0);
432  G4double betap = std::min(std::max(beta*(1.0+P2/beta),0.0),0.9999);
433 
434  G4double testf=0;
435 
436  if (G4UniformRand() < P1)
437  {
438  do{
439  cdt = 2.0*G4UniformRand()-1.0;
440  testf=2.0*G4UniformRand()-(1.0+cdt*cdt);
441  }while(testf>0);
442  }
443  else
444  {
445  do{
446  cdt = 2.0*G4UniformRand()-1.0;
447  testf=G4UniformRand()-(1.0-cdt*cdt);
448  }while(testf>0);
449  }
450  cdt = (cdt+betap)/(1.0+betap*cdt);
451 
452  //Get primary kinematics
453  sinTheta = std::sqrt(1. - cdt*cdt);
454  phi = twopi * G4UniformRand();
455  fLocalDirection.set(sinTheta* std::cos(phi),
456  sinTheta* std::sin(phi),
457  cdt);
458  //rotate
460  //return
461  return fLocalDirection;
462 
463 }
464 
465 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
466 
468  const G4double ,
469  const G4int )
470 {
471  G4cout << "WARNING: G4PenelopeBremsstrahlungAngular() does NOT support PolarAngle()" << G4endl;
472  G4cout << "Please use the alternative interface SampleDirection()" << G4endl;
473  G4Exception("G4PenelopeBremsstrahlungAngular::PolarAngle()",
474  "em0005",FatalException,"Unsupported interface");
475  return 0;
476 }
477 
478 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
479 
481 {
482  if (!theEffectiveZSq)
483  theEffectiveZSq = new std::map<const G4Material*,G4double>;
484 
485  //Already exists: return it
486  if (theEffectiveZSq->count(material))
487  return theEffectiveZSq->find(material)->second;
488 
489  //Helper for the calculation
490  std::vector<G4double> *StechiometricFactors = new std::vector<G4double>;
491  G4int nElements = material->GetNumberOfElements();
492  const G4ElementVector* elementVector = material->GetElementVector();
493  const G4double* fractionVector = material->GetFractionVector();
494  for (G4int i=0;i<nElements;i++)
495  {
496  G4double fraction = fractionVector[i];
497  G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
498  StechiometricFactors->push_back(fraction/atomicWeigth);
499  }
500  //Find max
501  G4double MaxStechiometricFactor = 0.;
502  for (G4int i=0;i<nElements;i++)
503  {
504  if ((*StechiometricFactors)[i] > MaxStechiometricFactor)
505  MaxStechiometricFactor = (*StechiometricFactors)[i];
506  }
507  //Normalize
508  for (G4int i=0;i<nElements;i++)
509  (*StechiometricFactors)[i] /= MaxStechiometricFactor;
510 
511  G4double sumz2 = 0;
512  G4double sums = 0;
513  for (G4int i=0;i<nElements;i++)
514  {
515  G4double Z = (*elementVector)[i]->GetZ();
516  sumz2 += (*StechiometricFactors)[i]*Z*Z;
517  sums += (*StechiometricFactors)[i];
518  }
519  delete StechiometricFactors;
520 
521  G4double ZBR = std::sqrt(sumz2/sums);
522  theEffectiveZSq->insert(std::make_pair(material,ZBR));
523 
524  return ZBR;
525 }
void set(double x, double y, double z)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
T max(const T t1, const T t2)
brief Return the largest of the two arguments
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).
G4double CalculateEffectiveZ(const G4Material *material)
const G4double * GetFractionVector() const
Definition: G4Material.hh:195
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
static constexpr double MeV
Definition: G4SIunits.hh:214
std::map< const G4Material *, G4double > * theEffectiveZSq
static constexpr double keV
Definition: G4SIunits.hh:216
#define G4endl
Definition: G4ios.hh:61
const G4ThreeVector & GetMomentumDirection() const
std::map< G4double, G4PhysicsTable * > * theLorentzTables1
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:191
static const G4double * P2[nN]
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:38
Double_t beta
G4double Value(G4double theEnergy, size_t &lastidx) const
G4double QQ2[NumberofZPoints][NumberofEPoints][NumberofKPoints]
const G4String & GetName() const
Definition: G4Material.hh:179
Float_t Z
static constexpr double g
Definition: G4SIunits.hh:183
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static constexpr double electron_mass_c2
static const G4double P11[nE]
#define G4UniformRand()
Definition: Randomize.hh:53
static constexpr double twopi
Definition: G4SIunits.hh:76
void SetSpline(G4bool)
static const G4double P10[nE]
G4double QQ1[NumberofZPoints][NumberofEPoints][NumberofKPoints]
std::ostream & tab(std::ostream &)
Definition: CCalutils.cc:89
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
void push_back(G4PhysicsVector *)
std::vector< G4Element * > G4ElementVector
int G4int
Definition: G4Types.hh:78
std::map< G4double, G4PhysicsTable * > * theLorentzTables2
TFile * file
static const G4double * P1[nN]
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4double PolarAngle(const G4double initial_energy, const G4double final_energy, const G4int Z)
static const G4double P20[nE]
static const G4double P21[nE]
void PutValue(size_t index, G4double energy, G4double dataValue)
static constexpr double mole
Definition: G4SIunits.hh:286
size_t GetNumberOfElements() const
Definition: G4Material.hh:187
T min(const T t1, const T t2)
brief Return the smallest of the two arguments