Geant4 9.6.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
 
const G4StringGetName () const
 

Protected Member Functions

G4ThreeVector PerpendicularVector (const G4ThreeVector &a) const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 52 of file G4PhotoElectricAngularGeneratorPolarized.hh.

Constructor & Destructor Documentation

◆ G4PhotoElectricAngularGeneratorPolarized()

G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized ( )

Definition at line 70 of file G4PhotoElectricAngularGeneratorPolarized.cc.

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

◆ ~G4PhotoElectricAngularGeneratorPolarized()

G4PhotoElectricAngularGeneratorPolarized::~G4PhotoElectricAngularGeneratorPolarized ( )

Definition at line 124 of file G4PhotoElectricAngularGeneratorPolarized.cc.

125{}

Member Function Documentation

◆ PerpendicularVector()

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

Definition at line 413 of file G4PhotoElectricAngularGeneratorPolarized.cc.

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

◆ PrintGeneratorInformation()

void G4PhotoElectricAngularGeneratorPolarized::PrintGeneratorInformation ( ) const

Definition at line 402 of file G4PhotoElectricAngularGeneratorPolarized.cc.

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

◆ SampleDirection()

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

Implements G4VEmAngularDistribution.

Definition at line 128 of file G4PhotoElectricAngularGeneratorPolarized.cc.

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

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