#include <G4PhotoElectricAngularGeneratorSauterGavrila.hh>
|
| G4PhotoElectricAngularGeneratorSauterGavrila () |
|
| G4PhotoElectricAngularGeneratorSauterGavrila (const G4PhotoElectricAngularGeneratorSauterGavrila &)=delete |
|
const G4String & | GetName () const |
|
G4PhotoElectricAngularGeneratorSauterGavrila & | operator= (const G4PhotoElectricAngularGeneratorSauterGavrila &right)=delete |
|
void | PrintGeneratorInformation () const override |
|
G4ThreeVector & | SampleDirection (const G4DynamicParticle *dp, G4double e=0.0, G4int shellId=0, const G4Material *mat=nullptr) override |
|
virtual G4ThreeVector & | SampleDirectionForShell (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) |
|
| ~G4PhotoElectricAngularGeneratorSauterGavrila () |
|
◆ G4PhotoElectricAngularGeneratorSauterGavrila() [1/2]
G4PhotoElectricAngularGeneratorSauterGavrila::G4PhotoElectricAngularGeneratorSauterGavrila |
( |
| ) |
|
|
explicit |
◆ ~G4PhotoElectricAngularGeneratorSauterGavrila()
G4PhotoElectricAngularGeneratorSauterGavrila::~G4PhotoElectricAngularGeneratorSauterGavrila |
( |
| ) |
|
◆ G4PhotoElectricAngularGeneratorSauterGavrila() [2/2]
◆ GetName()
const G4String & G4VEmAngularDistribution::GetName |
( |
| ) |
const |
|
inlineinherited |
◆ operator=()
◆ PrintGeneratorInformation()
void G4PhotoElectricAngularGeneratorSauterGavrila::PrintGeneratorInformation |
( |
| ) |
const |
|
overridevirtual |
◆ SampleDirection()
Implements G4VEmAngularDistribution.
Definition at line 67 of file G4PhotoElectricAngularGeneratorSauterGavrila.cc.
70{
71
72
73
74
81
82 if (gamma > 5.) {
85
86
87 }
88
89 G4double beta = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
90 G4double b = 0.5*gamma*(gamma - 1)*(gamma - 2);
91
93 if (gamma < 2.) grejsup = gamma*gamma*(1.+b-
beta*b);
94 else grejsup = gamma*gamma*(1.+b+
beta*b);
95
97 costeta = (rndm+
beta)/(rndm*
beta+1.);
98 term = 1.-
beta*costeta;
99 greject = (1.-costeta*costeta)*(1.+b*term)/(term*term);
101
102 sinteta = std::sqrt((1 - costeta)*(1 + costeta));
106}
static constexpr double twopi
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4ThreeVector fLocalDirection
References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, source.hepunit::electron_mass_c2, G4VEmAngularDistribution::fLocalDirection, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), and twopi.
◆ SampleDirectionForShell()
Reimplemented in G4DeltaAngle, G4DNABornAngle, and G4DNARuddAngle.
Definition at line 69 of file G4VEmAngularDistribution.cc.
74{
76}
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
References G4VEmAngularDistribution::SampleDirection(), and Z.
Referenced by G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), and G4MicroElecInelasticModel_new::SampleSecondaries().
◆ SamplePairDirections()
◆ fLocalDirection
Definition at line 103 of file G4VEmAngularDistribution.hh.
Referenced by G4VEmAngularDistribution::G4VEmAngularDistribution(), G4SauterGavrilaAngularDistribution::SampleDirection(), SampleDirection(), G4PhotoElectricAngularGeneratorPolarized::SampleDirection(), G4ModifiedMephi::SampleDirection(), G4DeltaAngle::SampleDirection(), G4DeltaAngleFreeScat::SampleDirection(), G4DipBustGenerator::SampleDirection(), G4ModifiedTsai::SampleDirection(), G4Generator2BN::SampleDirection(), G4Generator2BS::SampleDirection(), G4PenelopeBremsstrahlungAngular::SampleDirection(), G4RayleighAngularGenerator::SampleDirection(), G4AngleDirect::SampleDirection(), G4DNABornAngle::SampleDirectionForShell(), and G4DNARuddAngle::SampleDirectionForShell().
◆ fName
G4String G4VEmAngularDistribution::fName |
|
privateinherited |
◆ fPolarisation
G4bool G4VEmAngularDistribution::fPolarisation |
|
protectedinherited |
The documentation for this class was generated from the following files: