Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeBremsstrahlungAngular Class Reference

#include <G4PenelopeBremsstrahlungAngular.hh>

+ Inheritance diagram for G4PenelopeBremsstrahlungAngular:

Public Member Functions

 G4PenelopeBremsstrahlungAngular ()
 
 ~G4PenelopeBremsstrahlungAngular ()
 
G4double PolarAngle (const G4double initial_energy, const G4double final_energy, const G4int Z)
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double out_energy, G4int Z, const G4Material *mat=0)
 Samples the direction of the outgoing photon (in global coordinates).
 
void SetVerbosityLevel (G4int vl)
 Set/Get Verbosity level.
 
G4int GetVerbosityLevel ()
 
void Initialize ()
 
void PrepareTables (const G4Material *material, G4bool isMaster)
 Reserved for Master Model.
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
 
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)
 
const G4StringGetName () const
 
G4VEmAngularDistributionoperator= (const G4VEmAngularDistribution &right)=delete
 
 G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 55 of file G4PenelopeBremsstrahlungAngular.hh.

Constructor & Destructor Documentation

◆ G4PenelopeBremsstrahlungAngular()

G4PenelopeBremsstrahlungAngular::G4PenelopeBremsstrahlungAngular ( )

Definition at line 60 of file G4PenelopeBremsstrahlungAngular.cc.

60 :
61 G4VEmAngularDistribution("Penelope"), theEffectiveZSq(0),
62 theLorentzTables1(0),theLorentzTables2(0)
63
64{
65 dataRead = false;
66 verbosityLevel = 0;
67}

◆ ~G4PenelopeBremsstrahlungAngular()

G4PenelopeBremsstrahlungAngular::~G4PenelopeBremsstrahlungAngular ( )

Definition at line 71 of file G4PenelopeBremsstrahlungAngular.cc.

72{
73 ClearTables();
74}

Member Function Documentation

◆ GetVerbosityLevel()

G4int G4PenelopeBremsstrahlungAngular::GetVerbosityLevel ( )
inline

Definition at line 77 of file G4PenelopeBremsstrahlungAngular.hh.

77{return verbosityLevel;};

◆ 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{
80 ClearTables();
81}

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

◆ PolarAngle()

G4double G4PenelopeBremsstrahlungAngular::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().

Definition at line 466 of file G4PenelopeBremsstrahlungAngular.cc.

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}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ PrepareTables()

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

Reserved for Master Model.

Definition at line 172 of file G4PenelopeBremsstrahlungAngular.cc.

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}
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
int G4int
Definition: G4Types.hh:85
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)

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

◆ SampleDirection()

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

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

Implements G4VEmAngularDistribution.

Definition at line 327 of file G4PenelopeBremsstrahlungAngular.cc.

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}
#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

Referenced by G4PenelopeBremsstrahlungModel::SampleSecondaries().

◆ SetVerbosityLevel()

void G4PenelopeBremsstrahlungAngular::SetVerbosityLevel ( G4int  vl)
inline

Set/Get Verbosity level.

Definition at line 76 of file G4PenelopeBremsstrahlungAngular.hh.

76{verbosityLevel = vl;};

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