Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4PenelopeBremsstrahlungAngular Class Reference

#include <G4PenelopeBremsstrahlungAngular.hh>

Inheritance diagram for G4PenelopeBremsstrahlungAngular:
G4VEmAngularDistribution

Public Member Functions

 G4PenelopeBremsstrahlungAngular ()
 
const G4StringGetName () const
 
G4int GetVerbosityLevel ()
 
void Initialize ()
 
void PrepareTables (const G4Material *material, G4bool isMaster)
 Reserved for Master Model. More...
 
virtual void PrintGeneratorInformation () const
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=nullptr) override
 Samples the direction of the outgoing photon (in global coordinates). More...
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
 
void SetVerbosityLevel (G4int vl)
 Set/Get Verbosity level. More...
 
 ~G4PenelopeBremsstrahlungAngular ()
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Member Functions

G4double CalculateEffectiveZ (const G4Material *material)
 
void ClearTables ()
 
void ReadDataFile ()
 

Private Attributes

G4bool fDataRead
 
std::map< const G4Material *, G4double > * fEffectiveZSq
 
std::map< G4double, G4PhysicsTable * > * fLorentzTables1
 
std::map< G4double, G4PhysicsTable * > * fLorentzTables2
 
G4String fName
 
G4double fQQ1 [fNumberofZPoints][fNumberofEPoints][fNumberofKPoints]
 
G4double fQQ2 [fNumberofZPoints][fNumberofEPoints][fNumberofKPoints]
 
G4int fVerbosityLevel
 

Static Private Attributes

static const G4int fNumberofEPoints =6
 
static const G4int fNumberofKPoints =4
 
static const G4int fNumberofZPoints =6
 

Detailed Description

Definition at line 55 of file G4PenelopeBremsstrahlungAngular.hh.

Constructor & Destructor Documentation

◆ G4PenelopeBremsstrahlungAngular()

G4PenelopeBremsstrahlungAngular::G4PenelopeBremsstrahlungAngular ( )
explicit

Definition at line 60 of file G4PenelopeBremsstrahlungAngular.cc.

60 :
61 G4VEmAngularDistribution("Penelope"), fEffectiveZSq(nullptr),
62 fLorentzTables1(nullptr),fLorentzTables2(nullptr)
63
64{
65 fDataRead = false;
67}
std::map< G4double, G4PhysicsTable * > * fLorentzTables1
std::map< G4double, G4PhysicsTable * > * fLorentzTables2
std::map< const G4Material *, G4double > * fEffectiveZSq
G4VEmAngularDistribution(const G4String &name)

References fDataRead, and fVerbosityLevel.

◆ ~G4PenelopeBremsstrahlungAngular()

G4PenelopeBremsstrahlungAngular::~G4PenelopeBremsstrahlungAngular ( )

Member Function Documentation

◆ CalculateEffectiveZ()

G4double G4PenelopeBremsstrahlungAngular::CalculateEffectiveZ ( const G4Material material)
private

Definition at line 465 of file G4PenelopeBremsstrahlungAngular.cc.

