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