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

#include <G4PhotoElectricAngularGeneratorPolarized.hh>

+ Inheritance diagram for G4PhotoElectricAngularGeneratorPolarized:

Public Member Functions

 G4PhotoElectricAngularGeneratorPolarized ()
 
 ~G4PhotoElectricAngularGeneratorPolarized ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double eKinEnergy, G4int shellId, const G4Material *mat=0)
 
void PrintGeneratorInformation () const
 
- 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
 

Protected Member Functions

G4ThreeVector PerpendicularVector (const G4ThreeVector &a) const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 52 of file G4PhotoElectricAngularGeneratorPolarized.hh.

Constructor & Destructor Documentation

◆ G4PhotoElectricAngularGeneratorPolarized()

G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized ( )

Definition at line 71 of file G4PhotoElectricAngularGeneratorPolarized.cc.

72 :G4VEmAngularDistribution("AngularGenSauterGavrilaPolarized")
73{
74 const G4int arrayDim = 980;
75
76 //minimum electron beta parameter allowed
77 betaArray[0] = 0.02;
78 //beta step
79 betaArray[1] = 0.001;
80 //maximum index array for a and c tables
81 betaArray[2] = arrayDim - 1;
82
83 // read Majorant Surface Parameters. This are required in order to generate Gavrila angular photoelectron distribution
84 for(G4int level = 0; level < 2; level++){
85
86 char nameChar0[100] = "ftab0.dat"; // K-shell Majorant Surface Parameters
87 char nameChar1[100] = "ftab1.dat"; // L-shell Majorant Surface Parameters
88
89 G4String filename;
90 if(level == 0) filename = nameChar0;
91 if(level == 1) filename = nameChar1;
92
93 char* path = std::getenv("G4LEDATA");
94 if (!path)
95 {
96 G4String excep = "G4EMDataSet - G4LEDATA environment variable not set";
97 G4Exception("G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized",
98 "em0006",FatalException,"G4LEDATA environment variable not set");
99 return;
100 }
101
102 G4String pathString(path);
103 G4String dirFile = pathString + "/photoelectric_angular/" + filename;
104 std::ifstream infile(dirFile);
105 if (!infile.is_open())
106 {
107 G4String excep = "data file: " + dirFile + " not found";
108 G4Exception("G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized",
109 "em0003",FatalException,excep);
110 return;
111 }
112
113 // Read parameters into tables. The parameters are function of incident electron energy and shell level
114 G4float aRead=0,cRead=0, beta=0;
115 for(G4int i=0 ; i<arrayDim ;i++){
116 //fscanf(infile,"%f\t %e\t %e",&beta,&aRead,&cRead);
117 infile >> beta >> aRead >> cRead;
118 aMajorantSurfaceParameterTable[i][level] = aRead;
119 cMajorantSurfaceParameterTable[i][level] = cRead;
120 }
121 infile.close();
122 }
123}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
float G4float
Definition: G4Types.hh:84
int G4int
Definition: G4Types.hh:85

◆ ~G4PhotoElectricAngularGeneratorPolarized()

G4PhotoElectricAngularGeneratorPolarized::~G4PhotoElectricAngularGeneratorPolarized ( )

Definition at line 125 of file G4PhotoElectricAngularGeneratorPolarized.cc.

126{}

Member Function Documentation

◆ PerpendicularVector()

G4ThreeVector G4PhotoElectricAngularGeneratorPolarized::PerpendicularVector ( const G4ThreeVector a) const
protected

Definition at line 414 of file G4PhotoElectricAngularGeneratorPolarized.cc.

416{
417 G4double dx = a.x();
418 G4double dy = a.y();
419 G4double dz = a.z();
420 G4double x = dx < 0.0 ? -dx : dx;
421 G4double y = dy < 0.0 ? -dy : dy;
422 G4double z = dz < 0.0 ? -dz : dz;
423 if (x < y) {
424 return x < z ? G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy);
425 }else{
426 return y < z ? G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0);
427 }
428}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const

◆ PrintGeneratorInformation()

void G4PhotoElectricAngularGeneratorPolarized::PrintGeneratorInformation ( ) const

Definition at line 403 of file G4PhotoElectricAngularGeneratorPolarized.cc.

404{
405 G4cout << "\n" << G4endl;
406 G4cout << "Polarized Photoelectric Angular Generator" << G4endl;
407 G4cout << "PhotoElectric Electron Angular Generator based on the general Gavrila photoelectron angular distribution" << G4endl;
408 G4cout << "Includes polarization effects for K and L1 atomic shells, according to Gavrilla (1959, 1961)." << G4endl;
409 G4cout << "For higher shells the L1 cross-section is used." << G4endl;
410 G4cout << "(see Physics Reference Manual) \n" << G4endl;
411}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ SampleDirection()

G4ThreeVector & G4PhotoElectricAngularGeneratorPolarized::SampleDirection ( const G4DynamicParticle dp,
G4double  eKinEnergy,
G4int  shellId,
const G4Material mat = 0 
)
virtual

Implements G4VEmAngularDistribution.

Definition at line 129 of file G4PhotoElectricAngularGeneratorPolarized.cc.

134{
135 // (shellId == 0) - K-shell - Polarized model for K-shell
136 // (shellId > 0) - L1-shell - Polarized model for L1 and higher shells
137
138 // Calculate Lorentz term (gamma) and beta parameters
139 G4double gamma = 1 + eKinEnergy/electron_mass_c2;
140 G4double beta = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
141
142 const G4ThreeVector& direction = dp->GetMomentumDirection();
143 const G4ThreeVector& polarization = dp->GetPolarization();
144
145 G4double theta, phi = 0;
146 // Majorant surface parameters
147 // function of the outgoing electron kinetic energy
148 G4double aBeta = 0;
149 G4double cBeta = 0;
150
151 // For the outgoing kinetic energy
152 // find the current majorant surface parameters
153 PhotoElectronGetMajorantSurfaceAandCParameters(shellId, beta, &aBeta, &cBeta);
154
155 // Generate pho and theta according to the shell level
156 // and beta parameter of the electron
157 PhotoElectronGeneratePhiAndTheta(shellId, beta, aBeta, cBeta, &phi, &theta);
158
159 // Determine the rotation matrix
160 const G4RotationMatrix rotation =
161 PhotoElectronRotationMatrix(direction, polarization);
162
163 // Compute final direction of the outgoing electron
164 fLocalDirection = PhotoElectronComputeFinalDirection(rotation, theta, phi);
165
166 return fLocalDirection;
167}
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const

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