466{
467 if (!fEffectiveZSq)
468 fEffectiveZSq = new std::map<const G4Material*,G4double>;
469
470 //Already exists: return it
471 if (fEffectiveZSq->count(material))
472 return fEffectiveZSq->find(material)->second;
473
474 //Helper for the calculation
475 std::vector<G4double> *StechiometricFactors = new std::vector<G4double>;
476 G4int nElements = material->GetNumberOfElements();
477 const G4ElementVector* elementVector = material->GetElementVector();
478 const G4double* fractionVector = material->GetFractionVector();
479 for (G4int i=0;i<nElements;i++)
480 {
481 G4double fraction = fractionVector[i];
482 G4double atomicWeigth = (*elementVector)[i]->GetA()/(g/mole);
483 StechiometricFactors->push_back(fraction/atomicWeigth);
484 }
485 //Find max
486 G4double MaxStechiometricFactor = 0.;
487 for (G4int i=0;i<nElements;i++)
488 {
489 if ((*StechiometricFactors)[i] > MaxStechiometricFactor)
490 MaxStechiometricFactor = (*StechiometricFactors)[i];
491 }
492 //Normalize
493 for (G4int i=0;i<nElements;i++)
494 (*StechiometricFactors)[i] /= MaxStechiometricFactor;
495
496 G4double sumz2 = 0;
497 G4double sums = 0;
498 for (G4int i=0;i<nElements;i++)
499 {
500 G4double Z = (*elementVector)[i]->GetZ();
501 sumz2 += (*StechiometricFactors)[i]*Z*Z;
502 sums += (*StechiometricFactors)[i];
503 }
504 delete StechiometricFactors;
505
506 G4double ZBR = std::sqrt(sumz2/sums);
507 fEffectiveZSq->insert(std::make_pair(material,ZBR));
508
509 return ZBR;
510}
std::vector< const G4Element * > G4ElementVector
static constexpr double mole
Definition: G4SIunits.hh:279
static constexpr double g
Definition: G4SIunits.hh:168
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
string material
Definition: eplot.py:19

References fEffectiveZSq, g, eplot::material, mole, and Z.

Referenced by PrepareTables().

◆ ClearTables()

void G4PenelopeBremsstrahlungAngular::ClearTables ( )
private

Definition at line 85 of file G4PenelopeBremsstrahlungAngular.cc.

86{
88 {
89 for (auto j=fLorentzTables1->begin(); j != fLorentzTables1->end(); ++j)
90 {
91 G4PhysicsTable* tab = j->second;
92 tab->clearAndDestroy();
93 delete tab;
94 }
95 fLorentzTables1->clear();
96 delete fLorentzTables1;
97 fLorentzTables1 = nullptr;
98 }
99
100 if (fLorentzTables2)
101 {
102 for (auto j=fLorentzTables2->begin(); j != fLorentzTables2->end(); ++j)
103 {
104 G4PhysicsTable* tab = j->second;
105 tab->clearAndDestroy();
106 delete tab;
107 }
108 fLorentzTables2->clear();
109 delete fLorentzTables2;
110 fLorentzTables2 = nullptr;
111 }
112 if (fEffectiveZSq)
113 {
114 delete fEffectiveZSq;
115 fEffectiveZSq = nullptr;
116 }
117}
void clearAndDestroy()

References G4PhysicsTable::clearAndDestroy(), fEffectiveZSq, fLorentzTables1, and fLorentzTables2.

Referenced by Initialize(), and ~G4PenelopeBremsstrahlungAngular().

◆ GetName()

const G4String & G4VEmAngularDistribution::GetName ( ) const
inlineinherited

Definition at line 111 of file G4VEmAngularDistribution.hh.

112{
113 return fName;
114}

References G4VEmAngularDistribution::fName.

◆ GetVerbosityLevel()

G4int G4PenelopeBremsstrahlungAngular::GetVerbosityLevel ( )
inline

Definition at line 69 of file G4PenelopeBremsstrahlungAngular.hh.

69{return fVerbosityLevel;};

References fVerbosityLevel.

◆ Initialize()

void G4PenelopeBremsstrahlungAngular::Initialize ( )

Reserved for Master Model The Initialize() method forces the cleaning of tables

Definition at line 78 of file G4PenelopeBremsstrahlungAngular.cc.

79{
81}

References ClearTables().

Referenced by G4PenelopeBremsstrahlungModel::Initialise(), and G4PenelopeBremsstrahlungModel::InitialiseLocal().

◆ PrepareTables()

void G4PenelopeBremsstrahlungAngular::PrepareTables ( const G4Material material,
G4bool  isMaster 
)

Reserved for Master Model.

Definition at line 174 of file G4PenelopeBremsstrahlungAngular.cc.

