Geant4-11
|
#include <G4VEmAngularDistribution.hh>
Public Member Functions | |
G4VEmAngularDistribution (const G4String &name) | |
G4VEmAngularDistribution (const G4VEmAngularDistribution &)=delete | |
const G4String & | GetName () const |
G4VEmAngularDistribution & | operator= (const G4VEmAngularDistribution &right)=delete |
virtual void | PrintGeneratorInformation () const |
virtual G4ThreeVector & | SampleDirection (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0 |
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) |
virtual | ~G4VEmAngularDistribution () |
Protected Attributes | |
G4ThreeVector | fLocalDirection |
G4bool | fPolarisation |
Private Attributes | |
G4String | fName |
Definition at line 58 of file G4VEmAngularDistribution.hh.
|
explicit |
Definition at line 55 of file G4VEmAngularDistribution.cc.
References G4EmParameters::EnablePolarisation(), fLocalDirection, fPolarisation, G4EmParameters::Instance(), and CLHEP::Hep3Vector::set().
|
virtual |
Definition at line 64 of file G4VEmAngularDistribution.cc.
|
delete |
|
inline |
Definition at line 111 of file G4VEmAngularDistribution.hh.
References fName.
|
delete |
|
virtual |
Reimplemented in G4DeltaAngleFreeScat, G4DipBustGenerator, G4ModifiedTsai, G4DNABornAngle, G4DNARuddAngle, G4Generator2BN, G4Generator2BS, G4PhotoElectricAngularGeneratorPolarized, G4PhotoElectricAngularGeneratorSauterGavrila, G4ModifiedMephi, and G4SauterGavrilaAngularDistribution.
Definition at line 92 of file G4VEmAngularDistribution.cc.
|
pure virtual |
Implemented in G4SauterGavrilaAngularDistribution, G4PhotoElectricAngularGeneratorSauterGavrila, G4PhotoElectricAngularGeneratorPolarized, G4ModifiedMephi, G4DeltaAngle, G4DeltaAngleFreeScat, G4DNABornAngle, G4DNARuddAngle, G4DipBustGenerator, G4ModifiedTsai, G4Generator2BN, G4Generator2BS, G4PenelopeBremsstrahlungAngular, G4RayleighAngularGenerator, and G4AngleDirect.
Referenced by G4SynchrotronRadiation::PostStepDoIt(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), SampleDirectionForShell(), G4LivermoreBremsstrahlungModel::SampleSecondaries(), G4eBremParametrizedModel::SampleSecondaries(), G4eBremsstrahlungRelModel::SampleSecondaries(), G4SeltzerBergerModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4MuBremsstrahlungModel::SampleSecondaries(), G4AtimaEnergyLossModel::SampleSecondaries(), G4BetheBlochModel::SampleSecondaries(), G4BraggIonModel::SampleSecondaries(), G4BraggModel::SampleSecondaries(), G4ICRU73QOModel::SampleSecondaries(), G4LindhardSorensenIonModel::SampleSecondaries(), and G4MollerBhabhaModel::SampleSecondaries().
|
virtual |
Reimplemented in G4DeltaAngle, G4DNABornAngle, and G4DNARuddAngle.
Definition at line 69 of file G4VEmAngularDistribution.cc.
References SampleDirection(), and Z.
Referenced by G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), and G4MicroElecInelasticModel_new::SampleSecondaries().
|
virtual |
Reimplemented in G4ModifiedMephi, G4DipBustGenerator, and G4ModifiedTsai.
Definition at line 80 of file G4VEmAngularDistribution.cc.
References G4DynamicParticle::GetMomentumDirection().
Referenced by G4MuPairProductionModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), and G4PairProductionRelModel::SampleSecondaries().
|
protected |
Definition at line 103 of file G4VEmAngularDistribution.hh.
Referenced by G4VEmAngularDistribution(), G4SauterGavrilaAngularDistribution::SampleDirection(), G4PhotoElectricAngularGeneratorSauterGavrila::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().
|
private |
Definition at line 108 of file G4VEmAngularDistribution.hh.
Referenced by GetName().
|
protected |
Definition at line 104 of file G4VEmAngularDistribution.hh.
Referenced by G4VEmAngularDistribution().