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

#include <G4SimplePositronAtRestModel.hh>

+ Inheritance diagram for G4SimplePositronAtRestModel:

Public Member Functions

 G4SimplePositronAtRestModel ()
 
 ~G4SimplePositronAtRestModel () override=default
 
void SampleSecondaries (std::vector< G4DynamicParticle * > &secParticles, G4double &, const G4Material *) const override
 
void PrintGeneratorInformation () const override
 
G4SimplePositronAtRestModeloperator= (const G4SimplePositronAtRestModel &right)=delete
 
 G4SimplePositronAtRestModel (const G4SimplePositronAtRestModel &)=delete
 
- Public Member Functions inherited from G4VPositronAtRestModel
 G4VPositronAtRestModel (const G4String &name)
 
virtual ~G4VPositronAtRestModel ()=default
 
const G4StringGetName () const
 
G4VPositronAtRestModeloperator= (const G4VPositronAtRestModel &right)=delete
 
 G4VPositronAtRestModel (const G4VPositronAtRestModel &)=delete
 

Detailed Description

Definition at line 49 of file G4SimplePositronAtRestModel.hh.

Constructor & Destructor Documentation

◆ G4SimplePositronAtRestModel() [1/2]

G4SimplePositronAtRestModel::G4SimplePositronAtRestModel ( )

Definition at line 51 of file G4SimplePositronAtRestModel.cc.

52 : G4VPositronAtRestModel("Simple")
53{}
G4VPositronAtRestModel(const G4String &name)

Referenced by G4SimplePositronAtRestModel(), and operator=().

◆ ~G4SimplePositronAtRestModel()

G4SimplePositronAtRestModel::~G4SimplePositronAtRestModel ( )
overridedefault

◆ G4SimplePositronAtRestModel() [2/2]

G4SimplePositronAtRestModel::G4SimplePositronAtRestModel ( const G4SimplePositronAtRestModel & )
delete

Member Function Documentation

◆ operator=()

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

◆ PrintGeneratorInformation()

void G4SimplePositronAtRestModel::PrintGeneratorInformation ( ) const
overridevirtual

Implements G4VPositronAtRestModel.

Definition at line 83 of file G4SimplePositronAtRestModel.cc.

84{
85 G4cout << "\n" << G4endl;
86 G4cout << "Simple AtRest positron 2-gamma annihilation model" << G4endl;
87}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ SampleSecondaries()

void G4SimplePositronAtRestModel::SampleSecondaries ( std::vector< G4DynamicParticle * > & secParticles,
G4double & ,
const G4Material *  ) const
overridevirtual

Implements G4VPositronAtRestModel.

Definition at line 57 of file G4SimplePositronAtRestModel.cc.

60{
62 auto aGamma1 = new G4DynamicParticle(G4Gamma::Gamma(), dir1,
63 CLHEP::electron_mass_c2);
64 G4double phi = CLHEP::twopi * G4UniformRand();
65 G4double cosphi = std::cos(phi);
66 G4double sinphi = std::sin(phi);
67 G4ThreeVector pol1(cosphi, sinphi, 0.0);
68 pol1.rotateUz(dir1);
69 aGamma1->SetPolarization(pol1);
70 secParticles.push_back(aGamma1);
71
72 G4ThreeVector dir2 = -dir1;
73 auto aGamma2 = new G4DynamicParticle(G4Gamma::Gamma(), dir2,
74 CLHEP::electron_mass_c2);
75 G4ThreeVector pol2(-sinphi, cosphi, 0.0);
76 pol2.rotateUz(dir1);
77 aGamma2->SetPolarization(pol2);
78 secParticles.push_back(aGamma2);
79}
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52
static G4Gamma * Gamma()
Definition G4Gamma.cc:81

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