175{
176 //Unused at the moment: the G4PenelopeBremsstrahlungAngular is thread-local, so each worker
177 //builds its own version of the tables.
178
179 //Check if data file has already been read
180 if (!fDataRead)
181 {
182 ReadDataFile();
183 if (!fDataRead)
184 G4Exception("G4PenelopeBremsstrahlungAngular::PrepareInterpolationTables()",
185 "em2001",FatalException,"Unable to build interpolation table");
186 }
187
188 if (!fLorentzTables1)
189 fLorentzTables1 = new std::map<G4double,G4PhysicsTable*>;
190 if (!fLorentzTables2)
191 fLorentzTables2 = new std::map<G4double,G4PhysicsTable*>;
192
194
195 const G4int reducedEnergyGrid=21;
196 //Support arrays.
197 G4double betas[fNumberofEPoints]; //betas for interpolation
198 //tables for interpolation
201 //expanded tables for interpolation
202 G4double Q1E[fNumberofEPoints][reducedEnergyGrid];
203 G4double Q2E[fNumberofEPoints][reducedEnergyGrid];
204 G4double pZ[fNumberofZPoints] = {2.0,8.0,13.0,47.0,79.0,92.0};
205
206 G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
207 //Interpolation in Z
208 for (i=0;i<fNumberofEPoints;i++)
209 {
210 for (j=0;j<fNumberofKPoints;j++)
211 {
212 G4PhysicsFreeVector* QQ1vector =
213 new G4PhysicsFreeVector(fNumberofZPoints,/*spline =*/true);
214 G4PhysicsFreeVector* QQ2vector =
215 new G4PhysicsFreeVector(fNumberofZPoints,/*spline =*/true);
216
217 //fill vectors
218 for (k=0;k<fNumberofZPoints;k++)
219 {
220 QQ1vector->PutValues(k,pZ[k],G4Log(fQQ1[k][i][j]));
221 QQ2vector->PutValues(k,pZ[k],fQQ2[k][i][j]);
222 }
223 //Filled: now calculate derivatives
224 QQ1vector->FillSecondDerivatives();
225 QQ2vector->FillSecondDerivatives();
226
227 Q1[i][j]= G4Exp(QQ1vector->Value(Zmat));
228 Q2[i][j]=QQ2vector->Value(Zmat);
229 delete QQ1vector;
230 delete QQ2vector;
231 }
232 }
233 G4double pE[fNumberofEPoints] = {1.0e-03*MeV,5.0e-03*MeV,1.0e-02*MeV,5.0e-02*MeV,
234 1.0e-01*MeV,5.0e-01*MeV};
235 G4double pK[fNumberofKPoints] = {0.0,0.6,0.8,0.95};
236 G4double ppK[reducedEnergyGrid];
237
238 for(i=0;i<reducedEnergyGrid;i++)
239 ppK[i]=((G4double) i) * 0.05;
240
241
242 for(i=0;i<fNumberofEPoints;i++)
243 betas[i]=std::sqrt(pE[i]*(pE[i]+2*electron_mass_c2))/(pE[i]+electron_mass_c2);
244
245
246 for (i=0;i<fNumberofEPoints;i++)
247 {
248 for (j=0;j<fNumberofKPoints;j++)
249 Q1[i][j]=Q1[i][j]/Zmat;
250 }
251
252 //Expanded table of distribution parameters
253 for (i=0;i<fNumberofEPoints;i++)
254 {
257
258 for (j=0;j<fNumberofKPoints;j++)
259 {
260 Q1vector->PutValues(j,pK[j],G4Log(Q1[i][j])); //logarithmic
261 Q2vector->PutValues(j,pK[j],Q2[i][j]);
262 }
263
264 for (j=0;j<reducedEnergyGrid;j++)
265 {
266 Q1E[i][j]=Q1vector->Value(ppK[j]);
267 Q2E[i][j]=Q2vector->Value(ppK[j]);
268 }
269 delete Q1vector;
270 delete Q2vector;
271 }
272 //
273 //TABLES to be stored
274 //
275 G4PhysicsTable* theTable1 = new G4PhysicsTable();
276 G4PhysicsTable* theTable2 = new G4PhysicsTable();
277 // the table will contain reducedEnergyGrid G4PhysicsFreeVectors with different
278 // values of k,
279 // Each of the G4PhysicsFreeVectors has a profile of
280 // y vs. E
281 //
282 //reserve space of the vectors.
283 for (j=0;j<reducedEnergyGrid;j++)
284 {
285 G4PhysicsFreeVector* thevec = new G4PhysicsFreeVector(fNumberofEPoints,/*spline=*/true);
286 theTable1->push_back(thevec);
287 G4PhysicsFreeVector* thevec2 = new G4PhysicsFreeVector(fNumberofEPoints,/*spline=*/true);
288 theTable2->push_back(thevec2);
289 }
290
291 for (j=0;j<reducedEnergyGrid;j++)
292 {
293 G4PhysicsFreeVector* thevec = (G4PhysicsFreeVector*) (*theTable1)[j];
294 G4PhysicsFreeVector* thevec2 = (G4PhysicsFreeVector*) (*theTable2)[j];
295 for (i=0;i<fNumberofEPoints;i++)
296 {
297 thevec->PutValues(i,betas[i],Q1E[i][j]);
298 thevec2->PutValues(i,betas[i],Q2E[i][j]);
299 }
300 //Vectors are filled: calculate derivatives
301 thevec->FillSecondDerivatives();
302 thevec2->FillSecondDerivatives();
303 }
304
306 {
307 fLorentzTables1->insert(std::make_pair(Zmat,theTable1));
308 fLorentzTables2->insert(std::make_pair(Zmat,theTable2));
309 }
310 else
311 {
313 ed << "Unable to create tables of Lorentz coefficients for " << G4endl;
314 ed << "<Z>= " << Zmat << " in G4PenelopeBremsstrahlungAngular" << G4endl;
315 delete theTable1;
316 delete theTable2;
317 G4Exception("G4PenelopeBremsstrahlungAngular::PrepareInterpolationTables()",
318 "em2005",FatalException,ed);
319 }
320
321 return;
322}
@ 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
static constexpr double MeV
Definition: G4SIunits.hh:200
#define G4endl
Definition: G4ios.hh:57
G4double CalculateEffectiveZ(const G4Material *material)
G4double fQQ1[fNumberofZPoints][fNumberofEPoints][fNumberofKPoints]
G4double fQQ2[fNumberofZPoints][fNumberofEPoints][fNumberofKPoints]
void PutValues(const std::size_t index, const G4double energy, const G4double value)
void push_back(G4PhysicsVector *)
G4double Value(const G4double energy, std::size_t &lastidx) const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
float electron_mass_c2
Definition: hepunit.py:273

