Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4eeToTwoGammaModel.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 // $Id$
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4eeToTwoGammaModel
34 //
35 // Author: Vladimir Ivanchenko on base of Michel Maire code
36 //
37 // Creation date: 02.08.2004
38 //
39 // Modifications:
40 // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
41 // 18-04-05 Compute CrossSectionPerVolume (V.Ivanchenko)
42 // 06-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
43 // 29-06-06 Fix problem for zero energy incident positron (V.Ivanchenko)
44 // 20-10-06 Add theGamma as a member (V.Ivanchenko)
45 //
46 //
47 // Class Description:
48 //
49 // Implementation of e+ annihilation into 2 gamma
50 //
51 // The secondaries Gamma energies are sampled using the Heitler cross section.
52 //
53 // A modified version of the random number techniques of Butcher & Messel
54 // is used (Nuc Phys 20(1960),15).
55 //
56 // GEANT4 internal units.
57 //
58 // Note 1: The initial electron is assumed free and at rest.
59 //
60 // Note 2: The annihilation processes producing one or more than two photons are
61 // ignored, as negligible compared to the two photons process.
62 
63 
64 
65 //
66 // -------------------------------------------------------------------
67 //
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
71 #include "G4eeToTwoGammaModel.hh"
72 #include "G4PhysicalConstants.hh"
73 #include "G4SystemOfUnits.hh"
74 #include "G4TrackStatus.hh"
75 #include "G4Electron.hh"
76 #include "G4Positron.hh"
77 #include "G4Gamma.hh"
78 #include "Randomize.hh"
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82 
83 using namespace std;
84 
86  const G4String& nam)
87  : G4VEmModel(nam),
89  isInitialised(false)
90 {
91  theGamma = G4Gamma::Gamma();
92  fParticleChange = 0;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96 
98 {}
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101 
103  const G4DataVector&)
104 {
105  if(isInitialised) { return; }
106  fParticleChange = GetParticleChangeForGamma();
107  isInitialised = true;
108 }
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
111 
113  const G4ParticleDefinition*,
114  G4double kineticEnergy,
116 {
117  // Calculates the cross section per electron of annihilation into two photons
118  // from the Heilter formula.
119 
120  G4double ekin = std::max(eV,kineticEnergy);
121 
122  G4double tau = ekin/electron_mass_c2;
123  G4double gam = tau + 1.0;
124  G4double gamma2= gam*gam;
125  G4double bg2 = tau * (tau+2.0);
126  G4double bg = sqrt(bg2);
127 
128  G4double cross = pi_rcl2*((gamma2+4*gam+1.)*log(gam+bg) - (gam+3.)*bg)
129  / (bg2*(gam+1.));
130  return cross;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136  const G4ParticleDefinition* p,
137  G4double kineticEnergy, G4double Z,
139 {
140  // Calculates the cross section per atom of annihilation into two photons
141 
142  G4double cross = Z*ComputeCrossSectionPerElectron(p,kineticEnergy);
143  return cross;
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
149  const G4Material* material,
150  const G4ParticleDefinition* p,
151  G4double kineticEnergy,
153 {
154  // Calculates the cross section per volume of annihilation into two photons
155 
156  G4double eDensity = material->GetElectronDensity();
157  G4double cross = eDensity*ComputeCrossSectionPerElectron(p,kineticEnergy);
158  return cross;
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
163 void G4eeToTwoGammaModel::SampleSecondaries(vector<G4DynamicParticle*>* vdp,
164  const G4MaterialCutsCouple*,
165  const G4DynamicParticle* dp,
166  G4double,
167  G4double)
168 {
169  G4double PositKinEnergy = dp->GetKineticEnergy();
170  G4DynamicParticle *aGamma1, *aGamma2;
171 
172  // Case at rest
173  if(PositKinEnergy == 0.0) {
174  G4double cost = 2.*G4UniformRand()-1.;
175  G4double sint = sqrt((1. - cost)*(1. + cost));
176  G4double phi = twopi * G4UniformRand();
177  G4ThreeVector dir(sint*cos(phi), sint*sin(phi), cost);
178  phi = twopi * G4UniformRand();
179  G4ThreeVector pol(cos(phi), sin(phi), 0.0);
180  pol.rotateUz(dir);
181  aGamma1 = new G4DynamicParticle(theGamma, dir, electron_mass_c2);
182  aGamma1->SetPolarization(pol.x(),pol.y(),pol.z());
183  aGamma2 = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
184  aGamma1->SetPolarization(-pol.x(),-pol.y(),-pol.z());
185 
186  } else {
187 
188  G4ThreeVector PositDirection = dp->GetMomentumDirection();
189 
190  G4double tau = PositKinEnergy/electron_mass_c2;
191  G4double gam = tau + 1.0;
192  G4double tau2 = tau + 2.0;
193  G4double sqgrate = sqrt(tau/tau2)*0.5;
194  G4double sqg2m1 = sqrt(tau*tau2);
195 
196  // limits of the energy sampling
197  G4double epsilmin = 0.5 - sqgrate;
198  G4double epsilmax = 0.5 + sqgrate;
199  G4double epsilqot = epsilmax/epsilmin;
200 
201  //
202  // sample the energy rate of the created gammas
203  //
204  G4double epsil, greject;
205 
206  do {
207  epsil = epsilmin*pow(epsilqot,G4UniformRand());
208  greject = 1. - epsil + (2.*gam*epsil-1.)/(epsil*tau2*tau2);
209  } while( greject < G4UniformRand() );
210 
211  //
212  // scattered Gamma angles. ( Z - axis along the parent positron)
213  //
214 
215  G4double cost = (epsil*tau2-1.)/(epsil*sqg2m1);
216  if(std::abs(cost) > 1.0) {
217  G4cout << "### G4eeToTwoGammaModel WARNING cost= " << cost
218  << " positron Ekin(MeV)= " << PositKinEnergy
219  << " gamma epsil= " << epsil
220  << G4endl;
221  if(cost > 1.0) cost = 1.0;
222  else cost = -1.0;
223  }
224  G4double sint = sqrt((1.+cost)*(1.-cost));
225  G4double phi = twopi * G4UniformRand();
226 
227  //
228  // kinematic of the created pair
229  //
230 
231  G4double TotalAvailableEnergy = PositKinEnergy + 2.0*electron_mass_c2;
232  G4double Phot1Energy = epsil*TotalAvailableEnergy;
233 
234  G4ThreeVector Phot1Direction(sint*cos(phi), sint*sin(phi), cost);
235  Phot1Direction.rotateUz(PositDirection);
236  aGamma1 = new G4DynamicParticle (theGamma,Phot1Direction, Phot1Energy);
237  phi = twopi * G4UniformRand();
238  G4ThreeVector pol1(cos(phi), sin(phi), 0.0);
239  pol1.rotateUz(Phot1Direction);
240  aGamma1->SetPolarization(pol1.x(),pol1.y(),pol1.z());
241 
242  G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
243  G4double PositP= sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
244  G4ThreeVector dir = PositDirection*PositP - Phot1Direction*Phot1Energy;
245  G4ThreeVector Phot2Direction = dir.unit();
246 
247  // create G4DynamicParticle object for the particle2
248  aGamma2 = new G4DynamicParticle (theGamma,Phot2Direction, Phot2Energy);
249 
251  aGamma2->SetPolarization(-pol1.x(),-pol1.y(),-pol1.z());
252  }
253  /*
254  G4cout << "Annihilation in fly: e0= " << PositKinEnergy
255  << " m= " << electron_mass_c2
256  << " e1= " << Phot1Energy
257  << " e2= " << Phot2Energy << " dir= " << dir
258  << " -> " << Phot1Direction << " "
259  << Phot2Direction << G4endl;
260  */
261 
262  vdp->push_back(aGamma1);
263  vdp->push_back(aGamma2);
264  fParticleChange->SetProposedKineticEnergy(0.);
265  fParticleChange->ProposeTrackStatus(fStopAndKill);
266 }
267 
268 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....