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

#include <G4VEmAngularDistribution.hh>

+ Inheritance diagram for G4VEmAngularDistribution:

Public Member Functions

 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 Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 58 of file G4VEmAngularDistribution.hh.

Constructor & Destructor Documentation

◆ G4VEmAngularDistribution() [1/2]

G4VEmAngularDistribution::G4VEmAngularDistribution ( const G4String name)
explicit

Definition at line 55 of file G4VEmAngularDistribution.cc.

56 : fName(name)
57{
58 fLocalDirection.set(0.0,0.0,1.0);
60}
void set(double x, double y, double z)
static G4EmParameters * Instance()
G4bool EnablePolarisation() const

◆ ~G4VEmAngularDistribution()

G4VEmAngularDistribution::~G4VEmAngularDistribution ( )
virtual

Definition at line 64 of file G4VEmAngularDistribution.cc.

65{}

◆ G4VEmAngularDistribution() [2/2]

G4VEmAngularDistribution::G4VEmAngularDistribution ( const G4VEmAngularDistribution )
delete

Member Function Documentation

◆ GetName()

const G4String & G4VEmAngularDistribution::GetName ( ) const
inline

Definition at line 109 of file G4VEmAngularDistribution.hh.

110{
111 return fName;
112}

◆ operator=()

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

◆ SampleDirection()

◆ SampleDirectionForShell()

G4ThreeVector & G4VEmAngularDistribution::SampleDirectionForShell ( const G4DynamicParticle dp,
G4double  finalTotalEnergy,
G4int  Z,
G4int  shellID,
const G4Material mat 
)
virtual

◆ SamplePairDirections()

void G4VEmAngularDistribution::SamplePairDirections ( const G4DynamicParticle dp,
G4double  elecKinEnergy,
G4double  posiKinEnergy,
G4ThreeVector dirElectron,
G4ThreeVector dirPositron,
G4int  Z = 0,
const G4Material mat = nullptr 
)
virtual

Reimplemented in G4ModifiedMephi, G4DipBustGenerator, and G4ModifiedTsai.

Definition at line 80 of file G4VEmAngularDistribution.cc.

85{
86 dirElectron = dp->GetMomentumDirection();
87 dirPositron = dp->GetMomentumDirection();
88}
const G4ThreeVector & GetMomentumDirection() const

Referenced by G4MuPairProductionModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), and G4PairProductionRelModel::SampleSecondaries().

Member Data Documentation

◆ fLocalDirection

◆ fPolarisation

G4bool G4VEmAngularDistribution::fPolarisation
protected

Definition at line 102 of file G4VEmAngularDistribution.hh.

Referenced by G4VEmAngularDistribution().


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