References CalculateEffectiveZ(), source.hepunit::electron_mass_c2, FatalException, fDataRead, G4PhysicsVector::FillSecondDerivatives(), fLorentzTables1, fLorentzTables2, fNumberofEPoints, fNumberofKPoints, fNumberofZPoints, fQQ1, fQQ2, G4endl, G4Exception(), G4Exp(), G4Log(), eplot::material, MeV, G4PhysicsTable::push_back(), G4PhysicsFreeVector::PutValues(), ReadDataFile(), and G4PhysicsVector::Value().

Referenced by G4PenelopeBremsstrahlungModel::Initialise(), and G4PenelopeBremsstrahlungModel::InitialiseLocal().

◆ PrintGeneratorInformation()

void G4VEmAngularDistribution::PrintGeneratorInformation ( ) const
virtualinherited

◆ ReadDataFile()

void G4PenelopeBremsstrahlungAngular::ReadDataFile ( void  )
private

Definition at line 121 of file G4PenelopeBremsstrahlungAngular.cc.

122{
123 //Read information from DataBase file
124 char* path = std::getenv("G4LEDATA");
125 if (!path)
126 {
127 G4String excep =
128 "G4PenelopeBremsstrahlungAngular - G4LEDATA environment variable not set!";
129 G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
130 "em0006",FatalException,excep);
131 return;
132 }
133 G4String pathString(path);
134 G4String pathFile = pathString + "/penelope/bremsstrahlung/pdbrang.p08";
135 std::ifstream file(pathFile);
136
137 if (!file.is_open())
138 {
139 G4String excep = "G4PenelopeBremsstrahlungAngular - data file " + pathFile + " not found!";
140 G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
141 "em0003",FatalException,excep);
142 return;
143 }
144 G4int i=0,j=0,k=0; // i=index for Z, j=index for E, k=index for K
145
146 for (k=0;k<fNumberofKPoints;k++)
147 for (i=0;i<fNumberofZPoints;i++)
148 for (j=0;j<fNumberofEPoints;j++)
149 {
150 G4double a1,a2;
151 G4int ik1,iz1,ie1;
152 G4double zr,er,kr;
153 file >> iz1 >> ie1 >> ik1 >> zr >> er >> kr >> a1 >> a2;
154 //check the data are correct
155 if ((iz1-1 == i) && (ik1-1 == k) && (ie1-1 == j))
156 {
157 fQQ1[i][j][k]=a1;
158 fQQ2[i][j][k]=a2;
159 }
160 else
161 {
163 ed << "Corrupted data file " << pathFile << "?" << G4endl;
164 G4Exception("G4PenelopeBremsstrahlungAngular::ReadDataFile()",
165 "em0005",FatalException,ed);
166 }
167 }
168 file.close();
169 fDataRead = true;
170}

