Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OrePowellAtRestModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// GEANT4 Class file
28//
29//
30// File name: G4OrePowellAtRestModel
31//
32// Author: I.Semeniouk & D.Bernard
33//
34// Creation date: 04 Juin 2024
35//
36// -------------------------------------------------------------------
37//
38
40#include "G4DynamicParticle.hh"
41#include "G4Material.hh"
42#include "Randomize.hh"
43#include "G4Gamma.hh"
44#include "G4RandomDirection.hh"
45#include "G4ThreeVector.hh"
47#include "G4SystemOfUnits.hh"
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
56 std::vector<G4DynamicParticle*>& secParticles,
57 G4double&, const G4Material*) const
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}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144
146{
147 G4cout << "Orel Powell AtRest positron 3-gamma annihilation model" << G4endl;
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector & rotateUz(const Hep3Vector &)
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SampleSecondaries(std::vector< G4DynamicParticle * > &secParticles, G4double &, const G4Material *) const override
void PrintGeneratorInformation() const override
G4VPositronAtRestModel(const G4String &name)