source: trunk/source/processes/electromagnetic/standard/src/G4eeToTwoGammaModel.cc @ 1340

Last change on this file since 1340 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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