References FatalException, fDataRead, geant4_check_module_cycles::file, fNumberofEPoints, fNumberofKPoints, fNumberofZPoints, fQQ1, fQQ2, G4endl, and G4Exception().

Referenced by PrepareTables().

◆ SampleDirection()

G4ThreeVector & G4PenelopeBremsstrahlungAngular::SampleDirection ( const G4DynamicParticle dp,
G4double  out_energy,
G4int  Z,
const G4Material mat = nullptr 
)
overridevirtual

Samples the direction of the outgoing photon (in global coordinates).

Implements G4VEmAngularDistribution.

Definition at line 326 of file G4PenelopeBremsstrahlungAngular.cc.

330{
331 if (!material)
332 {
333 G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
334 "em2040",FatalException,"The pointer to G4Material* is nullptr");
335 return fLocalDirection;
336 }
337
338 //Retrieve the effective Z
339 G4double Zmat = 0;
340
341 if (!fEffectiveZSq)
342 {
343 G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
344 "em2040",FatalException,"EffectiveZ table not available");
345 return fLocalDirection;
346 }
347
348 //found in the table: return it
349 if (fEffectiveZSq->count(material))
350 Zmat = fEffectiveZSq->find(material)->second;
351 else
352 {
353 G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
354 "em2040",FatalException,"Material not found in the effectiveZ table");
355 return fLocalDirection;
356 }
357
358 if (fVerbosityLevel > 0)
359 {
360 G4cout << "Effective <Z> for material : " << material->GetName() <<
361 " = " << Zmat << G4endl;
362 }
363
364 G4double ePrimary = dp->GetKineticEnergy();
365
366 G4double beta = std::sqrt(ePrimary*(ePrimary+2*electron_mass_c2))/
367 (ePrimary+electron_mass_c2);
368 G4double cdt = 0;
369 G4double sinTheta = 0;
370 G4double phi = 0;
371
372 //Use a pure dipole distribution for energy above 500 keV
373 if (ePrimary > 500*keV)
374 {
375 cdt = 2.0*G4UniformRand() - 1.0;
376 if (G4UniformRand() > 0.75)
377 {
378 if (cdt<0)
379 cdt = -1.0*std::pow(-cdt,1./3.);
380 else
381 cdt = std::pow(cdt,1./3.);
382 }
383 cdt = (cdt+beta)/(1.0+beta*cdt);
384 //Get primary kinematics
385 sinTheta = std::sqrt(1. - cdt*cdt);
386 phi = twopi * G4UniformRand();
387 fLocalDirection.set(sinTheta* std::cos(phi),
388 sinTheta* std::sin(phi),
389 cdt);
390 //rotate
392 //return
393 return fLocalDirection;
394 }
395
396 if (!(fLorentzTables1->count(Zmat)) || !(fLorentzTables2->count(Zmat)))
397 {
399 ed << "Unable to retrieve Lorentz tables for Z= " << Zmat << G4endl;
400 G4Exception("G4PenelopeBremsstrahlungAngular::SampleDirection()",
401 "em2006",FatalException,ed);
402 }
403
404 //retrieve actual tables
405 const G4PhysicsTable* theTable1 = fLorentzTables1->find(Zmat)->second;
406 const G4PhysicsTable* theTable2 = fLorentzTables2->find(Zmat)->second;
407
408 G4double RK=20.0*eGamma/ePrimary;
409 G4int ik=std::min((G4int) RK,19);
410
411 G4double P10=0,P11=0,P1=0;
412 G4double P20=0,P21=0,P2=0;
413
414 //First coefficient
415 const G4PhysicsFreeVector* v1 = (G4PhysicsFreeVector*) (*theTable1)[ik];
416 const G4PhysicsFreeVector* v2 = (G4PhysicsFreeVector*) (*theTable1)[ik+1];
417 P10 = v1->Value(beta);
418 P11 = v2->Value(beta);
419 P1=P10+(RK-(G4double) ik)*(P11-P10);
420
421 //Second coefficient
422 const G4PhysicsFreeVector* v3 = (G4PhysicsFreeVector*) (*theTable2)[ik];
423 const G4PhysicsFreeVector* v4 = (G4PhysicsFreeVector*) (*theTable2)[ik+1];
424 P20=v3->Value(beta);
425 P21=v4->Value(beta);
426 P2=P20+(RK-(G4double) ik)*(P21-P20);
427
428 //Sampling from the Lorenz-trasformed dipole distributions
429 P1=std::min(G4Exp(P1)/beta,1.0);
430 G4double betap = std::min(std::max(beta*(1.0+P2/beta),0.0),0.9999);
431
432 G4double testf=0;
433
434 if (G4UniformRand() < P1)
435 {
436 do{
437 cdt = 2.0*G4UniformRand()-1.0;
438 testf=2.0*G4UniformRand()-(1.0+cdt*cdt);
439 }while(testf>0);
440 }
441 else
442 {
443 do{
444 cdt = 2.0*G4UniformRand()-1.0;
445 testf=G4UniformRand()-(1.0-cdt*cdt);
446 }while(testf>0);
447 }
448 cdt = (cdt+betap)/(1.0+betap*cdt);
449
450 //Get primary kinematics
451 sinTheta = std::sqrt(1. - cdt*cdt);
452 phi = twopi * G4UniformRand();
453 fLocalDirection.set(sinTheta* std::cos(phi),
454 sinTheta* std::sin(phi),
455 cdt);
456 //rotate
458 //return
459 return fLocalDirection;
460
461}
static const G4double P21[nE]
static const G4double P10[nE]
static const G4double * P1[nN]
static const G4double P11[nE]
static const G4double P20[nE]
static const G4double * P2[nN]
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double keV
Definition: G4SIunits.hh:202
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
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, source.hepunit::electron_mass_c2, FatalException, fEffectiveZSq, G4VEmAngularDistribution::fLocalDirection, fLorentzTables1, fLorentzTables2, fVerbosityLevel, G4cout, G4endl, G4Exception(), G4Exp(), G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), keV, eplot::material, G4INCL::Math::max(), G4INCL::Math::min(), P1, P10, P11, P2, P20, P21, CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), twopi, and G4PhysicsVector::Value().

