Geant4 11.1.1
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 ()
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double eKinEnergy, G4int shellId, const G4Material *mat=nullptr) override
 
void PrintGeneratorInformation () const override
 
G4PhotoElectricAngularGeneratorPolarizedoperator= (const G4PhotoElectricAngularGeneratorPolarized &right)=delete
 
 G4PhotoElectricAngularGeneratorPolarized (const G4PhotoElectricAngularGeneratorPolarized &)=delete
 
- 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)
 
virtual void PrintGeneratorInformation () const
 
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() [1/2]

G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized ( )
explicit

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
84 //Gavrila angular photoelectron distribution
85 for(G4int level = 0; level < 2; level++){
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 const char* path = G4FindDataDir("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
114 // energy and shell level
115 G4float aRead=0,cRead=0, beta=0;
116 for(G4int i=0 ; i<arrayDim ;++i){
117 infile >> beta >> aRead >> cRead;
118 aMajorantSurfaceParameterTable[i][level] = aRead;
119 cMajorantSurfaceParameterTable[i][level] = cRead;
120 }
121 infile.close();
122 }
123}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
float G4float
Definition: G4Types.hh:84
int G4int
Definition: G4Types.hh:85

◆ ~G4PhotoElectricAngularGeneratorPolarized()

G4PhotoElectricAngularGeneratorPolarized::~G4PhotoElectricAngularGeneratorPolarized ( )

Definition at line 127 of file G4PhotoElectricAngularGeneratorPolarized.cc.

128{}

◆ G4PhotoElectricAngularGeneratorPolarized() [2/2]

G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized ( const G4PhotoElectricAngularGeneratorPolarized )
delete

Member Function Documentation

◆ operator=()

G4PhotoElectricAngularGeneratorPolarized & G4PhotoElectricAngularGeneratorPolarized::operator= ( const G4PhotoElectricAngularGeneratorPolarized right)
delete

◆ PerpendicularVector()

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

Definition at line 424 of file G4PhotoElectricAngularGeneratorPolarized.cc.

426{
427 G4double dx = a.x();
428 G4double dy = a.y();
429 G4double dz = a.z();
430 G4double x = dx < 0.0 ? -dx : dx;
431 G4double y = dy < 0.0 ? -dy : dy;
432 G4double z = dz < 0.0 ? -dz : dz;
433 if (x < y) {
434 return x < z ? G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy);
435 }else{
436 return y < z ? G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0);
437 }
438}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
double z() const
double x() const
double y() const

◆ PrintGeneratorInformation()

void G4PhotoElectricAngularGeneratorPolarized::PrintGeneratorInformation ( ) const
overridevirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 411 of file G4PhotoElectricAngularGeneratorPolarized.cc.

412{
413 G4cout << "\n" << G4endl;
414 G4cout << "Polarized Photoelectric Angular Generator" << G4endl;
415 G4cout << "PhotoElectric Electron Angular Generator based on the general Gavrila photoelectron angular distribution" << G4endl;
416 G4cout << "Includes polarization effects for K and L1 atomic shells, according to Gavrilla (1959, 1961)." << G4endl;
417 G4cout << "For higher shells the L1 cross-section is used." << G4endl;
418 G4cout << "(see Physics Reference Manual) \n" << G4endl;
419}
#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 = nullptr 
)
overridevirtual

Implements G4VEmAngularDistribution.

Definition at line 133 of file G4PhotoElectricAngularGeneratorPolarized.cc.

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

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