Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DeltaAngleFreeScat Class Reference

#include <G4DeltaAngleFreeScat.hh>

+ Inheritance diagram for G4DeltaAngleFreeScat:

Public Member Functions

 G4DeltaAngleFreeScat (const G4String &name="")
 
 ~G4DeltaAngleFreeScat () override
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) final
 
void PrintGeneratorInformation () const final
 
G4DeltaAngleFreeScatoperator= (const G4DeltaAngleFreeScat &right)=delete
 
 G4DeltaAngleFreeScat (const G4DeltaAngleFreeScat &)=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
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Detailed Description

Definition at line 53 of file G4DeltaAngleFreeScat.hh.

Constructor & Destructor Documentation

◆ G4DeltaAngleFreeScat() [1/2]

G4DeltaAngleFreeScat::G4DeltaAngleFreeScat ( const G4String name = "")
explicit

◆ ~G4DeltaAngleFreeScat()

G4DeltaAngleFreeScat::~G4DeltaAngleFreeScat ( )
overridedefault

◆ G4DeltaAngleFreeScat() [2/2]

G4DeltaAngleFreeScat::G4DeltaAngleFreeScat ( const G4DeltaAngleFreeScat )
delete

Member Function Documentation

◆ operator=()

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

◆ PrintGeneratorInformation()

void G4DeltaAngleFreeScat::PrintGeneratorInformation ( ) const
finalvirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 81 of file G4DeltaAngleFreeScat.cc.

82{}

◆ SampleDirection()

G4ThreeVector & G4DeltaAngleFreeScat::SampleDirection ( const G4DynamicParticle dp,
G4double  kinEnergyFinal,
G4int  Z,
const G4Material mat = nullptr 
)
finalvirtual

Implements G4VEmAngularDistribution.

Definition at line 62 of file G4DeltaAngleFreeScat.cc.

65{
66 G4double deltaMomentum =
67 sqrt(kinEnergyFinal*(kinEnergyFinal + 2*electron_mass_c2));
68
69 G4double costet = kinEnergyFinal*(dp->GetTotalEnergy() + electron_mass_c2) /
70 (deltaMomentum * dp->GetTotalMomentum());
71
72 G4double phi = G4UniformRand()*twopi;
73 G4double sintet = sqrt((1 - costet)*(1 + costet));
74
75 fLocalDirection.set(sintet*cos(phi), sintet*sin(phi), costet);
77
78 return fLocalDirection;
79}
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const

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