Referenced by G4PenelopeBremsstrahlungModel::SampleSecondaries().

◆ SampleDirectionForShell()

G4ThreeVector & G4VEmAngularDistribution::SampleDirectionForShell ( const G4DynamicParticle dp,
G4double  finalTotalEnergy,
G4int  Z,
G4int  shellID,
const G4Material mat 
)
virtualinherited

◆ SamplePairDirections()

void G4VEmAngularDistribution::SamplePairDirections ( const G4DynamicParticle dp,
G4double  elecKinEnergy,
G4double  posiKinEnergy,
G4ThreeVector dirElectron,
G4ThreeVector dirPositron,
G4int  Z = 0,
const G4Material mat = nullptr 
)
virtualinherited

◆ SetVerbosityLevel()

void G4PenelopeBremsstrahlungAngular::SetVerbosityLevel ( G4int  vl)
inline

Set/Get Verbosity level.

Definition at line 68 of file G4PenelopeBremsstrahlungAngular.hh.

68{fVerbosityLevel = vl;};

References fVerbosityLevel.

Field Documentation

◆ fDataRead

G4bool G4PenelopeBremsstrahlungAngular::fDataRead
private

◆ fEffectiveZSq

std::map<const G4Material*,G4double>* G4PenelopeBremsstrahlungAngular::fEffectiveZSq
private

