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

#include <G4OrePowellAtRestModel.hh>

+ Inheritance diagram for G4OrePowellAtRestModel:

Public Member Functions

 G4OrePowellAtRestModel ()
 
 ~G4OrePowellAtRestModel () override=default
 
void SampleSecondaries (std::vector< G4DynamicParticle * > &secParticles, G4double &, const G4Material *) const override
 
void PrintGeneratorInformation () const override
 
G4OrePowellAtRestModeloperator= (const G4OrePowellAtRestModel &right)=delete
 
 G4OrePowellAtRestModel (const G4OrePowellAtRestModel &)=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 50 of file G4OrePowellAtRestModel.hh.

Constructor & Destructor Documentation

◆ G4OrePowellAtRestModel() [1/2]

G4OrePowellAtRestModel::G4OrePowellAtRestModel ( )

Definition at line 51 of file G4OrePowellAtRestModel.cc.

51: G4VPositronAtRestModel("OrePowell") {}
G4VPositronAtRestModel(const G4String &name)

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

◆ ~G4OrePowellAtRestModel()

G4OrePowellAtRestModel::~G4OrePowellAtRestModel ( )
overridedefault

◆ G4OrePowellAtRestModel() [2/2]

G4OrePowellAtRestModel::G4OrePowellAtRestModel ( const G4OrePowellAtRestModel & )
delete

Member Function Documentation

◆ operator=()

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

◆ PrintGeneratorInformation()

void G4OrePowellAtRestModel::PrintGeneratorInformation ( ) const
overridevirtual

Implements G4VPositronAtRestModel.

Definition at line 145 of file G4OrePowellAtRestModel.cc.

146{
147 G4cout << "Orel Powell AtRest positron 3-gamma annihilation model" << G4endl;
148}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout

◆ SampleSecondaries()

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

Implements G4VPositronAtRestModel.

Definition at line 55 of file G4OrePowellAtRestModel.cc.

58{
59 static const G4double PositronMass = CLHEP::electron_mass_c2;
60 const G4double ymax = 8.1;
61
62 CLHEP::HepRandomEngine* rndmEngine = G4Random::getTheEngine();
63
64 G4double cos12;
65 G4double cos13;
66 G4double sin12;
67 G4double sin13;
68 G4double r1;
69 G4double r2;
70 G4double r3;
71 G4double pdf;
72
73 G4double rndmv2[2];
74 G4double rndmv1;
75
76 do {
77 rndmv1 = rndmEngine->flat();
78
79 do {
80 rndmEngine->flatArray(2, rndmv2);
81 // energies of photon1 and photon2 normalized to electron rest mass
82 r1 = rndmv2[0];
83 r2 = rndmv2[1];
84 // energy conservation, with positronium assumed = 2 * electron rest mass
85 r3 = 2.0 - (r1+r2);
86 // cosine of angles between photons, from momentum conservation
87 cos12=(r3*r3 - r1*r1 -r2*r2)/(2*r1*r2);
88 cos13=(r2*r2 - r1*r1 -r3*r3)/(2*r1*r3);
89 // request both cosines < 1.
90 } while ( std::abs(cos12) > 1 || std::abs(cos13) > 1 );
91
92 sin12 = std::sqrt((1 + cos12)*(1 - cos12));
93 sin13 = -std::sqrt((1 + cos13)*(1 - cos13));
94
95 G4double cos23=cos12*cos13+sin12*sin13;
96
97 pdf = (1 - cos12)*(1 - cos12) + (1 - cos13)*(1 - cos13) + (1 - cos23)*(1 - cos23);
98 } while ( pdf < ymax * rndmv1 );
99 // END of Sampling
100
101 // photon directions in the decay plane, photon 1 along z, x perp to the plane.
102 G4ThreeVector PhotonMomentum1(0., 0., 1.);
103 G4ThreeVector PhotonMomentum2(0.,sin12,cos12);
104 G4ThreeVector PhotonMomentum3(0.,sin13,cos13);
105
106 // First Gamma direction
108
109 PhotonMomentum1.rotateUz(dir1);
110 PhotonMomentum2.rotateUz(dir1);
111 PhotonMomentum3.rotateUz(dir1);
112
113 auto aGamma1 = new G4DynamicParticle(G4Gamma::Gamma(), PhotonMomentum1,
114 r1 * PositronMass);
115
116 //Random polarization
117 G4double phi1 = CLHEP::twopi * G4UniformRand();
118 G4ThreeVector pol1(std::cos(phi1),std::sin(phi1),0.0);
119 pol1.rotateUz(PhotonMomentum1);
120 aGamma1->SetPolarization(pol1);
121 secParticles.push_back(aGamma1);
122
123 auto aGamma2 = new G4DynamicParticle(G4Gamma::Gamma(), PhotonMomentum2,
124 r2 * PositronMass);
125
126 G4double phi2 = CLHEP::twopi * G4UniformRand();
127 G4ThreeVector pol2(std::cos(phi2),std::sin(phi2),0.0);
128 pol2.rotateUz(PhotonMomentum2);
129 aGamma2->SetPolarization(pol2);
130 secParticles.push_back(aGamma2);
131
132 auto aGamma3 = new G4DynamicParticle(G4Gamma::Gamma(), PhotonMomentum3,
133 r3 * PositronMass);
134
135 G4double phi3 = CLHEP::twopi * G4UniformRand();
136 G4ThreeVector pol3(std::cos(phi3),std::sin(phi3),0.0);
137 pol3.rotateUz(PhotonMomentum3);
138 aGamma3->SetPolarization(pol3);
139 secParticles.push_back(aGamma3);
140
141}
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
#define G4UniformRand()
Definition Randomize.hh:52
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
static G4Gamma * Gamma()
Definition G4Gamma.cc:81

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