◆ fLocalDirection

G4ThreeVector G4VEmAngularDistribution::fLocalDirection
protectedinherited

◆ fLorentzTables1

std::map<G4double,G4PhysicsTable*>* G4PenelopeBremsstrahlungAngular::fLorentzTables1
private

Definition at line 87 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by ClearTables(), PrepareTables(), and SampleDirection().

◆ fLorentzTables2

std::map<G4double,G4PhysicsTable*>* G4PenelopeBremsstrahlungAngular::fLorentzTables2
private

Definition at line 88 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by ClearTables(), PrepareTables(), and SampleDirection().

◆ fName

G4String G4VEmAngularDistribution::fName
privateinherited

Definition at line 108 of file G4VEmAngularDistribution.hh.

Referenced by G4VEmAngularDistribution::GetName().

◆ fNumberofEPoints

const G4int G4PenelopeBremsstrahlungAngular::fNumberofEPoints =6
staticprivate

Definition at line 91 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by PrepareTables(), and ReadDataFile().

◆ fNumberofKPoints

const G4int G4PenelopeBremsstrahlungAngular::fNumberofKPoints =4
staticprivate

Definition at line 92 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by PrepareTables(), and ReadDataFile().

◆ fNumberofZPoints

const G4int G4PenelopeBremsstrahlungAngular::fNumberofZPoints =6
staticprivate

Definition at line 90 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by PrepareTables(), and ReadDataFile().

◆ fPolarisation

G4bool G4VEmAngularDistribution::fPolarisation
protectedinherited

◆ fQQ1

G4double G4PenelopeBremsstrahlungAngular::fQQ1[fNumberofZPoints][fNumberofEPoints][fNumberofKPoints]
private

Definition at line 94 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by PrepareTables(), and ReadDataFile().

◆ fQQ2

G4double G4PenelopeBremsstrahlungAngular::fQQ2[fNumberofZPoints][fNumberofEPoints][fNumberofKPoints]
private

Definition at line 95 of file G4PenelopeBremsstrahlungAngular.hh.

Referenced by PrepareTables(), and ReadDataFile().

◆ fVerbosityLevel

G4int G4PenelopeBremsstrahlungAngular::fVerbosityLevel
private

The documentation for this class was generated from the following files: